blockMeshCheck.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "blockMesh.H"
30 #include "IOmanip.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 bool Foam::blockMesh::checkDegenerate() const
35 {
36  const blockList& blocks = *this;
37 
38  for (const block& blk : blocks)
39  {
40  const cellShape& shape = blk.blockShape();
41 
42  if (shape.model().index() == cellModel::HEX)
43  {
44  // Check for collapsed edges
45  // - limit to HEX only for now.
46  for (label edgei = 0; edgei < shape.nEdges(); ++edgei)
47  {
48  edge e(shape.edge(edgei));
49 
50  if (!e.valid())
51  {
52  return true; // Looks like a collapsed edge
53  }
54  }
55  }
56  }
57 
58  return false;
59 }
60 
61 
62 void Foam::blockMesh::check
63 (
64  const polyMesh& bm,
65  const dictionary& dict
66 ) const
67 {
68  Info<< nl << "Check topology" << endl;
69 
70  bool ok = true;
71 
72  // Check for duplicate curved edge definitions
73  forAll(edges_, cei)
74  {
75  for (label cej=cei+1; cej<edges_.size(); cej++)
76  {
77  if (edges_[cei].compare(edges_[cej]) != 0)
78  {
79  Info<< " Curved edge ";
80  edges_[cej].write(Info, dict);
81  Info<< " is a duplicate of curved edge "
82  << edges_[cei] << endl;
83  ok = false;
84  break;
85  }
86  }
87  }
88 
89  // Check curved-edge/block-edge correspondence
90  //
91  // Loop over the edges of each block rather than the edgeList of the
92  // topological mesh due to problems with calcEdges for blocks with
93  // repeated point labels
94  const blockList& blocks = *this;
95 
96  for (const blockEdge& curvEdge : edges_)
97  {
98  bool found = false;
99 
100  for (const block& blk : blocks)
101  {
102  for (const edge& blkEdge : blk.blockShape().edges())
103  {
104  found = curvEdge.compare(blkEdge) != 0;
105  if (found) break;
106  }
107  if (found) break;
108  }
109 
110  if (!found)
111  {
112  Info<< " Curved edge ";
113  curvEdge.write(Info, dict);
114  Info<< " does not correspond to a block edge." << endl;
115  ok = false;
116  }
117  }
118 
119  const faceList& faces = bm.faces();
120 
121  // Check for duplicate curved face definitions
122  forAll(faces_, cfi)
123  {
124  for (label cfj=cfi+1; cfj<faces_.size(); cfj++)
125  {
126  if (faces_[cfi].compare(faces_[cfj]) != 0)
127  {
128  Info<< " Curved face ";
129  faces_[cfj].write(Info, dict);
130  Info<< " is a duplicate of curved face ";
131  faces_[cfi].write(Info, dict);
132  Info<< endl;
133  ok = false;
134  break;
135  }
136  }
137  }
138 
139  // Check curved-face/block-face correspondence
140  for (const blockFace& cface : faces_)
141  {
142  const face& cf = cface.vertices();
143 
144  bool found = false;
145 
146  if (cf.size() == 2)
147  {
148  const label bi = cf[0];
149  const label fi = cf[1];
150 
151  found =
152  (
153  bi >= 0 && bi < blocks.size()
154  && fi >= 0 && fi < blocks[bi].blockShape().nFaces()
155  );
156  }
157 
158  if (!found)
159  {
160  for (const face& bf : faces)
161  {
162  found = cface.compare(bf) != 0;
163  if (found) break;
164  }
165  }
166 
167  if (!found)
168  {
169  Info<< " Curved face ";
170  cface.write(Info, dict);
171  Info<< " does not correspond to a block face." << endl;
172  ok = false;
173  }
174  }
175 
176 
177  const pointField& points = bm.points();
178  const cellList& cells = bm.cells();
179  const polyPatchList& patches = bm.boundaryMesh();
180 
181  label nBoundaryFaces = 0;
182  for (const cell& c : cells)
183  {
184  nBoundaryFaces += c.nFaces();
185  }
186  nBoundaryFaces -= 2*bm.nInternalFaces();
187 
188  label nDefinedBoundaryFaces = 0;
189  for (const polyPatch& pp : patches)
190  {
191  nDefinedBoundaryFaces += pp.size();
192  }
193 
194 
195  if (verbose_)
196  {
197  Info<< nl << tab << "Basic statistics" << nl
198  << tab << tab << "Number of internal faces : "
199  << bm.nInternalFaces() << nl
200  << tab << tab << "Number of boundary faces : "
201  << nBoundaryFaces << nl
202  << tab << tab << "Number of defined boundary faces : "
203  << nDefinedBoundaryFaces << nl
204  << tab << tab << "Number of undefined boundary faces : "
205  << nBoundaryFaces - nDefinedBoundaryFaces << nl;
206 
207  if ((nBoundaryFaces - nDefinedBoundaryFaces) > 0)
208  {
209  Info<< tab << tab << tab
210  << "(Warning : only leave undefined the front and back planes "
211  << "of 2D planar geometries!)" << endl;
212  }
213 
214  Info<< tab << "Checking patch -> block consistency" << endl;
215  }
216 
217 
218  for (const polyPatch& pp : patches)
219  {
220  forAll(pp, patchFacei)
221  {
222  const face& patchFace = pp[patchFacei];
223 
224  bool patchFaceOK = false;
225 
226  for (const labelList& cellFaces : cells)
227  {
228  for (const label cellFacei : cellFaces)
229  {
230  const face& cellFace = faces[cellFacei];
231 
232  if (patchFace == cellFace)
233  {
234  patchFaceOK = true;
235 
236  if
237  (
238  (
239  patchFace.areaNormal(points)
240  & cellFace.areaNormal(points)
241  ) < 0
242  )
243  {
244  Info<< tab << tab
245  << "Face " << patchFacei
246  << " of patch " << pp.index()
247  << " (" << pp.name() << ')'
248  << " points inwards"
249  << endl;
250 
251  ok = false;
252  }
253  }
254  }
255  }
256 
257  if (!patchFaceOK)
258  {
259  Info<< tab << tab
260  << "Face " << patchFacei
261  << " of patch " << pp.index()
262  << " (" << pp.name() << ')'
263  << " does not match any block faces" << endl;
264 
265  ok = false;
266  }
267  }
268  }
269 
270  if (verbose_)
271  {
272  Info<< endl;
273  }
274 
275  // Report patch block/face correspondence
276  if (verbose_ > 1)
277  {
278  const labelList& own = bm.faceOwner();
279 
280  Info.stream().setf(ios_base::left);
281 
282  Info<< setw(20) << "patch" << "block/face" << nl
283  << setw(20) << "-----" << "----------" << nl;
284 
285  for (const polyPatch& pp : patches)
286  {
287  Info<< setw(20) << pp.name();
288 
289  label meshFacei = pp.start();
290 
291  forAll(pp, bfacei)
292  {
293  const label celli = own[meshFacei];
294  const label cellFacei = cells[celli].find(meshFacei);
295 
296  if (bfacei) Info<< token::SPACE;
297 
299  << celli << ' ' << cellFacei
300  << token::END_LIST;
301 
302  ++meshFacei;
303  }
304  Info<< nl;
305  }
306 
307  Info<< setw(20) << "-----" << "----------" << nl
308  << nl;
309  }
310 
311  if (!ok)
312  {
314  << "Block mesh topology incorrect, stopping mesh generation!"
315  << exit(FatalError);
316  }
317 }
318 
319 
320 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::cellModel::HEX
hex
Definition: cellModel.H:81
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
blockMesh.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::messageStream::stream
OSstream & stream(OSstream *alternative=nullptr)
Definition: messageStream.C:71
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatchList
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:47
Foam::cellList
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
IOmanip.H
Istream and Ostream manipulators taking arguments.
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::IOstream::setf
ios_base::fmtflags setf(const ios_base::fmtflags f)
Set flags of stream.
Definition: IOstream.H:378
Foam::blockMesh::blockList
PtrList< block > blockList
The list of blocks is stored as a PtrList.
Definition: blockMesh.H:174
found
bool found
Definition: TABSMDCalcMethod2.H:32
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::tab
constexpr char tab
Definition: Ostream.H:403
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::token::SPACE
Space [isspace].
Definition: token.H:125
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::token::END_LIST
End list [isseparator].
Definition: token.H:156
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::blockList
PtrList< block > blockList
A PtrList of blocks.
Definition: blockList.H:47
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::token::BEGIN_LIST
Begin list [isseparator].
Definition: token.H:155