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