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-2020 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 polyMesh_H
45 #define 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
91  enum readUpdateState
92  {
97  };
98 
99  //- Enumeration defining the decomposition of the cell for
100  // inside/outside test
101  enum cellDecomposition
102  {
103  FACE_PLANES, //- Faces considered as planes
104 
105  FACE_CENTRE_TRIS, //- Faces decomposed into triangles
106  // using face-centre
107 
108  FACE_DIAG_TRIS, //- Faces decomposed into triangles diagonally
109 
110  CELL_TETS //- Cell decomposed into tets
111  };
112 
113 
114 private:
115 
116  // 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  typedef polyMesh Mesh;
313 
314 
315  //- Runtime type information
316  TypeName("polyMesh");
317 
318  //- Return the default region name
319  static word defaultRegion;
320 
321  //- Return the mesh sub-directory name (usually "polyMesh")
322  static word meshSubDir;
323 
324 
325  // Constructors
326 
327  //- Read construct from IOobject
328  polyMesh(const IOobject& io, const bool doInit = true);
329 
330  //- Construct from IOobject or as zero-sized mesh
331  // Boundary is added using addPatches() member function
332  polyMesh(const IOobject& io, const zero, const bool syncPar=true);
333 
334  //- Construct from IOobject and components.
335  // Boundary is added using addPatches() member function
336  polyMesh
337  (
338  const IOobject& io,
339  pointField&& points,
340  faceList&& faces,
341  labelList&& owner,
342  labelList&& neighbour,
343  const bool syncPar = true
344  );
345 
346  //- Construct without boundary with cells rather than owner/neighbour.
347  // Boundary is added using addPatches() member function
348  polyMesh
349  (
350  const IOobject& io,
351  pointField&& points,
352  faceList&& faces,
353  cellList&& cells,
354  const bool syncPar = true
355  );
356 
357  //- Construct from cell shapes
358  polyMesh
359  (
360  const IOobject& io,
361  pointField&& points,
362  const cellShapeList& shapes,
363  const faceListList& boundaryFaces,
364  const wordList& boundaryPatchNames,
365  const wordList& boundaryPatchTypes,
366  const word& defaultBoundaryPatchName,
367  const word& defaultBoundaryPatchType,
368  const wordList& boundaryPatchPhysicalTypes,
369  const bool syncPar = true
370  );
371 
372  //- Construct from cell shapes, with patch information in dictionary
373  //- format.
374  polyMesh
375  (
376  const IOobject& io,
377  pointField&& points,
378  const cellShapeList& shapes,
379  const faceListList& boundaryFaces,
380  const wordList& boundaryPatchNames,
382  const word& defaultBoundaryPatchName,
383  const word& defaultBoundaryPatchType,
384  const bool syncPar = true
385  );
386 
387 
388  //- Destructor
389  virtual ~polyMesh();
390 
391 
392  // Member Functions
393 
394  // Database
395 
396  //- Override the objectRegistry dbDir for a single-region case
397  virtual const fileName& dbDir() const;
398 
399  //- Return the local mesh directory (dbDir()/meshSubDir)
400  fileName meshDir() const;
401 
402  //- Return the current instance directory for points
403  // Used in the construction of geometric mesh data dependent
404  // on points
405  const fileName& pointsInstance() const;
406 
407  //- Return the current instance directory for faces
408  const fileName& facesInstance() const;
409 
410  //- Set the instance for mesh files
411  void setInstance
412  (
413  const fileName& instance,
415  );
416 
417 
418  // Access
419 
420  //- Return raw points
421  virtual const pointField& points() const;
422 
423  //- Return true if io is up-to-date with points
424  virtual bool upToDatePoints(const regIOobject& io) const;
425 
426  //- Set io to be up-to-date with points
427  virtual void setUpToDatePoints(regIOobject& io) const;
428 
429  //- Return raw faces
430  virtual const faceList& faces() const;
431 
432  //- Return face owner
433  virtual const labelList& faceOwner() const;
434 
435  //- Return face neighbour
436  virtual const labelList& faceNeighbour() const;
437 
438  //- Return old points (mesh motion)
439  virtual const pointField& oldPoints() const;
440 
441  //- Return old cellCentres (mesh motion)
442  virtual const pointField& oldCellCentres() const;
443 
444  //- Return boundary mesh
445  const polyBoundaryMesh& boundaryMesh() const
446  {
447  return boundary_;
448  }
449 
450  //- Return mesh bounding box
451  const boundBox& bounds() const
452  {
453  return bounds_;
454  }
455 
456  //- Return the vector of geometric directions in mesh.
457  // Defined according to the presence of empty and wedge patches.
458  // 1 indicates unconstrained direction and -1 a constrained
459  // direction.
460  const Vector<label>& geometricD() const;
461 
462  //- Return the number of valid geometric dimensions in the mesh
463  label nGeometricD() const;
464 
465  //- Return the vector of solved-for directions in mesh.
466  // Differs from geometricD in that it includes for wedge cases
467  // the circumferential direction in case of swirl.
468  // 1 indicates valid direction and -1 an invalid direction.
469  const Vector<label>& solutionD() const;
470 
471  //- Return the number of valid solved-for dimensions in the mesh
472  label nSolutionD() const;
473 
474  //- Return the tetBasePtIs
475  const labelIOList& tetBasePtIs() const;
476 
477  //- Return the cell search tree
478  const indexedOctree<treeDataCell>& cellTree() const;
479 
480  //- Return point zone mesh
481  const pointZoneMesh& pointZones() const
482  {
483  return pointZones_;
484  }
485 
486  //- Return face zone mesh
487  const faceZoneMesh& faceZones() const
488  {
489  return faceZones_;
490  }
491 
492  //- Return cell zone mesh
493  const cellZoneMesh& cellZones() const
494  {
495  return cellZones_;
496  }
497 
498  //- Return parallel info
499  const globalMeshData& globalData() const;
500 
501  //- Return communicator used for parallel communication
502  label comm() const;
503 
504  //- Return communicator used for parallel communication
505  label& comm();
506 
507  //- Return the object registry
508  const objectRegistry& thisDb() const
509  {
510  return *this;
511  }
512 
513 
514  // Mesh motion
515 
516  //- Is mesh dynamic
517  virtual bool dynamic() const
518  {
519  return false;
520  }
521 
522  //- Is mesh moving
523  bool moving() const
524  {
525  return moving_;
526  }
527 
528  //- Set the mesh to be moving
529  bool moving(const bool m)
530  {
531  bool m0 = moving_;
532  moving_ = m;
533  return m0;
534  }
535 
536  //- Is mesh topology changing
537  bool topoChanging() const
538  {
539  return topoChanging_;
540  }
541 
542  //- Set the mesh topology to be changing
543  bool topoChanging(const bool c)
544  {
545  bool c0 = topoChanging_;
546  topoChanging_ = c;
547  return c0;
548  }
549 
550  //- Is mesh changing (topology changing and/or moving)
551  bool changing() const
552  {
553  return moving()||topoChanging();
554  }
555 
556  //- Move points, returns volumes swept by faces in motion
557  virtual tmp<scalarField> movePoints(const pointField&);
558 
559  //- Reset motion
560  void resetMotion() const;
561 
562 
563  // Topological change
564 
565  //- Return non-const access to the pointZones
567  {
568  return pointZones_;
569  }
570 
571  //- Return non-const access to the faceZones
573  {
574  return faceZones_;
575  }
576 
577  //- Return non-const access to the cellZones
579  {
580  return cellZones_;
581  }
582 
583  //- Add boundary patches
584  void addPatches
585  (
586  PtrList<polyPatch>& plist,
587  const bool validBoundary = true
588  );
589 
590  //- Add boundary patches
591  void addPatches
592  (
593  const List<polyPatch*>& p,
594  const bool validBoundary = true
595  );
596 
597  //- Add mesh zones
598  void addZones
599  (
600  const List<pointZone*>& pz,
601  const List<faceZone*>& fz,
602  const List<cellZone*>& cz
603  );
604 
605  //- Initialise all non-demand-driven data
606  virtual bool init(const bool doInit);
607 
608  //- Update the mesh based on the mesh files saved in
609  // time directories
610  virtual readUpdateState readUpdate();
611 
612  //- Update the mesh corresponding to given map
613  virtual void updateMesh(const mapPolyMesh& mpm);
614 
615  //- Remove boundary patches
616  void removeBoundary();
617 
618  //- Reset mesh primitive data. Assumes all patch info correct
619  // (so does e.g. parallel communication). If not use
620  // validBoundary=false
621  //
622  // \note Only autoPtr parameters that test as valid() are used
623  // for resetting, otherwise the existing entries are left
624  // untouched.
625  void resetPrimitives
626  (
629  autoPtr<labelList>&& owner,
630  autoPtr<labelList>&& neighbour,
631  const labelUList& patchSizes,
632  const labelUList& patchStarts,
633  const bool validBoundary = true
634  );
635 
636 
637  // Storage management
638 
639  //- Clear geometry
640  void clearGeom();
641 
642  //- Update geometry points; keep topology.
643  //- Optionally with new face decomposition
644  void updateGeomPoints
645  (
646  pointIOField&& newPoints,
647  autoPtr<labelIOList>& newTetBasePtIsPtr
648  );
649 
650  //- Clear addressing
651  void clearAddressing(const bool isMeshUpdate = false);
652 
653  //- Clear all geometry and addressing unnecessary for CFD
654  void clearOut();
655 
656  //- Clear primitive data (points, faces and cells)
657  void clearPrimitives();
658 
659  //- Clear tet base points
660  void clearTetBasePtIs();
661 
662  //- Clear cell tree data
663  void clearCellTree();
664 
665  //- Remove all files from mesh instance
666  void removeFiles(const fileName& instanceDir) const;
667 
668  //- Remove all files from mesh instance()
669  void removeFiles() const;
670 
671 
672  // Geometric checks. Selectively override primitiveMesh functionality.
673 
674  //- Check non-orthogonality
675  virtual bool checkFaceOrthogonality
676  (
677  const bool report = false,
678  labelHashSet* setPtr = nullptr
679  ) const;
680 
681  //- Check face skewness
682  virtual bool checkFaceSkewness
683  (
684  const bool report = false,
685  labelHashSet* setPtr = nullptr
686  ) const;
687 
688  //- Check edge alignment for 1D/2D cases
689  virtual bool checkEdgeAlignment
690  (
691  const bool report,
692  const Vector<label>& directions,
693  labelHashSet* setPtr
694  ) const;
695 
696  virtual bool checkCellDeterminant
697  (
698  const bool report,
699  labelHashSet* setPtr
700  ) const;
701 
702  //- Check mesh motion for correctness given motion points
703  virtual bool checkMeshMotion
704  (
705  const pointField& newPoints,
706  const bool report = false,
707  const bool detailedReport = false
708  ) const;
709 
710  //- Check for face weights
711  virtual bool checkFaceWeight
712  (
713  const bool report,
714  const scalar minWeight = 0.05,
715  labelHashSet* setPtr = nullptr
716  ) const;
717 
718  //- Check for neighbouring cell volumes
719  virtual bool checkVolRatio
720  (
721  const bool report,
722  const scalar minRatio = 0.01,
723  labelHashSet* setPtr = nullptr
724  ) const;
725 
726 
727  // Position search functions
728 
729  //- Find the cell, tetFacei and tetPti for point p
730  void findCellFacePt
731  (
732  const point& p,
733  label& celli,
734  label& tetFacei,
735  label& tetPti
736  ) const;
737 
738  //- Find the tetFacei and tetPti for point p in celli.
739  // tetFacei and tetPtI are set to -1 if not found
740  void findTetFacePt
741  (
742  const label celli,
743  const point& p,
744  label& tetFacei,
745  label& tetPti
746  ) const;
747 
748  //- Test if point p is in the celli
749  bool pointInCell
750  (
751  const point& p,
752  label celli,
754  ) const;
755 
756  //- Find cell enclosing this location and return index
757  // If not found -1 is returned
758  label findCell
759  (
760  const point& p,
762  ) const;
763 };
764 
765 
766 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
767 
768 } // End namespace Foam
769 
770 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
771 
772 #endif
773 
774 // ************************************************************************* //
Foam::polyMesh::addPatches
void addPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:961
Foam::polyMesh::changing
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:550
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
Foam::polyMesh::topoChanging
bool topoChanging() const
Is mesh topology changing.
Definition: polyMesh.H:536
Foam::polyMesh::oldCellCentres
virtual const pointField & oldCellCentres() const
Return old cellCentres (mesh motion)
Definition: polyMesh.C:1141
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::polyMesh::FACE_DIAG_TRIS
Definition: polyMesh.H:107
Foam::polyMesh::resetMotion
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1287
Foam::polyMesh::cellDecomposition
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:100
cellZoneMesh.H
Foam::cellZoneMesh.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::polyMesh::POINTS_MOVED
Definition: polyMesh.H:93
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:318
Foam::primitiveMesh::clearAddressing
void clearAddressing()
Clear topological data.
Definition: primitiveMeshClear.C:157
Foam::polyMesh::topoChanging
bool topoChanging(const bool c)
Set the mesh topology to be changing.
Definition: polyMesh.H:542
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::polyMesh::moving
bool moving(const bool m)
Set the mesh to be moving.
Definition: polyMesh.H:528
Foam::polyMesh::nGeometricD
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:869
Foam::polyMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1160
Foam::polyMesh::BoundaryMesh
polyBoundaryMesh BoundaryMesh
Definition: polyMesh.H:311
Foam::polyMesh::moving
bool moving() const
Is mesh moving.
Definition: polyMesh.H:522
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::polyMesh::dbDir
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:829
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobjectI.H:191
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:321
Foam::polyMesh::clearCellTree
void clearCellTree()
Clear cell tree data.
Definition: polyMeshClear.C:237
cellShapeList.H
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
objectRegistry.H
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:852
Foam::polyMesh::Mesh
polyMesh Mesh
Definition: polyMesh.H:310
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:492
pointIOField.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
primitiveMesh.H
Foam::polyMesh::clearOut
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
Definition: polyMeshClear.C:222
Foam::polyMesh::init
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: polyMesh.C:367
Foam::polyMesh::solutionD
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:875
Foam::polyMesh::removeFiles
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1349
faceList.H
Foam::directions
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:67
Foam::HashSet< label, Hash< label > >
Foam::polyMesh::geometricD
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:858
Foam::polyMesh::pointInCell
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1397
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:486
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::polyMesh::TOPO_PATCH_CHANGE
Definition: polyMesh.H:95
Foam::polyMesh::cellZones
cellZoneMesh & cellZones()
Return non-const access to the cellZones.
Definition: polyMesh.H:577
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:846
Foam::IOobject::writeOption
writeOption
Enumeration defining the write options.
Definition: IOobject.H:127
boundaryDicts
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())
Definition: createBlockMesh.H:64
Foam::Field< vector >
Foam::polyMesh::FACE_CENTRE_TRIS
Definition: polyMesh.H:104
Foam::polyMesh::clearPrimitives
void clearPrimitives()
Clear primitive data (points, faces and cells)
Definition: polyMeshClear.C:209
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:480
Foam::polyMesh::faceZones
faceZoneMesh & faceZones()
Return non-const access to the faceZones.
Definition: polyMesh.H:571
Foam::polyMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
Definition: polyMeshUpdate.C:42
Foam::polyMesh::clearTetBasePtIs
void clearTetBasePtIs()
Clear tet base points.
Definition: polyMeshClear.C:229
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::polyMesh::resetPrimitives
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:718
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:50
Foam::ZoneMesh< pointZone, polyMesh >
Foam::CompactIOList< face, label >
Foam::polyMesh::UNCHANGED
Definition: polyMesh.H:92
Foam::polyMesh::thisDb
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:507
cellVols
const scalarField & cellVols
Definition: temperatureAndPressureVariables.H:51
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
Foam::polyMesh::TOPO_CHANGE
Definition: polyMesh.H:94
cellList.H
Foam::polyMesh::comm
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1313
Foam::polyMesh::oldPoints
virtual const pointField & oldPoints() const
Return old points (mesh motion)
Definition: polyMesh.C:1119
pointZoneMesh.H
Foam::pointZoneMesh.
Foam::globalMeshData
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Definition: globalMeshData.H:107
Foam::polyMesh::findCellFacePt
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1356
faceZoneMesh.H
Foam::faceZoneMesh.
Foam::polyMesh::upToDatePoints
virtual bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1082
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
faceIOList.H
Foam::polyMesh::meshDir
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:840
Foam::polyMesh::removeBoundary
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:39
Foam::polyMesh::clearGeom
void clearGeom()
Clear geometry.
Definition: polyMeshClear.C:53
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
Foam::polyMesh::nSolutionD
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:886
Foam::polyMesh::bounds
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:450
pointField.H
boundBox.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::polyMesh::cellTree
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:926
Foam::polyMesh::~polyMesh
virtual ~polyMesh()
Destructor.
Definition: polyMesh.C:820
Foam::regIOobject
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:71
Foam::polyMesh::findCell
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1507
Foam::polyMesh::CELL_TETS
Definition: polyMesh.H:109
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
Foam::polyMesh::pointZones
pointZoneMesh & pointZones()
Return non-const access to the pointZones.
Definition: polyMesh.H:565
Foam::polyMesh::TypeName
TypeName("polyMesh")
Runtime type information.
Foam::Vector< label >
Foam::List< cell >
labelIOList.H
Foam::UList< label >
Foam::IOList< label >
polyBoundaryMesh.H
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::polyMesh::FACE_PLANES
Definition: polyMesh.H:102
Foam::pointCells
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition: pointCells.H:56
Foam::polyMesh::updateGeomPoints
void updateGeomPoints(pointIOField &&newPoints, autoPtr< labelIOList > &newTetBasePtIsPtr)
Definition: polyMeshClear.C:75
Foam::polyMesh::checkMeshMotion
virtual bool checkMeshMotion(const pointField &newPoints, const bool report=false, const bool detailedReport=false) const
Check mesh motion for correctness given motion points.
Definition: polyMeshCheck.C:740
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1295
Foam::polyMesh::findTetFacePt
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:1381
Foam::polyMesh::addZones
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:999
Foam::polyMesh::dynamic
virtual bool dynamic() const
Is mesh dynamic.
Definition: polyMesh.H:516
Foam::polyMesh::setUpToDatePoints
virtual void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1088
Foam::polyMesh::setInstance
void setInstance(const fileName &instance, const IOobject::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:36
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::polyMesh::readUpdate
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:75
Foam::zero
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:62
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:89
Foam::polyMesh::tetBasePtIs
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:892
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:78