FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

polyMeshGeometry.H

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Class
00025     Foam::polyMeshGeometry
00026 
00027 Description
00028     Updateable mesh geometry and checking routines.
00029 
00030     - non-ortho done across coupled faces.
00031     - faceWeight (delta factors) done across coupled faces.
00032 
00033 SourceFiles
00034     polyMeshGeometry.C
00035 
00036 \*---------------------------------------------------------------------------*/
00037 
00038 #ifndef polyMeshGeometry_H
00039 #define polyMeshGeometry_H
00040 
00041 #include <OpenFOAM/pointFields.H>
00042 #include <OpenFOAM/HashSet.H>
00043 
00044 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00045 
00046 namespace Foam
00047 {
00048 
00049 /*---------------------------------------------------------------------------*\
00050                            Class polyMeshGeometry Declaration
00051 \*---------------------------------------------------------------------------*/
00052 
00053 class polyMeshGeometry
00054 {
00055         //- Reference to polyMesh.
00056         const polyMesh& mesh_;
00057 
00058         //- Uptodate copy of face areas
00059         vectorField faceAreas_;
00060 
00061         //- Uptodate copy of face centres
00062         vectorField faceCentres_;
00063 
00064         //- Uptodate copy of cell centres
00065         vectorField cellCentres_;
00066 
00067         //- Uptodate copy of cell volumes
00068         scalarField cellVolumes_;
00069 
00070 
00071     // Private Member Functions
00072 
00073         //- Update face areas and centres on selected faces.
00074         void updateFaceCentresAndAreas
00075         (
00076             const pointField& p,
00077             const labelList& changedFaces
00078         );
00079 
00080         //- Update cell volumes and centres on selected cells. Requires
00081         //  cells and faces to be consistent set.
00082         void updateCellCentresAndVols
00083         (
00084             const labelList& changedCells,
00085             const labelList& changedFaces
00086         );
00087 
00088         //- Detect&report non-ortho error for single face.
00089         static scalar checkNonOrtho
00090         (
00091             const polyMesh& mesh,
00092             const bool report,
00093             const scalar severeNonorthogonalityThreshold,
00094             const label faceI,
00095             const vector& s,    // face area vector
00096             const vector& d,    // cc-cc vector
00097 
00098             label& severeNonOrth,
00099             label& errorNonOrth,
00100             labelHashSet* setPtr
00101         );
00102 
00103         //- Calculate skewness given two cell centres and one face centre.
00104         static scalar calcSkewness
00105         (
00106             const point& ownCc,
00107             const point& neiCc,
00108             const point& fc
00109         );
00110 
00111 public:
00112 
00113     ClassName("polyMeshGeometry");
00114 
00115     // Constructors
00116 
00117         //- Construct from mesh
00118         polyMeshGeometry(const polyMesh&);
00119 
00120 
00121     // Member Functions
00122 
00123         // Access
00124 
00125             const polyMesh& mesh() const
00126             {
00127                 return mesh_;
00128             }
00129 
00130             const vectorField& faceAreas() const
00131             {
00132                 return faceAreas_;
00133             }
00134             const vectorField& faceCentres() const
00135             {
00136                 return faceCentres_;
00137             }
00138             const vectorField& cellCentres() const
00139             {
00140                 return cellCentres_;
00141             }
00142             const scalarField& cellVolumes() const
00143             {
00144                 return cellVolumes_;
00145             }
00146 
00147         // Edit
00148 
00149             //- Take over properties from mesh
00150             void correct();
00151 
00152             //- Recalculate on selected faces. Recalculates cell properties
00153             //  on owner and neighbour of these cells.
00154             void correct
00155             (
00156                 const pointField& p,
00157                 const labelList& changedFaces
00158             );
00159 
00160             //- Helper function: get affected cells from faces
00161             static labelList affectedCells
00162             (
00163                 const polyMesh&,
00164                 const labelList& changedFaces
00165             );
00166 
00167 
00168         // Checking of selected faces with supplied geometry (mesh only used for
00169         // topology). Coupled aware (coupled faces treated as internal ones)
00170 
00171             //- See primitiveMesh
00172             static bool checkFaceDotProduct
00173             (
00174                 const bool report,
00175                 const scalar orthWarn,
00176                 const polyMesh&,
00177                 const vectorField& cellCentres,
00178                 const vectorField& faceAreas,
00179                 const labelList& checkFaces,
00180                 const List<labelPair>& baffles,
00181                 labelHashSet* setPtr
00182             );
00183 
00184             //- See primitiveMesh
00185             static bool checkFacePyramids
00186             (
00187                 const bool report,
00188                 const scalar minPyrVol,
00189                 const polyMesh&,
00190                 const vectorField& cellCentres,
00191                 const pointField& p,
00192                 const labelList& checkFaces,
00193                 const List<labelPair>& baffles,
00194                 labelHashSet*
00195             );
00196 
00197             //- See primitiveMesh
00198             static bool checkFaceSkewness
00199             (
00200                 const bool report,
00201                 const scalar internalSkew,
00202                 const scalar boundarySkew,
00203                 const polyMesh& mesh,
00204                 const vectorField& cellCentres,
00205                 const vectorField& faceCentres,
00206                 const vectorField& faceAreas,
00207                 const labelList& checkFaces,
00208                 const List<labelPair>& baffles,
00209                 labelHashSet* setPtr
00210             );
00211 
00212             //- Interpolation weights (0.5 for regular mesh)
00213             static bool checkFaceWeights
00214             (
00215                 const bool report,
00216                 const scalar warnWeight,
00217                 const polyMesh& mesh,
00218                 const vectorField& cellCentres,
00219                 const vectorField& faceCentres,
00220                 const vectorField& faceAreas,
00221                 const labelList& checkFaces,
00222                 const List<labelPair>& baffles,
00223                 labelHashSet* setPtr
00224             );
00225 
00226             //- Cell volume ratio of neighbouring cells (1 for regular mesh)
00227             static bool checkVolRatio
00228             (
00229                 const bool report,
00230                 const scalar warnRatio,
00231                 const polyMesh& mesh,
00232                 const scalarField& cellVolumes,
00233                 const labelList& checkFaces,
00234                 const List<labelPair>& baffles,
00235                 labelHashSet* setPtr
00236             );
00237 
00238             //- See primitiveMesh
00239             static bool checkFaceAngles
00240             (
00241                 const bool report,
00242                 const scalar maxDeg,
00243                 const polyMesh& mesh,
00244                 const vectorField& faceAreas,
00245                 const pointField& p,
00246                 const labelList& checkFaces,
00247                 labelHashSet* setPtr
00248             );
00249 
00250             //- Triangle (from face-centre decomposition) normal v.s.
00251             //  average face normal
00252             static bool checkFaceTwist
00253             (
00254                 const bool report,
00255                 const scalar minTwist,
00256                 const polyMesh&,
00257                 const vectorField& cellCentres,
00258                 const vectorField& faceAreas,
00259                 const vectorField& faceCentres,
00260                 const pointField& p,
00261                 const labelList& checkFaces,
00262                 labelHashSet* setPtr
00263             );
00264 
00265             //- Consecutive triangle (from face-centre decomposition) normals
00266             static bool checkTriangleTwist
00267             (
00268                 const bool report,
00269                 const scalar minTwist,
00270                 const polyMesh&,
00271                 const vectorField& faceAreas,
00272                 const vectorField& faceCentres,
00273                 const pointField& p,
00274                 const labelList& checkFaces,
00275                 labelHashSet* setPtr
00276             );
00277 
00278             //- Small faces
00279             static bool checkFaceArea
00280             (
00281                 const bool report,
00282                 const scalar minArea,
00283                 const polyMesh&,
00284                 const vectorField& faceAreas,
00285                 const labelList& checkFaces,
00286                 labelHashSet* setPtr
00287             );
00288 
00289             static bool checkCellDeterminant
00290             (
00291                 const bool report,
00292                 const scalar minDet,
00293                 const polyMesh&,
00294                 const vectorField& faceAreas,
00295                 const labelList& checkFaces,
00296                 const labelList& affectedCells,
00297                 labelHashSet* setPtr
00298             );
00299 
00300 
00301         // Checking of selected faces with local geometry. Uses above static
00302         // functions. Parallel aware.
00303 
00304             bool checkFaceDotProduct
00305             (
00306                 const bool report,
00307                 const scalar orthWarn,
00308                 const labelList& checkFaces,
00309                 const List<labelPair>& baffles,
00310                 labelHashSet* setPtr
00311             ) const;
00312 
00313             bool checkFacePyramids
00314             (
00315                 const bool report,
00316                 const scalar minPyrVol,
00317                 const pointField& p,
00318                 const labelList& checkFaces,
00319                 const List<labelPair>& baffles,
00320                 labelHashSet* setPtr
00321             ) const;
00322 
00323             bool checkFaceSkewness
00324             (
00325                 const bool report,
00326                 const scalar internalSkew,
00327                 const scalar boundarySkew,
00328                 const labelList& checkFaces,
00329                 const List<labelPair>& baffles,
00330                 labelHashSet* setPtr
00331             ) const;
00332 
00333             bool checkFaceWeights
00334             (
00335                 const bool report,
00336                 const scalar warnWeight,
00337                 const labelList& checkFaces,
00338                 const List<labelPair>& baffles,
00339                 labelHashSet* setPtr
00340             ) const;
00341 
00342             bool checkVolRatio
00343             (
00344                 const bool report,
00345                 const scalar warnRatio,
00346                 const labelList& checkFaces,
00347                 const List<labelPair>& baffles,
00348                 labelHashSet* setPtr
00349             ) const;
00350 
00351             bool checkFaceAngles
00352             (
00353                 const bool report,
00354                 const scalar maxDeg,
00355                 const pointField& p,
00356                 const labelList& checkFaces,
00357                 labelHashSet* setPtr
00358             ) const;
00359 
00360             bool checkFaceTwist
00361             (
00362                 const bool report,
00363                 const scalar minTwist,
00364                 const pointField& p,
00365                 const labelList& checkFaces,
00366                 labelHashSet* setPtr
00367             ) const;
00368 
00369             bool checkTriangleTwist
00370             (
00371                 const bool report,
00372                 const scalar minTwist,
00373                 const pointField& p,
00374                 const labelList& checkFaces,
00375                 labelHashSet* setPtr
00376             ) const;
00377 
00378             bool checkFaceArea
00379             (
00380                 const bool report,
00381                 const scalar minArea,
00382                 const labelList& checkFaces,
00383                 labelHashSet* setPtr
00384             ) const;
00385 
00386             bool checkCellDeterminant
00387             (
00388                 const bool report,
00389                 const scalar warnDet,
00390                 const labelList& checkFaces,
00391                 const labelList& affectedCells,
00392                 labelHashSet* setPtr
00393             ) const;
00394 };
00395 
00396 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00397 
00398 } // End namespace Foam
00399 
00400 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00401 
00402 #endif
00403 
00404 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines