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