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-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::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 "pointZoneMesh.H"
59 #include "faceZoneMesh.H"
60 #include "cellZoneMesh.H"
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 namespace Foam
65 {
66 
67 // Forward Declarations
68 class globalMeshData;
69 class mapPolyMesh;
70 class polyMeshTetDecomposition;
71 class treeDataCell;
72 template<class Type> class indexedOctree;
73 
74 /*---------------------------------------------------------------------------*\
75  Class polyMesh Declaration
76 \*---------------------------------------------------------------------------*/
77 
78 class polyMesh
79 :
80  public objectRegistry,
81  public primitiveMesh
82 {
83 public:
84 
85  // Public Data
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
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
105  FACE_CENTRE_TRIS, //- Faces decomposed into triangles
106  // using face-centre
107 
108  FACE_DIAG_TRIS, //- Faces decomposed into triangles diagonally
110  CELL_TETS //- Cell decomposed into tets
111  };
112 
113 
114 private:
115 
116  // Private 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  pointZoneMesh pointZones_;
165 
166  //- Face zones
167  faceZoneMesh faceZones_;
168 
169  //- Cell zones
170  cellZoneMesh cellZones_;
171 
172 
173  //- Parallel info
174  mutable autoPtr<globalMeshData> globalMeshDataPtr_;
175 
176 
177  // Mesh motion related data
178 
179  //- Is the mesh moving
180  bool moving_;
181 
182  //- Is the mesh topology changing
183  bool topoChanging_;
184 
185  //- Store old cell centres?
186  mutable bool storeOldCellCentres_;
187 
188  //- Current time index for mesh motion
189  mutable label curMotionTimeIndex_;
190 
191  //- Old points (for the last mesh motion)
192  mutable autoPtr<pointField> oldPointsPtr_;
193 
194  //- Old cell centres (for the last mesh motion)
195  mutable autoPtr<pointField> oldCellCentresPtr_;
196 
197 
198 
199  // Private Member Functions
200 
201  //- No copy construct
202  polyMesh(const polyMesh&) = delete;
203 
204  //- No copy assignment
205  void operator=(const polyMesh&) = delete;
206 
207  //- Initialise the polyMesh from the primitive data
208  void initMesh();
209 
210  //- Initialise the polyMesh from the given set of cells
211  void initMesh(cellList& c);
212 
213  //- Calculate the valid directions in the mesh from the boundaries
214  void calcDirections() const;
215 
216  //- Calculate the cell shapes from the primitive
217  // polyhedral information
218  void calcCellShapes() const;
219 
220  //- Read and return the tetBasePtIs
221  autoPtr<labelIOList> readTetBasePtIs() const;
222 
223 
224  // Helper functions for constructor from cell shapes
225 
226  labelListList cellShapePointCells(const cellShapeList&) const;
227 
228  labelList facePatchFaceCells
229  (
230  const faceList& patchFaces,
231  const labelListList& pointCells,
232  const faceListList& cellsFaceShapes,
233  const label patchID
234  ) const;
235 
236  void setTopology
237  (
238  const cellShapeList& cellsAsShapes,
239  const faceListList& boundaryFaces,
240  const wordList& boundaryPatchNames,
241  labelList& patchSizes,
242  labelList& patchStarts,
243  label& defaultPatchStart,
244  label& nFaces,
245  cellList& cells
246  );
247 
248 
249  // Geometry checks
250 
251  //- Check non-orthogonality
252  bool checkFaceOrthogonality
253  (
254  const vectorField& fAreas,
255  const vectorField& cellCtrs,
256  const bool report,
257  const bool detailedReport,
258  labelHashSet* setPtr
259  ) const;
260 
261  //- Check face skewness
262  bool checkFaceSkewness
263  (
264  const pointField& points,
265  const vectorField& fCtrs,
266  const vectorField& fAreas,
267  const vectorField& cellCtrs,
268  const bool report,
269  const bool detailedReport,
270  labelHashSet* setPtr
271  ) const;
272 
273  bool checkEdgeAlignment
274  (
275  const pointField& p,
276  const bool report,
277  const Vector<label>& directions,
278  labelHashSet* setPtr
279  ) const;
280 
281  bool checkCellDeterminant
282  (
283  const vectorField& faceAreas,
284  const bool report,
285  labelHashSet* setPtr,
286  const Vector<label>& meshD
287  ) const;
288 
289  bool checkFaceWeight
290  (
291  const vectorField& fCtrs,
292  const vectorField& fAreas,
293  const vectorField& cellCtrs,
294  const bool report,
295  const scalar minWeight,
296  labelHashSet* setPtr
297  ) const;
298 
299  bool checkVolRatio
300  (
301  const scalarField& cellVols,
302  const bool report,
303  const scalar minRatio,
304  labelHashSet* setPtr
305  ) const;
306 
307 public:
308 
309  // Public Typedefs
310 
311  //- The mesh type
312  typedef polyMesh Mesh;
313 
314  //- The boundary type associated with the mesh
316 
317 
318  //- Runtime type information
319  TypeName("polyMesh");
320 
321  //- Return the default region name
322  static word defaultRegion;
323 
324  //- Return the mesh sub-directory name (usually "polyMesh")
325  static word meshSubDir;
326 
327 
328  // Constructors
329 
330  //- Read construct from IOobject
331  explicit polyMesh(const IOobject& io, const bool doInit = true);
332 
333  //- Construct from IOobject or as zero-sized mesh
334  // Boundary is added using addPatches() member function
335  polyMesh(const IOobject& io, const zero, const bool syncPar=true);
336 
337  //- Construct from IOobject and components.
338  // Boundary is added using addPatches() member function
339  polyMesh
340  (
341  const IOobject& io,
342  pointField&& points,
343  faceList&& faces,
344  labelList&& owner,
345  labelList&& neighbour,
346  const bool syncPar = true
347  );
348 
349  //- Construct without boundary with cells rather than owner/neighbour.
350  // Boundary is added using addPatches() member function
351  polyMesh
352  (
353  const IOobject& io,
354  pointField&& points,
355  faceList&& faces,
356  cellList&& cells,
357  const bool syncPar = true
358  );
359 
360  //- Construct from cell shapes
361  polyMesh
362  (
363  const IOobject& io,
364  pointField&& points,
365  const cellShapeList& shapes,
366  const faceListList& boundaryFaces,
367  const wordList& boundaryPatchNames,
368  const wordList& boundaryPatchTypes,
369  const word& defaultBoundaryPatchName,
370  const word& defaultBoundaryPatchType,
371  const wordList& boundaryPatchPhysicalTypes,
372  const bool syncPar = true
373  );
374 
375  //- Construct from cell shapes, with patch information in dictionary
376  //- format.
377  polyMesh
378  (
379  const IOobject& io,
380  pointField&& points,
381  const cellShapeList& shapes,
382  const faceListList& boundaryFaces,
383  const wordList& boundaryPatchNames,
385  const word& defaultBoundaryPatchName,
386  const word& defaultBoundaryPatchType,
387  const bool syncPar = true
388  );
389 
390 
391  //- Destructor
392  virtual ~polyMesh();
393 
394 
395  // Member Functions
396 
397  // Database
398 
399  //- Override the objectRegistry dbDir for a single-region case
400  virtual const fileName& dbDir() const;
401 
402  //- Return the local mesh directory (dbDir()/meshSubDir)
403  fileName meshDir() const;
404 
405  //- Return the current instance directory for points
406  // Used in the construction of geometric mesh data dependent
407  // on points
408  const fileName& pointsInstance() const;
409 
410  //- Return the current instance directory for faces
411  const fileName& facesInstance() const;
412 
413  //- Set the instance for mesh files
414  void setInstance
415  (
416  const fileName& instance,
418  );
419 
420 
421  // Regions
422 
423  //- The mesh region name or word::null if polyMesh::defaultRegion
424  static const word& regionName(const word& region);
425 
426  //- The mesh region name or word::null if polyMesh::defaultRegion
427  const word& regionName() const;
428 
429 
430  // Access
431 
432  //- Return raw points
433  virtual const pointField& points() const;
434 
435  //- Return true if io is up-to-date with points
436  virtual bool upToDatePoints(const regIOobject& io) const;
437 
438  //- Set io to be up-to-date with points
439  virtual void setUpToDatePoints(regIOobject& io) const;
440 
441  //- Return raw faces
442  virtual const faceList& faces() const;
443 
444  //- Return face owner
445  virtual const labelList& faceOwner() const;
446 
447  //- Return face neighbour
448  virtual const labelList& faceNeighbour() const;
449 
450  //- Return old points (mesh motion)
451  virtual const pointField& oldPoints() const;
452 
453  //- Return old cellCentres (mesh motion)
454  virtual const pointField& oldCellCentres() const;
455 
456  //- Return boundary mesh
457  const polyBoundaryMesh& boundaryMesh() const noexcept
458  {
459  return boundary_;
460  }
461 
462  //- Return mesh bounding box
463  const boundBox& bounds() const noexcept
464  {
465  return bounds_;
466  }
467 
468  //- Return the vector of geometric directions in mesh.
469  // Defined according to the presence of empty and wedge patches.
470  // 1 indicates unconstrained direction and -1 a constrained
471  // direction.
472  const Vector<label>& geometricD() const;
473 
474  //- Return the number of valid geometric dimensions in the mesh
475  label nGeometricD() const;
476 
477  //- Return the vector of solved-for directions in mesh.
478  // Differs from geometricD in that it includes for wedge cases
479  // the circumferential direction in case of swirl.
480  // 1 indicates valid direction and -1 an invalid direction.
481  const Vector<label>& solutionD() const;
482 
483  //- Return the number of valid solved-for dimensions in the mesh
484  label nSolutionD() const;
485 
486  //- Return the tetBasePtIs
487  const labelIOList& tetBasePtIs() const;
488 
489  //- Return the cell search tree
490  const indexedOctree<treeDataCell>& cellTree() const;
491 
492  //- Return point zone mesh
493  const pointZoneMesh& pointZones() const noexcept
494  {
495  return pointZones_;
496  }
497 
498  //- Return face zone mesh
499  const faceZoneMesh& faceZones() const noexcept
500  {
501  return faceZones_;
502  }
503 
504  //- Return cell zone mesh
505  const cellZoneMesh& cellZones() const noexcept
506  {
507  return cellZones_;
508  }
509 
510  //- Return parallel info
511  const globalMeshData& globalData() const;
512 
513  //- Return communicator used for parallel communication
514  label comm() const noexcept;
515 
516  //- Return communicator used for parallel communication
517  label& comm() noexcept;
518 
519  //- Return the object registry
520  const objectRegistry& thisDb() const noexcept
521  {
522  return *this;
523  }
524 
525 
526  // Mesh motion
527 
528  //- Is mesh dynamic
529  virtual bool dynamic() const
530  {
531  return false;
532  }
533 
534  //- Is mesh moving
535  bool moving() const noexcept
536  {
537  return moving_;
538  }
539 
540  //- Set the mesh to be moving
541  bool moving(const bool on) noexcept
542  {
543  bool old(moving_);
544  moving_ = on;
545  return old;
546  }
547 
548  //- Is mesh topology changing
549  bool topoChanging() const noexcept
550  {
551  return topoChanging_;
552  }
553 
554  //- Set the mesh topology to be changing
555  bool topoChanging(const bool on) noexcept
556  {
557  bool old(topoChanging_);
558  topoChanging_ = on;
559  return old;
560  }
561 
562  //- Is mesh changing (topology changing and/or moving)
563  bool changing() const noexcept
564  {
565  return (moving() || topoChanging());
566  }
567 
568  //- Move points
569  virtual void movePoints(const pointField&);
570 
571  //- Reset motion
572  void resetMotion() const;
573 
574 
575  // Topological change
576 
577  //- Return non-const access to the pointZones
579  {
580  return pointZones_;
581  }
582 
583  //- Return non-const access to the faceZones
585  {
586  return faceZones_;
587  }
588 
589  //- Return non-const access to the cellZones
591  {
592  return cellZones_;
593  }
594 
595  //- Add boundary patches
596  void addPatches
597  (
598  polyPatchList& plist,
599  const bool validBoundary = true
600  );
601 
602  //- Add boundary patches
603  void addPatches
604  (
605  const List<polyPatch*>& p,
606  const bool validBoundary = true
607  );
608 
609  //- Add mesh zones
610  void addZones
611  (
612  PtrList<pointZone>&& pz,
613  PtrList<faceZone>&& fz,
614  PtrList<cellZone>&& cz
615  );
616 
617  //- Add mesh zones
618  void addZones
619  (
620  const List<pointZone*>& pz,
621  const List<faceZone*>& fz,
622  const List<cellZone*>& cz
623  );
624 
625  //- Initialise all non-demand-driven data
626  virtual bool init(const bool doInit);
627 
628  //- Update the mesh based on the mesh files saved in
629  // time directories
630  virtual readUpdateState readUpdate();
631 
632  //- Update the mesh corresponding to given map
633  virtual void updateMesh(const mapPolyMesh& mpm);
634 
635  //- Remove boundary patches
636  void removeBoundary();
637 
638  //- Reset mesh primitive data. Assumes all patch info correct
639  // (so does e.g. parallel communication). If not use
640  // validBoundary=false
641  //
642  // \note Only autoPtr parameters that test as valid() are used
643  // for resetting, otherwise the existing entries are left
644  // untouched.
645  void resetPrimitives
646  (
649  autoPtr<labelList>&& owner,
650  autoPtr<labelList>&& neighbour,
651  const labelUList& patchSizes,
652  const labelUList& patchStarts,
653  const bool validBoundary = true
654  );
655 
656 
657  // Storage management
658 
659  //- Clear geometry
660  void clearGeom();
661 
662  //- Update geometry points; keep topology.
663  //- Optionally with new face decomposition
664  void updateGeomPoints
665  (
666  pointIOField&& newPoints,
667  autoPtr<labelIOList>& newTetBasePtIsPtr
668  );
669 
670  //- Clear addressing
671  void clearAddressing(const bool isMeshUpdate = false);
672 
673  //- Clear all geometry and addressing unnecessary for CFD
674  void clearOut();
675 
676  //- Clear primitive data (points, faces and cells)
678 
679  //- Clear tet base points
680  void clearTetBasePtIs();
681 
682  //- Clear cell tree data
683  void clearCellTree();
684 
685  //- Remove all files from mesh instance
686  void removeFiles(const fileName& instanceDir) const;
687 
688  //- Remove all files from mesh instance()
689  void removeFiles() const;
690 
691 
692  bool hasTetBasePtIs() const { return bool(tetBasePtIsPtr_); }
693 
694 
695  // Geometric checks. Selectively override primitiveMesh functionality.
697  //- Check non-orthogonality
698  virtual bool checkFaceOrthogonality
699  (
700  const bool report = false,
701  labelHashSet* setPtr = nullptr
702  ) const;
703 
704  //- Check face skewness
705  virtual bool checkFaceSkewness
706  (
707  const bool report = false,
708  labelHashSet* setPtr = nullptr
709  ) const;
710 
711  //- Check edge alignment for 1D/2D cases
712  virtual bool checkEdgeAlignment
713  (
714  const bool report,
715  const Vector<label>& directions,
716  labelHashSet* setPtr
717  ) const;
718 
719  virtual bool checkCellDeterminant
720  (
721  const bool report,
722  labelHashSet* setPtr
723  ) const;
724 
725  //- Check mesh motion for correctness given motion points
726  virtual bool checkMeshMotion
727  (
728  const pointField& newPoints,
729  const bool report = false,
730  const bool detailedReport = false
731  ) const;
733  //- Check for face weights
734  virtual bool checkFaceWeight
735  (
736  const bool report,
737  const scalar minWeight = 0.05,
738  labelHashSet* setPtr = nullptr
739  ) const;
740 
741  //- Check for neighbouring cell volumes
742  virtual bool checkVolRatio
743  (
744  const bool report,
745  const scalar minRatio = 0.01,
746  labelHashSet* setPtr = nullptr
747  ) const;
748 
749 
750  // Position search functions
751 
752  //- Find the cell, tetFacei and tetPti for point p
754  (
755  const point& p,
756  label& celli,
757  label& tetFacei,
758  label& tetPti
759  ) const;
760 
761  //- Find the tetFacei and tetPti for point p in celli.
762  // tetFacei and tetPtI are set to -1 if not found
763  void findTetFacePt
764  (
765  const label celli,
766  const point& p,
767  label& tetFacei,
768  label& tetPti
769  ) const;
770 
771  //- Test if point p is in the celli
772  bool pointInCell
773  (
774  const point& p,
775  label celli,
777  ) const;
778 
779  //- Find cell enclosing this location and return index
780  // If not found -1 is returned
781  label findCell
782  (
783  const point& p,
785  ) const;
786 };
787 
788 
789 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
790 
791 } // End namespace Foam
792 
793 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
794 
795 #endif
796 
797 // ************************************************************************* //
writeOption
Enumeration defining write preferences.
polyMesh Mesh
The mesh type.
Definition: polyMesh.H:381
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:359
A class for handling file names.
Definition: fileName.H:71
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:1163
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:854
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1117
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:894
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:842
virtual void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1092
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1292
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:877
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:402
Foam::pointZoneMesh.
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:710
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
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:825
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:1402
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:102
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:871
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1073
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:860
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:848
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:1354
void clearPrimitives()
Clear primitive data (points, faces and cells)
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:584
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:1361
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:397
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
Definition: polyMesh.C:1004
virtual const pointField & oldPoints() const
Return old points (mesh motion)
Definition: polyMesh.C:1123
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1111
label comm() const noexcept
Return communicator used for parallel communication.
Definition: polyMesh.C:1318
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1300
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1098
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:714
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:646
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:732
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:233
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:888
virtual bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1086
bool moving() const noexcept
Is mesh moving.
Definition: polyMesh.H:696
const pointZoneMesh & pointZones() const noexcept
Return point zone mesh.
Definition: polyMesh.H:638
TypeName("polyMesh")
Runtime type information.
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:934
virtual const pointField & oldCellCentres() const
Return old cellCentres (mesh motion)
Definition: polyMesh.C:1145
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:836
Automatically write from objectRegistry::writeObject()
virtual bool dynamic() const
Is mesh dynamic.
Definition: polyMesh.H:688
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:1386
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:654
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1512
polyBoundaryMesh BoundaryMesh
The boundary type associated with the mesh.
Definition: polyMesh.H:386
virtual ~polyMesh()
Destructor.
Definition: polyMesh.C:808
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:65
const objectRegistry & thisDb() const noexcept
Return the object registry.
Definition: polyMesh.H:677
bool hasTetBasePtIs() const
Definition: polyMesh.H:909
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:592
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:73
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:89
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:966
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:171
Foam::faceZoneMesh.
A primitive field of type <T> with automated input and output.
void clearCellTree()
Clear cell tree data.
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.