immersedBoundaryFvPatch.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | foam-extend: Open Source CFD
4  \\ / O peration | Version: 3.2
5  \\ / A nd | Web: http://www.foam-extend.org
6  \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9  This file is part of foam-extend.
10 
11  foam-extend is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation, either version 3 of the License, or (at your
14  option) any later version.
15 
16  foam-extend is distributed in the hope that it will be useful, but
17  WITHOUT ANY WARRANTY; without even the implied warranty of
18  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19  General Public License for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::immersedBoundaryFvPatch
26 
27 Description
28  Immersed boundary FV patch
29 
30 Author
31  Zeljko Tukovic
32  Reorganisation by Hrvoje Jasak
33 
34 SourceFiles
35  immersedBoundaryFvPatch.C
36  immersedBoundaryFvPatchLeastSquaresFit.C
37  immersedBoundaryFvPatchTriAddressing.C
38  immersedBoundaryFvPatchSamplingWeights.C
39  immersedBoundaryFvPatchTemplates.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef immersedBoundaryFvPatch_H
44 #define immersedBoundaryFvPatch_H
45 
46 #include "fvPatch.H"
48 #include "volFieldsFwd.H"
49 #include "surfaceFieldsFwd.H"
50 #include "dynamicLabelList.H"
51 #include "labelPair.H"
52 #include "FieldFields.H"
53 #include "scalarMatrices.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 
60 // Forward declaration of classes
61 class fvMesh;
62 
63 /*---------------------------------------------------------------------------*\
64  Class immersedBoundaryFvPatch Declaration
65 \*---------------------------------------------------------------------------*/
66 
68 :
69  public fvPatch
70 {
71  // Private data
72 
73  //- Reference to processor patch
75 
76  //- Finite volume mesh reference
77  const fvMesh& mesh_;
78 
79  //- Time index for last update of mesh or moving boundary
80  mutable label ibUpdateTimeIndex_;
81 
82 
83  // Static data
84 
85  //- Fitting angle rejection factor (deg)
86  // Cells within the the radius will be used for the fitting function
87  static const debug::tolerancesSwitch angleFactor_;
88 
89  //- Fitting radius factor
90  // Cells within the the radius will be used for the fitting function
91  static const debug::tolerancesSwitch radiusFactor_;
92 
93  //- Maximum number of rows in cell-cell search
95 
96  //- Sampling point distance factor
97  // Sampling point is located distFactor further away from the wall
98  // from the immersed boundary cell centre
99  static const debug::tolerancesSwitch distFactor_;
100 
101 
102  // Demand-driven data
103 
104  // Immersed boundary data
105 
106  //- Fluid cells indicator, marking only live fluid cells
107  mutable volScalarField* gammaPtr_;
108 
109  //- Fluid cells indicator, marking live and IB cells
110  mutable volScalarField* gammaExtPtr_;
111 
112  //- Fluid faces indicator, marking faces between live cells
114 
115  //- List of fluid cells next to immersed boundary (IB cells)
116  mutable labelList* ibCellsPtr_;
117 
118  //- List of faces for which one neighbour is an IB cell
119  // and another neighbour is a live fluid cell (IB faces)
120  mutable labelList* ibFacesPtr_;
121 
122  //- List of IB cell index for each ibFace
123  mutable labelList* ibFaceCellsPtr_;
124 
125  //- List of IB face flip:
126  // false if IB face points into IB cell (out of the live cell)
127  // true if IB face points into a live cell
128  mutable boolList* ibFaceFlipsPtr_;
129 
130  //- List of fluid faces for which one neighbour is an IB cell
131  // and another neighbour is a dead cell (inside IB faces)
132  mutable labelList* ibInsideFacesPtr_;
133 
134  //- List of internal faces in the region bounded by IB faces
136 
137  //- Points at the immersed boundary (IB points)
138  // nearest to the IB cell centres
139  mutable vectorField* ibPointsPtr_;
140 
141  //- Normals at IB points
142  mutable vectorField* ibNormalsPtr_;
143 
144  //- List of faces (triangles) which are part of IB mesh
145  // nearest to the IB cell centres
146  mutable labelList* hitFacesPtr_;
147 
148  //- Sampling points for IB cells
150 
151  //- Interpolation weights for sampling weights
153 
154  //- Interpolation weights for sampling processor weights
156 
157  //- Interpolation addressing from IB points to tri faces
159 
160  //- Interpolation weights from IB points to tri faces
162 
163  //- Neighbour cells for immersed boundary cells
164  // (extended stencil)
166 
167  //- List of cells needed by neighbour processors
169 
170  //- Centres of cells from neighbour processors
172 
173  //- Gamma of cells from neighbour processors
175 
176  //- Cell-proc-cell addressing
178 
179  //- Dead cells list
180  mutable labelList* deadCellsPtr_;
181 
182  //- Extended dead cells list (dead cells + IB cells)
183  mutable labelList* deadCellsExtPtr_;
184 
185  //- Dead faces list
186  mutable labelList* deadFacesPtr_;
187 
188  //- List of live cells
189  mutable labelList* liveCellsPtr_;
190 
191  //- Average IB cell sizes
192  mutable scalarField* ibCellSizesPtr_;
193 
194  //- Inverse interpolation matrices for Dirichlet BC at the IB
196 
197  //- Inverse interpolation matrices for Neumann BC at the IB
199 
200  //- IB face area vectors
201  mutable vectorField* ibSfPtr_;
202 
203  //- IB face area vector magnitudess
204  mutable scalarField* ibMagSfPtr_;
205 
206  //- IB cell centre distances to IB
207  mutable scalarField* ibDeltaPtr_;
208 
209  //- IB cell centre distances to IB
211 
212 
213  // Tri surface data
214 
215  //- Tri surface face area vectors
216  mutable vectorField* triSfPtr_;
217 
218  //- Tri surface face labels contained in this mesh
219  mutable dynamicLabelList triFacesInMesh_;
220 
221 
222  // Private Member Functions
223 
224 
225  // Storage management
226 
227  //- Clear all demand-driven data
228  void clearOut();
229 
230 
231  // Make demand-driven data
232 
233  //- Make fluid cells indicator, marking only live fluid cells
234  void makeGamma() const;
235 
236  //- Make extended fluid cells indicator, marking live and IB cells
237  void makeGammaExt() const;
238 
239  //- Make fluid faces indicator
240  void makeSGamma() const;
241 
242  //- Make list of cells next to immersed boundary
243  void makeIbCells() const;
244 
245  //- Add corner points to IB cells list
246  void addIbCornerCells() const;
247 
248  //- Make IB faces
249  void makeIbFaces() const;
250 
251  //- Make tri addressing
252  void makeTriAddressing() const;
253 
254  //- Make inside IB faces
255  void makeIbInsideFaces() const;
256 
257  //- Make internal IB faces
258  void makeIbInternalFaces() const;
259 
260  //- Make immersed boundary points and normals
261  void makeIbPointsAndNormals() const;
262 
263  //- Make sampling point weights
264  void makeIbSamplingWeights() const;
265 
266  //- Make extended IB cells stencils
267  void makeIbCellCells() const;
268 
269  //- Make list of dead cells
270  void makeDeadCells() const;
271 
272  //- Make extended list of dead cells
273  void makeDeadCellsExt() const;
274 
275  //- Make list of dead faces
276  void makeDeadFaces() const;
277 
278  //- Make list of live cells
279  void makeLiveCells() const;
280 
281  //- Make immersed boundary cell sizes
282  void makeIbCellSizes() const;
283 
284  //- Make face area vectors and magnitudes
285  void makeIbSf() const;
286 
287  //- Make distance between IB cell centres
288  // and corresponding IB points
289  void makeIbDelta() const;
290 
291  //- Make distance between IB cell centres
292  // and corresponding sample IB points
293  void makeIbSamplingPointDelta() const;
294 
295  //- Make triangular surface face area vectors
296  void makeTriSf() const;
297 
298  //- Find nearest cell
299  label findNearestCell(const point& location) const;
300 
301  //- Return extended cell-cell addressing
302  void findCellCells
303  (
304  const vector& pt,
305  const label cellID,
306  labelList& cellCells
307  ) const;
308 
309  //- Calc cell size
310  scalar cellSize(label cellID) const;
311 
312  //- Calc cell projection area
313  scalar cellProjection(label cellID, const vector& dir) const;
314 
315 
316  // Boundary evaluation matrices
317 
318  //- Make inverse Dirichlet interpolation matrices
319  void makeInvDirichletMatrices() const;
320 
321  //- Make inverse Neumann interpolation matrices
322  void makeInvNeumannMatrices() const;
323 
324 
325 protected:
326 
327  // Protected Member Functions
328 
329  //- Initialise the patches for moving points
330  virtual void initMovePoints();
331 
332  //- Correct patches after moving points
333  virtual void movePoints();
334 
335 
336 public:
337 
338  //- Runtime type information
339  TypeName(immersedBoundaryPolyPatch::typeName_());
340 
341 
342  // Constructors
343 
344  //- Construct from polyPatch
346  (
347  const polyPatch& patch,
348  const fvBoundaryMesh& bm
349  );
350 
351 
352  //- Destructor
353  virtual ~immersedBoundaryFvPatch()
354  {}
355 
356 
357  // Member Functions
358 
359  // Access to immersed boundary components
360 
361  //- Return reference to immersed boundary polyPatch
363  {
364  return ibPolyPatch_;
365  }
366 
367  //- Return immersed boundary surface mesh
368  const triSurfaceMesh& ibMesh() const
369  {
370  return ibPolyPatch_.ibMesh();
371  }
372 
373  bool internalFlow() const
374  {
375  return ibPolyPatch_.internalFlow();
376  }
377 
378  //- Is the immersed boundary patch moving?
379  bool movingIb() const
380  {
381  return ibPolyPatch_.movingIb();
382  }
383 
384 
385  // Immersed boundary data access
386 
387  //- Get fluid cells indicator, marking only live fluid cells
388  const volScalarField& gamma() const;
389 
390  //- Get extended flud cells indicator, including live and IB cells
391  const volScalarField& gammaExt() const;
392 
393  //- Get fluid faces indicator, marking faces between live cells
394  const surfaceScalarField& sGamma() const;
395 
396 
397  //- Return list of fluid cells next to immersed boundary (IB cells)
398  const labelList& ibCells() const;
399 
400  //- Return list of faces for which one neighbour is an IB cell
401  // and another neighbour is a live fluid cell (IB faces)
402  const labelList& ibFaces() const;
403 
404  //- Return list of IB cell index for each ibFace
405  const labelList& ibFaceCells() const;
406 
407  //- Return list of IB face flip:
408  // false if IB face points into IB cell (out of the live cell)
409  // true if IB face points into a live cell
410  const boolList& ibFaceFlips() const;
411 
412  //- Return list of fluid faces for which one neighbour is an
413  // IB cell and another neighbour is a dead cell (inside IB faces)
414  const labelList& ibInsideFaces() const;
415 
416  //- Return list of internal faces in the region bounded by IB faces
417  const labelList& ibInternalFaces() const;
418 
419  //- Return IB points
420  const vectorField& ibPoints() const;
421 
422  //- Return IB normals
423  const vectorField& ibNormals() const;
424 
425  //- Return list of triangles in IB mesh nearest
426  // nearest to IB cell centres
427  const labelList& hitFaces() const;
428 
429  //- Return IB sampling points
430  const vectorField& ibSamplingPoints() const;
431 
432  //- Interpolation weights for sampling points
433  const scalarListList& ibSamplingWeights() const;
434 
435  //- Processor interpolation weights for sampling points
436  const scalarListList& ibSamplingProcWeights() const;
437 
438  //- Interpolation addressing from IB points to tri faces
439  const labelListList& cellsToTriAddr() const;
440 
441  //- Interpolation weights from IB points to tri faces
442  const scalarListList& cellsToTriWeights() const;
443 
444  //- Return IB cell extended stencil
445  const labelListList& ibCellCells() const;
446 
447  //- Return neighbour proc centres
448  // These are centres from neighbouring processors the
449  // local processor needs to receive
451 
452  //- Return neighbour proc gamma
453  // These are gamma values from neighbouring processors the
454  // local processor needs to receive
455  const FieldField<Field, scalar>& ibProcGamma() const;
456 
457  //- Return neighbour cell addressing
458  const List<List<labelPair> >& ibCellProcCells() const;
459 
460  //- Return neighbour proc cells
461  // These are cell data that other processors need from
462  // local processor
463  const labelListList& ibProcCells() const;
464 
465  //- Return dead cells
466  const labelList& deadCells() const;
467 
468  //- Return extended dead cells
469  const labelList& deadCellsExt() const;
470 
471  //- Return dead faces
472  const labelList& deadFaces() const;
473 
474  //- Return live cells
475  const labelList& liveCells() const;
476 
477  //- Return immersed boundary cell sizes
478  const scalarField& ibCellSizes() const;
479 
480  //- Get inverse Dirichlet interpolation matrix
482  invDirichletMatrices() const;
483 
484  //- Get inverse Neumann interpolation matrix
486  invNeumannMatrices() const;
487 
488 
489  //- Return IB face area vectors
490  const vectorField& ibSf() const;
491 
492  //- Return IB face area vector magnitudes
493  const scalarField& ibMagSf() const;
494 
495  //- Return distance to IB
496  const scalarField& ibDelta() const;
497 
498  //- Return distance to IB
499  const scalarField& ibSamplingPointDelta() const;
500 
501  //- Return triangular surface face area vectors
502  const vectorField& triSf() const;
503 
504  //- Return triangular surface face centres
505  const vectorField& triCf() const;
506 
507  //- Return labels of triangular faces which are inside the mesh
508  const dynamicLabelList& triFacesInMesh() const;
509 
510 
511  // Parallel communication
512 
513  //- Send and receive
514  template<class Type>
516  (
517  const Field<Type>& psi
518  ) const;
519 
520 
521  // Interpolation functions to and from triangular surface
522 
523  //- Collect ibPoint values: from tri face fields onto intersection
524  // points on the IB surface
525  template<class Type>
527  (
528  const Field<Type>& triValues
529  ) const;
530 
531  //- Collect ibPoint values: from tri face fields onto intersection
532  // points on the IB surface with tmp
533  template<class Type>
535  (
536  const tmp<Field<Type> >& ttriValues
537  ) const;
538 
539  //- triFace values: collect data from IB fields onto intersection
540  // points on the IB surface
541  template<class Type>
543  (
544  const Field<Type>& ibValues
545  ) const;
546 
547  //- triFace values: collect data from IB fields onto intersection
548  // points on the IB surface with tmp
549  template<class Type>
551  (
552  const tmp<Field<Type> >& tibValues
553  ) const;
554 
555 
556  //- Interpolation functions to sampling points from mesh cell centres
557  template<class Type>
559  (
560  const Field<Type>& cellValues
561  ) const;
562 
563 
564  // Helper functions
565 
566  //- Renumber Field to corespond to triangular faces contained
567  // inside the mesh
568  template<class Type>
569  const tmp<Field<Type> > renumberField(const Field<Type>& f) const;
570 
571 };
572 
573 
574 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
575 
576 } // End namespace Foam
577 
578 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
579 
580 #ifdef NoRepository
582 #endif
583 
584 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
585 
586 #endif
587 
588 // ************************************************************************* //
Foam::immersedBoundaryFvPatch::sGammaPtr_
surfaceScalarField * sGammaPtr_
Fluid faces indicator, marking faces between live cells.
Definition: immersedBoundaryFvPatch.H:112
Foam::immersedBoundaryFvPatch::ibProcGammaPtr_
FieldField< Field, scalar > * ibProcGammaPtr_
Gamma of cells from neighbour processors.
Definition: immersedBoundaryFvPatch.H:173
Foam::immersedBoundaryFvPatch::ibFaceFlipsPtr_
boolList * ibFaceFlipsPtr_
List of IB face flip:
Definition: immersedBoundaryFvPatch.H:127
Foam::immersedBoundaryFvPatch::maxCellCellRows_
static const debug::optimisationSwitch maxCellCellRows_
Maximum number of rows in cell-cell search.
Definition: immersedBoundaryFvPatch.H:93
volFieldsFwd.H
Foam::immersedBoundaryFvPatch::ibMesh
const triSurfaceMesh & ibMesh() const
Return immersed boundary surface mesh.
Definition: immersedBoundaryFvPatch.H:367
Foam::immersedBoundaryFvPatch::gammaExtPtr_
volScalarField * gammaExtPtr_
Fluid cells indicator, marking live and IB cells.
Definition: immersedBoundaryFvPatch.H:109
Foam::immersedBoundaryFvPatch::ibSamplingPointDelta
const scalarField & ibSamplingPointDelta() const
Return distance to IB.
Definition: immersedBoundaryFvPatch.C:2547
Foam::immersedBoundaryFvPatch::ibDeltaPtr_
scalarField * ibDeltaPtr_
IB cell centre distances to IB.
Definition: immersedBoundaryFvPatch.H:206
Foam::immersedBoundaryFvPatch::ibSamplingPointsPtr_
vectorField * ibSamplingPointsPtr_
Sampling points for IB cells.
Definition: immersedBoundaryFvPatch.H:148
Foam::immersedBoundaryFvPatch::triFacesInMesh_
dynamicLabelList triFacesInMesh_
Tri surface face labels contained in this mesh.
Definition: immersedBoundaryFvPatch.H:218
Foam::immersedBoundaryFvPatch::ibSf
const vectorField & ibSf() const
Return IB face area vectors.
Definition: immersedBoundaryFvPatch.C:2513
Foam::immersedBoundaryFvPatch::addIbCornerCells
void addIbCornerCells() const
Add corner points to IB cells list.
Definition: immersedBoundaryFvPatch.C:476
Foam::immersedBoundaryFvPatch::ibCellSizesPtr_
scalarField * ibCellSizesPtr_
Average IB cell sizes.
Definition: immersedBoundaryFvPatch.H:191
Foam::immersedBoundaryFvPatch::makeDeadFaces
void makeDeadFaces() const
Make list of dead faces.
Definition: immersedBoundaryFvPatch.C:1661
Foam::immersedBoundaryFvPatch::makeGamma
void makeGamma() const
Make fluid cells indicator, marking only live fluid cells.
Definition: immersedBoundaryFvPatch.C:83
Foam::immersedBoundaryFvPatch::ibCellSizes
const scalarField & ibCellSizes() const
Return immersed boundary cell sizes.
Definition: immersedBoundaryFvPatch.C:2502
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
Foam::immersedBoundaryFvPatch::movingIb
bool movingIb() const
Is the immersed boundary patch moving?
Definition: immersedBoundaryFvPatch.H:378
Foam::immersedBoundaryFvPatch::TypeName
TypeName(immersedBoundaryPolyPatch::typeName_())
Runtime type information.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::immersedBoundaryFvPatch::hitFaces
const labelList & hitFaces() const
Return list of triangles in IB mesh nearest.
Definition: immersedBoundaryFvPatch.C:2377
Foam::immersedBoundaryFvPatch::ibCellProcCells
const List< List< labelPair > > & ibCellProcCells() const
Return neighbour cell addressing.
Definition: immersedBoundaryFvPatch.C:2436
Foam::immersedBoundaryFvPatch::ibSamplingProcWeightsPtr_
scalarListList * ibSamplingProcWeightsPtr_
Interpolation weights for sampling processor weights.
Definition: immersedBoundaryFvPatch.H:154
Foam::immersedBoundaryFvPatch::ibProcGamma
const FieldField< Field, scalar > & ibProcGamma() const
Return neighbour proc gamma.
Definition: immersedBoundaryFvPatch.C:2424
Foam::immersedBoundaryFvPatch::ibFacesPtr_
labelList * ibFacesPtr_
List of faces for which one neighbour is an IB cell.
Definition: immersedBoundaryFvPatch.H:119
Foam::immersedBoundaryFvPatch::mesh_
const fvMesh & mesh_
Finite volume mesh reference.
Definition: immersedBoundaryFvPatch.H:76
FieldFields.H
Foam::immersedBoundaryFvPatch::makeIbSf
void makeIbSf() const
Make face area vectors and magnitudes.
Definition: immersedBoundaryFvPatch.C:1842
Foam::immersedBoundaryFvPatch::triSf
const vectorField & triSf() const
Return triangular surface face area vectors.
Definition: immersedBoundaryFvPatch.C:2558
Foam::immersedBoundaryFvPatch::makeIbSamplingWeights
void makeIbSamplingWeights() const
Make sampling point weights.
Definition: immersedBoundaryFvPatchSamplingWeights.C:32
Foam::immersedBoundaryFvPatch::liveCells
const labelList & liveCells() const
Return live cells.
Definition: immersedBoundaryFvPatch.C:2491
Foam::immersedBoundaryFvPatch::cellsToTriWeights
const scalarListList & cellsToTriWeights() const
Interpolation weights from IB points to tri faces.
Definition: immersedBoundaryFvPatchTriAddressing.C:223
Foam::immersedBoundaryFvPatch::ibFaceCellsPtr_
labelList * ibFaceCellsPtr_
List of IB cell index for each ibFace.
Definition: immersedBoundaryFvPatch.H:122
Foam::immersedBoundaryFvPatch::invNeumannMatrices
const PtrList< scalarRectangularMatrix > & invNeumannMatrices() const
Get inverse Neumann interpolation matrix.
Definition: immersedBoundaryFvPatchLeastSquaresFit.C:582
Foam::immersedBoundaryFvPatch::clearOut
void clearOut()
Clear all demand-driven data.
Definition: immersedBoundaryFvPatch.C:2149
Foam::immersedBoundaryFvPatch::movePoints
virtual void movePoints()
Correct patches after moving points.
Definition: immersedBoundaryFvPatch.C:2200
Foam::immersedBoundaryFvPatch::invDirichletMatricesPtr_
PtrList< scalarRectangularMatrix > * invDirichletMatricesPtr_
Inverse interpolation matrices for Dirichlet BC at the IB.
Definition: immersedBoundaryFvPatch.H:194
Foam::immersedBoundaryFvPatch::ibSamplingPoints
const vectorField & ibSamplingPoints() const
Return IB sampling points.
Definition: immersedBoundaryFvPatch.C:2389
Foam::immersedBoundaryFvPatch::ibFaceCells
const labelList & ibFaceCells() const
Return list of IB cell index for each ibFace.
Definition: immersedBoundaryFvPatch.C:2311
Foam::immersedBoundaryPolyPatch
Immersed boundary patch.
Definition: immersedBoundaryPolyPatch.H:56
Foam::immersedBoundaryFvPatch::sGamma
const surfaceScalarField & sGamma() const
Get fluid faces indicator, marking faces between live cells.
Definition: immersedBoundaryFvPatch.C:2278
Foam::immersedBoundaryFvPatch::makeDeadCells
void makeDeadCells() const
Make list of dead cells.
Definition: immersedBoundaryFvPatch.C:1592
Foam::triSurfaceMesh
IOoject and searching on triSurface.
Definition: triSurfaceMesh.H:63
Foam::immersedBoundaryFvPatch::ibCellCellsPtr_
labelListList * ibCellCellsPtr_
Neighbour cells for immersed boundary cells.
Definition: immersedBoundaryFvPatch.H:164
Foam::immersedBoundaryFvPatch::deadCells
const labelList & deadCells() const
Return dead cells.
Definition: immersedBoundaryFvPatch.C:2458
Foam::immersedBoundaryFvPatch::gammaPtr_
volScalarField * gammaPtr_
Fluid cells indicator, marking only live fluid cells.
Definition: immersedBoundaryFvPatch.H:106
Foam::immersedBoundaryFvPatch::ibProcCentresPtr_
FieldField< Field, vector > * ibProcCentresPtr_
Centres of cells from neighbour processors.
Definition: immersedBoundaryFvPatch.H:170
Foam::immersedBoundaryFvPatch::deadCellsPtr_
labelList * deadCellsPtr_
Dead cells list.
Definition: immersedBoundaryFvPatch.H:179
Foam::immersedBoundaryFvPatch::ibProcCellsPtr_
labelListList * ibProcCellsPtr_
List of cells needed by neighbour processors.
Definition: immersedBoundaryFvPatch.H:167
Foam::immersedBoundaryFvPatch::findNearestCell
label findNearestCell(const point &location) const
Find nearest cell.
Definition: immersedBoundaryFvPatch.C:1962
Foam::immersedBoundaryFvPatch::gamma
const volScalarField & gamma() const
Get fluid cells indicator, marking only live fluid cells.
Definition: immersedBoundaryFvPatch.C:2256
Foam::immersedBoundaryFvPatch::makeDeadCellsExt
void makeDeadCellsExt() const
Make extended list of dead cells.
Definition: immersedBoundaryFvPatch.C:1627
Foam::immersedBoundaryPolyPatch::internalFlow
bool internalFlow() const
Return true if solving for flow inside the immersed boundary.
Definition: immersedBoundaryPolyPatch.H:198
Foam::immersedBoundaryFvPatch::triFacesInMesh
const dynamicLabelList & triFacesInMesh() const
Return labels of triangular faces which are inside the mesh.
Definition: immersedBoundaryFvPatch.C:2576
Foam::immersedBoundaryFvPatch::ibInternalFacesPtr_
labelList * ibInternalFacesPtr_
List of internal faces in the region bounded by IB faces.
Definition: immersedBoundaryFvPatch.H:134
Foam::immersedBoundaryFvPatch::ibFaces
const labelList & ibFaces() const
Return list of faces for which one neighbour is an IB cell.
Definition: immersedBoundaryFvPatch.C:2300
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:52
Foam::debug::optimisationSwitch
int optimisationSwitch(const char *name, const int defaultValue=0)
Lookup optimisation switch or add default value.
Definition: debug.C:182
Foam::immersedBoundaryFvPatch::cellProjection
scalar cellProjection(label cellID, const vector &dir) const
Calc cell projection area.
Definition: immersedBoundaryFvPatch.C:2115
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::immersedBoundaryFvPatch::makeSGamma
void makeSGamma() const
Make fluid faces indicator.
Definition: immersedBoundaryFvPatch.C:217
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::immersedBoundaryFvPatch::ibCells
const labelList & ibCells() const
Return list of fluid cells next to immersed boundary (IB cells)
Definition: immersedBoundaryFvPatch.C:2289
Foam::immersedBoundaryFvPatch::ibInsideFacesPtr_
labelList * ibInsideFacesPtr_
List of fluid faces for which one neighbour is an IB cell.
Definition: immersedBoundaryFvPatch.H:131
Foam::immersedBoundaryFvPatch
Immersed boundary FV patch.
Definition: immersedBoundaryFvPatch.H:66
Foam::immersedBoundaryFvPatch::immersedBoundaryFvPatch
immersedBoundaryFvPatch(const polyPatch &patch, const fvBoundaryMesh &bm)
Construct from polyPatch.
Definition: immersedBoundaryFvPatch.C:2209
Foam::immersedBoundaryFvPatch::initMovePoints
virtual void initMovePoints()
Initialise the patches for moving points.
Definition: immersedBoundaryFvPatch.C:2196
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::immersedBoundaryFvPatch::ibCellProcCellsPtr_
List< List< labelPair > > * ibCellProcCellsPtr_
Cell-proc-cell addressing.
Definition: immersedBoundaryFvPatch.H:176
Foam::immersedBoundaryFvPatch::ibSamplingPointDeltaPtr_
scalarField * ibSamplingPointDeltaPtr_
IB cell centre distances to IB.
Definition: immersedBoundaryFvPatch.H:209
Foam::immersedBoundaryFvPatch::makeGammaExt
void makeGammaExt() const
Make extended fluid cells indicator, marking live and IB cells.
Definition: immersedBoundaryFvPatch.C:145
Foam::immersedBoundaryFvPatch::ibSamplingProcWeights
const scalarListList & ibSamplingProcWeights() const
Processor interpolation weights for sampling points.
Definition: immersedBoundaryFvPatchSamplingWeights.C:196
Foam::immersedBoundaryFvPatch::makeTriSf
void makeTriSf() const
Make triangular surface face area vectors.
Definition: immersedBoundaryFvPatch.C:1922
Foam::immersedBoundaryFvPatch::ibDelta
const scalarField & ibDelta() const
Return distance to IB.
Definition: immersedBoundaryFvPatch.C:2535
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::immersedBoundaryFvPatch::distFactor_
static const debug::tolerancesSwitch distFactor_
Sampling point distance factor.
Definition: immersedBoundaryFvPatch.H:98
Foam::immersedBoundaryFvPatch::triCf
const vectorField & triCf() const
Return triangular surface face centres.
Definition: immersedBoundaryFvPatch.C:2569
Foam::immersedBoundaryFvPatch::ibUpdateTimeIndex_
label ibUpdateTimeIndex_
Time index for last update of mesh or moving boundary.
Definition: immersedBoundaryFvPatch.H:79
Foam::immersedBoundaryFvPatch::ibMagSf
const scalarField & ibMagSf() const
Return IB face area vector magnitudes.
Definition: immersedBoundaryFvPatch.C:2524
immersedBoundaryPolyPatch.H
Foam::immersedBoundaryFvPatch::makeInvDirichletMatrices
void makeInvDirichletMatrices() const
Make inverse Dirichlet interpolation matrices.
Definition: immersedBoundaryFvPatchLeastSquaresFit.C:32
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::immersedBoundaryFvPatch::ibMagSfPtr_
scalarField * ibMagSfPtr_
IB face area vector magnitudess.
Definition: immersedBoundaryFvPatch.H:203
Foam::immersedBoundaryFvPatch::ibSfPtr_
vectorField * ibSfPtr_
IB face area vectors.
Definition: immersedBoundaryFvPatch.H:200
Foam::immersedBoundaryFvPatch::ibPointsPtr_
vectorField * ibPointsPtr_
Points at the immersed boundary (IB points)
Definition: immersedBoundaryFvPatch.H:138
Foam::immersedBoundaryFvPatch::makeIbPointsAndNormals
void makeIbPointsAndNormals() const
Make immersed boundary points and normals.
Definition: immersedBoundaryFvPatch.C:946
Foam::immersedBoundaryFvPatch::gammaExt
const volScalarField & gammaExt() const
Get extended flud cells indicator, including live and IB cells.
Definition: immersedBoundaryFvPatch.C:2267
Foam::immersedBoundaryFvPatch::makeIbCellCells
void makeIbCellCells() const
Make extended IB cells stencils.
Definition: immersedBoundaryFvPatch.C:1065
Foam::immersedBoundaryFvPatch::ibInternalFaces
const labelList & ibInternalFaces() const
Return list of internal faces in the region bounded by IB faces.
Definition: immersedBoundaryFvPatch.C:2344
Foam::immersedBoundaryFvPatch::makeInvNeumannMatrices
void makeInvNeumannMatrices() const
Make inverse Neumann interpolation matrices.
Definition: immersedBoundaryFvPatchLeastSquaresFit.C:268
Foam::immersedBoundaryFvPatch::ibFaceFlips
const boolList & ibFaceFlips() const
Return list of IB face flip:
Definition: immersedBoundaryFvPatch.C:2322
Foam::immersedBoundaryFvPatch::makeIbInternalFaces
void makeIbInternalFaces() const
Make internal IB faces.
Definition: immersedBoundaryFvPatch.C:852
Foam::immersedBoundaryFvPatch::cellsToTriAddrPtr_
labelListList * cellsToTriAddrPtr_
Interpolation addressing from IB points to tri faces.
Definition: immersedBoundaryFvPatch.H:157
Foam::immersedBoundaryFvPatch::ibNormalsPtr_
vectorField * ibNormalsPtr_
Normals at IB points.
Definition: immersedBoundaryFvPatch.H:141
Foam::immersedBoundaryFvPatch::ibCellCells
const labelListList & ibCellCells() const
Return IB cell extended stencil.
Definition: immersedBoundaryFvPatch.C:2400
Foam::immersedBoundaryPolyPatch::ibMesh
const triSurfaceMesh & ibMesh() const
Return immersed boundary surface mesh.
Definition: immersedBoundaryPolyPatch.H:192
Foam::immersedBoundaryFvPatch::triSfPtr_
vectorField * triSfPtr_
Tri surface face area vectors.
Definition: immersedBoundaryFvPatch.H:215
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::immersedBoundaryFvPatch::ibSamplingWeightsPtr_
scalarListList * ibSamplingWeightsPtr_
Interpolation weights for sampling weights.
Definition: immersedBoundaryFvPatch.H:151
Foam::immersedBoundaryFvPatch::ibPolyPatch_
const immersedBoundaryPolyPatch & ibPolyPatch_
Reference to processor patch.
Definition: immersedBoundaryFvPatch.H:73
Foam::immersedBoundaryFvPatch::ibNormals
const vectorField & ibNormals() const
Return IB normals.
Definition: immersedBoundaryFvPatch.C:2366
Foam::immersedBoundaryFvPatch::deadCellsExtPtr_
labelList * deadCellsExtPtr_
Extended dead cells list (dead cells + IB cells)
Definition: immersedBoundaryFvPatch.H:182
Foam::immersedBoundaryFvPatch::angleFactor_
static const debug::tolerancesSwitch angleFactor_
Fitting angle rejection factor (deg)
Definition: immersedBoundaryFvPatch.H:86
f
labelList f(nPoints)
Foam::immersedBoundaryFvPatch::ibCellsPtr_
labelList * ibCellsPtr_
List of fluid cells next to immersed boundary (IB cells)
Definition: immersedBoundaryFvPatch.H:115
Foam::immersedBoundaryFvPatch::toIbPoints
tmp< Field< Type > > toIbPoints(const Field< Type > &triValues) const
Collect ibPoint values: from tri face fields onto intersection.
Foam::Vector< scalar >
Foam::immersedBoundaryPolyPatch::movingIb
bool movingIb() const
Return true if immersed boundary is moving.
Definition: immersedBoundaryPolyPatch.H:207
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::immersedBoundaryFvPatch::deadFaces
const labelList & deadFaces() const
Return dead faces.
Definition: immersedBoundaryFvPatch.C:2480
Foam::immersedBoundaryFvPatch::ibSamplingWeights
const scalarListList & ibSamplingWeights() const
Interpolation weights for sampling points.
Definition: immersedBoundaryFvPatchSamplingWeights.C:184
Foam::immersedBoundaryFvPatch::makeIbCells
void makeIbCells() const
Make list of cells next to immersed boundary.
Definition: immersedBoundaryFvPatch.C:371
Foam::immersedBoundaryFvPatch::radiusFactor_
static const debug::tolerancesSwitch radiusFactor_
Fitting radius factor.
Definition: immersedBoundaryFvPatch.H:90
Foam::immersedBoundaryFvPatch::deadCellsExt
const labelList & deadCellsExt() const
Return extended dead cells.
Definition: immersedBoundaryFvPatch.C:2469
Foam::immersedBoundaryFvPatch::invNeumannMatricesPtr_
PtrList< scalarRectangularMatrix > * invNeumannMatricesPtr_
Inverse interpolation matrices for Neumann BC at the IB.
Definition: immersedBoundaryFvPatch.H:197
Foam::immersedBoundaryFvPatch::internalFlow
bool internalFlow() const
Definition: immersedBoundaryFvPatch.H:372
Foam::immersedBoundaryFvPatch::makeIbFaces
void makeIbFaces() const
Make IB faces.
Definition: immersedBoundaryFvPatch.C:642
immersedBoundaryFvPatchTemplates.C
surfaceFieldsFwd.H
Foam::immersedBoundaryFvPatch::ibPolyPatch
const immersedBoundaryPolyPatch & ibPolyPatch() const
Return reference to immersed boundary polyPatch.
Definition: immersedBoundaryFvPatch.H:361
Foam::immersedBoundaryFvPatch::ibInsideFaces
const labelList & ibInsideFaces() const
Return list of fluid faces for which one neighbour is an.
Definition: immersedBoundaryFvPatch.C:2333
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Foam::immersedBoundaryFvPatch::hitFacesPtr_
labelList * hitFacesPtr_
List of faces (triangles) which are part of IB mesh.
Definition: immersedBoundaryFvPatch.H:145
Foam::immersedBoundaryFvPatch::makeIbInsideFaces
void makeIbInsideFaces() const
Make inside IB faces.
Definition: immersedBoundaryFvPatch.C:765
Foam::immersedBoundaryFvPatch::makeLiveCells
void makeLiveCells() const
Make list of live cells.
Definition: immersedBoundaryFvPatch.C:1745
Foam::immersedBoundaryFvPatch::deadFacesPtr_
labelList * deadFacesPtr_
Dead faces list.
Definition: immersedBoundaryFvPatch.H:185
Foam::immersedBoundaryFvPatch::ibPoints
const vectorField & ibPoints() const
Return IB points.
Definition: immersedBoundaryFvPatch.C:2355
scalarMatrices.H
Foam::immersedBoundaryFvPatch::findCellCells
void findCellCells(const vector &pt, const label cellID, labelList &cellCells) const
Return extended cell-cell addressing.
Definition: immersedBoundaryFvPatch.C:1992
Foam::immersedBoundaryFvPatch::cellsToTriAddr
const labelListList & cellsToTriAddr() const
Interpolation addressing from IB points to tri faces.
Definition: immersedBoundaryFvPatchTriAddressing.C:211
fvPatch.H
Foam::immersedBoundaryFvPatch::toSamplingPoints
tmp< Field< Type > > toSamplingPoints(const Field< Type > &cellValues) const
Interpolation functions to sampling points from mesh cell centres.
Foam::immersedBoundaryFvPatch::cellSize
scalar cellSize(label cellID) const
Calc cell size.
Definition: immersedBoundaryFvPatch.C:2086
Foam::immersedBoundaryFvPatch::makeIbDelta
void makeIbDelta() const
Make distance between IB cell centres.
Definition: immersedBoundaryFvPatch.C:1868
Foam::immersedBoundaryFvPatch::ibProcCells
const labelListList & ibProcCells() const
Return neighbour proc cells.
Definition: immersedBoundaryFvPatch.C:2447
Foam::immersedBoundaryFvPatch::makeIbCellSizes
void makeIbCellSizes() const
Make immersed boundary cell sizes.
Definition: immersedBoundaryFvPatch.C:1779
Foam::immersedBoundaryFvPatch::renumberField
const tmp< Field< Type > > renumberField(const Field< Type > &f) const
Renumber Field to corespond to triangular faces contained.
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::immersedBoundaryFvPatch::sendAndReceive
tmp< FieldField< Field, Type > > sendAndReceive(const Field< Type > &psi) const
Send and receive.
Foam::immersedBoundaryFvPatch::makeIbSamplingPointDelta
void makeIbSamplingPointDelta() const
Make distance between IB cell centres.
Definition: immersedBoundaryFvPatch.C:1894
Foam::immersedBoundaryFvPatch::makeTriAddressing
void makeTriAddressing() const
Make tri addressing.
Definition: immersedBoundaryFvPatchTriAddressing.C:30
Foam::immersedBoundaryFvPatch::cellsToTriWeightsPtr_
scalarListList * cellsToTriWeightsPtr_
Interpolation weights from IB points to tri faces.
Definition: immersedBoundaryFvPatch.H:160
Foam::immersedBoundaryFvPatch::liveCellsPtr_
labelList * liveCellsPtr_
List of live cells.
Definition: immersedBoundaryFvPatch.H:188
Foam::immersedBoundaryFvPatch::toTriFaces
tmp< Field< Type > > toTriFaces(const Field< Type > &ibValues) const
triFace values: collect data from IB fields onto intersection
Foam::immersedBoundaryFvPatch::~immersedBoundaryFvPatch
virtual ~immersedBoundaryFvPatch()
Destructor.
Definition: immersedBoundaryFvPatch.H:352
Foam::immersedBoundaryFvPatch::ibProcCentres
const FieldField< Field, vector > & ibProcCentres() const
Return neighbour proc centres.
Definition: immersedBoundaryFvPatch.C:2412
Foam::immersedBoundaryFvPatch::invDirichletMatrices
const PtrList< scalarRectangularMatrix > & invDirichletMatrices() const
Get inverse Dirichlet interpolation matrix.
Definition: immersedBoundaryFvPatchLeastSquaresFit.C:570
labelPair.H