polyMesh.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::polyMesh
26 
27 Description
28  Mesh consisting of general polyhedral cells.
29 
30 SourceFiles
31  polyMesh.C
32  polyMeshInitMesh.C
33  polyMeshClear.C
34  polyMeshFromShapeMesh.C
35  polyMeshIO.C
36  polyMeshUpdate.C
37  polyMeshCheck.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef polyMesh_H
42 #define polyMesh_H
43 
44 #include "objectRegistry.H"
45 #include "primitiveMesh.H"
46 #include "pointField.H"
47 #include "faceList.H"
48 #include "cellList.H"
49 #include "cellShapeList.H"
50 #include "pointIOField.H"
51 #include "faceIOList.H"
52 #include "labelIOList.H"
53 #include "polyBoundaryMesh.H"
54 #include "boundBox.H"
55 #include "meshPointZones.H"
56 #include "meshFaceZones.H"
57 #include "meshCellZones.H"
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 namespace Foam
62 {
63 
64 // Forward declaration of classes
65 class globalMeshData;
66 class polyTopoChangeMap;
67 class polyMeshMap;
68 class polyDistributionMap;
69 class polyMeshTetDecomposition;
70 class treeDataCell;
71 template<class Type> class indexedOctree;
72 
73 /*---------------------------------------------------------------------------*\
74  Class polyMesh Declaration
75 \*---------------------------------------------------------------------------*/
76 
77 class polyMesh
78 :
79  public objectRegistry,
80  public primitiveMesh
81 {
82 
83 public:
84 
85  // Public data types
86 
87  //- Enumeration defining the state of the mesh after a read update.
88  // Used for post-processing applications, where the mesh
89  // needs to update based on the files written in time
90  // directories
91  enum readUpdateState
92  {
97  };
98 
99  //- Enumeration defining the decomposition of the cell for
100  // inside/outside test
101  enum cellDecomposition
102  {
103  FACE_PLANES, //- Faces considered as planes
104 
105  FACE_CENTRE_TRIS, //- Faces decomposed into triangles
106  // using face-centre
107 
108  FACE_DIAG_TRIS, //- Faces decomposed into triangles diagonally
109 
110  CELL_TETS //- Cell decomposed into tets
111  };
112 
113 
114 private:
115 
116  // Permanent data
117 
118  // Primitive mesh data
119 
120  //- Points
121  pointIOField points_;
122 
123  //- Faces
124  faceCompactIOList faces_;
125 
126  //- Face owner
127  labelIOList owner_;
128 
129  //- Face neighbour
130  labelIOList neighbour_;
131 
132  //- Have the primitives been cleared
133  bool clearedPrimitives_;
134 
135 
136  //- Boundary mesh
137  mutable polyBoundaryMesh boundary_;
138 
139  //- Mesh bounding-box.
140  // Created from points on construction, updated when the mesh moves
141  boundBox bounds_;
142 
143  //- Communicator used for parallel communication
144  label comm_;
145 
146  //- Vector of non-constrained directions in mesh
147  // defined according to the presence of empty and wedge patches
148  mutable Vector<label> geometricD_;
149 
150  //- Vector of valid directions in mesh
151  // defined according to the presence of empty patches
152  mutable Vector<label> solutionD_;
153 
154  //- Base point for face decomposition into tets
155  mutable autoPtr<labelIOList> tetBasePtIsPtr_;
156 
157  //- Search tree to allow spatial cell searching
158  mutable autoPtr<indexedOctree<treeDataCell>> cellTreePtr_;
159 
160 
161  // Zoning information
162 
163  //- Point zones
164  meshPointZones pointZones_;
165 
166  //- Face zones
167  meshFaceZones faceZones_;
168 
169  //- Cell zones
170  meshCellZones cellZones_;
171 
172 
173  //- Parallel info
174  mutable autoPtr<globalMeshData> globalMeshDataPtr_;
175 
176 
177  // Mesh motion related dat
178 
179  //- Current time index for mesh motion
180  mutable label curMotionTimeIndex_;
181 
182  //- Old points (for the last mesh motion)
183  mutable autoPtr<pointField> oldPointsPtr_;
184 
185  //- Old cell centres (for the last mesh motion)
186  mutable autoPtr<pointField> oldCellCentresPtr_;
187 
188  //- Whether or not to store the old cell centres
189  mutable bool storeOldCellCentres_;
190 
191 
192  // Private Member Functions
193 
194  //- Initialise the polyMesh from the primitive data
195  void initMesh();
196 
197  //- Initialise the polyMesh from the given set of cells
198  void initMesh(cellList& c);
199 
200  //- Calculate the valid directions in the mesh from the boundaries
201  void calcDirections() const;
202 
203  //- Calculate the cell shapes from the primitive
204  // polyhedral information
205  void calcCellShapes() const;
206 
207  //- Read and return the tetBasePtIs
208  autoPtr<labelIOList> readTetBasePtIs() const;
209 
210  //- Set the write option of the points
211  void setPointsWrite(const IOobject::writeOption wo);
212 
213  //- Set the write option of the topology
214  void setTopologyWrite(const IOobject::writeOption wo);
215 
216 
217  // Helper functions for constructor from cell shapes
218 
219  labelListList cellShapePointCells(const cellShapeList&) const;
220 
221  labelList facePatchFaceCells
222  (
223  const faceList& patchFaces,
224  const labelListList& pointCells,
225  const faceListList& cellsFaceShapes,
226  const label patchID
227  ) const;
228 
229  void setTopology
230  (
231  const cellShapeList& cellsAsShapes,
232  const faceListList& boundaryFaces,
233  const wordList& boundaryPatchNames,
234  labelList& patchSizes,
235  labelList& patchStarts,
236  label& defaultPatchStart,
237  label& nFaces,
238  cellList& cells
239  );
240 
241 
242 protected:
243 
244  // Protected Member Data
245 
246  //- Member data pending transfer to fvMesh
247 
248  //- Is the mesh moving
249  bool moving_;
250 
251  //- Has the mesh topology changed
252  bool topoChanged_;
253 
254 
255 public:
256 
257  // Public Typedefs
258 
259  typedef polyMesh Mesh;
261 
262 
263  //- Runtime type information
264  TypeName("polyMesh");
265 
266  //- Return the default region name
267  static word defaultRegion;
268 
269  //- Return the mesh sub-directory name (usually "polyMesh")
270  static word meshSubDir;
271 
272 
273  // Constructors
274 
275  //- Construct from IOobject
276  explicit polyMesh(const IOobject& io);
277 
278  //- Move construct from IOobject or from components.
279  // Boundary is added using addPatches() member function
280  polyMesh
281  (
282  const IOobject& io,
283  pointField&& points,
284  faceList&& faces,
285  labelList&& owner,
286  labelList&& neighbour,
287  const bool syncPar = true
288  );
289 
290  //- Move construct without boundary with cells rather than
291  // owner/neighbour.
292  // Boundary is added using addPatches() member function
293  polyMesh
294  (
295  const IOobject& io,
296  pointField&& points,
297  faceList&& faces,
298  cellList&& cells,
299  const bool syncPar = true
300  );
301 
302  //- Move construct from cell shapes
303  polyMesh
304  (
305  const IOobject& io,
306  pointField&& points,
307  const cellShapeList& shapes,
308  const faceListList& boundaryFaces,
309  const wordList& boundaryPatchNames,
310  const wordList& boundaryPatchTypes,
311  const word& defaultBoundaryPatchName,
312  const word& defaultBoundaryPatchType,
313  const wordList& boundaryPatchPhysicalTypes,
314  const bool syncPar = true
315  );
316 
317  //- Move construct from cell shapes with patch information in dictionary
318  // format.
319  polyMesh
320  (
321  const IOobject& io,
322  pointField&& points,
323  const cellShapeList& shapes,
324  const faceListList& boundaryFaces,
325  const wordList& boundaryPatchNames,
326  const PtrList<dictionary>& boundaryDicts,
327  const word& defaultBoundaryPatchName,
328  const word& defaultBoundaryPatchType,
329  const bool syncPar = true
330  );
331 
332  //- Move constructor
333  polyMesh(polyMesh&&);
334 
335  //- Disallow default bitwise copy construction
336  polyMesh(const polyMesh&) = delete;
337 
338 
339  //- Destructor
340  virtual ~polyMesh();
341 
342 
343  // Member Functions
344 
345  // Database
346 
347  //- Override the objectRegistry dbDir for a single-region case
348  virtual const fileName& dbDir() const;
349 
350  //- Return the local mesh directory (dbDir()/meshSubDir)
351  fileName meshDir() const;
352 
353  //- Return the current instance directory for points
354  // Used in the construction of geometric mesh data dependent
355  // on points
356  const fileName& pointsInstance() const;
357 
358  //- Return the current instance directory for faces
359  const fileName& facesInstance() const;
360 
361  //- Return the points write option
363 
364  //- Return the points write option
366 
367  //- Set the instance for the points files
368  void setPointsInstance(const fileName&);
369 
370  //- Set the instance for mesh files
371  void setInstance(const fileName&);
372 
373 
374  // Access
375 
376  // Primitive mesh data
377 
378  //- Return raw points
379  virtual const pointField& points() const;
380 
381  //- Return raw faces
382  virtual const faceList& faces() const;
383 
384  //- Return face owner
385  virtual const labelList& faceOwner() const;
386 
387  //- Return face neighbour
388  virtual const labelList& faceNeighbour() const;
389 
390  //- Return old points for mesh motion
391  virtual const pointField& oldPoints() const;
392 
393  //- Return old cell centres for mesh motion
394  virtual const pointField& oldCellCentres() const;
395 
396 
397  //- Return true if io is up-to-date with points
398  bool upToDatePoints(const regIOobject& io) const;
399 
400  //- Set io to be up-to-date with points
401  void setUpToDatePoints(regIOobject& io) const;
402 
403  //- Return boundary mesh
404  const polyBoundaryMesh& boundaryMesh() const
405  {
406  return boundary_;
407  }
408 
409  //- Return mesh bounding box
410  const boundBox& bounds() const
411  {
412  return bounds_;
413  }
414 
415  //- Return the vector of geometric directions in mesh.
416  // Defined according to the presence of empty and wedge patches.
417  // 1 indicates unconstrained direction and -1 a constrained
418  // direction.
419  const Vector<label>& geometricD() const;
420 
421  //- Return the number of valid geometric dimensions in the mesh
422  label nGeometricD() const;
423 
424  //- Return the vector of solved-for directions in mesh.
425  // Differs from geometricD in that it includes for wedge cases
426  // the circumferential direction in case of swirl.
427  // 1 indicates valid direction and -1 an invalid direction.
428  const Vector<label>& solutionD() const;
429 
430  //- Return the number of valid solved-for dimensions in the mesh
431  label nSolutionD() const;
432 
433  //- Return the tetBasePtIs
434  const labelIOList& tetBasePtIs() const;
435 
436  //- Return the cell search tree
437  const indexedOctree<treeDataCell>& cellTree() const;
438 
439  //- Return point zones
440  const meshPointZones& pointZones() const
441  {
442  return pointZones_;
443  }
444 
445  //- Return face zones
446  const meshFaceZones& faceZones() const
447  {
448  return faceZones_;
449  }
450 
451  //- Return cell zones
452  const meshCellZones& cellZones() const
453  {
454  return cellZones_;
455  }
456 
457  //- Return parallel info
458  const globalMeshData& globalData() const;
459 
460  //- Return communicator used for parallel communication
461  label comm() const;
462 
463  //- Return communicator used for parallel communication
464  label& comm();
465 
466  //- Return the object registry
467  const objectRegistry& thisDb() const
468  {
469  return *this;
470  }
471 
472 
473  // Mesh motion
474 
475  //- Is mesh moving
476  bool moving() const
477  {
478  return moving_;
479  }
480 
481  //- Has the mesh topology changed this time-step
482  bool topoChanged() const
483  {
484  return topoChanged_;
485  }
486 
487  //- Is mesh changing
488  // Moving or mesh topology changed this time-step)
489  bool changing() const
490  {
491  return moving() || topoChanged();
492  }
493 
494  //- Reset the points
495  // without storing old points or returning swept volumes
496  virtual void setPoints(const pointField&);
497 
498  //- Move points, returns volumes swept by faces in motion
499  virtual tmp<scalarField> movePoints(const pointField&);
500 
501  //- Reset motion
502  void resetMotion() const;
503 
504 
505  // Topological change
506 
507  //- Return non-const access to the pointZones
509  {
510  return pointZones_;
511  }
512 
513  //- Return non-const access to the faceZones
515  {
516  return faceZones_;
517  }
518 
519  //- Return non-const access to the cellZones
521  {
522  return cellZones_;
523  }
524 
525  //- Add boundary patches
526  void addPatches
527  (
528  const List<polyPatch*>&,
529  const bool validBoundary = true
530  );
531 
532  //- Add mesh zones
533  void addZones
534  (
535  const List<pointZone*>& pz,
536  const List<faceZone*>& fz,
537  const List<cellZone*>& cz
538  );
539 
540  //- Add/insert single patch. If validBoundary the new situation
541  // is consistent across processors.
542  virtual void addPatch
543  (
544  const label insertPatchi,
545  const polyPatch& patch,
546  const dictionary& patchFieldDict,
547  const word& defaultPatchFieldType,
548  const bool validBoundary
549  );
550 
551  //- Reorder and trim existing patches. If validBoundary the new
552  // situation is consistent across processors
553  virtual void reorderPatches
554  (
555  const labelUList& newToOld,
556  const bool validBoundary
557  );
558 
559  //- Update the mesh based on the mesh files saved in
560  // time directories
562 
563  //- Update topology using the given map
564  virtual void topoChange(const polyTopoChangeMap&);
565 
566  //- Update from another mesh using the given map
567  virtual void mapMesh(const polyMeshMap&);
568 
569  //- Redistribute or update using the given distribution map
570  virtual void distribute(const polyDistributionMap& map);
571 
572  //- Remove boundary patches
573  void removeBoundary();
574 
575  //- Reset mesh primitive data. Assumes all patch info correct
576  // (so does e.g. parallel communication). If not use
577  // validBoundary=false
578  void resetPrimitives
579  (
580  pointField&& points,
581  faceList&& faces,
582  labelList&& owner,
583  labelList&& neighbour,
584  const labelList& patchSizes,
585  const labelList& patchStarts,
586  const bool validBoundary = true
587  );
588 
589  //- Swap mesh
590  // For run-time mesh replacement and mesh to mesh mapping
591  void swap(polyMesh&);
592 
593 
594  // Storage management
595 
596  //- Clear geometry
597  void clearGeom();
598 
599  //- Clear addressing
600  void clearAddressing(const bool isMeshUpdate = false);
601 
602  //- Clear all geometry and addressing unnecessary for CFD
603  void clearOut();
604 
605  //- Clear primitive data (points, faces and cells)
606  void clearPrimitives();
607 
608  //- Clear tet base points
609  void clearTetBasePtIs();
610 
611  //- Clear cell tree data
612  void clearCellTree();
613 
614  //- Remove all files from mesh instance
615  void removeFiles(const fileName& instanceDir) const;
616 
617  //- Remove all files from mesh instance()
618  void removeFiles() const;
619 
620 
621  // Geometric checks. Selectively override primitiveMesh functionality.
622 
623  //- Check non-orthogonality
624  virtual bool checkFaceOrthogonality
625  (
626  const bool report = false,
627  labelHashSet* setPtr = nullptr
628  ) const;
629 
630  //- Check face skewness
631  virtual bool checkFaceSkewness
632  (
633  const bool report = false,
634  labelHashSet* setPtr = nullptr
635  ) const;
636 
637  //- Check edge alignment for 1D/2D cases
638  virtual bool checkEdgeAlignment
639  (
640  const bool report,
641  const Vector<label>& directions,
642  labelHashSet* setPtr
643  ) const;
644 
645  virtual bool checkCellDeterminant
646  (
647  const bool report,
648  labelHashSet* setPtr
649  ) const;
650 
651  //- Check for face weights
652  virtual bool checkFaceWeight
653  (
654  const bool report,
655  const scalar minWeight = 0.05,
656  labelHashSet* setPtr = nullptr
657  ) const;
658 
659  //- Check for neighbouring cell volumes
660  virtual bool checkVolRatio
661  (
662  const bool report,
663  const scalar minRatio = 0.01,
664  labelHashSet* setPtr = nullptr
665  ) const;
666 
667 
668  // Position search functions
669 
670  //- Find the cell, tetFacei and tetPti for point p
671  void findCellFacePt
672  (
673  const point& p,
674  label& celli,
675  label& tetFacei,
676  label& tetPti
677  ) const;
678 
679  //- Find the tetFacei and tetPti for point p in celli.
680  // tetFacei and tetPtI are set to -1 if not found
681  void findTetFacePt
682  (
683  const label celli,
684  const point& p,
685  label& tetFacei,
686  label& tetPti
687  ) const;
688 
689  //- Test if point p is in the celli
690  bool pointInCell
691  (
692  const point& p,
693  label celli,
695  ) const;
696 
697  //- Find cell enclosing this location and return index
698  // If not found -1 is returned
700  (
701  const point& p,
703  ) const;
704 
705 
706  // Write
707 
708  //- Write the underlying polyMesh
709  virtual bool writeObject
710  (
714  const bool write = true
715  ) const;
716 
717 
718  // Member Operators
719 
720  //- Disallow default bitwise assignment
721  void operator=(const polyMesh&) = delete;
722 };
723 
724 
725 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
726 
727 } // End namespace Foam
728 
729 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
730 
731 #endif
732 
733 // ************************************************************************* //
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
writeOption
Enumeration defining the write options.
Definition: IOobject.H:126
Version number type.
Definition: IOstream.H:97
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:87
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:194
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:67
A class for handling file names.
Definition: fileName.H:82
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Non-pointer based hierarchical recursive searching.
Definition: indexedOctree.H:72
Registry of regIOobjects.
Foam::polyBoundaryMesh.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual ~polyMesh()
Destructor.
Definition: polyMesh.C:992
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1624
bool topoChanged_
Has the mesh topology changed.
Definition: polyMesh.H:251
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:1026
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1055
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:269
polyMesh Mesh
Definition: polyMesh.H:258
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:445
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:1014
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1492
void clearCellTree()
Clear cell tree data.
virtual void addPatch(const label insertPatchi, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add/insert single patch. If validBoundary the new situation.
Definition: polyMesh.C:1279
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:101
bool moving_
Member data pending transfer to fvMesh.
Definition: polyMesh.H:248
IOobject::writeOption facesWriteOpt() const
Return the points write option.
Definition: polyMesh.C:1038
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1374
void clearGeom()
Clear geometry.
Definition: polyMeshClear.C:53
virtual bool checkFaceSkewness(const bool report=false, labelHashSet *setPtr=nullptr) const
Check face skewness.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
@ TOPO_PATCH_CHANGE
Definition: polyMesh.H:95
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:1001
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1775
void resetPrimitives(pointField &&points, faceList &&faces, labelList &&owner, labelList &&neighbour, const labelList &patchSizes, const labelList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:701
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1387
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:1078
void operator=(const polyMesh &)=delete
Disallow default bitwise assignment.
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:1649
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
TypeName("polyMesh")
Runtime type information.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1563
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:266
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:466
bool changing() const
Is mesh changing.
Definition: polyMesh.H:488
void clearPrimitives()
Clear primitive data (points, faces and cells)
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1665
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:1072
polyBoundaryMesh BoundaryMesh
Definition: polyMesh.H:259
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1399
void swap(polyMesh &)
Swap mesh.
Definition: polyMesh.C:802
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:1020
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1581
polyMesh(const IOobject &io)
Construct from IOobject.
Definition: polyMesh.C:163
const meshPointZones & pointZones() const
Return point zones.
Definition: polyMesh.H:439
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1555
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write the underlying polyMesh.
Definition: polyMeshIO.C:519
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:99
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1393
virtual bool checkFaceOrthogonality(const bool report=false, labelHashSet *setPtr=nullptr) const
Check non-orthogonality.
Definition: polyMeshCheck.C:34
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:1139
bool topoChanged() const
Has the mesh topology changed this time-step.
Definition: polyMesh.H:481
bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1437
void clearTetBasePtIs()
Clear tet base points.
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
void setPointsInstance(const fileName &)
Set the instance for the points files.
Definition: polyMeshIO.C:58
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:1183
virtual bool checkEdgeAlignment(const bool report, const Vector< label > &directions, labelHashSet *setPtr) const
Check edge alignment for 1D/2D cases.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1361
void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1443
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:1112
IOobject::writeOption pointsWriteOpt() const
Return the points write option.
Definition: polyMesh.C:1032
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1617
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
virtual const pointField & oldCellCentres() const
Return old cell centres for mesh motion.
Definition: polyMesh.C:1417
virtual bool checkFaceWeight(const bool report, const scalar minWeight=0.05, labelHashSet *setPtr=nullptr) const
Check for face weights.
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: polyMesh.C:1241
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
virtual bool checkCellDeterminant(const bool report, labelHashSet *setPtr) const
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:409
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:451
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
virtual bool checkVolRatio(const bool report, const scalar minRatio=0.01, labelHashSet *setPtr=nullptr) const
Check for neighbouring cell volumes.
virtual void setPoints(const pointField &)
Reset the points.
Definition: polyMesh.C:1449
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:1061
bool moving() const
Is mesh moving.
Definition: polyMesh.H:475
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:1044
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
const labelListList & pointCells() const
void clearAddressing()
Clear topological data.
label nFaces() const
const cellList & cells() const
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::meshCellZones.
Foam::meshFaceZones.
Foam::meshPointZones.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
volScalarField & p