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-2023 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 "MeshObject.H"
60 #include "polyMesh.H"
61 #include "lduMesh.H"
62 #include "faBoundaryMesh.H"
63 #include "edgeList.H"
64 #include "faceList.H"
65 #include "primitiveFieldsFwd.H"
66 #include "DimensionedField.H"
67 #include "areaFieldsFwd.H"
68 #include "edgeFieldsFwd.H"
69 #include "indirectPrimitivePatch.H"
70 #include "edgeInterpolation.H"
71 #include "labelIOList.H"
72 #include "FieldFields.H"
73 #include "faGlobalMeshData.H"
74 #include "faSchemes.H"
75 #include "faSolution.H"
76 
77 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79 namespace Foam
80 {
81 
82 // Forward Declarations
83 class faMeshBoundaryHalo;
84 class faMeshLduAddressing;
85 class faMeshMapper;
86 class faPatchData;
87 
88 /*---------------------------------------------------------------------------*\
89  Class faMesh Declaration
90 \*---------------------------------------------------------------------------*/
91 
92 class faMesh
93 :
94  public MeshObject<polyMesh, Foam::UpdateableMeshObject, faMesh>,
95  public lduMesh,
96  public faSchemes,
97  public edgeInterpolation, // may need input from faSchemes
98  public faSolution
99 {
100  // Private (internal) classes/structures
101 
102  //- A (proc, patchi, patchEdgei, meshFacei) tuple used internally
103  //- for managing patch/patch bookkeeping during construction.
104  // Finite-area patches are stored with negated indices (offset -2),
105  // which makes them readily identifiable and always sort before normal
106  // patches.
107  struct patchTuple
108  :
109  public FixedList<label, 4>
110  {
111 
112  // Constructors
113 
114  // Inherit constructors
116 
117  //- Default construct as 'invalid'
118  patchTuple()
119  {
120  clear();
121  }
122 
123 
124  // Static Member Functions
125 
126  // Globally consistent ordering
127  //
128  // 1. sort left/right as lower/higher processor connection
129  // 2. sort by proc/patch/patch index
130  static void sort(UList<Pair<patchTuple>>& list)
131  {
132  for (auto& tuples : list)
133  {
134  tuples.sort();
135  }
136  Foam::stableSort(list);
137  }
138 
139 
140  // Member Functions
141 
142  //- Reset to 'invalid'
143  void clear()
144  {
145  procNo(-1);
146  patchi(labelMax);
147  patchEdgei(-1);
148  meshFacei(-1);
149  }
150 
151  //- Valid if proc and patch (or patch edge) are non-negative
152  bool valid() const
153  {
154  return (procNo() >= 0 && patchi() != -1);
155  }
156 
157  // Processor is the first sort index
158  label procNo() const { return (*this)[0]; }
159  void procNo(label val) { (*this)[0] = val; }
160 
161  // PatchId is the second sort index (finiteArea patches are < -1)
162  label patchi() const { return (*this)[1]; }
163  void patchi(label val) { (*this)[1] = val; }
164 
165  // The patch edge index (on the finiteArea patch)
166  // is the third sort index
167  label patchEdgei() const { return (*this)[2]; }
168  void patchEdgei(label val) { (*this)[2] = val; }
169 
170  // The processor-local mesh face is the fourth sort index
171  label meshFacei() const { return (*this)[3]; }
172  void meshFacei(label val) { (*this)[3] = val; }
173 
174  //- Return the real patchId
175  label realPatchi() const
176  {
177  const label id = patchi();
178  return (id < -1 ? -(id + 2) : id);
179  }
180 
181  //- Set patchId as finiteArea
182  void faPatchi(label val)
183  {
184  patchi(-(val + 2));
185  }
186 
187  //- Considered to be finiteArea if (patchi < -1)
188  bool is_finiteArea() const noexcept
189  {
190  return (patchi() < -1);
191  }
192 
193  //- Considered to be processor local
194  bool is_localProc() const
195  {
196  return (procNo() == UPstream::myProcNo());
197  }
198  };
199 
200 
201  // Private Data
202 
203  //- Face labels
204  labelIOList faceLabels_;
205 
206  //- Boundary mesh
207  faBoundaryMesh boundary_;
208 
209 
210  // Primitive mesh data
211 
212  //- Edges, addressing into local point list
213  edgeList edges_;
214 
215  //- Edge owner
216  labelList edgeOwner_;
217 
218  //- Edge neighbour
219  labelList edgeNeighbour_;
220 
221 
222  // Primitive size data
223 
224  //- Total number of points
225  mutable label nPoints_;
226 
227  //- Total number of edges
228  mutable label nEdges_;
229 
230  //- Number of internal edges
231  mutable label nInternalEdges_;
232 
233  //- Number of faces
234  mutable label nFaces_;
235 
236 
237  // Communication support, updating
238 
239  //- Communicator used for parallel communication
240  label comm_;
241 
242  //- Current time index for motion
243  // Note. The whole mechanism will be replaced once the
244  // dimensionedField is created and the dimensionedField
245  // will take care of the old-time levels.
246  mutable label curTimeIndex_;
247 
248 
249  // Demand-driven data
250 
251  //- Primitive patch
252  mutable std::unique_ptr<uindirectPrimitivePatch> patchPtr_;
253 
254  //- List of poly patch/patch-face for faceLabels
255  mutable std::unique_ptr<List<labelPair>> polyPatchFacesPtr_;
256 
257  //- List of polyPatch ids used by the area mesh
258  mutable std::unique_ptr<labelList> polyPatchIdsPtr_;
259 
260  //- List of proc/mesh-face for boundary edge neighbours
261  mutable std::unique_ptr<List<labelPair>> bndConnectPtr_;
262 
263  //- Ldu addressing data
264  mutable faMeshLduAddressing* lduPtr_;
265 
266 
267  // Geometric Data
268 
269  //- Face areas
270  mutable DimensionedField<scalar, areaMesh>* SPtr_;
271 
272  //- Face areas old time level
273  mutable DimensionedField<scalar, areaMesh>* S0Ptr_;
274 
275  //- Face areas old-old time level
276  mutable DimensionedField<scalar, areaMesh>* S00Ptr_;
277 
278  //- Patch starts in the edge list
279  mutable labelList* patchStartsPtr_;
280 
281  //- Edge length vectors
282  mutable edgeVectorField* LePtr_;
283 
284  //- Mag edge length vectors
285  mutable edgeScalarField* magLePtr_;
286 
287  //- Face centres
288  mutable areaVectorField* faceCentresPtr_;
289 
290  //- Edge centres
291  mutable edgeVectorField* edgeCentresPtr_;
292 
293  //- Face area normals
294  mutable areaVectorField* faceAreaNormalsPtr_;
295 
296  //- Edge area normals
297  mutable edgeVectorField* edgeAreaNormalsPtr_;
298 
299  //- Point area normals
300  mutable std::unique_ptr<vectorField> pointAreaNormalsPtr_;
301 
302  //- Face curvatures
303  mutable areaScalarField* faceCurvaturesPtr_;
304 
305  //- Edge transformation tensors
306  mutable FieldField<Field, tensor>* edgeTransformTensorsPtr_;
307 
308  //- Whether point normals must be corrected for a patch
309  mutable boolList* correctPatchPointNormalsPtr_;
310 
311 
312  // Other mesh-related data
313 
314  //- Parallel info
315  mutable autoPtr<faGlobalMeshData> globalMeshDataPtr_;
316 
317  //- Mapping/swapping for boundary edge halo neighbours
318  mutable std::unique_ptr<faMeshBoundaryHalo> haloMapPtr_;
319 
320  //- Face centres for boundary edge halo neighbours
321  mutable std::unique_ptr<pointField> haloFaceCentresPtr_;
322 
323  //- Face normals for boundary edge halo neighbours
324  mutable std::unique_ptr<vectorField> haloFaceNormalsPtr_;
325 
326 
327  // Static Private Data
328 
329  //- Quadrics fit for pointAreaNormals (experimental)
330  static const int quadricsFit_;
331 
332 
333  // Private Member Functions
334 
335  //- No copy construct
336  faMesh(const faMesh&) = delete;
337 
338  //- No copy assignment
339  void operator=(const faMesh&) = delete;
340 
341  //- Set indirect patch, removing any old one.
342  // No communication
343  void initPatch() const;
344 
345  //- Set primitive mesh data.
346  // No communication
347  void setPrimitiveMeshData();
348 
349  //- Get list of (proc/patchi/patchEdgei/meshFacei) tuple pairs in an
350  //- globally consistent ordering
351  List<Pair<patchTuple>> getBoundaryEdgeConnections() const;
352 
353  //- Determine the boundary edge neighbour connections
354  void calcBoundaryConnections() const;
355 
356  //- Define boundary edge neighbours (proc/face) based on
357  //- gathered topology information
358  void setBoundaryConnections
359  (
360  const List<Pair<patchTuple>>& bndEdgeConnections
361  ) const;
362 
363 
364  // Private member functions to calculate demand driven data
365 
366  //- Calculate ldu addressing
367  void calcLduAddressing() const;
368 
369  //- Calculate patch starts in the edge list
370  void calcPatchStarts() const;
371 
372  //- Calculate which polyPatches, polyPatch/local-face
373  //- are related to the areaMesh
374  void calcWhichPatchFaces() const;
375 
376  //- Calculate the edge normals (direct calculation)
377  //- with flat boundary addressing
378  // Possible communication
379  tmp<vectorField> calcRawEdgeNormals(int calcType) const;
380 
381  //- Calculate edge lengths
382  // Triggers communication via calcEdgeAreaNormals
383  void calcLe() const;
384 
385  //- Calculate mag edge lengths
386  // No communication
387  void calcMagLe() const;
388 
389  //- Calculate face centres
390  // No communication
391  void calcFaceCentres() const;
392 
393  //- Calculate edge centres
394  // No communication
395  void calcEdgeCentres() const;
396 
397  //- Calculate face areas
398  // No communication
399  void calcS() const;
400 
401  //- Calculate face area normals
402  // Triggers communication via calcEdgeAreaNormals
403  void calcFaceAreaNormals() const;
404 
405  //- Calculate edge area normals
406  // Triggers communication via pointAreaNormals
407  void calcEdgeAreaNormals() const;
408 
409  //- Calculate point area normals
410  // Triggers communication for processor patches
411  void calcPointAreaNormals(vectorField& result) const;
412 
413  //- Calculate point area normals by quadrics fit
414  void calcPointAreaNormalsByQuadricsFit(vectorField& result) const;
415 
416  //- Calculate face curvatures
417  // Triggers communication via edgeLengthCorrection
418  void calcFaceCurvatures() const;
419 
420  //- Calculate edge transformation tensors
421  void calcEdgeTransformTensors() const;
422 
423  //- Clear geometry but not the face areas
424  void clearGeomNotAreas() const;
425 
426  //- Has halo face centers/normals
427  bool hasHaloFaceGeometry() const noexcept;
428 
429  //- Clear boundary halo information
430  void clearHalo() const;
431 
432  //- Clear geometry
433  void clearGeom() const;
434 
435  //- Clear addressing
436  void clearAddressing() const;
437 
438  //- Clear demand-driven data
439  void clearOut() const;
440 
441 
442  // Halo handling
443 
444  //- Calculate halo centres/normals
445  void calcHaloFaceGeometry() const;
446 
447 
448  // Helpers
449 
450  //- Create a single patch
451  faPatchList createOnePatch
452  (
453  const word& patchName,
454  const word& patchType = ""
455  ) const;
456 
457  //- Create list of patches from boundary definition
458  faPatchList createPatchList
459  (
460  const dictionary& bndDict,
461  const word& emptyPatchName = "",
462  const dictionary* defaultPatchDefinition = nullptr
463  ) const;
464 
465 
466  //- Fatal error if edge labels are out of range
467  void checkBoundaryEdgeLabelRange(const labelUList& edgeLabels) const;
468 
469  //- Extract list from contiguous (unordered) boundary data
470  //- to the locally sorted order.
471  template<class T>
472  List<T> boundarySubset
473  (
474  const UList<T>& bndField,
475  const labelUList& edgeLabels
476  ) const
477  {
478  #ifdef FULLDEBUG
479  checkBoundaryEdgeLabelRange(edgeLabels);
480  #endif
481  // Like UIndirectList but with an offset
482  List<T> result(edgeLabels.size());
483  forAll(edgeLabels, i)
484  {
485  result[i] = bndField[edgeLabels[i] - nInternalEdges_];
486  }
487  return result;
488  }
489 
490 
491  // Static Functions
492 
493  //- Test if faSchemes/faSolution files are available
494  static bool hasSystemFiles(const polyMesh& pMesh);
495 
496  //- Test if all files needed for read construction are available
497  static bool hasFiles(const polyMesh& pMesh);
498 
499 
500 public:
501 
502  // Public Typedefs
503 
504  //- The mesh type
505  typedef faMesh Mesh;
506 
507  //- The boundary type associated with the mesh
509 
510 
511  // Tuning switches
512 
513  //- Geometry treatment
514  static int geometryOrder_;
515 
516 
517  //- Runtime type information
518  TypeName("faMesh");
519 
520  //- The prefix to local: %finite-area
521  static const word prefix;
522 
523  //- The mesh sub-directory name (usually "faMesh")
524  static word meshSubDir;
525 
526 
527  // Constructors
528 
529  //- Read construct from polyMesh, using its IOobject properties
530  explicit faMesh(const polyMesh& pMesh, const bool doInit = true);
531 
532  //- Construct zero-sized from polyMesh
533  // Boundary is added using addFaPatches() member function
534  faMesh(const polyMesh& pMesh, const Foam::zero);
535 
536  //- Construct as copy (for dictionaries) and zero-sized
537  //- without boundary.
538  // Boundary is added using addFaPatches() member function
539  faMesh(const faMesh& baseMesh, const Foam::zero);
540 
541  //- Construct as copy (for dictionaries) and faceLabels
542  //- without boundary, using read properties from baseMesh.
543  // Boundary is added using addFaPatches() member function
544  faMesh(const faMesh& baseMesh, labelList&& faceLabels);
545 
546  //- Construct as copy (for dictionaries) and faceLabels
547  //- without boundary, using specified read properties.
548  // Boundary is added using addFaPatches() member function.
549  faMesh
550  (
551  const faMesh& baseMesh,
553  IOobjectOption ioOpt
554  );
555 
556  //- Construct from components (face labels) without boundary,
557  //- using IOobject properties from polyMesh.
558  // Boundary is added using addFaPatches() member function.
559  faMesh(const polyMesh& pMesh, labelList&& faceLabels);
560 
561  //- Construct from components (face labels) without boundary,
562  //- using specified read properties.
563  // Boundary is added using addFaPatches() member function.
564  faMesh
565  (
566  const polyMesh& pMesh,
568  IOobjectOption ioOpt
569  );
570 
571  //- Construct from single polyPatch
572  explicit faMesh(const polyPatch& pp, const bool doInit = true);
573 
574  //- Construct from definition
575  faMesh
576  (
577  const polyMesh& pMesh,
578  const dictionary& faMeshDefinition,
579  const bool doInit = true
580  );
581 
582 
583  //- Destructor
584  virtual ~faMesh();
585 
586 
587  // Static Functions
588 
589  //- Return the current geometry treatment
590  // A zero level or negative is with restricted neighbour information
591  static int geometryOrder() noexcept
592  {
593  return geometryOrder_;
594  }
595 
596  //- Set the preferred geometry treatment
597  // \return the previous value
598  static int geometryOrder(int order) noexcept
599  {
600  int old(geometryOrder_);
601  geometryOrder_ = order;
602  return old;
603  }
604 
605  //- Read construction from polyMesh if all files are available
606  static autoPtr<faMesh> TryNew(const polyMesh& pMesh);
607 
608 
609  // Member Functions
610 
611  // Topological Change
612 
613  //- Add boundary patches. Constructor helper
614  void addFaPatches
615  (
616  faPatchList& plist,
617  const bool validBoundary = true
618  );
619 
620  //- Add boundary patches. Constructor helper
621  void addFaPatches
622  (
623  const List<faPatch*>& p,
624  const bool validBoundary = true
625  );
626 
627  //- Initialise non-demand-driven data etc
628  bool init(const bool doInit);
629 
630  //- Processor/processor synchronisation for geometry fields.
631  // Largely internal use only (slightly hacky).
632  void syncGeom();
633 
634 
635  // Database
636 
637  //- The parent registry containing all finite-area meshes
638  //- on the polyMesh.
639  static const objectRegistry* registry(const polyMesh& pMesh);
640 
641  //- The single-region finite-area region on the polyMesh.
642  //- Uses lookupObject semantics - Fatal if non-existent
643  static const faMesh& mesh(const polyMesh& pMesh);
644 
645  //- Return access to polyMesh
646  inline const polyMesh& mesh() const;
647 
648  //- Return the local mesh directory (dbDir()/meshSubDir)
649  fileName meshDir() const;
650 
651  //- Return reference to time
652  const Time& time() const;
653 
654  //- Return the current instance directory for points
655  // Used in the construction of geometric mesh data dependent
656  // on points
657  const fileName& pointsInstance() const;
658 
659  //- Return the current instance directory for faces
660  const fileName& facesInstance() const;
661 
662  //- Const reference to the mesh and solver state data
663  const meshState& data() const { return mesh().data(); }
664 
665  //- Reference to the mesh and solver state data
666  meshState& data() { return const_cast<polyMesh&>(mesh()).data(); }
667 
668 
669  // Parallel
670 
671  //- Return communicator used for parallel communication
672  label comm() const noexcept { return comm_; }
673 
674  //- Return communicator used for parallel communication
675  label& comm() noexcept { return comm_; }
676 
677  //- Is demand-driven parallel info available?
678  bool hasGlobalData() const noexcept;
679 
680  //- Return parallel info (demand-driven)
681  const faGlobalMeshData& globalData() const;
682 
683 
684  // Solution Control
686  //- Non-null if faSchemes exists (can test as bool).
687  const faSchemes* hasSchemes() const;
688 
689  //- Non-null if faSolution exists (can test as bool).
690  const faSolution* hasSolution() const;
691 
692  //- Read-access to the faSchemes controls
693  const faSchemes& schemes() const;
694 
695  //- Read/write-access to the faSchemes controls
697 
698  //- Read-access to the faSolution controls
699  const faSolution& solution() const;
700 
701  //- Read/write-access to the faSolution controls
702  faSolution& solution();
703 
704 
705  // Access: Mesh size parameters
706 
707  //- Number of local mesh points
708  inline label nPoints() const noexcept;
709 
710  //- Number of local mesh edges
711  inline label nEdges() const noexcept;
712 
713  //- Number of internal faces
714  inline label nInternalEdges() const noexcept;
715 
716  //- Number of boundary edges (== nEdges - nInternalEdges)
717  inline label nBoundaryEdges() const noexcept;
718 
719  //- Number of patch faces
720  inline label nFaces() const noexcept;
721 
722 
723  // Access: Primitive mesh data
724 
725  //- Return local points
726  inline const pointField& points() const;
727 
728  //- Return local edges with reordered boundary
729  inline const edgeList& edges() const noexcept;
730 
731  //- Sub-list of local internal edges
732  inline const edgeList::subList internalEdges() const;
733 
734  //- Return local faces
735  inline const faceList& faces() const;
736 
737  //- Edge owner addressing
738  inline const labelList& edgeOwner() const noexcept;
739 
740  //- Edge neighbour addressing
741  inline const labelList& edgeNeighbour() const noexcept;
742 
743  //- True if the internalEdges use an ordering that does not
744  //- correspond 1-to-1 with the patch internalEdges.
745  inline bool hasInternalEdgeLabels() const noexcept;
746 
747 
748  // Database
749 
750  //- Return true if thisDb() is a valid DB
751  virtual bool hasDb() const;
752 
753  //- Return reference to the mesh database
754  virtual const objectRegistry& thisDb() const;
755 
756  //- Name function is needed to disambiguate those inherited
757  //- from base classes
758  const word& name() const
759  {
760  return thisDb().name();
761  }
762 
763 
764  // Regions
765 
766  //- The mesh region name or word::null if polyMesh::defaultRegion
767  const word& regionName() const;
768 
769 
770  // Access
771 
772  //- Return constant reference to boundary mesh
773  inline const faBoundaryMesh& boundary() const noexcept;
774 
775  //- Return the underlying polyMesh face labels
776  inline const labelList& faceLabels() const noexcept;
777 
778  //- The area-face corresponding to the mesh-face, -1 if not found
779  inline label whichFace(const label meshFacei) const;
780 
781  //- The polyPatches related to the areaMesh, in sorted order
782  inline const labelList& whichPolyPatches() const;
783 
784  //- The polyPatch/local-face for each faceLabels()
785  inline const List<labelPair>& whichPatchFaces() const;
786 
787  //- Return ldu addressing
788  virtual const lduAddressing& lduAddr() const;
789 
790  //- Return a list of pointers for each patch
791  // with only those pointing to interfaces being set
792  virtual lduInterfacePtrsList interfaces() const
793  {
794  return boundary().interfaces();
795  }
796 
797  //- Internal face owner
798  const labelUList& owner() const
799  {
800  return lduAddr().lowerAddr();
801  }
802 
803  //- Internal face neighbour
804  const labelUList& neighbour() const
805  {
806  return lduAddr().upperAddr();
807  }
808 
809  //- True if given edge label is internal to the mesh
810  bool isInternalEdge(const label edgeIndex) const noexcept
811  {
812  return (edgeIndex < nInternalEdges_);
813  }
814 
815  //- List of proc/face for the boundary edge neighbours
816  //- using primitive patch edge numbering.
817  inline const List<labelPair>& boundaryConnections() const;
818 
819  //- Boundary edge neighbour processors
820  //- (does not include own proc)
821  labelList boundaryProcs() const;
822 
823  //- List of proc/size for the boundary edge neighbour processors
824  //- (does not include own proc)
825  List<labelPair> boundaryProcSizes() const;
826 
827  //- Mapping/swapping for boundary halo neighbours
828  const faMeshBoundaryHalo& boundaryHaloMap() const;
829 
830  //- Face centres of boundary halo neighbours
831  const pointField& haloFaceCentres() const;
832 
833  //- Face unit-normals of boundary halo neighbours
834  const vectorField& haloFaceNormals() const;
835 
836  //- Face centres of boundary halo neighbours for specified patch
837  tmp<pointField> haloFaceCentres(const label patchi) const;
838 
839  //- Face unit-normals of boundary halo neighbours for specified patch
840  tmp<vectorField> haloFaceNormals(const label patchi) const;
841 
842 
843  // Interfacing
844 
845  //- The volume (owner) cells associated with the area-mesh
846  labelList faceCells() const;
847 
848 
849  // Storage Management
850 
851  //- Remove all files from mesh instance
852  void removeFiles(const fileName& instanceDir) const;
853 
854  //- Remove all files from mesh instance()
855  void removeFiles() const;
856 
857 
858  //- Has face areas: S()
859  bool hasFaceAreas() const noexcept { return bool(SPtr_); }
860 
861  //- Has face centres: areaCentres()
862  bool hasAreaCentres() const noexcept { return bool(faceCentresPtr_); }
863 
864  //- Has edge centres: edgeCentres()
865  bool hasEdgeCentres() const noexcept { return bool(edgeCentresPtr_); }
866 
867  //- Has edge length vectors: Le()
868  bool hasLe() const noexcept { return bool(LePtr_); }
869 
870  //- Has edge length magnitudes: magLe()
871  bool hasMagLe() const noexcept { return bool(magLePtr_); }
872 
873  //- Has face area normals: faceAreaNormals()
874  bool hasFaceAreaNormals() const noexcept
875  {
876  return bool(faceAreaNormalsPtr_);
877  }
878 
879  //- Has edge area normals: edgeAreaNormals()
880  bool hasEdgeAreaNormals() const noexcept
881  {
882  return bool(edgeAreaNormalsPtr_);
883  }
884 
885  //- Has point area normals: pointAreaNormals()
886  bool hasPointAreaNormals() const noexcept
887  {
888  return bool(pointAreaNormalsPtr_);
889  }
890 
891  //- Has face curvatures: faceCurvatures()
892  bool hasFaceCurvatures() const noexcept
893  {
894  return bool(faceCurvaturesPtr_);
895  }
896 
897  //- Has patch point normals corrections
899  {
900  return bool(correctPatchPointNormalsPtr_);
901  }
902 
903 
904  // Mesh motion and morphing
906  //- Is mesh moving
907  bool moving() const
908  {
909  return mesh().moving();
910  }
911 
912  //- Update after mesh motion
913  virtual bool movePoints();
914 
915  //- Update after topo change
916  virtual void updateMesh(const mapPolyMesh&);
917 
919  // Mapping
920 
921  //- Map all fields in time using given map.
922  virtual void mapFields(const faMeshMapper& mapper) const;
923 
924  //- Map face areas in time using given map.
925  virtual void mapOldAreas(const faMeshMapper& mapper) const;
926 
927 
928  // Demand-driven data
929 
930  //- Return constant reference to primitive patch
931  inline const uindirectPrimitivePatch& patch() const;
932 
933  //- Return reference to primitive patch
934  inline uindirectPrimitivePatch& patch();
935 
936  //- Return patch starts
937  const labelList& patchStarts() const;
938 
939  //- Return edge length vectors
940  const edgeVectorField& Le() const;
941 
942  //- Return edge length magnitudes
943  const edgeScalarField& magLe() const;
944 
945  //- Return normalised edge length vectors
947 
948  //- Return face centres as areaVectorField
949  const areaVectorField& areaCentres() const;
950 
951  //- Return edge centres as edgeVectorField
952  const edgeVectorField& edgeCentres() const;
953 
954  //- Return face areas
955  const DimensionedField<scalar, areaMesh>& S() const;
956 
957  //- Return old-time face areas
958  const DimensionedField<scalar, areaMesh>& S0() const;
959 
960  //- Return old-old-time face areas
962 
963  //- Return face area normals
964  const areaVectorField& faceAreaNormals() const;
965 
966  //- Return edge area normals
967  const edgeVectorField& edgeAreaNormals() const;
968 
969  //- Return point area normals
970  const vectorField& pointAreaNormals() const;
971 
972  //- Return face curvatures
973  const areaScalarField& faceCurvatures() const;
974 
975  //- Return edge transformation tensors
977 
978  //- Return internal point labels
979  labelList internalPoints() const;
980 
981  //- Return boundary point labels
982  labelList boundaryPoints() const;
983 
984  //- Return edge length correction
986 
987  //- Whether point normals should be corrected for a patch
988  bool correctPatchPointNormals(const label patchID) const;
989 
990  //- Set whether point normals should be corrected for a patch
992 
993 
994  // Storage management
995 
996 
997  //- Write mesh
998  virtual bool write(const bool writeOnProc = true) const;
999 
1000 
1001  // Member Operators
1002 
1003  bool operator!=(const faMesh& m) const;
1004 
1005  bool operator==(const faMesh& m) const;
1006 
1007 
1008  // Housekeeping
1009 
1010  //- No call operator. Prior to 2312 was used to obtain polyMesh
1011  void operator()() const = delete;
1012 };
1013 
1014 
1015 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1016 
1017 } // End namespace Foam
1018 
1019 
1020 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1021 
1022 #include "faMeshI.H"
1023 
1024 #ifdef NoRepository
1025  #include "faPatchFaMeshTemplates.C"
1026 #endif
1027 
1028 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1029 
1030 #endif
1031 
1032 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:87
labelList internalPoints() const
Return internal point labels.
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
tmp< edgeScalarField > edgeLengthCorrection() const
Return edge length correction.
const DimensionedField< scalar, areaMesh > & S00() const
Return old-old-time face areas.
Definition: faMesh.C:940
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...
const Time & time() const
Return reference to time.
Definition: faMesh.C:791
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
void syncGeom()
Processor/processor synchronisation for geometry fields.
Definition: faMesh.C:292
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:913
const edgeList::subList internalEdges() const
Sub-list of local internal edges.
Definition: faMeshI.H:79
const pointField & haloFaceCentres() const
Face centres of boundary halo neighbours.
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:914
Face to edge interpolation scheme. Included in faMesh.
bool operator!=(const faMesh &m) const
Definition: faMesh.C:1154
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:129
virtual const meshState & data() const noexcept
Const reference to the mesh and solver state data.
Definition: polyMesh.H:558
Various mesh related information for a parallel run.
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:135
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
const faGlobalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: faMesh.C:1033
const labelList & patchStarts() const
Return patch starts.
Definition: faMesh.C:858
bool isInternalEdge(const label edgeIndex) const noexcept
True if given edge label is internal to the mesh.
Definition: faMesh.H:1122
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:1074
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: faMesh.C:821
const List< labelPair > & boundaryConnections() const
List of proc/face for the boundary edge neighbours using primitive patch edge numbering.
Definition: faMeshI.H:162
const areaVectorField & areaCentres() const
Return face centres as areaVectorField.
Definition: faMesh.C:891
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:1215
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
virtual const objectRegistry & thisDb() const
Return reference to the mesh database.
Definition: faMesh.C:815
const edgeScalarField & magLe() const
Return edge length magnitudes.
Definition: faMesh.C:880
bool hasGlobalData() const noexcept
Is demand-driven parallel info available?
Definition: faMesh.C:1027
label nFaces() const noexcept
Number of patch faces.
Definition: faMeshI.H:61
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
label nPoints() const noexcept
Number of local mesh points.
Definition: faMeshI.H:37
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
#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:109
static int geometryOrder() noexcept
Return the current geometry treatment.
Definition: faMesh.H:797
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:85
const labelUList & neighbour() const
Internal face neighbour.
Definition: faMesh.H:1114
virtual bool hasDb() const
Return true if thisDb() is a valid DB.
Definition: faMesh.C:809
const faceList & faces() const
Return local faces.
Definition: faMeshI.H:85
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:1210
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
const labelList & faceLabels() const noexcept
Return the underlying polyMesh face labels.
Definition: faMeshI.H:103
virtual ~faMesh()
Destructor.
Definition: faMesh.C:741
bool hasEdgeAreaNormals() const noexcept
Has edge area normals: edgeAreaNormals()
Definition: faMesh.H:1228
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:785
bool operator==(const faMesh &m) const
Definition: faMesh.C:1160
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: faMesh.C:1044
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:869
bool hasFaceCurvatures() const noexcept
Has face curvatures: faceCurvatures()
Definition: faMesh.H:1244
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:900
bool hasPointAreaNormals() const noexcept
Has point area normals: pointAreaNormals()
Definition: faMesh.H:1236
bool hasAreaCentres() const noexcept
Has face centres: areaCentres()
Definition: faMesh.H:1200
static const word prefix
The prefix to local: finite-area.
Definition: faMesh.H:696
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:155
const pointField & points() const
Return local points.
Definition: faMeshI.H:67
const vectorField & haloFaceNormals() const
Face unit-normals of boundary halo neighbours.
patchWriters clear()
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
const FieldField< Field, tensor > & edgeTransformTensors() const
Return edge transformation tensors.
Definition: faMesh.C:1016
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:1055
virtual bool write(const bool writeOnProc=true) const
Write mesh.
Definition: faMesh.C:1143
const edgeList & edges() const noexcept
Return local edges with reordered boundary.
Definition: faMeshI.H:73
faBoundaryMesh BoundaryMesh
The boundary type associated with the mesh.
Definition: faMesh.H:677
const polyMesh & mesh() const
Return access to polyMesh.
Definition: faMeshI.H:24
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:31
const areaScalarField & faceCurvatures() const
Return face curvatures.
Definition: faMesh.C:1004
label nEdges() const noexcept
Number of local mesh edges.
Definition: faMeshI.H:43
boolList & correctPatchPointNormals() const
Set whether point normals should be corrected for a patch.
Definition: faMesh.C:1132
const labelUList & owner() const
Internal face owner.
Definition: faMesh.H:1106
bool moving() const noexcept
Is mesh moving.
Definition: polyMesh.H:731
labelList boundaryPoints() const
Return boundary point labels.
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:701
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: faMesh.C:803
label nInternalEdges() const noexcept
Number of internal faces.
Definition: faMeshI.H:49
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:902
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
bool moving() const
Is mesh moving.
Definition: faMesh.H:1263
bool hasEdgeCentres() const noexcept
Has edge centres: edgeCentres()
Definition: faMesh.H:1205
static autoPtr< faMesh > TryNew(const polyMesh &pMesh)
Read construction from polyMesh if all files are available.
Definition: faMeshNew.C:149
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
const faSchemes * hasSchemes() const
Non-null if faSchemes exists (can test as bool).
Definition: faMesh.C:749
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
bool hasFaceAreas() const noexcept
Has face areas: S()
Definition: faMesh.H:1195
const labelList & edgeNeighbour() const noexcept
Edge neighbour addressing.
Definition: faMeshI.H:97
const word & name() const
Name function is needed to disambiguate those inherited from base classes.
Definition: faMesh.H:1047
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: faMesh.C:797
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:975
bool init(const bool doInit)
Initialise non-demand-driven data etc.
Definition: faMesh.C:313
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:827
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:66
static int geometryOrder_
Geometry treatment.
Definition: faMesh.H:685
faMesh Mesh
The mesh type.
Definition: faMesh.H:672
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.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:92
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const List< labelPair > & whichPatchFaces() const
The polyPatch/local-face for each faceLabels()
Definition: faMeshI.H:145
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:852
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:55
bool hasFaceAreaNormals() const noexcept
Has face area normals: faceAreaNormals()
Definition: faMesh.H:1220
const vectorField & pointAreaNormals() const
Return point area normals.
Definition: faMesh.C:986
const areaVectorField & faceAreaNormals() const
Return face area normals.
Definition: faMesh.C:964
bool hasPatchPointNormalsCorrection() const noexcept
Has patch point normals corrections.
Definition: faMesh.H:1252
tmp< edgeVectorField > unitLe() const
Return normalised edge length vectors.
const faSolution * hasSolution() const
Non-null if faSolution exists (can test as bool).
Definition: faMesh.C:755
const labelList & edgeOwner() const noexcept
Edge owner addressing.
Definition: faMeshI.H:91
lduAddressing wrapper for faMesh
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:926
const faSchemes & schemes() const
Read-access to the faSchemes controls.
Definition: faMesh.C:761
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:1098