checkTopology.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 "checkTopology.H"
27 #include "polyMesh.H"
28 #include "Time.H"
29 #include "regionSplit.H"
30 #include "cellSet.H"
31 #include "faceSet.H"
32 #include "pointSet.H"
33 #include "IOmanip.H"
34 #include "emptyPolyPatch.H"
35 #include "processorPolyPatch.H"
36 #include "surfaceWriter.H"
37 #include "checkTools.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
42 (
43  const polyMesh& mesh,
44  const bool allTopology,
45  const bool allGeometry,
46  const autoPtr<surfaceWriter>& surfWriter,
47  const autoPtr<writer<scalar>>& setWriter
48 )
49 {
50  label noFailedChecks = 0;
51 
52  Info<< "Checking topology..." << endl;
53 
54  // Check if the boundary definition is unique
55  mesh.boundaryMesh().checkDefinition(true);
56 
57  // Check that empty patches cover all sides of the mesh
58  {
59  label nEmpty = 0;
60  forAll(mesh.boundaryMesh(), patchi)
61  {
62  if (isA<emptyPolyPatch>(mesh.boundaryMesh()[patchi]))
63  {
64  nEmpty += mesh.boundaryMesh()[patchi].size();
65  }
66  }
67  reduce(nEmpty, sumOp<label>());
68  label nTotCells = returnReduce(mesh.cells().size(), sumOp<label>());
69 
70  // These are actually warnings, not errors.
71  if (nTotCells && (nEmpty % nTotCells))
72  {
73  Info<< " ***Total number of faces on empty patches"
74  << " is not divisible by the number of cells in the mesh."
75  << " Hence this mesh is not 1D or 2D."
76  << endl;
77  }
78  }
79 
80  // Check if the boundary processor patches are correct
81  mesh.boundaryMesh().checkParallelSync(true);
82 
83  // Check names of zones are equal
84  mesh.cellZones().checkDefinition(true);
85  if (mesh.cellZones().checkParallelSync(true))
86  {
87  noFailedChecks++;
88  }
89  mesh.faceZones().checkDefinition(true);
90  if (mesh.faceZones().checkParallelSync(true))
91  {
92  noFailedChecks++;
93  }
94  mesh.pointZones().checkDefinition(true);
95  if (mesh.pointZones().checkParallelSync(true))
96  {
97  noFailedChecks++;
98  }
99 
100 
101  {
102  cellSet cells(mesh, "illegalCells", mesh.nCells()/100);
103 
104  forAll(mesh.cells(), celli)
105  {
106  const cell& cFaces = mesh.cells()[celli];
107 
108  if (cFaces.size() <= 3)
109  {
110  cells.insert(celli);
111  }
112  forAll(cFaces, i)
113  {
114  if (cFaces[i] < 0 || cFaces[i] >= mesh.nFaces())
115  {
116  cells.insert(celli);
117  break;
118  }
119  }
120  }
121  label nCells = returnReduce(cells.size(), sumOp<label>());
122 
123  if (nCells > 0)
124  {
125  Info<< " Illegal cells (less than 4 faces or out of range faces)"
126  << " found, number of cells: " << nCells << endl;
127  noFailedChecks++;
128 
129  Info<< " <<Writing " << nCells
130  << " illegal cells to set " << cells.name() << endl;
131  cells.instance() = mesh.pointsInstance();
132  cells.write();
133  if (surfWriter.valid())
134  {
135  mergeAndWrite(surfWriter(), cells);
136  }
137 
138  }
139  else
140  {
141  Info<< " Cell to face addressing OK." << endl;
142  }
143  }
144 
145 
146  {
147  pointSet points(mesh, "unusedPoints", mesh.nPoints()/100);
148  if (mesh.checkPoints(true, &points))
149  {
150  noFailedChecks++;
151 
152  label nPoints = returnReduce(points.size(), sumOp<label>());
153 
154  Info<< " <<Writing " << nPoints
155  << " unused points to set " << points.name() << endl;
156  points.instance() = mesh.pointsInstance();
157  points.write();
158  if (setWriter.valid())
159  {
160  mergeAndWrite(setWriter, points);
161  }
162  }
163  }
164 
165  {
166  faceSet faces(mesh, "upperTriangularFace", mesh.nFaces()/100);
167  if (mesh.checkUpperTriangular(true, &faces))
168  {
169  noFailedChecks++;
170  }
171 
172  label nFaces = returnReduce(faces.size(), sumOp<label>());
173 
174  if (nFaces > 0)
175  {
176  Info<< " <<Writing " << nFaces
177  << " unordered faces to set " << faces.name() << endl;
178  faces.instance() = mesh.pointsInstance();
179  faces.write();
180  if (surfWriter.valid())
181  {
182  mergeAndWrite(surfWriter(), faces);
183  }
184  }
185  }
186 
187  {
188  faceSet faces(mesh, "outOfRangeFaces", mesh.nFaces()/100);
189  if (mesh.checkFaceVertices(true, &faces))
190  {
191  noFailedChecks++;
192 
193  label nFaces = returnReduce(faces.size(), sumOp<label>());
194 
195  Info<< " <<Writing " << nFaces
196  << " faces with out-of-range or duplicate vertices to set "
197  << faces.name() << endl;
198  faces.instance() = mesh.pointsInstance();
199  faces.write();
200  if (surfWriter.valid())
201  {
202  mergeAndWrite(surfWriter(), faces);
203  }
204  }
205  }
206 
207  if (allTopology)
208  {
209  cellSet cells(mesh, "zipUpCells", mesh.nCells()/100);
210  if (mesh.checkCellsZipUp(true, &cells))
211  {
212  noFailedChecks++;
213 
214  label nCells = returnReduce(cells.size(), sumOp<label>());
215 
216  Info<< " <<Writing " << nCells
217  << " cells with over used edges to set " << cells.name()
218  << endl;
219  cells.instance() = mesh.pointsInstance();
220  cells.write();
221  if (surfWriter.valid())
222  {
223  mergeAndWrite(surfWriter(), cells);
224  }
225 
226  }
227  }
228 
229  if (allTopology)
230  {
231  faceSet faces(mesh, "edgeFaces", mesh.nFaces()/100);
232  if (mesh.checkFaceFaces(true, &faces))
233  {
234  noFailedChecks++;
235  }
236 
237  label nFaces = returnReduce(faces.size(), sumOp<label>());
238  if (nFaces > 0)
239  {
240  Info<< " <<Writing " << nFaces
241  << " faces with non-standard edge connectivity to set "
242  << faces.name() << endl;
243  faces.instance() = mesh.pointsInstance();
244  faces.write();
245  if (surfWriter.valid())
246  {
247  mergeAndWrite(surfWriter(), faces);
248  }
249  }
250  }
251 
252  if (allTopology)
253  {
254  labelList nInternalFaces(mesh.nCells(), 0);
255 
256  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
257  {
258  nInternalFaces[mesh.faceOwner()[facei]]++;
259  nInternalFaces[mesh.faceNeighbour()[facei]]++;
260  }
261  const polyBoundaryMesh& patches = mesh.boundaryMesh();
262  forAll(patches, patchi)
263  {
264  if (patches[patchi].coupled())
265  {
266  const labelUList& owners = patches[patchi].faceCells();
267 
268  forAll(owners, i)
269  {
270  nInternalFaces[owners[i]]++;
271  }
272  }
273  }
274 
275  cellSet oneCells(mesh, "oneInternalFaceCells", mesh.nCells()/100);
276  cellSet twoCells(mesh, "twoInternalFacesCells", mesh.nCells()/100);
277 
278  forAll(nInternalFaces, celli)
279  {
280  if (nInternalFaces[celli] <= 1)
281  {
282  oneCells.insert(celli);
283  }
284  else if (nInternalFaces[celli] == 2)
285  {
286  twoCells.insert(celli);
287  }
288  }
289 
290  label nOneCells = returnReduce(oneCells.size(), sumOp<label>());
291 
292  if (nOneCells > 0)
293  {
294  Info<< " <<Writing " << nOneCells
295  << " cells with zero or one non-boundary face to set "
296  << oneCells.name()
297  << endl;
298  oneCells.instance() = mesh.pointsInstance();
299  oneCells.write();
300  if (surfWriter.valid())
301  {
302  mergeAndWrite(surfWriter(), oneCells);
303  }
304  }
305 
306  label nTwoCells = returnReduce(twoCells.size(), sumOp<label>());
307 
308  if (nTwoCells > 0)
309  {
310  Info<< " <<Writing " << nTwoCells
311  << " cells with two non-boundary faces to set "
312  << twoCells.name()
313  << endl;
314  twoCells.instance() = mesh.pointsInstance();
315  twoCells.write();
316  if (surfWriter.valid())
317  {
318  mergeAndWrite(surfWriter(), twoCells);
319  }
320  }
321  }
322 
323  {
324  regionSplit rs(mesh);
325 
326  if (rs.nRegions() <= 1)
327  {
328  Info<< " Number of regions: " << rs.nRegions() << " (OK)."
329  << endl;
330 
331  }
332  else
333  {
334  Info<< " *Number of regions: " << rs.nRegions() << endl;
335 
336  Info<< " The mesh has multiple regions which are not connected "
337  "by any face." << endl
338  << " <<Writing region information to "
339  << mesh.time().timeName()/"cellToRegion"
340  << endl;
341 
342  labelIOList ctr
343  (
344  IOobject
345  (
346  "cellToRegion",
347  mesh.time().timeName(),
348  mesh,
349  IOobject::NO_READ,
350  IOobject::NO_WRITE
351  ),
352  rs
353  );
354  ctr.write();
355 
356 
357  // Points in multiple regions
358  pointSet points
359  (
360  mesh,
361  "multiRegionPoints",
362  mesh.nPoints()/1000
363  );
364 
365  // Is region disconnected
366  boolList regionDisconnected(rs.nRegions(), true);
367  if (allTopology)
368  {
369  // -1 : not assigned
370  // -2 : multiple regions
371  // >= 0 : single region
372  labelList pointToRegion(mesh.nPoints(), -1);
373 
374  for
375  (
376  label faceI = mesh.nInternalFaces();
377  faceI < mesh.nFaces();
378  faceI++
379  )
380  {
381  label regionI = rs[mesh.faceOwner()[faceI]];
382  const face& f = mesh.faces()[faceI];
383  forAll(f, fp)
384  {
385  label& pRegion = pointToRegion[f[fp]];
386  if (pRegion == -1)
387  {
388  pRegion = regionI;
389  }
390  else if (pRegion == -2)
391  {
392  // Already marked
393  regionDisconnected[regionI] = false;
394  }
395  else if (pRegion != regionI)
396  {
397  // Multiple regions
398  regionDisconnected[regionI] = false;
399  regionDisconnected[pRegion] = false;
400  pRegion = -2;
401  points.insert(f[fp]);
402  }
403  }
404  }
405 
406  Pstream::listCombineGather(regionDisconnected, andEqOp<bool>());
407  Pstream::listCombineScatter(regionDisconnected);
408  }
409 
410 
411 
412  // write cellSet for each region
413  PtrList<cellSet> cellRegions(rs.nRegions());
414  for (label i = 0; i < rs.nRegions(); i++)
415  {
416  cellRegions.set
417  (
418  i,
419  new cellSet
420  (
421  mesh,
422  "region" + Foam::name(i),
423  mesh.nCells()/100
424  )
425  );
426  }
427 
428  forAll(rs, i)
429  {
430  cellRegions[rs[i]].insert(i);
431  }
432 
433  for (label i = 0; i < rs.nRegions(); i++)
434  {
435  Info<< " <<Writing region " << i;
436  if (allTopology)
437  {
438  if (regionDisconnected[i])
439  {
440  Info<< " (fully disconnected)";
441  }
442  else
443  {
444  Info<< " (point connected)";
445  }
446  }
447  Info<< " with "
448  << returnReduce(cellRegions[i].size(), sumOp<scalar>())
449  << " cells to cellSet " << cellRegions[i].name() << endl;
450 
451  cellRegions[i].write();
452  }
453 
454  label nPoints = returnReduce(points.size(), sumOp<label>());
455  if (nPoints)
456  {
457  Info<< " <<Writing " << nPoints
458  << " points that are in multiple regions to set "
459  << points.name() << endl;
460  points.write();
461  if (setWriter.valid())
462  {
463  mergeAndWrite(setWriter, points);
464  }
465  }
466  }
467  }
468 
469 
470  {
471  if (!Pstream::parRun())
472  {
473  Info<< "\nChecking patch topology for multiply connected"
474  << " surfaces..." << endl;
475  }
476  else
477  {
478  Info<< "\nChecking basic patch addressing..." << endl;
479  }
480 
481 
482  const polyBoundaryMesh& patches = mesh.boundaryMesh();
483 
484  // Non-manifold points
485  pointSet points
486  (
487  mesh,
488  "nonManifoldPoints",
489  mesh.nPoints()/1000
490  );
491 
492  Pout.setf(ios_base::left);
493 
494  Info<< " "
495  << setw(20) << "Patch"
496  << setw(9) << "Faces"
497  << setw(9) << "Points";
498  if (!Pstream::parRun())
499  {
500  Info<< setw(34) << "Surface topology";
501  }
502  if (allGeometry)
503  {
504  Info<< " Bounding box";
505  }
506  Info<< endl;
507 
508  forAll(patches, patchi)
509  {
510  const polyPatch& pp = patches[patchi];
511 
512  if (!isA<processorPolyPatch>(pp))
513  {
514  Info<< " "
515  << setw(20) << pp.name()
516  << setw(9) << returnReduce(pp.size(), sumOp<label>())
517  << setw(9) << returnReduce(pp.nPoints(), sumOp<label>());
518 
519  if (!Pstream::parRun())
520  {
521  primitivePatch::surfaceTopo pTyp = pp.surfaceType();
522 
523  if (pp.empty())
524  {
525  Info<< setw(34) << "ok (empty)";
526  }
527  else if (pTyp == primitivePatch::MANIFOLD)
528  {
529  if (pp.checkPointManifold(true, &points))
530  {
531  Info<< setw(34)
532  << "multiply connected (shared point)";
533  }
534  else
535  {
536  Info<< setw(34) << "ok (closed singly connected)";
537  }
538 
539  // Add points on non-manifold edges to make set complete
540  pp.checkTopology(false, &points);
541  }
542  else
543  {
544  pp.checkTopology(false, &points);
545 
546  if (pTyp == primitivePatch::OPEN)
547  {
548  Info<< setw(34)
549  << "ok (non-closed singly connected)";
550  }
551  else
552  {
553  Info<< setw(34)
554  << "multiply connected (shared edge)";
555  }
556  }
557  }
558 
559  if (allGeometry)
560  {
561  const pointField& pts = pp.points();
562  const labelList& mp = pp.meshPoints();
563 
564  if (returnReduce(mp.size(), sumOp<label>()) > 0)
565  {
566  boundBox bb(point::max, point::min);
567  forAll(mp, i)
568  {
569  bb.min() = min(bb.min(), pts[mp[i]]);
570  bb.max() = max(bb.max(), pts[mp[i]]);
571  }
572  reduce(bb.min(), minOp<vector>());
573  reduce(bb.max(), maxOp<vector>());
574  Info<< ' ' << bb;
575  }
576  }
577  Info<< endl;
578  }
579  }
580 
581  if (points.size())
582  {
583  Info<< " <<Writing " << returnReduce(points.size(), sumOp<label>())
584  << " conflicting points to set "
585  << points.name() << endl;
586 
587  points.instance() = mesh.pointsInstance();
588  points.write();
589  }
590 
591  // Info.setf(ios_base::right);
592  }
593 
594  // Force creation of all addressing if requested.
595  // Errors will be reported as required
596  if (allTopology)
597  {
598  mesh.cells();
599  mesh.faces();
600  mesh.edges();
601  mesh.points();
602  mesh.faceOwner();
603  mesh.faceNeighbour();
604  mesh.cellCells();
605  mesh.edgeCells();
606  mesh.pointCells();
607  mesh.edgeFaces();
608  mesh.pointFaces();
609  mesh.cellEdges();
610  mesh.faceEdges();
611  mesh.pointEdges();
612  }
613 
614  return noFailedChecks;
615 }
616 
617 
618 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
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
ios_base::fmtflags setf(const ios_base::fmtflags f)
Set flags of stream.
Definition: IOstream.H:496
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
UList< label > labelUList
Definition: UList.H:65
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
const cellShapeList & cells
const pointField & points
label checkTopology(const polyMesh &, const bool, const bool, const autoPtr< surfaceWriter > &, const autoPtr< writer< scalar >> &)
List< label > labelList
A List of labels.
Definition: labelList.H:56
void mergeAndWrite(const polyMesh &mesh, const surfaceWriter &writer, const word &name, const indirectPrimitivePatch setPatch, const fileName &outputDir)
Generate merged surface on master and write. Needs input patch.
Istream and Ostream manipulators taking arguments.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42