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-2022 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
105  FACE_CENTRE_TRIS, //- Faces decomposed into triangles
106  // using face-centre
108  FACE_DIAG_TRIS, //- Faces decomposed into triangles diagonally
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  // Geometry checks
243 
244  //- Check non-orthogonality
245  bool checkFaceOrthogonality
246  (
247  const vectorField& fAreas,
248  const vectorField& cellCtrs,
249  const bool report,
250  const bool detailedReport,
251  labelHashSet* setPtr
252  ) const;
253 
254  //- Check face skewness
255  bool checkFaceSkewness
256  (
257  const pointField& points,
258  const vectorField& fCtrs,
259  const vectorField& fAreas,
260  const vectorField& cellCtrs,
261  const bool report,
262  const bool detailedReport,
263  labelHashSet* setPtr
264  ) const;
265 
266  bool checkEdgeAlignment
267  (
268  const pointField& p,
269  const bool report,
270  const Vector<label>& directions,
271  labelHashSet* setPtr
272  ) const;
273 
274  bool checkCellDeterminant
275  (
276  const vectorField& faceAreas,
277  const bool report,
278  labelHashSet* setPtr,
279  const Vector<label>& meshD
280  ) const;
281 
282  bool checkFaceWeight
283  (
284  const vectorField& fCtrs,
285  const vectorField& fAreas,
286  const vectorField& cellCtrs,
287  const bool report,
288  const scalar minWeight,
289  labelHashSet* setPtr
290  ) const;
291 
292  bool checkVolRatio
293  (
294  const scalarField& cellVols,
295  const bool report,
296  const scalar minRatio,
297  labelHashSet* setPtr
298  ) const;
299 
300 
301 protected:
302 
303  // Protected Member Data
304 
305  //- Member data pending transfer to fvMesh
306 
307  //- Is the mesh moving
308  bool moving_;
309 
310  //- Is the mesh topology changing
311  bool topoChanged_;
312 
313 
314 public:
315 
316  // Public Typedefs
318  typedef polyMesh Mesh;
320 
321 
322  //- Runtime type information
323  TypeName("polyMesh");
324 
325  //- Return the default region name
326  static word defaultRegion;
327 
328  //- Return the mesh sub-directory name (usually "polyMesh")
329  static word meshSubDir;
330 
331 
332  // Constructors
333 
334  //- Construct from IOobject
335  explicit polyMesh(const IOobject& io);
336 
337  //- Move construct from IOobject or from 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  //- Move construct without boundary with cells rather than
350  // owner/neighbour.
351  // Boundary is added using addPatches() member function
352  polyMesh
353  (
354  const IOobject& io,
355  pointField&& points,
356  faceList&& faces,
357  cellList&& cells,
358  const bool syncPar = true
359  );
360 
361  //- Move construct from cell shapes
362  polyMesh
363  (
364  const IOobject& io,
365  pointField&& points,
366  const cellShapeList& shapes,
367  const faceListList& boundaryFaces,
368  const wordList& boundaryPatchNames,
369  const wordList& boundaryPatchTypes,
370  const word& defaultBoundaryPatchName,
371  const word& defaultBoundaryPatchType,
372  const wordList& boundaryPatchPhysicalTypes,
373  const bool syncPar = true
374  );
375 
376  //- Move construct from cell shapes with patch information in dictionary
377  // format.
378  polyMesh
379  (
380  const IOobject& io,
381  pointField&& points,
382  const cellShapeList& shapes,
383  const faceListList& boundaryFaces,
384  const wordList& boundaryPatchNames,
385  const PtrList<dictionary>& boundaryDicts,
386  const word& defaultBoundaryPatchName,
387  const word& defaultBoundaryPatchType,
388  const bool syncPar = true
389  );
390 
391  //- Disallow default bitwise copy construction
392  polyMesh(const polyMesh&);
393 
394 
395  //- Destructor
396  virtual ~polyMesh();
397 
398 
399  // Member Functions
400 
401  // Database
402 
403  //- Override the objectRegistry dbDir for a single-region case
404  virtual const fileName& dbDir() const;
405 
406  //- Return the local mesh directory (dbDir()/meshSubDir)
407  fileName meshDir() const;
408 
409  //- Return the current instance directory for points
410  // Used in the construction of gemometric mesh data dependent
411  // on points
412  const fileName& pointsInstance() const;
413 
414  //- Return the current instance directory for faces
415  const fileName& facesInstance() const;
416 
417  //- Set the instance for mesh files
418  void setInstance(const fileName&);
419 
420  //- Set the instance for the points files
421  void setPointsInstance(const fileName&);
422 
423 
424  // Access
425 
426  // Primitive mesh data
427 
428  //- Return raw points
429  virtual const pointField& points() const;
430 
431  //- Return raw faces
432  virtual const faceList& faces() const;
433 
434  //- Return face owner
435  virtual const labelList& faceOwner() const;
436 
437  //- Return face neighbour
438  virtual const labelList& faceNeighbour() const;
439 
440  //- Return old points for mesh motion
441  virtual const pointField& oldPoints() const;
442 
443  //- Return old cell centres for mesh motion
444  virtual const pointField& oldCellCentres() const;
445 
446 
447  //- Return true if io is up-to-date with points
448  bool upToDatePoints(const regIOobject& io) const;
449 
450  //- Set io to be up-to-date with points
451  void setUpToDatePoints(regIOobject& io) const;
452 
453  //- Return boundary mesh
454  const polyBoundaryMesh& boundaryMesh() const
455  {
456  return boundary_;
457  }
458 
459  //- Return mesh bounding box
460  const boundBox& bounds() const
461  {
462  return bounds_;
463  }
464 
465  //- Return the vector of geometric directions in mesh.
466  // Defined according to the presence of empty and wedge patches.
467  // 1 indicates unconstrained direction and -1 a constrained
468  // direction.
469  const Vector<label>& geometricD() const;
470 
471  //- Return the number of valid geometric dimensions in the mesh
472  label nGeometricD() const;
473 
474  //- Return the vector of solved-for directions in mesh.
475  // Differs from geometricD in that it includes for wedge cases
476  // the circumferential direction in case of swirl.
477  // 1 indicates valid direction and -1 an invalid direction.
478  const Vector<label>& solutionD() const;
479 
480  //- Return the number of valid solved-for dimensions in the mesh
481  label nSolutionD() const;
482 
483  //- Return the tetBasePtIs
484  const labelIOList& tetBasePtIs() const;
485 
486  //- Return the cell search tree
487  const indexedOctree<treeDataCell>& cellTree() const;
488 
489  //- Return point zones
490  const meshPointZones& pointZones() const
491  {
492  return pointZones_;
493  }
494 
495  //- Return face zones
496  const meshFaceZones& faceZones() const
497  {
498  return faceZones_;
499  }
500 
501  //- Return cell zones
502  const meshCellZones& cellZones() const
503  {
504  return cellZones_;
505  }
506 
507  //- Return parallel info
508  const globalMeshData& globalData() const;
509 
510  //- Return communicator used for parallel communication
511  label comm() const;
512 
513  //- Return communicator used for parallel communication
514  label& comm();
515 
516  //- Return the object registry
517  const objectRegistry& thisDb() const
518  {
519  return *this;
520  }
521 
522 
523  // Mesh motion
524 
525  //- Is mesh moving
526  bool moving() const
527  {
528  return moving_;
529  }
530 
531  //- Is mesh topology changing
532  bool topoChanged() const
533  {
534  return topoChanged_;
535  }
536 
537  //- Is mesh changing (topology changing and/or moving)
538  bool changing() const
539  {
540  return moving() || topoChanged();
541  }
542 
543  //- Reset the points
544  // without storing old points or returning swept volumes
545  virtual void setPoints(const pointField&);
546 
547  //- Move points, returns volumes swept by faces in motion
548  virtual tmp<scalarField> movePoints(const pointField&);
549 
550  //- Reset motion
551  void resetMotion() const;
552 
553 
554  // Topological change
555 
556  //- Return non-const access to the pointZones
558  {
559  return pointZones_;
560  }
561 
562  //- Return non-const access to the faceZones
564  {
565  return faceZones_;
566  }
567 
568  //- Return non-const access to the cellZones
570  {
571  return cellZones_;
572  }
573 
574  //- Add boundary patches
575  void addPatches
576  (
577  const List<polyPatch*>&,
578  const bool validBoundary = true
579  );
580 
581  //- Add mesh zones
582  void addZones
583  (
584  const List<pointZone*>& pz,
585  const List<faceZone*>& fz,
586  const List<cellZone*>& cz
587  );
588 
589  //- Add/insert single patch. If validBoundary the new situation
590  // is consistent across processors.
591  virtual void addPatch
592  (
593  const label insertPatchi,
594  const polyPatch& patch,
595  const dictionary& patchFieldDict,
596  const word& defaultPatchFieldType,
597  const bool validBoundary
598  );
599 
600  //- Reorder and trim existing patches. If validBoundary the new
601  // situation is consistent across processors
602  virtual void reorderPatches
603  (
604  const labelUList& newToOld,
605  const bool validBoundary
606  );
607 
608  //- Update the mesh based on the mesh files saved in
609  // time directories
611 
612  //- Update topology using the given map
613  virtual void topoChange(const polyTopoChangeMap&);
614 
615  //- Update from another mesh using the given map
616  virtual void mapMesh(const polyMeshMap&);
617 
618  //- Redistribute or update using the given distribution map
619  virtual void distribute(const polyDistributionMap& map);
620 
621  //- Remove boundary patches
622  void removeBoundary();
623 
624  //- Reset mesh primitive data. Assumes all patch info correct
625  // (so does e.g. parallel communication). If not use
626  // validBoundary=false
627  void resetPrimitives
628  (
629  pointField&& points,
630  faceList&& faces,
631  labelList&& owner,
632  labelList&& neighbour,
633  const labelList& patchSizes,
634  const labelList& patchStarts,
635  const bool validBoundary = true
636  );
637 
638  //- Reset mesh
639  // For run-time mesh replacement and mesh to mesh mapping
640  void reset(const polyMesh&);
641 
642 
643  // Storage management
644 
645  //- Clear geometry
646  void clearGeom();
647 
648  //- Clear addressing
649  void clearAddressing(const bool isMeshUpdate = false);
650 
651  //- Clear all geometry and addressing unnecessary for CFD
652  void clearOut();
653 
654  //- Clear primitive data (points, faces and cells)
655  void clearPrimitives();
656 
657  //- Clear tet base points
658  void clearTetBasePtIs();
659 
660  //- Clear cell tree data
661  void clearCellTree();
662 
663  //- Remove all files from mesh instance
664  void removeFiles(const fileName& instanceDir) const;
665 
666  //- Remove all files from mesh instance()
667  void removeFiles() const;
668 
669 
670  // Geometric checks. Selectively override primitiveMesh functionality.
671 
672  //- Check non-orthogonality
673  virtual bool checkFaceOrthogonality
674  (
675  const bool report = false,
676  labelHashSet* setPtr = nullptr
677  ) const;
678 
679  //- Check face skewness
680  virtual bool checkFaceSkewness
681  (
682  const bool report = false,
683  labelHashSet* setPtr = nullptr
684  ) const;
685 
686  //- Check edge alignment for 1D/2D cases
687  virtual bool checkEdgeAlignment
688  (
689  const bool report,
690  const Vector<label>& directions,
691  labelHashSet* setPtr
692  ) const;
693 
694  virtual bool checkCellDeterminant
695  (
696  const bool report,
697  labelHashSet* setPtr
698  ) const;
699 
700  //- Check mesh motion for correctness given motion points
701  virtual bool checkMeshMotion
702  (
703  const pointField& newPoints,
704  const bool report = false,
705  const bool detailedReport = false
706  ) const;
707 
708  //- Check for face weights
709  virtual bool checkFaceWeight
710  (
711  const bool report,
712  const scalar minWeight = 0.05,
713  labelHashSet* setPtr = nullptr
714  ) const;
715 
716  //- Check for neighbouring cell volumes
717  virtual bool checkVolRatio
718  (
719  const bool report,
720  const scalar minRatio = 0.01,
721  labelHashSet* setPtr = nullptr
722  ) const;
723 
724 
725  // Position search functions
726 
727  //- Find the cell, tetFacei and tetPti for point p
728  void findCellFacePt
729  (
730  const point& p,
731  label& celli,
732  label& tetFacei,
733  label& tetPti
734  ) const;
735 
736  //- Find the tetFacei and tetPti for point p in celli.
737  // tetFacei and tetPtI are set to -1 if not found
738  void findTetFacePt
739  (
740  const label celli,
741  const point& p,
742  label& tetFacei,
743  label& tetPti
744  ) const;
745 
746  //- Test if point p is in the celli
747  bool pointInCell
748  (
749  const point& p,
750  label celli,
752  ) const;
753 
754  //- Find cell enclosing this location and return index
755  // If not found -1 is returned
757  (
758  const point& p,
760  ) const;
761 
762 
763  // Write
764 
765  //- Write the underlying polyMesh
766  virtual bool writeObject
767  (
771  const bool write = true
772  ) const;
773 
774 
775  // Member Operators
776 
777  //- Disallow default bitwise assignment
778  void operator=(const polyMesh&) = delete;
779 };
780 
781 
782 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
783 
784 } // End namespace Foam
785 
786 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
787 
788 #endif
789 
790 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
bool moving_
Member data pending transfer to fvMesh.
Definition: polyMesh.H:307
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: polyMesh.C:1091
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:537
polyMesh Mesh
Definition: polyMesh.H:317
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1342
writeOption
Enumeration defining the write options.
Definition: IOobject.H:125
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
void clearAddressing()
Clear topological data.
A class for handling file names.
Definition: fileName.H:79
virtual void distribute(const polyDistributionMap &map)
Redistribute or update using the given distribution map.
bool moving() const
Is mesh moving.
Definition: polyMesh.H:525
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:888
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:501
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1243
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:928
void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1293
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1427
label nFaces() const
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:911
Foam::meshFaceZones.
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:64
void reset(const polyMesh &)
Reset mesh.
Definition: polyMesh.C:772
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:325
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
void setPointsInstance(const fileName &)
Set the instance for the points files.
Definition: polyMeshIO.C:57
const cellList & cells() const
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:863
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write the underlying polyMesh.
Definition: polyMeshIO.C:505
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1537
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:100
const dimensionedScalar c
Speed of light in a vacuum.
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:905
const meshPointZones & pointZones() const
Return point zones.
Definition: polyMesh.H:489
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:894
bool topoChanged_
Is the mesh topology changing.
Definition: polyMesh.H:310
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:882
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1489
void clearPrimitives()
Clear primitive data (points, faces and cells)
Foam::meshPointZones.
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1496
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:1129
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1249
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1237
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
bool topoChanged() const
Is mesh topology changing.
Definition: polyMesh.H:531
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1435
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1453
void addZones(const List< pointZone *> &pz, const List< faceZone *> &fz, const List< cellZone *> &cz)
Add mesh zones.
Definition: polyMesh.C:1033
void addPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:989
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
const labelListList & pointCells() const
Foam::polyBoundaryMesh.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross product operators.
Definition: Vector.H:57
void clearGeom()
Clear geometry.
Definition: polyMeshClear.C:53
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
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:672
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
const scalarField & cellVols
void clearTetBasePtIs()
Clear tet base points.
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:77
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:922
bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1287
void operator=(const polyMesh &)=delete
Disallow default bitwise assignment.
TypeName("polyMesh")
Runtime type information.
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:962
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:459
virtual const pointField & oldCellCentres() const
Return old cell centres for mesh motion.
Definition: polyMesh.C:1267
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
Foam::meshCellZones.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:495
const vectorField & faceAreas() const
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:876
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:1521
Version number type.
Definition: IOstream.H:96
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1647
polyBoundaryMesh BoundaryMesh
Definition: polyMesh.H:318
virtual ~polyMesh()
Destructor.
Definition: polyMesh.C:854
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:52
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:516
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
virtual bool write(const bool write=true) const
Write using setting from DB.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:98
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual void setPoints(const pointField &)
Reset the points.
Definition: polyMesh.C:1299
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
void clearCellTree()
Clear cell tree data.
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.
polyMesh(const IOobject &io)
Construct from IOobject.
Definition: polyMesh.C:163