Loading...
Searching...
No Matches
splitMeshRegions.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-2017 OpenFOAM Foundation
9 Copyright (C) 2015-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
27Application
28 splitMeshRegions
29
30Group
31 grpMeshManipulationUtilities
32
33Description
34 Splits mesh into multiple regions.
35
36 Each region is defined as a domain whose cells can all be reached by
37 cell-face-cell walking without crossing
38 - boundary faces
39 - additional faces from faceset (-blockedFaces faceSet).
40 - any face between differing cellZones (-cellZones)
41
42 Output is:
43 - volScalarField with regions as different scalars (-detectOnly)
44 or
45 - mesh with multiple regions and mapped patches. These patches
46 either cover the whole interface between two region (default) or
47 only part according to faceZones (-useFaceZones)
48 or
49 - mesh with cells put into cellZones (-makeCellZones)
50
51 Note:
52
53 - multiple cellZones can be combined into a single region (cluster)
54 for further analysis using the 'addZones' or 'combineZones' option:
55 -addZones '((allSolids zoneA "zoneB.*")(allFluids none otherZone))'
56 or
57 -combineZones '((zoneA "zoneB.*")(none otherZone))
58 This can be combined with e.g. 'cellZones' or 'cellZonesOnly'. The
59 addZones option supplies the destination region name as first element in
60 the list. The combineZones option synthesises the region name e.g.
61 zoneA_zoneB0_zoneB1
62
63 - cellZonesOnly does not do a walk and uses the cellZones only. Use
64 this if you don't mind having disconnected domains in a single region.
65 This option requires all cells to be in one (and one only) cellZone.
66
67 - cellZonesFileOnly behaves like -cellZonesOnly but reads the cellZones
68 from the specified file. This allows one to explicitly specify the region
69 distribution and still have multiple cellZones per region.
70
71 - prefixRegion prefixes all normal patches with region name (interface
72 (patches already have region name prefix)
73
74 - Should work in parallel.
75 cellZones can differ on either side of processor boundaries in which case
76 the faces get moved from processor patch to mapped patch. Not very well
77 tested.
78
79 - If a cell zone gets split into more than one region it can detect
80 the largest matching region (-sloppyCellZones). This will accept any
81 region that covers more than 50% of the zone. It has to be a subset
82 so cannot have any cells in any other zone.
83
84 - If explicitly a single region has been selected (-largestOnly or
85 -insidePoint) its region name will be either
86 - name of a cellZone it matches to or
87 - "largestOnly" respectively "insidePoint" or
88 - polyMesh::defaultRegion if additionally -overwrite
89 (so it will overwrite the input mesh!)
90
91 - writes maps like decomposePar back to original mesh:
92 - pointRegionAddressing : for every point in this region the point in
93 the original mesh
94 - cellRegionAddressing : ,, cell ,, cell ,,
95 - faceRegionAddressing : ,, face ,, face in
96 the original mesh + 'turning index'. For a face in the same orientation
97 this is the original facelabel+1, for a turned face this is -facelabel-1
98 - boundaryRegionAddressing : for every patch in this region the
99 patch in the original mesh (or -1 if added patch)
100
101\*---------------------------------------------------------------------------*/
102
103#include "argList.H"
104#include "regionSplit.H"
105#include "fvMeshSubset.H"
106#include "IOobjectList.H"
107#include "volFields.H"
108#include "faceSet.H"
109#include "cellSet.H"
110#include "polyTopoChange.H"
111#include "removeCells.H"
112#include "edgeHashes.H"
113#include "syncTools.H"
114#include "ReadFields.H"
115#include "mappedWallPolyPatch.H"
116#include "fvMeshTools.H"
117#include "processorMeshes.H"
118
119using namespace Foam;
120
121// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122
123// Prepend prefix to selected patches.
124void renamePatches
125(
126 fvMesh& mesh,
127 const word& prefix,
128 const labelList& patchesToRename
129)
130{
131 polyBoundaryMesh& polyPatches =
132 const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
133 forAll(patchesToRename, i)
134 {
135 label patchi = patchesToRename[i];
136 polyPatch& pp = polyPatches[patchi];
137
139 {
141 << "Encountered coupled patch " << pp.name()
142 << ". Will only rename the patch itself,"
143 << " not any referred patches."
144 << " This might have to be done by hand."
145 << endl;
146 }
147
148 pp.name() = prefix + '_' + pp.name();
149 }
150 // Recalculate any demand driven data (e.g. group to name lookup)
151 polyPatches.updateMesh();
152}
153
154
155template<class GeoField>
156void subsetVolFields
157(
158 const fvMesh& mesh,
159 const fvMesh& subMesh,
160 const labelList& cellMap,
161 const labelList& faceMap,
162 const labelHashSet& addedPatches
163)
164{
165 const labelList patchMap(identity(mesh.boundaryMesh().size()));
166
167 for (const GeoField& fld : mesh.objectRegistry::csorted<GeoField>())
168 {
169 Info<< "Mapping field " << fld.name() << endl;
170
171 tmp<GeoField> tSubFld
172 (
174 (
175 fld,
176 subMesh,
177 patchMap,
178 cellMap,
179 faceMap
180 )
181 );
182
183 // Hack: set value to 0 for introduced patches (since don't
184 // get initialised.
185 forAll(tSubFld().boundaryField(), patchi)
186 {
187 if (addedPatches.found(patchi))
188 {
189 tSubFld.ref().boundaryFieldRef()[patchi] == Zero;
190 }
191 }
192
193 // Store on subMesh
194 GeoField* subFld = tSubFld.ptr();
195 subFld->rename(fld.name());
196 subFld->writeOpt(IOobject::AUTO_WRITE);
197 subFld->store();
198 }
199}
200
201
202template<class GeoField>
203void subsetSurfaceFields
204(
205 const fvMesh& mesh,
206 const fvMesh& subMesh,
207 const labelList& cellMap,
208 const labelList& faceMap,
209 const labelHashSet& addedPatches
210)
211{
212 const labelList patchMap(identity(mesh.boundaryMesh().size()));
213
214 for (const GeoField& fld : mesh.objectRegistry::csorted<GeoField>())
215 {
216 Info<< "Mapping field " << fld.name() << endl;
217
218 tmp<GeoField> tSubFld
219 (
221 (
222 fld,
223 subMesh,
224 patchMap,
225 cellMap,
226 faceMap
227 )
228 );
229
230 // Hack: set value to 0 for introduced patches (since don't
231 // get initialised.
232 forAll(tSubFld().boundaryField(), patchi)
233 {
234 if (addedPatches.found(patchi))
235 {
236 tSubFld.ref().boundaryFieldRef()[patchi] == Zero;
237 }
238 }
239
240 // Store on subMesh
241 GeoField* subFld = tSubFld.ptr();
242 subFld->rename(fld.name());
243 subFld->writeOpt(IOobject::AUTO_WRITE);
244 subFld->store();
245 }
246}
247
248// Select all cells not in the region
249labelList getNonRegionCells(const labelList& cellRegion, const label regionI)
250{
251 DynamicList<label> nonRegionCells(cellRegion.size());
252 forAll(cellRegion, celli)
253 {
254 if (cellRegion[celli] != regionI)
255 {
256 nonRegionCells.append(celli);
257 }
258 }
259 return labelList(std::move(nonRegionCells));
260}
261
262
263void addToInterface
264(
265 const polyMesh& mesh,
266 const label zoneID,
267 const label ownRegion,
268 const label neiRegion,
269 EdgeMap<Map<label>>& regionsToSize
270)
271{
272 const auto interface(edge::sorted(ownRegion, neiRegion));
273
274 auto iter = regionsToSize.find(interface);
275
276 if (iter.good())
277 {
278 // Check if zone present
279 auto zoneIter = iter().find(zoneID);
280 if (zoneIter.good())
281 {
282 // Found zone. Increment count.
283 ++(*zoneIter);
284 }
285 else
286 {
287 // New or no zone. Insert with count 1.
288 iter().insert(zoneID, 1);
289 }
290 }
291 else
292 {
293 // Create new interface of size 1.
294 regionsToSize(interface).insert(zoneID, 1);
295 }
296}
297
298
299// Get region-region interface name and sizes.
300// Returns interfaces as straight list for looping in identical order.
301void getInterfaceSizes
302(
303 const polyMesh& mesh,
304 const bool useFaceZones,
305 const labelList& cellRegion,
306 const wordList& regionNames,
307
308 edgeList& interfaces,
309 List<Pair<word>>& interfaceNames,
310 labelList& interfaceSizes,
311 labelList& faceToInterface
312)
313{
314 // From region-region to faceZone (or -1) to number of faces.
315
316 EdgeMap<Map<label>> regionsToSize;
317
318
319 // Internal faces
320 // ~~~~~~~~~~~~~~
321
322 forAll(mesh.faceNeighbour(), facei)
323 {
324 label ownRegion = cellRegion[mesh.faceOwner()[facei]];
325 label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
326
327 if (ownRegion != neiRegion)
328 {
329 addToInterface
330 (
331 mesh,
332 (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
333 ownRegion,
334 neiRegion,
335 regionsToSize
336 );
337 }
338 }
339
340 // Boundary faces
341 // ~~~~~~~~~~~~~~
342
343 // Neighbour cellRegion.
344 labelList coupledRegion(mesh.nBoundaryFaces());
345
346 forAll(coupledRegion, i)
347 {
348 label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
349 coupledRegion[i] = cellRegion[celli];
350 }
352
353 forAll(coupledRegion, i)
354 {
355 label facei = i+mesh.nInternalFaces();
356 label ownRegion = cellRegion[mesh.faceOwner()[facei]];
357 label neiRegion = coupledRegion[i];
358
359 if (ownRegion != neiRegion)
360 {
361 addToInterface
362 (
363 mesh,
364 (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
365 ownRegion,
366 neiRegion,
367 regionsToSize
368 );
369 }
370 }
371
372
373 if (UPstream::parRun())
374 {
375 if (UPstream::master())
376 {
377 // Receive and add to my sizes
378 for (const int proci : UPstream::subProcs())
379 {
380 EdgeMap<Map<label>> slaveSizes;
381 IPstream::recv(slaveSizes, proci);
382
383 forAllConstIters(slaveSizes, slaveIter)
384 {
385 const Map<label>& slaveInfo = *slaveIter;
386
387 auto masterIter = regionsToSize.find(slaveIter.key());
388
389 if (masterIter.good())
390 {
391 // Same inter-region
392 Map<label>& masterInfo = *masterIter;
393
394 forAllConstIters(slaveInfo, iter)
395 {
396 const label zoneID = iter.key();
397 const label slaveSize = iter.val();
398
399 auto zoneIter = masterInfo.find(zoneID);
400 if (zoneIter.good())
401 {
402 *zoneIter += slaveSize;
403 }
404 else
405 {
406 masterInfo.insert(zoneID, slaveSize);
407 }
408 }
409 }
410 else
411 {
412 regionsToSize.insert(slaveIter.key(), slaveInfo);
413 }
414 }
415 }
416 }
417 else
418 {
419 // send to master
420 OPstream::send(regionsToSize, UPstream::masterNo());
421 }
422 }
423
424 // Rework
425
426 Pstream::broadcast(regionsToSize);
427
428
429 // Now we have the global sizes of all inter-regions.
430 // Invert this on master and distribute.
431 label nInterfaces = 0;
432 forAllConstIters(regionsToSize, iter)
433 {
434 const Map<label>& info = iter.val();
435 nInterfaces += info.size();
436 }
437
438 interfaces.setSize(nInterfaces);
439 interfaceNames.setSize(nInterfaces);
440 interfaceSizes.setSize(nInterfaces);
441 EdgeMap<Map<label>> regionsToInterface(nInterfaces);
442
443 nInterfaces = 0;
444 forAllConstIters(regionsToSize, iter)
445 {
446 const edge& e = iter.key();
447 const Map<label>& info = iter.val();
448
449 const word& name0 = regionNames[e[0]];
450 const word& name1 = regionNames[e[1]];
451
452 forAllConstIters(info, infoIter)
453 {
454 interfaces[nInterfaces] = iter.key();
455 label zoneID = infoIter.key();
456 if (zoneID == -1)
457 {
458 interfaceNames[nInterfaces] = Pair<word>
459 (
460 name0 + "_to_" + name1,
461 name1 + "_to_" + name0
462 );
463 }
464 else
465 {
466 const word& zoneName = mesh.faceZones()[zoneID].name();
467 interfaceNames[nInterfaces] = Pair<word>
468 (
469 zoneName + "_" + name0 + "_to_" + name1,
470 zoneName + "_" + name1 + "_to_" + name0
471 );
472 }
473
474 interfaceSizes[nInterfaces] = infoIter();
475 regionsToInterface(e).insert(zoneID, nInterfaces);
476
477 nInterfaces++;
478 }
479 }
480
481
482 // Consistent interface information for all processors
484 (
486 interfaces,
487 interfaceNames,
488 interfaceSizes,
489 regionsToInterface
490 );
491
492 // Mark all inter-region faces.
493 faceToInterface.setSize(mesh.nFaces(), -1);
494
495 forAll(mesh.faceNeighbour(), facei)
496 {
497 label ownRegion = cellRegion[mesh.faceOwner()[facei]];
498 label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
499
500 if (ownRegion != neiRegion)
501 {
502 const auto interface(edge::sorted(ownRegion, neiRegion));
503
504 label zoneID = -1;
505 if (useFaceZones)
506 {
507 zoneID = mesh.faceZones().whichZone(facei);
508 }
509
510 faceToInterface[facei] = regionsToInterface[interface][zoneID];
511 }
512 }
513 forAll(coupledRegion, i)
514 {
515 label facei = i+mesh.nInternalFaces();
516 label ownRegion = cellRegion[mesh.faceOwner()[facei]];
517 label neiRegion = coupledRegion[i];
518
519 if (ownRegion != neiRegion)
520 {
521 const auto interface(edge::sorted(ownRegion, neiRegion));
522
523 label zoneID = -1;
524 if (useFaceZones)
525 {
526 zoneID = mesh.faceZones().whichZone(facei);
527 }
528
529 faceToInterface[facei] = regionsToInterface[interface][zoneID];
530 }
531 }
532}
533
534
535// Create mesh for region.
536autoPtr<mapPolyMesh> createRegionMesh
537(
538 const fvMesh& mesh,
539 // Region info
540 const labelList& cellRegion,
541 const label regionI,
542 const word& regionName,
543 // Interface info
544 const labelList& interfacePatches,
545 const labelList& faceToInterface,
546
547 autoPtr<fvMesh>& newMesh
548)
549{
550 // Create dummy system/fv*
552
553 // Neighbour cellRegion.
554 labelList coupledRegion(mesh.nBoundaryFaces());
555
556 forAll(coupledRegion, i)
557 {
558 label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
559 coupledRegion[i] = cellRegion[celli];
560 }
562
563
564 // Topology change container. Start off from existing mesh.
565 polyTopoChange meshMod(mesh);
566
567 // Cell remover engine
568 removeCells cellRemover(mesh);
569
570 // Select all but region cells
571 labelList cellsToRemove(getNonRegionCells(cellRegion, regionI));
572
573 // Find out which faces will get exposed. Note that this
574 // gets faces in mesh face order. So both regions will get same
575 // face in same order (important!)
576 labelList exposedFaces = cellRemover.getExposedFaces(cellsToRemove);
577
578 labelList exposedPatchIDs(exposedFaces.size());
579 forAll(exposedFaces, i)
580 {
581 label facei = exposedFaces[i];
582 label interfacei = faceToInterface[facei];
583
584 label ownRegion = cellRegion[mesh.faceOwner()[facei]];
585 label neiRegion = -1;
586
587 if (mesh.isInternalFace(facei))
588 {
589 neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
590 }
591 else
592 {
593 neiRegion = coupledRegion[facei-mesh.nInternalFaces()];
594 }
595
596
597 // Check which side is being kept - determines which of the two
598 // patches will be used.
599
600 label otherRegion = -1;
601
602 if (ownRegion == regionI && neiRegion != regionI)
603 {
604 otherRegion = neiRegion;
605 }
606 else if (ownRegion != regionI && neiRegion == regionI)
607 {
608 otherRegion = ownRegion;
609 }
610 else
611 {
613 << "Exposed face:" << facei
614 << " fc:" << mesh.faceCentres()[facei]
615 << " has owner region " << ownRegion
616 << " and neighbour region " << neiRegion
617 << " when handling region:" << regionI
618 << exit(FatalError);
619 }
620
621 // Find the patch.
622 if (regionI < otherRegion)
623 {
624 exposedPatchIDs[i] = interfacePatches[interfacei];
625 }
626 else
627 {
628 exposedPatchIDs[i] = interfacePatches[interfacei]+1;
629 }
630 }
631
632 // Remove faces
633 cellRemover.setRefinement
634 (
635 cellsToRemove,
636 exposedFaces,
637 exposedPatchIDs,
638 meshMod
639 );
640
641 autoPtr<mapPolyMesh> map = meshMod.makeMesh
642 (
643 newMesh,
645 (
647 mesh.time().timeName(),
648 mesh.time(),
649 IOobject::READ_IF_PRESENT, // read fv* if present
651 ),
652 mesh
653 );
654
655 return map;
656}
657
658
659void createAndWriteRegion
660(
661 const fvMesh& mesh,
662 const labelList& cellRegion,
663 const wordList& regionNames,
664 const bool prefixRegion,
665 const labelList& faceToInterface,
666 const labelList& interfacePatches,
667 const label regionI,
668 const word& newMeshInstance
669)
670{
671 Info<< "Creating mesh for region " << regionI
672 << ' ' << regionNames[regionI] << endl;
673
674 autoPtr<fvMesh> newMesh;
675 autoPtr<mapPolyMesh> map = createRegionMesh
676 (
677 mesh,
678 cellRegion,
679 regionI,
680 regionNames[regionI],
681 interfacePatches,
682 faceToInterface,
683 newMesh
684 );
685
686
687 // Make map of all added patches
688 labelHashSet addedPatches(2*interfacePatches.size());
689 forAll(interfacePatches, interfacei)
690 {
691 addedPatches.insert(interfacePatches[interfacei]);
692 addedPatches.insert(interfacePatches[interfacei]+1);
693 }
694
695
696 Info<< "Mapping fields" << endl;
697
698 // Map existing fields
699 newMesh().updateMesh(map());
700
701 // Add subsetted fields
702 subsetVolFields<volScalarField>
703 (
704 mesh,
705 newMesh(),
706 map().cellMap(),
707 map().faceMap(),
708 addedPatches
709 );
710 subsetVolFields<volVectorField>
711 (
712 mesh,
713 newMesh(),
714 map().cellMap(),
715 map().faceMap(),
716 addedPatches
717 );
718 subsetVolFields<volSphericalTensorField>
719 (
720 mesh,
721 newMesh(),
722 map().cellMap(),
723 map().faceMap(),
724 addedPatches
725 );
726 subsetVolFields<volSymmTensorField>
727 (
728 mesh,
729 newMesh(),
730 map().cellMap(),
731 map().faceMap(),
732 addedPatches
733 );
734 subsetVolFields<volTensorField>
735 (
736 mesh,
737 newMesh(),
738 map().cellMap(),
739 map().faceMap(),
740 addedPatches
741 );
742
743 subsetSurfaceFields<surfaceScalarField>
744 (
745 mesh,
746 newMesh(),
747 map().cellMap(),
748 map().faceMap(),
749 addedPatches
750 );
751 subsetSurfaceFields<surfaceVectorField>
752 (
753 mesh,
754 newMesh(),
755 map().cellMap(),
756 map().faceMap(),
757 addedPatches
758 );
759 subsetSurfaceFields<surfaceSphericalTensorField>
760 (
761 mesh,
762 newMesh(),
763 map().cellMap(),
764 map().faceMap(),
765 addedPatches
766 );
767 subsetSurfaceFields<surfaceSymmTensorField>
768 (
769 mesh,
770 newMesh(),
771 map().cellMap(),
772 map().faceMap(),
773 addedPatches
774 );
775 subsetSurfaceFields<surfaceTensorField>
776 (
777 mesh,
778 newMesh(),
779 map().cellMap(),
780 map().faceMap(),
781 addedPatches
782 );
783
784
785 const polyBoundaryMesh& newPatches = newMesh().boundaryMesh();
786 newPatches.checkParallelSync(true);
787
788 // Delete empty patches
789 // ~~~~~~~~~~~~~~~~~~~~
790
791 // Create reordering list to move patches-to-be-deleted to end
792 labelList oldToNew(newPatches.size(), -1);
793 DynamicList<label> sharedPatches(newPatches.size());
794 label newI = 0;
795
796 Info<< "Deleting empty patches" << endl;
797
798 // Assumes all non-proc boundaries are on all processors!
799 forAll(newPatches, patchi)
800 {
801 const polyPatch& pp = newPatches[patchi];
802
804 {
805 if (returnReduceOr(pp.size()))
806 {
807 oldToNew[patchi] = newI;
808 if (!addedPatches.found(patchi))
809 {
810 sharedPatches.append(newI);
811 }
812 newI++;
813 }
814 }
815 }
816
817 // Same for processor patches (but need no reduction)
818 forAll(newPatches, patchi)
819 {
820 const polyPatch& pp = newPatches[patchi];
821
822 if (isA<processorPolyPatch>(pp) && pp.size())
823 {
824 oldToNew[patchi] = newI++;
825 }
826 }
827
828 const label nNewPatches = newI;
829
830 // Move all deleteable patches to the end
831 forAll(oldToNew, patchi)
832 {
833 if (oldToNew[patchi] == -1)
834 {
835 oldToNew[patchi] = newI++;
836 }
837 }
838
839 //reorderPatches(newMesh(), oldToNew, nNewPatches);
840 fvMeshTools::reorderPatches(newMesh(), oldToNew, nNewPatches, true);
841
842 // Rename shared patches with region name
843 if (prefixRegion)
844 {
845 Info<< "Prefixing patches with region name" << endl;
846
847 renamePatches(newMesh(), regionNames[regionI], sharedPatches);
848 }
849
850
851 Info<< "Writing new mesh" << endl;
852
853 newMesh().setInstance(newMeshInstance);
854 newMesh().write();
855 topoSet::removeFiles(newMesh());
857
858 // Write addressing files like decomposePar
859 Info<< "Writing addressing to base mesh" << endl;
860
861 labelIOList pointProcAddressing
862 (
864 (
865 "pointRegionAddressing",
866 newMesh().facesInstance(),
868 newMesh(),
872 ),
873 map().pointMap()
874 );
875 Info<< "Writing map " << pointProcAddressing.name()
876 << " from region" << regionI
877 << " points back to base mesh." << endl;
878 pointProcAddressing.write();
879
880 labelIOList faceProcAddressing
881 (
883 (
884 "faceRegionAddressing",
885 newMesh().facesInstance(),
887 newMesh(),
891 ),
892 newMesh().nFaces()
893 );
894 forAll(faceProcAddressing, facei)
895 {
896 // face + turning index. (see decomposePar)
897 // Is the face pointing in the same direction?
898 label oldFacei = map().faceMap()[facei];
899
900 if
901 (
902 map().cellMap()[newMesh().faceOwner()[facei]]
903 == mesh.faceOwner()[oldFacei]
904 )
905 {
906 faceProcAddressing[facei] = oldFacei+1;
907 }
908 else
909 {
910 faceProcAddressing[facei] = -(oldFacei+1);
911 }
912 }
913 Info<< "Writing map " << faceProcAddressing.name()
914 << " from region" << regionI
915 << " faces back to base mesh." << endl;
916 faceProcAddressing.write();
917
918 labelIOList cellProcAddressing
919 (
921 (
922 "cellRegionAddressing",
923 newMesh().facesInstance(),
925 newMesh(),
929 ),
930 map().cellMap()
931 );
932 Info<< "Writing map " <<cellProcAddressing.name()
933 << " from region" << regionI
934 << " cells back to base mesh." << endl;
935 cellProcAddressing.write();
936
937 labelIOList boundaryProcAddressing
938 (
940 (
941 "boundaryRegionAddressing",
942 newMesh().facesInstance(),
944 newMesh(),
948 ),
949 labelList(nNewPatches, -1)
950 );
951 forAll(oldToNew, i)
952 {
953 if (!addedPatches.found(i))
954 {
955 label newI = oldToNew[i];
956 if (newI >= 0 && newI < nNewPatches)
957 {
958 boundaryProcAddressing[oldToNew[i]] = i;
959 }
960 }
961 }
962 Info<< "Writing map " << boundaryProcAddressing.name()
963 << " from region" << regionI
964 << " boundary back to base mesh." << endl;
965 boundaryProcAddressing.write();
966}
967
968
969// Create for every region-region interface with non-zero size two patches.
970// First one is for minimumregion to maximumregion.
971// Note that patches get created in same order on all processors (if parallel)
972// since looping over synchronised 'interfaces'.
973labelList addRegionPatches
974(
975 fvMesh& mesh,
976 const wordList& regionNames,
977 const edgeList& interfaces,
978 const List<Pair<word>>& interfaceNames
979)
980{
981 Info<< nl << "Adding patches" << nl << endl;
982
983 labelList interfacePatches(interfaces.size());
984
985 forAll(interfaces, interI)
986 {
987 const edge& e = interfaces[interI];
988 const Pair<word>& names = interfaceNames[interI];
989
990 //Info<< "For interface " << interI
991 // << " between regions " << e
992 // << " trying to add patches " << names << endl;
993
994
996 (
997 names[0],
998 0, // overridden
999 0, // overridden
1000 0, // overridden
1001 regionNames[e[1]], // sampleRegion
1003 names[1], // samplePatch
1004 point::zero, // offset
1005 mesh.boundaryMesh()
1006 );
1007
1008 interfacePatches[interI] = fvMeshTools::addPatch
1009 (
1010 mesh,
1011 patch1,
1012 dictionary(), //optional per field value
1014 true //validBoundary
1015 );
1016
1017 mappedWallPolyPatch patch2
1018 (
1019 names[1],
1020 0,
1021 0,
1022 0,
1023 regionNames[e[0]], // sampleRegion
1025 names[0],
1026 point::zero, // offset
1027 mesh.boundaryMesh()
1028 );
1030 (
1031 mesh,
1032 patch2,
1033 dictionary(), //optional per field value
1035 true //validBoundary
1036 );
1037
1038 Info<< "For interface between region " << regionNames[e[0]]
1039 << " and " << regionNames[e[1]] << " added patches" << endl
1040 << " " << interfacePatches[interI]
1041 << "\t" << mesh.boundaryMesh()[interfacePatches[interI]].name()
1042 << endl
1043 << " " << interfacePatches[interI]+1
1044 << "\t" << mesh.boundaryMesh()[interfacePatches[interI]+1].name()
1045 << endl;
1046 }
1047 return interfacePatches;
1048}
1049
1050
1051// Find region that covers most of cell zone
1052label findCorrespondingRegion
1053(
1054 const labelList& existingZoneID, // per cell the (unique) zoneID
1055 const labelList& cellRegion,
1056 const label nCellRegions,
1057 const label zoneI,
1058 const label minOverlapSize
1059)
1060{
1061 // Per region the number of cells in zoneI
1062 labelList cellsInZone(nCellRegions, Zero);
1063
1064 forAll(cellRegion, celli)
1065 {
1066 if (existingZoneID[celli] == zoneI)
1067 {
1068 cellsInZone[cellRegion[celli]]++;
1069 }
1070 }
1071
1072 Pstream::listReduce(cellsInZone, sumOp<label>());
1073
1074 // Pick region with largest overlap of zoneI
1075 label regionI = findMax(cellsInZone);
1076
1077
1078 if (cellsInZone[regionI] < minOverlapSize)
1079 {
1080 // Region covers too little of zone. Not good enough.
1081 regionI = -1;
1082 }
1083 else
1084 {
1085 // Check that region contains no cells that aren't in cellZone.
1086 forAll(cellRegion, celli)
1087 {
1088 if (cellRegion[celli] == regionI && existingZoneID[celli] != zoneI)
1089 {
1090 // celli in regionI but not in zoneI
1091 regionI = -1;
1092 break;
1093 }
1094 }
1095 // If one in error, all should be in error. Note that branch gets taken
1096 // on all procs.
1097 reduce(regionI, minOp<label>());
1098 }
1099
1100 return regionI;
1101}
1102
1103
1104void getClusterID
1105(
1106 const polyMesh& mesh,
1107 const cellZoneMesh& cellZones,
1108 const wordList& clusterNames,
1109 const labelListList& clusterToZones,
1110 labelList& clusterID,
1111 labelList& neiClusterID
1112)
1113{
1114 // Existing zoneID
1115 clusterID.setSize(mesh.nCells());
1116 clusterID = -1;
1117
1118 forAll(clusterToZones, clusterI)
1119 {
1120 for (const label zoneI : clusterToZones[clusterI])
1121 {
1122 const cellZone& cz = cellZones[zoneI];
1123
1124 forAll(cz, i)
1125 {
1126 label celli = cz[i];
1127 if (clusterID[celli] == -1)
1128 {
1129 clusterID[celli] = clusterI;
1130 }
1131 else
1132 {
1134 << "Cell " << celli << " with cell centre "
1135 << mesh.cellCentres()[celli]
1136 << " is multiple zones. This is not allowed." << endl
1137 << "It is in zone " << clusterNames[clusterID[celli]]
1138 << " and in zone " << clusterNames[clusterI]
1139 << exit(FatalError);
1140 }
1141 }
1142 }
1143 }
1144
1145 // Neighbour zoneID.
1146 syncTools::swapBoundaryCellList(mesh, clusterID, neiClusterID);
1147}
1148
1149
1150word makeRegionName
1151(
1152 const cellZoneMesh& czs,
1153 const label regioni,
1154 const labelList& zoneIDs
1155)
1156{
1157 // Synthesise region name. Equals the zone name if cluster consist of only
1158 // one zone
1159
1160 if (zoneIDs.empty())
1161 {
1162 return word("domain") + Foam::name(regioni);
1163 }
1164 else
1165 {
1166 // Zone indices are in cellZone order ...
1167 word regionName(czs[zoneIDs[0]].name());
1168
1169 // Synthesize name from appended cellZone names
1170 for (label i = 1; i < zoneIDs.size(); i++)
1171 {
1172 regionName += "_" + czs[zoneIDs[i]].name();
1173 }
1174 return regionName;
1175 }
1176}
1177
1178
1179void makeClusters
1180(
1181 const List<wordRes>& zoneClusters,
1182 const wordList& zoneClusterNames,
1183 const cellZoneMesh& cellZones,
1184 wordList& clusterNames,
1185 labelListList& clusterToZones,
1186 labelList& zoneToCluster
1187)
1188{
1189 // Check if there are clustering for zones. If none every zone goes into
1190 // its own cluster.
1191
1192 clusterNames.clear();
1193 clusterToZones.clear();
1194 zoneToCluster.setSize(cellZones.size());
1195 zoneToCluster = -1;
1196
1197 if (zoneClusters.size())
1198 {
1199 forAll(zoneClusters, clusteri)
1200 {
1201 const labelList zoneIDs(cellZones.indices(zoneClusters[clusteri]));
1202 UIndirectList<label>(zoneToCluster, zoneIDs) = clusteri;
1203 clusterNames.append
1204 (
1205 zoneClusterNames[clusteri].size()
1206 ? zoneClusterNames[clusteri]
1207 : makeRegionName
1208 (
1209 cellZones,
1210 clusteri,
1211 zoneIDs
1212 )
1213 );
1214 clusterToZones.append(std::move(zoneIDs));
1215 }
1216
1217 // Unclustered zone
1218 forAll(zoneToCluster, zonei)
1219 {
1220 if (zoneToCluster[zonei] == -1)
1221 {
1222 clusterNames.append(cellZones[zonei].name());
1223 clusterToZones.append(labelList(1, zonei));
1224 zoneToCluster[zonei] = clusterToZones.size();
1225 }
1226 }
1227 }
1228 else
1229 {
1230 for (const auto& cellZone : cellZones)
1231 {
1232 const label nClusters = clusterToZones.size();
1233 clusterNames.append(cellZone.name());
1234 clusterToZones.append(labelList(1, cellZone.index()));
1235 zoneToCluster[cellZone.index()] = nClusters;
1236 }
1237 }
1238}
1239
1240
1241void matchRegions
1242(
1243 const bool sloppyCellZones,
1244 const polyMesh& mesh,
1245
1246 const wordList& clusterNames,
1247 const labelListList& clusterToZones,
1248 const labelList& clusterID,
1249
1250 const label nCellRegions,
1251 const labelList& cellRegion,
1252
1253 labelListList& regionToZones,
1255 labelList& zoneToRegion
1256)
1257{
1258 const cellZoneMesh& cellZones = mesh.cellZones();
1259
1260 regionToZones.setSize(nCellRegions);
1261 regionToZones = labelList();
1262 regionNames.setSize(nCellRegions);
1264 zoneToRegion.setSize(cellZones.size(), -1);
1265
1266
1267 // Sizes per cluster
1268 labelList clusterSizes(clusterToZones.size(), Zero);
1269 forAll(clusterToZones, clusterI)
1270 {
1271 for (const label zoneI : clusterToZones[clusterI])
1272 {
1273 clusterSizes[clusterI] += cellZones[zoneI].size();
1274 }
1275 reduce(clusterSizes[clusterI], sumOp<label>());
1276 }
1277
1278 if (sloppyCellZones)
1279 {
1280 Info<< "Trying to match regions to existing cell zones;"
1281 << " region can be subset of cell zone." << nl << endl;
1282
1283 forAll(clusterToZones, clusterI)
1284 {
1285 label regionI = findCorrespondingRegion
1286 (
1287 clusterID,
1288 cellRegion,
1289 nCellRegions,
1290 clusterI,
1291 label(0.5*clusterSizes[clusterI]) // minimum overlap
1292 );
1293
1294 if (regionI != -1)
1295 {
1296 Info<< "Sloppily matched region " << regionI
1297 //<< " size " << regionSizes[regionI]
1298 << " to cluster " << clusterI
1299 << " size " << clusterSizes[clusterI]
1300 << endl;
1302 (
1303 zoneToRegion,
1304 clusterToZones[clusterI]
1305 ) = regionI;
1306 regionToZones[regionI] = clusterToZones[clusterI];
1307 regionNames[regionI] = clusterNames[clusterI];
1308 }
1309 }
1310 }
1311 else
1312 {
1313 Info<< "Trying to match regions to existing cell zones." << nl << endl;
1314
1315 forAll(clusterToZones, clusterI)
1316 {
1317 label regionI = findCorrespondingRegion
1318 (
1319 clusterID,
1320 cellRegion,
1321 nCellRegions,
1322 clusterI,
1323 clusterSizes[clusterI] // require exact match
1324 );
1325
1326 if (regionI != -1)
1327 {
1329 (
1330 zoneToRegion,
1331 clusterToZones[clusterI]
1332 ) = regionI;
1333 regionToZones[regionI] = clusterToZones[clusterI];
1334 regionNames[regionI] = clusterNames[clusterI];
1335 }
1336 }
1337 }
1338 // Allocate region names for unmatched regions.
1339 forAll(regionNames, regionI)
1340 {
1341 if (regionNames[regionI].empty())
1342 {
1343 regionNames[regionI] = makeRegionName
1344 (
1345 cellZones,
1346 regionI,
1347 regionToZones[regionI]
1348 );
1349 }
1350 }
1351}
1352
1353
1354void writeCellToRegion(const fvMesh& mesh, const labelList& cellRegion)
1355{
1356 // Write to manual decomposition option
1357 {
1358 labelIOList cellToRegion
1359 (
1360 IOobject
1361 (
1362 "cellToRegion",
1363 mesh.facesInstance(),
1364 mesh,
1368 ),
1369 cellRegion
1370 );
1371 cellToRegion.write();
1372
1373 Info<< "Writing region per cell file (for manual decomposition) to "
1374 << cellToRegion.objectPath() << nl << endl;
1375 }
1376 // Write for postprocessing
1377 {
1378 volScalarField cellToRegion
1379 (
1380 IOobject
1381 (
1382 "cellToRegion",
1383 mesh.time().timeName(),
1384 mesh,
1388 ),
1389 mesh,
1392 );
1393 forAll(cellRegion, celli)
1394 {
1395 cellToRegion[celli] = cellRegion[celli];
1396 }
1397 cellToRegion.write();
1398
1399 Info<< "Writing region per cell as volScalarField to "
1400 << cellToRegion.objectPath() << nl << endl;
1401 }
1402}
1403
1404
1405
1406
1407int main(int argc, char *argv[])
1408{
1410 (
1411 "Split mesh into multiple regions (detected by walking across faces)"
1412 );
1413 #include "addRegionOption.H"
1414 #include "addOverwriteOption.H"
1416 (
1417 "cellZones",
1418 "Additionally split cellZones off into separate regions"
1419 );
1421 (
1422 "cellZonesOnly",
1423 "Use cellZones only to split mesh into regions; do not use walking"
1424 );
1426 (
1427 "cellZonesFileOnly",
1428 "file",
1429 "Like -cellZonesOnly, but use specified file"
1430 );
1432 (
1433 "combineZones",
1434 "lists of zones",
1435 "Combine zones in follow-on analysis"
1436 );
1438 (
1439 "addZones",
1440 "lists of zones",
1441 "Combine zones in follow-on analysis"
1442 );
1444 (
1445 "blockedFaces",
1446 "faceSet",
1447 "Specify additional region boundaries that walking does not cross"
1448 );
1450 (
1451 "makeCellZones",
1452 "Place cells into cellZones instead of splitting mesh"
1453 );
1455 (
1456 "largestOnly",
1457 "Only write largest region"
1458 );
1460 (
1461 "insidePoint",
1462 "point",
1463 "Only write region containing point"
1464 );
1466 (
1467 "detectOnly",
1468 "Do not write mesh"
1469 );
1471 (
1472 "sloppyCellZones",
1473 "Try to match heuristically regions to existing cell zones"
1474 );
1476 (
1477 "useFaceZones",
1478 "Use faceZones to patch inter-region faces instead of single patch"
1479 );
1481 (
1482 "prefixRegion",
1483 "Prefix region name to all patches, not just coupling patches"
1484 );
1485
1486 argList::noFunctionObjects(); // Never use function objects
1487
1488 #include "setRootCase.H"
1489 #include "createTime.H"
1490 #include "createNamedMesh.H"
1491
1492 const word oldInstance = mesh.pointsInstance();
1493
1494 word blockedFacesName;
1495 if (args.readIfPresent("blockedFaces", blockedFacesName))
1496 {
1497 Info<< "Reading blocked internal faces from faceSet "
1498 << blockedFacesName << nl << endl;
1499 }
1500
1501 const bool makeCellZones = args.found("makeCellZones");
1502 const bool largestOnly = args.found("largestOnly");
1503 const bool insidePoint = args.found("insidePoint");
1504 const bool useCellZones = args.found("cellZones");
1505 const bool useCellZonesOnly = args.found("cellZonesOnly");
1506 const bool useCellZonesFile = args.found("cellZonesFileOnly");
1507 const bool combineZones = args.found("combineZones");
1508 const bool addZones = args.found("addZones");
1509 const bool overwrite = args.found("overwrite");
1510 const bool detectOnly = args.found("detectOnly");
1511 const bool sloppyCellZones = args.found("sloppyCellZones");
1512 const bool useFaceZones = args.found("useFaceZones");
1513 const bool prefixRegion = args.found("prefixRegion");
1514
1515
1516 if
1517 (
1518 (useCellZonesOnly || useCellZonesFile)
1519 && (useCellZones || blockedFacesName.size())
1520 )
1521 {
1523 << "You cannot specify both -cellZonesOnly or -cellZonesFileOnly"
1524 << " (which specify complete split)"
1525 << " in combination with -blockedFaces or -cellZones"
1526 << " (which imply a split based on topology)"
1527 << exit(FatalError);
1528 }
1529
1530
1531 if (useFaceZones)
1532 {
1533 Info<< "Using current faceZones to divide inter-region interfaces"
1534 << " into multiple patches."
1535 << nl << endl;
1536 }
1537 else
1538 {
1539 Info<< "Creating single patch per inter-region interface."
1540 << nl << endl;
1541 }
1542
1543
1544
1545 if (insidePoint && largestOnly)
1546 {
1548 << "You cannot specify both -largestOnly"
1549 << " (keep region with most cells)"
1550 << " and -insidePoint (keep region containing point)"
1551 << exit(FatalError);
1552 }
1553
1554
1555 // Make sure cellZone names consistent across processors
1556 mesh.cellZones().checkParallelSync(true);
1557
1558 List<wordRes> zoneClusters;
1559 wordList zoneClusterNames;
1560 if (combineZones)
1561 {
1562 if (addZones)
1563 {
1565 << "Cannot specify both combineZones and addZones"
1566 << exit(FatalError);
1567 }
1568 zoneClusters = args.get<List<wordRes>>("combineZones");
1569 zoneClusterNames.setSize(zoneClusters.size());
1570 }
1571 else if (addZones)
1572 {
1573 zoneClusters = args.get<List<wordRes>>("addZones");
1574 zoneClusterNames.setSize(zoneClusters.size());
1575 forAll(zoneClusters, clusteri)
1576 {
1577 // Pop of front - is name
1578
1579 wordRes& wrs = zoneClusters[clusteri];
1580
1581 zoneClusterNames[clusteri] = wrs[0];
1582
1583 for (label i = 1; i < wrs.size(); i++)
1584 {
1585 wrs[i-1] = wrs[i];
1586 }
1587 wrs.setSize(wrs.size()-1);
1588 }
1589 }
1590
1591
1592 // Determine per cell the region it belongs to
1593 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1594
1595 // cellRegion is the labelList with the region per cell.
1596 labelList cellRegion;
1597 // Region to zone(s)
1598 labelListList regionToZones;
1599 // Name of region
1601 // Zone to region
1602 labelList zoneToRegion;
1603
1604 label nCellRegions = 0;
1605 if (useCellZonesOnly)
1606 {
1607 Info<< "Using current cellZones to split mesh into regions."
1608 << " This requires all"
1609 << " cells to be in one and only one cellZone." << nl << endl;
1610
1611 // Collect sets of zones into clusters. If no cluster is just an identity
1612 // list (cluster 0 is cellZone 0 etc.)
1613 wordList clusterNames;
1614 labelListList clusterToZones;
1615 labelList zoneToCluster;
1616 makeClusters
1617 (
1618 zoneClusters,
1619 zoneClusterNames,
1620 mesh.cellZones(),
1621 clusterNames,
1622 clusterToZones,
1623 zoneToCluster
1624 );
1625
1626 // Existing clusterID
1627 labelList clusterID(mesh.nCells(), -1);
1628 // Neighbour clusterID.
1629 labelList neiClusterID(mesh.nBoundaryFaces());
1630 getClusterID
1631 (
1632 mesh,
1633 mesh.cellZones(),
1634 clusterNames,
1635 clusterToZones,
1636 clusterID,
1637 neiClusterID
1638 );
1639
1640 label unzonedCelli = clusterID.find(-1);
1641 if (unzonedCelli != -1)
1642 {
1644 << "For the cellZonesOnly option all cells "
1645 << "have to be in a cellZone." << endl
1646 << "Cell " << unzonedCelli
1647 << " at" << mesh.cellCentres()[unzonedCelli]
1648 << " is not in a cellZone. There might be more unzoned cells."
1649 << exit(FatalError);
1650 }
1651 cellRegion = clusterID;
1652 nCellRegions = gMax(cellRegion)+1;
1653 zoneToRegion = zoneToCluster;
1654 regionToZones = clusterToZones;
1655 regionNames = clusterNames;
1656 }
1657 else if (useCellZonesFile)
1658 {
1659 const word zoneFile(args["cellZonesFileOnly"]);
1660 Info<< "Reading split from cellZones file " << zoneFile << endl
1661 << "This requires all"
1662 << " cells to be in one and only one cellZone." << nl << endl;
1663
1664 cellZoneMesh newCellZones
1665 (
1666 IOobject
1667 (
1668 zoneFile,
1669 mesh.facesInstance(),
1671 mesh,
1675 ),
1676 mesh
1677 );
1678
1679 wordList clusterNames;
1680 labelListList clusterToZones;
1681 labelList zoneToCluster;
1682 makeClusters
1683 (
1684 zoneClusters,
1685 zoneClusterNames,
1686 newCellZones,
1687 clusterNames,
1688 clusterToZones,
1689 zoneToCluster
1690 );
1691
1692
1693 // Existing clusterID
1694 labelList clusterID(mesh.nCells(), -1);
1695 // Neighbour clusterID.
1696 labelList neiClusterID(mesh.nBoundaryFaces());
1697 getClusterID
1698 (
1699 mesh,
1700 newCellZones,
1701 clusterNames,
1702 clusterToZones,
1703 clusterID,
1704 neiClusterID
1705 );
1706
1707
1708 label unzonedCelli = clusterID.find(-1);
1709 if (unzonedCelli != -1)
1710 {
1712 << "For the cellZonesFileOnly option all cells "
1713 << "have to be in a cellZone." << endl
1714 << "Cell " << unzonedCelli
1715 << " at" << mesh.cellCentres()[unzonedCelli]
1716 << " is not in a cellZone. There might be more unzoned cells."
1717 << exit(FatalError);
1718 }
1719 cellRegion = clusterID;
1720 nCellRegions = gMax(cellRegion)+1;
1721 zoneToRegion = zoneToCluster;
1722 regionToZones = clusterToZones;
1723 regionNames = clusterNames;
1724 }
1725 else
1726 {
1727 // Determine connected regions
1728 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1729
1730 // Mark additional faces that are blocked
1731 boolList blockedFace;
1732
1733 // Read from faceSet
1734 if (blockedFacesName.size())
1735 {
1736 faceSet blockedFaceSet(mesh, blockedFacesName);
1737 Info<< "Read "
1738 << returnReduce(blockedFaceSet.size(), sumOp<label>())
1739 << " blocked faces from set " << blockedFacesName << nl << endl;
1740
1741 blockedFace.setSize(mesh.nFaces(), false);
1742
1743 for (const label facei : blockedFaceSet)
1744 {
1745 blockedFace[facei] = true;
1746 }
1747 }
1748
1749 // Collect sets of zones into clusters. If no cluster is just an
1750 // identity list (cluster 0 is cellZone 0 etc.)
1751 wordList clusterNames;
1752 labelListList clusterToZones;
1753 labelList zoneToCluster;
1754 makeClusters
1755 (
1756 zoneClusters,
1757 zoneClusterNames,
1758 mesh.cellZones(),
1759 clusterNames,
1760 clusterToZones,
1761 zoneToCluster
1762 );
1763
1764 // Existing clusterID
1765 labelList clusterID(mesh.nCells(), -1);
1766 // Neighbour clusterID.
1767 labelList neiClusterID(mesh.nBoundaryFaces());
1768 getClusterID
1769 (
1770 mesh,
1771 mesh.cellZones(),
1772 clusterNames,
1773 clusterToZones,
1774 clusterID,
1775 neiClusterID
1776 );
1777
1778
1779 // Imply from differing cellZones
1780 if (useCellZones)
1781 {
1782 blockedFace.setSize(mesh.nFaces(), false);
1783
1784 for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1785 {
1786 label ownCluster = clusterID[mesh.faceOwner()[facei]];
1787 label neiCluster = clusterID[mesh.faceNeighbour()[facei]];
1788
1789 if (ownCluster != neiCluster)
1790 {
1791 blockedFace[facei] = true;
1792 }
1793 }
1794
1795 // Different cellZones on either side of processor patch.
1796 forAll(neiClusterID, i)
1797 {
1798 label facei = i+mesh.nInternalFaces();
1799 label ownCluster = clusterID[mesh.faceOwner()[facei]];
1800 label neiCluster = neiClusterID[i];
1801
1802 if (ownCluster != neiCluster)
1803 {
1804 blockedFace[facei] = true;
1805 }
1806 }
1807 }
1808
1809 // Do a topological walk to determine regions
1810 regionSplit regions(mesh, blockedFace);
1811 nCellRegions = regions.nRegions();
1812 cellRegion.transfer(regions);
1813
1814 // Match regions to zones
1815 matchRegions
1816 (
1817 sloppyCellZones,
1818 mesh,
1819
1820 // cluster-to-zone and cluster-to-cell addressing
1821 clusterNames,
1822 clusterToZones,
1823 clusterID,
1824
1825 // cell-to-region addressing
1826 nCellRegions,
1827 cellRegion,
1828
1829 // matched zones
1830 regionToZones,
1832 zoneToRegion
1833 );
1834
1835 // Override any default region names if single region selected
1836 if (largestOnly || insidePoint)
1837 {
1838 forAll(regionToZones, regionI)
1839 {
1840 if (regionToZones[regionI].empty())
1841 {
1842 if (overwrite)
1843 {
1845 }
1846 else if (insidePoint)
1847 {
1848 regionNames[regionI] = "insidePoint";
1849 }
1850 else if (largestOnly)
1851 {
1852 regionNames[regionI] = "largestOnly";
1853 }
1854 }
1855 }
1856 }
1857 }
1858
1859 Info<< endl << "Number of regions:" << nCellRegions << nl << endl;
1860
1861
1862 // Write decomposition to file
1863 writeCellToRegion(mesh, cellRegion);
1864
1865
1866
1867 // Sizes per region
1868 // ~~~~~~~~~~~~~~~~
1869
1870 labelList regionSizes(nCellRegions, Zero);
1871
1872 forAll(cellRegion, celli)
1873 {
1874 regionSizes[cellRegion[celli]]++;
1875 }
1876 forAll(regionSizes, regionI)
1877 {
1878 reduce(regionSizes[regionI], sumOp<label>());
1879 }
1880
1881 Info<< "Region\tCells" << nl
1882 << "------\t-----" << endl;
1883
1884 forAll(regionSizes, regionI)
1885 {
1886 Info<< regionI << '\t' << regionSizes[regionI] << nl;
1887 }
1888 Info<< endl;
1889
1890
1891
1892 // Print region to zone
1893 Info<< "Region\tZone\tName" << nl
1894 << "------\t----\t----" << endl;
1895 forAll(regionToZones, regionI)
1896 {
1897 Info<< regionI << '\t' << flatOutput(regionToZones[regionI])
1898 << '\t'
1899 << regionNames[regionI] << nl;
1900 }
1901 Info<< endl;
1902
1903
1904
1905 // Since we're going to mess with patches and zones make sure all
1906 // is synchronised
1907 mesh.boundaryMesh().checkParallelSync(true);
1908 mesh.faceZones().checkParallelSync(true);
1909
1910
1911 // Interfaces
1912 // ----------
1913 // per interface:
1914 // - the two regions on either side
1915 // - the name
1916 // - the (global) size
1917 edgeList interfaces;
1918 List<Pair<word>> interfaceNames;
1919 labelList interfaceSizes;
1920 // per face the interface
1921 labelList faceToInterface;
1922
1923 getInterfaceSizes
1924 (
1925 mesh,
1926 useFaceZones,
1927 cellRegion,
1929
1930 interfaces,
1931 interfaceNames,
1932 interfaceSizes,
1933 faceToInterface
1934 );
1935
1936 Info<< "Sizes of interfaces between regions:" << nl << nl
1937 << "Interface\tRegion\tRegion\tFaces" << nl
1938 << "---------\t------\t------\t-----" << endl;
1939
1940 forAll(interfaces, interI)
1941 {
1942 const edge& e = interfaces[interI];
1943
1944 Info<< interI
1945 << "\t\t" << e[0] << "\t" << e[1]
1946 << "\t" << interfaceSizes[interI] << nl;
1947 }
1948 Info<< endl;
1949
1950
1951 if (detectOnly)
1952 {
1953 Info<< "End\n" << endl;
1954
1955 return 0;
1956 }
1957
1958
1959 // Read objects in time directory
1960 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1961
1962 IOobjectList objects(mesh, runTime.timeName());
1963
1964 // Read vol fields.
1965
1967 ReadFields(mesh, objects, vsFlds);
1968
1970 ReadFields(mesh, objects, vvFlds);
1971
1973 ReadFields(mesh, objects, vstFlds);
1974
1976 ReadFields(mesh, objects, vsymtFlds);
1977
1979 ReadFields(mesh, objects, vtFlds);
1980
1981 // Read surface fields.
1982
1984 ReadFields(mesh, objects, ssFlds);
1985
1987 ReadFields(mesh, objects, svFlds);
1988
1990 ReadFields(mesh, objects, sstFlds);
1991
1993 ReadFields(mesh, objects, ssymtFlds);
1994
1996 ReadFields(mesh, objects, stFlds);
1997
1998 Info<< endl;
1999
2000
2001 // Remove any demand-driven fields ('S', 'V' etc)
2002 mesh.clearOut();
2003
2004
2005 if (nCellRegions == 1)
2006 {
2007 Info<< "Only one region. Doing nothing." << endl;
2008 }
2009 else if (makeCellZones)
2010 {
2011 Info<< "Putting cells into cellZones instead of splitting mesh."
2012 << endl;
2013
2014 // Check if region overlaps with existing zone. If so keep.
2015
2016 for (label regionI = 0; regionI < nCellRegions; regionI++)
2017 {
2018 const labelList& zones = regionToZones[regionI];
2019
2020 if (zones.size() == 1 && zones[0] != -1)
2021 {
2022 // Existing zone
2023 const label zoneI = zones[0];
2024 Info<< " Region " << regionI << " : corresponds to existing"
2025 << " cellZone "
2026 << zoneI << ' ' << mesh.cellZones()[zoneI].name() << endl;
2027 }
2028 else
2029 {
2030 // Create new cellZone.
2031 const labelList regionCells(findIndices(cellRegion, regionI));
2032
2033 const word zoneName = "region" + Foam::name(regionI);
2034
2035 label zoneI = mesh.cellZones().findZoneID(zoneName);
2036
2037 if (zoneI == -1)
2038 {
2039 zoneI = mesh.cellZones().size();
2040 mesh.cellZones().setSize(zoneI+1);
2041 mesh.cellZones().set
2042 (
2043 zoneI,
2044 new cellZone
2045 (
2046 zoneName, //name
2047 std::move(regionCells), //addressing
2048 zoneI, //index
2049 mesh.cellZones() //cellZoneMesh
2050 )
2051 );
2052 }
2053 else
2054 {
2055 mesh.cellZones()[zoneI].clearAddressing();
2056 mesh.cellZones()[zoneI] = regionCells;
2057 }
2058 Info<< " Region " << regionI << " : created new cellZone "
2059 << zoneI << ' ' << mesh.cellZones()[zoneI].name() << endl;
2060 }
2061 }
2062 mesh.cellZones().writeOpt(IOobject::AUTO_WRITE);
2063
2064 if (!overwrite)
2065 {
2066 ++runTime;
2067 mesh.setInstance(runTime.timeName());
2068 }
2069 else
2070 {
2071 mesh.setInstance(oldInstance);
2072 }
2073
2074 Info<< "Writing cellZones as new mesh to time " << runTime.timeName()
2075 << nl << endl;
2076
2077 mesh.write();
2078
2079
2080 // Write cellSets for convenience
2081 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2082
2083 Info<< "Writing cellSets corresponding to cellZones." << nl << endl;
2084
2085 for (const auto& cz : mesh.cellZones())
2086 {
2087 cellSet(mesh, cz.name(), cz).write();
2088 }
2089 }
2090 else
2091 {
2092 // Add patches for interfaces
2093 // ~~~~~~~~~~~~~~~~~~~~~~~~~~
2094
2095 // Add all possible patches. Empty ones get filtered later on.
2096 Info<< nl << "Adding patches" << nl << endl;
2097
2098 labelList interfacePatches
2099 (
2100 addRegionPatches
2101 (
2102 mesh,
2104 interfaces,
2105 interfaceNames
2106 )
2107 );
2108
2109
2110 if (!overwrite)
2111 {
2112 ++runTime;
2113 }
2114
2115
2116 // Create regions
2117 // ~~~~~~~~~~~~~~
2118
2119 if (insidePoint)
2120 {
2121 const point insidePoint = args.get<point>("insidePoint");
2122
2123 label regionI = -1;
2124
2125 (void)mesh.tetBasePtIs();
2126
2127 label celli = mesh.findCell(insidePoint);
2128
2129 Info<< nl << "Found point " << insidePoint << " in cell " << celli
2130 << endl;
2131
2132 if (celli != -1)
2133 {
2134 regionI = cellRegion[celli];
2135 }
2136
2137 reduce(regionI, maxOp<label>());
2138
2139 Info<< nl
2140 << "Subsetting region " << regionI
2141 << " containing point " << insidePoint << endl;
2142
2143 if (regionI == -1)
2144 {
2146 << "Point " << insidePoint
2147 << " is not inside the mesh." << nl
2148 << "Bounding box of the mesh:" << mesh.bounds()
2149 << exit(FatalError);
2150 }
2151
2152 createAndWriteRegion
2153 (
2154 mesh,
2155 cellRegion,
2157 prefixRegion,
2158 faceToInterface,
2159 interfacePatches,
2160 regionI,
2161 (overwrite ? oldInstance : runTime.timeName())
2162 );
2163 }
2164 else if (largestOnly)
2165 {
2166 label regionI = findMax(regionSizes);
2167
2168 Info<< nl
2169 << "Subsetting region " << regionI
2170 << " of size " << regionSizes[regionI]
2171 << " as named region " << regionNames[regionI] << endl;
2172
2173 createAndWriteRegion
2174 (
2175 mesh,
2176 cellRegion,
2178 prefixRegion,
2179 faceToInterface,
2180 interfacePatches,
2181 regionI,
2182 (overwrite ? oldInstance : runTime.timeName())
2183 );
2184 }
2185 else
2186 {
2187 // Split all
2188 for (label regionI = 0; regionI < nCellRegions; regionI++)
2189 {
2190 Info<< nl
2191 << "Region " << regionI << nl
2192 << "-------- " << endl;
2193
2194 createAndWriteRegion
2195 (
2196 mesh,
2197 cellRegion,
2199 prefixRegion,
2200 faceToInterface,
2201 interfacePatches,
2202 regionI,
2203 (overwrite ? oldInstance : runTime.timeName())
2204 );
2205 }
2206 }
2207 }
2208
2209 Info<< "End\n" << endl;
2210
2211 return 0;
2212}
2213
2214
2215// ************************************************************************* //
Field reading functions for post-processing utilities.
Required Classes.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:59
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition HashTableI.H:86
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
static void recv(Type &value, const int fromProcNo, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, IOstreamOption::streamFormat fmt=IOstreamOption::BINARY)
Receive and deserialize a value. Uses operator>> for de-serialization.
Definition IPstream.H:80
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void setSize(label n)
Alias for resize().
Definition List.H:536
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
static void send(const Type &value, const UPstream::commsTypes commsType, const int toProcNo, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, IOstreamOption::streamFormat fmt=IOstreamOption::BINARY)
Serialize a value and send (buffered/blocking or standard mode). Uses operator<< for serialization.
Definition OPstream.H:80
static void broadcast(Type &value, const int communicator=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
static void broadcasts(const int communicator, Type &value, Args &&... values)
Broadcast multiple items to all communicator ranks. Does nothing in non-parallel.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label find(const T &val) const
Find index of the first occurrence of the value.
Definition UList.C:160
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition UPstream.H:1691
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition UPstream.H:1069
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition UPstream.H:1866
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static const Form zero
labelList indices(const wordRe &matcher, const bool useGroups=true) const
The (sorted) patch indices for all matches, optionally matching zone groups.
Definition ZoneMesh.C:562
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition argList.C:562
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition argList.C:389
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition argList.C:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A collection of cell labels.
Definition cellSet.H:50
A subset of mesh cells.
Definition cellZone.H:61
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
static edge sorted(label from, label to)
Create (in ascending order) from two vertex labels.
Definition edgeI.H:32
A list of face labels.
Definition faceSet.H:50
static tmp< DimensionedField< Type, volMesh > > interpolate(const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelUList &cellMap)
Map volume internal (dimensioned) field.
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches.
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Create additional fvSchemes/fvSolution files.
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add patch. Inserts patch before all processor patches.
Definition fvMeshTools.C:38
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & calculatedType() noexcept
The type name for calculated patch fields.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
@ NEARESTPATCHFACE
nearest face on selected patch
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order.
void updateMesh()
Correct polyBoundaryMesh after topology update.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Definition polyMesh.H:411
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Direct mesh changes based on v1.3 polyTopoChange syntax.
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Given list of cells to remove, insert all the topology changes.
Definition removeCells.H:60
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData, const bool parRun=UPstream::parRun())
Extract and swap to obtain neighbour cell values for all boundary faces.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled boundary face values. Uses eqOp.
Definition syncTools.H:524
A class for managing temporary objects.
Definition tmp.H:75
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition topoSet.C:693
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
label index() const noexcept
The index of this zone in the zone list.
const word & name() const noexcept
The zone name.
dynamicFvMesh & mesh
engineTime & runTime
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
Required Classes.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
auto & names
wordList regionNames
const labelIOList & zoneIDs
Definition correctPhi.H:59
return returnReduce(nRefine-oldNRefine, sumOp< label >())
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< word > wordList
List of word.
Definition fileName.H:60
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
const dimensionSet dimless
Dimensionless.
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
GeometricField< scalar, fvPatchField, volMesh > volScalarField
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
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).
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label findMax(const ListType &input, label start=0)
Linear search for the index of the max element, similar to std::max_element but for lists and returns...
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Foam::argList args(argc, argv)
volScalarField & e
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235