initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / dynamicMesh / motionSmoother / polyMeshGeometry / polyMeshGeometry.H
blob91f27e6a001137aceb5d532baa76e7a014594289
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::polyMeshGeometry
28 Description
29     Updateable mesh geometry and checking routines.
31     - non-ortho done across coupled faces.
32     - faceWeight (delta factors) done across coupled faces.
34 SourceFiles
35     polyMeshGeometry.C
37 \*---------------------------------------------------------------------------*/
39 #ifndef polyMeshGeometry_H
40 #define polyMeshGeometry_H
42 #include "pointFields.H"
43 #include "HashSet.H"
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 namespace Foam
50 /*---------------------------------------------------------------------------*\
51                            Class polyMeshGeometry Declaration
52 \*---------------------------------------------------------------------------*/
54 class polyMeshGeometry
56         //- Reference to polyMesh.
57         const polyMesh& mesh_;
59         //- Uptodate copy of face areas
60         vectorField faceAreas_;
62         //- Uptodate copy of face centres
63         vectorField faceCentres_;
65         //- Uptodate copy of cell centres
66         vectorField cellCentres_;
68         //- Uptodate copy of cell volumes
69         scalarField cellVolumes_;
72     // Private Member Functions
74         //- Update face areas and centres on selected faces.
75         void updateFaceCentresAndAreas
76         (
77             const pointField& p,
78             const labelList& changedFaces
79         );
81         //- Update cell volumes and centres on selected cells. Requires
82         //  cells and faces to be consistent set.
83         void updateCellCentresAndVols
84         (
85             const labelList& changedCells,
86             const labelList& changedFaces
87         );
89         //- Detect&report non-ortho error for single face.
90         static scalar checkNonOrtho
91         (
92             const polyMesh& mesh,
93             const bool report,
94             const scalar severeNonorthogonalityThreshold,
95             const label faceI,
96             const vector& s,    // face area vector
97             const vector& d,    // cc-cc vector
99             label& severeNonOrth,
100             label& errorNonOrth,
101             labelHashSet* setPtr
102         );
104         //- Calculate skewness given two cell centres and one face centre.
105         static scalar calcSkewness
106         (
107             const point& ownCc,
108             const point& neiCc,
109             const point& fc
110         );
112 public:
114     ClassName("polyMeshGeometry");
116     // Constructors
118         //- Construct from mesh
119         polyMeshGeometry(const polyMesh&);
122     // Member Functions
124         // Access
126             const polyMesh& mesh() const
127             {
128                 return mesh_;
129             }
131             const vectorField& faceAreas() const
132             {
133                 return faceAreas_;
134             }
135             const vectorField& faceCentres() const
136             {
137                 return faceCentres_;
138             }
139             const vectorField& cellCentres() const
140             {
141                 return cellCentres_;
142             }
143             const scalarField& cellVolumes() const
144             {
145                 return cellVolumes_;
146             }
148         // Edit
150             //- Take over properties from mesh
151             void correct();
153             //- Recalculate on selected faces. Recalculates cell properties
154             //  on owner and neighbour of these cells.
155             void correct
156             (
157                 const pointField& p,
158                 const labelList& changedFaces
159             );
161             //- Helper function: get affected cells from faces
162             static labelList affectedCells
163             (
164                 const polyMesh&,
165                 const labelList& changedFaces
166             );
169         // Checking of selected faces with supplied geometry (mesh only used for
170         // topology). Coupled aware (coupled faces treated as internal ones)
172             //- See primitiveMesh
173             static bool checkFaceDotProduct
174             (
175                 const bool report,
176                 const scalar orthWarn,
177                 const polyMesh&,
178                 const vectorField& cellCentres,
179                 const vectorField& faceAreas,
180                 const labelList& checkFaces,
181                 const List<labelPair>& baffles,
182                 labelHashSet* setPtr
183             );
185             //- See primitiveMesh
186             static bool checkFacePyramids
187             (
188                 const bool report,
189                 const scalar minPyrVol,
190                 const polyMesh&,
191                 const vectorField& cellCentres,
192                 const pointField& p,
193                 const labelList& checkFaces,
194                 const List<labelPair>& baffles,
195                 labelHashSet*
196             );
198             //- See primitiveMesh
199             static bool checkFaceSkewness
200             (
201                 const bool report,
202                 const scalar internalSkew,
203                 const scalar boundarySkew,
204                 const polyMesh& mesh,
205                 const vectorField& cellCentres,
206                 const vectorField& faceCentres,
207                 const vectorField& faceAreas,
208                 const labelList& checkFaces,
209                 const List<labelPair>& baffles,
210                 labelHashSet* setPtr
211             );
213             //- Interpolation weights (0.5 for regular mesh)
214             static bool checkFaceWeights
215             (
216                 const bool report,
217                 const scalar warnWeight,
218                 const polyMesh& mesh,
219                 const vectorField& cellCentres,
220                 const vectorField& faceCentres,
221                 const vectorField& faceAreas,
222                 const labelList& checkFaces,
223                 const List<labelPair>& baffles,
224                 labelHashSet* setPtr
225             );
227             //- Cell volume ratio of neighbouring cells (1 for regular mesh)
228             static bool checkVolRatio
229             (
230                 const bool report,
231                 const scalar warnRatio,
232                 const polyMesh& mesh,
233                 const scalarField& cellVolumes,
234                 const labelList& checkFaces,
235                 const List<labelPair>& baffles,
236                 labelHashSet* setPtr
237             );
239             //- See primitiveMesh
240             static bool checkFaceAngles
241             (
242                 const bool report,
243                 const scalar maxDeg,
244                 const polyMesh& mesh,
245                 const vectorField& faceAreas,
246                 const pointField& p,
247                 const labelList& checkFaces,
248                 labelHashSet* setPtr
249             );
251             //- Triangle (from face-centre decomposition) normal v.s.
252             //  average face normal
253             static bool checkFaceTwist
254             (
255                 const bool report,
256                 const scalar minTwist,
257                 const polyMesh&,
258                 const vectorField& cellCentres,
259                 const vectorField& faceAreas,
260                 const vectorField& faceCentres,
261                 const pointField& p,
262                 const labelList& checkFaces,
263                 labelHashSet* setPtr
264             );
266             //- Consecutive triangle (from face-centre decomposition) normals
267             static bool checkTriangleTwist
268             (
269                 const bool report,
270                 const scalar minTwist,
271                 const polyMesh&,
272                 const vectorField& faceAreas,
273                 const vectorField& faceCentres,
274                 const pointField& p,
275                 const labelList& checkFaces,
276                 labelHashSet* setPtr
277             );
279             //- Small faces
280             static bool checkFaceArea
281             (
282                 const bool report,
283                 const scalar minArea,
284                 const polyMesh&,
285                 const vectorField& faceAreas,
286                 const labelList& checkFaces,
287                 labelHashSet* setPtr
288             );
290             static bool checkCellDeterminant
291             (
292                 const bool report,
293                 const scalar minDet,
294                 const polyMesh&,
295                 const vectorField& faceAreas,
296                 const labelList& checkFaces,
297                 const labelList& affectedCells,
298                 labelHashSet* setPtr
299             );
302         // Checking of selected faces with local geometry. Uses above static
303         // functions. Parallel aware.
305             bool checkFaceDotProduct
306             (
307                 const bool report,
308                 const scalar orthWarn,
309                 const labelList& checkFaces,
310                 const List<labelPair>& baffles,
311                 labelHashSet* setPtr
312             ) const;
314             bool checkFacePyramids
315             (
316                 const bool report,
317                 const scalar minPyrVol,
318                 const pointField& p,
319                 const labelList& checkFaces,
320                 const List<labelPair>& baffles,
321                 labelHashSet* setPtr
322             ) const;
324             bool checkFaceSkewness
325             (
326                 const bool report,
327                 const scalar internalSkew,
328                 const scalar boundarySkew,
329                 const labelList& checkFaces,
330                 const List<labelPair>& baffles,
331                 labelHashSet* setPtr
332             ) const;
334             bool checkFaceWeights
335             (
336                 const bool report,
337                 const scalar warnWeight,
338                 const labelList& checkFaces,
339                 const List<labelPair>& baffles,
340                 labelHashSet* setPtr
341             ) const;
343             bool checkVolRatio
344             (
345                 const bool report,
346                 const scalar warnRatio,
347                 const labelList& checkFaces,
348                 const List<labelPair>& baffles,
349                 labelHashSet* setPtr
350             ) const;
352             bool checkFaceAngles
353             (
354                 const bool report,
355                 const scalar maxDeg,
356                 const pointField& p,
357                 const labelList& checkFaces,
358                 labelHashSet* setPtr
359             ) const;
361             bool checkFaceTwist
362             (
363                 const bool report,
364                 const scalar minTwist,
365                 const pointField& p,
366                 const labelList& checkFaces,
367                 labelHashSet* setPtr
368             ) const;
370             bool checkTriangleTwist
371             (
372                 const bool report,
373                 const scalar minTwist,
374                 const pointField& p,
375                 const labelList& checkFaces,
376                 labelHashSet* setPtr
377             ) const;
379             bool checkFaceArea
380             (
381                 const bool report,
382                 const scalar minArea,
383                 const labelList& checkFaces,
384                 labelHashSet* setPtr
385             ) const;
387             bool checkCellDeterminant
388             (
389                 const bool report,
390                 const scalar warnDet,
391                 const labelList& checkFaces,
392                 const labelList& affectedCells,
393                 labelHashSet* setPtr
394             ) const;
397 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
399 } // End namespace Foam
401 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403 #endif
405 // ************************************************************************* //