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-2022 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  // Store on subMesh
229  GeoField* subFld = tSubFld.ptr();
230  subFld->rename(fld.name());
231  subFld->writeOpt() = IOobject::AUTO_WRITE;
232  subFld->store();
233  }
234 }
235 
236 // Select all cells not in the region
237 labelList getNonRegionCells(const labelList& cellRegion, const label regioni)
238 {
239  DynamicList<label> nonRegionCells(cellRegion.size());
240  forAll(cellRegion, celli)
241  {
242  if (cellRegion[celli] != regioni)
243  {
244  nonRegionCells.append(celli);
245  }
246  }
247  return nonRegionCells.shrink();
248 }
249 
250 
251 void addToInterface
252 (
253  const polyMesh& mesh,
254  const label zoneID,
255  const label ownRegion,
256  const label neiRegion,
257  EdgeMap<Map<label>>& regionsToSize
258 )
259 {
260  edge interface
261  (
262  min(ownRegion, neiRegion),
263  max(ownRegion, neiRegion)
264  );
265 
266  EdgeMap<Map<label>>::iterator iter = regionsToSize.find
267  (
268  interface
269  );
270 
271  if (iter != regionsToSize.end())
272  {
273  // Check if zone present
274  Map<label>::iterator zoneFnd = iter().find(zoneID);
275  if (zoneFnd != iter().end())
276  {
277  // Found zone. Increment count.
278  zoneFnd()++;
279  }
280  else
281  {
282  // New or no zone. Insert with count 1.
283  iter().insert(zoneID, 1);
284  }
285  }
286  else
287  {
288  // Create new interface of size 1.
289  Map<label> zoneToSize;
290  zoneToSize.insert(zoneID, 1);
291  regionsToSize.insert(interface, zoneToSize);
292  }
293 }
294 
295 
296 // Get region-region interface name and sizes.
297 // Returns interfaces as straight list for looping in identical order.
298 void getInterfaceSizes
299 (
300  const polyMesh& mesh,
301  const bool useFaceZones,
302  const labelList& cellRegion,
303  const wordList& regionNames,
304 
305  edgeList& interfaces,
306  List<Pair<word>>& interfaceNames,
307  labelList& interfaceSizes,
308  labelList& faceToInterface
309 )
310 {
311  // From region-region to faceZone (or -1) to number of faces.
312 
313  EdgeMap<Map<label>> regionsToSize;
314 
315 
316  // Internal faces
317  // ~~~~~~~~~~~~~~
318 
319  forAll(mesh.faceNeighbour(), facei)
320  {
321  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
322  label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
323 
324  if (ownRegion != neiRegion)
325  {
326  addToInterface
327  (
328  mesh,
329  (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
330  ownRegion,
331  neiRegion,
332  regionsToSize
333  );
334  }
335  }
336 
337  // Boundary faces
338  // ~~~~~~~~~~~~~~
339 
340  // Neighbour cellRegion.
341  labelList coupledRegion(mesh.nFaces()-mesh.nInternalFaces());
342 
343  forAll(coupledRegion, i)
344  {
345  label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
346  coupledRegion[i] = cellRegion[celli];
347  }
348  syncTools::swapBoundaryFaceList(mesh, coupledRegion);
349 
350  forAll(coupledRegion, i)
351  {
352  label facei = i+mesh.nInternalFaces();
353  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
354  label neiRegion = coupledRegion[i];
355 
356  if (ownRegion != neiRegion)
357  {
358  addToInterface
359  (
360  mesh,
361  (useFaceZones ? mesh.faceZones().whichZone(facei) : -1),
362  ownRegion,
363  neiRegion,
364  regionsToSize
365  );
366  }
367  }
368 
369 
370  if (Pstream::parRun())
371  {
372  if (Pstream::master())
373  {
374  // Receive and add to my sizes
375  for
376  (
377  int slave=Pstream::firstSlave();
378  slave<=Pstream::lastSlave();
379  slave++
380  )
381  {
382  IPstream fromSlave(Pstream::commsTypes::blocking, slave);
383 
384  EdgeMap<Map<label>> slaveSizes(fromSlave);
385 
386  forAllConstIter(EdgeMap<Map<label>>, slaveSizes, slaveIter)
387  {
388  EdgeMap<Map<label>>::iterator masterIter =
389  regionsToSize.find(slaveIter.key());
390 
391  if (masterIter != regionsToSize.end())
392  {
393  // Same inter-region
394  const Map<label>& slaveInfo = slaveIter();
395  Map<label>& masterInfo = masterIter();
396 
397  forAllConstIter(Map<label>, slaveInfo, iter)
398  {
399  label zoneID = iter.key();
400  label slaveSize = iter();
401 
402  Map<label>::iterator zoneFnd = masterInfo.find
403  (
404  zoneID
405  );
406  if (zoneFnd != masterInfo.end())
407  {
408  zoneFnd() += slaveSize;
409  }
410  else
411  {
412  masterInfo.insert(zoneID, slaveSize);
413  }
414  }
415  }
416  else
417  {
418  regionsToSize.insert(slaveIter.key(), slaveIter());
419  }
420  }
421  }
422  }
423  else
424  {
425  // Send to master
426  {
427  OPstream toMaster
428  (
431  );
432  toMaster << regionsToSize;
433  }
434  }
435  }
436 
437  // Rework
438 
439  Pstream::scatter(regionsToSize);
440 
441 
442 
443  // Now we have the global sizes of all inter-regions.
444  // Invert this on master and distribute.
445  label nInterfaces = 0;
446  forAllConstIter(EdgeMap<Map<label>>, regionsToSize, iter)
447  {
448  const Map<label>& info = iter();
449  nInterfaces += info.size();
450  }
451 
452  interfaces.setSize(nInterfaces);
453  interfaceNames.setSize(nInterfaces);
454  interfaceSizes.setSize(nInterfaces);
455  EdgeMap<Map<label>> regionsToInterface(nInterfaces);
456 
457  nInterfaces = 0;
458  forAllConstIter(EdgeMap<Map<label>>, regionsToSize, iter)
459  {
460  const edge& e = iter.key();
461  const word& name0 = regionNames[e[0]];
462  const word& name1 = regionNames[e[1]];
463 
464  const Map<label>& info = iter();
465  forAllConstIter(Map<label>, info, infoIter)
466  {
467  interfaces[nInterfaces] = iter.key();
468  label zoneID = infoIter.key();
469  if (zoneID == -1)
470  {
471  interfaceNames[nInterfaces] = Pair<word>
472  (
473  name0 + "_to_" + name1,
474  name1 + "_to_" + name0
475  );
476  }
477  else
478  {
479  const word& zoneName = mesh.faceZones()[zoneID].name();
480  interfaceNames[nInterfaces] = Pair<word>
481  (
482  zoneName + "_" + name0 + "_to_" + name1,
483  zoneName + "_" + name1 + "_to_" + name0
484  );
485  }
486  interfaceSizes[nInterfaces] = infoIter();
487 
488  if (regionsToInterface.found(e))
489  {
490  regionsToInterface[e].insert(zoneID, nInterfaces);
491  }
492  else
493  {
494  Map<label> zoneAndInterface;
495  zoneAndInterface.insert(zoneID, nInterfaces);
496  regionsToInterface.insert(e, zoneAndInterface);
497  }
498  nInterfaces++;
499  }
500  }
501 
502 
503  // Now all processor have consistent interface information
504 
505  Pstream::scatter(interfaces);
506  Pstream::scatter(interfaceNames);
507  Pstream::scatter(interfaceSizes);
508  Pstream::scatter(regionsToInterface);
509 
510  // Mark all inter-region faces.
511  faceToInterface.setSize(mesh.nFaces(), -1);
512 
513  forAll(mesh.faceNeighbour(), facei)
514  {
515  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
516  label neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
517 
518  if (ownRegion != neiRegion)
519  {
520  label zoneID = -1;
521  if (useFaceZones)
522  {
523  zoneID = mesh.faceZones().whichZone(facei);
524  }
525 
526  edge interface
527  (
528  min(ownRegion, neiRegion),
529  max(ownRegion, neiRegion)
530  );
531 
532  faceToInterface[facei] = regionsToInterface[interface][zoneID];
533  }
534  }
535  forAll(coupledRegion, i)
536  {
537  label facei = i+mesh.nInternalFaces();
538  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
539  label neiRegion = coupledRegion[i];
540 
541  if (ownRegion != neiRegion)
542  {
543  label zoneID = -1;
544  if (useFaceZones)
545  {
546  zoneID = mesh.faceZones().whichZone(facei);
547  }
548 
549  edge interface
550  (
551  min(ownRegion, neiRegion),
552  max(ownRegion, neiRegion)
553  );
554 
555  faceToInterface[facei] = regionsToInterface[interface][zoneID];
556  }
557  }
558 }
559 
560 
561 // Create mesh for region.
562 autoPtr<polyTopoChangeMap> createRegionMesh
563 (
564  const fvMesh& mesh,
565  // Region info
566  const labelList& cellRegion,
567  const label regioni,
568  const word& regionName,
569  // Interface info
570  const labelList& interfacePatches,
571  const labelList& faceToInterface,
572 
573  autoPtr<fvMesh>& newMesh
574 )
575 {
576  // Create dummy system/fv*
577  {
579  (
580  "fvSchemes",
581  mesh.time().system(),
582  regionName,
583  mesh,
586  false
587  );
588 
589  Info<< "Testing:" << io.objectPath() << endl;
590 
591  if (!io.headerOk())
592  {
593  Info<< "Writing dummy " << regionName/io.name() << endl;
594  dictionary dummyDict;
595  dictionary divDict;
596  dummyDict.add("divSchemes", divDict);
597  dictionary gradDict;
598  dummyDict.add("gradSchemes", gradDict);
599  dictionary laplDict;
600  dummyDict.add("laplacianSchemes", laplDict);
601 
602  IOdictionary(io, dummyDict).regIOobject::write();
603  }
604  }
605  {
607  (
608  "fvSolution",
609  mesh.time().system(),
610  regionName,
611  mesh,
614  false
615  );
616 
617  if (!io.headerOk())
618  {
619  Info<< "Writing dummy " << regionName/io.name() << endl;
620  dictionary dummyDict;
621  IOdictionary(io, dummyDict).regIOobject::write();
622  }
623  }
624 
625 
626  // Neighbour cellRegion.
627  labelList coupledRegion(mesh.nFaces()-mesh.nInternalFaces());
628 
629  forAll(coupledRegion, i)
630  {
631  label celli = mesh.faceOwner()[i+mesh.nInternalFaces()];
632  coupledRegion[i] = cellRegion[celli];
633  }
634  syncTools::swapBoundaryFaceList(mesh, coupledRegion);
635 
636 
637  // Topology change container. Start off from existing mesh.
638  polyTopoChange meshMod(mesh);
639 
640  // Cell remover engine
641  removeCells cellRemover(mesh);
642 
643  // Select all but region cells
644  labelList cellsToRemove(getNonRegionCells(cellRegion, regioni));
645 
646  // Find out which faces will get exposed. Note that this
647  // gets faces in mesh face order. So both regions will get same
648  // face in same order (important!)
649  labelList exposedFaces = cellRemover.getExposedFaces(cellsToRemove);
650 
651  labelList exposedPatchIDs(exposedFaces.size());
652  forAll(exposedFaces, i)
653  {
654  label facei = exposedFaces[i];
655  label interfacei = faceToInterface[facei];
656 
657  label ownRegion = cellRegion[mesh.faceOwner()[facei]];
658  label neiRegion = -1;
659 
660  if (mesh.isInternalFace(facei))
661  {
662  neiRegion = cellRegion[mesh.faceNeighbour()[facei]];
663  }
664  else
665  {
666  neiRegion = coupledRegion[facei-mesh.nInternalFaces()];
667  }
668 
669 
670  // Check which side is being kept - determines which of the two
671  // patches will be used.
672 
673  label otherRegion = -1;
674 
675  if (ownRegion == regioni && neiRegion != regioni)
676  {
677  otherRegion = neiRegion;
678  }
679  else if (ownRegion != regioni && neiRegion == regioni)
680  {
681  otherRegion = ownRegion;
682  }
683  else
684  {
686  << "Exposed face:" << facei
687  << " fc:" << mesh.faceCentres()[facei]
688  << " has owner region " << ownRegion
689  << " and neighbour region " << neiRegion
690  << " when handling region:" << regioni
691  << exit(FatalError);
692  }
693 
694  // Find the patch.
695  if (regioni < otherRegion)
696  {
697  exposedPatchIDs[i] = interfacePatches[interfacei];
698  }
699  else
700  {
701  exposedPatchIDs[i] = interfacePatches[interfacei]+1;
702  }
703  }
704 
705  // Remove faces
706  cellRemover.setRefinement
707  (
708  cellsToRemove,
709  exposedFaces,
710  exposedPatchIDs,
711  meshMod
712  );
713 
714  autoPtr<polyTopoChangeMap> map = meshMod.makeMesh
715  (
716  newMesh,
717  IOobject
718  (
719  regionName,
720  mesh.time().timeName(),
721  mesh.time(),
724  ),
725  mesh
726  );
727 
728  return map;
729 }
730 
731 
732 void createAndWriteRegion
733 (
734  const fvMesh& mesh,
735  const labelList& cellRegion,
736  const wordList& regionNames,
737  const bool prefixRegion,
738  const labelList& faceToInterface,
739  const labelList& interfacePatches,
740  const label regioni,
741  const word& newMeshInstance
742 )
743 {
744  Info<< "Creating mesh for region " << regioni
745  << ' ' << regionNames[regioni] << endl;
746 
747  autoPtr<fvMesh> newMesh;
748  autoPtr<polyTopoChangeMap> map = createRegionMesh
749  (
750  mesh,
751  cellRegion,
752  regioni,
753  regionNames[regioni],
754  interfacePatches,
755  faceToInterface,
756  newMesh
757  );
758 
759 
760  // Make map of all added patches
761  labelHashSet addedPatches(2*interfacePatches.size());
762  forAll(interfacePatches, interfacei)
763  {
764  addedPatches.insert(interfacePatches[interfacei]);
765  addedPatches.insert(interfacePatches[interfacei]+1);
766  }
767 
768 
769  Info<< "Mapping fields" << endl;
770 
771  // Map existing fields
772  newMesh().topoChange(map());
773 
774  // Add subsetted fields
775  subsetVolFields<volScalarField>
776  (
777  mesh,
778  newMesh(),
779  map().cellMap(),
780  map().faceMap(),
781  addedPatches
782  );
783  subsetVolFields<volVectorField>
784  (
785  mesh,
786  newMesh(),
787  map().cellMap(),
788  map().faceMap(),
789  addedPatches
790  );
791  subsetVolFields<volSphericalTensorField>
792  (
793  mesh,
794  newMesh(),
795  map().cellMap(),
796  map().faceMap(),
797  addedPatches
798  );
799  subsetVolFields<volSymmTensorField>
800  (
801  mesh,
802  newMesh(),
803  map().cellMap(),
804  map().faceMap(),
805  addedPatches
806  );
807  subsetVolFields<volTensorField>
808  (
809  mesh,
810  newMesh(),
811  map().cellMap(),
812  map().faceMap(),
813  addedPatches
814  );
815 
816  subsetSurfaceFields<surfaceScalarField>
817  (
818  mesh,
819  newMesh(),
820  map().cellMap(),
821  map().faceMap(),
822  addedPatches
823  );
824  subsetSurfaceFields<surfaceVectorField>
825  (
826  mesh,
827  newMesh(),
828  map().cellMap(),
829  map().faceMap(),
830  addedPatches
831  );
832  subsetSurfaceFields<surfaceSphericalTensorField>
833  (
834  mesh,
835  newMesh(),
836  map().cellMap(),
837  map().faceMap(),
838  addedPatches
839  );
840  subsetSurfaceFields<surfaceSymmTensorField>
841  (
842  mesh,
843  newMesh(),
844  map().cellMap(),
845  map().faceMap(),
846  addedPatches
847  );
848  subsetSurfaceFields<surfaceTensorField>
849  (
850  mesh,
851  newMesh(),
852  map().cellMap(),
853  map().faceMap(),
854  addedPatches
855  );
856 
857 
858  const polyBoundaryMesh& newPatches = newMesh().boundaryMesh();
859  newPatches.checkParallelSync(true);
860 
861  // Delete empty patches
862  // ~~~~~~~~~~~~~~~~~~~~
863 
864  // Create reordering list to move patches-to-be-deleted to end
865  labelList oldToNew(newPatches.size(), -1);
866  DynamicList<label> sharedPatches(newPatches.size());
867  label newI = 0;
868 
869  Info<< "Deleting empty patches" << endl;
870 
871  // Assumes all non-proc boundaries are on all processors!
872  forAll(newPatches, patchi)
873  {
874  const polyPatch& pp = newPatches[patchi];
875 
876  if (!isA<processorPolyPatch>(pp))
877  {
878  if (returnReduce(pp.size(), sumOp<label>()) > 0)
879  {
880  oldToNew[patchi] = newI;
881  if (!addedPatches.found(patchi))
882  {
883  sharedPatches.append(newI);
884  }
885  newI++;
886  }
887  }
888  }
889 
890  // Same for processor patches (but need no reduction)
891  forAll(newPatches, patchi)
892  {
893  const polyPatch& pp = newPatches[patchi];
894 
895  if (isA<processorPolyPatch>(pp) && pp.size())
896  {
897  oldToNew[patchi] = newI++;
898  }
899  }
900 
901  const label nNewPatches = newI;
902 
903  // Move all deleteable patches to the end
904  forAll(oldToNew, patchi)
905  {
906  if (oldToNew[patchi] == -1)
907  {
908  oldToNew[patchi] = newI++;
909  }
910  }
911 
912  // reorderPatches(newMesh(), oldToNew, nNewPatches);
913  fvMeshTools::reorderPatches(newMesh(), oldToNew, nNewPatches, true);
914 
915  // Rename shared patches with region name
916  if (prefixRegion)
917  {
918  Info<< "Prefixing patches with region name" << endl;
919 
920  renamePatches(newMesh(), regionNames[regioni], sharedPatches);
921  }
922 
923 
924  Info<< "Writing new mesh" << endl;
925 
926  newMesh().setInstance(newMeshInstance);
927  newMesh().write();
928 
929  // Write addressing files like decomposePar
930  Info<< "Writing addressing to base mesh" << endl;
931 
932  labelIOList pointProcAddressing
933  (
934  IOobject
935  (
936  "pointRegionAddressing",
937  newMesh().facesInstance(),
938  newMesh().meshSubDir,
939  newMesh(),
940  IOobject::NO_READ,
941  IOobject::NO_WRITE,
942  false
943  ),
944  map().pointMap()
945  );
946  Info<< "Writing map " << pointProcAddressing.name()
947  << " from region" << regioni
948  << " points back to base mesh." << endl;
949  pointProcAddressing.write();
950 
951  labelIOList faceProcAddressing
952  (
953  IOobject
954  (
955  "faceRegionAddressing",
956  newMesh().facesInstance(),
957  newMesh().meshSubDir,
958  newMesh(),
959  IOobject::NO_READ,
960  IOobject::NO_WRITE,
961  false
962  ),
963  newMesh().nFaces()
964  );
965  forAll(faceProcAddressing, facei)
966  {
967  // face + turning index. (see decomposePar)
968  // Is the face pointing in the same direction?
969  label oldFacei = map().faceMap()[facei];
970 
971  if
972  (
973  map().cellMap()[newMesh().faceOwner()[facei]]
974  == mesh.faceOwner()[oldFacei]
975  )
976  {
977  faceProcAddressing[facei] = oldFacei+1;
978  }
979  else
980  {
981  faceProcAddressing[facei] = -(oldFacei+1);
982  }
983  }
984  Info<< "Writing map " << faceProcAddressing.name()
985  << " from region" << regioni
986  << " faces back to base mesh." << endl;
987  faceProcAddressing.write();
988 
989  labelIOList cellProcAddressing
990  (
991  IOobject
992  (
993  "cellRegionAddressing",
994  newMesh().facesInstance(),
995  newMesh().meshSubDir,
996  newMesh(),
997  IOobject::NO_READ,
998  IOobject::NO_WRITE,
999  false
1000  ),
1001  map().cellMap()
1002  );
1003  Info<< "Writing map " <<cellProcAddressing.name()
1004  << " from region" << regioni
1005  << " cells back to base mesh." << endl;
1006  cellProcAddressing.write();
1007 }
1008 
1009 
1010 // Create for every region-region interface with non-zero size two patches.
1011 // First one is for minimumregion to maximumregion.
1012 // Note that patches get created in same order on all processors (if parallel)
1013 // since looping over synchronised 'interfaces'.
1014 labelList addRegionPatches
1015 (
1016  fvMesh& mesh,
1017  const wordList& regionNames,
1018  const edgeList& interfaces,
1019  const List<Pair<word>>& interfaceNames
1020 )
1021 {
1022  Info<< nl << "Adding patches" << nl << endl;
1023 
1024  labelList interfacePatches(interfaces.size());
1025 
1026  forAll(interfaces, interI)
1027  {
1028  const edge& e = interfaces[interI];
1029  const Pair<word>& names = interfaceNames[interI];
1030 
1031  // Info<< "For interface " << interI
1032  // << " between regions " << e
1033  // << " trying to add patches " << names << endl;
1034 
1035 
1036  mappedWallPolyPatch patch1
1037  (
1038  names[0],
1039  0, // overridden
1040  0, // overridden
1041  0, // overridden
1042  regionNames[e[1]], // sampleRegion
1044  names[1], // samplePatch
1045  point::zero, // offset
1046  mesh.boundaryMesh()
1047  );
1048 
1049  interfacePatches[interI] = fvMeshTools::addPatch
1050  (
1051  mesh,
1052  patch1,
1053  dictionary(), // optional per field value
1055  true // validBoundary
1056  );
1057 
1058  mappedWallPolyPatch patch2
1059  (
1060  names[1],
1061  0,
1062  0,
1063  0,
1064  regionNames[e[0]], // sampleRegion
1066  names[0],
1067  point::zero, // offset
1068  mesh.boundaryMesh()
1069  );
1071  (
1072  mesh,
1073  patch2,
1074  dictionary(), // optional per field value
1076  true // validBoundary
1077  );
1078 
1079  Info<< "For interface between region " << regionNames[e[0]]
1080  << " and " << regionNames[e[1]] << " added patches" << endl
1081  << " " << interfacePatches[interI]
1082  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]].name()
1083  << endl
1084  << " " << interfacePatches[interI]+1
1085  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]+1].name()
1086  << endl;
1087  }
1088  return interfacePatches;
1089 }
1090 
1091 
1092 // Find region that covers most of cell zone
1093 label findCorrespondingRegion
1094 (
1095  const labelList& existingZoneID, // per cell the (unique) zoneID
1096  const labelList& cellRegion,
1097  const label nCellRegions,
1098  const label zoneI,
1099  const label minOverlapSize
1100 )
1101 {
1102  // Per region the number of cells in zoneI
1103  labelList cellsInZone(nCellRegions, 0);
1104 
1105  forAll(cellRegion, celli)
1106  {
1107  if (existingZoneID[celli] == zoneI)
1108  {
1109  cellsInZone[cellRegion[celli]]++;
1110  }
1111  }
1112 
1114  Pstream::listCombineScatter(cellsInZone);
1115 
1116  // Pick region with largest overlap of zoneI
1117  label regioni = findMax(cellsInZone);
1118 
1119 
1120  if (cellsInZone[regioni] < minOverlapSize)
1121  {
1122  // Region covers too little of zone. Not good enough.
1123  regioni = -1;
1124  }
1125  else
1126  {
1127  // Check that region contains no cells that aren't in cellZone.
1128  forAll(cellRegion, celli)
1129  {
1130  if (cellRegion[celli] == regioni && existingZoneID[celli] != zoneI)
1131  {
1132  // celli in regioni but not in zoneI
1133  regioni = -1;
1134  break;
1135  }
1136  }
1137  // If one in error, all should be in error. Note that branch gets taken
1138  // on all procs.
1139  reduce(regioni, minOp<label>());
1140  }
1141 
1142  return regioni;
1143 }
1144 
1145 
1146 // Get zone per cell
1147 // - non-unique zoning
1148 // - coupled zones
1149 void getZoneID
1150 (
1151  const polyMesh& mesh,
1152  const meshCellZones& cellZones,
1153  labelList& zoneID,
1154  labelList& neiZoneID
1155 )
1156 {
1157  // Existing zoneID
1158  zoneID.setSize(mesh.nCells());
1159  zoneID = -1;
1160 
1161  forAll(cellZones, zoneI)
1162  {
1163  const cellZone& cz = cellZones[zoneI];
1164 
1165  forAll(cz, i)
1166  {
1167  label celli = cz[i];
1168  if (zoneID[celli] == -1)
1169  {
1170  zoneID[celli] = zoneI;
1171  }
1172  else
1173  {
1175  << "Cell " << celli << " with cell centre "
1176  << mesh.cellCentres()[celli]
1177  << " is multiple zones. This is not allowed." << endl
1178  << "It is in zone " << cellZones[zoneID[celli]].name()
1179  << " and in zone " << cellZones[zoneI].name()
1180  << exit(FatalError);
1181  }
1182  }
1183  }
1184 
1185  // Neighbour zoneID.
1186  neiZoneID.setSize(mesh.nFaces()-mesh.nInternalFaces());
1187 
1188  forAll(neiZoneID, i)
1189  {
1190  neiZoneID[i] = zoneID[mesh.faceOwner()[i+mesh.nInternalFaces()]];
1191  }
1192  syncTools::swapBoundaryFaceList(mesh, neiZoneID);
1193 }
1194 
1195 
1196 void matchRegions
1197 (
1198  const bool sloppyCellZones,
1199  const polyMesh& mesh,
1200 
1201  const label nCellRegions,
1202  const labelList& cellRegion,
1203 
1204  const word& defaultRegionName,
1205 
1206  labelList& regionToZone,
1207  wordList& regionNames,
1208  labelList& zoneToRegion
1209 )
1210 {
1211  const meshCellZones& cellZones = mesh.cellZones();
1212 
1213  regionToZone.setSize(nCellRegions, -1);
1214  regionNames.setSize(nCellRegions);
1215  zoneToRegion.setSize(cellZones.size(), -1);
1216 
1217  // Get current per cell zoneID
1218  labelList zoneID(mesh.nCells(), -1);
1219  labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1220  getZoneID(mesh, cellZones, zoneID, neiZoneID);
1221 
1222  // Sizes per cellzone
1223  labelList zoneSizes(cellZones.size(), 0);
1224  {
1225  List<wordList> zoneNames(Pstream::nProcs());
1226  zoneNames[Pstream::myProcNo()] = cellZones.names();
1227  Pstream::gatherList(zoneNames);
1228  Pstream::scatterList(zoneNames);
1229 
1230  forAll(zoneNames, proci)
1231  {
1232  if (zoneNames[proci] != zoneNames[0])
1233  {
1235  << "cellZones not synchronised across processors." << endl
1236  << "Master has cellZones " << zoneNames[0] << endl
1237  << "Processor " << proci
1238  << " has cellZones " << zoneNames[proci]
1239  << exit(FatalError);
1240  }
1241  }
1242 
1243  forAll(cellZones, zoneI)
1244  {
1245  zoneSizes[zoneI] = returnReduce
1246  (
1247  cellZones[zoneI].size(),
1248  sumOp<label>()
1249  );
1250  }
1251  }
1252 
1253 
1254  if (sloppyCellZones)
1255  {
1256  Info<< "Trying to match regions to existing cell zones;"
1257  << " region can be subset of cell zone." << nl << endl;
1258 
1259  forAll(cellZones, zoneI)
1260  {
1261  label regioni = findCorrespondingRegion
1262  (
1263  zoneID,
1264  cellRegion,
1265  nCellRegions,
1266  zoneI,
1267  label(0.5*zoneSizes[zoneI]) // minimum overlap
1268  );
1269 
1270  if (regioni != -1)
1271  {
1272  Info<< "Sloppily matched region " << regioni
1273  //<< " size " << regionSizes[regioni]
1274  << " to zone " << zoneI << " size " << zoneSizes[zoneI]
1275  << endl;
1276  zoneToRegion[zoneI] = regioni;
1277  regionToZone[regioni] = zoneI;
1278  regionNames[regioni] = cellZones[zoneI].name();
1279  }
1280  }
1281  }
1282  else
1283  {
1284  Info<< "Trying to match regions to existing cell zones." << nl << endl;
1285 
1286  forAll(cellZones, zoneI)
1287  {
1288  label regioni = findCorrespondingRegion
1289  (
1290  zoneID,
1291  cellRegion,
1292  nCellRegions,
1293  zoneI,
1294  1 // minimum overlap
1295  );
1296 
1297  if (regioni != -1)
1298  {
1299  zoneToRegion[zoneI] = regioni;
1300  regionToZone[regioni] = zoneI;
1301  regionNames[regioni] = cellZones[zoneI].name();
1302  }
1303  }
1304  }
1305 
1306  label nUnmatchedRegions = 0;
1307 
1308  forAll(regionToZone, regioni)
1309  {
1310  if (regionToZone[regioni] == -1)
1311  {
1312  nUnmatchedRegions++;
1313  }
1314  }
1315 
1316  if (nUnmatchedRegions)
1317  {
1318  label nUnmatchedi = 1;
1319 
1320  // Allocate region names for unmatched regions.
1321  forAll(regionToZone, regioni)
1322  {
1323  if (regionToZone[regioni] == -1)
1324  {
1325  if
1326  (
1327  nUnmatchedRegions == 1
1328  && defaultRegionName != standardRegionName
1329  )
1330  {
1331  regionNames[regioni] = defaultRegionName;
1332  }
1333  else
1334  {
1335  regionNames[regioni] =
1336  defaultRegionName + Foam::name(nUnmatchedi++);
1337  }
1338  }
1339  }
1340  }
1341 }
1342 
1343 
1344 void writeCellToRegion(const fvMesh& mesh, const labelList& cellRegion)
1345 {
1346  // Write to manual decomposition option
1347  {
1348  labelIOList cellToRegion
1349  (
1350  IOobject
1351  (
1352  "cellToRegion",
1353  mesh.facesInstance(),
1354  mesh,
1357  false
1358  ),
1359  cellRegion
1360  );
1361  cellToRegion.write();
1362 
1363  Info<< "Writing region per cell file (for manual decomposition) to "
1364  << cellToRegion.relativeObjectPath() << nl << endl;
1365  }
1366  // Write for postprocessing
1367  {
1368  volScalarField cellToRegion
1369  (
1370  IOobject
1371  (
1372  "cellToRegion",
1373  mesh.time().timeName(),
1374  mesh,
1377  false
1378  ),
1379  mesh,
1381  zeroGradientFvPatchScalarField::typeName
1382  );
1383  forAll(cellRegion, celli)
1384  {
1385  cellToRegion[celli] = cellRegion[celli];
1386  }
1387  cellToRegion.write();
1388 
1389  Info<< "Writing region per cell as volScalarField to "
1390  << cellToRegion.relativeObjectPath() << nl << endl;
1391  }
1392 }
1393 
1394 
1395 
1396 
1397 int main(int argc, char *argv[])
1398 {
1400  (
1401  "splits mesh into multiple regions (detected by walking across faces)"
1402  );
1403  #include "addRegionOption.H"
1404  #include "addOverwriteOption.H"
1406  (
1407  "cellZones",
1408  "additionally split cellZones off into separate regions"
1409  );
1411  (
1412  "cellZonesOnly",
1413  "use cellZones only to split mesh into regions; do not use walking"
1414  );
1416  (
1417  "cellZonesFileOnly",
1418  "file",
1419  "like -cellZonesOnly, but use specified file"
1420  );
1422  (
1423  "blockedFaces",
1424  "faceSet",
1425  "specify additional region boundaries that walking does not cross"
1426  );
1428  (
1429  "makeCellZones",
1430  "place cells into cellZones instead of splitting mesh"
1431  );
1433  (
1434  "largestOnly",
1435  "only write largest region"
1436  );
1438  (
1439  "insidePoint",
1440  "point",
1441  "only write region containing point"
1442  );
1444  (
1445  "detectOnly",
1446  "do not write mesh"
1447  );
1449  (
1450  "sloppyCellZones",
1451  "try to match heuristically regions to existing cell zones"
1452  );
1454  (
1455  "useFaceZones",
1456  "use faceZones to patch inter-region faces instead of single patch"
1457  );
1459  (
1460  "prefixRegion",
1461  "prefix region name to all patches, not just coupling patches"
1462  );
1464  (
1465  "defaultRegionName",
1466  "name",
1467  "base name of the unspecified regions, defaults to \"region\""
1468  );
1470  (
1471  "noFields",
1472  "do not update fields"
1473  );
1474 
1475  #include "setRootCase.H"
1476  #include "createTime.H"
1477 
1478  runTime.functionObjects().off();
1479 
1480  #include "createNamedMesh.H"
1481 
1482  const word oldInstance = mesh.pointsInstance();
1483 
1484  word blockedFacesName;
1485  if (args.optionReadIfPresent("blockedFaces", blockedFacesName))
1486  {
1487  Info<< "Reading blocked internal faces from faceSet "
1488  << blockedFacesName << nl << endl;
1489  }
1490 
1491  const bool makeCellZones = args.optionFound("makeCellZones");
1492  const bool largestOnly = args.optionFound("largestOnly");
1493  const bool insidePoint = args.optionFound("insidePoint");
1494  const bool useCellZones = args.optionFound("cellZones");
1495  const bool useCellZonesOnly = args.optionFound("cellZonesOnly");
1496  const bool useCellZonesFile = args.optionFound("cellZonesFileOnly");
1497  const bool overwrite = args.optionFound("overwrite");
1498  const bool detectOnly = args.optionFound("detectOnly");
1499  const bool sloppyCellZones = args.optionFound("sloppyCellZones");
1500  const bool useFaceZones = args.optionFound("useFaceZones");
1501  const bool prefixRegion = args.optionFound("prefixRegion");
1502  const bool fields = !args.optionFound("noFields");
1503 
1504  if
1505  (
1506  (useCellZonesOnly || useCellZonesFile)
1507  && (useCellZones || blockedFacesName.size())
1508  )
1509  {
1511  << "You cannot specify both -cellZonesOnly or -cellZonesFileOnly"
1512  << " (which specify complete split)"
1513  << " in combination with -blockedFaces or -cellZones"
1514  << " (which imply a split based on topology)"
1515  << exit(FatalError);
1516  }
1517 
1518 
1519  if (useFaceZones)
1520  {
1521  Info<< "Using current faceZones to divide inter-region interfaces"
1522  << " into multiple patches."
1523  << nl << endl;
1524  }
1525  else
1526  {
1527  Info<< "Creating single patch per inter-region interface."
1528  << nl << endl;
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  const word defaultRegionName
1542  (
1543  args.optionLookupOrDefault("defaultRegionName", standardRegionName)
1544  );
1545 
1546  const meshCellZones& cellZones = mesh.cellZones();
1547 
1548  // Existing zoneID
1549  labelList zoneID(mesh.nCells(), -1);
1550  // Neighbour zoneID.
1551  labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1552  getZoneID(mesh, cellZones, zoneID, neiZoneID);
1553 
1554 
1555 
1556  // Determine per cell the region it belongs to
1557  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1558 
1559  // cellRegion is the labelList with the region per cell.
1560  labelList cellRegion;
1561  // Region per zone
1562  labelList regionToZone;
1563  // Name of region
1564  wordList regionNames;
1565  // Zone to region
1566  labelList zoneToRegion;
1567 
1568  label nCellRegions = 0;
1569  if (useCellZonesOnly)
1570  {
1571  Info<< "Using current cellZones to split mesh into regions."
1572  << " This requires all"
1573  << " cells to be in one and only one cellZone." << nl << endl;
1574 
1575  label unzonedCelli = findIndex(zoneID, -1);
1576  if (unzonedCelli != -1)
1577  {
1579  << "For the cellZonesOnly option all cells "
1580  << "have to be in a cellZone." << endl
1581  << "Cell " << unzonedCelli
1582  << " at" << mesh.cellCentres()[unzonedCelli]
1583  << " is not in a cellZone. There might be more unzoned cells."
1584  << exit(FatalError);
1585  }
1586  cellRegion = zoneID;
1587  nCellRegions = gMax(cellRegion)+1;
1588  regionToZone.setSize(nCellRegions);
1589  regionNames.setSize(nCellRegions);
1590  zoneToRegion.setSize(cellZones.size(), -1);
1591  for (label regioni = 0; regioni < nCellRegions; regioni++)
1592  {
1593  regionToZone[regioni] = regioni;
1594  zoneToRegion[regioni] = regioni;
1595  regionNames[regioni] = cellZones[regioni].name();
1596  }
1597  }
1598  else if (useCellZonesFile)
1599  {
1600  const word zoneFile = args.option("cellZonesFileOnly");
1601  Info<< "Reading split from cellZones file " << zoneFile << endl
1602  << "This requires all"
1603  << " cells to be in one and only one cellZone." << nl << endl;
1604 
1605  meshCellZones newCellZones
1606  (
1607  IOobject
1608  (
1609  zoneFile,
1610  mesh.facesInstance(),
1612  mesh,
1615  false
1616  ),
1617  mesh
1618  );
1619 
1620  labelList newZoneID(mesh.nCells(), -1);
1621  labelList newNeiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1622  getZoneID(mesh, newCellZones, newZoneID, newNeiZoneID);
1623 
1624  label unzonedCelli = findIndex(newZoneID, -1);
1625  if (unzonedCelli != -1)
1626  {
1628  << "For the cellZonesFileOnly option all cells "
1629  << "have to be in a cellZone." << endl
1630  << "Cell " << unzonedCelli
1631  << " at" << mesh.cellCentres()[unzonedCelli]
1632  << " is not in a cellZone. There might be more unzoned cells."
1633  << exit(FatalError);
1634  }
1635  cellRegion = newZoneID;
1636  nCellRegions = gMax(cellRegion)+1;
1637  zoneToRegion.setSize(newCellZones.size(), -1);
1638  regionToZone.setSize(nCellRegions);
1639  regionNames.setSize(nCellRegions);
1640  for (label regioni = 0; regioni < nCellRegions; regioni++)
1641  {
1642  regionToZone[regioni] = regioni;
1643  zoneToRegion[regioni] = regioni;
1644  regionNames[regioni] = newCellZones[regioni].name();
1645  }
1646  }
1647  else
1648  {
1649  // Determine connected regions
1650  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1651 
1652  // Mark additional faces that are blocked
1653  boolList blockedFace;
1654 
1655  // Read from faceSet
1656  if (blockedFacesName.size())
1657  {
1658  faceSet blockedFaceSet(mesh, blockedFacesName);
1659  Info<< "Read "
1660  << returnReduce(blockedFaceSet.size(), sumOp<label>())
1661  << " blocked faces from set " << blockedFacesName << nl << endl;
1662 
1663  blockedFace.setSize(mesh.nFaces(), false);
1664 
1665  forAllConstIter(faceSet, blockedFaceSet, iter)
1666  {
1667  blockedFace[iter.key()] = true;
1668  }
1669  }
1670 
1671  // Imply from differing cellZones
1672  if (useCellZones)
1673  {
1674  blockedFace.setSize(mesh.nFaces(), false);
1675 
1676  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1677  {
1678  label own = mesh.faceOwner()[facei];
1679  label nei = mesh.faceNeighbour()[facei];
1680 
1681  if (zoneID[own] != zoneID[nei])
1682  {
1683  blockedFace[facei] = true;
1684  }
1685  }
1686 
1687  // Different cellZones on either side of processor patch.
1688  forAll(neiZoneID, i)
1689  {
1690  label facei = i+mesh.nInternalFaces();
1691 
1692  if (zoneID[mesh.faceOwner()[facei]] != neiZoneID[i])
1693  {
1694  blockedFace[facei] = true;
1695  }
1696  }
1697  }
1698 
1699  // Do a topological walk to determine regions
1700  regionSplit regions(mesh, blockedFace);
1701  nCellRegions = regions.nRegions();
1702  cellRegion.transfer(regions);
1703 
1704  // Make up region names. If possible match them to existing zones.
1705  matchRegions
1706  (
1707  sloppyCellZones,
1708  mesh,
1709  nCellRegions,
1710  cellRegion,
1711  defaultRegionName,
1712 
1713  regionToZone,
1714  regionNames,
1715  zoneToRegion
1716  );
1717 
1718  // Override any default region names if single region selected
1719  if (largestOnly || insidePoint)
1720  {
1721  forAll(regionToZone, regioni)
1722  {
1723  if (regionToZone[regioni] == -1)
1724  {
1725  if (overwrite)
1726  {
1727  regionNames[regioni] = polyMesh::defaultRegion;
1728  }
1729  else if (insidePoint)
1730  {
1731  regionNames[regioni] = "insidePoint";
1732  }
1733  else if (largestOnly)
1734  {
1735  regionNames[regioni] = "largestOnly";
1736  }
1737  }
1738  }
1739  }
1740  }
1741 
1742  Info<< endl << "Number of regions:" << nCellRegions << nl << endl;
1743 
1744 
1745  // Write decomposition to file
1746  writeCellToRegion(mesh, cellRegion);
1747 
1748 
1749 
1750  // Sizes per region
1751  // ~~~~~~~~~~~~~~~~
1752 
1753  labelList regionSizes(nCellRegions, 0);
1754 
1755  forAll(cellRegion, celli)
1756  {
1757  regionSizes[cellRegion[celli]]++;
1758  }
1759  forAll(regionSizes, regioni)
1760  {
1761  reduce(regionSizes[regioni], sumOp<label>());
1762  }
1763 
1764  Info<< "Region\tCells" << nl
1765  << "------\t-----" << endl;
1766 
1767  forAll(regionSizes, regioni)
1768  {
1769  Info<< regioni << '\t' << regionSizes[regioni] << nl;
1770  }
1771  Info<< endl;
1772 
1773 
1774 
1775  // Print region to zone
1776  Info<< "Region\tZone\tName" << nl
1777  << "------\t----\t----" << endl;
1778  forAll(regionToZone, regioni)
1779  {
1780  Info<< regioni << '\t' << regionToZone[regioni] << '\t'
1781  << regionNames[regioni] << nl;
1782  }
1783  Info<< endl;
1784 
1785 
1786 
1787  // Since we're going to mess with patches and zones make sure all
1788  // is synchronised
1789  mesh.boundaryMesh().checkParallelSync(true);
1790  mesh.faceZones().checkParallelSync(true);
1791 
1792 
1793  // Interfaces
1794  // ----------
1795  // per interface:
1796  // - the two regions on either side
1797  // - the name
1798  // - the (global) size
1799  edgeList interfaces;
1800  List<Pair<word>> interfaceNames;
1801  labelList interfaceSizes;
1802  // per face the interface
1803  labelList faceToInterface;
1804 
1805  getInterfaceSizes
1806  (
1807  mesh,
1808  useFaceZones,
1809  cellRegion,
1810  regionNames,
1811 
1812  interfaces,
1813  interfaceNames,
1814  interfaceSizes,
1815  faceToInterface
1816  );
1817 
1818  Info<< "Sizes of interfaces between regions:" << nl << nl
1819  << "Interface\tRegion\tRegion\tFaces" << nl
1820  << "---------\t------\t------\t-----" << endl;
1821 
1822  forAll(interfaces, interI)
1823  {
1824  const edge& e = interfaces[interI];
1825 
1826  Info<< interI
1827  << "\t\t" << e[0] << '\t' << e[1]
1828  << '\t' << interfaceSizes[interI] << nl;
1829  }
1830  Info<< endl;
1831 
1832 
1833  if (detectOnly)
1834  {
1835  return 0;
1836  }
1837 
1838 
1839  // Read objects in time directory
1840  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1841 
1842  IOobjectList objects(mesh, runTime.timeName());
1843 
1844  if (fields) Info<< "Reading geometric fields" << nl << endl;
1845 
1846  #include "readVolFields.H"
1847  #include "readSurfaceFields.H"
1848  // #include "readPointFields.H"
1849 
1850  Info<< endl;
1851 
1852 
1853  // Remove any demand-driven fields ('S', 'V' etc)
1854  mesh.clearOut();
1855 
1856 
1857  if (nCellRegions == 1)
1858  {
1859  Info<< "Only one region. Doing nothing." << endl;
1860  }
1861  else if (makeCellZones)
1862  {
1863  Info<< "Putting cells into cellZones instead of splitting mesh."
1864  << endl;
1865 
1866  // Check if region overlaps with existing zone. If so keep.
1867 
1868  for (label regioni = 0; regioni < nCellRegions; regioni++)
1869  {
1870  label zoneI = regionToZone[regioni];
1871 
1872  if (zoneI != -1)
1873  {
1874  Info<< " Region " << regioni << " : corresponds to existing"
1875  << " cellZone "
1876  << zoneI << ' ' << cellZones[zoneI].name() << endl;
1877  }
1878  else
1879  {
1880  // Create new cellZone.
1881  labelList regionCells = findIndices(cellRegion, regioni);
1882 
1883  word zoneName = "region" + Foam::name(regioni);
1884 
1885  zoneI = cellZones.findZoneID(zoneName);
1886 
1887  if (zoneI == -1)
1888  {
1889  zoneI = cellZones.size();
1890  mesh.cellZones().setSize(zoneI+1);
1891  mesh.cellZones().set
1892  (
1893  zoneI,
1894  new cellZone
1895  (
1896  zoneName, // name
1897  regionCells, // addressing
1898  zoneI, // index
1899  cellZones // meshCellZones
1900  )
1901  );
1902  }
1903  else
1904  {
1905  mesh.cellZones()[zoneI].clearAddressing();
1906  mesh.cellZones()[zoneI] = regionCells;
1907  }
1908  Info<< " Region " << regioni << " : created new cellZone "
1909  << zoneI << ' ' << cellZones[zoneI].name() << endl;
1910  }
1911  }
1913 
1914  if (!overwrite)
1915  {
1916  runTime++;
1917  mesh.setInstance(runTime.timeName());
1918  }
1919  else
1920  {
1921  mesh.setInstance(oldInstance);
1922  }
1923 
1924  Info<< "Writing cellZones as new mesh to time " << runTime.timeName()
1925  << nl << endl;
1926 
1927  mesh.write();
1928 
1929 
1930  // Write cellSets for convenience
1931  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1932 
1933  Info<< "Writing cellSets corresponding to cellZones." << nl << endl;
1934 
1935  forAll(cellZones, zoneI)
1936  {
1937  const cellZone& cz = cellZones[zoneI];
1938 
1939  cellSet(mesh, cz.name(), cz).write();
1940  }
1941  }
1942  else
1943  {
1944  // Add patches for interfaces
1945  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
1946 
1947  // Add all possible patches. Empty ones get filtered later on.
1948  Info<< nl << "Adding patches" << nl << endl;
1949 
1950  labelList interfacePatches
1951  (
1952  addRegionPatches
1953  (
1954  mesh,
1955  regionNames,
1956  interfaces,
1957  interfaceNames
1958  )
1959  );
1960 
1961 
1962  if (!overwrite)
1963  {
1964  runTime++;
1965  }
1966 
1967 
1968  // Create regions
1969  // ~~~~~~~~~~~~~~
1970 
1971  if (insidePoint)
1972  {
1973  const point insidePoint = args.optionRead<point>("insidePoint");
1974 
1975  label regioni = -1;
1976 
1977  (void)mesh.tetBasePtIs();
1978 
1979  label celli = mesh.findCell(insidePoint);
1980 
1981  Info<< nl << "Found point " << insidePoint << " in cell " << celli
1982  << endl;
1983 
1984  if (celli != -1)
1985  {
1986  regioni = cellRegion[celli];
1987  }
1988 
1989  reduce(regioni, maxOp<label>());
1990 
1991  Info<< nl
1992  << "Subsetting region " << regioni
1993  << " containing point " << insidePoint << endl;
1994 
1995  if (regioni == -1)
1996  {
1998  << "Point " << insidePoint
1999  << " is not inside the mesh." << nl
2000  << "Bounding box of the mesh:" << mesh.bounds()
2001  << exit(FatalError);
2002  }
2003 
2004  createAndWriteRegion
2005  (
2006  mesh,
2007  cellRegion,
2008  regionNames,
2009  prefixRegion,
2010  faceToInterface,
2011  interfacePatches,
2012  regioni,
2013  (overwrite ? oldInstance : runTime.timeName())
2014  );
2015  }
2016  else if (largestOnly)
2017  {
2018  label regioni = findMax(regionSizes);
2019 
2020  Info<< nl
2021  << "Subsetting region " << regioni
2022  << " of size " << regionSizes[regioni]
2023  << " as named region " << regionNames[regioni] << endl;
2024 
2025  createAndWriteRegion
2026  (
2027  mesh,
2028  cellRegion,
2029  regionNames,
2030  prefixRegion,
2031  faceToInterface,
2032  interfacePatches,
2033  regioni,
2034  (overwrite ? oldInstance : runTime.timeName())
2035  );
2036  }
2037  else
2038  {
2039  // Split all
2040  for (label regioni = 0; regioni < nCellRegions; regioni++)
2041  {
2042  Info<< nl
2043  << "Region " << regioni << nl
2044  << "-------- " << endl;
2045 
2046  createAndWriteRegion
2047  (
2048  mesh,
2049  cellRegion,
2050  regionNames,
2051  prefixRegion,
2052  faceToInterface,
2053  interfacePatches,
2054  regioni,
2055  (overwrite ? oldInstance : runTime.timeName())
2056  );
2057  }
2058  }
2059  }
2060 
2061  Info<< "End\n" << endl;
2062 
2063  return 0;
2064 }
2065 
2066 
2067 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
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:258
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:141
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
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual Ostream & write(const char)=0
Write character.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
static int masterNo()
Process index of the master.
Definition: UPstream.H:417
const word & name() const
Return name.
Definition: IOobject.H:315
Foam::word regionName
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:888
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:501
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1243
void renamePatches(const wordUList &newNames, const bool validBoundary)
Rename the patches. If validBoundary is set this calls topoChange()
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:928
static int firstSlave()
Process index of first slave.
Definition: UPstream.H:446
label nFaces() const
Templated form of IOobject providing type information for file reading and header type checking...
Definition: IOobject.H:537
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
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
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:325
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.
const dimensionSet dimless
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
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:1658
fvMesh & mesh
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:1153
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:882
T optionRead(const word &opt) const
Read a value from the named option.
Definition: argListI.H:193
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
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
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
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:34
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:1237
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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
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,.
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:77
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:459
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:375
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
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, p_rgh, pimple.dict(), thermo.incompressible());mesh.schemes().setFluxRequired(p_rgh.name());hydrostaticInitialisation(p_rgh, p, rho, U, gh, ghf, pRef, thermo, pimple.dict());Info<< "Creating field dpdt\"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));dimensionedScalar initialMass=fvc::domainIntegrate(rho);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:131
label patchi
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:495
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
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:1647
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:76
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:436
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:98
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