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