polyMesh.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) 2011-2017, 2020 OpenFOAM Foundation
9  Copyright (C) 2018-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::polyMesh
29 
30 Description
31  Mesh consisting of general polyhedral cells.
32 
33 SourceFiles
34  polyMesh.C
35  polyMeshInitMesh.C
36  polyMeshClear.C
37  polyMeshFromShapeMesh.C
38  polyMeshIO.C
39  polyMeshUpdate.C
40  polyMeshCheck.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef Foam_polyMesh_H
45 #define Foam_polyMesh_H
46 
47 #include "objectRegistry.H"
48 #include "primitiveMesh.H"
49 #include "pointField.H"
50 #include "faceList.H"
51 #include "cellList.H"
52 #include "cellShapeList.H"
53 #include "pointIOField.H"
54 #include "faceIOList.H"
55 #include "labelIOList.H"
56 #include "polyBoundaryMesh.H"
57 #include "boundBox.H"
58 #include "globalIndex.H"
59 #include "pointZoneMesh.H"
60 #include "faceZoneMesh.H"
61 #include "cellZoneMesh.H"
62 #include "meshState.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
68 
69 // Forward Declarations
70 class globalMeshData;
71 class mapPolyMesh;
72 class polyMeshTetDecomposition;
73 class treeDataCell;
74 template<class Type> class indexedOctree;
75 
76 /*---------------------------------------------------------------------------*\
77  Class polyMesh Declaration
78 \*---------------------------------------------------------------------------*/
79 
80 class polyMesh
81 :
82  public objectRegistry,
83  public primitiveMesh
84 {
85 public:
86 
87  // Public Data
88 
89  //- Enumeration defining the state of the mesh after a read update.
90  // Used for post-processing applications, where the mesh
91  // needs to update based on the files written in time
92  // directories
94  {
99  };
100 
101  //- Enumeration defining the decomposition of the cell for
102  // inside/outside test
103  enum cellDecomposition
104  {
105  FACE_PLANES, //- Faces considered as planes
107  FACE_CENTRE_TRIS, //- Faces decomposed into triangles
108  // using face-centre
109 
110  FACE_DIAG_TRIS, //- Faces decomposed into triangles diagonally
112  CELL_TETS //- Cell decomposed into tets
113  };
114 
115 
116 private:
117 
118  // Private Data
119 
120  //- Mesh and solver state data
121  meshState data_;
122 
123 
124  // Primitive mesh data
125 
126  //- Points
127  pointIOField points_;
128 
129  //- Faces
130  faceCompactIOList faces_;
131 
132  //- Face owner
133  labelIOList owner_;
134 
135  //- Face neighbour
136  labelIOList neighbour_;
137 
138  //- Have the primitives been cleared
139  bool clearedPrimitives_;
140 
141 
142  //- Boundary mesh
143  mutable polyBoundaryMesh boundary_;
144 
145  //- Mesh bounding-box.
146  // Created from points on construction, updated when the mesh moves
147  boundBox bounds_;
148 
149  //- Communicator used for parallel communication
150  label comm_;
151 
152  //- Vector of non-constrained directions in mesh
153  // defined according to the presence of empty and wedge patches
154  mutable Vector<label> geometricD_;
155 
156  //- Vector of valid directions in mesh
157  // defined according to the presence of empty patches
158  mutable Vector<label> solutionD_;
159 
160  //- Base point for face decomposition into tets
161  mutable autoPtr<labelIOList> tetBasePtIsPtr_;
162 
163  //- Search tree to allow spatial cell searching
164  mutable autoPtr<indexedOctree<treeDataCell>> cellTreePtr_;
165 
166 
167  // Zoning information
168 
169  //- Point zones
170  pointZoneMesh pointZones_;
171 
172  //- Face zones
173  faceZoneMesh faceZones_;
174 
175  //- Cell zones
176  cellZoneMesh cellZones_;
177 
178 
179  //- Parallel info
180  mutable autoPtr<globalMeshData> globalMeshDataPtr_;
181 
182 
183  // Mesh motion related data
184 
185  //- Is the mesh moving
186  bool moving_;
187 
188  //- Is the mesh topology changing
189  bool topoChanging_;
190 
191  //- Store old cell centres?
192  mutable bool storeOldCellCentres_;
193 
194  //- Current time index for mesh motion
195  mutable label curMotionTimeIndex_;
196 
197  //- Old points (for the last mesh motion)
198  mutable autoPtr<pointField> oldPointsPtr_;
199 
200  //- Old cell centres (for the last mesh motion)
201  mutable autoPtr<pointField> oldCellCentresPtr_;
202 
203 
204 
205  // Private Member Functions
206 
207  //- No copy construct
208  polyMesh(const polyMesh&) = delete;
209 
210  //- No copy assignment
211  void operator=(const polyMesh&) = delete;
212 
213  //- Initialise the polyMesh from the primitive data
214  void initMesh();
215 
216  //- Initialise the polyMesh from the given set of cells
217  void initMesh(cellList& c);
218 
219  //- Calculate the valid directions in the mesh from the boundaries
220  void calcDirections() const;
221 
222  //- Calculate the cell shapes from the primitive
223  // polyhedral information
224  void calcCellShapes() const;
225 
226  //- Read and return the tetBasePtIs
227  autoPtr<labelIOList> readTetBasePtIs() const;
228 
229 
230  // Helper functions for constructor from cell shapes
231 
232  labelListList cellShapePointCells(const cellShapeList&) const;
233 
234  labelList facePatchFaceCells
235  (
236  const faceList& patchFaces,
237  const labelListList& pointCells,
238  const faceListList& cellsFaceShapes,
239  const label patchID
240  ) const;
241 
242  void setTopology
243  (
244  const cellShapeList& cellsAsShapes,
245  const faceListList& boundaryFaces,
246  const wordList& boundaryPatchNames,
247  labelList& patchSizes,
248  labelList& patchStarts,
249  label& defaultPatchStart,
250  label& nFaces,
251  cellList& cells
252  );
253 
254 
255  // Geometry checks
256 
257  //- Check non-orthogonality
258  bool checkFaceOrthogonality
259  (
260  const vectorField& fAreas,
261  const vectorField& cellCtrs,
262  const bool report,
263  const bool detailedReport,
264  labelHashSet* setPtr
265  ) const;
266 
267  //- Check face skewness
268  bool checkFaceSkewness
269  (
270  const pointField& points,
271  const vectorField& fCtrs,
272  const vectorField& fAreas,
273  const vectorField& cellCtrs,
274  const bool report,
275  const bool detailedReport,
276  labelHashSet* setPtr
277  ) const;
278 
279  bool checkEdgeAlignment
280  (
281  const pointField& p,
282  const bool report,
283  const Vector<label>& directions,
284  labelHashSet* setPtr
285  ) const;
286 
287  bool checkCellDeterminant
288  (
289  const vectorField& faceAreas,
290  const bool report,
291  labelHashSet* setPtr,
292  const Vector<label>& meshD
293  ) const;
294 
295  bool checkFaceWeight
296  (
297  const vectorField& fCtrs,
298  const vectorField& fAreas,
299  const vectorField& cellCtrs,
300  const bool report,
301  const scalar minWeight,
302  labelHashSet* setPtr
303  ) const;
304 
305  bool checkVolRatio
306  (
307  const scalarField& cellVols,
308  const bool report,
309  const scalar minRatio,
310  labelHashSet* setPtr
311  ) const;
312 
313 
314 public:
315 
316  // Public Typedefs
317 
318  //- The mesh type
319  typedef polyMesh Mesh;
320 
321  //- The boundary type associated with the mesh
323 
324 
325  //- Runtime type information
326  TypeName("polyMesh");
327 
328  //- Return the default region name
329  static word defaultRegion;
330 
331  //- Return the mesh sub-directory name (usually "polyMesh")
332  static word meshSubDir;
333 
334 
335  // Constructors
336 
337  //- Read construct from IOobject
338  explicit polyMesh(const IOobject& io, const bool doInit = true);
339 
340  //- Construct from IOobject or as zero-sized mesh
341  // Boundary is added using addPatches() member function
342  polyMesh(const IOobject& io, const Foam::zero, const bool syncPar=true);
343 
344  //- Construct from IOobject and components.
345  // Boundary is added using addPatches() member function
346  polyMesh
347  (
348  const IOobject& io,
349  pointField&& points,
350  faceList&& faces,
351  labelList&& owner,
352  labelList&& neighbour,
353  const bool syncPar = true
354  );
355 
356  //- Construct without boundary with cells rather than owner/neighbour.
357  // Boundary is added using addPatches() member function
358  polyMesh
359  (
360  const IOobject& io,
361  pointField&& points,
362  faceList&& faces,
363  cellList&& cells,
364  const bool syncPar = true
365  );
366 
367  //- Construct from cell shapes
368  polyMesh
369  (
370  const IOobject& io,
371  pointField&& points,
372  const cellShapeList& shapes,
373  const faceListList& boundaryFaces,
374  const wordList& boundaryPatchNames,
375  const wordList& boundaryPatchTypes,
376  const word& defaultBoundaryPatchName,
377  const word& defaultBoundaryPatchType,
378  const wordList& boundaryPatchPhysicalTypes,
379  const bool syncPar = true
380  );
381 
382  //- Construct from cell shapes, with patch information in dictionary
383  //- format.
384  polyMesh
385  (
386  const IOobject& io,
387  pointField&& points,
388  const cellShapeList& shapes,
389  const faceListList& boundaryFaces,
390  const wordList& boundaryPatchNames,
392  const word& defaultBoundaryPatchName,
393  const word& defaultBoundaryPatchType,
394  const bool syncPar = true
395  );
396 
397 
398  //- Destructor
399  virtual ~polyMesh();
400 
401 
402  // Member Functions
403 
404  // Database
405 
406  //- Override the objectRegistry dbDir for a single-region case
407  virtual const fileName& dbDir() const;
408 
409  //- Return the local mesh directory (dbDir()/meshSubDir)
410  fileName meshDir() const;
412  //- Return the current instance directory for points
413  // Used in the construction of geometric mesh data dependent
414  // on points
415  const fileName& pointsInstance() const;
416 
417  //- Return the current instance directory for faces
418  const fileName& facesInstance() const;
419 
420  //- Set the instance for mesh files
421  void setInstance
422  (
423  const fileName& instance,
425  );
426 
427 
428  // Regions
429 
430  //- Return the local mesh directory name (eg, "polyMesh")
431  //- after applying filter for defaultRegion
432  static fileName meshDir(const word& region);
433 
434  //- The mesh region name or word::null if polyMesh::defaultRegion
435  static const word& regionName(const word& region);
436 
437  //- The mesh region name or word::null if polyMesh::defaultRegion
438  const word& regionName() const;
439 
440 
441  // Access
442 
443  //- Const reference to the mesh and solver state data
444  virtual const meshState& data() const noexcept { return data_; }
445 
446  //- Reference to the mesh and solver state data
447  virtual meshState& data() noexcept { return data_; }
448 
449  //- Return raw points
450  virtual const pointField& points() const;
451 
452  //- Return true if io is up-to-date with points
453  virtual bool upToDatePoints(const regIOobject& io) const;
454 
455  //- Set io to be up-to-date with points
456  virtual void setUpToDatePoints(regIOobject& io) const;
457 
458  //- Return raw faces
459  virtual const faceList& faces() const;
460 
461  //- Return face owner
462  virtual const labelList& faceOwner() const;
463 
464  //- Return face neighbour
465  virtual const labelList& faceNeighbour() const;
466 
467  //- Return old points (mesh motion)
468  virtual const pointField& oldPoints() const;
469 
470  //- Return old cellCentres (mesh motion)
471  virtual const pointField& oldCellCentres() const;
472 
473  //- Return boundary mesh
474  const polyBoundaryMesh& boundaryMesh() const noexcept
475  {
476  return boundary_;
477  }
478 
479  //- Return mesh bounding box
480  const boundBox& bounds() const noexcept
481  {
482  return bounds_;
483  }
484 
485  //- Return the vector of geometric directions in mesh.
486  // Defined according to the presence of empty and wedge patches.
487  // 1 indicates unconstrained direction and -1 a constrained
488  // direction.
489  const Vector<label>& geometricD() const;
490 
491  //- Return the number of valid geometric dimensions in the mesh
492  label nGeometricD() const;
493 
494  //- Return the vector of solved-for directions in mesh.
495  // Differs from geometricD in that it includes for wedge cases
496  // the circumferential direction in case of swirl.
497  // 1 indicates valid direction and -1 an invalid direction.
498  const Vector<label>& solutionD() const;
499 
500  //- Return the number of valid solved-for dimensions in the mesh
501  label nSolutionD() const;
502 
503  //- Return the tetBasePtIs
504  const labelIOList& tetBasePtIs() const;
505 
506  //- Return the cell search tree
507  const indexedOctree<treeDataCell>& cellTree() const;
508 
509  //- Return point zone mesh
510  const pointZoneMesh& pointZones() const noexcept
511  {
512  return pointZones_;
513  }
514 
515  //- Return face zone mesh
516  const faceZoneMesh& faceZones() const noexcept
517  {
518  return faceZones_;
519  }
520 
521  //- Return cell zone mesh
522  const cellZoneMesh& cellZones() const noexcept
523  {
524  return cellZones_;
525  }
526 
527 
528  // Database
529 
530  //- Return the object registry
531  const objectRegistry& thisDb() const noexcept
532  {
533  return *this;
534  }
535 
536 
537  // Parallel
538 
539  //- The communicator used for parallel communication
540  label comm() const noexcept { return comm_; }
541 
542  //- The communicator used for parallel communication
543  label& comm() noexcept { return comm_; }
544 
545  //- Is demand-driven parallel info available?
546  bool hasGlobalData() const noexcept;
547 
548  //- Return parallel info (demand-driven)
549  const globalMeshData& globalData() const;
550 
551 
552  // Mesh motion
553 
554  //- Is mesh dynamic
555  virtual bool dynamic() const
556  {
557  return false;
558  }
560  //- Is mesh moving
561  bool moving() const noexcept
562  {
563  return moving_;
564  }
565 
566  //- Set the mesh to be moving
567  bool moving(const bool on) noexcept
568  {
569  bool old(moving_);
570  moving_ = on;
571  return old;
572  }
573 
574  //- Is mesh topology changing
575  bool topoChanging() const noexcept
576  {
577  return topoChanging_;
578  }
579 
580  //- Set the mesh topology to be changing
581  bool topoChanging(const bool on) noexcept
582  {
583  bool old(topoChanging_);
584  topoChanging_ = on;
585  return old;
586  }
587 
588  //- Is mesh changing (topology changing and/or moving)
589  bool changing() const noexcept
590  {
591  return (moving() || topoChanging());
592  }
593 
594  //- Move points
595  virtual void movePoints(const pointField&);
596 
597  //- Reset motion
598  void resetMotion() const;
599 
600 
601  // Topological change
602 
603  //- Return non-const access to the pointZones
605  {
606  return pointZones_;
607  }
608 
609  //- Return non-const access to the faceZones
611  {
612  return faceZones_;
613  }
614 
615  //- Return non-const access to the cellZones
617  {
618  return cellZones_;
619  }
620 
621  //- Add boundary patches
622  void addPatches
623  (
624  polyPatchList& plist,
625  const bool validBoundary = true
626  );
627 
628  //- Add boundary patches
629  void addPatches
630  (
631  const List<polyPatch*>& p,
632  const bool validBoundary = true
633  );
634 
635  //- Add mesh zones
636  void addZones
637  (
638  PtrList<pointZone>&& pz,
639  PtrList<faceZone>&& fz,
640  PtrList<cellZone>&& cz
641  );
642 
643  //- Add mesh zones
644  void addZones
645  (
646  const List<pointZone*>& pz,
647  const List<faceZone*>& fz,
648  const List<cellZone*>& cz
649  );
650 
651  //- Initialise all non-demand-driven data
652  virtual bool init(const bool doInit);
653 
654  //- Update the mesh based on the mesh files saved in
655  // time directories
656  virtual readUpdateState readUpdate();
657 
658  //- Update the mesh corresponding to given map
659  virtual void updateMesh(const mapPolyMesh& mpm);
660 
661  //- Remove boundary patches
662  void removeBoundary();
664  //- Reset mesh primitive data. Assumes all patch info correct
665  // (so does e.g. parallel communication). If not use
666  // validBoundary=false
667  //
668  // \note Only autoPtr parameters that test as good() are used
669  // for resetting, otherwise the existing entries are left
670  // untouched.
672  (
675  autoPtr<labelList>&& owner,
676  autoPtr<labelList>&& neighbour,
677  const labelUList& patchSizes,
678  const labelUList& patchStarts,
679  const bool validBoundary = true
680  );
681 
682 
683  // Storage management
684 
685  //- Clear geometry
686  void clearGeom();
687 
688  //- Update geometry points; keep topology.
689  //- Optionally with new face decomposition
691  (
692  pointIOField&& newPoints,
693  autoPtr<labelIOList>& newTetBasePtIsPtr
694  );
695 
696  //- Clear addressing
697  void clearAddressing(const bool isMeshUpdate = false);
698 
699  //- Clear all geometry and addressing
700  void clearOut(const bool isMeshUpdate = false);
702  //- Clear primitive data (points, faces and cells)
703  void clearPrimitives();
704 
705  //- Clear tet base points
707 
708  //- Clear cell tree data
709  void clearCellTree();
710 
711  //- Remove all files from mesh instance
712  void removeFiles(const fileName& instanceDir) const;
713 
714  //- Remove all files from mesh instance()
715  void removeFiles() const;
716 
717 
718  bool hasTetBasePtIs() const { return bool(tetBasePtIsPtr_); }
719 
720 
721  // Geometric checks. Selectively override primitiveMesh functionality.
722 
723  //- Check non-orthogonality
724  virtual bool checkFaceOrthogonality
725  (
726  const bool report = false,
727  labelHashSet* setPtr = nullptr
728  ) const;
729 
730  //- Check face skewness
731  virtual bool checkFaceSkewness
732  (
733  const bool report = false,
734  labelHashSet* setPtr = nullptr
735  ) const;
736 
737  //- Check edge alignment for 1D/2D cases
738  virtual bool checkEdgeAlignment
739  (
740  const bool report,
741  const Vector<label>& directions,
742  labelHashSet* setPtr
743  ) const;
744 
745  virtual bool checkCellDeterminant
746  (
747  const bool report,
748  labelHashSet* setPtr
749  ) const;
751  //- Check mesh motion for correctness given motion points
752  virtual bool checkMeshMotion
753  (
754  const pointField& newPoints,
755  const bool report = false,
756  const bool detailedReport = false
757  ) const;
759  //- Check for face weights
760  virtual bool checkFaceWeight
761  (
762  const bool report,
763  const scalar minWeight = 0.05,
764  labelHashSet* setPtr = nullptr
765  ) const;
766 
767  //- Check for neighbouring cell volumes
768  virtual bool checkVolRatio
769  (
770  const bool report,
771  const scalar minRatio = 0.01,
772  labelHashSet* setPtr = nullptr
773  ) const;
774 
775 
776  // Position search functions
777 
778  //- Find the cell, tetFacei and tetPti for point p
779  void findCellFacePt
780  (
781  const point& p,
782  label& celli,
783  label& tetFacei,
784  label& tetPti
785  ) const;
786 
787  //- Find the tetFacei and tetPti for point p in celli.
788  // tetFacei and tetPtI are set to -1 if not found
790  (
791  const label celli,
792  const point& p,
793  label& tetFacei,
794  label& tetPti
795  ) const;
796 
797  //- Test if point p is in the celli
798  bool pointInCell
799  (
800  const point& p,
801  label celli,
803  ) const;
804 
805  //- Find cell enclosing this location and return index
806  // If not found -1 is returned
807  label findCell
808  (
809  const point& p,
811  ) const;
812 
813 
814  // Write
815 
816  //- Write items held in the objectRegistry. Normally includes mesh
817  //- components (points, faces, etc) and any registered fields.
818  virtual bool writeObject
819  (
820  IOstreamOption streamOpt,
821  const bool writeOnProc = true
822  ) const;
823 };
824 
825 
826 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
827 
828 } // End namespace Foam
829 
830 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
831 
832 #endif
833 
834 // ************************************************************************* //
writeOption
Enumeration defining write preferences.
label comm() const noexcept
The communicator used for parallel communication.
Definition: polyMesh.H:701
polyMesh Mesh
The mesh type.
Definition: polyMesh.H:390
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
void clearAddressing()
Clear topological data.
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: polyMesh.C:346
A class for handling file names.
Definition: fileName.H:72
void updateGeomPoints(pointIOField &&newPoints, autoPtr< labelIOList > &newTetBasePtIsPtr)
Update geometry points; keep topology. Optionally with new face decomposition.
Definition: polyMeshClear.C:66
virtual void movePoints(const pointField &)
Move points.
Definition: polyMesh.C:1168
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:859
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:899
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:847
virtual void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1097
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1297
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:882
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
Definition: regIOobject.C:43
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:63
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:411
Database for mesh data, solution data, solver performance and other reduced data. ...
Definition: meshState.H:53
Foam::pointZoneMesh.
virtual const meshState & data() const noexcept
Const reference to the mesh and solver state data.
Definition: polyMesh.H:559
void resetPrimitives(autoPtr< pointField > &&points, autoPtr< faceList > &&faces, autoPtr< labelList > &&owner, autoPtr< labelList > &&neighbour, const labelUList &patchSizes, const labelUList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:704
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition: pointCells.H:52
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
const cellList & cells() const
A simple container for options an IOstream can normally have.
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:830
label nFaces() const noexcept
Number of mesh faces.
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1401
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:104
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:876
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:865
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
Foam::cellZoneMesh.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:853
IOList< label > labelIOList
IO for a List of label.
Definition: labelIOList.H:32
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1353
void clearPrimitives()
Clear primitive data (points, faces and cells)
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
bool hasGlobalData() const noexcept
Is demand-driven parallel info available?
Definition: polyMesh.C:1305
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1360
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:406
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
Definition: polyMesh.C:1009
virtual const pointField & oldPoints() const
Return old points (mesh motion)
Definition: polyMesh.C:1128
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
void clearGeom()
Clear geometry.
Definition: polyMeshClear.C:44
bool topoChanging() const noexcept
Is mesh topology changing.
Definition: polyMesh.H:750
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:32
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:671
Info<< "Creating cells"<< endl;cellShapes=b.shapes();Info<< "Creating boundary faces"<< endl;boundary.setSize(b.boundaryPatches().size());forAll(boundary, patchi) { faceList faces(b.boundaryPatches()[patchi].size());forAll(faces, facei) { faces[facei]=face(b.boundaryPatches()[patchi][facei]);} boundary[patchi].transfer(faces);} points.transfer(const_cast< pointField & >b.points()));}Info<< "Creating patch dictionaries"<< endl;wordList patchNames(boundary.size());forAll(patchNames, patchi){ patchNames[patchi]=polyPatch::defaultName(patchi);}PtrList< dictionary > boundaryDicts(boundary.size())
const scalarField & cellVols
void clearTetBasePtIs()
Clear tet base points.
bool changing() const noexcept
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:768
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:266
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:893
virtual bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1091
bool moving() const noexcept
Is mesh moving.
Definition: polyMesh.H:732
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:663
TypeName("polyMesh")
Runtime type information.
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:939
virtual const pointField & oldCellCentres() const
Return old cellCentres (mesh motion)
Definition: polyMesh.C:1150
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
const vectorField & faceAreas() const
const dimensionedScalar c
Speed of light in a vacuum.
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:841
Automatically write from objectRegistry::writeObject()
virtual bool dynamic() const
Is mesh dynamic.
Definition: polyMesh.H:724
void findTetFacePt(const label celli, const point &p, label &tetFacei, label &tetPti) const
Find the tetFacei and tetPti for point p in celli.
Definition: polyMesh.C:1385
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:679
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1511
polyBoundaryMesh BoundaryMesh
The boundary type associated with the mesh.
Definition: polyMesh.H:395
virtual ~polyMesh()
Destructor.
Definition: polyMesh.C:802
const objectRegistry & thisDb() const noexcept
Return the object registry.
Definition: polyMesh.H:690
bool hasTetBasePtIs() const
Definition: polyMesh.H:945
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition: polyMesh.H:617
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:68
Registry of regIOobjects.
void addPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:971
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Foam::faceZoneMesh.
A primitive field of type <T> with automated input and output.
void clearCellTree()
Clear cell tree data.
objectRegistry(const Time &db, const label initialCapacity=128)
Construct the time objectRegistry, with estimated table capacity (default: 128)
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
virtual bool checkMeshMotion(const pointField &newPoints, const bool report=false, const bool detailedReport=false) const
Check mesh motion for correctness given motion points.
Namespace for OpenFOAM.
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc=true) const
Write items held in the objectRegistry. Normally includes mesh components (points, faces, etc) and any registered fields.
Definition: polyMesh.C:1592