faMesh.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2017 Wikki Ltd
9  Copyright (C) 2021-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::faMesh
29 
30 Description
31  Finite area mesh (used for 2-D non-Euclidian finite area method)
32  defined using a \em patch of faces on a polyMesh
33  (ie, uindirectPrimitivePatch).
34 
35  The ordering of faces and points on the faMesh corresponds to
36  the localFaces and localPoints as per Foam::PrimitivePatch but
37  the edge addressing is handled slightly differently.
38  The internal edges of the faMesh will generally correspond identically
39  to the internalEdges of the PrimitivePatch (may change in the future)
40  but the boundaryEdges will be reordered compared to the PrimitivePatch
41  to allow edge boundary slices to be obtained.
42 
43 SourceFiles
44  faMesh.C
45  faMeshDemandDrivenData.C
46  faMeshPatches.C
47  faMeshTopology.C
48  faMeshUpdate.C
49 
50 Author
51  Zeljko Tukovic, FMENA
52  Hrvoje Jasak, Wikki Ltd.
53 
54 \*---------------------------------------------------------------------------*/
55 
56 #ifndef Foam_faMesh_H
57 #define Foam_faMesh_H
58 
59 #include "polyMesh.H"
60 #include "lduMesh.H"
61 #include "faBoundaryMesh.H"
62 #include "edgeList.H"
63 #include "faceList.H"
64 #include "primitiveFieldsFwd.H"
65 #include "DimensionedField.H"
66 #include "areaFieldsFwd.H"
67 #include "edgeFieldsFwd.H"
68 #include "indirectPrimitivePatch.H"
69 #include "edgeInterpolation.H"
70 #include "labelIOList.H"
71 #include "FieldFields.H"
72 #include "faGlobalMeshData.H"
73 #include "faSchemes.H"
74 #include "faSolution.H"
75 
76 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 
78 namespace Foam
79 {
80 
81 // Forward Declarations
82 class faMeshBoundaryHalo;
83 class faMeshLduAddressing;
84 class faMeshMapper;
85 class faPatchData;
86 
87 /*---------------------------------------------------------------------------*\
88  Class faMeshRegistry Declaration
89 \*---------------------------------------------------------------------------*/
90 
91 //- The objectRegistry for faMesh.
92 // This is a separate class to ensure it will be fully constructed before
93 // other data members use it, which ensures that its virtual methods are
94 // callable during construction (gcc).
95 class faMeshRegistry
96 :
97  public objectRegistry
98 {
99 public:
100 
101  // Constructors
102 
103  //- Construct an objectRegistry for given area region name.
104  // Uses faMeshesRegistry internally
105  faMeshRegistry(const word& areaRegion, const polyMesh& mesh);
106 
107 
108  // Member Functions
109 
110  //- True - thisDb() is a valid registry
111  virtual bool hasDb() const { return true; }
112 
113  //- Reference to the mesh database
114  virtual const objectRegistry& thisDb() const { return *this; }
115 
116  //- Local directory path of the objectRegistry relative to Time
117  //- with override for the single-region case
118  virtual const fileName& dbDir() const;
119 };
120 
121 
122 /*---------------------------------------------------------------------------*\
123  Class faMesh Declaration
124 \*---------------------------------------------------------------------------*/
125 
126 class faMesh
127 :
128  public faMeshRegistry,
129  public lduMesh,
130  public faSchemes,
131  public faSolution,
132  public edgeInterpolation // may need input from faSchemes
133 {
134  // Private (internal) classes/structures
135 
136  //- A (proc, patchi, patchEdgei, meshFacei) tuple used internally
137  //- for managing patch/patch bookkeeping during construction.
138  // Finite-area patches are stored with negated indices (offset -2),
139  // which makes them readily identifiable and always sort before normal
140  // patches.
141  struct patchTuple
142  :
143  public FixedList<label, 4>
144  {
145 
146  // Constructors
147 
148  // Inherit constructors
150 
151  //- Default construct as 'invalid'
152  patchTuple()
153  {
154  clear();
155  }
156 
157 
158  // Static Member Functions
159 
160  // Globally consistent ordering
161  //
162  // 1. sort left/right as lower/higher processor connection
163  // 2. sort by proc/patch/patch index
164  static void sort(UList<Pair<patchTuple>>& list)
165  {
166  for (auto& tuples : list)
167  {
168  tuples.sort();
169  }
170  Foam::stableSort(list);
171  }
172 
173 
174  // Member Functions
175 
176  //- Reset to 'invalid'
177  void clear()
178  {
179  procNo(-1);
180  patchi(labelMax);
181  patchEdgei(-1);
182  meshFacei(-1);
183  }
184 
185  //- Valid if proc and patch (or patch edge) are non-negative
186  bool valid() const
187  {
188  return (procNo() >= 0 && patchi() != -1);
189  }
190 
191  // Processor is the first sort index
192  label procNo() const { return (*this)[0]; }
193  void procNo(label val) { (*this)[0] = val; }
194 
195  // PatchId is the second sort index (finiteArea patches are < -1)
196  label patchi() const { return (*this)[1]; }
197  void patchi(label val) { (*this)[1] = val; }
198 
199  // The patch edge index (on the finiteArea patch)
200  // is the third sort index
201  label patchEdgei() const { return (*this)[2]; }
202  void patchEdgei(label val) { (*this)[2] = val; }
203 
204  // The processor-local mesh face is the fourth sort index
205  label meshFacei() const { return (*this)[3]; }
206  void meshFacei(label val) { (*this)[3] = val; }
207 
208  //- Return the real patchId
209  label realPatchi() const
210  {
211  const label id = patchi();
212  return (id < -1 ? -(id + 2) : id);
213  }
214 
215  //- Set patchId as finiteArea
216  void faPatchi(label val)
217  {
218  patchi(-(val + 2));
219  }
220 
221  //- Considered to be finiteArea if (patchi < -1)
222  bool is_finiteArea() const noexcept
223  {
224  return (patchi() < -1);
225  }
226 
227  //- Considered to be processor local
228  bool is_localProc() const
229  {
230  return (procNo() == UPstream::myProcNo());
231  }
232  };
233 
234 
235  // Private Data
236 
237  //- Face labels
238  labelIOList faceLabels_;
239 
240  //- Boundary mesh
241  faBoundaryMesh boundary_;
242 
243 
244  // Primitive mesh data
245 
246  //- Edges, addressing into local point list
247  edgeList edges_;
248 
249  //- Edge owner
250  labelList edgeOwner_;
251 
252  //- Edge neighbour
253  labelList edgeNeighbour_;
254 
255 
256  // Primitive size data
257 
258  //- Total number of points
259  mutable label nPoints_;
260 
261  //- Total number of edges
262  mutable label nEdges_;
263 
264  //- Number of internal edges
265  mutable label nInternalEdges_;
266 
267  //- Number of faces
268  mutable label nFaces_;
269 
270 
271  // Communication support, updating
272 
273  //- Communicator used for parallel communication
274  label comm_;
275 
276  //- Current time index for motion
277  // Note. The whole mechanism will be replaced once the
278  // dimensionedField is created and the dimensionedField
279  // will take care of the old-time levels.
280  mutable label curTimeIndex_;
281 
282 
283  // Demand-driven data
284 
285  //- Primitive patch
286  mutable std::unique_ptr<uindirectPrimitivePatch> patchPtr_;
287 
288  //- List of poly patch/patch-face for faceLabels
289  mutable std::unique_ptr<List<labelPair>> polyPatchFacesPtr_;
290 
291  //- List of polyPatch ids used by the area mesh
292  mutable std::unique_ptr<labelList> polyPatchIdsPtr_;
293 
294  //- List of proc/mesh-face for boundary edge neighbours
295  mutable std::unique_ptr<List<labelPair>> bndConnectPtr_;
296 
297  //- Ldu addressing data
298  mutable std::unique_ptr<faMeshLduAddressing> lduPtr_;
299 
300 
301  // Geometric Data
302 
303  //- Face areas
304  mutable std::unique_ptr<DimensionedField<scalar, areaMesh>> SPtr_;
305 
306  //- Face areas old time level
307  mutable std::unique_ptr<DimensionedField<scalar, areaMesh>> S0Ptr_;
308 
309  //- Face areas old-old time level
310  mutable std::unique_ptr<DimensionedField<scalar, areaMesh>> S00Ptr_;
311 
312  //- Patch starts in the edge list
313  mutable std::unique_ptr<labelList> patchStartsPtr_;
314 
315  //- Edge length vectors
316  mutable std::unique_ptr<edgeVectorField> LePtr_;
317 
318  //- Mag edge length vectors
319  mutable std::unique_ptr<edgeScalarField> magLePtr_;
320 
321  //- Face centres
322  mutable std::unique_ptr<areaVectorField> faceCentresPtr_;
323 
324  //- Edge centres
325  mutable std::unique_ptr<edgeVectorField> edgeCentresPtr_;
326 
327  //- Face area normals
328  mutable std::unique_ptr<areaVectorField> faceAreaNormalsPtr_;
329 
330  //- Edge area normals
331  mutable std::unique_ptr<edgeVectorField> edgeAreaNormalsPtr_;
332 
333  //- Point area normals
334  mutable std::unique_ptr<vectorField> pointAreaNormalsPtr_;
335 
336  //- Face curvatures
337  mutable std::unique_ptr<areaScalarField> faceCurvaturesPtr_;
338 
339  //- Edge transformation tensors
340  mutable std::unique_ptr<FieldField<Field, tensor>>
341  edgeTransformTensorsPtr_;
342 
343  //- Whether point normals must be corrected for a patch
344  mutable std::unique_ptr<boolList> correctPatchPointNormalsPtr_;
345 
346 
347  // Other mesh-related data
348 
349  //- Parallel info
350  mutable autoPtr<faGlobalMeshData> globalMeshDataPtr_;
351 
352  //- Mapping/swapping for boundary edge halo neighbours
353  mutable std::unique_ptr<faMeshBoundaryHalo> haloMapPtr_;
354 
355  //- Face centres for boundary edge halo neighbours
356  mutable std::unique_ptr<pointField> haloFaceCentresPtr_;
357 
358  //- Face normals for boundary edge halo neighbours
359  mutable std::unique_ptr<vectorField> haloFaceNormalsPtr_;
360 
361 
362  // Static Private Data
363 
364  //- The prefix to local: %finite-area
365  static const word prefix_;
366 
367  //- Quadrics fit for pointAreaNormals (experimental)
368  static const int quadricsFit_;
369 
370 
371  // Private Member Functions
372 
373  //- Set indirect patch, removing any old one.
374  // No communication
375  void initPatch() const;
376 
377  //- Set primitive mesh data.
378  // No communication
379  void setPrimitiveMeshData();
380 
381  //- Get list of (proc/patchi/patchEdgei/meshFacei) tuple pairs in an
382  //- globally consistent ordering
383  List<Pair<patchTuple>> getBoundaryEdgeConnections() const;
384 
385  //- Determine the boundary edge neighbour connections
386  void calcBoundaryConnections() const;
387 
388  //- Define boundary edge neighbours (proc/face) based on
389  //- gathered topology information
390  void setBoundaryConnections
391  (
392  const List<Pair<patchTuple>>& bndEdgeConnections
393  ) const;
394 
395 
396  // Private Member Functions to calculate demand driven data
397 
398  //- Calculate ldu addressing
399  void calcLduAddressing() const;
400 
401  //- Calculate patch starts in the edge list
402  void calcPatchStarts() const;
403 
404  //- Calculate which polyPatches, polyPatch/local-face
405  //- are related to the areaMesh
406  void calcWhichPatchFaces() const;
407 
408  //- Calculate the edge normals (direct calculation)
409  //- with flat boundary addressing
410  // Possible communication
411  tmp<vectorField> calcRawEdgeNormals(int calcType) const;
412 
413  //- Calculate edge lengths
414  // Triggers communication via calcEdgeAreaNormals
415  void calcLe() const;
416 
417  //- Calculate mag edge lengths
418  // No communication
419  void calcMagLe() const;
420 
421  //- Calculate face centres
422  // No communication
423  void calcFaceCentres() const;
424 
425  //- Calculate edge centres
426  // No communication
427  void calcEdgeCentres() const;
428 
429  //- Calculate face areas
430  // No communication
431  void calcS() const;
432 
433  //- Calculate face area normals
434  // Triggers communication via calcEdgeAreaNormals
435  void calcFaceAreaNormals() const;
436 
437  //- Calculate edge area normals
438  // Triggers communication via pointAreaNormals
439  void calcEdgeAreaNormals() const;
440 
441  //- Calculate point area normals
442  // Triggers communication for processor patches
443  void calcPointAreaNormals(vectorField& result) const;
444 
445  //- Calculate point area normals by quadrics fit
446  void calcPointAreaNormalsByQuadricsFit(vectorField& result) const;
447 
448  //- Calculate face curvatures
449  // Triggers communication via edgeLengthCorrection
450  void calcFaceCurvatures() const;
451 
452  //- Calculate edge transformation tensors
453  void calcEdgeTransformTensors() const;
454 
455  //- Clear geometry but not the face areas
456  void clearGeomNotAreas() const;
457 
458  //- Has halo face centers/normals
459  bool hasHaloFaceGeometry() const noexcept;
460 
461  //- Clear boundary halo information
462  void clearHalo() const;
463 
464  //- Clear geometry
465  void clearGeom() const;
466 
467  //- Clear addressing
468  void clearAddressing() const;
469 
470  //- Clear demand-driven data
471  void clearOut() const;
472 
473 
474  // Halo handling
475 
476  //- Calculate halo centres/normals
477  void calcHaloFaceGeometry() const;
478 
479 
480  // Helpers
481 
482  //- Create a single patch
483  faPatchList createOnePatch
484  (
485  const word& patchName,
486  const word& patchType = ""
487  ) const;
488 
489  //- Create list of patches from boundary definition
490  faPatchList createPatchList
491  (
492  const dictionary& bndDict,
493  const word& emptyPatchName = "",
494  const dictionary* defaultPatchDefinition = nullptr
495  ) const;
496 
497 
498  //- Fatal error if edge labels are out of range
499  void checkBoundaryEdgeLabelRange(const labelUList& edgeLabels) const;
500 
501  //- Extract list from contiguous (unordered) boundary data
502  //- to the locally sorted order.
503  template<class T>
504  List<T> boundarySubset
505  (
506  const UList<T>& bndField,
507  const labelUList& edgeLabels
508  ) const
509  {
510  #ifdef FULLDEBUG
511  checkBoundaryEdgeLabelRange(edgeLabels);
512  #endif
513  // Like UIndirectList but with an offset
514  List<T> result(edgeLabels.size());
515  forAll(edgeLabels, i)
516  {
517  result[i] = bndField[edgeLabels[i] - nInternalEdges_];
518  }
519  return result;
520  }
521 
522 
523  // Static Functions
524 
525  //- Test if faSchemes/faSolution files are available
526  static bool hasSystemFiles
527  (
528  const word& meshName,
529  const polyMesh& pMesh
530  );
531 
532  //- Test if mesh files needed for read construction are available
533  static bool hasMeshFiles
534  (
535  const word& meshName,
536  const polyMesh& pMesh
537  );
538 
539 public:
540 
541  // Public Typedefs
542 
543  //- The mesh type
544  typedef faMesh Mesh;
545 
546  //- The boundary type associated with the mesh
548 
549 
550  // Tuning switches
551 
552  //- Geometry treatment
553  static int geometryOrder_;
554 
555 
556  //- Runtime type information
557  TypeName("faMesh");
558 
559  //- The prefix to the parent registry name: %finite-area
560  static const word& prefix() noexcept;
561 
562  //- The mesh sub-directory name (usually "faMesh")
563  static word meshSubDir;
564 
565 
566  // Generated Methods
567 
568  //- No copy construct
569  faMesh(const faMesh&) = delete;
570 
571  //- No copy assignment
572  void operator=(const faMesh&) = delete;
573 
574 
575  // Constructors
576 
577  //- Construct zero-sized from polyMesh
578  // Boundary is added using addFaPatches() member function
579  faMesh(const word& meshName, const polyMesh&, Foam::zero);
580 
581  //- Construct zero-sized from polyMesh
582  // Boundary is added using addFaPatches() member function
583  faMesh(const polyMesh&, Foam::zero);
584 
585  //- Read construct from polyMesh, using its IOobject properties
586  faMesh(const word& meshName, const polyMesh&, const bool doInit = true);
587 
588  //- Read construct from polyMesh, using its IOobject properties
589  explicit faMesh(const polyMesh&, const bool doInit = true);
590 
591  //- Construct from components (face labels) without boundary,
592  //- using specified read properties.
593  // Boundary is added using addFaPatches() member function.
594  faMesh
595  (
596  const word& meshName,
597  const polyMesh& pMesh,
599  IOobjectOption ioOpt
600  );
601 
602  //- Construct from components (face labels) without boundary,
603  //- using IOobject properties from polyMesh.
604  // Boundary is added using addFaPatches() member function.
605  faMesh(const word& meshName, const polyMesh&, labelList&& faceLabels);
606 
607  //- Construct from components (face labels) without boundary,
608  //- using specified read properties.
609  // Boundary is added using addFaPatches() member function.
611 
612  //- Construct from components (face labels) without boundary,
613  //- using IOobject properties from polyMesh.
614  // Boundary is added using addFaPatches() member function.
615  faMesh(const polyMesh&, labelList&& faceLabels);
616 
617  //- Construct as copy (for dictionaries) and zero-sized
618  //- without boundary.
619  // Boundary is added using addFaPatches() member function
620  faMesh(const word& meshName, const faMesh& baseMesh, Foam::zero);
621 
622  //- Construct as copy (for dictionaries) and zero-sized
623  //- without boundary, using IOobject properties from polyMesh.
624  // Boundary is added using addFaPatches() member function
625  faMesh(const faMesh& baseMesh, Foam::zero);
626 
627  //- Construct as copy (for dictionaries) and faceLabels
628  //- without boundary, using specified read properties.
629  // Boundary is added using addFaPatches() member function
630  faMesh
631  (
632  const word& meshName,
633  const faMesh& baseMesh,
635  IOobjectOption ioOpt
636  );
637 
638  //- Construct as copy (for dictionaries) and faceLabels
639  //- without boundary, using specified read properties.
640  // Boundary is added using addFaPatches() member function
641  faMesh
642  (
643  const faMesh& baseMesh,
645  IOobjectOption ioOpt
646  );
647 
648  //- Construct as copy (for dictionaries) and faceLabels
649  //- without boundary, using IOobject properties from polyMesh.
650  // Boundary is added using addFaPatches() member function
651  faMesh
652  (
653  const word& meshName,
654  const faMesh& baseMesh,
656  );
657 
658  //- Construct as copy (for dictionaries) and faceLabels
659  //- without boundary, using read properties from baseMesh.
660  // Boundary is added using addFaPatches() member function
661  faMesh(const faMesh& baseMesh, labelList&& faceLabels);
662 
663  //- Construct from single polyPatch
664  faMesh
665  (
666  const word& meshName,
667  const polyPatch& pp,
668  const bool doInit = true
669  );
670 
671  //- Construct from single polyPatch
672  explicit faMesh(const polyPatch& pp, const bool doInit = true);
673 
674  //- Construct from definition
675  faMesh
676  (
677  const word& meshName,
678  const polyMesh& pMesh,
679  const dictionary& faMeshDefinition,
680  const bool doInit = true
681  );
682 
683  //- Construct from definition
684  faMesh
685  (
686  const polyMesh& pMesh,
687  const dictionary& faMeshDefinition,
688  const bool doInit = true
689  );
690 
691 
692  //- Destructor
693  virtual ~faMesh();
694 
695 
696  // Static Functions
697 
698  //- Return the current geometry treatment
699  // A zero level or negative is with restricted neighbour information
700  static int geometryOrder() noexcept
701  {
702  return geometryOrder_;
703  }
704 
705  //- Set the preferred geometry treatment
706  // \return the previous value
707  static int geometryOrder(int order) noexcept
708  {
709  int old(geometryOrder_);
710  geometryOrder_ = order;
711  return old;
712  }
713 
714  //- Read construction from polyMesh if all files are available
715  static autoPtr<faMesh> TryNew
716  (
717  const word& meshName,
718  const polyMesh& pMesh
719  );
720 
721  //- Read construction from polyMesh if all files are available
722  static autoPtr<faMesh> TryNew(const polyMesh& pMesh);
723 
724 
725  // Member Functions
727  // Topological Change
728 
729  //- Add boundary patches. Constructor helper
730  void addFaPatches
731  (
732  faPatchList& plist,
733  const bool validBoundary = true
734  );
735 
736  //- Add boundary patches. Constructor helper
737  void addFaPatches
738  (
739  const List<faPatch*>& p,
740  const bool validBoundary = true
741  );
742 
743  //- Initialise non-demand-driven data etc
744  bool init(const bool doInit);
745 
746  //- Processor/processor synchronisation for geometry fields.
747  // Largely internal use only (slightly hacky).
748  void syncGeom();
749 
751  // Database
752 
753  //- The parent registry containing all finite-area meshes
754  //- on the polyMesh.
755  static const objectRegistry* registry(const polyMesh& pMesh);
756 
757  //- The single-region finite-area region on the polyMesh.
758  //- Uses lookupObject semantics - Fatal if non-existent
759  static const faMesh& mesh(const polyMesh& pMesh);
760 
761  //- The single-region or specified finite-area region on the polyMesh.
762  //- Uses lookupObject semantics - Fatal if non-existent
763  static const faMesh& mesh(const polyMesh&, const word& areaRegion);
764 
765  //- Return access to polyMesh
766  const polyMesh& mesh() const;
767 
768  //- Return the local mesh directory (dbDir()/meshSubDir)
769  fileName meshDir() const;
770 
771  //- Return reference to time
772  const Time& time() const;
773 
774  //- Return the current instance directory for points
775  // Used in the construction of geometric mesh data dependent
776  // on points
777  const fileName& pointsInstance() const;
778 
779  //- Return the current instance directory for faces
780  const fileName& facesInstance() const;
781 
782  //- Const reference to the mesh and solver state data
783  const meshState& data() const { return mesh().data(); }
784 
785  //- Reference to the mesh and solver state data
786  meshState& data() { return const_cast<polyMesh&>(mesh()).data(); }
787 
788 
789  // Parallel
790 
791  //- Return communicator used for parallel communication
792  label comm() const noexcept { return comm_; }
793 
794  //- Return communicator used for parallel communication
795  label& comm() noexcept { return comm_; }
796 
797  //- Is demand-driven parallel info available?
798  bool hasGlobalData() const noexcept;
799 
800  //- Return parallel info (demand-driven)
801  const faGlobalMeshData& globalData() const;
802 
803 
804  // Solution Control
805 
806  //- Non-null if faSchemes exists (can test as bool).
807  const faSchemes* hasSchemes() const;
808 
809  //- Non-null if faSolution exists (can test as bool).
810  const faSolution* hasSolution() const;
811 
812  //- Read-access to the faSchemes controls
813  const faSchemes& schemes() const;
814 
815  //- Read/write-access to the faSchemes controls
816  faSchemes& schemes();
817 
818  //- Read-access to the faSolution controls
819  const faSolution& solution() const;
820 
821  //- Read/write-access to the faSolution controls
822  faSolution& solution();
823 
824 
825  // Access: Mesh size parameters
826 
827  //- Number of local mesh points
828  inline label nPoints() const noexcept;
829 
830  //- Number of local mesh edges
831  inline label nEdges() const noexcept;
832 
833  //- Number of internal faces
834  inline label nInternalEdges() const noexcept;
835 
836  //- Number of boundary edges (== nEdges - nInternalEdges)
837  inline label nBoundaryEdges() const noexcept;
838 
839  //- Number of patch faces
840  inline label nFaces() const noexcept;
841 
842 
843  // Access: Primitive mesh data
844 
845  //- Return local points
846  inline const pointField& points() const;
847 
848  //- Return local edges with reordered boundary
849  inline const edgeList& edges() const noexcept;
850 
851  //- Sub-list of local internal edges
852  inline const edgeList::subList internalEdges() const;
853 
854  //- Return local faces
855  inline const faceList& faces() const;
856 
857  //- Edge owner addressing
858  inline const labelList& edgeOwner() const noexcept;
859 
860  //- Edge neighbour addressing
861  inline const labelList& edgeNeighbour() const noexcept;
862 
863  //- True if the internalEdges use an ordering that does not
864  //- correspond 1-to-1 with the patch internalEdges.
865  inline bool hasInternalEdgeLabels() const noexcept;
866 
867 
868  // Database
869 
870  //- True - thisDb() is a valid registry
871  virtual bool hasDb() const { return true; }
872 
873  //- Reference to the mesh database
874  virtual const objectRegistry& thisDb() const
875  {
876  return faMeshRegistry::thisDb();
877  }
878 
879  //- Local directory path of the objectRegistry relative to Time
880  //- with override for the single-region case
881  virtual const fileName& dbDir() const
882  {
883  return faMeshRegistry::dbDir();
884  }
885 
886  //- Name function is needed to disambiguate those inherited
887  //- from base classes
888  const word& name() const
889  {
890  return faMeshRegistry::thisDb().name();
891  }
892 
893 
894  // Regions
895 
896  //- Local registry directory path (relative to Time) for specified
897  //- area mesh (of a single-region volume mesh)
898  static fileName dbDir(const word& areaRegion);
899 
900  //- Local registry directory path (relative to Time) for specified
901  //- volume mesh and area mesh combination
902  static fileName dbDir(const word& volRegion, const word& areaRegion);
903 
904  //- Local registry directory path (relative to Time) for specified
905  //- volume mesh and area mesh combination
906  static fileName dbDir
907  (
908  const polyMesh& pMesh,
909  const word& areaRegion = word::null
910  );
911 
912  //- The local mesh directory name (eg, "faMesh") for specified
913  //- area mesh (of a single-region volume mesh)
914  static fileName meshDir(const word& areaRegion);
915 
916  //- The local mesh directory name (eg, "faMesh") for specified
917  //- volume mesh and area mesh combination
918  static fileName meshDir(const word& volRegion, const word& areaRegion);
919 
920  //- The local mesh directory name (eg, "faMesh") for specified
921  //- volume mesh and area mesh combination
922  static fileName meshDir
923  (
924  const polyMesh& pMesh,
925  const word& areaRegion = word::null
926  );
927 
928  //- The mesh region name or word::null if polyMesh::defaultRegion
929  const word& regionName() const;
930 
931 
932  // Access
933 
934  //- Return constant reference to boundary mesh
935  inline const faBoundaryMesh& boundary() const noexcept;
936 
937  //- Return the underlying polyMesh face labels
938  inline const labelList& faceLabels() const noexcept;
939 
940  //- The area-face corresponding to the mesh-face, -1 if not found
941  inline label whichFace(const label meshFacei) const;
942 
943  //- The polyPatches related to the areaMesh, in sorted order
944  inline const labelList& whichPolyPatches() const;
945 
946  //- The polyPatch/local-face for each faceLabels()
947  inline const List<labelPair>& whichPatchFaces() const;
948 
949  //- Return ldu addressing
950  virtual const lduAddressing& lduAddr() const;
951 
952  //- Return a list of pointers for each patch
953  // with only those pointing to interfaces being set
955  {
956  return boundary().interfaces();
957  }
958 
959  //- Internal face owner
960  const labelUList& owner() const
961  {
962  return lduAddr().lowerAddr();
963  }
964 
965  //- Internal face neighbour
966  const labelUList& neighbour() const
967  {
968  return lduAddr().upperAddr();
969  }
970 
971  //- True if given edge label is internal to the mesh
972  bool isInternalEdge(const label edgeIndex) const noexcept
973  {
974  return (edgeIndex < nInternalEdges_);
975  }
976 
977  //- List of proc/face for the boundary edge neighbours
978  //- using primitive patch edge numbering.
979  inline const List<labelPair>& boundaryConnections() const;
980 
981  //- Boundary edge neighbour processors
982  //- (does not include own proc)
983  labelList boundaryProcs() const;
984 
985  //- List of proc/size for the boundary edge neighbour processors
986  //- (does not include own proc)
987  List<labelPair> boundaryProcSizes() const;
988 
989  //- Mapping/swapping for boundary halo neighbours
990  const faMeshBoundaryHalo& boundaryHaloMap() const;
991 
992  //- Face centres of boundary halo neighbours
993  const pointField& haloFaceCentres() const;
994 
995  //- Face unit-normals of boundary halo neighbours
996  const vectorField& haloFaceNormals() const;
997 
998  //- Face centres of boundary halo neighbours for specified patch
999  tmp<pointField> haloFaceCentres(const label patchi) const;
1000 
1001  //- Face unit-normals of boundary halo neighbours for specified patch
1002  tmp<vectorField> haloFaceNormals(const label patchi) const;
1003 
1004 
1005  // Interfacing
1006 
1007  //- The volume (owner) cells associated with the area-mesh
1008  labelList faceCells() const;
1009 
1010 
1011  // Storage Management
1012 
1013  //- Remove all files from mesh instance
1014  void removeFiles(const fileName& instanceDir) const;
1015 
1016  //- Remove all files from mesh instance()
1017  void removeFiles() const;
1018 
1019 
1020  //- Has face areas: S()
1021  bool hasFaceAreas() const noexcept { return bool(SPtr_); }
1022 
1023  //- Has face centres: areaCentres()
1024  bool hasAreaCentres() const noexcept { return bool(faceCentresPtr_); }
1025 
1026  //- Has edge centres: edgeCentres()
1027  bool hasEdgeCentres() const noexcept { return bool(edgeCentresPtr_); }
1028 
1029  //- Has edge length vectors: Le()
1030  bool hasLe() const noexcept { return bool(LePtr_); }
1031 
1032  //- Has edge length magnitudes: magLe()
1033  bool hasMagLe() const noexcept { return bool(magLePtr_); }
1034 
1035  //- Has face area normals: faceAreaNormals()
1036  bool hasFaceAreaNormals() const noexcept
1037  {
1038  return bool(faceAreaNormalsPtr_);
1039  }
1040 
1041  //- Has edge area normals: edgeAreaNormals()
1042  bool hasEdgeAreaNormals() const noexcept
1043  {
1044  return bool(edgeAreaNormalsPtr_);
1045  }
1046 
1047  //- Has point area normals: pointAreaNormals()
1048  bool hasPointAreaNormals() const noexcept
1049  {
1050  return bool(pointAreaNormalsPtr_);
1051  }
1052 
1053  //- Has face curvatures: faceCurvatures()
1054  bool hasFaceCurvatures() const noexcept
1055  {
1056  return bool(faceCurvaturesPtr_);
1057  }
1058 
1059  //- Has patch point normals corrections
1061  {
1062  return bool(correctPatchPointNormalsPtr_);
1063  }
1064 
1065 
1066  // Mesh motion and morphing
1068  //- Is mesh moving
1069  bool moving() const
1070  {
1071  return mesh().moving();
1072  }
1073 
1074  //- Update after mesh motion
1075  virtual bool movePoints();
1076 
1077  //- Update after topo change
1078  virtual void updateMesh(const mapPolyMesh&);
1079 
1081  // Mapping
1082 
1083  //- Map all fields in time using given map.
1084  virtual void mapFields(const faMeshMapper& mapper) const;
1085 
1086  //- Map face areas in time using given map.
1087  virtual void mapOldAreas(const faMeshMapper& mapper) const;
1088 
1089 
1090  // Demand-driven data
1091 
1092  //- Return constant reference to primitive patch
1093  inline const uindirectPrimitivePatch& patch() const;
1094 
1095  //- Return reference to primitive patch
1096  inline uindirectPrimitivePatch& patch();
1097 
1098  //- Return patch starts
1099  const labelList& patchStarts() const;
1100 
1101  //- Return edge length vectors
1102  const edgeVectorField& Le() const;
1103 
1104  //- Return edge length magnitudes
1105  const edgeScalarField& magLe() const;
1106 
1107  //- Return normalised edge length vectors
1108  tmp<edgeVectorField> unitLe() const;
1109 
1110  //- Return face centres as areaVectorField
1111  const areaVectorField& areaCentres() const;
1112 
1113  //- Return edge centres as edgeVectorField
1114  const edgeVectorField& edgeCentres() const;
1115 
1116  //- Return face areas
1117  const DimensionedField<scalar, areaMesh>& S() const;
1118 
1119  //- Return old-time face areas
1120  const DimensionedField<scalar, areaMesh>& S0() const;
1121 
1122  //- Return old-old-time face areas
1123  const DimensionedField<scalar, areaMesh>& S00() const;
1124 
1125  //- Return face area normals
1126  const areaVectorField& faceAreaNormals() const;
1127 
1128  //- Return edge area normals
1129  const edgeVectorField& edgeAreaNormals() const;
1130 
1131  //- Return point area normals
1132  const vectorField& pointAreaNormals() const;
1133 
1134  //- Return face curvatures
1135  const areaScalarField& faceCurvatures() const;
1136 
1137  //- Return edge transformation tensors
1139 
1140  //- Return internal point labels
1141  labelList internalPoints() const;
1142 
1143  //- Return boundary point labels
1144  labelList boundaryPoints() const;
1145 
1146  //- Return edge length correction
1148 
1149  //- Whether point normals should be corrected for a patch
1150  bool correctPatchPointNormals(const label patchID) const;
1151 
1152  //- Set whether point normals should be corrected for a patch
1154 
1155 
1156  // Storage management
1157 
1158 
1159  //- Write mesh
1160  virtual bool write(const bool writeOnProc = true) const;
1161 
1162 
1163  // Member Operators
1164 
1165  bool operator!=(const faMesh& m) const;
1166 
1167  bool operator==(const faMesh& m) const;
1168 
1169 
1170  // Housekeeping
1171 
1172  //- No call operator. Prior to 2312 was used to obtain polyMesh
1173  void operator()() const = delete;
1174 };
1175 
1176 
1177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1178 
1179 } // End namespace Foam
1180 
1181 
1182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1183 
1184 #include "faMeshI.H"
1185 
1186 #ifdef NoRepository
1187  #include "faPatchFaMeshTemplates.C"
1188 #endif
1189 
1190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1191 
1192 #endif
1193 
1194 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:133
labelList internalPoints() const
Return internal point labels.
faMesh(const faMesh &)=delete
No copy construct.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual void updateMesh(const mapPolyMesh &)
Update after topo change.
Definition: faMeshUpdate.C:31
void addFaPatches(faPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: faMeshPatches.C:68
A class for handling file names.
Definition: fileName.H:72
faSolution(const faSolution &)=delete
No copy construct.
tmp< edgeScalarField > edgeLengthCorrection() const
Return edge length correction.
const DimensionedField< scalar, areaMesh > & S00() const
Return old-old-time face areas.
Definition: faMesh.C:1159
virtual bool hasDb() const
True - thisDb() is a valid registry.
Definition: faMesh.H:114
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch with only those pointing to interfaces being set...
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
const Time & time() const
Return reference to time.
Definition: faMesh.C:1022
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
void syncGeom()
Processor/processor synchronisation for geometry fields.
Definition: faMesh.C:396
virtual const fileName & dbDir() const
Local directory path of the objectRegistry relative to Time with override for the single-region case...
virtual void mapFields(const faMeshMapper &mapper) const
Map all fields in time using given map.
Definition: faMeshUpdate.C:144
label comm() const noexcept
Return communicator used for parallel communication.
Definition: faMesh.H:1075
const edgeList::subList internalEdges() const
Sub-list of local internal edges.
Definition: faMeshI.H:72
const pointField & haloFaceCentres() const
Face centres of boundary halo neighbours.
The objectRegistry for faMesh.
Definition: faMesh.H:93
Database for mesh data, solution data, solver performance and other reduced data. ...
Definition: meshState.H:53
const DimensionedField< scalar, areaMesh > & S() const
Return face areas.
Definition: faMesh.C:1133
Face to edge interpolation scheme. Included in faMesh.
bool operator!=(const faMesh &m) const
Definition: faMesh.C:1374
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:53
label whichFace(const label meshFacei) const
The area-face corresponding to the mesh-face, -1 if not found.
Definition: faMeshI.H:122
virtual const meshState & data() const noexcept
Const reference to the mesh and solver state data.
Definition: polyMesh.H:559
void operator()() const =delete
No call operator. Prior to 2312 was used to obtain polyMesh.
void stableSort(UList< T > &list)
Stable sort the list.
Definition: UList.C:312
const labelList & whichPolyPatches() const
The polyPatches related to the areaMesh, in sorted order.
Definition: faMeshI.H:128
Generic GeometricField class.
const faGlobalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: faMesh.C:1252
const labelList & patchStarts() const
Return patch starts.
Definition: faMesh.C:1077
bool isInternalEdge(const label edgeIndex) const noexcept
True if given edge label is internal to the mesh.
Definition: faMesh.H:1340
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1086
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: faMesh.C:1040
const List< labelPair > & boundaryConnections() const
List of proc/face for the boundary edge neighbours using primitive patch edge numbering.
Definition: faMeshI.H:155
const areaVectorField & areaCentres() const
Return face centres as areaVectorField.
Definition: faMesh.C:1110
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
labelList boundaryProcs() const
Boundary edge neighbour processors (does not include own proc)
bool hasMagLe() const noexcept
Has edge length magnitudes: magLe()
Definition: faMesh.H:1433
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
const edgeScalarField & magLe() const
Return edge length magnitudes.
Definition: faMesh.C:1099
bool hasGlobalData() const noexcept
Is demand-driven parallel info available?
Definition: faMesh.C:1246
label nFaces() const noexcept
Number of patch faces.
Definition: faMeshI.H:54
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
label nPoints() const noexcept
Number of local mesh points.
Definition: faMeshI.H:30
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
virtual const objectRegistry & thisDb() const
Reference to the mesh database.
Definition: faMesh.H:1203
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const uindirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition: faMeshI.H:102
void operator=(const faMesh &)=delete
No copy assignment.
static int geometryOrder() noexcept
Return the current geometry treatment.
Definition: faMesh.H:944
virtual const fileName & dbDir() const
Local directory path of the objectRegistry relative to Time with override for the single-region case...
Definition: faMesh.H:1212
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
const labelUList & neighbour() const
Internal face neighbour.
Definition: faMesh.H:1332
const faceList & faces() const
Return local faces.
Definition: faMeshI.H:78
A list of faces which address into the list of points.
List< labelPair > boundaryProcSizes() const
List of proc/size for the boundary edge neighbour processors (does not include own proc) ...
bool hasLe() const noexcept
Has edge length vectors: Le()
Definition: faMesh.H:1428
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
Forwards for edge field types.
Forward declarations of the specialisations of Field<T> for scalar, vector and tensor.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: instant.H:46
dynamicFvMesh & mesh
const labelList & faceLabels() const noexcept
Return the underlying polyMesh face labels.
Definition: faMeshI.H:96
virtual ~faMesh()
Destructor.
Definition: faMesh.C:966
bool hasEdgeAreaNormals() const noexcept
Has edge area normals: edgeAreaNormals()
Definition: faMesh.H:1446
A class for handling words, derived from Foam::string.
Definition: word.H:63
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: faMesh.C:1016
bool operator==(const faMesh &m) const
Definition: faMesh.C:1380
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: faMesh.C:1263
virtual const labelUList & upperAddr() const =0
Return upper addressing.
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:296
const edgeVectorField & Le() const
Return edge length vectors.
Definition: faMesh.C:1088
bool hasFaceCurvatures() const noexcept
Has face curvatures: faceCurvatures()
Definition: faMesh.H:1462
const faMeshBoundaryHalo & boundaryHaloMap() const
Mapping/swapping for boundary halo neighbours.
const meshState & data() const
Const reference to the mesh and solver state data.
Definition: faMesh.H:1062
faMeshRegistry(const word &areaRegion, const polyMesh &mesh)
Construct an objectRegistry for given area region name.
bool hasPointAreaNormals() const noexcept
Has point area normals: pointAreaNormals()
Definition: faMesh.H:1454
bool hasAreaCentres() const noexcept
Has face centres: areaCentres()
Definition: faMesh.H:1418
bool hasInternalEdgeLabels() const noexcept
True if the internalEdges use an ordering that does not correspond 1-to-1 with the patch internalEdge...
Definition: faMeshI.H:148
static const word null
An empty word.
Definition: word.H:84
const pointField & points() const
Return local points.
Definition: faMeshI.H:60
const vectorField & haloFaceNormals() const
Face unit-normals of boundary halo neighbours.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
virtual bool hasDb() const
True - thisDb() is a valid registry.
Definition: faMesh.H:1198
const FieldField< Field, tensor > & edgeTransformTensors() const
Return edge transformation tensors.
Definition: faMesh.C:1235
virtual void mapOldAreas(const faMeshMapper &mapper) const
Map face areas in time using given map.
Definition: faMeshUpdate.C:166
const direction noexcept
Definition: Scalar.H:258
virtual bool movePoints()
Update after mesh motion.
Definition: faMesh.C:1274
virtual bool write(const bool writeOnProc=true) const
Write mesh.
Definition: faMesh.C:1363
const edgeList & edges() const noexcept
Return local edges with reordered boundary.
Definition: faMeshI.H:66
faBoundaryMesh BoundaryMesh
The boundary type associated with the mesh.
Definition: faMesh.H:726
const polyMesh & mesh() const
Return access to polyMesh.
Definition: faMesh.C:1010
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:24
const areaScalarField & faceCurvatures() const
Return face curvatures.
Definition: faMesh.C:1223
faSchemes(const faSchemes &)=delete
No copy construct.
label nEdges() const noexcept
Number of local mesh edges.
Definition: faMeshI.H:36
boolList & correctPatchPointNormals() const
Set whether point normals should be corrected for a patch.
Definition: faMesh.C:1351
const labelUList & owner() const
Internal face owner.
Definition: faMesh.H:1324
bool moving() const noexcept
Is mesh moving.
Definition: polyMesh.H:732
labelList boundaryPoints() const
Return boundary point labels.
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:750
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: faMesh.C:1034
label nInternalEdges() const noexcept
Number of internal faces.
Definition: faMeshI.H:42
Selector class for finite area differencing schemes. faMesh is derived from faSchemes so that all fie...
Definition: faSchemes.H:51
const edgeVectorField & edgeCentres() const
Return edge centres as edgeVectorField.
Definition: faMesh.C:1121
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
static autoPtr< faMesh > TryNew(const word &meshName, const polyMesh &pMesh)
Read construction from polyMesh if all files are available.
Definition: faMeshNew.C:188
bool moving() const
Is mesh moving.
Definition: faMesh.H:1481
bool hasEdgeCentres() const noexcept
Has edge centres: edgeCentres()
Definition: faMesh.H:1423
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const faSchemes * hasSchemes() const
Non-null if faSchemes exists (can test as bool).
Definition: faMesh.C:974
bool hasFaceAreas() const noexcept
Has face areas: S()
Definition: faMesh.H:1413
const labelList & edgeNeighbour() const noexcept
Edge neighbour addressing.
Definition: faMeshI.H:90
void clear()
Clear all entries from the registry.
const word & name() const
Name function is needed to disambiguate those inherited from base classes.
Definition: faMesh.H:1221
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: faMesh.C:1028
Finite area boundary mesh.
Selector class for finite area solution. faMesh is derived from faSolution so that all fields have ac...
Definition: faSolution.H:51
const edgeVectorField & edgeAreaNormals() const
Return edge area normals.
Definition: faMesh.C:1194
bool init(const bool doInit)
Initialise non-demand-driven data etc.
Definition: faMesh.C:417
constexpr label labelMax
Definition: label.H:55
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
labelList faceCells() const
The volume (owner) cells associated with the area-mesh.
Definition: faMesh.C:1046
Field< vector > vectorField
Specialisation of Field<T> for vector.
static const objectRegistry * registry(const polyMesh &pMesh)
The parent registry containing all finite-area meshes on the polyMesh.
Definition: faMesh.C:145
static int geometryOrder_
Geometry treatment.
Definition: faMesh.H:734
faMesh Mesh
The mesh type.
Definition: faMesh.H:721
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
const List< labelPair > & whichPatchFaces() const
The polyPatch/local-face for each faceLabels()
Definition: faMeshI.H:138
virtual const objectRegistry & thisDb() const
Reference to the mesh database.
Definition: faMesh.H:119
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
TypeName("faMesh")
Runtime type information.
void removeFiles() const
Remove all files from mesh instance()
Definition: faMesh.C:1071
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
label nBoundaryEdges() const noexcept
Number of boundary edges (== nEdges - nInternalEdges)
Definition: faMeshI.H:48
bool hasFaceAreaNormals() const noexcept
Has face area normals: faceAreaNormals()
Definition: faMesh.H:1438
static const word & prefix() noexcept
The prefix to the parent registry name: finite-area.
Definition: faMesh.C:66
const vectorField & pointAreaNormals() const
Return point area normals.
Definition: faMesh.C:1205
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition: faMesh.C:1183
bool hasPatchPointNormalsCorrection() const noexcept
Has patch point normals corrections.
Definition: faMesh.H:1470
tmp< edgeVectorField > unitLe() const
Return normalised edge length vectors.
const faSolution & solution() const
Read-access to the faSolution controls.
Definition: faMesh.C:998
const faSolution * hasSolution() const
Non-null if faSolution exists (can test as bool).
Definition: faMesh.C:980
objectRegistry(const Time &db, const label initialCapacity=128)
Construct the time objectRegistry, with estimated table capacity (default: 128)
const labelList & edgeOwner() const noexcept
Edge owner addressing.
Definition: faMeshI.H:84
Forwards and collection of common area field types.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
const DimensionedField< scalar, areaMesh > & S0() const
Return old-time face areas.
Definition: faMesh.C:1145
const faSchemes & schemes() const
Read-access to the faSchemes controls.
Definition: faMesh.C:986
Class holds all the necessary information for mapping fields associated with faMesh.
Definition: faMeshMapper.H:64
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: faMesh.H:1316