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