polyMeshGeometry.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-2016 OpenFOAM Foundation
9  Copyright (C) 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::polyMeshGeometry
29 
30 Description
31  Updateable mesh geometry and checking routines.
32 
33  - non-ortho done across coupled faces.
34  - faceWeight (delta factors) done across coupled faces.
35 
36 SourceFiles
37  polyMeshGeometry.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef polyMeshGeometry_H
42 #define polyMeshGeometry_H
43 
44 #include "pointFields.H"
45 #include "HashSet.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class polyMeshGeometry Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 class polyMeshGeometry
57 {
58  //- Reference to polyMesh.
59  const polyMesh& mesh_;
60 
61  //- Uptodate copy of face areas
62  vectorField faceAreas_;
63 
64  //- Uptodate copy of face centres
65  vectorField faceCentres_;
66 
67  //- Uptodate copy of cell centres
68  vectorField cellCentres_;
69 
70  //- Uptodate copy of cell volumes
71  scalarField cellVolumes_;
72 
73 
74  // Private Member Functions
75 
76  //- Update face areas and centres on selected faces.
77  void updateFaceCentresAndAreas
78  (
79  const pointField& p,
80  const labelList& changedFaces
81  );
82 
83  //- Update cell volumes and centres on selected cells. Requires
84  // cells and faces to be consistent set.
85  void updateCellCentresAndVols
86  (
87  const labelList& changedCells,
88  const labelList& changedFaces
89  );
90 
91  //- Detect&report non-ortho error for single face.
92  static scalar checkNonOrtho
93  (
94  const polyMesh& mesh,
95  const bool report,
96  const scalar severeNonorthogonalityThreshold,
97  const label facei,
98  const vector& s, // face area vector
99  const vector& d, // cc-cc vector
100 
101  label& severeNonOrth,
102  label& errorNonOrth,
103  labelHashSet* setPtr
104  );
105 
106  //- Calculate skewness given two cell centres and one face centre.
107  static scalar calcSkewness
108  (
109  const point& ownCc,
110  const point& neiCc,
111  const point& fc
112  );
113 
114  //- Detect&report incorrect face-triangle orientation
115  static bool checkFaceTet
116  (
117  const polyMesh&,
118  const bool report,
119  const scalar minTetQuality,
120  const pointField& p,
121  const label facei,
122  const point& fc, // face centre
123  const point& cc, // cell centre
124  labelHashSet* setPtr
125  );
126 
127 
128 public:
129 
130  ClassName("polyMeshGeometry");
131 
132  // Constructors
133 
134  //- Construct from mesh
135  explicit polyMeshGeometry(const polyMesh&);
136 
137 
138  // Member Functions
139 
140  // Access
141 
142  const polyMesh& mesh() const
143  {
144  return mesh_;
145  }
146 
147  const vectorField& faceAreas() const
148  {
149  return faceAreas_;
150  }
151  const vectorField& faceCentres() const
152  {
153  return faceCentres_;
154  }
155  const vectorField& cellCentres() const
156  {
157  return cellCentres_;
158  }
159  const scalarField& cellVolumes() const
160  {
161  return cellVolumes_;
162  }
163 
164  // Edit
165 
166  //- Take over properties from mesh
167  void correct();
168 
169  //- Recalculate on selected faces. Recalculates cell properties
170  // on owner and neighbour of these cells.
171  void correct
172  (
173  const pointField& p,
174  const labelList& changedFaces
175  );
176 
177  //- Helper function: get affected cells from faces
178  static labelList affectedCells
179  (
180  const polyMesh&,
181  const labelList& changedFaces
182  );
183 
184 
185  // Checking of selected faces with supplied geometry (mesh only used for
186  // topology). Coupled aware (coupled faces treated as internal ones)
187 
188  //- See primitiveMesh
189  static bool checkFaceDotProduct
190  (
191  const bool report,
192  const scalar orthWarn,
193  const polyMesh&,
194  const vectorField& cellCentres,
195  const vectorField& faceAreas,
196  const labelList& checkFaces,
197  const List<labelPair>& baffles,
198  labelHashSet* setPtr
199  );
200 
201  //- See primitiveMesh
202  static bool checkFacePyramids
203  (
204  const bool report,
205  const scalar minPyrVol,
206  const polyMesh&,
207  const vectorField& cellCentres,
208  const pointField& p,
209  const labelList& checkFaces,
210  const List<labelPair>& baffles,
211  labelHashSet*
212  );
213 
214  //- See primitiveMesh
215  static bool checkFaceTets
216  (
217  const bool report,
218  const scalar minPyrVol,
219  const polyMesh&,
220  const vectorField& cellCentres,
221  const vectorField& faceCentres,
222  const pointField& p,
223  const labelList& checkFaces,
224  const List<labelPair>& baffles,
225  labelHashSet*
226  );
227 
228  //- See primitiveMesh
229  static bool checkFaceSkewness
230  (
231  const bool report,
232  const scalar internalSkew,
233  const scalar boundarySkew,
234  const polyMesh& mesh,
235  const pointField& points,
236  const vectorField& cellCentres,
237  const vectorField& faceCentres,
238  const vectorField& faceAreas,
239  const labelList& checkFaces,
240  const List<labelPair>& baffles,
241  labelHashSet* setPtr
242  );
243 
244  //- Interpolation weights (0.5 for regular mesh)
245  static bool checkFaceWeights
246  (
247  const bool report,
248  const scalar warnWeight,
249  const polyMesh& mesh,
250  const vectorField& cellCentres,
251  const vectorField& faceCentres,
252  const vectorField& faceAreas,
253  const labelList& checkFaces,
254  const List<labelPair>& baffles,
255  labelHashSet* setPtr
256  );
257 
258  //- Cell volume ratio of neighbouring cells (1 for regular mesh)
259  static bool checkVolRatio
260  (
261  const bool report,
262  const scalar warnRatio,
263  const polyMesh& mesh,
264  const scalarField& cellVolumes,
265  const labelList& checkFaces,
266  const List<labelPair>& baffles,
267  labelHashSet* setPtr
268  );
269 
270  //- See primitiveMesh
271  static bool checkFaceAngles
272  (
273  const bool report,
274  const scalar maxDeg,
275  const polyMesh& mesh,
276  const vectorField& faceAreas,
277  const pointField& p,
278  const labelList& checkFaces,
279  labelHashSet* setPtr
280  );
281 
282  //- Triangle (from face-centre decomposition) normal v.s.
283  // average face normal
284  static bool checkFaceTwist
285  (
286  const bool report,
287  const scalar minTwist,
288  const polyMesh&,
289  const vectorField& cellCentres,
290  const vectorField& faceAreas,
291  const vectorField& faceCentres,
292  const pointField& p,
293  const labelList& checkFaces,
294  labelHashSet* setPtr
295  );
296 
297  //- Consecutive triangle (from face-centre decomposition) normals
298  static bool checkTriangleTwist
299  (
300  const bool report,
301  const scalar minTwist,
302  const polyMesh&,
303  const vectorField& faceAreas,
304  const vectorField& faceCentres,
305  const pointField& p,
306  const labelList& checkFaces,
307  labelHashSet* setPtr
308  );
309 
310  //- Area of faces v.s. sum of triangle areas
311  static bool checkFaceFlatness
312  (
313  const bool report,
314  const scalar minFlatness,
315  const polyMesh&,
316  const vectorField& faceAreas,
317  const vectorField& faceCentres,
318  const pointField& p,
319  const labelList& checkFaces,
320  labelHashSet* setPtr
321  );
322 
323  //- Small faces
324  static bool checkFaceArea
325  (
326  const bool report,
327  const scalar minArea,
328  const polyMesh&,
329  const vectorField& faceAreas,
330  const labelList& checkFaces,
331  labelHashSet* setPtr
332  );
333 
334  //- Area of internal faces v.s. boundary faces
335  static bool checkCellDeterminant
336  (
337  const bool report,
338  const scalar minDet,
339  const polyMesh&,
340  const vectorField& faceAreas,
341  const labelList& checkFaces,
342  const labelList& affectedCells,
343  labelHashSet* setPtr
344  );
345 
346 
347  // Checking of selected faces with local geometry. Uses above static
348  // functions. Parallel aware.
349 
351  (
352  const bool report,
353  const scalar orthWarn,
354  const labelList& checkFaces,
355  const List<labelPair>& baffles,
356  labelHashSet* setPtr
357  ) const;
358 
359  bool checkFacePyramids
360  (
361  const bool report,
362  const scalar minPyrVol,
363  const pointField& p,
364  const labelList& checkFaces,
365  const List<labelPair>& baffles,
366  labelHashSet* setPtr
367  ) const;
368 
369  bool checkFaceTets
370  (
371  const bool report,
372  const scalar minTetQuality,
373  const pointField& p,
374  const labelList& checkFaces,
375  const List<labelPair>& baffles,
376  labelHashSet* setPtr
377  ) const;
378 
379  bool checkFaceWeights
380  (
381  const bool report,
382  const scalar warnWeight,
383  const labelList& checkFaces,
384  const List<labelPair>& baffles,
385  labelHashSet* setPtr
386  ) const;
387 
388  bool checkVolRatio
389  (
390  const bool report,
391  const scalar warnRatio,
392  const labelList& checkFaces,
393  const List<labelPair>& baffles,
394  labelHashSet* setPtr
395  ) const;
396 
397  bool checkFaceAngles
398  (
399  const bool report,
400  const scalar maxDeg,
401  const pointField& p,
402  const labelList& checkFaces,
403  labelHashSet* setPtr
404  ) const;
405 
406  bool checkFaceTwist
407  (
408  const bool report,
409  const scalar minTwist,
410  const pointField& p,
411  const labelList& checkFaces,
412  labelHashSet* setPtr
413  ) const;
414 
415  bool checkTriangleTwist
416  (
417  const bool report,
418  const scalar minTwist,
419  const pointField& p,
420  const labelList& checkFaces,
421  labelHashSet* setPtr
422  ) const;
423 
424  bool checkFaceFlatness
425  (
426  const bool report,
427  const scalar minFlatness,
428  const pointField& p,
429  const labelList& checkFaces,
430  labelHashSet* setPtr
431  ) const;
432 
433  bool checkFaceArea
434  (
435  const bool report,
436  const scalar minArea,
437  const labelList& checkFaces,
438  labelHashSet* setPtr
439  ) const;
440 
442  (
443  const bool report,
444  const scalar warnDet,
445  const labelList& checkFaces,
446  const labelList& affectedCells,
447  labelHashSet* setPtr
448  ) const;
449 };
450 
451 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
452 
453 } // End namespace Foam
454 
455 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
456 
457 #endif
458 
459 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMeshGeometry::checkFaceDotProduct
static bool checkFaceDotProduct(const bool report, const scalar orthWarn, const polyMesh &, const vectorField &cellCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
See primitiveMesh.
Definition: polyMeshGeometry.C:360
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::polyMeshGeometry::ClassName
ClassName("polyMeshGeometry")
Foam::polyMeshGeometry::faceAreas
const vectorField & faceAreas() const
Definition: polyMeshGeometry.H:146
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::polyMeshGeometry::cellVolumes
const scalarField & cellVolumes() const
Definition: polyMeshGeometry.H:158
Foam::polyMeshGeometry::checkFaceWeights
static bool checkFaceWeights(const bool report, const scalar warnWeight, const polyMesh &mesh, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
Interpolation weights (0.5 for regular mesh)
Definition: polyMeshGeometry.C:1135
Foam::polyMeshGeometry::correct
void correct()
Take over properties from mesh.
Definition: polyMeshGeometry.C:337
Foam::polyMeshGeometry::checkFaceTets
static bool checkFaceTets(const bool report, const scalar minPyrVol, const polyMesh &, const vectorField &cellCentres, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *)
See primitiveMesh.
Definition: polyMeshGeometry.C:709
Foam::polyMeshGeometry::polyMeshGeometry
polyMeshGeometry(const polyMesh &)
Construct from mesh.
Definition: polyMeshGeometry.C:327
Foam::HashSet< label, Hash< label > >
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::polyMeshGeometry::faceCentres
const vectorField & faceCentres() const
Definition: polyMeshGeometry.H:150
Foam::polyMeshGeometry::checkFaceSkewness
static bool checkFaceSkewness(const bool report, const scalar internalSkew, const scalar boundarySkew, const polyMesh &mesh, const pointField &points, const vectorField &cellCentres, const vectorField &faceCentres, const vectorField &faceAreas, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
See primitiveMesh.
Definition: polyMeshGeometry.C:930
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::polyMeshGeometry::affectedCells
static labelList affectedCells(const polyMesh &, const labelList &changedFaces)
Helper function: get affected cells from faces.
Definition: polyMeshGeometry.C:183
Foam::polyMeshGeometry::checkFaceArea
static bool checkFaceArea(const bool report, const scalar minArea, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Small faces.
Definition: polyMeshGeometry.C:1932
Foam::Field< vector >
Foam::polyMeshGeometry::mesh
const polyMesh & mesh() const
Definition: polyMeshGeometry.H:141
HashSet.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::polyMeshGeometry::checkFaceTwist
static bool checkFaceTwist(const bool report, const scalar minTwist, const polyMesh &, const vectorField &cellCentres, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Triangle (from face-centre decomposition) normal v.s.
Definition: polyMeshGeometry.C:1547
Foam::polyMeshGeometry::checkFaceFlatness
static bool checkFaceFlatness(const bool report, const scalar minFlatness, const polyMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Area of faces v.s. sum of triangle areas.
Definition: polyMeshGeometry.C:1839
Foam::Vector< scalar >
Foam::List< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::polyMeshGeometry
Updateable mesh geometry and checking routines.
Definition: polyMeshGeometry.H:55
Foam::polyMeshGeometry::checkFacePyramids
static bool checkFacePyramids(const bool report, const scalar minPyrVol, const polyMesh &, const vectorField &cellCentres, const pointField &p, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *)
See primitiveMesh.
Definition: polyMeshGeometry.C:538
Foam::polyMeshGeometry::checkCellDeterminant
static bool checkCellDeterminant(const bool report, const scalar minDet, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, const labelList &affectedCells, labelHashSet *setPtr)
Area of internal faces v.s. boundary faces.
Definition: polyMeshGeometry.C:1989
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::polyMeshGeometry::checkTriangleTwist
static bool checkTriangleTwist(const bool report, const scalar minTwist, const polyMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Consecutive triangle (from face-centre decomposition) normals.
Definition: polyMeshGeometry.C:1693
Foam::polyMeshGeometry::cellCentres
const vectorField & cellCentres() const
Definition: polyMeshGeometry.H:154
pointFields.H
Foam::polyMeshGeometry::checkVolRatio
static bool checkVolRatio(const bool report, const scalar warnRatio, const polyMesh &mesh, const scalarField &cellVolumes, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *setPtr)
Cell volume ratio of neighbouring cells (1 for regular mesh)
Definition: polyMeshGeometry.C:1283
Foam::polyMeshGeometry::checkFaceAngles
static bool checkFaceAngles(const bool report, const scalar maxDeg, const polyMesh &mesh, const vectorField &faceAreas, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
See primitiveMesh.
Definition: polyMeshGeometry.C:1419