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 | Copyright (C) 2011-2015 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 "pointZoneMesh.H"
56 #include "faceZoneMesh.H"
57 #include "cellZoneMesh.H"
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 namespace Foam
62 {
63 
64 // Forward declaration of classes
65 class globalMeshData;
66 class mapPolyMesh;
67 class polyMeshTetDecomposition;
68 class treeDataCell;
69 template<class Type> class indexedOctree;
70 
71 /*---------------------------------------------------------------------------*\
72  Class polyMesh Declaration
73 \*---------------------------------------------------------------------------*/
74 
75 class polyMesh
76 :
77  public objectRegistry,
78  public primitiveMesh
79 {
80 
81 public:
82 
83  // Public data types
84 
85  //- Enumeration defining the state of the mesh after a read update.
86  // Used for post-processing applications, where the mesh
87  // needs to update based on the files written in time
88  // directores
89  enum readUpdateState
90  {
95  };
96 
97  //- Enumeration defining the decomposition of the cell for
98  // inside/outside test
100  {
101  FACE_PLANES, //- Faces considered as planes
102 
103  FACE_CENTRE_TRIS, //- Faces decomposed into triangles
104  // using face-centre
105 
106  FACE_DIAG_TRIS, //- Faces decomposed into triangles diagonally
107 
108  CELL_TETS //- Cell decomposed into tets
109  };
110 
111 
112 private:
113 
114  // Permanent data
115 
116  // Primitive mesh data
117 
118  //- Points
120 
121  //- Faces
123 
124  //- Face owner
126 
127  //- Face neighbour
129 
130  //- Have the primitives been cleared
131  bool clearedPrimitives_;
132 
133 
134  //- Boundary mesh
135  mutable polyBoundaryMesh boundary_;
136 
137  //- Mesh bounding-box.
138  // Created from points on construction, updated when the mesh moves
140 
141  //- Communicator used for parallel communication
142  label comm_;
143 
144  //- Vector of non-constrained directions in mesh
145  // defined according to the presence of empty and wedge patches
146  mutable Vector<label> geometricD_;
147 
148  //- Vector of valid directions in mesh
149  // defined according to the presence of empty patches
150  mutable Vector<label> solutionD_;
151 
152  //- Base point for face decomposition into tets
154 
155  //- Search tree to allow spatial cell searching
157 
158 
159  // Zoning information
160 
161  //- Point zones
163 
164  //- Face zones
166 
167  //- Cell zones
169 
170 
171  //- Parallel info
173 
174 
175  // Mesh motion related data
176 
177  //- Is the mesh moving
178  bool moving_;
179 
180  //- Is the mesh topology changing
181  bool topoChanging_;
182 
183  //- Current time index for mesh motion
184  mutable label curMotionTimeIndex_;
185 
186  //- Old points (for the last mesh motion)
188 
189 
190  // Private Member Functions
191 
192  //- Disallow construct as copy
193  polyMesh(const polyMesh&);
194 
195  //- Disallow default bitwise assignment
196  void operator=(const polyMesh&);
197 
198  //- Initialise the polyMesh from the primitive data
199  void initMesh();
200 
201  //- Initialise the polyMesh from the given set of cells
202  void initMesh(cellList& c);
203 
204  //- Calculate the valid directions in the mesh from the boundaries
205  void calcDirections() const;
206 
207  //- Calculate the cell shapes from the primitive
208  // polyhedral information
209  void calcCellShapes() const;
210 
211 
212  // Helper functions for constructor from cell shapes
213 
215 
217  (
218  const faceList& patchFaces,
219  const labelListList& pointCells,
220  const faceListList& cellsFaceShapes,
221  const label patchID
222  ) const;
223 
224  void setTopology
225  (
226  const cellShapeList& cellsAsShapes,
227  const faceListList& boundaryFaces,
228  const wordList& boundaryPatchNames,
229  labelList& patchSizes,
230  labelList& patchStarts,
231  label& defaultPatchStart,
232  label& nFaces,
233  cellList& cells
234  );
235 
236 
237  // Geometry checks
238 
239  //- Check non-orthogonality
241  (
242  const vectorField& fAreas,
243  const vectorField& cellCtrs,
244  const bool report,
245  const bool detailedReport,
246  labelHashSet* setPtr
247  ) const;
248 
249  //- Check face skewness
250  bool checkFaceSkewness
251  (
252  const pointField& points,
253  const vectorField& fCtrs,
254  const vectorField& fAreas,
255  const vectorField& cellCtrs,
256  const bool report,
257  const bool detailedReport,
258  labelHashSet* setPtr
259  ) const;
260 
261  bool checkEdgeAlignment
262  (
263  const pointField& p,
264  const bool report,
265  const Vector<label>& directions,
266  labelHashSet* setPtr
267  ) const;
268 
270  (
271  const vectorField& faceAreas,
272  const bool report,
273  labelHashSet* setPtr,
274  const Vector<label>& meshD
275  ) const;
276 
277  bool checkFaceWeight
278  (
279  const vectorField& fCtrs,
280  const vectorField& fAreas,
281  const vectorField& cellCtrs,
282  const bool report,
283  const scalar minWeight,
284  labelHashSet* setPtr
285  ) const;
286 
287  bool checkVolRatio
288  (
289  const scalarField& cellVols,
290  const bool report,
291  const scalar minRatio,
292  labelHashSet* setPtr
293  ) const;
294 
295 public:
296 
297  // Public typedefs
298 
299  typedef polyMesh Mesh;
301 
302 
303  //- Runtime type information
304  TypeName("polyMesh");
305 
306  //- Return the default region name
307  static word defaultRegion;
308 
309  //- Return the mesh sub-directory name (usually "polyMesh")
310  static word meshSubDir;
311 
312 
313  // Constructors
314 
315  //- Construct from IOobject
316  explicit polyMesh(const IOobject& io);
317 
318  //- Construct from IOobject or from components.
319  // Boundary is added using addPatches() member function
320  polyMesh
321  (
322  const IOobject& io,
323  const Xfer<pointField>& points,
324  const Xfer<faceList>& faces,
325  const Xfer<labelList>& owner,
326  const Xfer<labelList>& neighbour,
327  const bool syncPar = true
328  );
329 
330  //- Construct without boundary with cells rather than owner/neighbour.
331  // Boundary is added using addPatches() member function
332  polyMesh
333  (
334  const IOobject& io,
335  const Xfer<pointField>& points,
336  const Xfer<faceList>& faces,
337  const Xfer<cellList>& cells,
338  const bool syncPar = true
339  );
340 
341  //- Construct from cell shapes
342  polyMesh
343  (
344  const IOobject& io,
345  const Xfer<pointField>& points,
346  const cellShapeList& shapes,
347  const faceListList& boundaryFaces,
348  const wordList& boundaryPatchNames,
349  const wordList& boundaryPatchTypes,
350  const word& defaultBoundaryPatchName,
351  const word& defaultBoundaryPatchType,
352  const wordList& boundaryPatchPhysicalTypes,
353  const bool syncPar = true
354  );
355 
356  //- Construct from cell shapes with patch information in dictionary
357  // format.
358  polyMesh
359  (
360  const IOobject& io,
361  const Xfer<pointField>& points,
362  const cellShapeList& shapes,
363  const faceListList& boundaryFaces,
364  const wordList& boundaryPatchNames,
365  const PtrList<dictionary>& boundaryDicts,
366  const word& defaultBoundaryPatchName,
367  const word& defaultBoundaryPatchType,
368  const bool syncPar = true
369  );
370 
371 
372  //- Destructor
373  virtual ~polyMesh();
374 
375 
376  // Member Functions
377 
378  // Database
379 
380  //- Override the objectRegistry dbDir for a single-region case
381  virtual const fileName& dbDir() const;
382 
383  //- Return the local mesh directory (dbDir()/meshSubDir)
384  fileName meshDir() const;
385 
386  //- Return the current instance directory for points
387  // Used in the consruction of gemometric mesh data dependent
388  // on points
389  const fileName& pointsInstance() const;
390 
391  //- Return the current instance directory for faces
392  const fileName& facesInstance() const;
393 
394  //- Set the instance for mesh files
395  void setInstance(const fileName&);
396 
397 
398  // Access
399 
400  //- Return raw points
401  virtual const pointField& points() const;
402 
403  //- Return true if io is up-to-date with points
404  virtual bool upToDatePoints(const regIOobject& io) const;
405 
406  //- Set io to be up-to-date with points
407  virtual void setUpToDatePoints(regIOobject& io) const;
408 
409  //- Return raw faces
410  virtual const faceList& faces() const;
411 
412  //- Return face owner
413  virtual const labelList& faceOwner() const;
414 
415  //- Return face neighbour
416  virtual const labelList& faceNeighbour() const;
417 
418  //- Return old points for mesh motion
419  virtual const pointField& oldPoints() const;
420 
421  //- Return boundary mesh
422  const polyBoundaryMesh& boundaryMesh() const
423  {
424  return boundary_;
425  }
426 
427  //- Return mesh bounding box
428  const boundBox& bounds() const
429  {
430  return bounds_;
431  }
432 
433  //- Return the vector of geometric directions in mesh.
434  // Defined according to the presence of empty and wedge patches.
435  // 1 indicates unconstrained direction and -1 a constrained
436  // direction.
437  const Vector<label>& geometricD() const;
438 
439  //- Return the number of valid geometric dimensions in the mesh
440  label nGeometricD() const;
441 
442  //- Return the vector of solved-for directions in mesh.
443  // Differs from geometricD in that it includes for wedge cases
444  // the circumferential direction in case of swirl.
445  // 1 indicates valid direction and -1 an invalid direction.
446  const Vector<label>& solutionD() const;
447 
448  //- Return the number of valid solved-for dimensions in the mesh
449  label nSolutionD() const;
450 
451  //- Return the tetBasePtIs
452  const labelList& tetBasePtIs() const;
453 
454  //- Return the cell search tree
455  const indexedOctree<treeDataCell>& cellTree() const;
456 
457  //- Return point zone mesh
458  const pointZoneMesh& pointZones() const
459  {
460  return pointZones_;
461  }
462 
463  //- Return face zone mesh
464  const faceZoneMesh& faceZones() const
465  {
466  return faceZones_;
467  }
468 
469  //- Return cell zone mesh
470  const cellZoneMesh& cellZones() const
471  {
472  return cellZones_;
473  }
474 
475  //- Return parallel info
476  const globalMeshData& globalData() const;
477 
478  //- Return communicator used for parallel communication
479  label comm() const;
480 
481  //- Return communicator used for parallel communication
482  label& comm();
483 
484  //- Return the object registry
485  const objectRegistry& thisDb() const
486  {
487  return *this;
488  }
489 
490 
491  // Mesh motion
492 
493  //- Is mesh moving
494  bool moving() const
495  {
496  return moving_;
497  }
498 
499  //- Set the mesh to be moving
500  bool moving(const bool m)
501  {
502  bool m0 = moving_;
503  moving_ = m;
504  return m0;
505  }
506 
507  //- Is mesh topology changing
508  bool topoChanging() const
509  {
510  return topoChanging_;
511  }
512 
513  //- Set the mesh topology to be changing
514  bool topoChanging(const bool c)
515  {
516  bool c0 = topoChanging_;
517  topoChanging_ = c;
518  return c0;
519  }
520 
521  //- Is mesh changing (topology changing and/or moving)
522  bool changing() const
523  {
524  return moving()||topoChanging();
525  }
526 
527  //- Move points, returns volumes swept by faces in motion
528  virtual tmp<scalarField> movePoints(const pointField&);
529 
530  //- Reset motion
531  void resetMotion() const;
532 
533 
534  // Topological change
535 
536  //- Return non-const access to the pointZones
538  {
539  return pointZones_;
540  }
541 
542  //- Return non-const access to the faceZones
544  {
545  return faceZones_;
546  }
547 
548  //- Return non-const access to the cellZones
550  {
551  return cellZones_;
552  }
553 
554  //- Add boundary patches
555  void addPatches
556  (
557  const List<polyPatch*>&,
558  const bool validBoundary = true
559  );
560 
561  //- Add mesh zones
562  void addZones
563  (
564  const List<pointZone*>& pz,
565  const List<faceZone*>& fz,
566  const List<cellZone*>& cz
567  );
568 
569  //- Update the mesh based on the mesh files saved in
570  // time directories
571  virtual readUpdateState readUpdate();
572 
573  //- Update the mesh corresponding to given map
574  virtual void updateMesh(const mapPolyMesh& mpm);
575 
576  //- Remove boundary patches
577  void removeBoundary();
578 
579  //- Reset mesh primitive data. Assumes all patch info correct
580  // (so does e.g. parallel communication). If not use
581  // validBoundary=false
582  void resetPrimitives
583  (
584  const Xfer<pointField>& points,
585  const Xfer<faceList>& faces,
586  const Xfer<labelList>& owner,
587  const Xfer<labelList>& neighbour,
588  const labelList& patchSizes,
589  const labelList& patchStarts,
590  const bool validBoundary = true
591  );
592 
593 
594  // Storage management
595 
596  //- Clear geometry
597  void clearGeom();
598 
599  //- Clear addressing
600  void clearAddressing(const bool isMeshUpdate = false);
601 
602  //- Clear all geometry and addressing unnecessary for CFD
603  void clearOut();
604 
605  //- Clear primitive data (points, faces and cells)
606  void clearPrimitives();
607 
608  //- Clear geometry not used for CFD (cellTree, tetBasePtIs)
609  void clearAdditionalGeom();
610 
611  //- Clear cell tree data
612  void clearCellTree();
613 
614  //- Remove all files from mesh instance
615  void removeFiles(const fileName& instanceDir) const;
616 
617  //- Remove all files from mesh instance()
618  void removeFiles() const;
619 
620 
621  // Geometric checks. Selectively override primitiveMesh functionality.
622 
623  //- Check non-orthogonality
624  virtual bool checkFaceOrthogonality
625  (
626  const bool report = false,
627  labelHashSet* setPtr = NULL
628  ) const;
629 
630  //- Check face skewness
631  virtual bool checkFaceSkewness
632  (
633  const bool report = false,
634  labelHashSet* setPtr = NULL
635  ) const;
636 
637  //- Check edge alignment for 1D/2D cases
638  virtual bool checkEdgeAlignment
639  (
640  const bool report,
641  const Vector<label>& directions,
642  labelHashSet* setPtr
643  ) const;
644 
645  virtual bool checkCellDeterminant
646  (
647  const bool report,
648  labelHashSet* setPtr
649  ) const;
650 
651  //- Check mesh motion for correctness given motion points
652  virtual bool checkMeshMotion
653  (
654  const pointField& newPoints,
655  const bool report = false,
656  const bool detailedReport = false
657  ) const;
658 
659  //- Check for face weights
660  virtual bool checkFaceWeight
661  (
662  const bool report,
663  const scalar minWeight = 0.05,
664  labelHashSet* setPtr = NULL
665  ) const;
666 
667  //- Check for neighbouring cell volumes
668  virtual bool checkVolRatio
669  (
670  const bool report,
671  const scalar minRatio = 0.01,
672  labelHashSet* setPtr = NULL
673  ) const;
674 
675 
676  // Position search functions
677 
678  //- Find the cell, tetFacei and tetPti for point p
679  void findCellFacePt
680  (
681  const point& p,
682  label& celli,
683  label& tetFacei,
684  label& tetPti
685  ) const;
686 
687  //- Find the tetFacei and tetPti for point p in celli.
688  // tetFaceI and tetPtI are set to -1 if not found
689  void findTetFacePt
690  (
691  const label celli,
692  const point& p,
693  label& tetFacei,
694  label& tetPti
695  ) const;
696 
697  //- Test if point p is in the celli
698  bool pointInCell
699  (
700  const point& p,
701  label celli,
703  ) const;
704 
705  //- Find cell enclosing this location and return index
706  // If not found -1 is returned
708  (
709  const point& p,
711  ) const;
712 };
713 
714 
715 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
716 
717 } // End namespace Foam
718 
719 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
720 
721 #endif
722 
723 // ************************************************************************* //
Foam::polyMesh::checkCellDeterminant
bool checkCellDeterminant(const vectorField &faceAreas, const bool report, labelHashSet *setPtr, const Vector< label > &meshD) const
Definition: polyMeshCheck.C:414
Foam::polyMesh::changing
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:521
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::polyMesh::topoChanging
bool topoChanging() const
Is mesh topology changing.
Definition: polyMesh.H:507
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::polyMesh::calcCellShapes
void calcCellShapes() const
Calculate the cell shapes from the primitive.
Foam::polyMesh::FACE_DIAG_TRIS
@ FACE_DIAG_TRIS
Definition: polyMesh.H:105
Foam::polyMesh::checkVolRatio
bool checkVolRatio(const scalarField &cellVols, const bool report, const scalar minRatio, labelHashSet *setPtr) const
Definition: polyMeshCheck.C:595
Foam::polyMesh::resetMotion
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1135
Foam::polyMesh::cellDecomposition
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:98
cellZoneMesh.H
Foam::cellZoneMesh.
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::IOField< vector >
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::polyMesh::POINTS_MOVED
@ POINTS_MOVED
Definition: polyMesh.H:91
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Foam::primitiveMesh::clearAddressing
void clearAddressing()
Clear topological data.
Definition: primitiveMeshClear.C:142
Foam::polyMesh::topoChanging
bool topoChanging(const bool c)
Set the mesh topology to be changing.
Definition: polyMesh.H:513
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::polyMesh::faces_
faceCompactIOList faces_
Faces.
Definition: polyMesh.H:121
Foam::polyMesh::moving
bool moving(const bool m)
Set the mesh to be moving.
Definition: polyMesh.H:499
Foam::polyMesh::nGeometricD
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:795
Foam::polyMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1048
Foam::polyMesh::BoundaryMesh
polyBoundaryMesh BoundaryMesh
Definition: polyMesh.H:299
Foam::polyMesh::tetBasePtIs
const labelList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:818
Foam::polyMesh::facePatchFaceCells
labelList facePatchFaceCells(const faceList &patchFaces, const labelListList &pointCells, const faceListList &cellsFaceShapes, const label patchID) const
Definition: polyMeshFromShapeMesh.C:74
Foam::polyMesh::moving
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
Foam::polyMesh::dbDir
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:753
Foam::polyMesh::cellShapePointCells
labelListList cellShapePointCells(const cellShapeList &) const
Definition: polyMeshFromShapeMesh.C:37
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
Foam::polyMesh::clearCellTree
void clearCellTree()
Clear cell tree data.
Definition: polyMeshClear.C:180
Foam::polyMesh::faceZones_
faceZoneMesh faceZones_
Face zones.
Definition: polyMesh.H:164
Foam::polyMesh::cellTreePtr_
autoPtr< indexedOctree< treeDataCell > > cellTreePtr_
Search tree to allow spatial cell searching.
Definition: polyMesh.H:155
cellShapeList.H
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:136
objectRegistry.H
Foam::polyMesh::neighbour_
labelIOList neighbour_
Face neighbour.
Definition: polyMesh.H:127
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
Foam::polyMesh::Mesh
polyMesh Mesh
Definition: polyMesh.H:298
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
pointIOField.H
Foam::polyMesh::comm
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1161
Foam::polyMesh::bounds_
boundBox bounds_
Mesh bounding-box.
Definition: polyMesh.H:138
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::polyMesh::clearOut
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
Definition: polyMeshClear.C:173
Foam::polyMesh::operator=
void operator=(const polyMesh &)
Disallow default bitwise assignment.
Foam::polyMesh::solutionD
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:801
Foam::polyMesh::removeFiles
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1197
Foam::polyMesh::geometricD_
Vector< label > geometricD_
Vector of non-constrained directions in mesh.
Definition: polyMesh.H:145
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:63
Foam::HashSet< label, Hash< label > >
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
Foam::polyMesh::geometricD
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:784
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:1292
Foam::polyMesh::globalMeshDataPtr_
autoPtr< globalMeshData > globalMeshDataPtr_
Parallel info.
Definition: polyMesh.H:171
Foam::Xfer
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::polyMesh::setTopology
void setTopology(const cellShapeList &cellsAsShapes, const faceListList &boundaryFaces, const wordList &boundaryPatchNames, labelList &patchSizes, labelList &patchStarts, label &defaultPatchStart, label &nFaces, cellList &cells)
Set faces_, calculate cells and patchStarts.
Definition: polyMeshFromShapeMesh.C:132
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
Foam::polyMesh::topoChanging_
bool topoChanging_
Is the mesh topology changing.
Definition: polyMesh.H:180
Foam::polyMesh::clearAdditionalGeom
void clearAdditionalGeom()
Clear geometry not used for CFD (cellTree, tetBasePtIs)
Definition: polyMeshClear.C:83
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::polyMesh::TOPO_PATCH_CHANGE
@ TOPO_PATCH_CHANGE
Definition: polyMesh.H:93
Foam::polyMesh::cellZones
cellZoneMesh & cellZones()
Return non-const access to the cellZones.
Definition: polyMesh.H:548
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
Foam::polyMesh::owner_
labelIOList owner_
Face owner.
Definition: polyMesh.H:124
Foam::polyMesh::comm_
label comm_
Communicator used for parallel communication.
Definition: polyMesh.H:141
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::polyMesh::FACE_CENTRE_TRIS
@ FACE_CENTRE_TRIS
Definition: polyMesh.H:102
Foam::polyMesh::clearPrimitives
void clearPrimitives()
Clear primitive data (points, faces and cells)
Definition: polyMeshClear.C:160
Foam::polyMesh::pointZones
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:457
Foam::polyMesh::faceZones
faceZoneMesh & faceZones()
Return non-const access to the faceZones.
Definition: polyMesh.H:542
Foam::polyMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
Definition: polyMeshUpdate.C:39
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
Foam::polyMesh::boundary_
polyBoundaryMesh boundary_
Boundary mesh.
Definition: polyMesh.H:134
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataTriSurface.H:47
Foam::ZoneMesh
A list of mesh zones.
Definition: cellZoneMeshFwd.H:39
Foam::CompactIOList
A List of objects of type <T> with automated input and output using a compact storage....
Definition: CompactIOList.H:53
Foam::polyMesh::UNCHANGED
@ UNCHANGED
Definition: polyMesh.H:90
Foam::polyMesh::curMotionTimeIndex_
label curMotionTimeIndex_
Current time index for mesh motion.
Definition: polyMesh.H:183
Foam::polyMesh::thisDb
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:484
cellVols
const scalarField & cellVols
Definition: temperatureAndPressureVariables.H:49
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::polyMesh::TOPO_CHANGE
@ TOPO_CHANGE
Definition: polyMesh.H:92
cellList.H
Foam::polyMesh::oldPoints
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1029
Foam::polyMesh::pointZones_
pointZoneMesh pointZones_
Point zones.
Definition: polyMesh.H:161
Foam::polyMesh::resetPrimitives
void resetPrimitives(const Xfer< pointField > &points, const Xfer< faceList > &faces, const Xfer< labelList > &owner, const Xfer< labelList > &neighbour, const labelList &patchSizes, const labelList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:642
Foam::polyMesh::points_
pointIOField points_
Points.
Definition: polyMesh.H:118
pointZoneMesh.H
Foam::pointZoneMesh.
Foam::globalMeshData
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Definition: globalMeshData.H:106
Foam::polyMesh::checkFaceSkewness
bool checkFaceSkewness(const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check face skewness.
Definition: polyMeshCheck.C:173
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:1204
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:992
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
faceIOList.H
Foam::polyMesh::oldPointsPtr_
autoPtr< pointField > oldPointsPtr_
Old points (for the last mesh motion)
Definition: polyMesh.H:186
Foam::polyMesh::checkFaceOrthogonality
bool checkFaceOrthogonality(const vectorField &fAreas, const vectorField &cellCtrs, const bool report, const bool detailedReport, labelHashSet *setPtr) const
Check non-orthogonality.
Definition: polyMeshCheck.C:34
Foam::polyMesh::meshDir
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:766
Foam::polyMesh::removeBoundary
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
Foam::polyMesh::clearGeom
void clearGeom()
Clear geometry.
Definition: polyMeshClear.C:55
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
Foam::polyMesh::nSolutionD
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:812
Foam::polyMesh::setInstance
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
Foam::polyMesh::bounds
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:427
pointField.H
boundBox.H
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::polyMesh::cellTree
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:843
Foam::polyMesh::~polyMesh
virtual ~polyMesh()
Destructor.
Definition: polyMesh.C:744
Foam::regIOobject
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:60
Foam::polyMesh::findCell
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1411
Foam::polyMesh::CELL_TETS
@ CELL_TETS
Definition: polyMesh.H:107
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::polyMesh::tetBasePtIsPtr_
autoPtr< labelList > tetBasePtIsPtr_
Base point for face decomposition into tets.
Definition: polyMesh.H:152
Foam::polyMesh::pointZones
pointZoneMesh & pointZones()
Return non-const access to the pointZones.
Definition: polyMesh.H:536
Foam::polyMesh::calcDirections
void calcDirections() const
Calculate the valid directions in the mesh from the boundaries.
Definition: polyMesh.C:52
Foam::polyMesh::TypeName
TypeName("polyMesh")
Runtime type information.
Foam::Vector< label >
Foam::polyMesh::addPatches
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:878
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
labelIOList.H
Foam::IOList< label >
Foam::primitiveMesh::pointCells
const labelListList & pointCells() const
Definition: primitiveMeshPointCells.C:108
Foam::polyMesh::moving_
bool moving_
Is the mesh moving.
Definition: polyMesh.H:177
polyBoundaryMesh.H
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::polyMesh::checkEdgeAlignment
bool checkEdgeAlignment(const pointField &p, const bool report, const Vector< label > &directions, labelHashSet *setPtr) const
Definition: polyMeshCheck.C:283
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::polyMesh::solutionD_
Vector< label > solutionD_
Vector of valid directions in mesh.
Definition: polyMesh.H:149
Foam::polyMesh::FACE_PLANES
@ FACE_PLANES
Definition: polyMesh.H:100
Foam::polyMesh::cellZones_
cellZoneMesh cellZones_
Cell zones.
Definition: polyMesh.H:167
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:796
Foam::polyMesh::checkFaceWeight
bool checkFaceWeight(const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs, const bool report, const scalar minWeight, labelHashSet *setPtr) const
Definition: polyMeshCheck.C:498
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
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:1276
Foam::polyMesh::addZones
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:922
Foam::polyMesh::setUpToDatePoints
virtual void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:998
Foam::polyMesh::polyMesh
polyMesh(const polyMesh &)
Disallow construct as copy.
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
Foam::polyMesh::readUpdate
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:66
Foam::polyMesh::initMesh
void initMesh()
Initialise the polyMesh from the primitive data.
Definition: polyMeshInitMesh.C:30
Foam::polyMesh::clearedPrimitives_
bool clearedPrimitives_
Have the primitives been cleared.
Definition: polyMesh.H:130
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:141
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79