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-2022 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 "vtkCoordSetWriter.H"
40 #include "vtkSurfaceWriter.H"
41 #include "checkTools.H"
42 #include "treeBoundBox.H"
43 #include "syncTools.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 template<class PatchType>
49 (
50  const bool allGeometry,
51  const word& name,
52  const polyMesh& mesh,
53  const PatchType& pp,
54  const labelList& meshFaces,
55  const labelList& meshEdges,
56  pointSet& points
57 )
58 {
59  typedef typename PatchType::surfaceTopo TopoType;
60 
61  const label globalSize = returnReduce(pp.size(), sumOp<label>());
62 
63  Info<< " "
64  << setw(20) << name
65  << setw(9) << globalSize
66  << setw(9) << returnReduce(pp.nPoints(), sumOp<label>());
67 
68  if (globalSize == 0)
69  {
70  Info<< setw(34) << "ok (empty)";
71  }
72  else if (Pstream::parRun())
73  {
74  // Parallel - use mesh edges
75  // - no check for point-pinch
76  // - no check for consistent orientation (if that is posible to
77  // check?)
78 
79  // (see addPatchCellLayer::globalEdgeFaces)
80  // From mesh edge to global face labels. Non-empty sublists only for
81  // pp edges.
82  labelListList globalEdgeFaces(mesh.nEdges());
83 
84  const labelListList& edgeFaces = pp.edgeFaces();
85 
86  // Global numbering
87  const globalIndex globalFaces(mesh.nFaces());
88 
89  forAll(edgeFaces, edgei)
90  {
91  label meshEdgei = meshEdges[edgei];
92  const labelList& eFaces = edgeFaces[edgei];
93 
94  // Store face and processor as unique tag.
95  labelList& globalEFaces = globalEdgeFaces[meshEdgei];
96  globalEFaces.setSize(eFaces.size());
97  forAll(eFaces, i)
98  {
99  globalEFaces[i] = globalFaces.toGlobal(meshFaces[eFaces[i]]);
100  }
101  //Pout<< "At edge:" << meshEdgei
102  // << " ctr:" << mesh.edges()[meshEdgei].centre(mesh.points())
103  // << " have eFaces:" << globalEdgeFaces[meshEdgei]
104  // << endl;
105  }
106 
107  // Synchronise across coupled edges.
108  syncTools::syncEdgeList
109  (
110  mesh,
111  globalEdgeFaces,
112  ListOps::uniqueEqOp<label>(),
113  labelList() // null value
114  );
115 
116 
117  label labelTyp = TopoType::MANIFOLD;
118  forAll(meshEdges, edgei)
119  {
120  const label meshEdgei = meshEdges[edgei];
121  const labelList& globalEFaces = globalEdgeFaces[meshEdgei];
122  if (globalEFaces.size() == 1)
123  {
124  //points.insert(mesh.edges()[meshEdgei]);
125  labelTyp = max(labelTyp, TopoType::OPEN);
126  }
127  else if (globalEFaces.size() == 0 || globalEFaces.size() > 2)
128  {
129  points.insert(mesh.edges()[meshEdgei]);
130  labelTyp = max(labelTyp, TopoType::ILLEGAL);
131  }
132  }
133  reduce(labelTyp, maxOp<label>());
134 
135  if (labelTyp == TopoType::MANIFOLD)
136  {
137  Info<< setw(34) << "ok (closed singly connected)";
138  }
139  else if (labelTyp == TopoType::OPEN)
140  {
141  Info<< setw(34)
142  << "ok (non-closed singly connected)";
143  }
144  else
145  {
146  Info<< setw(34)
147  << "multiply connected (shared edge)";
148  }
149  }
150  else
151  {
152  TopoType pTyp = pp.surfaceType();
153 
154  if (pTyp == TopoType::MANIFOLD)
155  {
156  if (pp.checkPointManifold(true, &points))
157  {
158  Info<< setw(34)
159  << "multiply connected (shared point)";
160  }
161  else
162  {
163  Info<< setw(34) << "ok (closed singly connected)";
164  }
165 
166  // Add points on non-manifold edges to make set complete
167  pp.checkTopology(false, &points);
168  }
169  else
170  {
171  pp.checkTopology(false, &points);
172 
173  if (pTyp == TopoType::OPEN)
174  {
175  Info<< setw(34)
176  << "ok (non-closed singly connected)";
177  }
178  else
179  {
180  Info<< setw(34)
181  << "multiply connected (shared edge)";
182  }
183  }
184  }
185 
186  if (allGeometry)
187  {
188  const labelList& mp = pp.meshPoints();
189 
190  if (returnReduceOr(mp.size()))
191  {
192  boundBox bb(pp.points(), mp, true); // reduce
193  Info<< ' ' << bb;
194  }
195  }
196 }
197 
198 
199 template<class Zone>
200 Foam::label Foam::checkZones
201 (
202  const polyMesh& mesh,
203  const ZoneMesh<Zone, polyMesh>& zones,
204  topoSet& set
205 )
206 {
207  labelList zoneID(set.maxSize(mesh), -1);
208  for (const auto& zone : zones)
209  {
210  for (const label elem : zone)
211  {
212  if
213  (
214  zoneID[elem] != -1
215  && zoneID[elem] != zone.index()
216  )
217  {
218  set.insert(elem);
219  }
220  zoneID[elem] = zone.index();
221  }
222  }
223 
224  return returnReduce(set.size(), sumOp<label>());
225 }
226 
227 
228 Foam::label Foam::checkTopology
229 (
230  const polyMesh& mesh,
231  const bool allTopology,
232  const bool allGeometry,
233  autoPtr<surfaceWriter>& surfWriter,
234  autoPtr<coordSetWriter>& setWriter
235 )
236 {
237  label noFailedChecks = 0;
238 
239  Info<< "Checking topology..." << endl;
240 
241  // Check if the boundary definition is unique
242  mesh.boundaryMesh().checkDefinition(true);
243 
244  // Check that empty patches cover all sides of the mesh
245  {
246  label nEmpty = 0;
247  forAll(mesh.boundaryMesh(), patchi)
248  {
249  if (isA<emptyPolyPatch>(mesh.boundaryMesh()[patchi]))
250  {
251  nEmpty += mesh.boundaryMesh()[patchi].size();
252  }
253  }
254  reduce(nEmpty, sumOp<label>());
255  const label nCells = returnReduce(mesh.cells().size(), sumOp<label>());
256 
257  // These are actually warnings, not errors.
258  if (nCells && (nEmpty % nCells))
259  {
260  Info<< " ***Total number of faces on empty patches"
261  << " is not divisible by the number of cells in the mesh."
262  << " Hence this mesh is not 1D or 2D."
263  << endl;
264  }
265  }
266 
267  // Check if the boundary processor patches are correct
268  mesh.boundaryMesh().checkParallelSync(true);
269 
270  // Check names of zones are equal
271  mesh.cellZones().checkDefinition(true);
272  if (mesh.cellZones().checkParallelSync(true))
273  {
274  noFailedChecks++;
275  }
276  mesh.faceZones().checkDefinition(true);
277  if (mesh.faceZones().checkParallelSync(true))
278  {
279  noFailedChecks++;
280  }
281  mesh.pointZones().checkDefinition(true);
282  if (mesh.pointZones().checkParallelSync(true))
283  {
284  noFailedChecks++;
285  }
286 
287 
288  {
289  cellSet cells(mesh, "illegalCells", mesh.nCells()/100);
290 
291  forAll(mesh.cells(), celli)
292  {
293  const cell& cFaces = mesh.cells()[celli];
294 
295  if (cFaces.size() <= 3)
296  {
297  cells.insert(celli);
298  }
299  for (const label facei : cFaces)
300  {
301  if (facei < 0 || facei >= mesh.nFaces())
302  {
303  cells.insert(celli);
304  break;
305  }
306  }
307  }
308  const label nCells = returnReduce(cells.size(), sumOp<label>());
309 
310  if (nCells > 0)
311  {
312  Info<< " Illegal cells (less than 4 faces or out of range faces)"
313  << " found, number of cells: " << nCells << endl;
314  noFailedChecks++;
315 
316  Info<< " <<Writing " << nCells
317  << " illegal cells to set " << cells.name() << endl;
318  cells.instance() = mesh.pointsInstance();
319  cells.write();
320  if (surfWriter && surfWriter->enabled())
321  {
322  mergeAndWrite(*surfWriter, cells);
323  }
324  }
325  else
326  {
327  Info<< " Cell to face addressing OK." << endl;
328  }
329  }
330 
331 
332  {
333  pointSet points(mesh, "unusedPoints", mesh.nPoints()/100);
334  if (mesh.checkPoints(true, &points))
335  {
336  noFailedChecks++;
337 
338  const label nPoints = returnReduce(points.size(), sumOp<label>());
339 
340  Info<< " <<Writing " << nPoints
341  << " unused points to set " << points.name() << endl;
342  points.instance() = mesh.pointsInstance();
343  points.write();
344  if (setWriter && setWriter->enabled())
345  {
346  mergeAndWrite(*setWriter, points);
347  }
348  }
349  }
350 
351  {
352  faceSet faces(mesh, "upperTriangularFace", mesh.nFaces()/100);
353  if (mesh.checkUpperTriangular(true, &faces))
354  {
355  noFailedChecks++;
356  }
357 
358  const label nFaces = returnReduce(faces.size(), sumOp<label>());
359 
360  if (nFaces > 0)
361  {
362  Info<< " <<Writing " << nFaces
363  << " unordered faces to set " << faces.name() << endl;
364  faces.instance() = mesh.pointsInstance();
365  faces.write();
366  if (surfWriter && surfWriter->enabled())
367  {
368  mergeAndWrite(*surfWriter, faces);
369  }
370  }
371  }
372 
373  {
374  faceSet faces(mesh, "outOfRangeFaces", mesh.nFaces()/100);
375  if (mesh.checkFaceVertices(true, &faces))
376  {
377  noFailedChecks++;
378 
379  const label nFaces = returnReduce(faces.size(), sumOp<label>());
380 
381  Info<< " <<Writing " << nFaces
382  << " faces with out-of-range or duplicate vertices to set "
383  << faces.name() << endl;
384  faces.instance() = mesh.pointsInstance();
385  faces.write();
386  if (surfWriter && surfWriter->enabled())
387  {
388  mergeAndWrite(*surfWriter, faces);
389  }
390  }
391  }
392 
393  if (allTopology)
394  {
395  cellSet cells(mesh, "zipUpCells", mesh.nCells()/100);
396  if (mesh.checkCellsZipUp(true, &cells))
397  {
398  noFailedChecks++;
399 
400  const label nCells = returnReduce(cells.size(), sumOp<label>());
401 
402  Info<< " <<Writing " << nCells
403  << " cells with over used edges to set " << cells.name()
404  << endl;
405  cells.instance() = mesh.pointsInstance();
406  cells.write();
407  if (surfWriter && surfWriter->enabled())
408  {
409  mergeAndWrite(*surfWriter, cells);
410  }
411  }
412  }
413 
414  if (allTopology)
415  {
416  faceSet faces(mesh, "edgeFaces", mesh.nFaces()/100);
417  if (mesh.checkFaceFaces(true, &faces))
418  {
419  noFailedChecks++;
420  }
421 
422  const label nFaces = returnReduce(faces.size(), sumOp<label>());
423  if (nFaces > 0)
424  {
425  Info<< " <<Writing " << nFaces
426  << " faces with non-standard edge connectivity to set "
427  << faces.name() << endl;
428  faces.instance() = mesh.pointsInstance();
429  faces.write();
430  if (surfWriter && surfWriter->enabled())
431  {
432  mergeAndWrite(*surfWriter, faces);
433  }
434  }
435  }
436 
437  if (allTopology)
438  {
439  labelList nInternalFaces(mesh.nCells(), Zero);
440 
441  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
442  {
443  nInternalFaces[mesh.faceOwner()[facei]]++;
444  nInternalFaces[mesh.faceNeighbour()[facei]]++;
445  }
446  const polyBoundaryMesh& patches = mesh.boundaryMesh();
447  forAll(patches, patchi)
448  {
449  if (patches[patchi].coupled())
450  {
451  const labelUList& owners = patches[patchi].faceCells();
452 
453  for (const label facei : owners)
454  {
455  nInternalFaces[facei]++;
456  }
457  }
458  }
459 
460  cellSet oneCells(mesh, "oneInternalFaceCells", mesh.nCells()/100);
461  cellSet twoCells(mesh, "twoInternalFacesCells", mesh.nCells()/100);
462 
463  forAll(nInternalFaces, celli)
464  {
465  if (nInternalFaces[celli] <= 1)
466  {
467  oneCells.insert(celli);
468  }
469  else if (nInternalFaces[celli] == 2)
470  {
471  twoCells.insert(celli);
472  }
473  }
474 
475  const label nOneCells = returnReduce(oneCells.size(), sumOp<label>());
476 
477  if (nOneCells > 0)
478  {
479  Info<< " <<Writing " << nOneCells
480  << " cells with zero or one non-boundary face to set "
481  << oneCells.name()
482  << endl;
483  oneCells.instance() = mesh.pointsInstance();
484  oneCells.write();
485  if (surfWriter && surfWriter->enabled())
486  {
487  mergeAndWrite(*surfWriter, oneCells);
488  }
489  }
490 
491  const label nTwoCells = returnReduce(twoCells.size(), sumOp<label>());
492 
493  if (nTwoCells > 0)
494  {
495  Info<< " <<Writing " << nTwoCells
496  << " cells with two non-boundary faces to set "
497  << twoCells.name()
498  << endl;
499  twoCells.instance() = mesh.pointsInstance();
500  twoCells.write();
501  if (surfWriter && surfWriter->enabled())
502  {
503  mergeAndWrite(*surfWriter, twoCells);
504  }
505  }
506  }
507 
508  {
509  regionSplit rs(mesh);
510 
511  if (rs.nRegions() <= 1)
512  {
513  Info<< " Number of regions: " << rs.nRegions() << " (OK)."
514  << endl;
515 
516  }
517  else
518  {
519  Info<< " *Number of regions: " << rs.nRegions() << endl;
520 
521  Info<< " The mesh has multiple regions which are not connected "
522  "by any face." << endl
523  << " <<Writing region information to "
524  << mesh.time().timeName()/"cellToRegion"
525  << endl;
526 
527  labelIOList ctr
528  (
529  IOobject
530  (
531  "cellToRegion",
532  mesh.time().timeName(),
533  mesh,
534  IOobject::NO_READ,
535  IOobject::NO_WRITE
536  ),
537  rs
538  );
539  ctr.write();
540 
541 
542  // Points in multiple regions
543  pointSet points
544  (
545  mesh,
546  "multiRegionPoints",
547  mesh.nPoints()/1000
548  );
549 
550  // Is region disconnected
551  boolList regionDisconnected(rs.nRegions(), true);
552  if (allTopology)
553  {
554  // -1 : not assigned
555  // -2 : multiple regions
556  // >= 0 : single region
557  labelList pointToRegion(mesh.nPoints(), -1);
558 
559  for
560  (
561  label facei = mesh.nInternalFaces();
562  facei < mesh.nFaces();
563  ++facei
564  )
565  {
566  const label regioni = rs[mesh.faceOwner()[facei]];
567  const face& f = mesh.faces()[facei];
568  for (const label verti : f)
569  {
570  label& pRegion = pointToRegion[verti];
571  if (pRegion == -1)
572  {
573  pRegion = regioni;
574  }
575  else if (pRegion == -2)
576  {
577  // Already marked
578  regionDisconnected[regioni] = false;
579  }
580  else if (pRegion != regioni)
581  {
582  // Multiple regions
583  regionDisconnected[regioni] = false;
584  regionDisconnected[pRegion] = false;
585  pRegion = -2;
586  points.insert(verti);
587  }
588  }
589  }
590 
591  Pstream::listCombineReduce(regionDisconnected, andEqOp<bool>());
592  }
593 
594 
595 
596  // write cellSet for each region
597  PtrList<cellSet> cellRegions(rs.nRegions());
598  for (label i = 0; i < rs.nRegions(); i++)
599  {
600  cellRegions.set
601  (
602  i,
603  new cellSet
604  (
605  mesh,
606  "region" + Foam::name(i),
607  mesh.nCells()/100
608  )
609  );
610  }
611 
612  forAll(rs, i)
613  {
614  cellRegions[rs[i]].insert(i);
615  }
616 
617  for (label i = 0; i < rs.nRegions(); i++)
618  {
619  Info<< " <<Writing region " << i;
620  if (allTopology)
621  {
622  if (regionDisconnected[i])
623  {
624  Info<< " (fully disconnected)";
625  }
626  else
627  {
628  Info<< " (point connected)";
629  }
630  }
631  Info<< " with "
632  << returnReduce(cellRegions[i].size(), sumOp<scalar>())
633  << " cells to cellSet " << cellRegions[i].name() << endl;
634 
635  cellRegions[i].write();
636  }
637 
638  const label nPoints = returnReduce(points.size(), sumOp<label>());
639  if (nPoints)
640  {
641  Info<< " <<Writing " << nPoints
642  << " points that are in multiple regions to set "
643  << points.name() << endl;
644  points.write();
645  if (setWriter && setWriter->enabled())
646  {
647  mergeAndWrite(*setWriter, points);
648  }
649  }
650  }
651  }
652 
653  // Non-manifold points
654  pointSet points
655  (
656  mesh,
657  "nonManifoldPoints",
658  mesh.nPoints()/1000
659  );
660 
661  {
662  Info<< "\nChecking patch topology for multiply connected"
663  << " surfaces..." << endl;
664 
665  const polyBoundaryMesh& patches = mesh.boundaryMesh();
666 
667  Pout.setf(ios_base::left);
668 
669  Info<< " "
670  << setw(20) << "Patch"
671  << setw(9) << "Faces"
672  << setw(9) << "Points"
673  << "Surface topology";
674  if (allGeometry)
675  {
676  Info<< " Bounding box";
677  }
678  Info<< endl;
679 
680  forAll(patches, patchi)
681  {
682  const polyPatch& pp = patches[patchi];
683 
684  if (!isA<processorPolyPatch>(pp))
685  {
686  checkPatch
687  (
688  allGeometry,
689  pp.name(),
690  mesh,
691  pp,
692  identity(pp.size(), pp.start()),
693  pp.meshEdges(),
694  points
695  );
696  Info<< endl;
697  }
698  }
699 
700  //Info.setf(ios_base::right);
701  }
702 
703  {
704  Info<< "\nChecking faceZone topology for multiply connected"
705  << " surfaces..." << endl;
706 
707  Pout.setf(ios_base::left);
708 
709  const faceZoneMesh& faceZones = mesh.faceZones();
710 
711  if (faceZones.size())
712  {
713  Info<< " "
714  << setw(20) << "FaceZone"
715  << setw(9) << "Faces"
716  << setw(9) << "Points"
717  << setw(34) << "Surface topology";
718  if (allGeometry)
719  {
720  Info<< " Bounding box";
721  }
722  Info<< endl;
723 
724  for (const faceZone& fz : faceZones)
725  {
726  checkPatch
727  (
728  allGeometry,
729  fz.name(),
730  mesh,
731  fz(), // patch
732  fz, // mesh face labels
733  fz.meshEdges(), // mesh edge labels
734  points
735  );
736  Info<< endl;
737  }
738 
739  // Check for duplicates
740  if (allTopology)
741  {
742  faceSet mzFaces(mesh, "multiZoneFaces", mesh.nFaces()/100);
743  const label nMulti = checkZones(mesh, faceZones, mzFaces);
744  if (nMulti)
745  {
746  Info<< " <<Writing " << nMulti
747  << " faces that are in multiple zones"
748  << " to set " << mzFaces.name() << endl;
749  mzFaces.instance() = mesh.pointsInstance();
750  mzFaces.write();
751  if (surfWriter && surfWriter->enabled())
752  {
753  mergeAndWrite(*surfWriter, mzFaces);
754  }
755  }
756  }
757  }
758  else
759  {
760  Info<< " No faceZones found."<<endl;
761  }
762  }
763 
764  const label nPoints = returnReduce(points.size(), sumOp<label>());
765 
766  if (nPoints)
767  {
768  Info<< " <<Writing " << nPoints
769  << " conflicting points to set " << points.name() << endl;
770  points.instance() = mesh.pointsInstance();
771  points.write();
772  if (setWriter && setWriter->enabled())
773  {
774  mergeAndWrite(*setWriter, points);
775  }
776  }
777 
778  {
779  Info<< "\nChecking basic cellZone addressing..." << endl;
780 
781  Pout.setf(ios_base::left);
782 
783  const cellZoneMesh& cellZones = mesh.cellZones();
784 
785  if (cellZones.size())
786  {
787  Info<< " "
788  << setw(20) << "CellZone"
789  << setw(13) << "Cells"
790  << setw(13) << "Points"
791  << setw(13) << "Volume"
792  << "BoundingBox" << endl;
793 
794  const cellList& cells = mesh.cells();
795  const faceList& faces = mesh.faces();
796  const scalarField& cellVolumes = mesh.cellVolumes();
797 
798  bitSet isZonePoint(mesh.nPoints());
799 
800  for (const cellZone& cZone : cellZones)
801  {
802  boundBox bb;
803  isZonePoint.reset(); // clears all bits (reset count)
804  scalar v = 0.0;
805 
806  for (const label celli : cZone)
807  {
808  v += cellVolumes[celli];
809  for (const label facei : cells[celli])
810  {
811  const face& f = faces[facei];
812  for (const label verti : f)
813  {
814  if (isZonePoint.set(verti))
815  {
816  bb.add(mesh.points()[verti]);
817  }
818  }
819  }
820  }
821 
822  bb.reduce(); // Global min/max
823 
824  Info<< " "
825  << setw(19) << cZone.name()
826  << ' ' << setw(12)
827  << returnReduce(cZone.size(), sumOp<label>())
828  << ' ' << setw(12)
829  << returnReduce(isZonePoint.count(), sumOp<label>())
830  << ' ' << setw(12)
831  << returnReduce(v, sumOp<scalar>())
832  << ' ' << bb << endl;
833  }
834 
835 
836  // Check for duplicates
837  if (allTopology)
838  {
839  cellSet mzCells(mesh, "multiZoneCells", mesh.nCells()/100);
840  const label nMulti = checkZones(mesh, cellZones, mzCells);
841  if (nMulti)
842  {
843  Info<< " <<Writing " << nMulti
844  << " cells that are in multiple zones"
845  << " to set " << mzCells.name() << endl;
846  mzCells.instance() = mesh.pointsInstance();
847  mzCells.write();
848  if (surfWriter && surfWriter->enabled())
849  {
850  mergeAndWrite(*surfWriter, mzCells);
851  }
852  }
853  }
854  }
855  else
856  {
857  Info<< " No cellZones found."<<endl;
858  }
859  }
860 
861 
862  {
863  Info<< "\nChecking basic pointZone addressing..." << endl;
864 
865  Pout.setf(ios_base::left);
866 
867  const pointZoneMesh& pointZones = mesh.pointZones();
868 
869  if (pointZones.size())
870  {
871  Info<< " "
872  << setw(20) << "PointZone"
873  << setw(8) << "Points"
874  << "BoundingBox" << nl;
875 
876  for (const auto& zone : pointZones)
877  {
878  boundBox bb
879  (
880  mesh.points(),
881  static_cast<const labelUList&>(zone),
882  true // Reduce (global min/max)
883  );
884 
885  Info<< " "
886  << setw(20) << zone.name()
887  << setw(8)
888  << returnReduce(zone.size(), sumOp<label>())
889  << bb << endl;
890  }
891 
892 
893  // Check for duplicates
894  if (allTopology)
895  {
896  pointSet mzPoints(mesh, "multiZonePoints", mesh.nPoints()/100);
897  const label nMulti = checkZones(mesh, pointZones, mzPoints);
898  if (nMulti)
899  {
900  Info<< " <<Writing " << nMulti
901  << " points that are in multiple zones"
902  << " to set " << mzPoints.name() << endl;
903  mzPoints.instance() = mesh.pointsInstance();
904  mzPoints.write();
905  if (setWriter && setWriter->enabled())
906  {
907  mergeAndWrite(*setWriter, mzPoints);
908  }
909  }
910  }
911  }
912  else
913  {
914  Info<< " No pointZones found."<<endl;
915  }
916  }
917 
918 
919  // Force creation of all addressing if requested.
920  // Errors will be reported as required
921  if (allTopology)
922  {
923  mesh.cells();
924  mesh.faces();
925  mesh.edges();
926  mesh.points();
927  mesh.faceOwner();
928  mesh.faceNeighbour();
929  mesh.cellCells();
930  mesh.edgeCells();
931  mesh.pointCells();
932  mesh.edgeFaces();
933  mesh.pointFaces();
934  mesh.cellEdges();
935  mesh.faceEdges();
936  mesh.pointEdges();
937  }
938 
939  return noFailedChecks;
940 }
941 
942 
943 // ************************************************************************* //
label checkTopology(const polyMesh &mesh, const bool allTopology, const bool allGeometry, autoPtr< surfaceWriter > &surfWriter, autoPtr< coordSetWriter > &setWriter)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:51
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
ios_base::fmtflags setf(const ios_base::fmtflags f)
Set flags of stream.
Definition: IOstream.H:465
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
void checkPatch(const bool allGeometry, const word &name, const polyMesh &mesh, const PatchType &pp, const labelList &meshFaces, const labelList &meshEdges, pointSet &points)
List< face > faceList
A List of faces.
Definition: faceListFwd.H:41
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
void mergeAndWrite(const polyMesh &mesh, surfaceWriter &writer, const word &name, const indirectPrimitivePatch &setPatch, const fileName &outputDir)
Generate merged surface on master and write. Needs input patch.
label checkZones(const polyMesh &mesh, const ZoneMesh< Zone, polyMesh > &zones, topoSet &set)
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
UList< label > labelUList
A UList of labels.
Definition: UList.H:80
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
void setSize(const label n)
Alias for resize()
Definition: List.H:289
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
const cellShapeList & cells
const pointField & points
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:31
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
Istream and Ostream manipulators taking arguments.
labelList f(nPoints)
const polyBoundaryMesh & patches
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
List< label > labelList
A List of labels.
Definition: List.H:62
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
List< bool > boolList
A List of bools.
Definition: List.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:41
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
const dimensionedScalar mp
Proton mass.
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:38
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157