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