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-2023 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  Map<label> zoneToSize;
299  zoneToSize.insert(zoneID, 1);
300  regionsToSize.insert(interface, zoneToSize);
301  }
302 }
303 
304 
305 // Get region-region interface name and sizes.
306 // Returns interfaces as straight list for looping in identical order.
307 void getInterfaceSizes
308 (
309  const polyMesh& mesh,
310  const bool useFaceZones,
311  const labelList& cellRegion,
312  const wordList& regionNames,
313 
314  edgeList& interfaces,
315  List<Pair<word>>& interfaceNames,
316  labelList& interfaceSizes,
317  labelList& faceToInterface
318 )
319 {
320  // From region-region to faceZone (or -1) to number of faces.
321 
322  EdgeMap<Map<label>> regionsToSize;
323 
324 
325  // Internal faces
326  // ~~~~~~~~~~~~~~
327 
328  forAll(mesh.faceNeighbour(), facei)
329  {
330  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
331  label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
332 
333  if (ownRegion != neiRegion)
334  {
335  addToInterface
336  (
337  mesh,
338  (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
339  ownRegion,
340  neiRegion,
341  regionsToSize
342  );
343  }
344  }
345 
346  // Boundary faces
347  // ~~~~~~~~~~~~~~
348 
349  // Neighbour cellRegion.
350  labelList coupledRegion(mesh.nBoundaryFaces());
351 
352  forAll(coupledRegion, i)
353  {
354  label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
355  coupledRegion[i] = cellRegion[celli];
356  }
357  syncTools::swapBoundaryFaceList(mesh, coupledRegion);
358 
359  forAll(coupledRegion, i)
360  {
361  label facei = i+mesh.nInternalFaces();
362  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
363  label neiRegion = coupledRegion[i];
364 
365  if (ownRegion != neiRegion)
366  {
367  addToInterface
368  (
369  mesh,
370  (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
371  ownRegion,
372  neiRegion,
373  regionsToSize
374  );
375  }
376  }
377 
378 
379  if (Pstream::parRun())
380  {
381  if (Pstream::master())
382  {
383  // Receive and add to my sizes
384  for (const int slave : Pstream::subProcs())
385  {
386  IPstream fromSlave(Pstream::commsTypes::blocking, slave);
387 
388  EdgeMap<Map<label>> slaveSizes(fromSlave);
389 
390  forAllConstIters(slaveSizes, slaveIter)
391  {
392  const Map<label>& slaveInfo = *slaveIter;
393 
394  auto masterIter = regionsToSize.find(slaveIter.key());
395 
396  if (masterIter.good())
397  {
398  // Same inter-region
399  Map<label>& masterInfo = *masterIter;
400 
401  forAllConstIters(slaveInfo, iter)
402  {
403  const label zoneID = iter.key();
404  const label slaveSize = iter.val();
405 
406  auto zoneIter = masterInfo.find(zoneID);
407  if (zoneIter.good())
408  {
409  *zoneIter += slaveSize;
410  }
411  else
412  {
413  masterInfo.insert(zoneID, slaveSize);
414  }
415  }
416  }
417  else
418  {
419  regionsToSize.insert(slaveIter.key(), slaveInfo);
420  }
421  }
422  }
423  }
424  else
425  {
426  // Send to master
427  {
428  OPstream toMaster
429  (
432  );
433  toMaster << regionsToSize;
434  }
435  }
436  }
437 
438  // Rework
439 
440  Pstream::broadcast(regionsToSize);
441 
442 
443  // Now we have the global sizes of all inter-regions.
444  // Invert this on master and distribute.
445  label nInterfaces = 0;
446  forAllConstIters(regionsToSize, iter)
447  {
448  const Map<label>& info = iter.val();
449  nInterfaces += info.size();
450  }
451 
452  interfaces.setSize(nInterfaces);
453  interfaceNames.setSize(nInterfaces);
454  interfaceSizes.setSize(nInterfaces);
455  EdgeMap<Map<label>> regionsToInterface(nInterfaces);
456 
457  nInterfaces = 0;
458  forAllConstIters(regionsToSize, iter)
459  {
460  const edge& e = iter.key();
461  const Map<label>& info = iter.val();
462 
463  const word& name0 = regionNames[e[0]];
464  const word& name1 = regionNames[e[1]];
465 
466  forAllConstIters(info, infoIter)
467  {
468  interfaces[nInterfaces] = iter.key();
469  label zoneID = infoIter.key();
470  if (zoneID == -1)
471  {
472  interfaceNames[nInterfaces] = Pair<word>
473  (
474  name0 + "_to_" + name1,
475  name1 + "_to_" + name0
476  );
477  }
478  else
479  {
480  const word& zoneName = mesh.faceZones()[zoneID].name();
481  interfaceNames[nInterfaces] = Pair<word>
482  (
483  zoneName + "_" + name0 + "_to_" + name1,
484  zoneName + "_" + name1 + "_to_" + name0
485  );
486  }
487  interfaceSizes[nInterfaces] = infoIter();
488 
489  if (regionsToInterface.found(e))
490  {
491  regionsToInterface[e].insert(zoneID, nInterfaces);
492  }
493  else
494  {
495  Map<label> zoneAndInterface;
496  zoneAndInterface.insert(zoneID, nInterfaces);
497  regionsToInterface.insert(e, zoneAndInterface);
498  }
499  nInterfaces++;
500  }
501  }
502 
503 
504  // Consistent interface information for all processors
506  (
508  interfaces,
509  interfaceNames,
510  interfaceSizes,
511  regionsToInterface
512  );
513 
514  // Mark all inter-region faces.
515  faceToInterface.setSize(mesh.nFaces(), -1);
516 
517  forAll(mesh.faceNeighbour(), facei)
518  {
519  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
520  label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
521 
522  if (ownRegion != neiRegion)
523  {
524  label zoneID = -1;
525  if (useFaceZones)
526  {
527  zoneID = mesh.faceZones().whichZone(facei);
528  }
529 
531  (
532  min(ownRegion, neiRegion),
533  max(ownRegion, neiRegion)
534  );
535 
536  faceToInterface[facei] = regionsToInterface[interface][zoneID];
537  }
538  }
539  forAll(coupledRegion, i)
540  {
541  label facei = i+mesh.nInternalFaces();
542  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
543  label neiRegion = coupledRegion[i];
544 
545  if (ownRegion != neiRegion)
546  {
547  label zoneID = -1;
548  if (useFaceZones)
549  {
550  zoneID = mesh.faceZones().whichZone(facei);
551  }
552 
554  (
555  min(ownRegion, neiRegion),
556  max(ownRegion, neiRegion)
557  );
558 
559  faceToInterface[facei] = regionsToInterface[interface][zoneID];
560  }
561  }
562 }
563 
564 
565 // Create mesh for region.
566 autoPtr<mapPolyMesh> createRegionMesh
567 (
568  const fvMesh& mesh,
569  // Region info
570  const labelList& cellRegion,
571  const label regionI,
572  const word& regionName,
573  // Interface info
574  const labelList& interfacePatches,
575  const labelList& faceToInterface,
576 
577  autoPtr<fvMesh>& newMesh
578 )
579 {
580  // Create dummy system/fv*
582 
583  // Neighbour cellRegion.
584  labelList coupledRegion(mesh.nBoundaryFaces());
585 
586  forAll(coupledRegion, i)
587  {
588  label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
589  coupledRegion[i] = cellRegion[celli];
590  }
591  syncTools::swapBoundaryFaceList(mesh, coupledRegion);
592 
593 
594  // Topology change container. Start off from existing mesh.
595  polyTopoChange meshMod(mesh);
596 
597  // Cell remover engine
598  removeCells cellRemover(mesh);
599 
600  // Select all but region cells
601  labelList cellsToRemove(getNonRegionCells(cellRegion, regionI));
602 
603  // Find out which faces will get exposed. Note that this
604  // gets faces in mesh face order. So both regions will get same
605  // face in same order (important!)
606  labelList exposedFaces = cellRemover.getExposedFaces(cellsToRemove);
607 
608  labelList exposedPatchIDs(exposedFaces.size());
609  forAll(exposedFaces, i)
610  {
611  label facei = exposedFaces[i];
612  label interfacei = faceToInterface[facei];
613 
614  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
615  label neiRegion = -1;
616 
617  if (mesh.isInternalFace(facei))
618  {
619  neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
620  }
621  else
622  {
623  neiRegion = coupledRegion[facei-mesh.nInternalFaces()];
624  }
625 
626 
627  // Check which side is being kept - determines which of the two
628  // patches will be used.
629 
630  label otherRegion = -1;
631 
632  if (ownRegion == regionI && neiRegion != regionI)
633  {
634  otherRegion = neiRegion;
635  }
636  else if (ownRegion != regionI && neiRegion == regionI)
637  {
638  otherRegion = ownRegion;
639  }
640  else
641  {
643  << "Exposed face:" << facei
644  << " fc:" << mesh.faceCentres()[facei]
645  << " has owner region " << ownRegion
646  << " and neighbour region " << neiRegion
647  << " when handling region:" << regionI
648  << exit(FatalError);
649  }
650 
651  // Find the patch.
652  if (regionI < otherRegion)
653  {
654  exposedPatchIDs[i] = interfacePatches[interfacei];
655  }
656  else
657  {
658  exposedPatchIDs[i] = interfacePatches[interfacei]+1;
659  }
660  }
661 
662  // Remove faces
663  cellRemover.setRefinement
664  (
665  cellsToRemove,
666  exposedFaces,
667  exposedPatchIDs,
668  meshMod
669  );
670 
671  autoPtr<mapPolyMesh> map = meshMod.makeMesh
672  (
673  newMesh,
674  IOobject
675  (
676  regionName,
677  mesh.time().timeName(),
678  mesh.time(),
679  IOobject::READ_IF_PRESENT, // read fv* if present
681  ),
682  mesh
683  );
684 
685  return map;
686 }
687 
688 
689 void createAndWriteRegion
690 (
691  const fvMesh& mesh,
692  const labelList& cellRegion,
693  const wordList& regionNames,
694  const bool prefixRegion,
695  const labelList& faceToInterface,
696  const labelList& interfacePatches,
697  const label regionI,
698  const word& newMeshInstance
699 )
700 {
701  Info<< "Creating mesh for region " << regionI
702  << ' ' << regionNames[regionI] << endl;
703 
704  autoPtr<fvMesh> newMesh;
705  autoPtr<mapPolyMesh> map = createRegionMesh
706  (
707  mesh,
708  cellRegion,
709  regionI,
710  regionNames[regionI],
711  interfacePatches,
712  faceToInterface,
713  newMesh
714  );
715 
716 
717  // Make map of all added patches
718  labelHashSet addedPatches(2*interfacePatches.size());
719  forAll(interfacePatches, interfacei)
720  {
721  addedPatches.insert(interfacePatches[interfacei]);
722  addedPatches.insert(interfacePatches[interfacei]+1);
723  }
724 
725 
726  Info<< "Mapping fields" << endl;
727 
728  // Map existing fields
729  newMesh().updateMesh(map());
730 
731  // Add subsetted fields
732  subsetVolFields<volScalarField>
733  (
734  mesh,
735  newMesh(),
736  map().cellMap(),
737  map().faceMap(),
738  addedPatches
739  );
740  subsetVolFields<volVectorField>
741  (
742  mesh,
743  newMesh(),
744  map().cellMap(),
745  map().faceMap(),
746  addedPatches
747  );
748  subsetVolFields<volSphericalTensorField>
749  (
750  mesh,
751  newMesh(),
752  map().cellMap(),
753  map().faceMap(),
754  addedPatches
755  );
756  subsetVolFields<volSymmTensorField>
757  (
758  mesh,
759  newMesh(),
760  map().cellMap(),
761  map().faceMap(),
762  addedPatches
763  );
764  subsetVolFields<volTensorField>
765  (
766  mesh,
767  newMesh(),
768  map().cellMap(),
769  map().faceMap(),
770  addedPatches
771  );
772 
773  subsetSurfaceFields<surfaceScalarField>
774  (
775  mesh,
776  newMesh(),
777  map().cellMap(),
778  map().faceMap(),
779  addedPatches
780  );
781  subsetSurfaceFields<surfaceVectorField>
782  (
783  mesh,
784  newMesh(),
785  map().cellMap(),
786  map().faceMap(),
787  addedPatches
788  );
789  subsetSurfaceFields<surfaceSphericalTensorField>
790  (
791  mesh,
792  newMesh(),
793  map().cellMap(),
794  map().faceMap(),
795  addedPatches
796  );
797  subsetSurfaceFields<surfaceSymmTensorField>
798  (
799  mesh,
800  newMesh(),
801  map().cellMap(),
802  map().faceMap(),
803  addedPatches
804  );
805  subsetSurfaceFields<surfaceTensorField>
806  (
807  mesh,
808  newMesh(),
809  map().cellMap(),
810  map().faceMap(),
811  addedPatches
812  );
813 
814 
815  const polyBoundaryMesh& newPatches = newMesh().boundaryMesh();
816  newPatches.checkParallelSync(true);
817 
818  // Delete empty patches
819  // ~~~~~~~~~~~~~~~~~~~~
820 
821  // Create reordering list to move patches-to-be-deleted to end
822  labelList oldToNew(newPatches.size(), -1);
823  DynamicList<label> sharedPatches(newPatches.size());
824  label newI = 0;
825 
826  Info<< "Deleting empty patches" << endl;
827 
828  // Assumes all non-proc boundaries are on all processors!
829  forAll(newPatches, patchi)
830  {
831  const polyPatch& pp = newPatches[patchi];
832 
833  if (!isA<processorPolyPatch>(pp))
834  {
835  if (returnReduceOr(pp.size()))
836  {
837  oldToNew[patchi] = newI;
838  if (!addedPatches.found(patchi))
839  {
840  sharedPatches.append(newI);
841  }
842  newI++;
843  }
844  }
845  }
846 
847  // Same for processor patches (but need no reduction)
848  forAll(newPatches, patchi)
849  {
850  const polyPatch& pp = newPatches[patchi];
851 
852  if (isA<processorPolyPatch>(pp) && pp.size())
853  {
854  oldToNew[patchi] = newI++;
855  }
856  }
857 
858  const label nNewPatches = newI;
859 
860  // Move all deleteable patches to the end
861  forAll(oldToNew, patchi)
862  {
863  if (oldToNew[patchi] == -1)
864  {
865  oldToNew[patchi] = newI++;
866  }
867  }
868 
869  //reorderPatches(newMesh(), oldToNew, nNewPatches);
870  fvMeshTools::reorderPatches(newMesh(), oldToNew, nNewPatches, true);
871 
872  // Rename shared patches with region name
873  if (prefixRegion)
874  {
875  Info<< "Prefixing patches with region name" << endl;
876 
877  renamePatches(newMesh(), regionNames[regionI], sharedPatches);
878  }
879 
880 
881  Info<< "Writing new mesh" << endl;
882 
883  newMesh().setInstance(newMeshInstance);
884  newMesh().write();
885  topoSet::removeFiles(newMesh());
886  processorMeshes::removeFiles(newMesh());
887 
888  // Write addressing files like decomposePar
889  Info<< "Writing addressing to base mesh" << endl;
890 
891  labelIOList pointProcAddressing
892  (
893  IOobject
894  (
895  "pointRegionAddressing",
896  newMesh().facesInstance(),
898  newMesh(),
902  ),
903  map().pointMap()
904  );
905  Info<< "Writing map " << pointProcAddressing.name()
906  << " from region" << regionI
907  << " points back to base mesh." << endl;
908  pointProcAddressing.write();
909 
910  labelIOList faceProcAddressing
911  (
912  IOobject
913  (
914  "faceRegionAddressing",
915  newMesh().facesInstance(),
917  newMesh(),
921  ),
922  newMesh().nFaces()
923  );
924  forAll(faceProcAddressing, facei)
925  {
926  // face + turning index. (see decomposePar)
927  // Is the face pointing in the same direction?
928  label oldFacei = map().faceMap()[facei];
929 
930  if
931  (
932  map().cellMap()[newMesh().faceOwner()[facei]]
933  == mesh.faceOwner()[oldFacei]
934  )
935  {
936  faceProcAddressing[facei] = oldFacei+1;
937  }
938  else
939  {
940  faceProcAddressing[facei] = -(oldFacei+1);
941  }
942  }
943  Info<< "Writing map " << faceProcAddressing.name()
944  << " from region" << regionI
945  << " faces back to base mesh." << endl;
946  faceProcAddressing.write();
947 
948  labelIOList cellProcAddressing
949  (
950  IOobject
951  (
952  "cellRegionAddressing",
953  newMesh().facesInstance(),
955  newMesh(),
959  ),
960  map().cellMap()
961  );
962  Info<< "Writing map " <<cellProcAddressing.name()
963  << " from region" << regionI
964  << " cells back to base mesh." << endl;
965  cellProcAddressing.write();
966 
967  labelIOList boundaryProcAddressing
968  (
969  IOobject
970  (
971  "boundaryRegionAddressing",
972  newMesh().facesInstance(),
974  newMesh(),
978  ),
979  labelList(nNewPatches, -1)
980  );
981  forAll(oldToNew, i)
982  {
983  if (!addedPatches.found(i))
984  {
985  label newI = oldToNew[i];
986  if (newI >= 0 && newI < nNewPatches)
987  {
988  boundaryProcAddressing[oldToNew[i]] = i;
989  }
990  }
991  }
992  Info<< "Writing map " << boundaryProcAddressing.name()
993  << " from region" << regionI
994  << " boundary back to base mesh." << endl;
995  boundaryProcAddressing.write();
996 }
997 
998 
999 // Create for every region-region interface with non-zero size two patches.
1000 // First one is for minimumregion to maximumregion.
1001 // Note that patches get created in same order on all processors (if parallel)
1002 // since looping over synchronised 'interfaces'.
1003 labelList addRegionPatches
1004 (
1005  fvMesh& mesh,
1006  const wordList& regionNames,
1007  const edgeList& interfaces,
1008  const List<Pair<word>>& interfaceNames
1009 )
1010 {
1011  Info<< nl << "Adding patches" << nl << endl;
1012 
1013  labelList interfacePatches(interfaces.size());
1014 
1015  forAll(interfaces, interI)
1016  {
1017  const edge& e = interfaces[interI];
1018  const Pair<word>& names = interfaceNames[interI];
1019 
1020  //Info<< "For interface " << interI
1021  // << " between regions " << e
1022  // << " trying to add patches " << names << endl;
1023 
1024 
1025  mappedWallPolyPatch patch1
1026  (
1027  names[0],
1028  0, // overridden
1029  0, // overridden
1030  0, // overridden
1031  regionNames[e[1]], // sampleRegion
1033  names[1], // samplePatch
1034  point::zero, // offset
1035  mesh.boundaryMesh()
1036  );
1037 
1038  interfacePatches[interI] = fvMeshTools::addPatch
1039  (
1040  mesh,
1041  patch1,
1042  dictionary(), //optional per field value
1044  true //validBoundary
1045  );
1046 
1047  mappedWallPolyPatch patch2
1048  (
1049  names[1],
1050  0,
1051  0,
1052  0,
1053  regionNames[e[0]], // sampleRegion
1055  names[0],
1056  point::zero, // offset
1057  mesh.boundaryMesh()
1058  );
1060  (
1061  mesh,
1062  patch2,
1063  dictionary(), //optional per field value
1065  true //validBoundary
1066  );
1067 
1068  Info<< "For interface between region " << regionNames[e[0]]
1069  << " and " << regionNames[e[1]] << " added patches" << endl
1070  << " " << interfacePatches[interI]
1071  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]].name()
1072  << endl
1073  << " " << interfacePatches[interI]+1
1074  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]+1].name()
1075  << endl;
1076  }
1077  return interfacePatches;
1078 }
1079 
1080 
1081 // Find region that covers most of cell zone
1082 label findCorrespondingRegion
1083 (
1084  const labelList& existingZoneID, // per cell the (unique) zoneID
1085  const labelList& cellRegion,
1086  const label nCellRegions,
1087  const label zoneI,
1088  const label minOverlapSize
1089 )
1090 {
1091  // Per region the number of cells in zoneI
1092  labelList cellsInZone(nCellRegions, Zero);
1093 
1094  forAll(cellRegion, celli)
1095  {
1096  if (existingZoneID[celli] == zoneI)
1097  {
1098  cellsInZone[cellRegion[celli]]++;
1099  }
1100  }
1101 
1103 
1104  // Pick region with largest overlap of zoneI
1105  label regionI = findMax(cellsInZone);
1106 
1107 
1108  if (cellsInZone[regionI] < minOverlapSize)
1109  {
1110  // Region covers too little of zone. Not good enough.
1111  regionI = -1;
1112  }
1113  else
1114  {
1115  // Check that region contains no cells that aren't in cellZone.
1116  forAll(cellRegion, celli)
1117  {
1118  if (cellRegion[celli] == regionI && existingZoneID[celli] != zoneI)
1119  {
1120  // celli in regionI but not in zoneI
1121  regionI = -1;
1122  break;
1123  }
1124  }
1125  // If one in error, all should be in error. Note that branch gets taken
1126  // on all procs.
1127  reduce(regionI, minOp<label>());
1128  }
1129 
1130  return regionI;
1131 }
1132 
1133 
1134 void getClusterID
1135 (
1136  const polyMesh& mesh,
1137  const cellZoneMesh& cellZones,
1138  const wordList& clusterNames,
1139  const labelListList& clusterToZones,
1140  labelList& clusterID,
1141  labelList& neiClusterID
1142 )
1143 {
1144  // Existing zoneID
1145  clusterID.setSize(mesh.nCells());
1146  clusterID = -1;
1147 
1148  forAll(clusterToZones, clusterI)
1149  {
1150  for (const label zoneI : clusterToZones[clusterI])
1151  {
1152  const cellZone& cz = cellZones[zoneI];
1153 
1154  forAll(cz, i)
1155  {
1156  label celli = cz[i];
1157  if (clusterID[celli] == -1)
1158  {
1159  clusterID[celli] = clusterI;
1160  }
1161  else
1162  {
1164  << "Cell " << celli << " with cell centre "
1165  << mesh.cellCentres()[celli]
1166  << " is multiple zones. This is not allowed." << endl
1167  << "It is in zone " << clusterNames[clusterID[celli]]
1168  << " and in zone " << clusterNames[clusterI]
1169  << exit(FatalError);
1170  }
1171  }
1172  }
1173  }
1174 
1175  // Neighbour zoneID.
1176  syncTools::swapBoundaryCellList(mesh, clusterID, neiClusterID);
1177 }
1178 
1179 
1180 word makeRegionName
1181 (
1182  const cellZoneMesh& czs,
1183  const label regioni,
1184  const labelList& zoneIDs
1185 )
1186 {
1187  // Synthesise region name. Equals the zone name if cluster consist of only
1188  // one zone
1189 
1190  if (zoneIDs.empty())
1191  {
1192  return word("domain") + Foam::name(regioni);
1193  }
1194  else
1195  {
1196  // Zone indices are in cellZone order ...
1197  word regionName(czs[zoneIDs[0]].name());
1198 
1199  // Synthesize name from appended cellZone names
1200  for (label i = 1; i < zoneIDs.size(); i++)
1201  {
1202  regionName += "_" + czs[zoneIDs[i]].name();
1203  }
1204  return regionName;
1205  }
1206 }
1207 
1208 
1209 void makeClusters
1210 (
1211  const List<wordRes>& zoneClusters,
1212  const wordList& zoneClusterNames,
1213  const cellZoneMesh& cellZones,
1214  wordList& clusterNames,
1215  labelListList& clusterToZones,
1216  labelList& zoneToCluster
1217 )
1218 {
1219  // Check if there are clustering for zones. If none every zone goes into
1220  // its own cluster.
1221 
1222  clusterNames.clear();
1223  clusterToZones.clear();
1224  zoneToCluster.setSize(cellZones.size());
1225  zoneToCluster = -1;
1226 
1227  if (zoneClusters.size())
1228  {
1229  forAll(zoneClusters, clusteri)
1230  {
1231  const labelList zoneIDs(cellZones.indices(zoneClusters[clusteri]));
1232  UIndirectList<label>(zoneToCluster, zoneIDs) = clusteri;
1233  clusterNames.append
1234  (
1235  zoneClusterNames[clusteri].size()
1236  ? zoneClusterNames[clusteri]
1237  : makeRegionName
1238  (
1239  cellZones,
1240  clusteri,
1241  zoneIDs
1242  )
1243  );
1244  clusterToZones.append(std::move(zoneIDs));
1245  }
1246 
1247  // Unclustered zone
1248  forAll(zoneToCluster, zonei)
1249  {
1250  if (zoneToCluster[zonei] == -1)
1251  {
1252  clusterNames.append(cellZones[zonei].name());
1253  clusterToZones.append(labelList(1, zonei));
1254  zoneToCluster[zonei] = clusterToZones.size();
1255  }
1256  }
1257  }
1258  else
1259  {
1260  for (const auto& cellZone : cellZones)
1261  {
1262  const label nClusters = clusterToZones.size();
1263  clusterNames.append(cellZone.name());
1264  clusterToZones.append(labelList(1, cellZone.index()));
1265  zoneToCluster[cellZone.index()] = nClusters;
1266  }
1267  }
1268 }
1269 
1270 
1271 void matchRegions
1272 (
1273  const bool sloppyCellZones,
1274  const polyMesh& mesh,
1275 
1276  const wordList& clusterNames,
1277  const labelListList& clusterToZones,
1278  const labelList& clusterID,
1279 
1280  const label nCellRegions,
1281  const labelList& cellRegion,
1282 
1283  labelListList& regionToZones,
1285  labelList& zoneToRegion
1286 )
1287 {
1288  const cellZoneMesh& cellZones = mesh.cellZones();
1289 
1290  regionToZones.setSize(nCellRegions);
1291  regionToZones = labelList();
1292  regionNames.setSize(nCellRegions);
1294  zoneToRegion.setSize(cellZones.size(), -1);
1295 
1296 
1297  // Sizes per cluster
1298  labelList clusterSizes(clusterToZones.size(), Zero);
1299  forAll(clusterToZones, clusterI)
1300  {
1301  for (const label zoneI : clusterToZones[clusterI])
1302  {
1303  clusterSizes[clusterI] += cellZones[zoneI].size();
1304  }
1305  reduce(clusterSizes[clusterI], sumOp<label>());
1306  }
1307 
1308  if (sloppyCellZones)
1309  {
1310  Info<< "Trying to match regions to existing cell zones;"
1311  << " region can be subset of cell zone." << nl << endl;
1312 
1313  forAll(clusterToZones, clusterI)
1314  {
1315  label regionI = findCorrespondingRegion
1316  (
1317  clusterID,
1318  cellRegion,
1319  nCellRegions,
1320  clusterI,
1321  label(0.5*clusterSizes[clusterI]) // minimum overlap
1322  );
1323 
1324  if (regionI != -1)
1325  {
1326  Info<< "Sloppily matched region " << regionI
1327  //<< " size " << regionSizes[regionI]
1328  << " to cluster " << clusterI
1329  << " size " << clusterSizes[clusterI]
1330  << endl;
1332  (
1333  zoneToRegion,
1334  clusterToZones[clusterI]
1335  ) = regionI;
1336  regionToZones[regionI] = clusterToZones[clusterI];
1337  regionNames[regionI] = clusterNames[clusterI];
1338  }
1339  }
1340  }
1341  else
1342  {
1343  Info<< "Trying to match regions to existing cell zones." << nl << endl;
1344 
1345  forAll(clusterToZones, clusterI)
1346  {
1347  label regionI = findCorrespondingRegion
1348  (
1349  clusterID,
1350  cellRegion,
1351  nCellRegions,
1352  clusterI,
1353  clusterSizes[clusterI] // require exact match
1354  );
1355 
1356  if (regionI != -1)
1357  {
1359  (
1360  zoneToRegion,
1361  clusterToZones[clusterI]
1362  ) = regionI;
1363  regionToZones[regionI] = clusterToZones[clusterI];
1364  regionNames[regionI] = clusterNames[clusterI];
1365  }
1366  }
1367  }
1368  // Allocate region names for unmatched regions.
1369  forAll(regionNames, regionI)
1370  {
1371  if (regionNames[regionI].empty())
1372  {
1373  regionNames[regionI] = makeRegionName
1374  (
1375  cellZones,
1376  regionI,
1377  regionToZones[regionI]
1378  );
1379  }
1380  }
1381 }
1382 
1383 
1384 void writeCellToRegion(const fvMesh& mesh, const labelList& cellRegion)
1385 {
1386  // Write to manual decomposition option
1387  {
1388  labelIOList cellToRegion
1389  (
1390  IOobject
1391  (
1392  "cellToRegion",
1393  mesh.facesInstance(),
1394  mesh,
1398  ),
1399  cellRegion
1400  );
1401  cellToRegion.write();
1402 
1403  Info<< "Writing region per cell file (for manual decomposition) to "
1404  << cellToRegion.objectPath() << nl << endl;
1405  }
1406  // Write for postprocessing
1407  {
1408  volScalarField cellToRegion
1409  (
1410  IOobject
1411  (
1412  "cellToRegion",
1413  mesh.time().timeName(),
1414  mesh,
1418  ),
1419  mesh,
1422  );
1423  forAll(cellRegion, celli)
1424  {
1425  cellToRegion[celli] = cellRegion[celli];
1426  }
1427  cellToRegion.write();
1428 
1429  Info<< "Writing region per cell as volScalarField to "
1430  << cellToRegion.objectPath() << nl << endl;
1431  }
1432 }
1433 
1434 
1435 
1436 
1437 int main(int argc, char *argv[])
1438 {
1440  (
1441  "Split mesh into multiple regions (detected by walking across faces)"
1442  );
1443  #include "addRegionOption.H"
1444  #include "addOverwriteOption.H"
1446  (
1447  "cellZones",
1448  "Additionally split cellZones off into separate regions"
1449  );
1451  (
1452  "cellZonesOnly",
1453  "Use cellZones only to split mesh into regions; do not use walking"
1454  );
1456  (
1457  "cellZonesFileOnly",
1458  "file",
1459  "Like -cellZonesOnly, but use specified file"
1460  );
1462  (
1463  "combineZones",
1464  "lists of zones",
1465  "Combine zones in follow-on analysis"
1466  );
1468  (
1469  "addZones",
1470  "lists of zones",
1471  "Combine zones in follow-on analysis"
1472  );
1474  (
1475  "blockedFaces",
1476  "faceSet",
1477  "Specify additional region boundaries that walking does not cross"
1478  );
1480  (
1481  "makeCellZones",
1482  "Place cells into cellZones instead of splitting mesh"
1483  );
1485  (
1486  "largestOnly",
1487  "Only write largest region"
1488  );
1490  (
1491  "insidePoint",
1492  "point",
1493  "Only write region containing point"
1494  );
1496  (
1497  "detectOnly",
1498  "Do not write mesh"
1499  );
1501  (
1502  "sloppyCellZones",
1503  "Try to match heuristically regions to existing cell zones"
1504  );
1506  (
1507  "useFaceZones",
1508  "Use faceZones to patch inter-region faces instead of single patch"
1509  );
1511  (
1512  "prefixRegion",
1513  "Prefix region name to all patches, not just coupling patches"
1514  );
1515 
1516  argList::noFunctionObjects(); // Never use function objects
1517 
1518  #include "setRootCase.H"
1519  #include "createTime.H"
1520  #include "createNamedMesh.H"
1521 
1522  const word oldInstance = mesh.pointsInstance();
1523 
1524  word blockedFacesName;
1525  if (args.readIfPresent("blockedFaces", blockedFacesName))
1526  {
1527  Info<< "Reading blocked internal faces from faceSet "
1528  << blockedFacesName << nl << endl;
1529  }
1530 
1531  const bool makeCellZones = args.found("makeCellZones");
1532  const bool largestOnly = args.found("largestOnly");
1533  const bool insidePoint = args.found("insidePoint");
1534  const bool useCellZones = args.found("cellZones");
1535  const bool useCellZonesOnly = args.found("cellZonesOnly");
1536  const bool useCellZonesFile = args.found("cellZonesFileOnly");
1537  const bool combineZones = args.found("combineZones");
1538  const bool addZones = args.found("addZones");
1539  const bool overwrite = args.found("overwrite");
1540  const bool detectOnly = args.found("detectOnly");
1541  const bool sloppyCellZones = args.found("sloppyCellZones");
1542  const bool useFaceZones = args.found("useFaceZones");
1543  const bool prefixRegion = args.found("prefixRegion");
1544 
1545 
1546  if
1547  (
1548  (useCellZonesOnly || useCellZonesFile)
1549  && (useCellZones || blockedFacesName.size())
1550  )
1551  {
1553  << "You cannot specify both -cellZonesOnly or -cellZonesFileOnly"
1554  << " (which specify complete split)"
1555  << " in combination with -blockedFaces or -cellZones"
1556  << " (which imply a split based on topology)"
1557  << exit(FatalError);
1558  }
1559 
1560 
1561  if (useFaceZones)
1562  {
1563  Info<< "Using current faceZones to divide inter-region interfaces"
1564  << " into multiple patches."
1565  << nl << endl;
1566  }
1567  else
1568  {
1569  Info<< "Creating single patch per inter-region interface."
1570  << nl << endl;
1571  }
1572 
1573 
1574 
1575  if (insidePoint && largestOnly)
1576  {
1578  << "You cannot specify both -largestOnly"
1579  << " (keep region with most cells)"
1580  << " and -insidePoint (keep region containing point)"
1581  << exit(FatalError);
1582  }
1583 
1584 
1585  // Make sure cellZone names consistent across processors
1587 
1588  List<wordRes> zoneClusters;
1589  wordList zoneClusterNames;
1590  if (combineZones)
1591  {
1592  if (addZones)
1593  {
1595  << "Cannot specify both combineZones and addZones"
1596  << exit(FatalError);
1597  }
1598  zoneClusters = args.get<List<wordRes>>("combineZones");
1599  zoneClusterNames.setSize(zoneClusters.size());
1600  }
1601  else if (addZones)
1602  {
1603  zoneClusters = args.get<List<wordRes>>("addZones");
1604  zoneClusterNames.setSize(zoneClusters.size());
1605  forAll(zoneClusters, clusteri)
1606  {
1607  // Pop of front - is name
1608 
1609  wordRes& wrs = zoneClusters[clusteri];
1610 
1611  zoneClusterNames[clusteri] = wrs[0];
1612 
1613  for (label i = 1; i < wrs.size(); i++)
1614  {
1615  wrs[i-1] = wrs[i];
1616  }
1617  wrs.setSize(wrs.size()-1);
1618  }
1619  }
1620 
1621 
1622  // Determine per cell the region it belongs to
1623  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1624 
1625  // cellRegion is the labelList with the region per cell.
1626  labelList cellRegion;
1627  // Region to zone(s)
1628  labelListList regionToZones;
1629  // Name of region
1631  // Zone to region
1632  labelList zoneToRegion;
1633 
1634  label nCellRegions = 0;
1635  if (useCellZonesOnly)
1636  {
1637  Info<< "Using current cellZones to split mesh into regions."
1638  << " This requires all"
1639  << " cells to be in one and only one cellZone." << nl << endl;
1640 
1641  // Collect sets of zones into clusters. If no cluster is just an identity
1642  // list (cluster 0 is cellZone 0 etc.)
1643  wordList clusterNames;
1644  labelListList clusterToZones;
1645  labelList zoneToCluster;
1646  makeClusters
1647  (
1648  zoneClusters,
1649  zoneClusterNames,
1650  mesh.cellZones(),
1651  clusterNames,
1652  clusterToZones,
1653  zoneToCluster
1654  );
1655 
1656  // Existing clusterID
1657  labelList clusterID(mesh.nCells(), -1);
1658  // Neighbour clusterID.
1659  labelList neiClusterID(mesh.nBoundaryFaces());
1660  getClusterID
1661  (
1662  mesh,
1663  mesh.cellZones(),
1664  clusterNames,
1665  clusterToZones,
1666  clusterID,
1667  neiClusterID
1668  );
1669 
1670  label unzonedCelli = clusterID.find(-1);
1671  if (unzonedCelli != -1)
1672  {
1674  << "For the cellZonesOnly option all cells "
1675  << "have to be in a cellZone." << endl
1676  << "Cell " << unzonedCelli
1677  << " at" << mesh.cellCentres()[unzonedCelli]
1678  << " is not in a cellZone. There might be more unzoned cells."
1679  << exit(FatalError);
1680  }
1681  cellRegion = clusterID;
1682  nCellRegions = gMax(cellRegion)+1;
1683  zoneToRegion = zoneToCluster;
1684  regionToZones = clusterToZones;
1685  regionNames = clusterNames;
1686  }
1687  else if (useCellZonesFile)
1688  {
1689  const word zoneFile(args["cellZonesFileOnly"]);
1690  Info<< "Reading split from cellZones file " << zoneFile << endl
1691  << "This requires all"
1692  << " cells to be in one and only one cellZone." << nl << endl;
1693 
1694  cellZoneMesh newCellZones
1695  (
1696  IOobject
1697  (
1698  zoneFile,
1699  mesh.facesInstance(),
1701  mesh,
1705  ),
1706  mesh
1707  );
1708 
1709  wordList clusterNames;
1710  labelListList clusterToZones;
1711  labelList zoneToCluster;
1712  makeClusters
1713  (
1714  zoneClusters,
1715  zoneClusterNames,
1716  newCellZones,
1717  clusterNames,
1718  clusterToZones,
1719  zoneToCluster
1720  );
1721 
1722 
1723  // Existing clusterID
1724  labelList clusterID(mesh.nCells(), -1);
1725  // Neighbour clusterID.
1726  labelList neiClusterID(mesh.nBoundaryFaces());
1727  getClusterID
1728  (
1729  mesh,
1730  newCellZones,
1731  clusterNames,
1732  clusterToZones,
1733  clusterID,
1734  neiClusterID
1735  );
1736 
1737 
1738  label unzonedCelli = clusterID.find(-1);
1739  if (unzonedCelli != -1)
1740  {
1742  << "For the cellZonesFileOnly option all cells "
1743  << "have to be in a cellZone." << endl
1744  << "Cell " << unzonedCelli
1745  << " at" << mesh.cellCentres()[unzonedCelli]
1746  << " is not in a cellZone. There might be more unzoned cells."
1747  << exit(FatalError);
1748  }
1749  cellRegion = clusterID;
1750  nCellRegions = gMax(cellRegion)+1;
1751  zoneToRegion = zoneToCluster;
1752  regionToZones = clusterToZones;
1753  regionNames = clusterNames;
1754  }
1755  else
1756  {
1757  // Determine connected regions
1758  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1759 
1760  // Mark additional faces that are blocked
1761  boolList blockedFace;
1762 
1763  // Read from faceSet
1764  if (blockedFacesName.size())
1765  {
1766  faceSet blockedFaceSet(mesh, blockedFacesName);
1767  Info<< "Read "
1768  << returnReduce(blockedFaceSet.size(), sumOp<label>())
1769  << " blocked faces from set " << blockedFacesName << nl << endl;
1770 
1771  blockedFace.setSize(mesh.nFaces(), false);
1772 
1773  for (const label facei : blockedFaceSet)
1774  {
1775  blockedFace[facei] = true;
1776  }
1777  }
1778 
1779  // Collect sets of zones into clusters. If no cluster is just an
1780  // identity list (cluster 0 is cellZone 0 etc.)
1781  wordList clusterNames;
1782  labelListList clusterToZones;
1783  labelList zoneToCluster;
1784  makeClusters
1785  (
1786  zoneClusters,
1787  zoneClusterNames,
1788  mesh.cellZones(),
1789  clusterNames,
1790  clusterToZones,
1791  zoneToCluster
1792  );
1793 
1794  // Existing clusterID
1795  labelList clusterID(mesh.nCells(), -1);
1796  // Neighbour clusterID.
1797  labelList neiClusterID(mesh.nBoundaryFaces());
1798  getClusterID
1799  (
1800  mesh,
1801  mesh.cellZones(),
1802  clusterNames,
1803  clusterToZones,
1804  clusterID,
1805  neiClusterID
1806  );
1807 
1808 
1809  // Imply from differing cellZones
1810  if (useCellZones)
1811  {
1812  blockedFace.setSize(mesh.nFaces(), false);
1813 
1814  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1815  {
1816  label ownCluster = clusterID[mesh.faceOwner()[facei]];
1817  label neiCluster = clusterID[mesh.faceNeighbour()[facei]];
1818 
1819  if (ownCluster != neiCluster)
1820  {
1821  blockedFace[facei] = true;
1822  }
1823  }
1824 
1825  // Different cellZones on either side of processor patch.
1826  forAll(neiClusterID, i)
1827  {
1828  label facei = i+mesh.nInternalFaces();
1829  label ownCluster = clusterID[mesh.faceOwner()[facei]];
1830  label neiCluster = neiClusterID[i];
1831 
1832  if (ownCluster != neiCluster)
1833  {
1834  blockedFace[facei] = true;
1835  }
1836  }
1837  }
1838 
1839  // Do a topological walk to determine regions
1840  regionSplit regions(mesh, blockedFace);
1841  nCellRegions = regions.nRegions();
1842  cellRegion.transfer(regions);
1843 
1844  // Match regions to zones
1845  matchRegions
1846  (
1847  sloppyCellZones,
1848  mesh,
1849 
1850  // cluster-to-zone and cluster-to-cell addressing
1851  clusterNames,
1852  clusterToZones,
1853  clusterID,
1854 
1855  // cell-to-region addressing
1856  nCellRegions,
1857  cellRegion,
1858 
1859  // matched zones
1860  regionToZones,
1861  regionNames,
1862  zoneToRegion
1863  );
1864 
1865  // Override any default region names if single region selected
1866  if (largestOnly || insidePoint)
1867  {
1868  forAll(regionToZones, regionI)
1869  {
1870  if (regionToZones[regionI].empty())
1871  {
1872  if (overwrite)
1873  {
1875  }
1876  else if (insidePoint)
1877  {
1878  regionNames[regionI] = "insidePoint";
1879  }
1880  else if (largestOnly)
1881  {
1882  regionNames[regionI] = "largestOnly";
1883  }
1884  }
1885  }
1886  }
1887  }
1888 
1889  Info<< endl << "Number of regions:" << nCellRegions << nl << endl;
1890 
1891 
1892  // Write decomposition to file
1893  writeCellToRegion(mesh, cellRegion);
1894 
1895 
1896 
1897  // Sizes per region
1898  // ~~~~~~~~~~~~~~~~
1899 
1900  labelList regionSizes(nCellRegions, Zero);
1901 
1902  forAll(cellRegion, celli)
1903  {
1904  regionSizes[cellRegion[celli]]++;
1905  }
1906  forAll(regionSizes, regionI)
1907  {
1908  reduce(regionSizes[regionI], sumOp<label>());
1909  }
1910 
1911  Info<< "Region\tCells" << nl
1912  << "------\t-----" << endl;
1913 
1914  forAll(regionSizes, regionI)
1915  {
1916  Info<< regionI << '\t' << regionSizes[regionI] << nl;
1917  }
1918  Info<< endl;
1919 
1920 
1921 
1922  // Print region to zone
1923  Info<< "Region\tZone\tName" << nl
1924  << "------\t----\t----" << endl;
1925  forAll(regionToZones, regionI)
1926  {
1927  Info<< regionI << '\t' << flatOutput(regionToZones[regionI])
1928  << '\t'
1929  << regionNames[regionI] << nl;
1930  }
1931  Info<< endl;
1932 
1933 
1934 
1935  // Since we're going to mess with patches and zones make sure all
1936  // is synchronised
1939 
1940 
1941  // Interfaces
1942  // ----------
1943  // per interface:
1944  // - the two regions on either side
1945  // - the name
1946  // - the (global) size
1947  edgeList interfaces;
1948  List<Pair<word>> interfaceNames;
1949  labelList interfaceSizes;
1950  // per face the interface
1951  labelList faceToInterface;
1952 
1953  getInterfaceSizes
1954  (
1955  mesh,
1956  useFaceZones,
1957  cellRegion,
1958  regionNames,
1959 
1960  interfaces,
1961  interfaceNames,
1962  interfaceSizes,
1963  faceToInterface
1964  );
1965 
1966  Info<< "Sizes of interfaces between regions:" << nl << nl
1967  << "Interface\tRegion\tRegion\tFaces" << nl
1968  << "---------\t------\t------\t-----" << endl;
1969 
1970  forAll(interfaces, interI)
1971  {
1972  const edge& e = interfaces[interI];
1973 
1974  Info<< interI
1975  << "\t\t" << e[0] << "\t" << e[1]
1976  << "\t" << interfaceSizes[interI] << nl;
1977  }
1978  Info<< endl;
1979 
1980 
1981  if (detectOnly)
1982  {
1983  Info<< "End\n" << endl;
1984 
1985  return 0;
1986  }
1987 
1988 
1989  // Read objects in time directory
1990  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1991 
1992  IOobjectList objects(mesh, runTime.timeName());
1993 
1994  // Read vol fields.
1995 
1996  PtrList<volScalarField> vsFlds;
1997  ReadFields(mesh, objects, vsFlds);
1998 
1999  PtrList<volVectorField> vvFlds;
2000  ReadFields(mesh, objects, vvFlds);
2001 
2003  ReadFields(mesh, objects, vstFlds);
2004 
2005  PtrList<volSymmTensorField> vsymtFlds;
2006  ReadFields(mesh, objects, vsymtFlds);
2007 
2008  PtrList<volTensorField> vtFlds;
2009  ReadFields(mesh, objects, vtFlds);
2010 
2011  // Read surface fields.
2012 
2014  ReadFields(mesh, objects, ssFlds);
2015 
2017  ReadFields(mesh, objects, svFlds);
2018 
2020  ReadFields(mesh, objects, sstFlds);
2021 
2023  ReadFields(mesh, objects, ssymtFlds);
2024 
2026  ReadFields(mesh, objects, stFlds);
2027 
2028  Info<< endl;
2029 
2030 
2031  // Remove any demand-driven fields ('S', 'V' etc)
2032  mesh.clearOut();
2033 
2034 
2035  if (nCellRegions == 1)
2036  {
2037  Info<< "Only one region. Doing nothing." << endl;
2038  }
2039  else if (makeCellZones)
2040  {
2041  Info<< "Putting cells into cellZones instead of splitting mesh."
2042  << endl;
2043 
2044  // Check if region overlaps with existing zone. If so keep.
2045 
2046  for (label regionI = 0; regionI < nCellRegions; regionI++)
2047  {
2048  const labelList& zones = regionToZones[regionI];
2049 
2050  if (zones.size() == 1 && zones[0] != -1)
2051  {
2052  // Existing zone
2053  const label zoneI = zones[0];
2054  Info<< " Region " << regionI << " : corresponds to existing"
2055  << " cellZone "
2056  << zoneI << ' ' << mesh.cellZones()[zoneI].name() << endl;
2057  }
2058  else
2059  {
2060  // Create new cellZone.
2061  const labelList regionCells(findIndices(cellRegion, regionI));
2062 
2063  const word zoneName = "region" + Foam::name(regionI);
2064 
2065  label zoneI = mesh.cellZones().findZoneID(zoneName);
2066 
2067  if (zoneI == -1)
2068  {
2069  zoneI = mesh.cellZones().size();
2070  mesh.cellZones().setSize(zoneI+1);
2071  mesh.cellZones().set
2072  (
2073  zoneI,
2074  new cellZone
2075  (
2076  zoneName, //name
2077  std::move(regionCells), //addressing
2078  zoneI, //index
2079  mesh.cellZones() //cellZoneMesh
2080  )
2081  );
2082  }
2083  else
2084  {
2085  mesh.cellZones()[zoneI].clearAddressing();
2086  mesh.cellZones()[zoneI] = regionCells;
2087  }
2088  Info<< " Region " << regionI << " : created new cellZone "
2089  << zoneI << ' ' << mesh.cellZones()[zoneI].name() << endl;
2090  }
2091  }
2093 
2094  if (!overwrite)
2095  {
2096  ++runTime;
2098  }
2099  else
2100  {
2101  mesh.setInstance(oldInstance);
2102  }
2103 
2104  Info<< "Writing cellZones as new mesh to time " << runTime.timeName()
2105  << nl << endl;
2106 
2107  mesh.write();
2108 
2109 
2110  // Write cellSets for convenience
2111  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2112 
2113  Info<< "Writing cellSets corresponding to cellZones." << nl << endl;
2114 
2115  for (const auto& cz : mesh.cellZones())
2116  {
2117  cellSet(mesh, cz.name(), cz).write();
2118  }
2119  }
2120  else
2121  {
2122  // Add patches for interfaces
2123  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
2124 
2125  // Add all possible patches. Empty ones get filtered later on.
2126  Info<< nl << "Adding patches" << nl << endl;
2127 
2128  labelList interfacePatches
2129  (
2130  addRegionPatches
2131  (
2132  mesh,
2133  regionNames,
2134  interfaces,
2135  interfaceNames
2136  )
2137  );
2138 
2139 
2140  if (!overwrite)
2141  {
2142  ++runTime;
2143  }
2144 
2145 
2146  // Create regions
2147  // ~~~~~~~~~~~~~~
2148 
2149  if (insidePoint)
2150  {
2151  const point insidePoint = args.get<point>("insidePoint");
2152 
2153  label regionI = -1;
2154 
2155  (void)mesh.tetBasePtIs();
2156 
2157  label celli = mesh.findCell(insidePoint);
2158 
2159  Info<< nl << "Found point " << insidePoint << " in cell " << celli
2160  << endl;
2161 
2162  if (celli != -1)
2163  {
2164  regionI = cellRegion[celli];
2165  }
2166 
2167  reduce(regionI, maxOp<label>());
2168 
2169  Info<< nl
2170  << "Subsetting region " << regionI
2171  << " containing point " << insidePoint << endl;
2172 
2173  if (regionI == -1)
2174  {
2176  << "Point " << insidePoint
2177  << " is not inside the mesh." << nl
2178  << "Bounding box of the mesh:" << mesh.bounds()
2179  << exit(FatalError);
2180  }
2181 
2182  createAndWriteRegion
2183  (
2184  mesh,
2185  cellRegion,
2186  regionNames,
2187  prefixRegion,
2188  faceToInterface,
2189  interfacePatches,
2190  regionI,
2191  (overwrite ? oldInstance : runTime.timeName())
2192  );
2193  }
2194  else if (largestOnly)
2195  {
2196  label regionI = findMax(regionSizes);
2197 
2198  Info<< nl
2199  << "Subsetting region " << regionI
2200  << " of size " << regionSizes[regionI]
2201  << " as named region " << regionNames[regionI] << endl;
2202 
2203  createAndWriteRegion
2204  (
2205  mesh,
2206  cellRegion,
2207  regionNames,
2208  prefixRegion,
2209  faceToInterface,
2210  interfacePatches,
2211  regionI,
2212  (overwrite ? oldInstance : runTime.timeName())
2213  );
2214  }
2215  else
2216  {
2217  // Split all
2218  for (label regionI = 0; regionI < nCellRegions; regionI++)
2219  {
2220  Info<< nl
2221  << "Region " << regionI << nl
2222  << "-------- " << endl;
2223 
2224  createAndWriteRegion
2225  (
2226  mesh,
2227  cellRegion,
2228  regionNames,
2229  prefixRegion,
2230  faceToInterface,
2231  interfacePatches,
2232  regionI,
2233  (overwrite ? oldInstance : runTime.timeName())
2234  );
2235  }
2236  }
2237  }
2238 
2239  Info<< "End\n" << endl;
2240 
2241  return 0;
2242 }
2243 
2244 
2245 // ************************************************************************* //
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:869
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:136
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:227
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:815
"blocking" : (MPI_Bsend, MPI_Recv)
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:598
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
bool found(const Key &key) const
Same as contains()
Definition: HashTable.H:1354
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:899
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:666
Given list of cells to remove, insert all the topology changes.
Definition: removeCells.H:59
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:410
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:1049
Required Classes.
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:409
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
Input inter-processor communications stream.
Definition: IPstream.H:49
bool insert(const Key &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:342
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:618
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:316
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:435
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:629
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:137
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
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:405
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition: UPstream.H:1059
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:1102
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:670
Output inter-processor communications stream.
Definition: OPstream.H:49
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:346
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
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:56
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:1082
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:678
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
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:616
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:74
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:1185
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.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Do not request registration (bool: false)
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 swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:485
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
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