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