Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "foamVtkLineWriter.H"
40#include "vtkCoordSetWriter.H"
41#include "vtkSurfaceWriter.H"
42#include "checkTools.H"
43#include "treeBoundBox.H"
44#include "syncTools.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48// NEWER CODE
49// ~~~~~~~~~~
50// Return true on error
51template<class PatchType>
53(
54 const bool allGeometry,
55 const std::string& name,
56 const polyMesh& mesh,
57 const PatchType& pp,
58 const labelUList& meshEdges,
59 labelHashSet* pointSetPtr,
60 labelHashSet* badEdgesPtr
61)
62{
63 if (badEdgesPtr)
64 {
65 badEdgesPtr->clear();
66 }
67
68 typedef typename PatchType::surfaceTopo TopoType;
69
70 bool foundError = false;
71
72 const label globalSize = returnReduce(pp.size(), sumOp<label>());
73
74 Info<< " "
75 << setw(20) << name.c_str()
76 << setw(9) << globalSize
77 << setw(9) << returnReduce(pp.nPoints(), sumOp<label>());
78
79 if (globalSize == 0)
80 {
81 Info<< setw(34) << "ok (empty)";
82 }
83 else if (UPstream::parRun())
84 {
85 // Parallel - use mesh edges
86 // - no check for point-pinch
87 // - no check for consistent orientation (if that is posible to
88 // check?)
89
90 // Count number of edge/face connections (globally)
91 labelList nEdgeConnections(mesh.nEdges(), Zero);
92
93 const labelListList& edgeFaces = pp.edgeFaces();
94
95 forAll(edgeFaces, edgei)
96 {
97 nEdgeConnections[meshEdges[edgei]] = edgeFaces[edgei].size();
98 }
99
100 // Synchronise across coupled edges
101 syncTools::syncEdgeList
102 (
103 mesh,
104 nEdgeConnections,
105 plusEqOp<label>(),
106 label(0) // null value
107 );
108
109 label labelTyp = TopoType::MANIFOLD;
110 forAll(meshEdges, edgei)
111 {
112 const label meshEdgei = meshEdges[edgei];
113 const label numNbrs = nEdgeConnections[meshEdgei];
114 if (numNbrs == 1)
115 {
116 //if (pointSetPtr) pointSetPtr->insert(mesh.edges()[meshEdgei]);
117 labelTyp = max(labelTyp, TopoType::OPEN);
118 }
119 else if (numNbrs == 0 || numNbrs > 2)
120 {
121 if (pointSetPtr) pointSetPtr->insert(mesh.edges()[meshEdgei]);
122 if (badEdgesPtr) badEdgesPtr->insert(edgei);
123 labelTyp = max(labelTyp, TopoType::ILLEGAL);
124 }
125 }
126 reduce(labelTyp, maxOp<label>());
127
128 if (labelTyp == TopoType::MANIFOLD)
129 {
130 Info<< setw(34) << "ok (closed singly connected)";
131 }
132 else if (labelTyp == TopoType::OPEN)
133 {
134 Info<< setw(34) << "ok (non-closed singly connected)";
135 }
136 else
137 {
138 Info<< setw(34) << "multiply connected (shared edge)";
139 }
140
141 foundError = (labelTyp == TopoType::ILLEGAL);
142 }
143 else
144 {
145 const TopoType pTyp = pp.surfaceType(badEdgesPtr);
146
147 if (pTyp == TopoType::MANIFOLD)
148 {
149 if (pp.checkPointManifold(true, pointSetPtr))
150 {
151 Info<< setw(34) << "multiply connected (shared point)";
152 }
153 else
154 {
155 Info<< setw(34) << "ok (closed singly connected)";
156 }
157
158 if (pointSetPtr)
159 {
160 // Add points on non-manifold edges to make set complete
161 pp.checkTopology(false, pointSetPtr);
162 }
163 }
164 else
165 {
166 if (pointSetPtr)
167 {
168 pp.checkTopology(false, pointSetPtr);
169 }
170
171 if (pTyp == TopoType::OPEN)
172 {
173 Info<< setw(34) << "ok (non-closed singly connected)";
174 }
175 else
176 {
177 Info<< setw(34) << "multiply connected (shared edge)";
178 }
179 }
180
181 foundError = (pTyp == TopoType::ILLEGAL);
182 }
183
184 if (allGeometry)
185 {
186 boundBox bb(pp.box());
187 bb.reduce();
188
189 if (bb.good())
190 {
191 Info<< ' ' << bb;
192 }
193 }
194
195 return foundError;
196}
197
198
199// OLDER CODE
200// ~~~~~~~~~~
201// Return true on error
202template<class PatchType>
204(
205 const bool allGeometry,
206 const std::string& name,
207 const polyMesh& mesh,
208 const PatchType& pp,
209 const labelUList& meshFaces,
210 const labelUList& meshEdges,
211 labelHashSet* pointSetPtr,
212 labelHashSet* badEdgesPtr
213)
214{
215 if (badEdgesPtr)
216 {
217 badEdgesPtr->clear();
218 }
219
220 typedef typename PatchType::surfaceTopo TopoType;
221
222 bool foundError = false;
223
224 const label globalSize = returnReduce(pp.size(), sumOp<label>());
225
226 Info<< " "
227 << setw(20) << name.c_str()
228 << setw(9) << globalSize
229 << setw(9) << returnReduce(pp.nPoints(), sumOp<label>());
230
231 if (globalSize == 0)
232 {
233 Info<< setw(34) << "ok (empty)";
234 }
235 else if (UPstream::parRun())
236 {
237 // Parallel - use mesh edges
238 // - no check for point-pinch
239 // - no check for consistent orientation (if that is posible to
240 // check?)
241
242 // OLDER CODE
243 // ~~~~~~~~~~
244 // Synchronise connected faces using global face numbering
245 //
246 // (see addPatchCellLayer::globalEdgeFaces)
247 // From mesh edge to global face labels. Non-empty sublists only for
248 // pp edges.
249 labelListList globalEdgeFaces(mesh.nEdges());
250
251 const labelListList& edgeFaces = pp.edgeFaces();
252
253 // Global numbering
254 const globalIndex globalFaces(mesh.nFaces());
255
256 forAll(edgeFaces, edgei)
257 {
258 const label meshEdgei = meshEdges[edgei];
259 const labelList& eFaces = edgeFaces[edgei];
260
261 // Store face and processor as unique tag.
262 labelList& globalEFaces = globalEdgeFaces[meshEdgei];
263 globalEFaces.resize(eFaces.size());
264 forAll(eFaces, i)
265 {
266 globalEFaces[i] = globalFaces.toGlobal(meshFaces[eFaces[i]]);
267 }
268 //Pout<< "At edge:" << meshEdgei
269 // << " ctr:" << mesh.edges()[meshEdgei].centre(mesh.points())
270 // << " have eFaces:" << globalEdgeFaces[meshEdgei]
271 // << endl;
272 }
273
274 // Synchronise across coupled edges
275 syncTools::syncEdgeList
276 (
277 mesh,
278 globalEdgeFaces,
279 ListOps::uniqueEqOp<label>(),
280 labelList() // null value
281 );
282
283 label labelTyp = TopoType::MANIFOLD;
284 forAll(meshEdges, edgei)
285 {
286 const label meshEdgei = meshEdges[edgei];
287 const labelList& globalEFaces = globalEdgeFaces[meshEdgei];
288 const label numNbrs = globalEFaces.size();
289 if (numNbrs == 1)
290 {
291 //if (pointSetPtr) pointSetPtr->insert(mesh.edges()[meshEdgei]);
292 labelTyp = max(labelTyp, TopoType::OPEN);
293 }
294 else if (numNbrs == 0 || numNbrs > 2)
295 {
296 if (pointSetPtr) pointSetPtr->insert(mesh.edges()[meshEdgei]);
297 if (badEdgesPtr) badEdgesPtr->insert(edgei);
298 labelTyp = max(labelTyp, TopoType::ILLEGAL);
299 }
300 }
301 reduce(labelTyp, maxOp<label>());
302
303 if (labelTyp == TopoType::MANIFOLD)
304 {
305 Info<< setw(34) << "ok (closed singly connected)";
306 }
307 else if (labelTyp == TopoType::OPEN)
308 {
309 Info<< setw(34) << "ok (non-closed singly connected)";
310 }
311 else
312 {
313 Info<< setw(34) << "multiply connected (shared edge)";
314 }
315
316 foundError = (labelTyp == TopoType::ILLEGAL);
317 }
318 else
319 {
320 const TopoType pTyp = pp.surfaceType(badEdgesPtr);
321
322 if (pTyp == TopoType::MANIFOLD)
323 {
324 if (pp.checkPointManifold(true, pointSetPtr))
325 {
326 Info<< setw(34) << "multiply connected (shared point)";
327 }
328 else
329 {
330 Info<< setw(34) << "ok (closed singly connected)";
331 }
332
333 if (pointSetPtr)
334 {
335 // Add points on non-manifold edges to make set complete
336 pp.checkTopology(false, pointSetPtr);
337 }
338 }
339 else
340 {
341 if (pointSetPtr)
342 {
343 pp.checkTopology(false, pointSetPtr);
344 }
345
346 if (pTyp == TopoType::OPEN)
347 {
348 Info<< setw(34) << "ok (non-closed singly connected)";
349 }
350 else
351 {
352 Info<< setw(34) << "multiply connected (shared edge)";
353 }
354 }
355
356 foundError = (pTyp == TopoType::ILLEGAL);
357 }
358
359 if (allGeometry)
360 {
361 boundBox bb(pp.box());
362 bb.reduce();
363
364 if (bb.good())
365 {
366 Info<< ' ' << bb;
367 }
368 }
369
370 return foundError;
371}
372
373
374template<class Zone>
375Foam::label Foam::checkZones
376(
377 const polyMesh& mesh,
378 const ZoneMesh<Zone, polyMesh>& zones,
379 topoSet& set
380)
381{
382 labelList zoneID(set.maxSize(mesh), -1);
383 for (const auto& zone : zones)
384 {
385 for (const label elem : zone)
386 {
387 if
388 (
389 zoneID[elem] != -1
390 && zoneID[elem] != zone.index()
391 )
392 {
393 set.insert(elem);
394 }
395 zoneID[elem] = zone.index();
396 }
397 }
398
399 return returnReduce(set.size(), sumOp<label>());
400}
401
402
403Foam::label Foam::checkTopology
404(
405 const polyMesh& mesh,
406 const bool allTopology,
407 const bool allGeometry,
408 autoPtr<surfaceWriter>& surfWriter,
409 autoPtr<coordSetWriter>& setWriter,
410 const bool writeBadEdges
411)
412{
413 label noFailedChecks = 0;
414
415 Info<< "Checking topology..." << endl;
416
417 // Check if the boundary definition is unique
418 mesh.boundaryMesh().checkDefinition(true);
419
420 // Check that empty patches cover all sides of the mesh
421 {
422 label nEmpty = 0;
423 forAll(mesh.boundaryMesh(), patchi)
424 {
425 if (isA<emptyPolyPatch>(mesh.boundaryMesh()[patchi]))
426 {
427 nEmpty += mesh.boundaryMesh()[patchi].size();
428 }
429 }
430 reduce(nEmpty, sumOp<label>());
431 const label nCells = returnReduce(mesh.cells().size(), sumOp<label>());
432
433 // These are actually warnings, not errors.
434 if (nCells && (nEmpty % nCells))
435 {
436 Info<< " ***Total number of faces on empty patches"
437 << " is not divisible by the number of cells in the mesh."
438 << " Hence this mesh is not 1D or 2D."
439 << endl;
440 }
441 }
442
443 // Check if the boundary processor patches are correct
444 mesh.boundaryMesh().checkParallelSync(true);
445
446 // Check names of zones are equal
447 mesh.cellZones().checkDefinition(true);
448 if (mesh.cellZones().checkParallelSync(true))
449 {
450 noFailedChecks++;
451 }
452 mesh.faceZones().checkDefinition(true);
453 if (mesh.faceZones().checkParallelSync(true))
454 {
455 noFailedChecks++;
456 }
457 mesh.pointZones().checkDefinition(true);
458 if (mesh.pointZones().checkParallelSync(true))
459 {
460 noFailedChecks++;
461 }
462
463
464 {
465 cellSet cells(mesh, "illegalCells", mesh.nCells()/100);
466
467 forAll(mesh.cells(), celli)
468 {
469 const cell& cFaces = mesh.cells()[celli];
470
471 if (cFaces.size() <= 3)
472 {
473 cells.insert(celli);
474 }
475 for (const label facei : cFaces)
476 {
477 if (facei < 0 || facei >= mesh.nFaces())
478 {
479 cells.insert(celli);
480 break;
481 }
482 }
483 }
484 const label nCells = returnReduce(cells.size(), sumOp<label>());
485
486 if (nCells > 0)
487 {
488 Info<< " Illegal cells (less than 4 faces or out of range faces)"
489 << " found, number of cells: " << nCells << endl;
490 noFailedChecks++;
491
492 Info<< " <<Writing " << nCells
493 << " illegal cells to set " << cells.name() << endl;
494 cells.instance() = mesh.pointsInstance();
495 cells.write();
496 if (surfWriter && surfWriter->enabled())
497 {
498 mergeAndWrite(*surfWriter, cells);
499 }
500 }
501 else
502 {
503 Info<< " Cell to face addressing OK." << endl;
504 }
505 }
506
507
508 {
509 pointSet points(mesh, "unusedPoints", mesh.nPoints()/100);
510 if (mesh.checkPoints(true, &points))
511 {
512 noFailedChecks++;
513
514 const label nPoints = returnReduce(points.size(), sumOp<label>());
515
516 Info<< " <<Writing " << nPoints
517 << " unused points to set " << points.name() << endl;
518 points.instance() = mesh.pointsInstance();
519 points.write();
520 if (setWriter && setWriter->enabled())
521 {
522 mergeAndWrite(*setWriter, points);
523 }
524 }
525 }
526
527 {
528 faceSet faces(mesh, "upperTriangularFace", mesh.nFaces()/100);
529 if (mesh.checkUpperTriangular(true, &faces))
530 {
531 noFailedChecks++;
532 }
533
534 const label nFaces = returnReduce(faces.size(), sumOp<label>());
535
536 if (nFaces > 0)
537 {
538 Info<< " <<Writing " << nFaces
539 << " unordered faces to set " << faces.name() << endl;
540 faces.instance() = mesh.pointsInstance();
541 faces.write();
542 if (surfWriter && surfWriter->enabled())
543 {
544 mergeAndWrite(*surfWriter, faces);
545 }
546 }
547 }
548
549 {
550 faceSet faces(mesh, "outOfRangeFaces", mesh.nFaces()/100);
551 if (mesh.checkFaceVertices(true, &faces))
552 {
553 noFailedChecks++;
554
555 const label nFaces = returnReduce(faces.size(), sumOp<label>());
556
557 Info<< " <<Writing " << nFaces
558 << " faces with out-of-range or duplicate vertices to set "
559 << faces.name() << endl;
560 faces.instance() = mesh.pointsInstance();
561 faces.write();
562 if (surfWriter && surfWriter->enabled())
563 {
564 mergeAndWrite(*surfWriter, faces);
565 }
566 }
567 }
568
569 if (allTopology)
570 {
571 cellSet cells(mesh, "zipUpCells", mesh.nCells()/100);
572 if (mesh.checkCellsZipUp(true, &cells))
573 {
574 noFailedChecks++;
575
576 const label nCells = returnReduce(cells.size(), sumOp<label>());
577
578 Info<< " <<Writing " << nCells
579 << " cells with over used edges to set " << cells.name()
580 << endl;
581 cells.instance() = mesh.pointsInstance();
582 cells.write();
583 if (surfWriter && surfWriter->enabled())
584 {
585 mergeAndWrite(*surfWriter, cells);
586 }
587 }
588 }
589
590 if (allTopology)
591 {
592 faceSet faces(mesh, "edgeFaces", mesh.nFaces()/100);
593 if (mesh.checkFaceFaces(true, &faces))
594 {
595 noFailedChecks++;
596 }
597
598 const label nFaces = returnReduce(faces.size(), sumOp<label>());
599 if (nFaces > 0)
600 {
601 Info<< " <<Writing " << nFaces
602 << " faces with non-standard edge connectivity to set "
603 << faces.name() << endl;
604 faces.instance() = mesh.pointsInstance();
605 faces.write();
606 if (surfWriter && surfWriter->enabled())
607 {
608 mergeAndWrite(*surfWriter, faces);
609 }
610 }
611 }
612
613 if (allTopology)
614 {
615 labelList nInternalFaces(mesh.nCells(), Zero);
616
617 for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
618 {
619 nInternalFaces[mesh.faceOwner()[facei]]++;
620 nInternalFaces[mesh.faceNeighbour()[facei]]++;
621 }
622 const polyBoundaryMesh& patches = mesh.boundaryMesh();
623 forAll(patches, patchi)
624 {
625 if (patches[patchi].coupled())
626 {
627 const labelUList& owners = patches[patchi].faceCells();
628
629 for (const label facei : owners)
630 {
631 nInternalFaces[facei]++;
632 }
633 }
634 }
635
636 cellSet oneCells(mesh, "oneInternalFaceCells", mesh.nCells()/100);
637 cellSet twoCells(mesh, "twoInternalFacesCells", mesh.nCells()/100);
638
639 forAll(nInternalFaces, celli)
640 {
641 if (nInternalFaces[celli] <= 1)
642 {
643 oneCells.insert(celli);
644 }
645 else if (nInternalFaces[celli] == 2)
646 {
647 twoCells.insert(celli);
648 }
649 }
650
651 const label nOneCells = returnReduce(oneCells.size(), sumOp<label>());
652
653 if (nOneCells > 0)
654 {
655 Info<< " <<Writing " << nOneCells
656 << " cells with zero or one non-boundary face to set "
657 << oneCells.name()
658 << endl;
659 oneCells.instance() = mesh.pointsInstance();
660 oneCells.write();
661 if (surfWriter && surfWriter->enabled())
662 {
663 mergeAndWrite(*surfWriter, oneCells);
664 }
665 }
666
667 const label nTwoCells = returnReduce(twoCells.size(), sumOp<label>());
668
669 if (nTwoCells > 0)
670 {
671 Info<< " <<Writing " << nTwoCells
672 << " cells with two non-boundary faces to set "
673 << twoCells.name()
674 << endl;
675 twoCells.instance() = mesh.pointsInstance();
676 twoCells.write();
677 if (surfWriter && surfWriter->enabled())
678 {
679 mergeAndWrite(*surfWriter, twoCells);
680 }
681 }
682 }
683
684 {
685 regionSplit rs(mesh);
686
687 if (rs.nRegions() <= 1)
688 {
689 Info<< " Number of regions: " << rs.nRegions() << " (OK)."
690 << endl;
691
692 }
693 else
694 {
695 Info<< " *Number of regions: " << rs.nRegions() << endl;
696
697 Info<< " The mesh has multiple regions which are not connected "
698 "by any face." << endl
699 << " <<Writing region information to "
700 << mesh.time().timeName()/"cellToRegion"
701 << endl;
702
703 IOList<label>::writeContents
704 (
705 IOobject("cellToRegion", mesh.time().timeName(), mesh),
706 rs
707 );
708
709
710 // Points in multiple regions
711 pointSet points
712 (
713 mesh,
714 "multiRegionPoints",
715 mesh.nPoints()/1000
716 );
717
718 // Is region disconnected
719 boolList regionDisconnected(rs.nRegions(), true);
720 if (allTopology)
721 {
722 // -1 : not assigned
723 // -2 : multiple regions
724 // >= 0 : single region
725 labelList pointToRegion(mesh.nPoints(), -1);
726
727 for
728 (
729 label facei = mesh.nInternalFaces();
730 facei < mesh.nFaces();
731 ++facei
732 )
733 {
734 const label regioni = rs[mesh.faceOwner()[facei]];
735 const face& f = mesh.faces()[facei];
736 for (const label verti : f)
737 {
738 label& pRegion = pointToRegion[verti];
739 if (pRegion == -1)
740 {
741 pRegion = regioni;
742 }
743 else if (pRegion == -2)
744 {
745 // Already marked
746 regionDisconnected[regioni] = false;
747 }
748 else if (pRegion != regioni)
749 {
750 // Multiple regions
751 regionDisconnected[regioni] = false;
752 regionDisconnected[pRegion] = false;
753 pRegion = -2;
754 points.insert(verti);
755 }
756 }
757 }
758
759 Pstream::listCombineReduce(regionDisconnected, andEqOp<bool>());
760 }
761
762
763
764 // write cellSet for each region
765 PtrList<cellSet> cellRegions(rs.nRegions());
766 for (label i = 0; i < rs.nRegions(); i++)
767 {
768 cellRegions.set
769 (
770 i,
771 new cellSet
772 (
773 mesh,
774 "region" + Foam::name(i),
775 mesh.nCells()/100
776 )
777 );
778 }
779
780 forAll(rs, i)
781 {
782 cellRegions[rs[i]].insert(i);
783 }
784
785 for (label i = 0; i < rs.nRegions(); i++)
786 {
787 Info<< " <<Writing region " << i;
788 if (allTopology)
789 {
790 if (regionDisconnected[i])
791 {
792 Info<< " (fully disconnected)";
793 }
794 else
795 {
796 Info<< " (point connected)";
797 }
798 }
799 Info<< " with "
800 << returnReduce(cellRegions[i].size(), sumOp<scalar>())
801 << " cells to cellSet " << cellRegions[i].name() << endl;
802
803 cellRegions[i].write();
804 }
805
806 const label nPoints = returnReduce(points.size(), sumOp<label>());
807 if (nPoints)
808 {
809 Info<< " <<Writing " << nPoints
810 << " points that are in multiple regions to set "
811 << points.name() << endl;
812 points.write();
813 if (setWriter && setWriter->enabled())
814 {
815 mergeAndWrite(*setWriter, points);
816 }
817 }
818 }
819 }
820
821 // Non-manifold points
822 pointSet nonManifoldPoints
823 (
824 mesh,
825 "nonManifoldPoints",
826 mesh.nPoints()/1000
827 );
828
829 {
830 Info<< "\nChecking patch topology for multiply connected"
831 << " surfaces..." << endl;
832
833 const polyBoundaryMesh& patches = mesh.boundaryMesh();
834
835 Pout.setf(ios_base::left);
836
837 Info<< " "
838 << setw(20) << "Patch"
839 << setw(9) << "Faces"
840 << setw(9) << "Points"
841 << " Surface topology";
842 if (allGeometry)
843 {
844 Info<< " Bounding box";
845 }
846 Info<< endl;
847
848 for (const polyPatch& pp : patches)
849 {
850 if (!UPstream::parRun() || !isA<processorPolyPatch>(pp))
851 {
853 (
854 allGeometry,
855 pp.name(),
856 mesh,
857 pp,
858 // identity(pp.size(), pp.start()),
859 pp.meshEdges(),
860 &nonManifoldPoints
861 );
862 Info<< endl;
863 }
864 }
865
866 // All non-processor boundary patches
867 {
868 const label nGlobalPatches
869 (
870 UPstream::parRun()
871 ? patches.nNonProcessor()
872 : patches.size()
873 );
874
876 (
877 identity
878 (
879 (
880 patches.range(nGlobalPatches-1).end_value()
881 - patches.start()
882 ),
883 patches.start()
884 )
885 );
886
888 (
889 UIndirectList<face>(mesh.faces(), faceLabels),
890 mesh.points()
891 );
892
893 // Non-manifold
894 labelHashSet badEdges(pp.nEdges()/20);
895
896 bool hadBadEdges = checkPatch
897 (
898 allGeometry,
899 "\".*\"",
900 mesh,
901 pp,
902 // faceLabels,
903 pp.meshEdges(mesh.edges(), mesh.pointEdges()),
904 nullptr, // No point set
905 &badEdges
906 );
907 Info<< nl << endl;
908
909 if (writeBadEdges && hadBadEdges)
910 {
911 edgeList dumpEdges(pp.edges(), badEdges.sortedToc());
912
913 vtk::lineWriter writer
914 (
915 pp.localPoints(),
916 dumpEdges,
917 fileName
918 (
919 mesh.time().globalPath()
920 / ("checkMesh-illegal-edges")
921 )
922 );
923
924 writer.writeGeometry();
925
926 // CellData
927 writer.beginCellData();
928 writer.writeProcIDs();
929
930 Info<< "Wrote "
931 << returnReduce(dumpEdges.size(), sumOp<label>())
932 << " bad edges: " << writer.output().name() << nl;
933 writer.close();
934 }
935 else if (hadBadEdges)
936 {
937 Info<< "Detected "
938 << returnReduce(badEdges.size(), sumOp<label>())
939 << " bad edges (possibly relevant for finite-area)" << nl;
940 }
941 }
942
943 //Info.setf(ios_base::right);
944 }
945
946 {
947 Info<< "\nChecking faceZone topology for multiply connected"
948 << " surfaces..." << endl;
949
950 Pout.setf(ios_base::left);
951
952 const faceZoneMesh& faceZones = mesh.faceZones();
953
954 if (faceZones.size())
955 {
956 Info<< " "
957 << setw(20) << "FaceZone"
958 << setw(9) << "Faces"
959 << setw(9) << "Points"
960 << setw(34) << "Surface topology";
961 if (allGeometry)
962 {
963 Info<< " Bounding box";
964 }
965 Info<< endl;
966
967 for (const faceZone& fz : faceZones)
968 {
970 (
971 allGeometry,
972 fz.name(),
973 mesh,
974 fz(), // patch
975 // fz, // mesh face labels
976 fz.meshEdges(), // mesh edge labels
977 &nonManifoldPoints
978 );
979 Info<< endl;
980 }
981
982 // Check for duplicates
983 if (allTopology)
984 {
985 faceSet mzFaces(mesh, "multiZoneFaces", mesh.nFaces()/100);
986 const label nMulti = checkZones(mesh, faceZones, mzFaces);
987 if (nMulti)
988 {
989 Info<< " <<Writing " << nMulti
990 << " faces that are in multiple zones"
991 << " to set " << mzFaces.name() << endl;
992 mzFaces.instance() = mesh.pointsInstance();
993 mzFaces.write();
994 if (surfWriter && surfWriter->enabled())
995 {
996 mergeAndWrite(*surfWriter, mzFaces);
997 }
998 }
999 }
1000 }
1001 else
1002 {
1003 Info<< " No faceZones found."<<endl;
1004 }
1005 }
1006
1007
1008 const label nPoints =
1009 returnReduce(nonManifoldPoints.size(), sumOp<label>());
1010
1011 if (nPoints)
1012 {
1013 Info<< " <<Writing " << nPoints
1014 << " conflicting points to set " << nonManifoldPoints.name()
1015 << endl;
1016 nonManifoldPoints.instance() = mesh.pointsInstance();
1017 nonManifoldPoints.write();
1018 if (setWriter && setWriter->enabled())
1019 {
1020 mergeAndWrite(*setWriter, nonManifoldPoints);
1021 }
1022 }
1023
1024 {
1025 Info<< "\nChecking basic cellZone addressing..." << endl;
1026
1027 Pout.setf(ios_base::left);
1028
1029 const cellZoneMesh& cellZones = mesh.cellZones();
1030
1031 if (cellZones.size())
1032 {
1033 Info<< " "
1034 << setw(20) << "CellZone"
1035 << setw(13) << "Cells"
1036 << setw(13) << "Points"
1037 << setw(13) << "Volume"
1038 << "BoundingBox" << endl;
1039
1040 const cellList& cells = mesh.cells();
1041 const faceList& faces = mesh.faces();
1042 const scalarField& cellVolumes = mesh.cellVolumes();
1043
1044 bitSet isZonePoint(mesh.nPoints());
1045
1046 for (const cellZone& cZone : cellZones)
1047 {
1048 boundBox bb;
1049 isZonePoint.reset(); // clears all bits (reset count)
1050 scalar v = 0.0;
1051
1052 for (const label celli : cZone)
1053 {
1054 v += cellVolumes[celli];
1055 for (const label facei : cells[celli])
1056 {
1057 const face& f = faces[facei];
1058 for (const label verti : f)
1059 {
1060 if (isZonePoint.set(verti))
1061 {
1062 bb.add(mesh.points()[verti]);
1063 }
1064 }
1065 }
1066 }
1067
1068 bb.reduce(); // Global min/max
1069
1070 Info<< " "
1071 << setw(19) << cZone.name()
1072 << ' ' << setw(12)
1073 << returnReduce(cZone.size(), sumOp<label>())
1074 << ' ' << setw(12)
1075 << returnReduce(isZonePoint.count(), sumOp<label>())
1076 << ' ' << setw(12)
1077 << returnReduce(v, sumOp<scalar>())
1078 << ' ' << bb << endl;
1079 }
1080
1081
1082 // Check for duplicates
1083 if (allTopology)
1084 {
1085 cellSet mzCells(mesh, "multiZoneCells", mesh.nCells()/100);
1086 const label nMulti = checkZones(mesh, cellZones, mzCells);
1087 if (nMulti)
1088 {
1089 Info<< " <<Writing " << nMulti
1090 << " cells that are in multiple zones"
1091 << " to set " << mzCells.name() << endl;
1092 mzCells.instance() = mesh.pointsInstance();
1093 mzCells.write();
1094 if (surfWriter && surfWriter->enabled())
1095 {
1096 mergeAndWrite(*surfWriter, mzCells);
1097 }
1098 }
1099 }
1100 }
1101 else
1102 {
1103 Info<< " No cellZones found."<<endl;
1104 }
1105 }
1106
1107
1108 {
1109 Info<< "\nChecking basic pointZone addressing..." << endl;
1110
1111 Pout.setf(ios_base::left);
1112
1113 const pointZoneMesh& pointZones = mesh.pointZones();
1114
1115 if (pointZones.size())
1116 {
1117 Info<< " "
1118 << setw(20) << "PointZone"
1119 << setw(8) << "Points"
1120 << "BoundingBox" << nl;
1121
1122 for (const auto& zone : pointZones)
1123 {
1124 boundBox bb
1125 (
1126 mesh.points(),
1127 static_cast<const labelUList&>(zone),
1128 true // Reduce (global min/max)
1129 );
1130
1131 Info<< " "
1132 << setw(20) << zone.name()
1133 << setw(8)
1134 << returnReduce(zone.size(), sumOp<label>())
1135 << bb << endl;
1136 }
1137
1138
1139 // Check for duplicates
1140 if (allTopology)
1141 {
1142 pointSet mzPoints(mesh, "multiZonePoints", mesh.nPoints()/100);
1143 const label nMulti = checkZones(mesh, pointZones, mzPoints);
1144 if (nMulti)
1145 {
1146 Info<< " <<Writing " << nMulti
1147 << " points that are in multiple zones"
1148 << " to set " << mzPoints.name() << endl;
1149 mzPoints.instance() = mesh.pointsInstance();
1150 mzPoints.write();
1151 if (setWriter && setWriter->enabled())
1152 {
1153 mergeAndWrite(*setWriter, mzPoints);
1154 }
1155 }
1156 }
1157 }
1158 else
1159 {
1160 Info<< " No pointZones found."<<endl;
1161 }
1162 }
1163
1164
1165 // Force creation of all addressing if requested.
1166 // Errors will be reported as required
1167 if (allTopology)
1168 {
1169 mesh.cells();
1170 mesh.faces();
1171 mesh.edges();
1172 mesh.points();
1173 mesh.faceOwner();
1174 mesh.faceNeighbour();
1175 mesh.cellCells();
1176 mesh.edgeCells();
1177 mesh.pointCells();
1178 mesh.edgeFaces();
1179 mesh.pointFaces();
1180 mesh.cellEdges();
1181 mesh.faceEdges();
1182 mesh.pointEdges();
1183 }
1184
1185 return noFailedChecks;
1186}
1187
1188
1189// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
labelHashSet * pointSetPtr
bool foundError
labelList faceLabels(nFaceLabels)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
labelHashSet badEdges(pp.nEdges()/20)
labelHashSet * badEdgesPtr
std::ios_base::fmtflags setf(std::ios_base::fmtflags f)
Set stream flag(s), return old stream flags.
Definition IOstream.H:506
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const polyBoundaryMesh & patches
bool coupled
dynamicFvMesh & mesh
auto & name
const pointField & points
label nPoints
const cellShapeList & cells
label nMulti
return returnReduce(nRefine-oldNRefine, sumOp< label >())
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Definition BitOps.C:35
List< edge > edgeList
List of edge.
Definition edgeList.H:32
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with pointZone content on a polyMesh.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
label checkTopology(const polyMesh &mesh, const bool allTopology, const bool allGeometry, autoPtr< surfaceWriter > &surfWriter, autoPtr< coordSetWriter > &setWriter, const bool writeBadEdges=false)
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Omanip< int > setw(const int i)
Definition IOmanip.H:199
label checkZones(const polyMesh &mesh, const ZoneMesh< Zone, polyMesh > &zones, topoSet &set)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
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.
List< bool > boolList
A List of bools.
Definition List.H:60
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
bool checkPatch(const bool allGeometry, const std::string &name, const polyMesh &mesh, const PatchType &pp, const labelUList &meshEdges, labelHashSet *pointSetPtr=nullptr, labelHashSet *badEdgesPtr=nullptr)
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299