polyMeshCheck.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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 Description
25  PolyMesh checking routines. Checks various criteria for a mesh and supplied
26  geometry, with the mesh only used for topology. Coupled patch aware (i.e.,
27  coupled faces are treated as internal).
28 
29 SourceFiles
30  polyMeshCheck.C
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef polyMeshCheck_H
35 #define polyMeshCheck_H
36 
37 #include "pointFields.H"
38 #include "HashSet.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace polyMeshCheck
45 {
46 
47 //- Check face non-orthogonality
49 (
50  const bool report,
51  const scalar orthWarn,
52  const polyMesh&,
53  const vectorField& cellCentres,
54  const vectorField& faceAreas,
55  const labelList& checkFaces,
56  const List<labelPair>& baffles,
57  labelHashSet* setPtr
58 );
59 
60 //- Check face pyramid volumes
62 (
63  const bool report,
64  const scalar minPyrVol,
65  const polyMesh&,
66  const vectorField& cellCentres,
67  const pointField& p,
68  const labelList& checkFaces,
69  const List<labelPair>& baffles,
71 );
72 
73 //- Check face tetrahedra volumes
74 bool checkFaceTets
75 (
76  const bool report,
77  const scalar minPyrVol,
78  const polyMesh&,
79  const vectorField& cellCentres,
80  const vectorField& faceCentres,
81  const pointField& p,
82  const labelList& checkFaces,
83  const List<labelPair>& baffles,
85 );
86 
87 //- Check face skewness
89 (
90  const bool report,
91  const scalar internalSkew,
92  const scalar boundarySkew,
93  const polyMesh& mesh,
94  const pointField& points,
95  const vectorField& cellCentres,
96  const vectorField& faceCentres,
97  const vectorField& faceAreas,
98  const labelList& checkFaces,
99  const List<labelPair>& baffles,
100  labelHashSet* setPtr
101 );
102 
103 //- Check interpolation weights (0.5 for regular mesh)
104 bool checkFaceWeights
105 (
106  const bool report,
107  const scalar warnWeight,
108  const polyMesh& mesh,
109  const vectorField& cellCentres,
110  const vectorField& faceCentres,
111  const vectorField& faceAreas,
112  const labelList& checkFaces,
113  const List<labelPair>& baffles,
114  labelHashSet* setPtr
115 );
116 
117 //- Cell volume ratio of neighbouring cells (1 for regular mesh)
118 bool checkVolRatio
119 (
120  const bool report,
121  const scalar warnRatio,
122  const polyMesh& mesh,
123  const scalarField& cellVolumes,
124  const labelList& checkFaces,
125  const List<labelPair>& baffles,
126  labelHashSet* setPtr
127 );
128 
129 //- Check convexity of angles in a face. See primitiveMesh for explanation.
130 bool checkFaceAngles
131 (
132  const bool report,
133  const scalar maxDeg,
134  const polyMesh& mesh,
135  const vectorField& faceAreas,
136  const pointField& p,
137  const labelList& checkFaces,
138  labelHashSet* setPtr
139 );
140 
141 // Check the difference between normals of individual face-triangles (from
142 // face-centre decomposition) and the cell-cell centre vector
143 bool checkFaceTwist
144 (
145  const bool report,
146  const scalar minTwist,
147  const polyMesh&,
148  const vectorField& cellCentres,
149  const vectorField& faceAreas,
150  const vectorField& faceCentres,
151  const pointField& p,
152  const labelList& checkFaces,
153  labelHashSet* setPtr
154 );
155 
156 //- Check consecutive face-triangle (from face-centre decomposition) normals
158 (
159  const bool report,
160  const scalar minTwist,
161  const polyMesh&,
162  const vectorField& faceAreas,
163  const vectorField& faceCentres,
164  const pointField& p,
165  const labelList& checkFaces,
166  labelHashSet* setPtr
167 );
168 
169 //- Check for face areas v.s. sum of face-triangle (from face-centre
170 // decomposition) areas
172 (
173  const bool report,
174  const scalar minFlatness,
175  const polyMesh&,
176  const vectorField& faceAreas,
177  const vectorField& faceCentres,
178  const pointField& p,
179  const labelList& checkFaces,
180  labelHashSet* setPtr
181 );
182 
183 //- Check for small faces
184 bool checkFaceArea
185 (
186  const bool report,
187  const scalar minArea,
188  const polyMesh&,
189  const vectorField& faceAreas,
190  const labelList& checkFaces,
191  labelHashSet* setPtr
192 );
193 
194 //- Check the area of internal faces v.s. boundary faces
196 (
197  const bool report,
198  const scalar minDet,
199  const polyMesh&,
200  const vectorField& faceAreas,
201  const labelList& checkFaces,
202  labelHashSet* setPtr
203 );
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace polyMeshCheck
208 } // End namespace Foam
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 #endif
213 
214 // ************************************************************************* //
bool checkTriangleTwist(const bool report, const scalar minTwist, const polyMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Check consecutive face-triangle (from face-centre decomposition) normals.
bool checkFaceArea(const bool report, const scalar minArea, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Check for small faces.
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)
Check face non-orthogonality.
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)
Check face skewness.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
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 *)
Check face tetrahedra volumes.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
const pointField & points
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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)
bool checkFaceFlatness(const bool report, const scalar minFlatness, const polyMesh &, const vectorField &faceAreas, const vectorField &faceCentres, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Check for face areas v.s. sum of face-triangle (from face-centre.
bool checkFaceAngles(const bool report, const scalar maxDeg, const polyMesh &mesh, const vectorField &faceAreas, const pointField &p, const labelList &checkFaces, labelHashSet *setPtr)
Check convexity of angles in a face. See primitiveMesh for explanation.
bool checkCellDeterminant(const bool report, const scalar minDet, const polyMesh &, const vectorField &faceAreas, const labelList &checkFaces, labelHashSet *setPtr)
Check the area of internal faces v.s. boundary faces.
bool checkFacePyramids(const bool report, const scalar minPyrVol, const polyMesh &, const vectorField &cellCentres, const pointField &p, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *)
Check face pyramid volumes.
Field< vector > vectorField
Specialisation of Field<T> for vector.
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)
Check interpolation weights (0.5 for regular mesh)
volScalarField & p
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)
Namespace for OpenFOAM.