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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "blockMesh.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
32 void Foam::blockMesh::check(const polyMesh& bm, const dictionary& dict) const
33 {
34  Info<< nl << "Check topology" << endl;
35 
36  bool ok = true;
37 
38  // Check for duplicate curved edge definitions
39  forAll(edges_, cei)
40  {
41  for (label cej=cei+1; cej<edges_.size(); cej++)
42  {
43  if (edges_[cei].compare(edges_[cej]) != 0)
44  {
45  Info<< " Curved edge ";
46  edges_[cej].write(Info, dict);
47  Info<< " is a duplicate of curved edge " << edges_[cei]
48  << endl;
49  ok = false;
50  break;
51  }
52  }
53  }
54 
55  // Check curved-edge/block-edge correspondence
56  //
57  // Loop over the edges of each block rather than the edgeList of the
58  // topological mesh due to problems with calcEdges for blocks with
59  // repeated point labels
60  const blockList& blocks = *this;
61 
62  forAll(edges_, cei)
63  {
64  bool found = false;
65 
66  forAll(blocks, blocki)
67  {
68  edgeList edges = blocks[blocki].blockShape().edges();
69 
70  forAll(edges, ei)
71  {
72  found = edges_[cei].compare(edges[ei][0], edges[ei][1]) != 0;
73  if (found) break;
74  }
75  if (found) break;
76  }
77 
78  if (!found)
79  {
80  Info<< " Curved edge ";
81  edges_[cei].write(Info, dict);
82  Info<< " does not correspond to a block edge." << endl;
83  ok = false;
84  }
85  }
86 
87  const faceList& faces = bm.faces();
88 
89  // Check for duplicate curved face definitions
90  forAll(faces_, cfi)
91  {
92  for (label cfj=cfi+1; cfj<faces_.size(); cfj++)
93  {
94  if (faces_[cfi].compare(faces_[cfj]) != 0)
95  {
96  Info<< " Curved face ";
97  faces_[cfj].write(Info, dict);
98  Info<< " is a duplicate of curved face ";
99  faces_[cfi].write(Info, dict);
100  Info<< endl;
101  ok = false;
102  break;
103  }
104  }
105  }
106 
107  // Check curved-face/block-face correspondence
108  forAll(faces_, cfi)
109  {
110  bool found = false;
111 
112  forAll(faces, fi)
113  {
114  found = faces_[cfi].compare(faces[fi]) != 0;
115  if (found) break;
116  }
117 
118  if (!found)
119  {
120  Info<< " Curved face ";
121  faces_[cfi].write(Info, dict);
122  Info<< " does not correspond to a block face." << endl;
123  ok = false;
124  }
125  }
126 
127  const pointField& points = bm.points();
128  const cellList& cells = bm.cells();
129  const polyPatchList& patches = bm.boundaryMesh();
130 
131  label nBoundaryFaces = 0;
132  forAll(cells, celli)
133  {
134  nBoundaryFaces += cells[celli].nFaces();
135  }
136 
137  nBoundaryFaces -= 2*bm.nInternalFaces();
138 
139  label nDefinedBoundaryFaces = 0;
140  forAll(patches, patchi)
141  {
142  nDefinedBoundaryFaces += patches[patchi].size();
143  }
144 
145 
146  if (verbose_)
147  {
148  Info<< nl << tab << "Basic statistics" << nl
149  << tab << tab << "Number of internal faces : "
150  << bm.nInternalFaces() << nl
151  << tab << tab << "Number of boundary faces : "
152  << nBoundaryFaces << nl
153  << tab << tab << "Number of defined boundary faces : "
154  << nDefinedBoundaryFaces << nl
155  << tab << tab << "Number of undefined boundary faces : "
156  << nBoundaryFaces - nDefinedBoundaryFaces << nl;
157 
158  if ((nBoundaryFaces - nDefinedBoundaryFaces) > 0)
159  {
160  Info<< tab << tab << tab
161  << "(Warning : only leave undefined the front and back planes "
162  << "of 2D planar geometries!)" << endl;
163  }
164 
165  Info<< tab << "Checking patch -> block consistency" << endl;
166  }
167 
168 
169  forAll(patches, patchi)
170  {
171  const faceList& Patch = patches[patchi];
172 
173  forAll(Patch, patchFacei)
174  {
175  const face& patchFace = Patch[patchFacei];
176  bool patchFaceOK = false;
177 
178  forAll(cells, celli)
179  {
180  const labelList& cellFaces = cells[celli];
181 
182  forAll(cellFaces, cellFacei)
183  {
184  if (patchFace == faces[cellFaces[cellFacei]])
185  {
186  patchFaceOK = true;
187 
188  if
189  (
190  (
191  patchFace.areaNormal(points)
192  & faces[cellFaces[cellFacei]].areaNormal(points)
193  ) < 0
194  )
195  {
196  Info<< tab << tab
197  << "Face " << patchFacei
198  << " of patch " << patchi
199  << " (" << patches[patchi].name() << ")"
200  << " points inwards"
201  << endl;
202 
203  ok = false;
204  }
205  }
206  }
207  }
208 
209  if (!patchFaceOK)
210  {
211  Info<< tab << tab
212  << "Face " << patchFacei
213  << " of patch " << patchi
214  << " (" << patches[patchi].name() << ")"
215  << " does not match any block faces" << endl;
216 
217  ok = false;
218  }
219  }
220  }
221 
222  if (verbose_)
223  {
224  Info<< endl;
225  }
226 
227  if (!ok)
228  {
230  << "Block mesh topology incorrect, stopping mesh generation!"
231  << exit(FatalError);
232  }
233 }
234 
235 
236 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::blockMesh::patches
const faceListList & patches() const
Return the patch face lists.
Definition: blockMesh.C:198
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
blockMesh.H
Foam::blockMesh::faces
const blockFaceList & faces() const
Return the curved faces.
Definition: blockMesh.H:351
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
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
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::blockMesh::blockList
PtrList< block > blockList
The list of blocks is stored as a PtrList.
Definition: blockMesh.H:157
found
bool found
Definition: TABSMDCalcMethod2.H:32
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::tab
constexpr char tab
Definition: Ostream.H:384
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::blockMesh::points
const pointField & points() const
The points for the entire mesh.
Definition: blockMesh.C:176
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::blockMesh::edges
const blockEdgeList & edges() const
Return the curved edges.
Definition: blockMesh.H:345
Foam::blockMesh::cells
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:187