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-2016 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::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  {
443  toMaster << regionsToSize;
444  }
445  }
446  }
447 
448  // Rework
449 
450  Pstream::scatter(regionsToSize);
451 
452 
453 
454  // Now we have the global sizes of all inter-regions.
455  // Invert this on master and distribute.
456  label nInterfaces = 0;
457  forAllConstIter(EdgeMap<Map<label>>, regionsToSize, iter)
458  {
459  const Map<label>& info = iter();
460  nInterfaces += info.size();
461  }
462 
463  interfaces.setSize(nInterfaces);
464  interfaceNames.setSize(nInterfaces);
465  interfaceSizes.setSize(nInterfaces);
466  EdgeMap<Map<label>> regionsToInterface(nInterfaces);
467 
468  nInterfaces = 0;
469  forAllConstIter(EdgeMap<Map<label>>, regionsToSize, iter)
470  {
471  const edge& e = iter.key();
472  const word& name0 = regionNames[e[0]];
473  const word& name1 = regionNames[e[1]];
474 
475  const Map<label>& info = iter();
476  forAllConstIter(Map<label>, info, infoIter)
477  {
478  interfaces[nInterfaces] = iter.key();
479  label zoneID = infoIter.key();
480  if (zoneID == -1)
481  {
482  interfaceNames[nInterfaces] = Pair<word>
483  (
484  name0 + "_to_" + name1,
485  name1 + "_to_" + name0
486  );
487  }
488  else
489  {
490  const word& zoneName = mesh.faceZones()[zoneID].name();
491  interfaceNames[nInterfaces] = Pair<word>
492  (
493  zoneName + "_" + name0 + "_to_" + name1,
494  zoneName + "_" + name1 + "_to_" + name0
495  );
496  }
497  interfaceSizes[nInterfaces] = infoIter();
498 
499  if (regionsToInterface.found(e))
500  {
501  regionsToInterface[e].insert(zoneID, nInterfaces);
502  }
503  else
504  {
505  Map<label> zoneAndInterface;
506  zoneAndInterface.insert(zoneID, nInterfaces);
507  regionsToInterface.insert(e, zoneAndInterface);
508  }
509  nInterfaces++;
510  }
511  }
512 
513 
514  // Now all processor have consistent interface information
515 
516  Pstream::scatter(interfaces);
517  Pstream::scatter(interfaceNames);
518  Pstream::scatter(interfaceSizes);
519  Pstream::scatter(regionsToInterface);
520 
521  // Mark all inter-region faces.
522  faceToInterface.setSize(mesh.nFaces(), -1);
523 
524  forAll(mesh.faceNeighbour(), facei)
525  {
526  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
527  label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
528 
529  if (ownRegion != neiRegion)
530  {
531  label zoneID = -1;
532  if (useFaceZones)
533  {
534  zoneID = mesh.faceZones().whichZone(facei);
535  }
536 
538  (
539  min(ownRegion, neiRegion),
540  max(ownRegion, neiRegion)
541  );
542 
543  faceToInterface[facei] = regionsToInterface[interface][zoneID];
544  }
545  }
546  forAll(coupledRegion, i)
547  {
548  label facei = i+mesh.nInternalFaces();
549  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
550  label neiRegion = coupledRegion[i];
551 
552  if (ownRegion != neiRegion)
553  {
554  label zoneID = -1;
555  if (useFaceZones)
556  {
557  zoneID = mesh.faceZones().whichZone(facei);
558  }
559 
561  (
562  min(ownRegion, neiRegion),
563  max(ownRegion, neiRegion)
564  );
565 
566  faceToInterface[facei] = regionsToInterface[interface][zoneID];
567  }
568  }
569 }
570 
571 
572 // Create mesh for region.
573 autoPtr<mapPolyMesh> createRegionMesh
574 (
575  const fvMesh& mesh,
576  // Region info
577  const labelList& cellRegion,
578  const label regionI,
579  const word& regionName,
580  // Interface info
581  const labelList& interfacePatches,
582  const labelList& faceToInterface,
583 
584  autoPtr<fvMesh>& newMesh
585 )
586 {
587  // Create dummy system/fv*
588  {
589  IOobject io
590  (
591  "fvSchemes",
592  mesh.time().system(),
593  regionName,
594  mesh,
597  false
598  );
599 
600  Info<< "Testing:" << io.objectPath() << endl;
601 
602  if (!io.headerOk())
603  // if (!exists(io.objectPath()))
604  {
605  Info<< "Writing dummy " << regionName/io.name() << endl;
606  dictionary dummyDict;
607  dictionary divDict;
608  dummyDict.add("divSchemes", divDict);
609  dictionary gradDict;
610  dummyDict.add("gradSchemes", gradDict);
611  dictionary laplDict;
612  dummyDict.add("laplacianSchemes", laplDict);
613 
614  IOdictionary(io, dummyDict).regIOobject::write();
615  }
616  }
617  {
618  IOobject io
619  (
620  "fvSolution",
621  mesh.time().system(),
622  regionName,
623  mesh,
626  false
627  );
628 
629  if (!io.headerOk())
630  //if (!exists(io.objectPath()))
631  {
632  Info<< "Writing dummy " << regionName/io.name() << endl;
633  dictionary dummyDict;
634  IOdictionary(io, dummyDict).regIOobject::write();
635  }
636  }
637 
638 
639  // Neighbour cellRegion.
640  labelList coupledRegion(mesh.nFaces()-mesh.nInternalFaces());
641 
642  forAll(coupledRegion, i)
643  {
644  label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
645  coupledRegion[i] = cellRegion[celli];
646  }
647  syncTools::swapBoundaryFaceList(mesh, coupledRegion);
648 
649 
650  // Topology change container. Start off from existing mesh.
651  polyTopoChange meshMod(mesh);
652 
653  // Cell remover engine
654  removeCells cellRemover(mesh);
655 
656  // Select all but region cells
657  labelList cellsToRemove(getNonRegionCells(cellRegion, regionI));
658 
659  // Find out which faces will get exposed. Note that this
660  // gets faces in mesh face order. So both regions will get same
661  // face in same order (important!)
662  labelList exposedFaces = cellRemover.getExposedFaces(cellsToRemove);
663 
664  labelList exposedPatchIDs(exposedFaces.size());
665  forAll(exposedFaces, i)
666  {
667  label facei = exposedFaces[i];
668  label interfacei = faceToInterface[facei];
669 
670  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
671  label neiRegion = -1;
672 
673  if (mesh.isInternalFace(facei))
674  {
675  neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
676  }
677  else
678  {
679  neiRegion = coupledRegion[facei-mesh.nInternalFaces()];
680  }
681 
682 
683  // Check which side is being kept - determines which of the two
684  // patches will be used.
685 
686  label otherRegion = -1;
687 
688  if (ownRegion == regionI && neiRegion != regionI)
689  {
690  otherRegion = neiRegion;
691  }
692  else if (ownRegion != regionI && neiRegion == regionI)
693  {
694  otherRegion = ownRegion;
695  }
696  else
697  {
699  << "Exposed face:" << facei
700  << " fc:" << mesh.faceCentres()[facei]
701  << " has owner region " << ownRegion
702  << " and neighbour region " << neiRegion
703  << " when handling region:" << regionI
704  << exit(FatalError);
705  }
706 
707  // Find the patch.
708  if (regionI < otherRegion)
709  {
710  exposedPatchIDs[i] = interfacePatches[interfacei];
711  }
712  else
713  {
714  exposedPatchIDs[i] = interfacePatches[interfacei]+1;
715  }
716  }
717 
718  // Remove faces
719  cellRemover.setRefinement
720  (
721  cellsToRemove,
722  exposedFaces,
723  exposedPatchIDs,
724  meshMod
725  );
726 
727  autoPtr<mapPolyMesh> map = meshMod.makeMesh
728  (
729  newMesh,
730  IOobject
731  (
732  regionName,
733  mesh.time().timeName(),
734  mesh.time(),
737  ),
738  mesh
739  );
740 
741  return map;
742 }
743 
744 
745 void createAndWriteRegion
746 (
747  const fvMesh& mesh,
748  const labelList& cellRegion,
749  const wordList& regionNames,
750  const bool prefixRegion,
751  const labelList& faceToInterface,
752  const labelList& interfacePatches,
753  const label regionI,
754  const word& newMeshInstance
755 )
756 {
757  Info<< "Creating mesh for region " << regionI
758  << ' ' << regionNames[regionI] << endl;
759 
760  autoPtr<fvMesh> newMesh;
761  autoPtr<mapPolyMesh> map = createRegionMesh
762  (
763  mesh,
764  cellRegion,
765  regionI,
766  regionNames[regionI],
767  interfacePatches,
768  faceToInterface,
769  newMesh
770  );
771 
772 
773  // Make map of all added patches
774  labelHashSet addedPatches(2*interfacePatches.size());
775  forAll(interfacePatches, interfacei)
776  {
777  addedPatches.insert(interfacePatches[interfacei]);
778  addedPatches.insert(interfacePatches[interfacei]+1);
779  }
780 
781 
782  Info<< "Mapping fields" << endl;
783 
784  // Map existing fields
785  newMesh().updateMesh(map());
786 
787  // Add subsetted fields
788  subsetVolFields<volScalarField>
789  (
790  mesh,
791  newMesh(),
792  map().cellMap(),
793  map().faceMap(),
794  addedPatches
795  );
796  subsetVolFields<volVectorField>
797  (
798  mesh,
799  newMesh(),
800  map().cellMap(),
801  map().faceMap(),
802  addedPatches
803  );
804  subsetVolFields<volSphericalTensorField>
805  (
806  mesh,
807  newMesh(),
808  map().cellMap(),
809  map().faceMap(),
810  addedPatches
811  );
812  subsetVolFields<volSymmTensorField>
813  (
814  mesh,
815  newMesh(),
816  map().cellMap(),
817  map().faceMap(),
818  addedPatches
819  );
820  subsetVolFields<volTensorField>
821  (
822  mesh,
823  newMesh(),
824  map().cellMap(),
825  map().faceMap(),
826  addedPatches
827  );
828 
829  subsetSurfaceFields<surfaceScalarField>
830  (
831  mesh,
832  newMesh(),
833  map().cellMap(),
834  map().faceMap(),
835  addedPatches
836  );
837  subsetSurfaceFields<surfaceVectorField>
838  (
839  mesh,
840  newMesh(),
841  map().cellMap(),
842  map().faceMap(),
843  addedPatches
844  );
845  subsetSurfaceFields<surfaceSphericalTensorField>
846  (
847  mesh,
848  newMesh(),
849  map().cellMap(),
850  map().faceMap(),
851  addedPatches
852  );
853  subsetSurfaceFields<surfaceSymmTensorField>
854  (
855  mesh,
856  newMesh(),
857  map().cellMap(),
858  map().faceMap(),
859  addedPatches
860  );
861  subsetSurfaceFields<surfaceTensorField>
862  (
863  mesh,
864  newMesh(),
865  map().cellMap(),
866  map().faceMap(),
867  addedPatches
868  );
869 
870 
871  const polyBoundaryMesh& newPatches = newMesh().boundaryMesh();
872  newPatches.checkParallelSync(true);
873 
874  // Delete empty patches
875  // ~~~~~~~~~~~~~~~~~~~~
876 
877  // Create reordering list to move patches-to-be-deleted to end
878  labelList oldToNew(newPatches.size(), -1);
879  DynamicList<label> sharedPatches(newPatches.size());
880  label newI = 0;
881 
882  Info<< "Deleting empty patches" << endl;
883 
884  // Assumes all non-proc boundaries are on all processors!
885  forAll(newPatches, patchi)
886  {
887  const polyPatch& pp = newPatches[patchi];
888 
889  if (!isA<processorPolyPatch>(pp))
890  {
891  if (returnReduce(pp.size(), sumOp<label>()) > 0)
892  {
893  oldToNew[patchi] = newI;
894  if (!addedPatches.found(patchi))
895  {
896  sharedPatches.append(newI);
897  }
898  newI++;
899  }
900  }
901  }
902 
903  // Same for processor patches (but need no reduction)
904  forAll(newPatches, patchi)
905  {
906  const polyPatch& pp = newPatches[patchi];
907 
908  if (isA<processorPolyPatch>(pp) && pp.size())
909  {
910  oldToNew[patchi] = newI++;
911  }
912  }
913 
914  const label nNewPatches = newI;
915 
916  // Move all deleteable patches to the end
917  forAll(oldToNew, patchi)
918  {
919  if (oldToNew[patchi] == -1)
920  {
921  oldToNew[patchi] = newI++;
922  }
923  }
924 
925  //reorderPatches(newMesh(), oldToNew, nNewPatches);
926  fvMeshTools::reorderPatches(newMesh(), oldToNew, nNewPatches, true);
927 
928  // Rename shared patches with region name
929  if (prefixRegion)
930  {
931  Info<< "Prefixing patches with region name" << endl;
932 
933  renamePatches(newMesh(), regionNames[regionI], sharedPatches);
934  }
935 
936 
937  Info<< "Writing new mesh" << endl;
938 
939  newMesh().setInstance(newMeshInstance);
940  newMesh().write();
941 
942  // Write addressing files like decomposePar
943  Info<< "Writing addressing to base mesh" << endl;
944 
945  labelIOList pointProcAddressing
946  (
947  IOobject
948  (
949  "pointRegionAddressing",
950  newMesh().facesInstance(),
951  newMesh().meshSubDir,
952  newMesh(),
953  IOobject::NO_READ,
954  IOobject::NO_WRITE,
955  false
956  ),
957  map().pointMap()
958  );
959  Info<< "Writing map " << pointProcAddressing.name()
960  << " from region" << regionI
961  << " points back to base mesh." << endl;
962  pointProcAddressing.write();
963 
965  (
966  IOobject
967  (
968  "faceRegionAddressing",
969  newMesh().facesInstance(),
970  newMesh().meshSubDir,
971  newMesh(),
972  IOobject::NO_READ,
973  IOobject::NO_WRITE,
974  false
975  ),
976  newMesh().nFaces()
977  );
978  forAll(faceProcAddressing, facei)
979  {
980  // face + turning index. (see decomposePar)
981  // Is the face pointing in the same direction?
982  label oldFacei = map().faceMap()[facei];
983 
984  if
985  (
986  map().cellMap()[newMesh().faceOwner()[facei]]
987  == mesh.faceOwner()[oldFacei]
988  )
989  {
990  faceProcAddressing[facei] = oldFacei+1;
991  }
992  else
993  {
994  faceProcAddressing[facei] = -(oldFacei+1);
995  }
996  }
997  Info<< "Writing map " << faceProcAddressing.name()
998  << " from region" << regionI
999  << " faces back to base mesh." << endl;
1000  faceProcAddressing.write();
1001 
1002  labelIOList cellProcAddressing
1003  (
1004  IOobject
1005  (
1006  "cellRegionAddressing",
1007  newMesh().facesInstance(),
1008  newMesh().meshSubDir,
1009  newMesh(),
1010  IOobject::NO_READ,
1011  IOobject::NO_WRITE,
1012  false
1013  ),
1014  map().cellMap()
1015  );
1016  Info<< "Writing map " <<cellProcAddressing.name()
1017  << " from region" << regionI
1018  << " cells back to base mesh." << endl;
1019  cellProcAddressing.write();
1020 
1021  labelIOList boundaryProcAddressing
1022  (
1023  IOobject
1024  (
1025  "boundaryRegionAddressing",
1026  newMesh().facesInstance(),
1027  newMesh().meshSubDir,
1028  newMesh(),
1029  IOobject::NO_READ,
1030  IOobject::NO_WRITE,
1031  false
1032  ),
1033  labelList(nNewPatches, -1)
1034  );
1035  forAll(oldToNew, i)
1036  {
1037  if (!addedPatches.found(i))
1038  {
1039  label newI = oldToNew[i];
1040  if (newI >= 0 && newI < nNewPatches)
1041  {
1042  boundaryProcAddressing[oldToNew[i]] = i;
1043  }
1044  }
1045  }
1046  Info<< "Writing map " << boundaryProcAddressing.name()
1047  << " from region" << regionI
1048  << " boundary back to base mesh." << endl;
1049  boundaryProcAddressing.write();
1050 }
1051 
1052 
1053 // Create for every region-region interface with non-zero size two patches.
1054 // First one is for minimumregion to maximumregion.
1055 // Note that patches get created in same order on all processors (if parallel)
1056 // since looping over synchronised 'interfaces'.
1057 labelList addRegionPatches
1058 (
1059  fvMesh& mesh,
1060  const wordList& regionNames,
1061  const edgeList& interfaces,
1062  const List<Pair<word>>& interfaceNames
1063 )
1064 {
1065  Info<< nl << "Adding patches" << nl << endl;
1066 
1067  labelList interfacePatches(interfaces.size());
1068 
1069  forAll(interfaces, interI)
1070  {
1071  const edge& e = interfaces[interI];
1072  const Pair<word>& names = interfaceNames[interI];
1073 
1074  //Info<< "For interface " << interI
1075  // << " between regions " << e
1076  // << " trying to add patches " << names << endl;
1077 
1078 
1079  mappedWallPolyPatch patch1
1080  (
1081  names[0],
1082  0, // overridden
1083  0, // overridden
1084  0, // overridden
1085  regionNames[e[1]], // sampleRegion
1087  names[1], // samplePatch
1088  point::zero, // offset
1089  mesh.boundaryMesh()
1090  );
1091 
1092  interfacePatches[interI] = fvMeshTools::addPatch
1093  (
1094  mesh,
1095  patch1,
1096  dictionary(), //optional per field value
1098  true //validBoundary
1099  );
1100 
1101  mappedWallPolyPatch patch2
1102  (
1103  names[1],
1104  0,
1105  0,
1106  0,
1107  regionNames[e[0]], // sampleRegion
1109  names[0],
1110  point::zero, // offset
1111  mesh.boundaryMesh()
1112  );
1114  (
1115  mesh,
1116  patch2,
1117  dictionary(), //optional per field value
1119  true //validBoundary
1120  );
1121 
1122  Info<< "For interface between region " << regionNames[e[0]]
1123  << " and " << regionNames[e[1]] << " added patches" << endl
1124  << " " << interfacePatches[interI]
1125  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]].name()
1126  << endl
1127  << " " << interfacePatches[interI]+1
1128  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]+1].name()
1129  << endl;
1130  }
1131  return interfacePatches;
1132 }
1133 
1134 
1135 // Find region that covers most of cell zone
1136 label findCorrespondingRegion
1137 (
1138  const labelList& existingZoneID, // per cell the (unique) zoneID
1139  const labelList& cellRegion,
1140  const label nCellRegions,
1141  const label zoneI,
1142  const label minOverlapSize
1143 )
1144 {
1145  // Per region the number of cells in zoneI
1146  labelList cellsInZone(nCellRegions, 0);
1147 
1148  forAll(cellRegion, celli)
1149  {
1150  if (existingZoneID[celli] == zoneI)
1151  {
1152  cellsInZone[cellRegion[celli]]++;
1153  }
1154  }
1155 
1157  Pstream::listCombineScatter(cellsInZone);
1158 
1159  // Pick region with largest overlap of zoneI
1160  label regionI = findMax(cellsInZone);
1161 
1162 
1163  if (cellsInZone[regionI] < minOverlapSize)
1164  {
1165  // Region covers too little of zone. Not good enough.
1166  regionI = -1;
1167  }
1168  else
1169  {
1170  // Check that region contains no cells that aren't in cellZone.
1171  forAll(cellRegion, celli)
1172  {
1173  if (cellRegion[celli] == regionI && existingZoneID[celli] != zoneI)
1174  {
1175  // celli in regionI but not in zoneI
1176  regionI = -1;
1177  break;
1178  }
1179  }
1180  // If one in error, all should be in error. Note that branch gets taken
1181  // on all procs.
1182  reduce(regionI, minOp<label>());
1183  }
1184 
1185  return regionI;
1186 }
1187 
1188 
1189 // Get zone per cell
1190 // - non-unique zoning
1191 // - coupled zones
1192 void getZoneID
1193 (
1194  const polyMesh& mesh,
1195  const cellZoneMesh& cellZones,
1196  labelList& zoneID,
1197  labelList& neiZoneID
1198 )
1199 {
1200  // Existing zoneID
1201  zoneID.setSize(mesh.nCells());
1202  zoneID = -1;
1203 
1204  forAll(cellZones, zoneI)
1205  {
1206  const cellZone& cz = cellZones[zoneI];
1207 
1208  forAll(cz, i)
1209  {
1210  label celli = cz[i];
1211  if (zoneID[celli] == -1)
1212  {
1213  zoneID[celli] = zoneI;
1214  }
1215  else
1216  {
1218  << "Cell " << celli << " with cell centre "
1219  << mesh.cellCentres()[celli]
1220  << " is multiple zones. This is not allowed." << endl
1221  << "It is in zone " << cellZones[zoneID[celli]].name()
1222  << " and in zone " << cellZones[zoneI].name()
1223  << exit(FatalError);
1224  }
1225  }
1226  }
1227 
1228  // Neighbour zoneID.
1229  neiZoneID.setSize(mesh.nFaces()-mesh.nInternalFaces());
1230 
1231  forAll(neiZoneID, i)
1232  {
1233  neiZoneID[i] = zoneID[mesh.faceOwner()[i+mesh.nInternalFaces()]];
1234  }
1235  syncTools::swapBoundaryFaceList(mesh, neiZoneID);
1236 }
1237 
1238 
1239 void matchRegions
1240 (
1241  const bool sloppyCellZones,
1242  const polyMesh& mesh,
1243 
1244  const label nCellRegions,
1245  const labelList& cellRegion,
1246 
1247  labelList& regionToZone,
1248  wordList& regionNames,
1249  labelList& zoneToRegion
1250 )
1251 {
1252  const cellZoneMesh& cellZones = mesh.cellZones();
1253 
1254  regionToZone.setSize(nCellRegions, -1);
1255  regionNames.setSize(nCellRegions);
1256  zoneToRegion.setSize(cellZones.size(), -1);
1257 
1258  // Get current per cell zoneID
1259  labelList zoneID(mesh.nCells(), -1);
1260  labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1261  getZoneID(mesh, cellZones, zoneID, neiZoneID);
1262 
1263  // Sizes per cellzone
1264  labelList zoneSizes(cellZones.size(), 0);
1265  {
1266  List<wordList> zoneNames(Pstream::nProcs());
1267  zoneNames[Pstream::myProcNo()] = cellZones.names();
1268  Pstream::gatherList(zoneNames);
1269  Pstream::scatterList(zoneNames);
1270 
1271  forAll(zoneNames, proci)
1272  {
1273  if (zoneNames[proci] != zoneNames[0])
1274  {
1276  << "cellZones not synchronised across processors." << endl
1277  << "Master has cellZones " << zoneNames[0] << endl
1278  << "Processor " << proci
1279  << " has cellZones " << zoneNames[proci]
1280  << exit(FatalError);
1281  }
1282  }
1283 
1284  forAll(cellZones, zoneI)
1285  {
1286  zoneSizes[zoneI] = returnReduce
1287  (
1288  cellZones[zoneI].size(),
1289  sumOp<label>()
1290  );
1291  }
1292  }
1293 
1294 
1295  if (sloppyCellZones)
1296  {
1297  Info<< "Trying to match regions to existing cell zones;"
1298  << " region can be subset of cell zone." << nl << endl;
1299 
1300  forAll(cellZones, zoneI)
1301  {
1302  label regionI = findCorrespondingRegion
1303  (
1304  zoneID,
1305  cellRegion,
1306  nCellRegions,
1307  zoneI,
1308  label(0.5*zoneSizes[zoneI]) // minimum overlap
1309  );
1310 
1311  if (regionI != -1)
1312  {
1313  Info<< "Sloppily matched region " << regionI
1314  //<< " size " << regionSizes[regionI]
1315  << " to zone " << zoneI << " size " << zoneSizes[zoneI]
1316  << endl;
1317  zoneToRegion[zoneI] = regionI;
1318  regionToZone[regionI] = zoneI;
1319  regionNames[regionI] = cellZones[zoneI].name();
1320  }
1321  }
1322  }
1323  else
1324  {
1325  Info<< "Trying to match regions to existing cell zones." << nl << endl;
1326 
1327  forAll(cellZones, zoneI)
1328  {
1329  label regionI = findCorrespondingRegion
1330  (
1331  zoneID,
1332  cellRegion,
1333  nCellRegions,
1334  zoneI,
1335  1 // minimum overlap
1336  );
1337 
1338  if (regionI != -1)
1339  {
1340  zoneToRegion[zoneI] = regionI;
1341  regionToZone[regionI] = zoneI;
1342  regionNames[regionI] = cellZones[zoneI].name();
1343  }
1344  }
1345  }
1346  // Allocate region names for unmatched regions.
1347  forAll(regionToZone, regionI)
1348  {
1349  if (regionToZone[regionI] == -1)
1350  {
1351  regionNames[regionI] = "domain" + Foam::name(regionI);
1352  }
1353  }
1354 }
1355 
1356 
1357 void writeCellToRegion(const fvMesh& mesh, const labelList& cellRegion)
1358 {
1359  // Write to manual decomposition option
1360  {
1361  labelIOList cellToRegion
1362  (
1363  IOobject
1364  (
1365  "cellToRegion",
1366  mesh.facesInstance(),
1367  mesh,
1370  false
1371  ),
1372  cellRegion
1373  );
1374  cellToRegion.write();
1375 
1376  Info<< "Writing region per cell file (for manual decomposition) to "
1377  << cellToRegion.objectPath() << nl << endl;
1378  }
1379  // Write for postprocessing
1380  {
1381  volScalarField cellToRegion
1382  (
1383  IOobject
1384  (
1385  "cellToRegion",
1386  mesh.time().timeName(),
1387  mesh,
1390  false
1391  ),
1392  mesh,
1393  dimensionedScalar("zero", dimless, 0),
1394  zeroGradientFvPatchScalarField::typeName
1395  );
1396  forAll(cellRegion, celli)
1397  {
1398  cellToRegion[celli] = cellRegion[celli];
1399  }
1400  cellToRegion.write();
1401 
1402  Info<< "Writing region per cell as volScalarField to "
1403  << cellToRegion.objectPath() << nl << endl;
1404  }
1405 }
1406 
1407 
1408 
1409 
1410 int main(int argc, char *argv[])
1411 {
1413  (
1414  "splits mesh into multiple regions (detected by walking across faces)"
1415  );
1416  #include "addRegionOption.H"
1417  #include "addOverwriteOption.H"
1419  (
1420  "cellZones",
1421  "additionally split cellZones off into separate regions"
1422  );
1424  (
1425  "cellZonesOnly",
1426  "use cellZones only to split mesh into regions; do not use walking"
1427  );
1429  (
1430  "cellZonesFileOnly",
1431  "file",
1432  "like -cellZonesOnly, but use specified file"
1433  );
1435  (
1436  "blockedFaces",
1437  "faceSet",
1438  "specify additional region boundaries that walking does not cross"
1439  );
1441  (
1442  "makeCellZones",
1443  "place cells into cellZones instead of splitting mesh"
1444  );
1446  (
1447  "largestOnly",
1448  "only write largest region"
1449  );
1451  (
1452  "insidePoint",
1453  "point",
1454  "only write region containing point"
1455  );
1457  (
1458  "detectOnly",
1459  "do not write mesh"
1460  );
1462  (
1463  "sloppyCellZones",
1464  "try to match heuristically regions to existing cell zones"
1465  );
1467  (
1468  "useFaceZones",
1469  "use faceZones to patch inter-region faces instead of single patch"
1470  );
1472  (
1473  "prefixRegion",
1474  "prefix region name to all patches, not just coupling patches"
1475  );
1476 
1477  #include "setRootCase.H"
1478  #include "createTime.H"
1479  runTime.functionObjects().off();
1480  #include "createNamedMesh.H"
1481  const word oldInstance = mesh.pointsInstance();
1482 
1483  word blockedFacesName;
1484  if (args.optionReadIfPresent("blockedFaces", blockedFacesName))
1485  {
1486  Info<< "Reading blocked internal faces from faceSet "
1487  << blockedFacesName << nl << endl;
1488  }
1489 
1490  const bool makeCellZones = args.optionFound("makeCellZones");
1491  const bool largestOnly = args.optionFound("largestOnly");
1492  const bool insidePoint = args.optionFound("insidePoint");
1493  const bool useCellZones = args.optionFound("cellZones");
1494  const bool useCellZonesOnly = args.optionFound("cellZonesOnly");
1495  const bool useCellZonesFile = args.optionFound("cellZonesFileOnly");
1496  const bool overwrite = args.optionFound("overwrite");
1497  const bool detectOnly = args.optionFound("detectOnly");
1498  const bool sloppyCellZones = args.optionFound("sloppyCellZones");
1499  const bool useFaceZones = args.optionFound("useFaceZones");
1500  const bool prefixRegion = args.optionFound("prefixRegion");
1501 
1502 
1503  if
1504  (
1505  (useCellZonesOnly || useCellZonesFile)
1506  && (useCellZones || blockedFacesName.size())
1507  )
1508  {
1510  << "You cannot specify both -cellZonesOnly or -cellZonesFileOnly"
1511  << " (which specify complete split)"
1512  << " in combination with -blockedFaces or -cellZones"
1513  << " (which imply a split based on topology)"
1514  << exit(FatalError);
1515  }
1516 
1517 
1518  if (useFaceZones)
1519  {
1520  Info<< "Using current faceZones to divide inter-region interfaces"
1521  << " into multiple patches."
1522  << nl << endl;
1523  }
1524  else
1525  {
1526  Info<< "Creating single patch per inter-region interface."
1527  << nl << endl;
1528  }
1529 
1530 
1531 
1532  if (insidePoint && largestOnly)
1533  {
1535  << "You cannot specify both -largestOnly"
1536  << " (keep region with most cells)"
1537  << " and -insidePoint (keep region containing point)"
1538  << exit(FatalError);
1539  }
1540 
1541 
1542  const cellZoneMesh& cellZones = mesh.cellZones();
1543 
1544  // Existing zoneID
1545  labelList zoneID(mesh.nCells(), -1);
1546  // Neighbour zoneID.
1547  labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1548  getZoneID(mesh, cellZones, zoneID, neiZoneID);
1549 
1550 
1551 
1552  // Determine per cell the region it belongs to
1553  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1554 
1555  // cellRegion is the labelList with the region per cell.
1556  labelList cellRegion;
1557  // Region per zone
1558  labelList regionToZone;
1559  // Name of region
1560  wordList regionNames;
1561  // Zone to region
1562  labelList zoneToRegion;
1563 
1564  label nCellRegions = 0;
1565  if (useCellZonesOnly)
1566  {
1567  Info<< "Using current cellZones to split mesh into regions."
1568  << " This requires all"
1569  << " cells to be in one and only one cellZone." << nl << endl;
1570 
1571  label unzonedCelli = findIndex(zoneID, -1);
1572  if (unzonedCelli != -1)
1573  {
1575  << "For the cellZonesOnly option all cells "
1576  << "have to be in a cellZone." << endl
1577  << "Cell " << unzonedCelli
1578  << " at" << mesh.cellCentres()[unzonedCelli]
1579  << " is not in a cellZone. There might be more unzoned cells."
1580  << exit(FatalError);
1581  }
1582  cellRegion = zoneID;
1583  nCellRegions = gMax(cellRegion)+1;
1584  regionToZone.setSize(nCellRegions);
1585  regionNames.setSize(nCellRegions);
1586  zoneToRegion.setSize(cellZones.size(), -1);
1587  for (label regionI = 0; regionI < nCellRegions; regionI++)
1588  {
1589  regionToZone[regionI] = regionI;
1590  zoneToRegion[regionI] = regionI;
1591  regionNames[regionI] = cellZones[regionI].name();
1592  }
1593  }
1594  else if (useCellZonesFile)
1595  {
1596  const word zoneFile = args.option("cellZonesFileOnly");
1597  Info<< "Reading split from cellZones file " << zoneFile << endl
1598  << "This requires all"
1599  << " cells to be in one and only one cellZone." << nl << endl;
1600 
1601  cellZoneMesh newCellZones
1602  (
1603  IOobject
1604  (
1605  zoneFile,
1606  mesh.facesInstance(),
1608  mesh,
1611  false
1612  ),
1613  mesh
1614  );
1615 
1616  labelList newZoneID(mesh.nCells(), -1);
1617  labelList newNeiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1618  getZoneID(mesh, newCellZones, newZoneID, newNeiZoneID);
1619 
1620  label unzonedCelli = findIndex(newZoneID, -1);
1621  if (unzonedCelli != -1)
1622  {
1624  << "For the cellZonesFileOnly option all cells "
1625  << "have to be in a cellZone." << endl
1626  << "Cell " << unzonedCelli
1627  << " at" << mesh.cellCentres()[unzonedCelli]
1628  << " is not in a cellZone. There might be more unzoned cells."
1629  << exit(FatalError);
1630  }
1631  cellRegion = newZoneID;
1632  nCellRegions = gMax(cellRegion)+1;
1633  zoneToRegion.setSize(newCellZones.size(), -1);
1634  regionToZone.setSize(nCellRegions);
1635  regionNames.setSize(nCellRegions);
1636  for (label regionI = 0; regionI < nCellRegions; regionI++)
1637  {
1638  regionToZone[regionI] = regionI;
1639  zoneToRegion[regionI] = regionI;
1640  regionNames[regionI] = newCellZones[regionI].name();
1641  }
1642  }
1643  else
1644  {
1645  // Determine connected regions
1646  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1647 
1648  // Mark additional faces that are blocked
1649  boolList blockedFace;
1650 
1651  // Read from faceSet
1652  if (blockedFacesName.size())
1653  {
1654  faceSet blockedFaceSet(mesh, blockedFacesName);
1655  Info<< "Read "
1656  << returnReduce(blockedFaceSet.size(), sumOp<label>())
1657  << " blocked faces from set " << blockedFacesName << nl << endl;
1658 
1659  blockedFace.setSize(mesh.nFaces(), false);
1660 
1661  forAllConstIter(faceSet, blockedFaceSet, iter)
1662  {
1663  blockedFace[iter.key()] = true;
1664  }
1665  }
1666 
1667  // Imply from differing cellZones
1668  if (useCellZones)
1669  {
1670  blockedFace.setSize(mesh.nFaces(), false);
1671 
1672  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1673  {
1674  label own = mesh.faceOwner()[facei];
1675  label nei = mesh.faceNeighbour()[facei];
1676 
1677  if (zoneID[own] != zoneID[nei])
1678  {
1679  blockedFace[facei] = true;
1680  }
1681  }
1682 
1683  // Different cellZones on either side of processor patch.
1684  forAll(neiZoneID, i)
1685  {
1686  label facei = i+mesh.nInternalFaces();
1687 
1688  if (zoneID[mesh.faceOwner()[facei]] != neiZoneID[i])
1689  {
1690  blockedFace[facei] = true;
1691  }
1692  }
1693  }
1694 
1695  // Do a topological walk to determine regions
1696  regionSplit regions(mesh, blockedFace);
1697  nCellRegions = regions.nRegions();
1698  cellRegion.transfer(regions);
1699 
1700  // Make up region names. If possible match them to existing zones.
1701  matchRegions
1702  (
1703  sloppyCellZones,
1704  mesh,
1705  nCellRegions,
1706  cellRegion,
1707 
1708  regionToZone,
1709  regionNames,
1710  zoneToRegion
1711  );
1712 
1713  // Override any default region names if single region selected
1714  if (largestOnly || insidePoint)
1715  {
1716  forAll(regionToZone, regionI)
1717  {
1718  if (regionToZone[regionI] == -1)
1719  {
1720  if (overwrite)
1721  {
1722  regionNames[regionI] = polyMesh::defaultRegion;
1723  }
1724  else if (insidePoint)
1725  {
1726  regionNames[regionI] = "insidePoint";
1727  }
1728  else if (largestOnly)
1729  {
1730  regionNames[regionI] = "largestOnly";
1731  }
1732  }
1733  }
1734  }
1735  }
1736 
1737  Info<< endl << "Number of regions:" << nCellRegions << nl << endl;
1738 
1739 
1740  // Write decomposition to file
1741  writeCellToRegion(mesh, cellRegion);
1742 
1743 
1744 
1745  // Sizes per region
1746  // ~~~~~~~~~~~~~~~~
1747 
1748  labelList regionSizes(nCellRegions, 0);
1749 
1750  forAll(cellRegion, celli)
1751  {
1752  regionSizes[cellRegion[celli]]++;
1753  }
1754  forAll(regionSizes, regionI)
1755  {
1756  reduce(regionSizes[regionI], sumOp<label>());
1757  }
1758 
1759  Info<< "Region\tCells" << nl
1760  << "------\t-----" << endl;
1761 
1762  forAll(regionSizes, regionI)
1763  {
1764  Info<< regionI << '\t' << regionSizes[regionI] << nl;
1765  }
1766  Info<< endl;
1767 
1768 
1769 
1770  // Print region to zone
1771  Info<< "Region\tZone\tName" << nl
1772  << "------\t----\t----" << endl;
1773  forAll(regionToZone, regionI)
1774  {
1775  Info<< regionI << '\t' << regionToZone[regionI] << '\t'
1776  << regionNames[regionI] << nl;
1777  }
1778  Info<< endl;
1779 
1780 
1781 
1782  // Since we're going to mess with patches and zones make sure all
1783  // is synchronised
1784  mesh.boundaryMesh().checkParallelSync(true);
1785  mesh.faceZones().checkParallelSync(true);
1786 
1787 
1788  // Interfaces
1789  // ----------
1790  // per interface:
1791  // - the two regions on either side
1792  // - the name
1793  // - the (global) size
1794  edgeList interfaces;
1795  List<Pair<word>> interfaceNames;
1796  labelList interfaceSizes;
1797  // per face the interface
1798  labelList faceToInterface;
1799 
1800  getInterfaceSizes
1801  (
1802  mesh,
1803  useFaceZones,
1804  cellRegion,
1805  regionNames,
1806 
1807  interfaces,
1808  interfaceNames,
1809  interfaceSizes,
1810  faceToInterface
1811  );
1812 
1813  Info<< "Sizes of interfaces between regions:" << nl << nl
1814  << "Interface\tRegion\tRegion\tFaces" << nl
1815  << "---------\t------\t------\t-----" << endl;
1816 
1817  forAll(interfaces, interI)
1818  {
1819  const edge& e = interfaces[interI];
1820 
1821  Info<< interI
1822  << "\t\t" << e[0] << '\t' << e[1]
1823  << '\t' << interfaceSizes[interI] << nl;
1824  }
1825  Info<< endl;
1826 
1827 
1828  if (detectOnly)
1829  {
1830  return 0;
1831  }
1832 
1833 
1834  // Read objects in time directory
1835  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1836 
1837  IOobjectList objects(mesh, runTime.timeName());
1838 
1839  // Read vol fields.
1840 
1841  PtrList<volScalarField> vsFlds;
1842  ReadFields(mesh, objects, vsFlds);
1843 
1844  PtrList<volVectorField> vvFlds;
1845  ReadFields(mesh, objects, vvFlds);
1846 
1848  ReadFields(mesh, objects, vstFlds);
1849 
1850  PtrList<volSymmTensorField> vsymtFlds;
1851  ReadFields(mesh, objects, vsymtFlds);
1852 
1853  PtrList<volTensorField> vtFlds;
1854  ReadFields(mesh, objects, vtFlds);
1855 
1856  // Read surface fields.
1857 
1859  ReadFields(mesh, objects, ssFlds);
1860 
1862  ReadFields(mesh, objects, svFlds);
1863 
1865  ReadFields(mesh, objects, sstFlds);
1866 
1868  ReadFields(mesh, objects, ssymtFlds);
1869 
1871  ReadFields(mesh, objects, stFlds);
1872 
1873  Info<< endl;
1874 
1875 
1876  // Remove any demand-driven fields ('S', 'V' etc)
1877  mesh.clearOut();
1878 
1879 
1880  if (nCellRegions == 1)
1881  {
1882  Info<< "Only one region. Doing nothing." << endl;
1883  }
1884  else if (makeCellZones)
1885  {
1886  Info<< "Putting cells into cellZones instead of splitting mesh."
1887  << endl;
1888 
1889  // Check if region overlaps with existing zone. If so keep.
1890 
1891  for (label regionI = 0; regionI < nCellRegions; regionI++)
1892  {
1893  label zoneI = regionToZone[regionI];
1894 
1895  if (zoneI != -1)
1896  {
1897  Info<< " Region " << regionI << " : corresponds to existing"
1898  << " cellZone "
1899  << zoneI << ' ' << cellZones[zoneI].name() << endl;
1900  }
1901  else
1902  {
1903  // Create new cellZone.
1904  labelList regionCells = findIndices(cellRegion, regionI);
1905 
1906  word zoneName = "region" + Foam::name(regionI);
1907 
1908  zoneI = cellZones.findZoneID(zoneName);
1909 
1910  if (zoneI == -1)
1911  {
1912  zoneI = cellZones.size();
1913  mesh.cellZones().setSize(zoneI+1);
1914  mesh.cellZones().set
1915  (
1916  zoneI,
1917  new cellZone
1918  (
1919  zoneName, //name
1920  regionCells, //addressing
1921  zoneI, //index
1922  cellZones //cellZoneMesh
1923  )
1924  );
1925  }
1926  else
1927  {
1928  mesh.cellZones()[zoneI].clearAddressing();
1929  mesh.cellZones()[zoneI] = regionCells;
1930  }
1931  Info<< " Region " << regionI << " : created new cellZone "
1932  << zoneI << ' ' << cellZones[zoneI].name() << endl;
1933  }
1934  }
1936 
1937  if (!overwrite)
1938  {
1939  runTime++;
1940  mesh.setInstance(runTime.timeName());
1941  }
1942  else
1943  {
1944  mesh.setInstance(oldInstance);
1945  }
1946 
1947  Info<< "Writing cellZones as new mesh to time " << runTime.timeName()
1948  << nl << endl;
1949 
1950  mesh.write();
1951 
1952 
1953  // Write cellSets for convenience
1954  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1955 
1956  Info<< "Writing cellSets corresponding to cellZones." << nl << endl;
1957 
1958  forAll(cellZones, zoneI)
1959  {
1960  const cellZone& cz = cellZones[zoneI];
1961 
1962  cellSet(mesh, cz.name(), cz).write();
1963  }
1964  }
1965  else
1966  {
1967  // Add patches for interfaces
1968  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
1969 
1970  // Add all possible patches. Empty ones get filtered later on.
1971  Info<< nl << "Adding patches" << nl << endl;
1972 
1973  labelList interfacePatches
1974  (
1975  addRegionPatches
1976  (
1977  mesh,
1978  regionNames,
1979  interfaces,
1980  interfaceNames
1981  )
1982  );
1983 
1984 
1985  if (!overwrite)
1986  {
1987  runTime++;
1988  }
1989 
1990 
1991  // Create regions
1992  // ~~~~~~~~~~~~~~
1993 
1994  if (insidePoint)
1995  {
1996  const point insidePoint = args.optionRead<point>("insidePoint");
1997 
1998  label regionI = -1;
1999 
2000  (void)mesh.tetBasePtIs();
2001 
2002  label celli = mesh.findCell(insidePoint);
2003 
2004  Info<< nl << "Found point " << insidePoint << " in cell " << celli
2005  << endl;
2006 
2007  if (celli != -1)
2008  {
2009  regionI = cellRegion[celli];
2010  }
2011 
2012  reduce(regionI, maxOp<label>());
2013 
2014  Info<< nl
2015  << "Subsetting region " << regionI
2016  << " containing point " << insidePoint << endl;
2017 
2018  if (regionI == -1)
2019  {
2021  << "Point " << insidePoint
2022  << " is not inside the mesh." << nl
2023  << "Bounding box of the mesh:" << mesh.bounds()
2024  << exit(FatalError);
2025  }
2026 
2027  createAndWriteRegion
2028  (
2029  mesh,
2030  cellRegion,
2031  regionNames,
2032  prefixRegion,
2033  faceToInterface,
2034  interfacePatches,
2035  regionI,
2036  (overwrite ? oldInstance : runTime.timeName())
2037  );
2038  }
2039  else if (largestOnly)
2040  {
2041  label regionI = findMax(regionSizes);
2042 
2043  Info<< nl
2044  << "Subsetting region " << regionI
2045  << " of size " << regionSizes[regionI]
2046  << " as named region " << regionNames[regionI] << endl;
2047 
2048  createAndWriteRegion
2049  (
2050  mesh,
2051  cellRegion,
2052  regionNames,
2053  prefixRegion,
2054  faceToInterface,
2055  interfacePatches,
2056  regionI,
2057  (overwrite ? oldInstance : runTime.timeName())
2058  );
2059  }
2060  else
2061  {
2062  // Split all
2063  for (label regionI = 0; regionI < nCellRegions; regionI++)
2064  {
2065  Info<< nl
2066  << "Region " << regionI << nl
2067  << "-------- " << endl;
2068 
2069  createAndWriteRegion
2070  (
2071  mesh,
2072  cellRegion,
2073  regionNames,
2074  prefixRegion,
2075  faceToInterface,
2076  interfacePatches,
2077  regionI,
2078  (overwrite ? oldInstance : runTime.timeName())
2079  );
2080  }
2081  }
2082  }
2083 
2084  Info<< "End\n" << endl;
2085 
2086  return 0;
2087 }
2088 
2089 
2090 // ************************************************************************* //
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:116
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
static int masterNo()
Process index of the master.
Definition: UPstream.H:405
T optionRead(const word &opt) const
Read a value from the named option.
Definition: argListI.H:187
A list of face labels.
Definition: faceSet.H:48
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:427
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:106
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
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
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())
static int firstSlave()
Process index of first slave.
Definition: UPstream.H:434
PtrList< labelIOList > & faceProcAddressing
Foam::word regionName
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:417
Given list of cells to remove insert all the topology changes.
Definition: removeCells.H:59
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:411
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
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:715
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
writeOption writeOpt() const
Definition: IOobject.H:314
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order. Return.
Definition: ZoneMesh.C:428
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:737
label nCells() const
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Input inter-processor communications stream.
Definition: IPstream.H:50
bool insert(const edge &, const T &newElmt)
Insert a new hashedEntry.
iterator find(const edge &)
Find and return an iterator set at the hashedEntry.
const word & name() const
Return name.
Definition: zone.H:150
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:341
dynamicFvMesh & mesh
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
A class for handling words, derived from string.
Definition: word.H:59
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 ...
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:256
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
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:93
const vectorField & cellCentres() const
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.
An STL-conforming hash table.
Definition: HashTable.H:61
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
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
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).
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)
const labelList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:818
const string & option(const word &opt) const
Return the argument string associated with the named option.
Definition: argListI.H:102
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 nFaces() const
void setSize(const label)
Reset size of List.
Definition: List.C:295
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:393
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:399
label patchi
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
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:62
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
A collection of cell labels.
Definition: cellSet.H:48
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:221
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:1408
messageStream Info
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
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:53
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:83
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
const word & system() const
Return system name.
Definition: TimePaths.H:114
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:124
A class for managing temporary objects.
Definition: PtrList.H:54
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
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:91
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.
label nInternalFaces() const
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:365
runTime write()
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:870
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static int lastSlave(const label communicator=0)
Process index of last slave.
Definition: UPstream.H:440