initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / utilities / mesh / generation / blockMesh / checkBlockMesh.C
blob3d093a46dc825f08bb8d605b0faa1149b8230b6e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 Description
27 \*---------------------------------------------------------------------------*/
29 #include "blockMesh.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace Foam
36 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
38 // Check the blockMesh topology
40 void blockMesh::checkBlockMesh(const polyMesh& bm)
42     Info<< nl << "Check block mesh topology" << endl;
44     bool blockMeshOK = true;
46     const pointField& points = bm.points();
47     const faceList& faces = bm.faces();
48     const cellList& cells = bm.cells();
49     const polyPatchList& patches = bm.boundaryMesh();
51     label nBoundaryFaces=0;
52     forAll(cells, celli)
53     {
54         nBoundaryFaces += cells[celli].nFaces();
55     }
57     nBoundaryFaces -= 2*bm.nInternalFaces();
59     label nDefinedBoundaryFaces=0;
60     forAll(patches, patchi)
61     {
62         nDefinedBoundaryFaces += patches[patchi].size();
63     }
66     Info<< nl << tab << "Basic statistics" << endl;
68     Info<< tab << tab << "Number of internal faces : "
69         << bm.nInternalFaces() << endl;
71     Info<< tab << tab << "Number of boundary faces : "
72         << nBoundaryFaces << endl;
74     Info<< tab << tab << "Number of defined boundary faces : "
75         << nDefinedBoundaryFaces << endl;
77     Info<< tab << tab << "Number of undefined boundary faces : "
78         << nBoundaryFaces - nDefinedBoundaryFaces << endl;
80     if ((nBoundaryFaces - nDefinedBoundaryFaces) > 0)
81     {
82         Info<< tab << tab << tab
83             << "(Warning : only leave undefined the front and back planes "
84             << "of 2D planar geometries!)" << endl;
85     }
87     Info<< nl << tab << "Checking patch -> block consistency" << endl;
90     forAll(patches, patchi)
91     {
92         const faceList& Patch = patches[patchi];
94         forAll(Patch, patchFacei)
95         {
96             const face& patchFace = Patch[patchFacei];
97             bool patchFaceOK = false;
99             forAll(cells, celli)
100             {
101                 const labelList& cellFaces = cells[celli];
103                 forAll(cellFaces, cellFacei)
104                 {
105                     if (patchFace == faces[cellFaces[cellFacei]])
106                     {
107                         patchFaceOK = true;
109                         if
110                         (
111                             (
112                                 patchFace.normal(points)
113                               & faces[cellFaces[cellFacei]].normal(points)
114                             ) < 0.0
115                         )
116                         {
117                             Info<< tab << tab
118                                 << "Face " << patchFacei
119                                 << " of patch " << patchi
120                                 << " (" << patches[patchi].name() << ")"
121                                 << " points inwards"
122                                 << endl;
124                             blockMeshOK = false;
125                         }
126                     }
127                 }
128             }
130             if (!patchFaceOK)
131             {
132                 Info<< tab << tab
133                     << "Face " << patchFacei
134                     << " of patch " << patchi
135                     << " (" << patches[patchi].name() << ")"
136                     << " does not match any block faces" << endl;
138                 blockMeshOK = false;
139             }
140         }
141     }
143     if (!blockMeshOK)
144     {
145         FatalErrorIn("blockMesh::checkBlockMesh(const polyMesh& bm)")
146             << "Block mesh topology incorrect, stopping mesh generation!"
147             << exit(FatalError);
148     }
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 } // End namespace Foam
156 // ************************************************************************* //