checkTopology.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) 2017-2020 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 "checkTopology.H"
30 #include "polyMesh.H"
31 #include "Time.H"
32 #include "regionSplit.H"
33 #include "cellSet.H"
34 #include "faceSet.H"
35 #include "pointSet.H"
36 #include "IOmanip.H"
37 #include "emptyPolyPatch.H"
38 #include "processorPolyPatch.H"
39 #include "vtkSurfaceWriter.H"
40 #include "checkTools.H"
41 #include "treeBoundBox.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 template<class PatchType>
47 (
48  const bool allGeometry,
49  const word& name,
50  const PatchType& pp,
51  pointSet& points
52 )
53 {
54  Info<< " "
55  << setw(20) << name
56  << setw(9) << returnReduce(pp.size(), sumOp<label>())
57  << setw(9) << returnReduce(pp.nPoints(), sumOp<label>());
58 
59  if (!Pstream::parRun())
60  {
61  typedef typename PatchType::surfaceTopo TopoType;
62  TopoType pTyp = pp.surfaceType();
63 
64  if (pp.empty())
65  {
66  Info<< setw(34) << "ok (empty)";
67  }
68  else if (pTyp == TopoType::MANIFOLD)
69  {
70  if (pp.checkPointManifold(true, &points))
71  {
72  Info<< setw(34)
73  << "multiply connected (shared point)";
74  }
75  else
76  {
77  Info<< setw(34) << "ok (closed singly connected)";
78  }
79 
80  // Add points on non-manifold edges to make set complete
81  pp.checkTopology(false, &points);
82  }
83  else
84  {
85  pp.checkTopology(false, &points);
86 
87  if (pTyp == TopoType::OPEN)
88  {
89  Info<< setw(34)
90  << "ok (non-closed singly connected)";
91  }
92  else
93  {
94  Info<< setw(34)
95  << "multiply connected (shared edge)";
96  }
97  }
98  }
99 
100  if (allGeometry)
101  {
102  const labelList& mp = pp.meshPoints();
103 
104  if (returnReduce(mp.size(), sumOp<label>()) > 0)
105  {
106  boundBox bb(pp.points(), mp, true); // reduce
107  Info<< ' ' << bb;
108  }
109  }
110 }
111 
112 
113 Foam::label Foam::checkTopology
114 (
115  const polyMesh& mesh,
116  const bool allTopology,
117  const bool allGeometry,
118  autoPtr<surfaceWriter>& surfWriter,
119  const autoPtr<writer<scalar>>& setWriter
120 )
121 {
122  label noFailedChecks = 0;
123 
124  Info<< "Checking topology..." << endl;
125 
126  // Check if the boundary definition is unique
127  mesh.boundaryMesh().checkDefinition(true);
128 
129  // Check that empty patches cover all sides of the mesh
130  {
131  label nEmpty = 0;
132  forAll(mesh.boundaryMesh(), patchi)
133  {
134  if (isA<emptyPolyPatch>(mesh.boundaryMesh()[patchi]))
135  {
136  nEmpty += mesh.boundaryMesh()[patchi].size();
137  }
138  }
139  reduce(nEmpty, sumOp<label>());
140  label nTotCells = returnReduce(mesh.cells().size(), sumOp<label>());
141 
142  // These are actually warnings, not errors.
143  if (nTotCells && (nEmpty % nTotCells))
144  {
145  Info<< " ***Total number of faces on empty patches"
146  << " is not divisible by the number of cells in the mesh."
147  << " Hence this mesh is not 1D or 2D."
148  << endl;
149  }
150  }
151 
152  // Check if the boundary processor patches are correct
153  mesh.boundaryMesh().checkParallelSync(true);
154 
155  // Check names of zones are equal
156  mesh.cellZones().checkDefinition(true);
157  if (mesh.cellZones().checkParallelSync(true))
158  {
159  noFailedChecks++;
160  }
161  mesh.faceZones().checkDefinition(true);
162  if (mesh.faceZones().checkParallelSync(true))
163  {
164  noFailedChecks++;
165  }
166  mesh.pointZones().checkDefinition(true);
167  if (mesh.pointZones().checkParallelSync(true))
168  {
169  noFailedChecks++;
170  }
171 
172 
173  {
174  cellSet cells(mesh, "illegalCells", mesh.nCells()/100);
175 
176  forAll(mesh.cells(), celli)
177  {
178  const cell& cFaces = mesh.cells()[celli];
179 
180  if (cFaces.size() <= 3)
181  {
182  cells.insert(celli);
183  }
184  for (const label facei : cFaces)
185  {
186  if (facei < 0 || facei >= mesh.nFaces())
187  {
188  cells.insert(celli);
189  break;
190  }
191  }
192  }
193  const label nCells = returnReduce(cells.size(), sumOp<label>());
194 
195  if (nCells > 0)
196  {
197  Info<< " Illegal cells (less than 4 faces or out of range faces)"
198  << " found, number of cells: " << nCells << endl;
199  noFailedChecks++;
200 
201  Info<< " <<Writing " << nCells
202  << " illegal cells to set " << cells.name() << endl;
203  cells.instance() = mesh.pointsInstance();
204  cells.write();
205  if (surfWriter)
206  {
207  mergeAndWrite(*surfWriter, cells);
208  }
209  }
210  else
211  {
212  Info<< " Cell to face addressing OK." << endl;
213  }
214  }
215 
216 
217  {
218  pointSet points(mesh, "unusedPoints", mesh.nPoints()/100);
219  if (mesh.checkPoints(true, &points))
220  {
221  noFailedChecks++;
222 
223  label nPoints = returnReduce(points.size(), sumOp<label>());
224 
225  Info<< " <<Writing " << nPoints
226  << " unused points to set " << points.name() << endl;
227  points.instance() = mesh.pointsInstance();
228  points.write();
229  if (setWriter)
230  {
231  mergeAndWrite(*setWriter, points);
232  }
233  }
234  }
235 
236  {
237  faceSet faces(mesh, "upperTriangularFace", mesh.nFaces()/100);
238  if (mesh.checkUpperTriangular(true, &faces))
239  {
240  noFailedChecks++;
241  }
242 
243  const label nFaces = returnReduce(faces.size(), sumOp<label>());
244 
245  if (nFaces > 0)
246  {
247  Info<< " <<Writing " << nFaces
248  << " unordered faces to set " << faces.name() << endl;
249  faces.instance() = mesh.pointsInstance();
250  faces.write();
251  if (surfWriter)
252  {
253  mergeAndWrite(*surfWriter, faces);
254  }
255  }
256  }
257 
258  {
259  faceSet faces(mesh, "outOfRangeFaces", mesh.nFaces()/100);
260  if (mesh.checkFaceVertices(true, &faces))
261  {
262  noFailedChecks++;
263 
264  const label nFaces = returnReduce(faces.size(), sumOp<label>());
265 
266  Info<< " <<Writing " << nFaces
267  << " faces with out-of-range or duplicate vertices to set "
268  << faces.name() << endl;
269  faces.instance() = mesh.pointsInstance();
270  faces.write();
271  if (surfWriter)
272  {
273  mergeAndWrite(*surfWriter, faces);
274  }
275  }
276  }
277 
278  if (allTopology)
279  {
280  cellSet cells(mesh, "zipUpCells", mesh.nCells()/100);
281  if (mesh.checkCellsZipUp(true, &cells))
282  {
283  noFailedChecks++;
284 
285  const label nCells = returnReduce(cells.size(), sumOp<label>());
286 
287  Info<< " <<Writing " << nCells
288  << " cells with over used edges to set " << cells.name()
289  << endl;
290  cells.instance() = mesh.pointsInstance();
291  cells.write();
292  if (surfWriter)
293  {
294  mergeAndWrite(*surfWriter, cells);
295  }
296 
297  }
298  }
299 
300  if (allTopology)
301  {
302  faceSet faces(mesh, "edgeFaces", mesh.nFaces()/100);
303  if (mesh.checkFaceFaces(true, &faces))
304  {
305  noFailedChecks++;
306  }
307 
308  const label nFaces = returnReduce(faces.size(), sumOp<label>());
309  if (nFaces > 0)
310  {
311  Info<< " <<Writing " << nFaces
312  << " faces with non-standard edge connectivity to set "
313  << faces.name() << endl;
314  faces.instance() = mesh.pointsInstance();
315  faces.write();
316  if (surfWriter)
317  {
318  mergeAndWrite(*surfWriter, faces);
319  }
320  }
321  }
322 
323  if (allTopology)
324  {
325  labelList nInternalFaces(mesh.nCells(), Zero);
326 
327  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
328  {
329  nInternalFaces[mesh.faceOwner()[facei]]++;
330  nInternalFaces[mesh.faceNeighbour()[facei]]++;
331  }
332  const polyBoundaryMesh& patches = mesh.boundaryMesh();
333  forAll(patches, patchi)
334  {
335  if (patches[patchi].coupled())
336  {
337  const labelUList& owners = patches[patchi].faceCells();
338 
339  for (const label facei : owners)
340  {
341  nInternalFaces[facei]++;
342  }
343  }
344  }
345 
346  cellSet oneCells(mesh, "oneInternalFaceCells", mesh.nCells()/100);
347  cellSet twoCells(mesh, "twoInternalFacesCells", mesh.nCells()/100);
348 
349  forAll(nInternalFaces, celli)
350  {
351  if (nInternalFaces[celli] <= 1)
352  {
353  oneCells.insert(celli);
354  }
355  else if (nInternalFaces[celli] == 2)
356  {
357  twoCells.insert(celli);
358  }
359  }
360 
361  label nOneCells = returnReduce(oneCells.size(), sumOp<label>());
362 
363  if (nOneCells > 0)
364  {
365  Info<< " <<Writing " << nOneCells
366  << " cells with zero or one non-boundary face to set "
367  << oneCells.name()
368  << endl;
369  oneCells.instance() = mesh.pointsInstance();
370  oneCells.write();
371  if (surfWriter)
372  {
373  mergeAndWrite(*surfWriter, oneCells);
374  }
375  }
376 
377  label nTwoCells = returnReduce(twoCells.size(), sumOp<label>());
378 
379  if (nTwoCells > 0)
380  {
381  Info<< " <<Writing " << nTwoCells
382  << " cells with two non-boundary faces to set "
383  << twoCells.name()
384  << endl;
385  twoCells.instance() = mesh.pointsInstance();
386  twoCells.write();
387  if (surfWriter)
388  {
389  mergeAndWrite(*surfWriter, twoCells);
390  }
391  }
392  }
393 
394  {
395  regionSplit rs(mesh);
396 
397  if (rs.nRegions() <= 1)
398  {
399  Info<< " Number of regions: " << rs.nRegions() << " (OK)."
400  << endl;
401 
402  }
403  else
404  {
405  Info<< " *Number of regions: " << rs.nRegions() << endl;
406 
407  Info<< " The mesh has multiple regions which are not connected "
408  "by any face." << endl
409  << " <<Writing region information to "
410  << mesh.time().timeName()/"cellToRegion"
411  << endl;
412 
413  labelIOList ctr
414  (
415  IOobject
416  (
417  "cellToRegion",
418  mesh.time().timeName(),
419  mesh,
420  IOobject::NO_READ,
421  IOobject::NO_WRITE
422  ),
423  rs
424  );
425  ctr.write();
426 
427 
428  // Points in multiple regions
429  pointSet points
430  (
431  mesh,
432  "multiRegionPoints",
433  mesh.nPoints()/1000
434  );
435 
436  // Is region disconnected
437  boolList regionDisconnected(rs.nRegions(), true);
438  if (allTopology)
439  {
440  // -1 : not assigned
441  // -2 : multiple regions
442  // >= 0 : single region
443  labelList pointToRegion(mesh.nPoints(), -1);
444 
445  for
446  (
447  label facei = mesh.nInternalFaces();
448  facei < mesh.nFaces();
449  ++facei
450  )
451  {
452  const label regioni = rs[mesh.faceOwner()[facei]];
453  const face& f = mesh.faces()[facei];
454  for (const label verti : f)
455  {
456  label& pRegion = pointToRegion[verti];
457  if (pRegion == -1)
458  {
459  pRegion = regioni;
460  }
461  else if (pRegion == -2)
462  {
463  // Already marked
464  regionDisconnected[regioni] = false;
465  }
466  else if (pRegion != regioni)
467  {
468  // Multiple regions
469  regionDisconnected[regioni] = false;
470  regionDisconnected[pRegion] = false;
471  pRegion = -2;
472  points.insert(verti);
473  }
474  }
475  }
476 
477  Pstream::listCombineGather(regionDisconnected, andEqOp<bool>());
478  Pstream::listCombineScatter(regionDisconnected);
479  }
480 
481 
482 
483  // write cellSet for each region
484  PtrList<cellSet> cellRegions(rs.nRegions());
485  for (label i = 0; i < rs.nRegions(); i++)
486  {
487  cellRegions.set
488  (
489  i,
490  new cellSet
491  (
492  mesh,
493  "region" + Foam::name(i),
494  mesh.nCells()/100
495  )
496  );
497  }
498 
499  forAll(rs, i)
500  {
501  cellRegions[rs[i]].insert(i);
502  }
503 
504  for (label i = 0; i < rs.nRegions(); i++)
505  {
506  Info<< " <<Writing region " << i;
507  if (allTopology)
508  {
509  if (regionDisconnected[i])
510  {
511  Info<< " (fully disconnected)";
512  }
513  else
514  {
515  Info<< " (point connected)";
516  }
517  }
518  Info<< " with "
519  << returnReduce(cellRegions[i].size(), sumOp<scalar>())
520  << " cells to cellSet " << cellRegions[i].name() << endl;
521 
522  cellRegions[i].write();
523  }
524 
525  label nPoints = returnReduce(points.size(), sumOp<label>());
526  if (nPoints)
527  {
528  Info<< " <<Writing " << nPoints
529  << " points that are in multiple regions to set "
530  << points.name() << endl;
531  points.write();
532  if (setWriter)
533  {
534  mergeAndWrite(*setWriter, points);
535  }
536  }
537  }
538  }
539 
540  // Non-manifold points
541  pointSet points
542  (
543  mesh,
544  "nonManifoldPoints",
545  mesh.nPoints()/1000
546  );
547 
548  {
549  if (!Pstream::parRun())
550  {
551  Info<< "\nChecking patch topology for multiply connected"
552  << " surfaces..." << endl;
553  }
554  else
555  {
556  Info<< "\nChecking basic patch addressing..." << endl;
557  }
558 
559  const polyBoundaryMesh& patches = mesh.boundaryMesh();
560 
561  Pout.setf(ios_base::left);
562 
563  Info<< " "
564  << setw(20) << "Patch"
565  << setw(9) << "Faces"
566  << setw(9) << "Points";
567  if (!Pstream::parRun())
568  {
569  Info<< setw(34) << "Surface topology";
570  }
571  if (allGeometry)
572  {
573  Info<< " Bounding box";
574  }
575  Info<< endl;
576 
577  forAll(patches, patchi)
578  {
579  const polyPatch& pp = patches[patchi];
580 
581  if (!isA<processorPolyPatch>(pp))
582  {
583  checkPatch(allGeometry, pp.name(), pp, points);
584  Info<< endl;
585  }
586  }
587 
588  //Info.setf(ios_base::right);
589  }
590 
591  {
592  if (!Pstream::parRun())
593  {
594  Info<< "\nChecking faceZone topology for multiply connected"
595  << " surfaces..." << endl;
596  }
597  else
598  {
599  Info<< "\nChecking basic faceZone addressing..." << endl;
600  }
601 
602  Pout.setf(ios_base::left);
603 
604  const faceZoneMesh& faceZones = mesh.faceZones();
605 
606  if (faceZones.size())
607  {
608  Info<< " "
609  << setw(20) << "FaceZone"
610  << setw(9) << "Faces"
611  << setw(9) << "Points";
612 
613  if (!Pstream::parRun())
614  {
615  Info<< setw(34) << "Surface topology";
616  }
617  if (allGeometry)
618  {
619  Info<< " Bounding box";
620  }
621  Info<< endl;
622 
623  for (const faceZone& fz : faceZones)
624  {
625  checkPatch(allGeometry, fz.name(), fz(), points);
626  Info<< endl;
627  }
628  }
629  else
630  {
631  Info<< " No faceZones found."<<endl;
632  }
633  }
634 
635  const label nPoints = returnReduce(points.size(), sumOp<label>());
636 
637  if (nPoints)
638  {
639  Info<< " <<Writing " << nPoints
640  << " conflicting points to set " << points.name() << endl;
641  points.instance() = mesh.pointsInstance();
642  points.write();
643  if (setWriter)
644  {
645  mergeAndWrite(*setWriter, points);
646  }
647  }
648 
649  {
650  Info<< "\nChecking basic cellZone addressing..." << endl;
651 
652  Pout.setf(ios_base::left);
653 
654  const cellZoneMesh& cellZones = mesh.cellZones();
655 
656  if (cellZones.size())
657  {
658  Info<< " "
659  << setw(20) << "CellZone"
660  << setw(13) << "Cells"
661  << setw(13) << "Points"
662  << setw(13) << "Volume"
663  << "BoundingBox" << endl;
664 
665  const cellList& cells = mesh.cells();
666  const faceList& faces = mesh.faces();
667  const scalarField& cellVolumes = mesh.cellVolumes();
668 
669  bitSet isZonePoint(mesh.nPoints());
670 
671  for (const cellZone& cZone : cellZones)
672  {
673  boundBox bb;
674  isZonePoint.reset(); // clears all bits (reset count)
675  scalar v = 0.0;
676 
677  for (const label celli : cZone)
678  {
679  v += cellVolumes[celli];
680  for (const label facei : cells[celli])
681  {
682  const face& f = faces[facei];
683  for (const label verti : f)
684  {
685  if (isZonePoint.set(verti))
686  {
687  bb.add(mesh.points()[verti]);
688  }
689  }
690  }
691  }
692 
693  bb.reduce(); // Global min/max
694 
695  Info<< " "
696  << setw(19) << cZone.name()
697  << ' ' << setw(12)
698  << returnReduce(cZone.size(), sumOp<label>())
699  << ' ' << setw(12)
700  << returnReduce(isZonePoint.count(), sumOp<label>())
701  << ' ' << setw(12)
702  << returnReduce(v, sumOp<scalar>())
703  << ' ' << bb << endl;
704  }
705  }
706  else
707  {
708  Info<< " No cellZones found."<<endl;
709  }
710  }
711 
712  // Force creation of all addressing if requested.
713  // Errors will be reported as required
714  if (allTopology)
715  {
716  mesh.cells();
717  mesh.faces();
718  mesh.edges();
719  mesh.points();
720  mesh.faceOwner();
721  mesh.faceNeighbour();
722  mesh.cellCells();
723  mesh.edgeCells();
724  mesh.pointCells();
725  mesh.edgeFaces();
726  mesh.pointFaces();
727  mesh.cellEdges();
728  mesh.faceEdges();
729  mesh.pointEdges();
730  }
731 
732  return noFailedChecks;
733 }
734 
735 
736 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:63
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
Foam::constant::atomic::mp
const dimensionedScalar mp
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:88
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:61
checkTopology.H
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::checkTopology
label checkTopology(const polyMesh &mesh, const bool allTopology, const bool allGeometry, autoPtr< surfaceWriter > &surfWriter, const autoPtr< writer< scalar >> &setWriter)
Foam::Pout
prefixOSstream Pout
Foam::cellZoneMesh
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
Definition: cellZoneMeshFwd.H:38
polyMesh.H
Foam::labelIOList
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:38
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::checkPatch
void checkPatch(const bool allGeometry, const word &name, const PatchType &pp, pointSet &points)
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::faceZoneMesh
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
Definition: faceZoneMeshFwd.H:38
treeBoundBox.H
Foam::Info
messageStream Info
faceSet.H
regionSplit.H
Foam::cellList
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:41
IOmanip.H
Istream and Ostream manipulators taking arguments.
processorPolyPatch.H
reduce
reduce(hasMovingMesh, orOp< bool >())
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::Ostream::write
virtual bool write(const token &tok)=0
Foam::setw
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
emptyPolyPatch.H
Foam::IOstream::setf
ios_base::fmtflags setf(const ios_base::fmtflags f)
Definition: IOstream.H:374
Time.H
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:41
vtkSurfaceWriter.H
points
const pointField & points
Definition: gmvOutputHeader.H:1
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:59
Foam::mergeAndWrite
void mergeAndWrite(const polyMesh &mesh, surfaceWriter &writer, const word &name, const indirectPrimitivePatch &setPatch, const fileName &outputDir)
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
cellSet.H
coupled
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:81
checkTools.H
pointSet.H