All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
splitMeshRegions.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 in between 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 // Standard default region base name
113 const word standardRegionName("region");
114 
115 
116 // Prepend prefix to selected patches.
117 void renamePatches
118 (
119  fvMesh& mesh,
120  const word& prefix,
121  const labelList& patchesToRename
122 )
123 {
124  polyBoundaryMesh& patches =
125  const_cast<polyBoundaryMesh&>(mesh.boundaryMesh());
126 
127  // Create a list of all the new names
128  wordList newNames = patches.names();
129  forAll(patchesToRename, i)
130  {
131  const label patchi = patchesToRename[i];
132  newNames[patchi] = prefix + '_' + newNames[patchi];
133  }
134 
135  // Apply to the patches
136  patches.renamePatches(newNames, true);
137 }
138 
139 
140 template<class GeoField>
141 void subsetVolFields
142 (
143  const fvMesh& mesh,
144  const fvMesh& subMesh,
145  const labelList& cellMap,
146  const labelList& faceMap,
147  const labelHashSet& addedPatches
148 )
149 {
150  const labelList patchMap(identity(mesh.boundaryMesh().size()));
151 
153  (
154  mesh.objectRegistry::lookupClass<GeoField>()
155  );
156  forAllConstIter(typename HashTable<const GeoField*>, fields, iter)
157  {
158  const GeoField& fld = *iter();
159 
160  Info<< "Mapping field " << fld.name() << endl;
161 
162  tmp<GeoField> tSubFld
163  (
165  (
166  fld,
167  subMesh,
168  patchMap,
169  cellMap,
170  faceMap
171  )
172  );
173 
174  // Hack: set value to 0 for introduced patches (since don't
175  // get initialised.
176  forAll(tSubFld().boundaryField(), patchi)
177  {
178  if (addedPatches.found(patchi))
179  {
180  tSubFld.ref().boundaryFieldRef()[patchi] ==
181  typename GeoField::value_type(Zero);
182  }
183  }
184 
185  // Store on subMesh
186  GeoField* subFld = tSubFld.ptr();
187  subFld->rename(fld.name());
188  subFld->writeOpt() = IOobject::AUTO_WRITE;
189  subFld->store();
190  }
191 }
192 
193 
194 template<class GeoField>
195 void subsetSurfaceFields
196 (
197  const fvMesh& mesh,
198  const fvMesh& subMesh,
199  const labelList& cellMap,
200  const labelList& faceMap,
201  const labelHashSet& addedPatches
202 )
203 {
204  const labelList patchMap(identity(mesh.boundaryMesh().size()));
205 
207  (
208  mesh.objectRegistry::lookupClass<GeoField>()
209  );
210  forAllConstIter(typename HashTable<const GeoField*>, fields, iter)
211  {
212  const GeoField& fld = *iter();
213 
214  Info<< "Mapping field " << fld.name() << endl;
215 
216  tmp<GeoField> tSubFld
217  (
219  (
220  fld,
221  subMesh,
222  patchMap,
223  cellMap,
224  faceMap
225  )
226  );
227 
228  // Hack: set value to 0 for introduced patches (since don't
229  // get initialised.
230  forAll(tSubFld().boundaryField(), patchi)
231  {
232  if (addedPatches.found(patchi))
233  {
234  tSubFld.ref().boundaryFieldRef()[patchi] ==
235  typename GeoField::value_type(Zero);
236  }
237  }
238 
239  // Store on subMesh
240  GeoField* subFld = tSubFld.ptr();
241  subFld->rename(fld.name());
242  subFld->writeOpt() = IOobject::AUTO_WRITE;
243  subFld->store();
244  }
245 }
246 
247 // Select all cells not in the region
248 labelList getNonRegionCells(const labelList& cellRegion, const label regioni)
249 {
250  DynamicList<label> nonRegionCells(cellRegion.size());
251  forAll(cellRegion, celli)
252  {
253  if (cellRegion[celli] != regioni)
254  {
255  nonRegionCells.append(celli);
256  }
257  }
258  return nonRegionCells.shrink();
259 }
260 
261 
262 void addToInterface
263 (
264  const polyMesh& mesh,
265  const label zoneID,
266  const label ownRegion,
267  const label neiRegion,
268  EdgeMap<Map<label>>& regionsToSize
269 )
270 {
271  edge interface
272  (
273  min(ownRegion, neiRegion),
274  max(ownRegion, neiRegion)
275  );
276 
277  EdgeMap<Map<label>>::iterator iter = regionsToSize.find
278  (
279  interface
280  );
281 
282  if (iter != regionsToSize.end())
283  {
284  // Check if zone present
285  Map<label>::iterator zoneFnd = iter().find(zoneID);
286  if (zoneFnd != iter().end())
287  {
288  // Found zone. Increment count.
289  zoneFnd()++;
290  }
291  else
292  {
293  // New or no zone. Insert with count 1.
294  iter().insert(zoneID, 1);
295  }
296  }
297  else
298  {
299  // Create new interface of size 1.
300  Map<label> zoneToSize;
301  zoneToSize.insert(zoneID, 1);
302  regionsToSize.insert(interface, zoneToSize);
303  }
304 }
305 
306 
307 // Get region-region interface name and sizes.
308 // Returns interfaces as straight list for looping in identical order.
309 void getInterfaceSizes
310 (
311  const polyMesh& mesh,
312  const bool useFaceZones,
313  const labelList& cellRegion,
314  const wordList& regionNames,
315 
316  edgeList& interfaces,
317  List<Pair<word>>& interfaceNames,
318  labelList& interfaceSizes,
319  labelList& faceToInterface
320 )
321 {
322  // From region-region to faceZone (or -1) to number of faces.
323 
324  EdgeMap<Map<label>> regionsToSize;
325 
326 
327  // Internal faces
328  // ~~~~~~~~~~~~~~
329 
330  forAll(mesh.faceNeighbour(), facei)
331  {
332  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
333  label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
334 
335  if (ownRegion != neiRegion)
336  {
337  addToInterface
338  (
339  mesh,
340  (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
341  ownRegion,
342  neiRegion,
343  regionsToSize
344  );
345  }
346  }
347 
348  // Boundary faces
349  // ~~~~~~~~~~~~~~
350 
351  // Neighbour cellRegion.
352  labelList coupledRegion(mesh.nFaces()-mesh.nInternalFaces());
353 
354  forAll(coupledRegion, i)
355  {
356  label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
357  coupledRegion[i] = cellRegion[celli];
358  }
359  syncTools::swapBoundaryFaceList(mesh, coupledRegion);
360 
361  forAll(coupledRegion, i)
362  {
363  label facei = i+mesh.nInternalFaces();
364  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
365  label neiRegion = coupledRegion[i];
366 
367  if (ownRegion != neiRegion)
368  {
369  addToInterface
370  (
371  mesh,
372  (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
373  ownRegion,
374  neiRegion,
375  regionsToSize
376  );
377  }
378  }
379 
380 
381  if (Pstream::parRun())
382  {
383  if (Pstream::master())
384  {
385  // Receive and add to my sizes
386  for
387  (
388  int slave=Pstream::firstSlave();
389  slave<=Pstream::lastSlave();
390  slave++
391  )
392  {
393  IPstream fromSlave(Pstream::commsTypes::blocking, slave);
394 
395  EdgeMap<Map<label>> slaveSizes(fromSlave);
396 
397  forAllConstIter(EdgeMap<Map<label>>, slaveSizes, slaveIter)
398  {
399  EdgeMap<Map<label>>::iterator masterIter =
400  regionsToSize.find(slaveIter.key());
401 
402  if (masterIter != regionsToSize.end())
403  {
404  // Same inter-region
405  const Map<label>& slaveInfo = slaveIter();
406  Map<label>& masterInfo = masterIter();
407 
408  forAllConstIter(Map<label>, slaveInfo, iter)
409  {
410  label zoneID = iter.key();
411  label slaveSize = iter();
412 
413  Map<label>::iterator zoneFnd = masterInfo.find
414  (
415  zoneID
416  );
417  if (zoneFnd != masterInfo.end())
418  {
419  zoneFnd() += slaveSize;
420  }
421  else
422  {
423  masterInfo.insert(zoneID, slaveSize);
424  }
425  }
426  }
427  else
428  {
429  regionsToSize.insert(slaveIter.key(), slaveIter());
430  }
431  }
432  }
433  }
434  else
435  {
436  // Send to master
437  {
438  OPstream toMaster
439  (
442  );
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 
537  edge interface
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 
560  edge interface
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.typeHeaderOk<IOdictionary>(true))
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.typeHeaderOk<IOdictionary>(true))
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 meshCellZones& 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  const word& defaultRegionName,
1248 
1249  labelList& regionToZone,
1250  wordList& regionNames,
1251  labelList& zoneToRegion
1252 )
1253 {
1254  const meshCellZones& cellZones = mesh.cellZones();
1255 
1256  regionToZone.setSize(nCellRegions, -1);
1257  regionNames.setSize(nCellRegions);
1258  zoneToRegion.setSize(cellZones.size(), -1);
1259 
1260  // Get current per cell zoneID
1261  labelList zoneID(mesh.nCells(), -1);
1262  labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1263  getZoneID(mesh, cellZones, zoneID, neiZoneID);
1264 
1265  // Sizes per cellzone
1266  labelList zoneSizes(cellZones.size(), 0);
1267  {
1268  List<wordList> zoneNames(Pstream::nProcs());
1269  zoneNames[Pstream::myProcNo()] = cellZones.names();
1270  Pstream::gatherList(zoneNames);
1271  Pstream::scatterList(zoneNames);
1272 
1273  forAll(zoneNames, proci)
1274  {
1275  if (zoneNames[proci] != zoneNames[0])
1276  {
1278  << "cellZones not synchronised across processors." << endl
1279  << "Master has cellZones " << zoneNames[0] << endl
1280  << "Processor " << proci
1281  << " has cellZones " << zoneNames[proci]
1282  << exit(FatalError);
1283  }
1284  }
1285 
1286  forAll(cellZones, zoneI)
1287  {
1288  zoneSizes[zoneI] = returnReduce
1289  (
1290  cellZones[zoneI].size(),
1291  sumOp<label>()
1292  );
1293  }
1294  }
1295 
1296 
1297  if (sloppyCellZones)
1298  {
1299  Info<< "Trying to match regions to existing cell zones;"
1300  << " region can be subset of cell zone." << nl << endl;
1301 
1302  forAll(cellZones, zoneI)
1303  {
1304  label regioni = findCorrespondingRegion
1305  (
1306  zoneID,
1307  cellRegion,
1308  nCellRegions,
1309  zoneI,
1310  label(0.5*zoneSizes[zoneI]) // minimum overlap
1311  );
1312 
1313  if (regioni != -1)
1314  {
1315  Info<< "Sloppily matched region " << regioni
1316  //<< " size " << regionSizes[regioni]
1317  << " to zone " << zoneI << " size " << zoneSizes[zoneI]
1318  << endl;
1319  zoneToRegion[zoneI] = regioni;
1320  regionToZone[regioni] = zoneI;
1321  regionNames[regioni] = cellZones[zoneI].name();
1322  }
1323  }
1324  }
1325  else
1326  {
1327  Info<< "Trying to match regions to existing cell zones." << nl << endl;
1328 
1329  forAll(cellZones, zoneI)
1330  {
1331  label regioni = findCorrespondingRegion
1332  (
1333  zoneID,
1334  cellRegion,
1335  nCellRegions,
1336  zoneI,
1337  1 // minimum overlap
1338  );
1339 
1340  if (regioni != -1)
1341  {
1342  zoneToRegion[zoneI] = regioni;
1343  regionToZone[regioni] = zoneI;
1344  regionNames[regioni] = cellZones[zoneI].name();
1345  }
1346  }
1347  }
1348 
1349  label nUnmatchedRegions = 0;
1350 
1351  forAll(regionToZone, regioni)
1352  {
1353  if (regionToZone[regioni] == -1)
1354  {
1355  nUnmatchedRegions++;
1356  }
1357  }
1358 
1359  if (nUnmatchedRegions)
1360  {
1361  label nUnmatchedi = 1;
1362 
1363  // Allocate region names for unmatched regions.
1364  forAll(regionToZone, regioni)
1365  {
1366  if (regionToZone[regioni] == -1)
1367  {
1368  if
1369  (
1370  nUnmatchedRegions == 1
1371  && defaultRegionName != standardRegionName
1372  )
1373  {
1374  regionNames[regioni] = defaultRegionName;
1375  }
1376  else
1377  {
1378  regionNames[regioni] =
1379  defaultRegionName + Foam::name(nUnmatchedi++);
1380  }
1381  }
1382  }
1383  }
1384 }
1385 
1386 
1387 void writeCellToRegion(const fvMesh& mesh, const labelList& cellRegion)
1388 {
1389  // Write to manual decomposition option
1390  {
1391  labelIOList cellToRegion
1392  (
1393  IOobject
1394  (
1395  "cellToRegion",
1396  mesh.facesInstance(),
1397  mesh,
1400  false
1401  ),
1402  cellRegion
1403  );
1404  cellToRegion.write();
1405 
1406  Info<< "Writing region per cell file (for manual decomposition) to "
1407  << cellToRegion.localObjectPath() << nl << endl;
1408  }
1409  // Write for postprocessing
1410  {
1411  volScalarField cellToRegion
1412  (
1413  IOobject
1414  (
1415  "cellToRegion",
1416  mesh.time().timeName(),
1417  mesh,
1420  false
1421  ),
1422  mesh,
1424  zeroGradientFvPatchScalarField::typeName
1425  );
1426  forAll(cellRegion, celli)
1427  {
1428  cellToRegion[celli] = cellRegion[celli];
1429  }
1430  cellToRegion.write();
1431 
1432  Info<< "Writing region per cell as volScalarField to "
1433  << cellToRegion.localObjectPath() << nl << endl;
1434  }
1435 }
1436 
1437 
1438 
1439 
1440 int main(int argc, char *argv[])
1441 {
1443  (
1444  "splits mesh into multiple regions (detected by walking across faces)"
1445  );
1446  #include "addRegionOption.H"
1447  #include "addOverwriteOption.H"
1449  (
1450  "cellZones",
1451  "additionally split cellZones off into separate regions"
1452  );
1454  (
1455  "cellZonesOnly",
1456  "use cellZones only to split mesh into regions; do not use walking"
1457  );
1459  (
1460  "cellZonesFileOnly",
1461  "file",
1462  "like -cellZonesOnly, but use specified file"
1463  );
1465  (
1466  "blockedFaces",
1467  "faceSet",
1468  "specify additional region boundaries that walking does not cross"
1469  );
1471  (
1472  "makeCellZones",
1473  "place cells into cellZones instead of splitting mesh"
1474  );
1476  (
1477  "largestOnly",
1478  "only write largest region"
1479  );
1481  (
1482  "insidePoint",
1483  "point",
1484  "only write region containing point"
1485  );
1487  (
1488  "detectOnly",
1489  "do not write mesh"
1490  );
1492  (
1493  "sloppyCellZones",
1494  "try to match heuristically regions to existing cell zones"
1495  );
1497  (
1498  "useFaceZones",
1499  "use faceZones to patch inter-region faces instead of single patch"
1500  );
1502  (
1503  "prefixRegion",
1504  "prefix region name to all patches, not just coupling patches"
1505  );
1507  (
1508  "defaultRegionName",
1509  "name",
1510  "base name of the unspecified regions, defaults to \"region\""
1511  );
1513  (
1514  "noFields",
1515  "do not update fields"
1516  );
1517 
1518  #include "setRootCase.H"
1519  #include "createTime.H"
1520 
1522 
1523  #include "createNamedMesh.H"
1524 
1525  const word oldInstance = mesh.pointsInstance();
1526 
1527  word blockedFacesName;
1528  if (args.optionReadIfPresent("blockedFaces", blockedFacesName))
1529  {
1530  Info<< "Reading blocked internal faces from faceSet "
1531  << blockedFacesName << nl << endl;
1532  }
1533 
1534  const bool makeCellZones = args.optionFound("makeCellZones");
1535  const bool largestOnly = args.optionFound("largestOnly");
1536  const bool insidePoint = args.optionFound("insidePoint");
1537  const bool useCellZones = args.optionFound("cellZones");
1538  const bool useCellZonesOnly = args.optionFound("cellZonesOnly");
1539  const bool useCellZonesFile = args.optionFound("cellZonesFileOnly");
1540  const bool overwrite = args.optionFound("overwrite");
1541  const bool detectOnly = args.optionFound("detectOnly");
1542  const bool sloppyCellZones = args.optionFound("sloppyCellZones");
1543  const bool useFaceZones = args.optionFound("useFaceZones");
1544  const bool prefixRegion = args.optionFound("prefixRegion");
1545  const bool fields = !args.optionFound("noFields");
1546 
1547  if
1548  (
1549  (useCellZonesOnly || useCellZonesFile)
1550  && (useCellZones || blockedFacesName.size())
1551  )
1552  {
1554  << "You cannot specify both -cellZonesOnly or -cellZonesFileOnly"
1555  << " (which specify complete split)"
1556  << " in combination with -blockedFaces or -cellZones"
1557  << " (which imply a split based on topology)"
1558  << exit(FatalError);
1559  }
1560 
1561 
1562  if (useFaceZones)
1563  {
1564  Info<< "Using current faceZones to divide inter-region interfaces"
1565  << " into multiple patches."
1566  << nl << endl;
1567  }
1568  else
1569  {
1570  Info<< "Creating single patch per inter-region interface."
1571  << nl << endl;
1572  }
1573 
1574 
1575  if (insidePoint && largestOnly)
1576  {
1578  << "You cannot specify both -largestOnly"
1579  << " (keep region with most cells)"
1580  << " and -insidePoint (keep region containing point)"
1581  << exit(FatalError);
1582  }
1583 
1584  const word defaultRegionName
1585  (
1586  args.optionLookupOrDefault("defaultRegionName", standardRegionName)
1587  );
1588 
1589  const meshCellZones& cellZones = mesh.cellZones();
1590 
1591  // Existing zoneID
1592  labelList zoneID(mesh.nCells(), -1);
1593  // Neighbour zoneID.
1594  labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1595  getZoneID(mesh, cellZones, zoneID, neiZoneID);
1596 
1597 
1598 
1599  // Determine per cell the region it belongs to
1600  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1601 
1602  // cellRegion is the labelList with the region per cell.
1603  labelList cellRegion;
1604  // Region per zone
1605  labelList regionToZone;
1606  // Name of region
1607  wordList regionNames;
1608  // Zone to region
1609  labelList zoneToRegion;
1610 
1611  label nCellRegions = 0;
1612  if (useCellZonesOnly)
1613  {
1614  Info<< "Using current cellZones to split mesh into regions."
1615  << " This requires all"
1616  << " cells to be in one and only one cellZone." << nl << endl;
1617 
1618  label unzonedCelli = findIndex(zoneID, -1);
1619  if (unzonedCelli != -1)
1620  {
1622  << "For the cellZonesOnly option all cells "
1623  << "have to be in a cellZone." << endl
1624  << "Cell " << unzonedCelli
1625  << " at" << mesh.cellCentres()[unzonedCelli]
1626  << " is not in a cellZone. There might be more unzoned cells."
1627  << exit(FatalError);
1628  }
1629  cellRegion = zoneID;
1630  nCellRegions = gMax(cellRegion)+1;
1631  regionToZone.setSize(nCellRegions);
1632  regionNames.setSize(nCellRegions);
1633  zoneToRegion.setSize(cellZones.size(), -1);
1634  for (label regioni = 0; regioni < nCellRegions; regioni++)
1635  {
1636  regionToZone[regioni] = regioni;
1637  zoneToRegion[regioni] = regioni;
1638  regionNames[regioni] = cellZones[regioni].name();
1639  }
1640  }
1641  else if (useCellZonesFile)
1642  {
1643  const word zoneFile = args.option("cellZonesFileOnly");
1644  Info<< "Reading split from cellZones file " << zoneFile << endl
1645  << "This requires all"
1646  << " cells to be in one and only one cellZone." << nl << endl;
1647 
1648  meshCellZones newCellZones
1649  (
1650  IOobject
1651  (
1652  zoneFile,
1653  mesh.facesInstance(),
1655  mesh,
1658  false
1659  ),
1660  mesh
1661  );
1662 
1663  labelList newZoneID(mesh.nCells(), -1);
1664  labelList newNeiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1665  getZoneID(mesh, newCellZones, newZoneID, newNeiZoneID);
1666 
1667  label unzonedCelli = findIndex(newZoneID, -1);
1668  if (unzonedCelli != -1)
1669  {
1671  << "For the cellZonesFileOnly option all cells "
1672  << "have to be in a cellZone." << endl
1673  << "Cell " << unzonedCelli
1674  << " at" << mesh.cellCentres()[unzonedCelli]
1675  << " is not in a cellZone. There might be more unzoned cells."
1676  << exit(FatalError);
1677  }
1678  cellRegion = newZoneID;
1679  nCellRegions = gMax(cellRegion)+1;
1680  zoneToRegion.setSize(newCellZones.size(), -1);
1681  regionToZone.setSize(nCellRegions);
1682  regionNames.setSize(nCellRegions);
1683  for (label regioni = 0; regioni < nCellRegions; regioni++)
1684  {
1685  regionToZone[regioni] = regioni;
1686  zoneToRegion[regioni] = regioni;
1687  regionNames[regioni] = newCellZones[regioni].name();
1688  }
1689  }
1690  else
1691  {
1692  // Determine connected regions
1693  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1694 
1695  // Mark additional faces that are blocked
1696  boolList blockedFace;
1697 
1698  // Read from faceSet
1699  if (blockedFacesName.size())
1700  {
1701  faceSet blockedFaceSet(mesh, blockedFacesName);
1702  Info<< "Read "
1703  << returnReduce(blockedFaceSet.size(), sumOp<label>())
1704  << " blocked faces from set " << blockedFacesName << nl << endl;
1705 
1706  blockedFace.setSize(mesh.nFaces(), false);
1707 
1708  forAllConstIter(faceSet, blockedFaceSet, iter)
1709  {
1710  blockedFace[iter.key()] = true;
1711  }
1712  }
1713 
1714  // Imply from differing cellZones
1715  if (useCellZones)
1716  {
1717  blockedFace.setSize(mesh.nFaces(), false);
1718 
1719  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1720  {
1721  label own = mesh.faceOwner()[facei];
1722  label nei = mesh.faceNeighbour()[facei];
1723 
1724  if (zoneID[own] != zoneID[nei])
1725  {
1726  blockedFace[facei] = true;
1727  }
1728  }
1729 
1730  // Different cellZones on either side of processor patch.
1731  forAll(neiZoneID, i)
1732  {
1733  label facei = i+mesh.nInternalFaces();
1734 
1735  if (zoneID[mesh.faceOwner()[facei]] != neiZoneID[i])
1736  {
1737  blockedFace[facei] = true;
1738  }
1739  }
1740  }
1741 
1742  // Do a topological walk to determine regions
1743  regionSplit regions(mesh, blockedFace);
1744  nCellRegions = regions.nRegions();
1745  cellRegion.transfer(regions);
1746 
1747  // Make up region names. If possible match them to existing zones.
1748  matchRegions
1749  (
1750  sloppyCellZones,
1751  mesh,
1752  nCellRegions,
1753  cellRegion,
1754  defaultRegionName,
1755 
1756  regionToZone,
1757  regionNames,
1758  zoneToRegion
1759  );
1760 
1761  // Override any default region names if single region selected
1762  if (largestOnly || insidePoint)
1763  {
1764  forAll(regionToZone, regioni)
1765  {
1766  if (regionToZone[regioni] == -1)
1767  {
1768  if (overwrite)
1769  {
1770  regionNames[regioni] = polyMesh::defaultRegion;
1771  }
1772  else if (insidePoint)
1773  {
1774  regionNames[regioni] = "insidePoint";
1775  }
1776  else if (largestOnly)
1777  {
1778  regionNames[regioni] = "largestOnly";
1779  }
1780  }
1781  }
1782  }
1783  }
1784 
1785  Info<< endl << "Number of regions:" << nCellRegions << nl << endl;
1786 
1787 
1788  // Write decomposition to file
1789  writeCellToRegion(mesh, cellRegion);
1790 
1791 
1792 
1793  // Sizes per region
1794  // ~~~~~~~~~~~~~~~~
1795 
1796  labelList regionSizes(nCellRegions, 0);
1797 
1798  forAll(cellRegion, celli)
1799  {
1800  regionSizes[cellRegion[celli]]++;
1801  }
1802  forAll(regionSizes, regioni)
1803  {
1804  reduce(regionSizes[regioni], sumOp<label>());
1805  }
1806 
1807  Info<< "Region\tCells" << nl
1808  << "------\t-----" << endl;
1809 
1810  forAll(regionSizes, regioni)
1811  {
1812  Info<< regioni << '\t' << regionSizes[regioni] << nl;
1813  }
1814  Info<< endl;
1815 
1816 
1817 
1818  // Print region to zone
1819  Info<< "Region\tZone\tName" << nl
1820  << "------\t----\t----" << endl;
1821  forAll(regionToZone, regioni)
1822  {
1823  Info<< regioni << '\t' << regionToZone[regioni] << '\t'
1824  << regionNames[regioni] << nl;
1825  }
1826  Info<< endl;
1827 
1828 
1829 
1830  // Since we're going to mess with patches and zones make sure all
1831  // is synchronised
1832  mesh.boundaryMesh().checkParallelSync(true);
1833  mesh.faceZones().checkParallelSync(true);
1834 
1835 
1836  // Interfaces
1837  // ----------
1838  // per interface:
1839  // - the two regions on either side
1840  // - the name
1841  // - the (global) size
1842  edgeList interfaces;
1843  List<Pair<word>> interfaceNames;
1844  labelList interfaceSizes;
1845  // per face the interface
1846  labelList faceToInterface;
1847 
1848  getInterfaceSizes
1849  (
1850  mesh,
1851  useFaceZones,
1852  cellRegion,
1853  regionNames,
1854 
1855  interfaces,
1856  interfaceNames,
1857  interfaceSizes,
1858  faceToInterface
1859  );
1860 
1861  Info<< "Sizes of interfaces between regions:" << nl << nl
1862  << "Interface\tRegion\tRegion\tFaces" << nl
1863  << "---------\t------\t------\t-----" << endl;
1864 
1865  forAll(interfaces, interI)
1866  {
1867  const edge& e = interfaces[interI];
1868 
1869  Info<< interI
1870  << "\t\t" << e[0] << '\t' << e[1]
1871  << '\t' << interfaceSizes[interI] << nl;
1872  }
1873  Info<< endl;
1874 
1875 
1876  if (detectOnly)
1877  {
1878  return 0;
1879  }
1880 
1881 
1882  // Read objects in time directory
1883  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1884 
1886 
1887  if (fields) Info<< "Reading geometric fields" << nl << endl;
1888 
1889  #include "readVolFields.H"
1890  #include "readSurfaceFields.H"
1891  // #include "readPointFields.H"
1892 
1893  Info<< endl;
1894 
1895 
1896  // Remove any demand-driven fields ('S', 'V' etc)
1897  mesh.clearOut();
1898 
1899 
1900  if (nCellRegions == 1)
1901  {
1902  Info<< "Only one region. Doing nothing." << endl;
1903  }
1904  else if (makeCellZones)
1905  {
1906  Info<< "Putting cells into cellZones instead of splitting mesh."
1907  << endl;
1908 
1909  // Check if region overlaps with existing zone. If so keep.
1910 
1911  for (label regioni = 0; regioni < nCellRegions; regioni++)
1912  {
1913  label zoneI = regionToZone[regioni];
1914 
1915  if (zoneI != -1)
1916  {
1917  Info<< " Region " << regioni << " : corresponds to existing"
1918  << " cellZone "
1919  << zoneI << ' ' << cellZones[zoneI].name() << endl;
1920  }
1921  else
1922  {
1923  // Create new cellZone.
1924  labelList regionCells = findIndices(cellRegion, regioni);
1925 
1926  word zoneName = "region" + Foam::name(regioni);
1927 
1928  zoneI = cellZones.findZoneID(zoneName);
1929 
1930  if (zoneI == -1)
1931  {
1932  zoneI = cellZones.size();
1933  mesh.cellZones().setSize(zoneI+1);
1934  mesh.cellZones().set
1935  (
1936  zoneI,
1937  new cellZone
1938  (
1939  zoneName, // name
1940  regionCells, // addressing
1941  zoneI, // index
1942  cellZones // meshCellZones
1943  )
1944  );
1945  }
1946  else
1947  {
1948  mesh.cellZones()[zoneI].clearAddressing();
1949  mesh.cellZones()[zoneI] = regionCells;
1950  }
1951  Info<< " Region " << regioni << " : created new cellZone "
1952  << zoneI << ' ' << cellZones[zoneI].name() << endl;
1953  }
1954  }
1956 
1957  if (!overwrite)
1958  {
1959  runTime++;
1960  mesh.setInstance(runTime.timeName());
1961  }
1962  else
1963  {
1964  mesh.setInstance(oldInstance);
1965  }
1966 
1967  Info<< "Writing cellZones as new mesh to time " << runTime.timeName()
1968  << nl << endl;
1969 
1970  mesh.write();
1971 
1972 
1973  // Write cellSets for convenience
1974  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1975 
1976  Info<< "Writing cellSets corresponding to cellZones." << nl << endl;
1977 
1978  forAll(cellZones, zoneI)
1979  {
1980  const cellZone& cz = cellZones[zoneI];
1981 
1982  cellSet(mesh, cz.name(), cz).write();
1983  }
1984  }
1985  else
1986  {
1987  // Add patches for interfaces
1988  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
1989 
1990  // Add all possible patches. Empty ones get filtered later on.
1991  Info<< nl << "Adding patches" << nl << endl;
1992 
1993  labelList interfacePatches
1994  (
1995  addRegionPatches
1996  (
1997  mesh,
1998  regionNames,
1999  interfaces,
2000  interfaceNames
2001  )
2002  );
2003 
2004 
2005  if (!overwrite)
2006  {
2007  runTime++;
2008  }
2009 
2010 
2011  // Create regions
2012  // ~~~~~~~~~~~~~~
2013 
2014  if (insidePoint)
2015  {
2016  const point insidePoint = args.optionRead<point>("insidePoint");
2017 
2018  label regioni = -1;
2019 
2020  (void)mesh.tetBasePtIs();
2021 
2022  label celli = mesh.findCell(insidePoint);
2023 
2024  Info<< nl << "Found point " << insidePoint << " in cell " << celli
2025  << endl;
2026 
2027  if (celli != -1)
2028  {
2029  regioni = cellRegion[celli];
2030  }
2031 
2032  reduce(regioni, maxOp<label>());
2033 
2034  Info<< nl
2035  << "Subsetting region " << regioni
2036  << " containing point " << insidePoint << endl;
2037 
2038  if (regioni == -1)
2039  {
2041  << "Point " << insidePoint
2042  << " is not inside the mesh." << nl
2043  << "Bounding box of the mesh:" << mesh.bounds()
2044  << exit(FatalError);
2045  }
2046 
2047  createAndWriteRegion
2048  (
2049  mesh,
2050  cellRegion,
2051  regionNames,
2052  prefixRegion,
2053  faceToInterface,
2054  interfacePatches,
2055  regioni,
2056  (overwrite ? oldInstance : runTime.timeName())
2057  );
2058  }
2059  else if (largestOnly)
2060  {
2061  label regioni = findMax(regionSizes);
2062 
2063  Info<< nl
2064  << "Subsetting region " << regioni
2065  << " of size " << regionSizes[regioni]
2066  << " as named region " << regionNames[regioni] << endl;
2067 
2068  createAndWriteRegion
2069  (
2070  mesh,
2071  cellRegion,
2072  regionNames,
2073  prefixRegion,
2074  faceToInterface,
2075  interfacePatches,
2076  regioni,
2077  (overwrite ? oldInstance : runTime.timeName())
2078  );
2079  }
2080  else
2081  {
2082  // Split all
2083  for (label regioni = 0; regioni < nCellRegions; regioni++)
2084  {
2085  Info<< nl
2086  << "Region " << regioni << nl
2087  << "-------- " << endl;
2088 
2089  createAndWriteRegion
2090  (
2091  mesh,
2092  cellRegion,
2093  regionNames,
2094  prefixRegion,
2095  faceToInterface,
2096  interfacePatches,
2097  regioni,
2098  (overwrite ? oldInstance : runTime.timeName())
2099  );
2100  }
2101  }
2102  }
2103 
2104  Info<< "End\n" << endl;
2105 
2106  return 0;
2107 }
2108 
2109 
2110 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
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:240
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:140
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
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:417
const word & name() const
Return name.
Definition: IOobject.H:303
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
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:808
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
void off()
Switch the function objects off.
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:482
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1175
void renamePatches(const wordUList &newNames, const bool validBoundary)
Rename the patches. If validBoundary is set this calls updateMesh()
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:848
static int firstSlave()
Process index of first slave.
Definition: UPstream.H:446
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:164
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Given list of cells to remove insert all the topology changes.
Definition: removeCells.H:59
label nCells() const
engineTime & runTime
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order. Return.
Definition: MeshZones.C:428
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
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:636
const dimensionSet dimless
Info<< "Calculating turbulent flame speed field St\"<< endl;volScalarField St(IOobject("St", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:228
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurrences of given element. Linear search.
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1093
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
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:204
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:1133
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
wordList names() const
Return a list of zone names.
Definition: MeshZones.C:256
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
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:243
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:802
dynamicFvMesh & mesh
T optionRead(const word &opt) const
Read a value from the named option.
Definition: argListI.H:193
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
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:33
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 patch names.
This boundary condition is not designed to be evaluated; it is assumed 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:128
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
const word & name() const
Return name.
Definition: zone.H:147
An STL-conforming hash table.
Definition: HashTable.H:61
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:113
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:260
Map from edge (expressed as its endpoints) to value.
Definition: EdgeMap.H:47
Type gMax(const FieldField< Field, Type > &f)
objects
Output inter-processor communications stream.
Definition: OPstream.H:50
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence 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 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
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: MeshZones.C:221
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:440
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.
writeOption writeOpt() const
Definition: IOobject.H:363
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
label patchi
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:410
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A collection of cell labels.
Definition: cellSet.H:48
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:476
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:1606
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
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
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:118
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:159
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
void clearAddressing()
Clear addressing.
Definition: MeshZones.C:387
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
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:108
Namespace for OpenFOAM.
static int lastSlave(const label communicator=0)
Process index of last slave.
Definition: UPstream.H:452