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 | Copyright (C) 2011-2017 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  splitMeshRegions
26 
27 Description
28  Splits mesh into multiple regions.
29 
30  Each region is defined as a domain whose cells can all be reached by
31  cell-face-cell walking without crossing
32  - boundary faces
33  - additional faces from faceset (-blockedFaces faceSet).
34  - any face inbetween differing cellZones (-cellZones)
35 
36  Output is:
37  - volScalarField with regions as different scalars (-detectOnly)
38  or
39  - mesh with multiple regions and mapped patches. These patches
40  either cover the whole interface between two region (default) or
41  only part according to faceZones (-useFaceZones)
42  or
43  - mesh with cells put into cellZones (-makeCellZones)
44 
45  Note:
46  - cellZonesOnly does not do a walk and uses the cellZones only. Use
47  this if you don't mind having disconnected domains in a single region.
48  This option requires all cells to be in one (and one only) cellZone.
49 
50  - cellZonesFileOnly behaves like -cellZonesOnly but reads the cellZones
51  from the specified file. This allows one to explicitly specify the region
52  distribution and still have multiple cellZones per region.
53 
54  - useCellZonesOnly does not do a walk and uses the cellZones only. Use
55  this if you don't mind having disconnected domains in a single region.
56  This option requires all cells to be in one (and one only) cellZone.
57 
58  - prefixRegion prefixes all normal patches with region name (interface
59  (patches already have region name prefix)
60 
61  - Should work in parallel.
62  cellZones can differ on either side of processor boundaries in which case
63  the faces get moved from processor patch to directMapped patch. Not
64  the faces get moved from processor patch to mapped patch. Not
65  very well tested.
66 
67  - If a cell zone gets split into more than one region it can detect
68  the largest matching region (-sloppyCellZones). This will accept any
69  region that covers more than 50% of the zone. It has to be a subset
70  so cannot have any cells in any other zone.
71 
72  - If explicitly a single region has been selected (-largestOnly or
73  -insidePoint) its region name will be either
74  - name of a cellZone it matches to or
75  - "largestOnly" respectively "insidePoint" or
76  - polyMesh::defaultRegion if additionally -overwrite
77  (so it will overwrite the input mesh!)
78 
79  - writes maps like decomposePar back to original mesh:
80  - pointRegionAddressing : for every point in this region the point in
81  the original mesh
82  - cellRegionAddressing : ,, cell ,, cell ,,
83  - faceRegionAddressing : ,, face ,, face in
84  the original mesh + 'turning index'. For a face in the same orientation
85  this is the original facelabel+1, for a turned face this is -facelabel-1
86  - boundaryRegionAddressing : for every patch in this region the
87  patch in the original mesh (or -1 if added patch)
88 
89 \*---------------------------------------------------------------------------*/
90 
91 #include "SortableList.H"
92 #include "argList.H"
93 #include "regionSplit.H"
94 #include "fvMeshSubset.H"
95 #include "IOobjectList.H"
96 #include "volFields.H"
97 #include "faceSet.H"
98 #include "cellSet.H"
99 #include "polyTopoChange.H"
100 #include "removeCells.H"
101 #include "EdgeMap.H"
102 #include "syncTools.H"
103 #include "ReadFields.H"
104 #include "mappedWallPolyPatch.H"
105 #include "fvMeshTools.H"
107 
108 using namespace Foam;
109 
110 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
111 
112 // Prepend prefix to selected patches.
113 void renamePatches
114 (
115  fvMesh& mesh,
116  const word& prefix,
117  const labelList& patchesToRename
118 )
119 {
120  polyBoundaryMesh& polyPatches =
121  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
122  forAll(patchesToRename, i)
123  {
124  label patchi = patchesToRename[i];
125  polyPatch& pp = polyPatches[patchi];
126 
127  if (isA<coupledPolyPatch>(pp))
128  {
130  << "Encountered coupled patch " << pp.name()
131  << ". Will only rename the patch itself,"
132  << " not any referred patches."
133  << " This might have to be done by hand."
134  << endl;
135  }
136 
137  pp.name() = prefix + '_' + pp.name();
138  }
139  // Recalculate any demand driven data (e.g. group to name lookup)
140  polyPatches.updateMesh();
141 }
142 
143 
144 template<class GeoField>
145 void subsetVolFields
146 (
147  const fvMesh& mesh,
148  const fvMesh& subMesh,
149  const labelList& cellMap,
150  const labelList& faceMap,
151  const labelHashSet& addedPatches
152 )
153 {
154  const labelList patchMap(identity(mesh.boundaryMesh().size()));
155 
157  (
158  mesh.objectRegistry::lookupClass<GeoField>()
159  );
161  {
162  const GeoField& fld = *iter();
163 
164  Info<< "Mapping field " << fld.name() << endl;
165 
166  tmp<GeoField> tSubFld
167  (
169  (
170  fld,
171  subMesh,
172  patchMap,
173  cellMap,
174  faceMap
175  )
176  );
177 
178  // Hack: set value to 0 for introduced patches (since don't
179  // get initialised.
180  forAll(tSubFld().boundaryField(), patchi)
181  {
182  if (addedPatches.found(patchi))
183  {
184  tSubFld.ref().boundaryFieldRef()[patchi] ==
185  typename GeoField::value_type(Zero);
186  }
187  }
188 
189  // Store on subMesh
190  GeoField* subFld = tSubFld.ptr();
191  subFld->rename(fld.name());
192  subFld->writeOpt() = IOobject::AUTO_WRITE;
193  subFld->store();
194  }
195 }
196 
197 
198 template<class GeoField>
199 void subsetSurfaceFields
200 (
201  const fvMesh& mesh,
202  const fvMesh& subMesh,
203  const labelList& cellMap,
204  const labelList& faceMap,
205  const labelHashSet& addedPatches
206 )
207 {
208  const labelList patchMap(identity(mesh.boundaryMesh().size()));
209 
211  (
212  mesh.objectRegistry::lookupClass<GeoField>()
213  );
215  {
216  const GeoField& fld = *iter();
217 
218  Info<< "Mapping field " << fld.name() << endl;
219 
220  tmp<GeoField> tSubFld
221  (
223  (
224  fld,
225  subMesh,
226  patchMap,
227  cellMap,
228  faceMap
229  )
230  );
231 
232  // Hack: set value to 0 for introduced patches (since don't
233  // get initialised.
234  forAll(tSubFld().boundaryField(), patchi)
235  {
236  if (addedPatches.found(patchi))
237  {
238  tSubFld.ref().boundaryFieldRef()[patchi] ==
239  typename GeoField::value_type(Zero);
240  }
241  }
242 
243  // Store on subMesh
244  GeoField* subFld = tSubFld.ptr();
245  subFld->rename(fld.name());
246  subFld->writeOpt() = IOobject::AUTO_WRITE;
247  subFld->store();
248  }
249 }
250 
251 // Select all cells not in the region
252 labelList getNonRegionCells(const labelList& cellRegion, const label regionI)
253 {
254  DynamicList<label> nonRegionCells(cellRegion.size());
255  forAll(cellRegion, celli)
256  {
257  if (cellRegion[celli] != regionI)
258  {
259  nonRegionCells.append(celli);
260  }
261  }
262  return nonRegionCells.shrink();
263 }
264 
265 
266 void addToInterface
267 (
268  const polyMesh& mesh,
269  const label zoneID,
270  const label ownRegion,
271  const label neiRegion,
272  EdgeMap<Map<label>>& regionsToSize
273 )
274 {
276  (
277  min(ownRegion, neiRegion),
278  max(ownRegion, neiRegion)
279  );
280 
281  EdgeMap<Map<label>>::iterator iter = regionsToSize.find
282  (
283  interface
284  );
285 
286  if (iter != regionsToSize.end())
287  {
288  // Check if zone present
289  Map<label>::iterator zoneFnd = iter().find(zoneID);
290  if (zoneFnd != iter().end())
291  {
292  // Found zone. Increment count.
293  zoneFnd()++;
294  }
295  else
296  {
297  // New or no zone. Insert with count 1.
298  iter().insert(zoneID, 1);
299  }
300  }
301  else
302  {
303  // Create new interface of size 1.
304  Map<label> zoneToSize;
305  zoneToSize.insert(zoneID, 1);
306  regionsToSize.insert(interface, zoneToSize);
307  }
308 }
309 
310 
311 // Get region-region interface name and sizes.
312 // Returns interfaces as straight list for looping in identical order.
313 void getInterfaceSizes
314 (
315  const polyMesh& mesh,
316  const bool useFaceZones,
317  const labelList& cellRegion,
318  const wordList& regionNames,
319 
320  edgeList& interfaces,
321  List<Pair<word>>& interfaceNames,
322  labelList& interfaceSizes,
323  labelList& faceToInterface
324 )
325 {
326  // From region-region to faceZone (or -1) to number of faces.
327 
328  EdgeMap<Map<label>> regionsToSize;
329 
330 
331  // Internal faces
332  // ~~~~~~~~~~~~~~
333 
334  forAll(mesh.faceNeighbour(), facei)
335  {
336  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
337  label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
338 
339  if (ownRegion != neiRegion)
340  {
341  addToInterface
342  (
343  mesh,
344  (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
345  ownRegion,
346  neiRegion,
347  regionsToSize
348  );
349  }
350  }
351 
352  // Boundary faces
353  // ~~~~~~~~~~~~~~
354 
355  // Neighbour cellRegion.
356  labelList coupledRegion(mesh.nFaces()-mesh.nInternalFaces());
357 
358  forAll(coupledRegion, i)
359  {
360  label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
361  coupledRegion[i] = cellRegion[celli];
362  }
363  syncTools::swapBoundaryFaceList(mesh, coupledRegion);
364 
365  forAll(coupledRegion, i)
366  {
367  label facei = i+mesh.nInternalFaces();
368  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
369  label neiRegion = coupledRegion[i];
370 
371  if (ownRegion != neiRegion)
372  {
373  addToInterface
374  (
375  mesh,
376  (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
377  ownRegion,
378  neiRegion,
379  regionsToSize
380  );
381  }
382  }
383 
384 
385  if (Pstream::parRun())
386  {
387  if (Pstream::master())
388  {
389  // Receive and add to my sizes
390  for
391  (
392  int slave=Pstream::firstSlave();
393  slave<=Pstream::lastSlave();
394  slave++
395  )
396  {
397  IPstream fromSlave(Pstream::commsTypes::blocking, slave);
398 
399  EdgeMap<Map<label>> slaveSizes(fromSlave);
400 
401  forAllConstIter(EdgeMap<Map<label>>, slaveSizes, slaveIter)
402  {
403  EdgeMap<Map<label>>::iterator masterIter =
404  regionsToSize.find(slaveIter.key());
405 
406  if (masterIter != regionsToSize.end())
407  {
408  // Same inter-region
409  const Map<label>& slaveInfo = slaveIter();
410  Map<label>& masterInfo = masterIter();
411 
412  forAllConstIter(Map<label>, slaveInfo, iter)
413  {
414  label zoneID = iter.key();
415  label slaveSize = iter();
416 
417  Map<label>::iterator zoneFnd = masterInfo.find
418  (
419  zoneID
420  );
421  if (zoneFnd != masterInfo.end())
422  {
423  zoneFnd() += slaveSize;
424  }
425  else
426  {
427  masterInfo.insert(zoneID, slaveSize);
428  }
429  }
430  }
431  else
432  {
433  regionsToSize.insert(slaveIter.key(), slaveIter());
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::scatter(regionsToSize);
455 
456 
457 
458  // Now we have the global sizes of all inter-regions.
459  // Invert this on master and distribute.
460  label nInterfaces = 0;
461  forAllConstIter(EdgeMap<Map<label>>, regionsToSize, iter)
462  {
463  const Map<label>& info = iter();
464  nInterfaces += info.size();
465  }
466 
467  interfaces.setSize(nInterfaces);
468  interfaceNames.setSize(nInterfaces);
469  interfaceSizes.setSize(nInterfaces);
470  EdgeMap<Map<label>> regionsToInterface(nInterfaces);
471 
472  nInterfaces = 0;
473  forAllConstIter(EdgeMap<Map<label>>, regionsToSize, iter)
474  {
475  const edge& e = iter.key();
476  const word& name0 = regionNames[e[0]];
477  const word& name1 = regionNames[e[1]];
478 
479  const Map<label>& info = iter();
480  forAllConstIter(Map<label>, 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  // Now all processor have consistent interface information
519 
520  Pstream::scatter(interfaces);
521  Pstream::scatter(interfaceNames);
522  Pstream::scatter(interfaceSizes);
523  Pstream::scatter(regionsToInterface);
524 
525  // Mark all inter-region faces.
526  faceToInterface.setSize(mesh.nFaces(), -1);
527 
528  forAll(mesh.faceNeighbour(), facei)
529  {
530  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
531  label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
532 
533  if (ownRegion != neiRegion)
534  {
535  label zoneID = -1;
536  if (useFaceZones)
537  {
538  zoneID = mesh.faceZones().whichZone(facei);
539  }
540 
542  (
543  min(ownRegion, neiRegion),
544  max(ownRegion, neiRegion)
545  );
546 
547  faceToInterface[facei] = regionsToInterface[interface][zoneID];
548  }
549  }
550  forAll(coupledRegion, i)
551  {
552  label facei = i+mesh.nInternalFaces();
553  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
554  label neiRegion = coupledRegion[i];
555 
556  if (ownRegion != neiRegion)
557  {
558  label zoneID = -1;
559  if (useFaceZones)
560  {
561  zoneID = mesh.faceZones().whichZone(facei);
562  }
563 
565  (
566  min(ownRegion, neiRegion),
567  max(ownRegion, neiRegion)
568  );
569 
570  faceToInterface[facei] = regionsToInterface[interface][zoneID];
571  }
572  }
573 }
574 
575 
576 // Create mesh for region.
577 autoPtr<mapPolyMesh> createRegionMesh
578 (
579  const fvMesh& mesh,
580  // Region info
581  const labelList& cellRegion,
582  const label regionI,
583  const word& regionName,
584  // Interface info
585  const labelList& interfacePatches,
586  const labelList& faceToInterface,
587 
588  autoPtr<fvMesh>& newMesh
589 )
590 {
591  // Create dummy system/fv*
592  {
593  IOobject io
594  (
595  "fvSchemes",
596  mesh.time().system(),
597  regionName,
598  mesh,
601  false
602  );
603 
604  Info<< "Testing:" << io.objectPath() << endl;
605 
606  if (!io.typeHeaderOk<IOdictionary>(true))
607  // if (!exists(io.objectPath()))
608  {
609  Info<< "Writing dummy " << regionName/io.name() << endl;
610  dictionary dummyDict;
611  dictionary divDict;
612  dummyDict.add("divSchemes", divDict);
613  dictionary gradDict;
614  dummyDict.add("gradSchemes", gradDict);
615  dictionary laplDict;
616  dummyDict.add("laplacianSchemes", laplDict);
617 
618  IOdictionary(io, dummyDict).regIOobject::write();
619  }
620  }
621  {
622  IOobject io
623  (
624  "fvSolution",
625  mesh.time().system(),
626  regionName,
627  mesh,
630  false
631  );
632 
633  if (!io.typeHeaderOk<IOdictionary>(true))
634  //if (!exists(io.objectPath()))
635  {
636  Info<< "Writing dummy " << regionName/io.name() << endl;
637  dictionary dummyDict;
638  IOdictionary(io, dummyDict).regIOobject::write();
639  }
640  }
641 
642 
643  // Neighbour cellRegion.
644  labelList coupledRegion(mesh.nFaces()-mesh.nInternalFaces());
645 
646  forAll(coupledRegion, i)
647  {
648  label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
649  coupledRegion[i] = cellRegion[celli];
650  }
651  syncTools::swapBoundaryFaceList(mesh, coupledRegion);
652 
653 
654  // Topology change container. Start off from existing mesh.
655  polyTopoChange meshMod(mesh);
656 
657  // Cell remover engine
658  removeCells cellRemover(mesh);
659 
660  // Select all but region cells
661  labelList cellsToRemove(getNonRegionCells(cellRegion, regionI));
662 
663  // Find out which faces will get exposed. Note that this
664  // gets faces in mesh face order. So both regions will get same
665  // face in same order (important!)
666  labelList exposedFaces = cellRemover.getExposedFaces(cellsToRemove);
667 
668  labelList exposedPatchIDs(exposedFaces.size());
669  forAll(exposedFaces, i)
670  {
671  label facei = exposedFaces[i];
672  label interfacei = faceToInterface[facei];
673 
674  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
675  label neiRegion = -1;
676 
677  if (mesh.isInternalFace(facei))
678  {
679  neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
680  }
681  else
682  {
683  neiRegion = coupledRegion[facei-mesh.nInternalFaces()];
684  }
685 
686 
687  // Check which side is being kept - determines which of the two
688  // patches will be used.
689 
690  label otherRegion = -1;
691 
692  if (ownRegion == regionI && neiRegion != regionI)
693  {
694  otherRegion = neiRegion;
695  }
696  else if (ownRegion != regionI && neiRegion == regionI)
697  {
698  otherRegion = ownRegion;
699  }
700  else
701  {
703  << "Exposed face:" << facei
704  << " fc:" << mesh.faceCentres()[facei]
705  << " has owner region " << ownRegion
706  << " and neighbour region " << neiRegion
707  << " when handling region:" << regionI
708  << exit(FatalError);
709  }
710 
711  // Find the patch.
712  if (regionI < otherRegion)
713  {
714  exposedPatchIDs[i] = interfacePatches[interfacei];
715  }
716  else
717  {
718  exposedPatchIDs[i] = interfacePatches[interfacei]+1;
719  }
720  }
721 
722  // Remove faces
723  cellRemover.setRefinement
724  (
725  cellsToRemove,
726  exposedFaces,
727  exposedPatchIDs,
728  meshMod
729  );
730 
731  autoPtr<mapPolyMesh> map = meshMod.makeMesh
732  (
733  newMesh,
734  IOobject
735  (
736  regionName,
737  mesh.time().timeName(),
738  mesh.time(),
741  ),
742  mesh
743  );
744 
745  return map;
746 }
747 
748 
749 void createAndWriteRegion
750 (
751  const fvMesh& mesh,
752  const labelList& cellRegion,
753  const wordList& regionNames,
754  const bool prefixRegion,
755  const labelList& faceToInterface,
756  const labelList& interfacePatches,
757  const label regionI,
758  const word& newMeshInstance
759 )
760 {
761  Info<< "Creating mesh for region " << regionI
762  << ' ' << regionNames[regionI] << endl;
763 
764  autoPtr<fvMesh> newMesh;
765  autoPtr<mapPolyMesh> map = createRegionMesh
766  (
767  mesh,
768  cellRegion,
769  regionI,
770  regionNames[regionI],
771  interfacePatches,
772  faceToInterface,
773  newMesh
774  );
775 
776 
777  // Make map of all added patches
778  labelHashSet addedPatches(2*interfacePatches.size());
779  forAll(interfacePatches, interfacei)
780  {
781  addedPatches.insert(interfacePatches[interfacei]);
782  addedPatches.insert(interfacePatches[interfacei]+1);
783  }
784 
785 
786  Info<< "Mapping fields" << endl;
787 
788  // Map existing fields
789  newMesh().updateMesh(map());
790 
791  // Add subsetted fields
792  subsetVolFields<volScalarField>
793  (
794  mesh,
795  newMesh(),
796  map().cellMap(),
797  map().faceMap(),
798  addedPatches
799  );
800  subsetVolFields<volVectorField>
801  (
802  mesh,
803  newMesh(),
804  map().cellMap(),
805  map().faceMap(),
806  addedPatches
807  );
808  subsetVolFields<volSphericalTensorField>
809  (
810  mesh,
811  newMesh(),
812  map().cellMap(),
813  map().faceMap(),
814  addedPatches
815  );
816  subsetVolFields<volSymmTensorField>
817  (
818  mesh,
819  newMesh(),
820  map().cellMap(),
821  map().faceMap(),
822  addedPatches
823  );
824  subsetVolFields<volTensorField>
825  (
826  mesh,
827  newMesh(),
828  map().cellMap(),
829  map().faceMap(),
830  addedPatches
831  );
832 
833  subsetSurfaceFields<surfaceScalarField>
834  (
835  mesh,
836  newMesh(),
837  map().cellMap(),
838  map().faceMap(),
839  addedPatches
840  );
841  subsetSurfaceFields<surfaceVectorField>
842  (
843  mesh,
844  newMesh(),
845  map().cellMap(),
846  map().faceMap(),
847  addedPatches
848  );
849  subsetSurfaceFields<surfaceSphericalTensorField>
850  (
851  mesh,
852  newMesh(),
853  map().cellMap(),
854  map().faceMap(),
855  addedPatches
856  );
857  subsetSurfaceFields<surfaceSymmTensorField>
858  (
859  mesh,
860  newMesh(),
861  map().cellMap(),
862  map().faceMap(),
863  addedPatches
864  );
865  subsetSurfaceFields<surfaceTensorField>
866  (
867  mesh,
868  newMesh(),
869  map().cellMap(),
870  map().faceMap(),
871  addedPatches
872  );
873 
874 
875  const polyBoundaryMesh& newPatches = newMesh().boundaryMesh();
876  newPatches.checkParallelSync(true);
877 
878  // Delete empty patches
879  // ~~~~~~~~~~~~~~~~~~~~
880 
881  // Create reordering list to move patches-to-be-deleted to end
882  labelList oldToNew(newPatches.size(), -1);
883  DynamicList<label> sharedPatches(newPatches.size());
884  label newI = 0;
885 
886  Info<< "Deleting empty patches" << endl;
887 
888  // Assumes all non-proc boundaries are on all processors!
889  forAll(newPatches, patchi)
890  {
891  const polyPatch& pp = newPatches[patchi];
892 
893  if (!isA<processorPolyPatch>(pp))
894  {
895  if (returnReduce(pp.size(), sumOp<label>()) > 0)
896  {
897  oldToNew[patchi] = newI;
898  if (!addedPatches.found(patchi))
899  {
900  sharedPatches.append(newI);
901  }
902  newI++;
903  }
904  }
905  }
906 
907  // Same for processor patches (but need no reduction)
908  forAll(newPatches, patchi)
909  {
910  const polyPatch& pp = newPatches[patchi];
911 
912  if (isA<processorPolyPatch>(pp) && pp.size())
913  {
914  oldToNew[patchi] = newI++;
915  }
916  }
917 
918  const label nNewPatches = newI;
919 
920  // Move all deleteable patches to the end
921  forAll(oldToNew, patchi)
922  {
923  if (oldToNew[patchi] == -1)
924  {
925  oldToNew[patchi] = newI++;
926  }
927  }
928 
929  //reorderPatches(newMesh(), oldToNew, nNewPatches);
930  fvMeshTools::reorderPatches(newMesh(), oldToNew, nNewPatches, true);
931 
932  // Rename shared patches with region name
933  if (prefixRegion)
934  {
935  Info<< "Prefixing patches with region name" << endl;
936 
937  renamePatches(newMesh(), regionNames[regionI], sharedPatches);
938  }
939 
940 
941  Info<< "Writing new mesh" << endl;
942 
943  newMesh().setInstance(newMeshInstance);
944  newMesh().write();
945 
946  // Write addressing files like decomposePar
947  Info<< "Writing addressing to base mesh" << endl;
948 
949  labelIOList pointProcAddressing
950  (
951  IOobject
952  (
953  "pointRegionAddressing",
954  newMesh().facesInstance(),
955  newMesh().meshSubDir,
956  newMesh(),
957  IOobject::NO_READ,
958  IOobject::NO_WRITE,
959  false
960  ),
961  map().pointMap()
962  );
963  Info<< "Writing map " << pointProcAddressing.name()
964  << " from region" << regionI
965  << " points back to base mesh." << endl;
966  pointProcAddressing.write();
967 
969  (
970  IOobject
971  (
972  "faceRegionAddressing",
973  newMesh().facesInstance(),
974  newMesh().meshSubDir,
975  newMesh(),
976  IOobject::NO_READ,
977  IOobject::NO_WRITE,
978  false
979  ),
980  newMesh().nFaces()
981  );
982  forAll(faceProcAddressing, facei)
983  {
984  // face + turning index. (see decomposePar)
985  // Is the face pointing in the same direction?
986  label oldFacei = map().faceMap()[facei];
987 
988  if
989  (
990  map().cellMap()[newMesh().faceOwner()[facei]]
991  == mesh.faceOwner()[oldFacei]
992  )
993  {
994  faceProcAddressing[facei] = oldFacei+1;
995  }
996  else
997  {
998  faceProcAddressing[facei] = -(oldFacei+1);
999  }
1000  }
1001  Info<< "Writing map " << faceProcAddressing.name()
1002  << " from region" << regionI
1003  << " faces back to base mesh." << endl;
1004  faceProcAddressing.write();
1005 
1006  labelIOList cellProcAddressing
1007  (
1008  IOobject
1009  (
1010  "cellRegionAddressing",
1011  newMesh().facesInstance(),
1012  newMesh().meshSubDir,
1013  newMesh(),
1014  IOobject::NO_READ,
1015  IOobject::NO_WRITE,
1016  false
1017  ),
1018  map().cellMap()
1019  );
1020  Info<< "Writing map " <<cellProcAddressing.name()
1021  << " from region" << regionI
1022  << " cells back to base mesh." << endl;
1023  cellProcAddressing.write();
1024 
1025  labelIOList boundaryProcAddressing
1026  (
1027  IOobject
1028  (
1029  "boundaryRegionAddressing",
1030  newMesh().facesInstance(),
1031  newMesh().meshSubDir,
1032  newMesh(),
1033  IOobject::NO_READ,
1034  IOobject::NO_WRITE,
1035  false
1036  ),
1037  labelList(nNewPatches, -1)
1038  );
1039  forAll(oldToNew, i)
1040  {
1041  if (!addedPatches.found(i))
1042  {
1043  label newI = oldToNew[i];
1044  if (newI >= 0 && newI < nNewPatches)
1045  {
1046  boundaryProcAddressing[oldToNew[i]] = i;
1047  }
1048  }
1049  }
1050  Info<< "Writing map " << boundaryProcAddressing.name()
1051  << " from region" << regionI
1052  << " boundary back to base mesh." << endl;
1053  boundaryProcAddressing.write();
1054 }
1055 
1056 
1057 // Create for every region-region interface with non-zero size two patches.
1058 // First one is for minimumregion to maximumregion.
1059 // Note that patches get created in same order on all processors (if parallel)
1060 // since looping over synchronised 'interfaces'.
1061 labelList addRegionPatches
1062 (
1063  fvMesh& mesh,
1064  const wordList& regionNames,
1065  const edgeList& interfaces,
1066  const List<Pair<word>>& interfaceNames
1067 )
1068 {
1069  Info<< nl << "Adding patches" << nl << endl;
1070 
1071  labelList interfacePatches(interfaces.size());
1072 
1073  forAll(interfaces, interI)
1074  {
1075  const edge& e = interfaces[interI];
1076  const Pair<word>& names = interfaceNames[interI];
1077 
1078  //Info<< "For interface " << interI
1079  // << " between regions " << e
1080  // << " trying to add patches " << names << endl;
1081 
1082 
1083  mappedWallPolyPatch patch1
1084  (
1085  names[0],
1086  0, // overridden
1087  0, // overridden
1088  0, // overridden
1089  regionNames[e[1]], // sampleRegion
1091  names[1], // samplePatch
1092  point::zero, // offset
1093  mesh.boundaryMesh()
1094  );
1095 
1096  interfacePatches[interI] = fvMeshTools::addPatch
1097  (
1098  mesh,
1099  patch1,
1100  dictionary(), //optional per field value
1102  true //validBoundary
1103  );
1104 
1105  mappedWallPolyPatch patch2
1106  (
1107  names[1],
1108  0,
1109  0,
1110  0,
1111  regionNames[e[0]], // sampleRegion
1113  names[0],
1114  point::zero, // offset
1115  mesh.boundaryMesh()
1116  );
1118  (
1119  mesh,
1120  patch2,
1121  dictionary(), //optional per field value
1123  true //validBoundary
1124  );
1125 
1126  Info<< "For interface between region " << regionNames[e[0]]
1127  << " and " << regionNames[e[1]] << " added patches" << endl
1128  << " " << interfacePatches[interI]
1129  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]].name()
1130  << endl
1131  << " " << interfacePatches[interI]+1
1132  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]+1].name()
1133  << endl;
1134  }
1135  return interfacePatches;
1136 }
1137 
1138 
1139 // Find region that covers most of cell zone
1140 label findCorrespondingRegion
1141 (
1142  const labelList& existingZoneID, // per cell the (unique) zoneID
1143  const labelList& cellRegion,
1144  const label nCellRegions,
1145  const label zoneI,
1146  const label minOverlapSize
1147 )
1148 {
1149  // Per region the number of cells in zoneI
1150  labelList cellsInZone(nCellRegions, 0);
1151 
1152  forAll(cellRegion, celli)
1153  {
1154  if (existingZoneID[celli] == zoneI)
1155  {
1156  cellsInZone[cellRegion[celli]]++;
1157  }
1158  }
1159 
1161  Pstream::listCombineScatter(cellsInZone);
1162 
1163  // Pick region with largest overlap of zoneI
1164  label regionI = findMax(cellsInZone);
1165 
1166 
1167  if (cellsInZone[regionI] < minOverlapSize)
1168  {
1169  // Region covers too little of zone. Not good enough.
1170  regionI = -1;
1171  }
1172  else
1173  {
1174  // Check that region contains no cells that aren't in cellZone.
1175  forAll(cellRegion, celli)
1176  {
1177  if (cellRegion[celli] == regionI && existingZoneID[celli] != zoneI)
1178  {
1179  // celli in regionI but not in zoneI
1180  regionI = -1;
1181  break;
1182  }
1183  }
1184  // If one in error, all should be in error. Note that branch gets taken
1185  // on all procs.
1186  reduce(regionI, minOp<label>());
1187  }
1188 
1189  return regionI;
1190 }
1191 
1192 
1193 // Get zone per cell
1194 // - non-unique zoning
1195 // - coupled zones
1196 void getZoneID
1197 (
1198  const polyMesh& mesh,
1199  const cellZoneMesh& cellZones,
1200  labelList& zoneID,
1201  labelList& neiZoneID
1202 )
1203 {
1204  // Existing zoneID
1205  zoneID.setSize(mesh.nCells());
1206  zoneID = -1;
1207 
1208  forAll(cellZones, zoneI)
1209  {
1210  const cellZone& cz = cellZones[zoneI];
1211 
1212  forAll(cz, i)
1213  {
1214  label celli = cz[i];
1215  if (zoneID[celli] == -1)
1216  {
1217  zoneID[celli] = zoneI;
1218  }
1219  else
1220  {
1222  << "Cell " << celli << " with cell centre "
1223  << mesh.cellCentres()[celli]
1224  << " is multiple zones. This is not allowed." << endl
1225  << "It is in zone " << cellZones[zoneID[celli]].name()
1226  << " and in zone " << cellZones[zoneI].name()
1227  << exit(FatalError);
1228  }
1229  }
1230  }
1231 
1232  // Neighbour zoneID.
1233  neiZoneID.setSize(mesh.nFaces()-mesh.nInternalFaces());
1234 
1235  forAll(neiZoneID, i)
1236  {
1237  neiZoneID[i] = zoneID[mesh.faceOwner()[i+mesh.nInternalFaces()]];
1238  }
1239  syncTools::swapBoundaryFaceList(mesh, neiZoneID);
1240 }
1241 
1242 
1243 void matchRegions
1244 (
1245  const bool sloppyCellZones,
1246  const polyMesh& mesh,
1247 
1248  const label nCellRegions,
1249  const labelList& cellRegion,
1250 
1251  labelList& regionToZone,
1252  wordList& regionNames,
1253  labelList& zoneToRegion
1254 )
1255 {
1256  const cellZoneMesh& cellZones = mesh.cellZones();
1257 
1258  regionToZone.setSize(nCellRegions, -1);
1259  regionNames.setSize(nCellRegions);
1260  zoneToRegion.setSize(cellZones.size(), -1);
1261 
1262  // Get current per cell zoneID
1263  labelList zoneID(mesh.nCells(), -1);
1264  labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1265  getZoneID(mesh, cellZones, zoneID, neiZoneID);
1266 
1267  // Sizes per cellzone
1268  labelList zoneSizes(cellZones.size(), 0);
1269  {
1270  List<wordList> zoneNames(Pstream::nProcs());
1271  zoneNames[Pstream::myProcNo()] = cellZones.names();
1272  Pstream::gatherList(zoneNames);
1273  Pstream::scatterList(zoneNames);
1274 
1275  forAll(zoneNames, proci)
1276  {
1277  if (zoneNames[proci] != zoneNames[0])
1278  {
1280  << "cellZones not synchronised across processors." << endl
1281  << "Master has cellZones " << zoneNames[0] << endl
1282  << "Processor " << proci
1283  << " has cellZones " << zoneNames[proci]
1284  << exit(FatalError);
1285  }
1286  }
1287 
1288  forAll(cellZones, zoneI)
1289  {
1290  zoneSizes[zoneI] = returnReduce
1291  (
1292  cellZones[zoneI].size(),
1293  sumOp<label>()
1294  );
1295  }
1296  }
1297 
1298 
1299  if (sloppyCellZones)
1300  {
1301  Info<< "Trying to match regions to existing cell zones;"
1302  << " region can be subset of cell zone." << nl << endl;
1303 
1304  forAll(cellZones, zoneI)
1305  {
1306  label regionI = findCorrespondingRegion
1307  (
1308  zoneID,
1309  cellRegion,
1310  nCellRegions,
1311  zoneI,
1312  label(0.5*zoneSizes[zoneI]) // minimum overlap
1313  );
1314 
1315  if (regionI != -1)
1316  {
1317  Info<< "Sloppily matched region " << regionI
1318  //<< " size " << regionSizes[regionI]
1319  << " to zone " << zoneI << " size " << zoneSizes[zoneI]
1320  << endl;
1321  zoneToRegion[zoneI] = regionI;
1322  regionToZone[regionI] = zoneI;
1323  regionNames[regionI] = cellZones[zoneI].name();
1324  }
1325  }
1326  }
1327  else
1328  {
1329  Info<< "Trying to match regions to existing cell zones." << nl << endl;
1330 
1331  forAll(cellZones, zoneI)
1332  {
1333  label regionI = findCorrespondingRegion
1334  (
1335  zoneID,
1336  cellRegion,
1337  nCellRegions,
1338  zoneI,
1339  1 // minimum overlap
1340  );
1341 
1342  if (regionI != -1)
1343  {
1344  zoneToRegion[zoneI] = regionI;
1345  regionToZone[regionI] = zoneI;
1346  regionNames[regionI] = cellZones[zoneI].name();
1347  }
1348  }
1349  }
1350  // Allocate region names for unmatched regions.
1351  forAll(regionToZone, regionI)
1352  {
1353  if (regionToZone[regionI] == -1)
1354  {
1355  regionNames[regionI] = "domain" + Foam::name(regionI);
1356  }
1357  }
1358 }
1359 
1360 
1361 void writeCellToRegion(const fvMesh& mesh, const labelList& cellRegion)
1362 {
1363  // Write to manual decomposition option
1364  {
1365  labelIOList cellToRegion
1366  (
1367  IOobject
1368  (
1369  "cellToRegion",
1370  mesh.facesInstance(),
1371  mesh,
1374  false
1375  ),
1376  cellRegion
1377  );
1378  cellToRegion.write();
1379 
1380  Info<< "Writing region per cell file (for manual decomposition) to "
1381  << cellToRegion.objectPath() << nl << endl;
1382  }
1383  // Write for postprocessing
1384  {
1385  volScalarField cellToRegion
1386  (
1387  IOobject
1388  (
1389  "cellToRegion",
1390  mesh.time().timeName(),
1391  mesh,
1394  false
1395  ),
1396  mesh,
1397  dimensionedScalar("zero", dimless, 0),
1398  zeroGradientFvPatchScalarField::typeName
1399  );
1400  forAll(cellRegion, celli)
1401  {
1402  cellToRegion[celli] = cellRegion[celli];
1403  }
1404  cellToRegion.write();
1405 
1406  Info<< "Writing region per cell as volScalarField to "
1407  << cellToRegion.objectPath() << nl << endl;
1408  }
1409 }
1410 
1411 
1412 
1413 
1414 int main(int argc, char *argv[])
1415 {
1417  (
1418  "splits mesh into multiple regions (detected by walking across faces)"
1419  );
1420  #include "addRegionOption.H"
1421  #include "addOverwriteOption.H"
1423  (
1424  "cellZones",
1425  "additionally split cellZones off into separate regions"
1426  );
1428  (
1429  "cellZonesOnly",
1430  "use cellZones only to split mesh into regions; do not use walking"
1431  );
1433  (
1434  "cellZonesFileOnly",
1435  "file",
1436  "like -cellZonesOnly, but use specified file"
1437  );
1439  (
1440  "blockedFaces",
1441  "faceSet",
1442  "specify additional region boundaries that walking does not cross"
1443  );
1445  (
1446  "makeCellZones",
1447  "place cells into cellZones instead of splitting mesh"
1448  );
1450  (
1451  "largestOnly",
1452  "only write largest region"
1453  );
1455  (
1456  "insidePoint",
1457  "point",
1458  "only write region containing point"
1459  );
1461  (
1462  "detectOnly",
1463  "do not write mesh"
1464  );
1466  (
1467  "sloppyCellZones",
1468  "try to match heuristically regions to existing cell zones"
1469  );
1471  (
1472  "useFaceZones",
1473  "use faceZones to patch inter-region faces instead of single patch"
1474  );
1476  (
1477  "prefixRegion",
1478  "prefix region name to all patches, not just coupling patches"
1479  );
1480 
1481  #include "setRootCase.H"
1482  #include "createTime.H"
1483  runTime.functionObjects().off();
1484  #include "createNamedMesh.H"
1485  const word oldInstance = mesh.pointsInstance();
1486 
1487  word blockedFacesName;
1488  if (args.optionReadIfPresent("blockedFaces", blockedFacesName))
1489  {
1490  Info<< "Reading blocked internal faces from faceSet "
1491  << blockedFacesName << nl << endl;
1492  }
1493 
1494  const bool makeCellZones = args.optionFound("makeCellZones");
1495  const bool largestOnly = args.optionFound("largestOnly");
1496  const bool insidePoint = args.optionFound("insidePoint");
1497  const bool useCellZones = args.optionFound("cellZones");
1498  const bool useCellZonesOnly = args.optionFound("cellZonesOnly");
1499  const bool useCellZonesFile = args.optionFound("cellZonesFileOnly");
1500  const bool overwrite = args.optionFound("overwrite");
1501  const bool detectOnly = args.optionFound("detectOnly");
1502  const bool sloppyCellZones = args.optionFound("sloppyCellZones");
1503  const bool useFaceZones = args.optionFound("useFaceZones");
1504  const bool prefixRegion = args.optionFound("prefixRegion");
1505 
1506 
1507  if
1508  (
1509  (useCellZonesOnly || useCellZonesFile)
1510  && (useCellZones || blockedFacesName.size())
1511  )
1512  {
1514  << "You cannot specify both -cellZonesOnly or -cellZonesFileOnly"
1515  << " (which specify complete split)"
1516  << " in combination with -blockedFaces or -cellZones"
1517  << " (which imply a split based on topology)"
1518  << exit(FatalError);
1519  }
1520 
1521 
1522  if (useFaceZones)
1523  {
1524  Info<< "Using current faceZones to divide inter-region interfaces"
1525  << " into multiple patches."
1526  << nl << endl;
1527  }
1528  else
1529  {
1530  Info<< "Creating single patch per inter-region interface."
1531  << nl << endl;
1532  }
1533 
1534 
1535 
1536  if (insidePoint && largestOnly)
1537  {
1539  << "You cannot specify both -largestOnly"
1540  << " (keep region with most cells)"
1541  << " and -insidePoint (keep region containing point)"
1542  << exit(FatalError);
1543  }
1544 
1545 
1546  const cellZoneMesh& cellZones = mesh.cellZones();
1547 
1548  // Existing zoneID
1549  labelList zoneID(mesh.nCells(), -1);
1550  // Neighbour zoneID.
1551  labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1552  getZoneID(mesh, cellZones, zoneID, neiZoneID);
1553 
1554 
1555 
1556  // Determine per cell the region it belongs to
1557  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1558 
1559  // cellRegion is the labelList with the region per cell.
1560  labelList cellRegion;
1561  // Region per zone
1562  labelList regionToZone;
1563  // Name of region
1564  wordList regionNames;
1565  // Zone to region
1566  labelList zoneToRegion;
1567 
1568  label nCellRegions = 0;
1569  if (useCellZonesOnly)
1570  {
1571  Info<< "Using current cellZones to split mesh into regions."
1572  << " This requires all"
1573  << " cells to be in one and only one cellZone." << nl << endl;
1574 
1575  label unzonedCelli = findIndex(zoneID, -1);
1576  if (unzonedCelli != -1)
1577  {
1579  << "For the cellZonesOnly option all cells "
1580  << "have to be in a cellZone." << endl
1581  << "Cell " << unzonedCelli
1582  << " at" << mesh.cellCentres()[unzonedCelli]
1583  << " is not in a cellZone. There might be more unzoned cells."
1584  << exit(FatalError);
1585  }
1586  cellRegion = zoneID;
1587  nCellRegions = gMax(cellRegion)+1;
1588  regionToZone.setSize(nCellRegions);
1589  regionNames.setSize(nCellRegions);
1590  zoneToRegion.setSize(cellZones.size(), -1);
1591  for (label regionI = 0; regionI < nCellRegions; regionI++)
1592  {
1593  regionToZone[regionI] = regionI;
1594  zoneToRegion[regionI] = regionI;
1595  regionNames[regionI] = cellZones[regionI].name();
1596  }
1597  }
1598  else if (useCellZonesFile)
1599  {
1600  const word zoneFile = args.option("cellZonesFileOnly");
1601  Info<< "Reading split from cellZones file " << zoneFile << endl
1602  << "This requires all"
1603  << " cells to be in one and only one cellZone." << nl << endl;
1604 
1605  cellZoneMesh newCellZones
1606  (
1607  IOobject
1608  (
1609  zoneFile,
1610  mesh.facesInstance(),
1612  mesh,
1615  false
1616  ),
1617  mesh
1618  );
1619 
1620  labelList newZoneID(mesh.nCells(), -1);
1621  labelList newNeiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1622  getZoneID(mesh, newCellZones, newZoneID, newNeiZoneID);
1623 
1624  label unzonedCelli = findIndex(newZoneID, -1);
1625  if (unzonedCelli != -1)
1626  {
1628  << "For the cellZonesFileOnly option all cells "
1629  << "have to be in a cellZone." << endl
1630  << "Cell " << unzonedCelli
1631  << " at" << mesh.cellCentres()[unzonedCelli]
1632  << " is not in a cellZone. There might be more unzoned cells."
1633  << exit(FatalError);
1634  }
1635  cellRegion = newZoneID;
1636  nCellRegions = gMax(cellRegion)+1;
1637  zoneToRegion.setSize(newCellZones.size(), -1);
1638  regionToZone.setSize(nCellRegions);
1639  regionNames.setSize(nCellRegions);
1640  for (label regionI = 0; regionI < nCellRegions; regionI++)
1641  {
1642  regionToZone[regionI] = regionI;
1643  zoneToRegion[regionI] = regionI;
1644  regionNames[regionI] = newCellZones[regionI].name();
1645  }
1646  }
1647  else
1648  {
1649  // Determine connected regions
1650  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1651 
1652  // Mark additional faces that are blocked
1653  boolList blockedFace;
1654 
1655  // Read from faceSet
1656  if (blockedFacesName.size())
1657  {
1658  faceSet blockedFaceSet(mesh, blockedFacesName);
1659  Info<< "Read "
1660  << returnReduce(blockedFaceSet.size(), sumOp<label>())
1661  << " blocked faces from set " << blockedFacesName << nl << endl;
1662 
1663  blockedFace.setSize(mesh.nFaces(), false);
1664 
1665  forAllConstIter(faceSet, blockedFaceSet, iter)
1666  {
1667  blockedFace[iter.key()] = true;
1668  }
1669  }
1670 
1671  // Imply from differing cellZones
1672  if (useCellZones)
1673  {
1674  blockedFace.setSize(mesh.nFaces(), false);
1675 
1676  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1677  {
1678  label own = mesh.faceOwner()[facei];
1679  label nei = mesh.faceNeighbour()[facei];
1680 
1681  if (zoneID[own] != zoneID[nei])
1682  {
1683  blockedFace[facei] = true;
1684  }
1685  }
1686 
1687  // Different cellZones on either side of processor patch.
1688  forAll(neiZoneID, i)
1689  {
1690  label facei = i+mesh.nInternalFaces();
1691 
1692  if (zoneID[mesh.faceOwner()[facei]] != neiZoneID[i])
1693  {
1694  blockedFace[facei] = true;
1695  }
1696  }
1697  }
1698 
1699  // Do a topological walk to determine regions
1700  regionSplit regions(mesh, blockedFace);
1701  nCellRegions = regions.nRegions();
1702  cellRegion.transfer(regions);
1703 
1704  // Make up region names. If possible match them to existing zones.
1705  matchRegions
1706  (
1707  sloppyCellZones,
1708  mesh,
1709  nCellRegions,
1710  cellRegion,
1711 
1712  regionToZone,
1713  regionNames,
1714  zoneToRegion
1715  );
1716 
1717  // Override any default region names if single region selected
1718  if (largestOnly || insidePoint)
1719  {
1720  forAll(regionToZone, regionI)
1721  {
1722  if (regionToZone[regionI] == -1)
1723  {
1724  if (overwrite)
1725  {
1726  regionNames[regionI] = polyMesh::defaultRegion;
1727  }
1728  else if (insidePoint)
1729  {
1730  regionNames[regionI] = "insidePoint";
1731  }
1732  else if (largestOnly)
1733  {
1734  regionNames[regionI] = "largestOnly";
1735  }
1736  }
1737  }
1738  }
1739  }
1740 
1741  Info<< endl << "Number of regions:" << nCellRegions << nl << endl;
1742 
1743 
1744  // Write decomposition to file
1745  writeCellToRegion(mesh, cellRegion);
1746 
1747 
1748 
1749  // Sizes per region
1750  // ~~~~~~~~~~~~~~~~
1751 
1752  labelList regionSizes(nCellRegions, 0);
1753 
1754  forAll(cellRegion, celli)
1755  {
1756  regionSizes[cellRegion[celli]]++;
1757  }
1758  forAll(regionSizes, regionI)
1759  {
1760  reduce(regionSizes[regionI], sumOp<label>());
1761  }
1762 
1763  Info<< "Region\tCells" << nl
1764  << "------\t-----" << endl;
1765 
1766  forAll(regionSizes, regionI)
1767  {
1768  Info<< regionI << '\t' << regionSizes[regionI] << nl;
1769  }
1770  Info<< endl;
1771 
1772 
1773 
1774  // Print region to zone
1775  Info<< "Region\tZone\tName" << nl
1776  << "------\t----\t----" << endl;
1777  forAll(regionToZone, regionI)
1778  {
1779  Info<< regionI << '\t' << regionToZone[regionI] << '\t'
1780  << regionNames[regionI] << nl;
1781  }
1782  Info<< endl;
1783 
1784 
1785 
1786  // Since we're going to mess with patches and zones make sure all
1787  // is synchronised
1788  mesh.boundaryMesh().checkParallelSync(true);
1789  mesh.faceZones().checkParallelSync(true);
1790 
1791 
1792  // Interfaces
1793  // ----------
1794  // per interface:
1795  // - the two regions on either side
1796  // - the name
1797  // - the (global) size
1798  edgeList interfaces;
1799  List<Pair<word>> interfaceNames;
1800  labelList interfaceSizes;
1801  // per face the interface
1802  labelList faceToInterface;
1803 
1804  getInterfaceSizes
1805  (
1806  mesh,
1807  useFaceZones,
1808  cellRegion,
1809  regionNames,
1810 
1811  interfaces,
1812  interfaceNames,
1813  interfaceSizes,
1814  faceToInterface
1815  );
1816 
1817  Info<< "Sizes of interfaces between regions:" << nl << nl
1818  << "Interface\tRegion\tRegion\tFaces" << nl
1819  << "---------\t------\t------\t-----" << endl;
1820 
1821  forAll(interfaces, interI)
1822  {
1823  const edge& e = interfaces[interI];
1824 
1825  Info<< interI
1826  << "\t\t" << e[0] << '\t' << e[1]
1827  << '\t' << interfaceSizes[interI] << nl;
1828  }
1829  Info<< endl;
1830 
1831 
1832  if (detectOnly)
1833  {
1834  return 0;
1835  }
1836 
1837 
1838  // Read objects in time directory
1839  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1840 
1841  IOobjectList objects(mesh, runTime.timeName());
1842 
1843  // Read vol fields.
1844 
1845  PtrList<volScalarField> vsFlds;
1846  ReadFields(mesh, objects, vsFlds);
1847 
1848  PtrList<volVectorField> vvFlds;
1849  ReadFields(mesh, objects, vvFlds);
1850 
1852  ReadFields(mesh, objects, vstFlds);
1853 
1854  PtrList<volSymmTensorField> vsymtFlds;
1855  ReadFields(mesh, objects, vsymtFlds);
1856 
1857  PtrList<volTensorField> vtFlds;
1858  ReadFields(mesh, objects, vtFlds);
1859 
1860  // Read surface fields.
1861 
1863  ReadFields(mesh, objects, ssFlds);
1864 
1866  ReadFields(mesh, objects, svFlds);
1867 
1869  ReadFields(mesh, objects, sstFlds);
1870 
1872  ReadFields(mesh, objects, ssymtFlds);
1873 
1875  ReadFields(mesh, objects, stFlds);
1876 
1877  Info<< endl;
1878 
1879 
1880  // Remove any demand-driven fields ('S', 'V' etc)
1881  mesh.clearOut();
1882 
1883 
1884  if (nCellRegions == 1)
1885  {
1886  Info<< "Only one region. Doing nothing." << endl;
1887  }
1888  else if (makeCellZones)
1889  {
1890  Info<< "Putting cells into cellZones instead of splitting mesh."
1891  << endl;
1892 
1893  // Check if region overlaps with existing zone. If so keep.
1894 
1895  for (label regionI = 0; regionI < nCellRegions; regionI++)
1896  {
1897  label zoneI = regionToZone[regionI];
1898 
1899  if (zoneI != -1)
1900  {
1901  Info<< " Region " << regionI << " : corresponds to existing"
1902  << " cellZone "
1903  << zoneI << ' ' << cellZones[zoneI].name() << endl;
1904  }
1905  else
1906  {
1907  // Create new cellZone.
1908  labelList regionCells = findIndices(cellRegion, regionI);
1909 
1910  word zoneName = "region" + Foam::name(regionI);
1911 
1912  zoneI = cellZones.findZoneID(zoneName);
1913 
1914  if (zoneI == -1)
1915  {
1916  zoneI = cellZones.size();
1917  mesh.cellZones().setSize(zoneI+1);
1918  mesh.cellZones().set
1919  (
1920  zoneI,
1921  new cellZone
1922  (
1923  zoneName, //name
1924  regionCells, //addressing
1925  zoneI, //index
1926  cellZones //cellZoneMesh
1927  )
1928  );
1929  }
1930  else
1931  {
1932  mesh.cellZones()[zoneI].clearAddressing();
1933  mesh.cellZones()[zoneI] = regionCells;
1934  }
1935  Info<< " Region " << regionI << " : created new cellZone "
1936  << zoneI << ' ' << cellZones[zoneI].name() << endl;
1937  }
1938  }
1940 
1941  if (!overwrite)
1942  {
1943  runTime++;
1944  mesh.setInstance(runTime.timeName());
1945  }
1946  else
1947  {
1948  mesh.setInstance(oldInstance);
1949  }
1950 
1951  Info<< "Writing cellZones as new mesh to time " << runTime.timeName()
1952  << nl << endl;
1953 
1954  mesh.write();
1955 
1956 
1957  // Write cellSets for convenience
1958  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1959 
1960  Info<< "Writing cellSets corresponding to cellZones." << nl << endl;
1961 
1962  forAll(cellZones, zoneI)
1963  {
1964  const cellZone& cz = cellZones[zoneI];
1965 
1966  cellSet(mesh, cz.name(), cz).write();
1967  }
1968  }
1969  else
1970  {
1971  // Add patches for interfaces
1972  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
1973 
1974  // Add all possible patches. Empty ones get filtered later on.
1975  Info<< nl << "Adding patches" << nl << endl;
1976 
1977  labelList interfacePatches
1978  (
1979  addRegionPatches
1980  (
1981  mesh,
1982  regionNames,
1983  interfaces,
1984  interfaceNames
1985  )
1986  );
1987 
1988 
1989  if (!overwrite)
1990  {
1991  runTime++;
1992  }
1993 
1994 
1995  // Create regions
1996  // ~~~~~~~~~~~~~~
1997 
1998  if (insidePoint)
1999  {
2000  const point insidePoint = args.optionRead<point>("insidePoint");
2001 
2002  label regionI = -1;
2003 
2004  (void)mesh.tetBasePtIs();
2005 
2006  label celli = mesh.findCell(insidePoint);
2007 
2008  Info<< nl << "Found point " << insidePoint << " in cell " << celli
2009  << endl;
2010 
2011  if (celli != -1)
2012  {
2013  regionI = cellRegion[celli];
2014  }
2015 
2016  reduce(regionI, maxOp<label>());
2017 
2018  Info<< nl
2019  << "Subsetting region " << regionI
2020  << " containing point " << insidePoint << endl;
2021 
2022  if (regionI == -1)
2023  {
2025  << "Point " << insidePoint
2026  << " is not inside the mesh." << nl
2027  << "Bounding box of the mesh:" << mesh.bounds()
2028  << exit(FatalError);
2029  }
2030 
2031  createAndWriteRegion
2032  (
2033  mesh,
2034  cellRegion,
2035  regionNames,
2036  prefixRegion,
2037  faceToInterface,
2038  interfacePatches,
2039  regionI,
2040  (overwrite ? oldInstance : runTime.timeName())
2041  );
2042  }
2043  else if (largestOnly)
2044  {
2045  label regionI = findMax(regionSizes);
2046 
2047  Info<< nl
2048  << "Subsetting region " << regionI
2049  << " of size " << regionSizes[regionI]
2050  << " as named region " << regionNames[regionI] << endl;
2051 
2052  createAndWriteRegion
2053  (
2054  mesh,
2055  cellRegion,
2056  regionNames,
2057  prefixRegion,
2058  faceToInterface,
2059  interfacePatches,
2060  regionI,
2061  (overwrite ? oldInstance : runTime.timeName())
2062  );
2063  }
2064  else
2065  {
2066  // Split all
2067  for (label regionI = 0; regionI < nCellRegions; regionI++)
2068  {
2069  Info<< nl
2070  << "Region " << regionI << nl
2071  << "-------- " << endl;
2072 
2073  createAndWriteRegion
2074  (
2075  mesh,
2076  cellRegion,
2077  regionNames,
2078  prefixRegion,
2079  faceToInterface,
2080  interfacePatches,
2081  regionI,
2082  (overwrite ? oldInstance : runTime.timeName())
2083  );
2084  }
2085  }
2086  }
2087 
2088  Info<< "End\n" << endl;
2089 
2090  return 0;
2091 }
2092 
2093 
2094 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order. Return.
Definition: ZoneMesh.C:428
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:119
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:230
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches. If validBoundary call is parallel.
Definition: fvMeshTools.C:324
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Read all fields of the specified type.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:387
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
static int masterNo()
Process index of the master.
Definition: UPstream.H:406
const word & name() const
Return name.
Definition: IOobject.H:291
A list of face labels.
Definition: faceSet.H:48
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:466
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:801
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:110
const double e
Elementary charge.
Definition: doubleFloat.H:78
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
interfaceProperties interface(alpha1, U, mixture())
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1055
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:841
static int firstSlave()
Process index of first slave.
Definition: UPstream.H:435
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:877
PtrList< labelIOList > & faceProcAddressing
label nFaces() const
Foam::word regionName
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:418
Given list of cells to remove insert all the topology changes.
Definition: removeCells.H:59
label nCells() const
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:412
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Field reading functions for post-processing utilities.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurences of given element. Linear search.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:814
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:104
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
Input inter-processor communications stream.
Definition: IPstream.H:50
bool insert(const edge &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const edge &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:795
dynamicFvMesh & mesh
T optionRead(const word &opt) const
Read a value from the named option.
Definition: argListI.H:187
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:341
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
A class for handling words, derived from string.
Definition: word.H:59
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:472
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:32
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:95
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1049
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:91
const word & name() const
Return name.
Definition: zone.H:147
An STL-conforming hash table.
Definition: HashTable.H:62
const vectorField & cellCentres() const
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
const word & system() const
Return system name.
Definition: TimePaths.H:114
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Foam::polyBoundaryMesh.
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
static const char nl
Definition: Ostream.H:262
Map from edge (expressed as its endpoints) to value.
Definition: EdgeMap.H:47
void updateMesh()
Correct polyBoundaryMesh after topology update.
Type gMax(const FieldField< Field, Type > &f)
Output inter-processor communications stream.
Definition: OPstream.H:50
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
A subset of mesh cells.
Definition: cellZone.H:61
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:221
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:430
writeOption writeOpt() const
Definition: IOobject.H:363
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:394
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:400
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
A collection of cell labels.
Definition: cellSet.H:48
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Direct mesh changes based on v1.3 polyTopoChange syntax.
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1394
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:256
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual Ostream & write(const token &)=0
Write next token to stream.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:85
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:126
A class for managing temporary objects.
Definition: PtrList.H:53
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
const string & option(const word &opt) const
Return the argument string associated with the named option.
Definition: argListI.H:102
Namespace for OpenFOAM.
static int lastSlave(const label communicator=0)
Process index of last slave.
Definition: UPstream.H:441