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  //- Boundary mesh
136  mutable polyBoundaryMesh boundary_;
137 
138  //- Mesh bounding-box.
139  // Created from points on construction, updated when the mesh moves
140  boundBox bounds_;
141 
142  //- Communicator used for parallel communication
143  label comm_;
144 
145  //- Vector of non-constrained directions in mesh
146  // defined according to the presence of empty and wedge patches
147  mutable Vector<label> geometricD_;
148 
149  //- Vector of valid directions in mesh
150  // defined according to the presence of empty patches
151  mutable Vector<label> solutionD_;
152 
153  //- Base point for face decomposition into tets
154  mutable autoPtr<labelIOList> tetBasePtIsPtr_;
155 
156  //- Search tree to allow spatial cell searching
157  mutable autoPtr<indexedOctree<treeDataCell>> cellTreePtr_;
158 
159 
160  // Zoning information
161 
162  //- Point zones
163  meshPointZones pointZones_;
164 
165  //- Face zones
166  meshFaceZones faceZones_;
167 
168  //- Cell zones
169  meshCellZones cellZones_;
170 
171 
172  //- Parallel info
173  mutable autoPtr<globalMeshData> globalMeshDataPtr_;
174 
175 
176  // Mesh motion related dat
177 
178  //- Current time index for mesh motion
179  mutable label curMotionTimeIndex_;
180 
181  //- Old points (for the last mesh motion)
182  mutable autoPtr<pointField> oldPointsPtr_;
183 
184  //- Old cell centres (for the last mesh motion)
185  mutable autoPtr<pointField> oldCellCentresPtr_;
186 
187  //- Whether or not to store the old cell centres
188  mutable bool storeOldCellCentres_;
189 
190 
191  // Private Member Functions
192 
193  //- Initialise the polyMesh from the primitive data
194  void initMesh();
195 
196  //- Initialise the polyMesh from the given set of cells
197  void initMesh(cellList& c);
198 
199  //- Calculate the valid directions in the mesh from the boundaries
200  void calcDirections() const;
201 
202  //- Calculate the cell shapes from the primitive
203  // polyhedral information
204  void calcCellShapes() const;
205 
206  //- Read and return the tetBasePtIs
207  autoPtr<labelIOList> readTetBasePtIs() const;
208 
209  //- Set the write option of the points
210  void setPointsWrite(const IOobject::writeOption wo);
211 
212  //- Set the write option of the topology
213  void setTopologyWrite(const IOobject::writeOption wo);
214 
215 
216  // Helper functions for constructor from cell shapes
217 
218  labelListList cellShapePointCells(const cellShapeList&) const;
219 
220  labelList facePatchFaceCells
221  (
222  const faceList& patchFaces,
223  const labelListList& pointCells,
224  const faceListList& cellsFaceShapes,
225  const label patchID
226  ) const;
227 
228  void setTopology
229  (
230  const cellShapeList& cellsAsShapes,
231  const faceListList& boundaryFaces,
232  const wordList& boundaryPatchNames,
233  labelList& patchSizes,
234  labelList& patchStarts,
235  label& defaultPatchStart,
236  label& nFaces,
237  cellList& cells
238  );
239 
240 
241 protected:
242 
243  // Protected Member Data
244 
245  //- Member data pending transfer to fvMesh
246 
247  //- Is the mesh moving
248  bool moving_;
249 
250  //- Has the mesh topology changed
251  bool topoChanged_;
252 
253 
254 public:
255 
256  // Public Typedefs
257 
258  typedef polyMesh Mesh;
260 
261 
262  //- Runtime type information
263  TypeName("polyMesh");
264 
265 
266  // Static data
267 
268  //- Return the default region name
269  static word defaultRegion;
270 
271  //- Return the mesh sub-directory name (usually "polyMesh")
272  static word meshSubDir;
273 
274 
275  // Constructors
276 
277  //- Construct from IOobject
278  explicit polyMesh(const IOobject& io);
279 
280  //- Move construct from IOobject or from components.
281  // Boundary is added using addPatches() member function
282  polyMesh
283  (
284  const IOobject& io,
285  pointField&& points,
286  faceList&& faces,
287  labelList&& owner,
288  labelList&& neighbour,
289  const bool syncPar = true
290  );
291 
292  //- Move construct without boundary with cells rather than
293  // owner/neighbour.
294  // Boundary is added using addPatches() member function
295  polyMesh
296  (
297  const IOobject& io,
298  pointField&& points,
299  faceList&& faces,
300  cellList&& cells,
301  const bool syncPar = true
302  );
303 
304  //- Move construct from cell shapes
305  polyMesh
306  (
307  const IOobject& io,
308  pointField&& points,
309  const cellShapeList& shapes,
310  const faceListList& boundaryFaces,
311  const wordList& boundaryPatchNames,
312  const wordList& boundaryPatchTypes,
313  const word& defaultBoundaryPatchName,
314  const word& defaultBoundaryPatchType,
315  const wordList& boundaryPatchPhysicalTypes,
316  const bool syncPar = true
317  );
318 
319  //- Move construct from cell shapes with patch information in dictionary
320  // format.
321  polyMesh
322  (
323  const IOobject& io,
324  pointField&& points,
325  const cellShapeList& shapes,
326  const faceListList& boundaryFaces,
327  const wordList& boundaryPatchNames,
328  const PtrList<dictionary>& boundaryDicts,
329  const word& defaultBoundaryPatchName,
330  const word& defaultBoundaryPatchType,
331  const bool syncPar = true
332  );
333 
334  //- Move constructor
335  polyMesh(polyMesh&&);
336 
337  //- Disallow default bitwise copy construction
338  polyMesh(const polyMesh&) = delete;
339 
340 
341  //- Destructor
342  virtual ~polyMesh();
343 
344 
345  // Member Functions
346 
347  // Database
348 
349  //- Override the objectRegistry dbDir for a single-region case
350  virtual const fileName& dbDir() const;
351 
352  //- Return the local mesh directory (dbDir()/meshSubDir)
353  fileName meshDir() const;
354 
355  //- Return the current instance directory for points
356  // Used in the construction of geometric mesh data dependent
357  // on points
358  const fileName& pointsInstance() const;
359 
360  //- Return the current instance directory for faces
361  const fileName& facesInstance() const;
362 
363  //- Return the points write option
365 
366  //- Return the points write option
368 
369  //- Set the instance for the points files
370  void setPointsInstance(const fileName&);
371 
372  //- Set the instance for mesh files
373  void setInstance(const fileName&);
374 
375 
376  // Access
377 
378  // Primitive mesh data
379 
380  //- Return raw points
381  virtual const pointField& points() const;
382 
383  //- Return raw faces
384  virtual const faceList& faces() const;
385 
386  //- Return face owner
387  virtual const labelList& faceOwner() const;
388 
389  //- Return face neighbour
390  virtual const labelList& faceNeighbour() const;
391 
392  //- Return old points for mesh motion
393  virtual const pointField& oldPoints() const;
394 
395  //- Return old cell centres for mesh motion
396  virtual const pointField& oldCellCentres() const;
397 
398 
399  //- Return true if io is up-to-date with points
400  bool upToDatePoints(const regIOobject& io) const;
401 
402  //- Set io to be up-to-date with points
403  void setUpToDatePoints(regIOobject& io) const;
404 
405  //- Return boundary mesh
406  const polyBoundaryMesh& boundaryMesh() const
407  {
408  return boundary_;
409  }
410 
411  //- Return mesh bounding box
412  const boundBox& bounds() const
413  {
414  return bounds_;
415  }
416 
417  //- Return the vector of geometric directions in mesh.
418  // Defined according to the presence of empty and wedge patches.
419  // 1 indicates unconstrained direction and -1 a constrained
420  // direction.
421  const Vector<label>& geometricD() const;
422 
423  //- Return the number of valid geometric dimensions in the mesh
424  label nGeometricD() const;
425 
426  //- Return the vector of solved-for directions in mesh.
427  // Differs from geometricD in that it includes for wedge cases
428  // the circumferential direction in case of swirl.
429  // 1 indicates valid direction and -1 an invalid direction.
430  const Vector<label>& solutionD() const;
431 
432  //- Return the number of valid solved-for dimensions in the mesh
433  label nSolutionD() const;
434 
435  //- Return the tetBasePtIs
436  const labelIOList& tetBasePtIs() const;
437 
438  //- Return the cell search tree
439  const indexedOctree<treeDataCell>& cellTree() const;
440 
441  //- Return point zones
442  const meshPointZones& pointZones() const
443  {
444  return pointZones_;
445  }
446 
447  //- Return face zones
448  const meshFaceZones& faceZones() const
449  {
450  return faceZones_;
451  }
452 
453  //- Return cell zones
454  const meshCellZones& cellZones() const
455  {
456  return cellZones_;
457  }
458 
459  //- Return parallel info
460  const globalMeshData& globalData() const;
461 
462  //- Return communicator used for parallel communication
463  label comm() const;
464 
465  //- Return communicator used for parallel communication
466  label& comm();
467 
468  //- Return the object registry
469  const objectRegistry& thisDb() const
470  {
471  return *this;
472  }
473 
474 
475  // Mesh motion
476 
477  //- Is mesh moving
478  bool moving() const
479  {
480  return moving_;
481  }
482 
483  //- Has the mesh topology changed this time-step
484  bool topoChanged() const
485  {
486  return topoChanged_;
487  }
488 
489  //- Is mesh changing
490  // Moving or mesh topology changed this time-step)
491  bool changing() const
492  {
493  return moving() || topoChanged();
494  }
495 
496  //- Reset the points
497  // without storing old points or returning swept volumes
498  virtual void setPoints(const pointField&);
499 
500  //- Move points, returns volumes swept by faces in motion
501  virtual tmp<scalarField> movePoints(const pointField&);
502 
503  //- Reset motion
504  void resetMotion() const;
505 
506 
507  // Topological change
508 
509  //- Return non-const access to the pointZones
511  {
512  return pointZones_;
513  }
514 
515  //- Return non-const access to the faceZones
517  {
518  return faceZones_;
519  }
520 
521  //- Return non-const access to the cellZones
523  {
524  return cellZones_;
525  }
526 
527  //- Add boundary patches
528  void addPatches
529  (
530  const List<polyPatch*>&,
531  const bool validBoundary = true
532  );
533 
534  //- Add mesh zones
535  void addZones
536  (
537  const List<pointZone*>& pz,
538  const List<faceZone*>& fz,
539  const List<cellZone*>& cz
540  );
541 
542  //- Add/insert single patch. If validBoundary the new situation
543  // is consistent across processors.
544  virtual void addPatch
545  (
546  const label insertPatchi,
547  const polyPatch& patch,
548  const dictionary& patchFieldDict,
549  const word& defaultPatchFieldType,
550  const bool validBoundary
551  );
552 
553  //- Reorder and trim existing patches. If validBoundary the new
554  // situation is consistent across processors
555  virtual void reorderPatches
556  (
557  const labelUList& newToOld,
558  const bool validBoundary
559  );
560 
561  //- Update the mesh based on the mesh files saved in
562  // time directories
564 
565  //- Update topology using the given map
566  virtual void topoChange(const polyTopoChangeMap&);
567 
568  //- Update from another mesh using the given map
569  virtual void mapMesh(const polyMeshMap&);
570 
571  //- Redistribute or update using the given distribution map
572  virtual void distribute(const polyDistributionMap& map);
573 
574  //- Remove boundary patches
575  void removeBoundary();
576 
577  //- Reset mesh primitive data. Assumes all patch info correct
578  // (so does e.g. parallel communication). If not use
579  // validBoundary=false
580  void resetPrimitives
581  (
582  pointField&& points,
583  faceList&& faces,
584  labelList&& owner,
585  labelList&& neighbour,
586  const labelList& patchSizes,
587  const labelList& patchStarts,
588  const bool validBoundary = true
589  );
590 
591  //- Swap mesh
592  // For run-time mesh replacement and mesh to mesh mapping
593  void swap(polyMesh&);
594 
595 
596  // Storage management
597 
598  //- Clear geometry
599  void clearGeom();
600 
601  //- Clear addressing
602  void clearAddressing(const bool isMeshUpdate = false);
603 
604  //- Clear all geometry and addressing unnecessary for CFD
605  void clearOut();
606 
607  //- Clear primitive data (points, faces and cells)
608  void clearPrimitives();
609 
610  //- Clear tet base points
611  void clearTetBasePtIs();
612 
613  //- Clear cell tree data
614  void clearCellTree();
615 
616  //- Remove all files from mesh instance
617  void removeFiles(const fileName& instanceDir) const;
618 
619  //- Remove all files from mesh instance()
620  void removeFiles() const;
621 
622 
623  // Geometric checks. Selectively override primitiveMesh functionality.
624 
625  //- Check non-orthogonality
626  virtual bool checkFaceOrthogonality
627  (
628  const bool report = false,
629  labelHashSet* setPtr = nullptr
630  ) const;
631 
632  //- Check face skewness
633  virtual bool checkFaceSkewness
634  (
635  const bool report = false,
636  labelHashSet* setPtr = nullptr
637  ) const;
638 
639  //- Check edge alignment for 1D/2D cases
640  virtual bool checkEdgeAlignment
641  (
642  const bool report,
643  const Vector<label>& directions,
644  labelHashSet* setPtr
645  ) const;
646 
647  virtual bool checkCellDeterminant
648  (
649  const bool report,
650  labelHashSet* setPtr
651  ) const;
652 
653  //- Check for face weights
654  virtual bool checkFaceWeight
655  (
656  const bool report,
657  const scalar minWeight = 0.05,
658  labelHashSet* setPtr = nullptr
659  ) const;
660 
661  //- Check for neighbouring cell volumes
662  virtual bool checkVolRatio
663  (
664  const bool report,
665  const scalar minRatio = 0.01,
666  labelHashSet* setPtr = nullptr
667  ) const;
668 
669 
670  // Position search functions
671 
672  //- Find the cell, tetFacei and tetPti for point p
673  void findCellFacePt
674  (
675  const point& p,
676  label& celli,
677  label& tetFacei,
678  label& tetPti
679  ) const;
680 
681  //- Find the tetFacei and tetPti for point p in celli.
682  // tetFacei and tetPtI are set to -1 if not found
683  void findTetFacePt
684  (
685  const label celli,
686  const point& p,
687  label& tetFacei,
688  label& tetPti
689  ) const;
690 
691  //- Test if point p is in the celli
692  bool pointInCell
693  (
694  const point& p,
695  label celli,
697  ) const;
698 
699  //- Find cell enclosing this location and return index
700  // If not found -1 is returned
702  (
703  const point& p,
705  ) const;
706 
707 
708  // Write
709 
710  //- Write the underlying polyMesh
711  virtual bool writeObject
712  (
716  const bool write = true
717  ) const;
718 
719 
720  // Member Operators
721 
722  //- Disallow default bitwise assignment
723  void operator=(const polyMesh&) = delete;
724 };
725 
726 
727 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
728 
729 } // End namespace Foam
730 
731 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
732 
733 #endif
734 
735 // ************************************************************************* //
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:991
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1623
bool topoChanged_
Has the mesh topology changed.
Definition: polyMesh.H:250
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:1025
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1054
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:268
polyMesh Mesh
Definition: polyMesh.H:257
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:447
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:1013
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1491
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:1278
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:101
bool moving_
Member data pending transfer to fvMesh.
Definition: polyMesh.H:247
IOobject::writeOption facesWriteOpt() const
Return the points write option.
Definition: polyMesh.C:1037
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1373
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:1000
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1774
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:700
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1386
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:1077
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:1648
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:1562
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:468
bool changing() const
Is mesh changing.
Definition: polyMesh.H:490
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:1664
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:1071
polyBoundaryMesh BoundaryMesh
Definition: polyMesh.H:258
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1398
void swap(polyMesh &)
Swap mesh.
Definition: polyMesh.C:801
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:1019
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1580
polyMesh(const IOobject &io)
Construct from IOobject.
Definition: polyMesh.C:162
const meshPointZones & pointZones() const
Return point zones.
Definition: polyMesh.H:441
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1554
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:405
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1392
virtual bool checkFaceOrthogonality(const bool report=false, labelHashSet *setPtr=nullptr) const
Check non-orthogonality.
Definition: polyMeshCheck.C:42
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:1138
bool topoChanged() const
Has the mesh topology changed this time-step.
Definition: polyMesh.H:483
bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1436
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:1182
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:1360
void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1442
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:1111
IOobject::writeOption pointsWriteOpt() const
Return the points write option.
Definition: polyMesh.C:1031
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1616
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
virtual const pointField & oldCellCentres() const
Return old cell centres for mesh motion.
Definition: polyMesh.C:1416
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:1240
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:271
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:411
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:453
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:1448
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:1060
bool moving() const
Is mesh moving.
Definition: polyMesh.H:477
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:1043
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