blockMeshCheck.C
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-2018 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 \*---------------------------------------------------------------------------*/
25 
26 #include "blockMesh.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 void Foam::blockMesh::check(const polyMesh& bm, const dictionary& dict) const
31 {
32  Info<< nl << "Check topology" << endl;
33 
34  bool ok = true;
35 
36  // Check for duplicate curved edge definitions
37  forAll(edges_, cei)
38  {
39  for (label cej=cei+1; cej<edges_.size(); cej++)
40  {
41  if (edges_[cei].compare(edges_[cej]) != 0)
42  {
43  Info<< " Curved edge ";
44  edges_[cej].write(Info, dict);
45  Info<< " is a duplicate of curved edge " << edges_[cei]
46  << endl;
47  ok = false;
48  break;
49  }
50  }
51  }
52 
53  // Check curved-edge/block-edge correspondence
54  //
55  // Loop over the edges of each block rather than the edgeList of the
56  // topological mesh due to problems with calcEdges for blocks with
57  // repeated point labels
58  const blockList& blocks = *this;
59 
60  forAll(edges_, cei)
61  {
62  bool found = false;
63 
64  forAll(blocks, blocki)
65  {
66  edgeList edges = blocks[blocki].blockShape().edges();
67 
68  forAll(edges, ei)
69  {
70  found = edges_[cei].compare(edges[ei][0], edges[ei][1]) != 0;
71  if (found) break;
72  }
73  if (found) break;
74  }
75 
76  if (!found)
77  {
78  Info<< " Curved edge ";
79  edges_[cei].write(Info, dict);
80  Info<< " does not correspond to a block edge."
81  << endl;
82  ok = false;
83  }
84  }
85 
86  const faceList& faces = bm.faces();
87 
88  // Check for duplicate curved face definitions
89  forAll(faces_, cfi)
90  {
91  for (label cfj=cfi+1; cfj<faces_.size(); cfj++)
92  {
93  if (faces_[cfi].compare(faces_[cfj]) != 0)
94  {
95  Info<< " Curved face ";
96  faces_[cfj].write(Info, dict);
97  Info<< " is a duplicate of curved face ";
98  faces_[cfi].write(Info, dict);
99  Info<< endl;
100  ok = false;
101  break;
102  }
103  }
104  }
105 
106  // Check curved-face/block-face correspondence
107  forAll(faces_, cfi)
108  {
109  bool found = false;
110 
111  forAll(faces, fi)
112  {
113  found = faces_[cfi].compare(faces[fi]) != 0;
114  if (found) break;
115  }
116 
117  if (!found)
118  {
119  Info<< " Curved face ";
120  faces_[cfi].write(Info, dict);
121  Info<< " does not correspond to a block face." << endl;
122  ok = false;
123  }
124  }
125 
126  const pointField& points = bm.points();
127  const cellList& cells = bm.cells();
128  const polyPatchList& patches = bm.boundaryMesh();
129 
130  label nBoundaryFaces = 0;
131  forAll(cells, celli)
132  {
133  nBoundaryFaces += cells[celli].nFaces();
134  }
135 
136  nBoundaryFaces -= 2*bm.nInternalFaces();
137 
138  label nDefinedBoundaryFaces = 0;
139  forAll(patches, patchi)
140  {
141  nDefinedBoundaryFaces += patches[patchi].size();
142  }
143 
144 
145  if (verboseOutput)
146  {
147  Info<< nl << tab << "Basic statistics" << nl
148  << tab << tab << "Number of internal faces : "
149  << bm.nInternalFaces() << nl
150  << tab << tab << "Number of boundary faces : "
151  << nBoundaryFaces << nl
152  << tab << tab << "Number of defined boundary faces : "
153  << nDefinedBoundaryFaces << nl
154  << tab << tab << "Number of undefined boundary faces : "
155  << nBoundaryFaces - nDefinedBoundaryFaces << nl;
156 
157  if ((nBoundaryFaces - nDefinedBoundaryFaces) > 0)
158  {
159  Info<< tab << tab << tab
160  << "(Warning : only leave undefined the front and back planes "
161  << "of 2D planar geometries!)" << endl;
162  }
163 
164  Info<< tab << "Checking patch -> block consistency" << endl;
165  }
166 
167 
168  forAll(patches, patchi)
169  {
170  const faceList& Patch = patches[patchi];
171 
172  forAll(Patch, patchFacei)
173  {
174  const face& patchFace = Patch[patchFacei];
175  bool patchFaceOK = false;
176 
177  forAll(cells, celli)
178  {
179  const labelList& cellFaces = cells[celli];
180 
181  forAll(cellFaces, cellFacei)
182  {
183  if (patchFace == faces[cellFaces[cellFacei]])
184  {
185  patchFaceOK = true;
186 
187  if
188  (
189  (
190  patchFace.area(points)
191  & faces[cellFaces[cellFacei]].area(points)
192  ) < 0.0
193  )
194  {
195  Info<< tab << tab
196  << "Face " << patchFacei
197  << " of patch " << patchi
198  << " (" << patches[patchi].name() << ")"
199  << " points inwards"
200  << endl;
201 
202  ok = false;
203  }
204  }
205  }
206  }
207 
208  if (!patchFaceOK)
209  {
210  Info<< tab << tab
211  << "Face " << patchFacei
212  << " of patch " << patchi
213  << " (" << patches[patchi].name() << ")"
214  << " does not match any block faces" << endl;
215 
216  ok = false;
217  }
218  }
219  }
220 
221  if (verboseOutput)
222  {
223  Info<< endl;
224  }
225 
226  if (!ok)
227  {
229  << "Block mesh topology incorrect, stopping mesh generation!"
230  << exit(FatalError);
231  }
232 }
233 
234 
235 // ************************************************************************* //
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const char tab
Definition: Ostream.H:264
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
List< face > faceList
Definition: faceListFwd.H:43
PtrList< block > blockList
A PtrList of blocks.
Definition: blockList.H:45
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< edge > edgeList
Definition: edgeList.H:38
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const char nl
Definition: Ostream.H:265
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label patchi
messageStream Info
List< cell > cellList
list of cells
Definition: cellList.H:42