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-2021 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"
118 #include "processorMeshes.H"
119 
120 using namespace Foam;
121 
122 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
123 
124 // Prepend prefix to selected patches.
125 void renamePatches
126 (
127  fvMesh& mesh,
128  const word& prefix,
129  const labelList& patchesToRename
130 )
131 {
132  polyBoundaryMesh& polyPatches =
133  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
134  forAll(patchesToRename, i)
135  {
136  label patchi = patchesToRename[i];
137  polyPatch& pp = polyPatches[patchi];
138 
139  if (isA<coupledPolyPatch>(pp))
140  {
142  << "Encountered coupled patch " << pp.name()
143  << ". Will only rename the patch itself,"
144  << " not any referred patches."
145  << " This might have to be done by hand."
146  << endl;
147  }
148 
149  pp.name() = prefix + '_' + pp.name();
150  }
151  // Recalculate any demand driven data (e.g. group to name lookup)
152  polyPatches.updateMesh();
153 }
154 
155 
156 template<class GeoField>
157 void subsetVolFields
158 (
159  const fvMesh& mesh,
160  const fvMesh& subMesh,
161  const labelList& cellMap,
162  const labelList& faceMap,
163  const labelHashSet& addedPatches
164 )
165 {
166  const labelList patchMap(identity(mesh.boundaryMesh().size()));
167 
169  (
170  mesh.objectRegistry::lookupClass<GeoField>()
171  );
172  forAllConstIters(fields, iter)
173  {
174  const GeoField& fld = *iter.val();
175 
176  Info<< "Mapping field " << fld.name() << endl;
177 
178  tmp<GeoField> tSubFld
179  (
181  (
182  fld,
183  subMesh,
184  patchMap,
185  cellMap,
186  faceMap
187  )
188  );
189 
190  // Hack: set value to 0 for introduced patches (since don't
191  // get initialised.
192  forAll(tSubFld().boundaryField(), patchi)
193  {
194  if (addedPatches.found(patchi))
195  {
196  tSubFld.ref().boundaryFieldRef()[patchi] ==
197  typename GeoField::value_type(Zero);
198  }
199  }
200 
201  // Store on subMesh
202  GeoField* subFld = tSubFld.ptr();
203  subFld->rename(fld.name());
204  subFld->writeOpt(IOobject::AUTO_WRITE);
205  subFld->store();
206  }
207 }
208 
209 
210 template<class GeoField>
211 void subsetSurfaceFields
212 (
213  const fvMesh& mesh,
214  const fvMesh& subMesh,
215  const labelList& cellMap,
216  const labelList& faceMap,
217  const labelHashSet& addedPatches
218 )
219 {
220  const labelList patchMap(identity(mesh.boundaryMesh().size()));
221 
223  (
224  mesh.objectRegistry::lookupClass<GeoField>()
225  );
226  forAllConstIters(fields, iter)
227  {
228  const GeoField& fld = *iter.val();
229 
230  Info<< "Mapping field " << fld.name() << endl;
231 
232  tmp<GeoField> tSubFld
233  (
235  (
236  fld,
237  subMesh,
238  patchMap,
239  cellMap,
240  faceMap
241  )
242  );
243 
244  // Hack: set value to 0 for introduced patches (since don't
245  // get initialised.
246  forAll(tSubFld().boundaryField(), patchi)
247  {
248  if (addedPatches.found(patchi))
249  {
250  tSubFld.ref().boundaryFieldRef()[patchi] ==
251  typename GeoField::value_type(Zero);
252  }
253  }
254 
255  // Store on subMesh
256  GeoField* subFld = tSubFld.ptr();
257  subFld->rename(fld.name());
258  subFld->writeOpt(IOobject::AUTO_WRITE);
259  subFld->store();
260  }
261 }
262 
263 // Select all cells not in the region
264 labelList getNonRegionCells(const labelList& cellRegion, const label regionI)
265 {
266  DynamicList<label> nonRegionCells(cellRegion.size());
267  forAll(cellRegion, celli)
268  {
269  if (cellRegion[celli] != regionI)
270  {
271  nonRegionCells.append(celli);
272  }
273  }
274  return nonRegionCells.shrink();
275 }
276 
277 
278 void addToInterface
279 (
280  const polyMesh& mesh,
281  const label zoneID,
282  const label ownRegion,
283  const label neiRegion,
284  EdgeMap<Map<label>>& regionsToSize
285 )
286 {
288  (
289  min(ownRegion, neiRegion),
290  max(ownRegion, neiRegion)
291  );
292 
293  auto iter = regionsToSize.find(interface);
294 
295  if (iter.found())
296  {
297  // Check if zone present
298  auto zoneIter = iter().find(zoneID);
299  if (zoneIter.found())
300  {
301  // Found zone. Increment count.
302  ++(*zoneIter);
303  }
304  else
305  {
306  // New or no zone. Insert with count 1.
307  iter().insert(zoneID, 1);
308  }
309  }
310  else
311  {
312  // Create new interface of size 1.
313  Map<label> zoneToSize;
314  zoneToSize.insert(zoneID, 1);
315  regionsToSize.insert(interface, zoneToSize);
316  }
317 }
318 
319 
320 // Get region-region interface name and sizes.
321 // Returns interfaces as straight list for looping in identical order.
322 void getInterfaceSizes
323 (
324  const polyMesh& mesh,
325  const bool useFaceZones,
326  const labelList& cellRegion,
327  const wordList& regionNames,
328 
329  edgeList& interfaces,
330  List<Pair<word>>& interfaceNames,
331  labelList& interfaceSizes,
332  labelList& faceToInterface
333 )
334 {
335  // From region-region to faceZone (or -1) to number of faces.
336 
337  EdgeMap<Map<label>> regionsToSize;
338 
339 
340  // Internal faces
341  // ~~~~~~~~~~~~~~
342 
343  forAll(mesh.faceNeighbour(), facei)
344  {
345  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
346  label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
347 
348  if (ownRegion != neiRegion)
349  {
350  addToInterface
351  (
352  mesh,
353  (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
354  ownRegion,
355  neiRegion,
356  regionsToSize
357  );
358  }
359  }
360 
361  // Boundary faces
362  // ~~~~~~~~~~~~~~
363 
364  // Neighbour cellRegion.
365  labelList coupledRegion(mesh.nBoundaryFaces());
366 
367  forAll(coupledRegion, i)
368  {
369  label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
370  coupledRegion[i] = cellRegion[celli];
371  }
372  syncTools::swapBoundaryFaceList(mesh, coupledRegion);
373 
374  forAll(coupledRegion, i)
375  {
376  label facei = i+mesh.nInternalFaces();
377  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
378  label neiRegion = coupledRegion[i];
379 
380  if (ownRegion != neiRegion)
381  {
382  addToInterface
383  (
384  mesh,
385  (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
386  ownRegion,
387  neiRegion,
388  regionsToSize
389  );
390  }
391  }
392 
393 
394  if (Pstream::parRun())
395  {
396  if (Pstream::master())
397  {
398  // Receive and add to my sizes
399  for (const int slave : Pstream::subProcs())
400  {
401  IPstream fromSlave(Pstream::commsTypes::blocking, slave);
402 
403  EdgeMap<Map<label>> slaveSizes(fromSlave);
404 
405  forAllConstIters(slaveSizes, slaveIter)
406  {
407  const Map<label>& slaveInfo = *slaveIter;
408 
409  auto masterIter = regionsToSize.find(slaveIter.key());
410 
411  if (masterIter.found())
412  {
413  // Same inter-region
414  Map<label>& masterInfo = *masterIter;
415 
416  forAllConstIters(slaveInfo, iter)
417  {
418  const label zoneID = iter.key();
419  const label slaveSize = iter.val();
420 
421  auto zoneIter = masterInfo.find(zoneID);
422  if (zoneIter.found())
423  {
424  *zoneIter += slaveSize;
425  }
426  else
427  {
428  masterInfo.insert(zoneID, slaveSize);
429  }
430  }
431  }
432  else
433  {
434  regionsToSize.insert(slaveIter.key(), slaveInfo);
435  }
436  }
437  }
438  }
439  else
440  {
441  // Send to master
442  {
443  OPstream toMaster
444  (
447  );
448  toMaster << regionsToSize;
449  }
450  }
451  }
452 
453  // Rework
454 
455  Pstream::scatter(regionsToSize);
456 
457 
458 
459  // Now we have the global sizes of all inter-regions.
460  // Invert this on master and distribute.
461  label nInterfaces = 0;
462  forAllConstIters(regionsToSize, iter)
463  {
464  const Map<label>& info = iter.val();
465  nInterfaces += info.size();
466  }
467 
468  interfaces.setSize(nInterfaces);
469  interfaceNames.setSize(nInterfaces);
470  interfaceSizes.setSize(nInterfaces);
471  EdgeMap<Map<label>> regionsToInterface(nInterfaces);
472 
473  nInterfaces = 0;
474  forAllConstIters(regionsToSize, iter)
475  {
476  const edge& e = iter.key();
477  const Map<label>& info = iter.val();
478 
479  const word& name0 = regionNames[e[0]];
480  const word& name1 = regionNames[e[1]];
481 
482  forAllConstIters(info, infoIter)
483  {
484  interfaces[nInterfaces] = iter.key();
485  label zoneID = infoIter.key();
486  if (zoneID == -1)
487  {
488  interfaceNames[nInterfaces] = Pair<word>
489  (
490  name0 + "_to_" + name1,
491  name1 + "_to_" + name0
492  );
493  }
494  else
495  {
496  const word& zoneName = mesh.faceZones()[zoneID].name();
497  interfaceNames[nInterfaces] = Pair<word>
498  (
499  zoneName + "_" + name0 + "_to_" + name1,
500  zoneName + "_" + name1 + "_to_" + name0
501  );
502  }
503  interfaceSizes[nInterfaces] = infoIter();
504 
505  if (regionsToInterface.found(e))
506  {
507  regionsToInterface[e].insert(zoneID, nInterfaces);
508  }
509  else
510  {
511  Map<label> zoneAndInterface;
512  zoneAndInterface.insert(zoneID, nInterfaces);
513  regionsToInterface.insert(e, zoneAndInterface);
514  }
515  nInterfaces++;
516  }
517  }
518 
519 
520  // Now all processor have consistent interface information
521 
522  Pstream::scatter(interfaces);
523  Pstream::scatter(interfaceNames);
524  Pstream::scatter(interfaceSizes);
525  Pstream::scatter(regionsToInterface);
526 
527  // Mark all inter-region faces.
528  faceToInterface.setSize(mesh.nFaces(), -1);
529 
530  forAll(mesh.faceNeighbour(), facei)
531  {
532  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
533  label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
534 
535  if (ownRegion != neiRegion)
536  {
537  label zoneID = -1;
538  if (useFaceZones)
539  {
540  zoneID = mesh.faceZones().whichZone(facei);
541  }
542 
544  (
545  min(ownRegion, neiRegion),
546  max(ownRegion, neiRegion)
547  );
548 
549  faceToInterface[facei] = regionsToInterface[interface][zoneID];
550  }
551  }
552  forAll(coupledRegion, i)
553  {
554  label facei = i+mesh.nInternalFaces();
555  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
556  label neiRegion = coupledRegion[i];
557 
558  if (ownRegion != neiRegion)
559  {
560  label zoneID = -1;
561  if (useFaceZones)
562  {
563  zoneID = mesh.faceZones().whichZone(facei);
564  }
565 
567  (
568  min(ownRegion, neiRegion),
569  max(ownRegion, neiRegion)
570  );
571 
572  faceToInterface[facei] = regionsToInterface[interface][zoneID];
573  }
574  }
575 }
576 
577 
578 // Create mesh for region.
579 autoPtr<mapPolyMesh> createRegionMesh
580 (
581  const fvMesh& mesh,
582  // Region info
583  const labelList& cellRegion,
584  const label regionI,
585  const word& regionName,
586  // Interface info
587  const labelList& interfacePatches,
588  const labelList& faceToInterface,
589 
590  autoPtr<fvMesh>& newMesh
591 )
592 {
593  // Create dummy system/fv*
595 
596  // Neighbour cellRegion.
597  labelList coupledRegion(mesh.nBoundaryFaces());
598 
599  forAll(coupledRegion, i)
600  {
601  label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
602  coupledRegion[i] = cellRegion[celli];
603  }
604  syncTools::swapBoundaryFaceList(mesh, coupledRegion);
605 
606 
607  // Topology change container. Start off from existing mesh.
608  polyTopoChange meshMod(mesh);
609 
610  // Cell remover engine
611  removeCells cellRemover(mesh);
612 
613  // Select all but region cells
614  labelList cellsToRemove(getNonRegionCells(cellRegion, regionI));
615 
616  // Find out which faces will get exposed. Note that this
617  // gets faces in mesh face order. So both regions will get same
618  // face in same order (important!)
619  labelList exposedFaces = cellRemover.getExposedFaces(cellsToRemove);
620 
621  labelList exposedPatchIDs(exposedFaces.size());
622  forAll(exposedFaces, i)
623  {
624  label facei = exposedFaces[i];
625  label interfacei = faceToInterface[facei];
626 
627  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
628  label neiRegion = -1;
629 
630  if (mesh.isInternalFace(facei))
631  {
632  neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
633  }
634  else
635  {
636  neiRegion = coupledRegion[facei-mesh.nInternalFaces()];
637  }
638 
639 
640  // Check which side is being kept - determines which of the two
641  // patches will be used.
642 
643  label otherRegion = -1;
644 
645  if (ownRegion == regionI && neiRegion != regionI)
646  {
647  otherRegion = neiRegion;
648  }
649  else if (ownRegion != regionI && neiRegion == regionI)
650  {
651  otherRegion = ownRegion;
652  }
653  else
654  {
656  << "Exposed face:" << facei
657  << " fc:" << mesh.faceCentres()[facei]
658  << " has owner region " << ownRegion
659  << " and neighbour region " << neiRegion
660  << " when handling region:" << regionI
661  << exit(FatalError);
662  }
663 
664  // Find the patch.
665  if (regionI < otherRegion)
666  {
667  exposedPatchIDs[i] = interfacePatches[interfacei];
668  }
669  else
670  {
671  exposedPatchIDs[i] = interfacePatches[interfacei]+1;
672  }
673  }
674 
675  // Remove faces
676  cellRemover.setRefinement
677  (
678  cellsToRemove,
679  exposedFaces,
680  exposedPatchIDs,
681  meshMod
682  );
683 
684  autoPtr<mapPolyMesh> map = meshMod.makeMesh
685  (
686  newMesh,
687  IOobject
688  (
689  regionName,
690  mesh.time().timeName(),
691  mesh.time(),
692  IOobject::READ_IF_PRESENT, // read fv* if present
694  ),
695  mesh
696  );
697 
698  return map;
699 }
700 
701 
702 void createAndWriteRegion
703 (
704  const fvMesh& mesh,
705  const labelList& cellRegion,
706  const wordList& regionNames,
707  const bool prefixRegion,
708  const labelList& faceToInterface,
709  const labelList& interfacePatches,
710  const label regionI,
711  const word& newMeshInstance
712 )
713 {
714  Info<< "Creating mesh for region " << regionI
715  << ' ' << regionNames[regionI] << endl;
716 
717  autoPtr<fvMesh> newMesh;
718  autoPtr<mapPolyMesh> map = createRegionMesh
719  (
720  mesh,
721  cellRegion,
722  regionI,
723  regionNames[regionI],
724  interfacePatches,
725  faceToInterface,
726  newMesh
727  );
728 
729 
730  // Make map of all added patches
731  labelHashSet addedPatches(2*interfacePatches.size());
732  forAll(interfacePatches, interfacei)
733  {
734  addedPatches.insert(interfacePatches[interfacei]);
735  addedPatches.insert(interfacePatches[interfacei]+1);
736  }
737 
738 
739  Info<< "Mapping fields" << endl;
740 
741  // Map existing fields
742  newMesh().updateMesh(map());
743 
744  // Add subsetted fields
745  subsetVolFields<volScalarField>
746  (
747  mesh,
748  newMesh(),
749  map().cellMap(),
750  map().faceMap(),
751  addedPatches
752  );
753  subsetVolFields<volVectorField>
754  (
755  mesh,
756  newMesh(),
757  map().cellMap(),
758  map().faceMap(),
759  addedPatches
760  );
761  subsetVolFields<volSphericalTensorField>
762  (
763  mesh,
764  newMesh(),
765  map().cellMap(),
766  map().faceMap(),
767  addedPatches
768  );
769  subsetVolFields<volSymmTensorField>
770  (
771  mesh,
772  newMesh(),
773  map().cellMap(),
774  map().faceMap(),
775  addedPatches
776  );
777  subsetVolFields<volTensorField>
778  (
779  mesh,
780  newMesh(),
781  map().cellMap(),
782  map().faceMap(),
783  addedPatches
784  );
785 
786  subsetSurfaceFields<surfaceScalarField>
787  (
788  mesh,
789  newMesh(),
790  map().cellMap(),
791  map().faceMap(),
792  addedPatches
793  );
794  subsetSurfaceFields<surfaceVectorField>
795  (
796  mesh,
797  newMesh(),
798  map().cellMap(),
799  map().faceMap(),
800  addedPatches
801  );
802  subsetSurfaceFields<surfaceSphericalTensorField>
803  (
804  mesh,
805  newMesh(),
806  map().cellMap(),
807  map().faceMap(),
808  addedPatches
809  );
810  subsetSurfaceFields<surfaceSymmTensorField>
811  (
812  mesh,
813  newMesh(),
814  map().cellMap(),
815  map().faceMap(),
816  addedPatches
817  );
818  subsetSurfaceFields<surfaceTensorField>
819  (
820  mesh,
821  newMesh(),
822  map().cellMap(),
823  map().faceMap(),
824  addedPatches
825  );
826 
827 
828  const polyBoundaryMesh& newPatches = newMesh().boundaryMesh();
829  newPatches.checkParallelSync(true);
830 
831  // Delete empty patches
832  // ~~~~~~~~~~~~~~~~~~~~
833 
834  // Create reordering list to move patches-to-be-deleted to end
835  labelList oldToNew(newPatches.size(), -1);
836  DynamicList<label> sharedPatches(newPatches.size());
837  label newI = 0;
838 
839  Info<< "Deleting empty patches" << endl;
840 
841  // Assumes all non-proc boundaries are on all processors!
842  forAll(newPatches, patchi)
843  {
844  const polyPatch& pp = newPatches[patchi];
845 
846  if (!isA<processorPolyPatch>(pp))
847  {
848  if (returnReduce(pp.size(), sumOp<label>()) > 0)
849  {
850  oldToNew[patchi] = newI;
851  if (!addedPatches.found(patchi))
852  {
853  sharedPatches.append(newI);
854  }
855  newI++;
856  }
857  }
858  }
859 
860  // Same for processor patches (but need no reduction)
861  forAll(newPatches, patchi)
862  {
863  const polyPatch& pp = newPatches[patchi];
864 
865  if (isA<processorPolyPatch>(pp) && pp.size())
866  {
867  oldToNew[patchi] = newI++;
868  }
869  }
870 
871  const label nNewPatches = newI;
872 
873  // Move all deleteable patches to the end
874  forAll(oldToNew, patchi)
875  {
876  if (oldToNew[patchi] == -1)
877  {
878  oldToNew[patchi] = newI++;
879  }
880  }
881 
882  //reorderPatches(newMesh(), oldToNew, nNewPatches);
883  fvMeshTools::reorderPatches(newMesh(), oldToNew, nNewPatches, true);
884 
885  // Rename shared patches with region name
886  if (prefixRegion)
887  {
888  Info<< "Prefixing patches with region name" << endl;
889 
890  renamePatches(newMesh(), regionNames[regionI], sharedPatches);
891  }
892 
893 
894  Info<< "Writing new mesh" << endl;
895 
896  newMesh().setInstance(newMeshInstance);
897  newMesh().write();
898  topoSet::removeFiles(newMesh());
899  processorMeshes::removeFiles(newMesh());
900 
901  // Write addressing files like decomposePar
902  Info<< "Writing addressing to base mesh" << endl;
903 
904  labelIOList pointProcAddressing
905  (
906  IOobject
907  (
908  "pointRegionAddressing",
909  newMesh().facesInstance(),
910  newMesh().meshSubDir,
911  newMesh(),
914  false
915  ),
916  map().pointMap()
917  );
918  Info<< "Writing map " << pointProcAddressing.name()
919  << " from region" << regionI
920  << " points back to base mesh." << endl;
921  pointProcAddressing.write();
922 
924  (
925  IOobject
926  (
927  "faceRegionAddressing",
928  newMesh().facesInstance(),
929  newMesh().meshSubDir,
930  newMesh(),
933  false
934  ),
935  newMesh().nFaces()
936  );
937  forAll(faceProcAddressing, facei)
938  {
939  // face + turning index. (see decomposePar)
940  // Is the face pointing in the same direction?
941  label oldFacei = map().faceMap()[facei];
942 
943  if
944  (
945  map().cellMap()[newMesh().faceOwner()[facei]]
946  == mesh.faceOwner()[oldFacei]
947  )
948  {
949  faceProcAddressing[facei] = oldFacei+1;
950  }
951  else
952  {
953  faceProcAddressing[facei] = -(oldFacei+1);
954  }
955  }
956  Info<< "Writing map " << faceProcAddressing.name()
957  << " from region" << regionI
958  << " faces back to base mesh." << endl;
959  faceProcAddressing.write();
960 
961  labelIOList cellProcAddressing
962  (
963  IOobject
964  (
965  "cellRegionAddressing",
966  newMesh().facesInstance(),
967  newMesh().meshSubDir,
968  newMesh(),
971  false
972  ),
973  map().cellMap()
974  );
975  Info<< "Writing map " <<cellProcAddressing.name()
976  << " from region" << regionI
977  << " cells back to base mesh." << endl;
978  cellProcAddressing.write();
979 
980  labelIOList boundaryProcAddressing
981  (
982  IOobject
983  (
984  "boundaryRegionAddressing",
985  newMesh().facesInstance(),
986  newMesh().meshSubDir,
987  newMesh(),
990  false
991  ),
992  labelList(nNewPatches, -1)
993  );
994  forAll(oldToNew, i)
995  {
996  if (!addedPatches.found(i))
997  {
998  label newI = oldToNew[i];
999  if (newI >= 0 && newI < nNewPatches)
1000  {
1001  boundaryProcAddressing[oldToNew[i]] = i;
1002  }
1003  }
1004  }
1005  Info<< "Writing map " << boundaryProcAddressing.name()
1006  << " from region" << regionI
1007  << " boundary back to base mesh." << endl;
1008  boundaryProcAddressing.write();
1009 }
1010 
1011 
1012 // Create for every region-region interface with non-zero size two patches.
1013 // First one is for minimumregion to maximumregion.
1014 // Note that patches get created in same order on all processors (if parallel)
1015 // since looping over synchronised 'interfaces'.
1016 labelList addRegionPatches
1017 (
1018  fvMesh& mesh,
1019  const wordList& regionNames,
1020  const edgeList& interfaces,
1021  const List<Pair<word>>& interfaceNames
1022 )
1023 {
1024  Info<< nl << "Adding patches" << nl << endl;
1025 
1026  labelList interfacePatches(interfaces.size());
1027 
1028  forAll(interfaces, interI)
1029  {
1030  const edge& e = interfaces[interI];
1031  const Pair<word>& names = interfaceNames[interI];
1032 
1033  //Info<< "For interface " << interI
1034  // << " between regions " << e
1035  // << " trying to add patches " << names << endl;
1036 
1037 
1038  mappedWallPolyPatch patch1
1039  (
1040  names[0],
1041  0, // overridden
1042  0, // overridden
1043  0, // overridden
1044  regionNames[e[1]], // sampleRegion
1046  names[1], // samplePatch
1047  point::zero, // offset
1048  mesh.boundaryMesh()
1049  );
1050 
1051  interfacePatches[interI] = fvMeshTools::addPatch
1052  (
1053  mesh,
1054  patch1,
1055  dictionary(), //optional per field value
1057  true //validBoundary
1058  );
1059 
1060  mappedWallPolyPatch patch2
1061  (
1062  names[1],
1063  0,
1064  0,
1065  0,
1066  regionNames[e[0]], // sampleRegion
1068  names[0],
1069  point::zero, // offset
1070  mesh.boundaryMesh()
1071  );
1073  (
1074  mesh,
1075  patch2,
1076  dictionary(), //optional per field value
1078  true //validBoundary
1079  );
1080 
1081  Info<< "For interface between region " << regionNames[e[0]]
1082  << " and " << regionNames[e[1]] << " added patches" << endl
1083  << " " << interfacePatches[interI]
1084  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]].name()
1085  << endl
1086  << " " << interfacePatches[interI]+1
1087  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]+1].name()
1088  << endl;
1089  }
1090  return interfacePatches;
1091 }
1092 
1093 
1094 // Find region that covers most of cell zone
1095 label findCorrespondingRegion
1096 (
1097  const labelList& existingZoneID, // per cell the (unique) zoneID
1098  const labelList& cellRegion,
1099  const label nCellRegions,
1100  const label zoneI,
1101  const label minOverlapSize
1102 )
1103 {
1104  // Per region the number of cells in zoneI
1105  labelList cellsInZone(nCellRegions, Zero);
1106 
1107  forAll(cellRegion, celli)
1108  {
1109  if (existingZoneID[celli] == zoneI)
1110  {
1111  cellsInZone[cellRegion[celli]]++;
1112  }
1113  }
1114 
1116  Pstream::listCombineScatter(cellsInZone);
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,
1411  false
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,
1431  false
1432  ),
1433  mesh,
1435  zeroGradientFvPatchScalarField::typeName
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,
1718  false
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\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\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\t" << e[0] << "\t\t" << e[1]
1990  << "\t\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  }
2106  mesh.cellZones().writeOpt(IOobject::AUTO_WRITE);
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 // ************************************************************************* //
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:63
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::maxOp
Definition: ops.H:217
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
faceProcAddressing
PtrList< labelIOList > & faceProcAddressing
Definition: checkFaceAddressingComp.H:9
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:190
mappedWallPolyPatch.H
Foam::UPstream::masterNo
static constexpr int masterNo() noexcept
Definition: UPstream.H:447
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:87
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Definition: fvMesh.C:1034
Foam::fvMeshTools::createDummyFvMeshFiles
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Definition: fvMeshTools.C:753
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:88
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:59
Foam::polyMesh::defaultRegion
static word defaultRegion
Definition: polyMesh.H:314
Foam::ReadFields
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fields, const bool syncPar=true, const bool readOldTime=false)
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:57
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::removeCells
Given list of cells to remove, insert all the topology changes.
Definition: removeCells.H:59
Foam::DynamicList< label >
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:49
Foam::minOp
Definition: ops.H:218
fvMeshSubset.H
interface
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
addOverwriteOption.H
Foam::polyMesh::meshSubDir
static word meshSubDir
Definition: polyMesh.H:317
polyTopoChange.H
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:95
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:59
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Definition: polyMesh.C:845
Foam::HashTable< T, edge, Hash< edge > >::insert
bool insert(const Key &key, const T &obj)
Definition: HashTableI.H:173
Foam::Map< label >
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Definition: UPstream.H:453
Foam::syncTools::swapBoundaryFaceList
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Definition: syncTools.H:441
Foam::faceSet
A list of face labels.
Definition: faceSet.H:47
Foam::Pstream::scatter
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Definition: gatherScatter.C:144
Foam::List::append
void append(const T &val)
IOobjectList.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Definition: polyMesh.H:440
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
removeCells.H
Foam::HashSet
A HashTable with keys but without contents that is similar to std::unordered_set.
Definition: HashSet.H:73
processorMeshes.H
Foam::argList::get
T get(const label index) const
Definition: argListI.H:271
syncTools.H
Foam::argList::readIfPresent
bool readIfPresent(const word &optName, T &val) const
Definition: argListI.H:316
zoneIDs
const labelIOList & zoneIDs
Definition: correctPhi.H:59
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Definition: hashSets.C:26
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
regionNames
wordList regionNames
Definition: getAllRegionOptions.H:31
Foam::sumOp
Definition: ops.H:207
Foam::ZoneMesh::indices
labelList indices(const wordRe &matcher, const bool useGroups=true) const
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:58
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Definition: polyMesh.C:839
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::polyBoundaryMesh::checkParallelSync
bool checkParallelSync(const bool report=false) const
Definition: polyBoundaryMesh.C:965
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:45
Foam::argList::noFunctionObjects
static void noFunctionObjects(bool addWithOption=false)
Definition: argList.C:466
Foam::primitiveMesh::nCells
label nCells() const noexcept
Definition: primitiveMeshI.H:89
Foam::UPstream::subProcs
static rangeType subProcs(const label communicator=worldComm)
Definition: UPstream.H:511
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Definition: regIOobjectWrite.C:125
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:64
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Definition: polyMesh.H:488
Foam::List::setSize
void setSize(const label n)
Definition: List.H:218
argList.H
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Definition: polyMesh.C:1100
faceSet.H
Foam::List::transfer
void transfer(List< T > &list)
regionSplit.H
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:183
addRegionOption.H
Foam::primitiveMesh::nBoundaryFaces
label nBoundaryFaces() const noexcept
Definition: primitiveMeshI.H:77
Foam::ZoneMesh< cellZone, polyMesh >
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Definition: polyMesh.H:482
Foam::regionSplit
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:136
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
Foam::fvMeshTools::reorderPatches
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Definition: fvMeshTools.C:322
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
Foam::PtrList::append
void append(T *ptr)
Definition: PtrListI.H:106
fld
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))
Definition: gmvOutputLagrangian.H:23
Foam::syncTools::swapBoundaryCellList
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Definition: syncToolsTemplates.C:1385
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Definition: ZoneMesh.C:282
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
createNamedMesh.H
Required Variables.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:36
zoneID
const labelIOList & zoneID
Definition: interpolatedFaces.H:22
Foam::Ostream::write
virtual bool write(const token &tok)=0
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
Foam::calculatedFvPatchField
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Definition: calculatedFvPatchField.H:62
Foam
Definition: atmBoundaryLayer.C:26
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Definition: FlatOutput.H:217
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:47
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Definition: ZoneMesh.C:512
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:51
Foam::HashTable< T, edge, Hash< edge > >::find
iterator find(const Key &key)
Definition: HashTableI.H:107
Foam::mappedWallPolyPatch
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedWallPolyPatch.H:55
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::zoneIdentifier::index
label index() const noexcept
Definition: zoneIdentifier.H:131
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:101
Foam::polyBoundaryMesh::updateMesh
void updateMesh()
Definition: polyBoundaryMesh.C:1156
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:284
Foam::polyMesh::bounds
const boundBox & bounds() const
Definition: polyMesh.H:446
Foam::IOobject::name
const word & name() const noexcept
Definition: IOobjectI.H:58
Foam::argList::addBoolOption
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Definition: argList.C:317
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:49
setRootCase.H
Foam::EdgeMap
Map from edge (expressed as its endpoints) to value. For easier forward declaration it is currently i...
Definition: EdgeMap.H:45
Foam::polyMesh::findCell
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Definition: polyMesh.C:1500
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:78
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const noexcept
Definition: primitiveMeshI.H:96
Foam::fvMeshSubset::interpolate
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelUList &patchMap, const labelUList &cellMap, const labelUList &faceMap, const bool allowUnmapped=false)
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
edgeHashes.H
ReadFields.H
Field reading functions for post-processing utilities.
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::Pair< word >
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::zoneIdentifier::name
const word & name() const noexcept
Definition: zoneIdentifier.H:119
fvMeshTools.H
Foam::mappedPatchBase::NEARESTPATCHFACE
@ NEARESTPATCHFACE
nearest face on selected patch
Definition: mappedPatchBase.H:118
Foam::Vector< scalar >
Foam::findMax
label findMax(const ListType &input, label start=0)
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const noexcept
Definition: primitiveMeshI.H:71
Foam::UPstream::parRun
static bool & parRun() noexcept
Definition: UPstream.H:429
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::fvMeshTools::addPatch
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Definition: fvMeshTools.C:30
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:71
Foam::processorMeshes::removeFiles
static void removeFiles(const polyMesh &mesh)
Definition: processorMeshes.C:267
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::identity
labelList identity(const label len, label start=0)
Definition: labelList.C:31
Foam::IOList< label >
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:47
Foam::List::clear
void clear()
Definition: ListI.H:109
Foam::HashSet::insert
bool insert(const Key &key)
Definition: HashSet.H:191
Foam::word::null
static const word null
Definition: word.H:78
createTime.H
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
Foam::fvMesh::time
const Time & time() const
Definition: fvMesh.H:276
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Form zero
Definition: VectorSpace.H:111
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:56
Foam::findIndices
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:49
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Definition: primitiveMeshI.H:83
Foam::patchIdentifier::name
const word & name() const noexcept
Definition: patchIdentifier.H:131
Foam::plusEqOp
Definition: ops.H:66
cellSet.H
Foam::ZoneMesh::clearAddressing
void clearAddressing()
Definition: ZoneMesh.C:702
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Definition: combineGatherScatter.C:426
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::argList::addOption
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Definition: argList.C:328
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
args
Foam::argList args(argc, argv)
Foam::PtrListOps::names
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Foam::topoSet::removeFiles
static void removeFiles(const polyMesh &)
Definition: topoSet.C:628
zeroGradientFvPatchFields.H
Foam::ZoneMesh::checkParallelSync
bool checkParallelSync(const bool report=false) const
Definition: ZoneMesh.C:745
WarningInFunction
#define WarningInFunction
Definition: messageStream.H:353
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:585
Foam::polyMesh::setInstance
void setInstance(const fileName &instance, const IOobject::writeOption wOpt=IOobject::AUTO_WRITE)
Definition: polyMeshIO.C:29
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Definition: polyMesh.C:1106
Foam::dimless
const dimensionSet dimless
Definition: dimensionSets.C:182
Foam::fvMesh::clearOut
void clearOut()
Definition: fvMesh.C:232
Foam::argList::found
bool found(const word &optName) const
Definition: argListI.H:171
Foam::fvMesh::name
const word & name() const
Definition: fvMesh.H:296
Foam::UPstream::commsTypes::blocking
@ blocking
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:181
Foam::polyMesh::tetBasePtIs
const labelIOList & tetBasePtIs() const
Definition: polyMesh.C:885