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

PatchToolsCheck.C

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 \*---------------------------------------------------------------------------*/
00025 
00026 #include "PatchTools.H"
00027 
00028 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00029 
00030 template
00031 <
00032     class Face,
00033     template<class> class FaceList,
00034     class PointField,
00035     class PointType
00036 >
00037 
00038 bool
00039 Foam::PatchTools::checkOrientation
00040 (
00041     const PrimitivePatch<Face, FaceList, PointField, PointType>& p,
00042     const bool report,
00043     labelHashSet* setPtr
00044 )
00045 {
00046     bool foundError = false;
00047 
00048     // Check edge normals, face normals, point normals.
00049     forAll(p.faceEdges(), faceI)
00050     {
00051         const labelList& edgeLabels = p.faceEdges()[faceI];
00052         bool valid = true;
00053 
00054         if (edgeLabels.size() < 3)
00055         {
00056             if (report)
00057             {
00058                 Info<< "Face[" << faceI << "] " << p[faceI]
00059                     << " has fewer than 3 edges. Edges: " << edgeLabels
00060                     << endl;
00061             }
00062             valid = false;
00063         }
00064         else
00065         {
00066             forAll(edgeLabels, i)
00067             {
00068                 if (edgeLabels[i] < 0 || edgeLabels[i] >= p.nEdges())
00069                 {
00070                     if (report)
00071                     {
00072                         Info<< "edge number " << edgeLabels[i]
00073                             << " on face " << faceI
00074                             << " out-of-range\n"
00075                             << "This usually means the input surface has "
00076                             << "edges with more than 2 faces connected."
00077                             << endl;
00078                     }
00079                     valid = false;
00080                 }
00081             }
00082         }
00083 
00084         if (!valid)
00085         {
00086             foundError = true;
00087             continue;
00088         }
00089 
00090 
00091         //
00092         //- Compute normal from 3 points, use the first as the origin
00093         // minor warpage should not be a problem
00094         const Face& f = p[faceI];
00095         const point& p0 = p.points()[f[0]];
00096         const point& p1 = p.points()[f[1]];
00097         const point& p2 = p.points()[f[f.size()-1]];
00098 
00099         const vector pointNormal((p1 - p0) ^ (p2 - p0));
00100         if ((pointNormal & p.faceNormals()[faceI]) < 0)
00101         {
00102             foundError = false;
00103 
00104             if (report)
00105             {
00106                 Info
00107                     << "Normal calculated from points inconsistent"
00108                     << " with faceNormal" << nl
00109                     << "face: " << f << nl
00110                     << "points: " << p0 << ' ' << p1 << ' ' << p2 << nl
00111                     << "pointNormal:" << pointNormal << nl
00112                     << "faceNormal:" << p.faceNormals()[faceI] << endl;
00113             }
00114         }
00115     }
00116 
00117 
00118     forAll(p.edges(), edgeI)
00119     {
00120         const edge& e = p.edges()[edgeI];
00121         const labelList& neighbouringFaces = p.edgeFaces()[edgeI];
00122 
00123         if (neighbouringFaces.size() == 2)
00124         {
00125             // we use localFaces() since edges() are LOCAL
00126             // these are both already available
00127             const Face& faceA = p.localFaces()[neighbouringFaces[0]];
00128             const Face& faceB = p.localFaces()[neighbouringFaces[1]];
00129 
00130             // If the faces are correctly oriented, the edges must go in
00131             // different directions on connected faces.
00132             if (faceA.edgeDirection(e) == faceB.edgeDirection(e))
00133             {
00134                 if (report)
00135                 {
00136                     Info<< "face orientation incorrect." << nl
00137                         << "localEdge[" << edgeI << "] " << e
00138                         << " between faces:" << nl
00139                         << "  face[" << neighbouringFaces[0] << "] "
00140                         << p[neighbouringFaces[0]]
00141                         << "   localFace: " << faceA
00142                         << nl
00143                         << "  face[" << neighbouringFaces[1] << "] "
00144                         << p[neighbouringFaces[1]]
00145                         << "   localFace: " << faceB
00146                         << endl;
00147                 }
00148                 if (setPtr)
00149                 {
00150                     setPtr->insert(edgeI);
00151                 }
00152 
00153                 foundError = true;
00154             }
00155         }
00156         else if (neighbouringFaces.size() != 1)
00157         {
00158             if (report)
00159             {
00160                 Info
00161                     << "Wrong number of edge neighbours." << nl
00162                     << "edge[" << edgeI << "] " << e
00163                     << " with points:" << p.localPoints()[e.start()]
00164                     << ' ' << p.localPoints()[e.end()]
00165                     << " has neighbouringFaces:" << neighbouringFaces << endl;
00166             }
00167 
00168             if (setPtr)
00169             {
00170                 setPtr->insert(edgeI);
00171             }
00172 
00173             foundError = true;
00174         }
00175     }
00176 
00177     return foundError;
00178 }
00179 
00180 
00181 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines