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-2023 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 {
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(identityMap(mesh.boundaryMesh().size()));
151 
153  (
154  mesh.objectRegistry::lookupClass<GeoField>()
155  );
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(identityMap(mesh.boundaryMesh().size()));
205 
207  (
208  mesh.objectRegistry::lookupClass<GeoField>()
209  );
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().name(),
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 delete-able 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(),
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(),
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(),
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  mappedWallPolyPatch patch1
1036  (
1037  names[0],
1038  0, // overridden
1039  0, // overridden
1040  0, // overridden
1041  regionNames[e[1]], // neighbourRegion
1042  names[1], // neighbourPatch
1043  mesh.boundaryMesh()
1044  );
1045 
1046  interfacePatches[interI] = fvMeshTools::addPatch
1047  (
1048  mesh,
1049  patch1,
1050  dictionary(), // optional per field value
1052  true // validBoundary
1053  );
1054 
1055  mappedWallPolyPatch patch2
1056  (
1057  names[1],
1058  0,
1059  0,
1060  0,
1061  regionNames[e[0]], // neighbourRegion
1062  names[0],
1063  mesh.boundaryMesh()
1064  );
1065 
1067  (
1068  mesh,
1069  patch2,
1070  dictionary(), // optional per field value
1072  true // validBoundary
1073  );
1074 
1075  Info<< "For interface between region " << regionNames[e[0]]
1076  << " and " << regionNames[e[1]] << " added patches" << endl
1077  << " " << interfacePatches[interI]
1078  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]].name()
1079  << endl
1080  << " " << interfacePatches[interI]+1
1081  << "\t" << mesh.boundaryMesh()[interfacePatches[interI]+1].name()
1082  << endl;
1083  }
1084  return interfacePatches;
1085 }
1086 
1087 
1088 // Find region that covers most of cell zone
1089 label findCorrespondingRegion
1090 (
1091  const labelList& existingZoneID, // per cell the (unique) zoneID
1092  const labelList& cellRegion,
1093  const label nCellRegions,
1094  const label zoneI,
1095  const label minOverlapSize
1096 )
1097 {
1098  // Per region the number of cells in zoneI
1099  labelList cellsInZone(nCellRegions, 0);
1100 
1101  forAll(cellRegion, celli)
1102  {
1103  if (existingZoneID[celli] == zoneI)
1104  {
1105  cellsInZone[cellRegion[celli]]++;
1106  }
1107  }
1108 
1110  Pstream::listCombineScatter(cellsInZone);
1111 
1112  // Pick region with largest overlap of zoneI
1113  label regioni = findMax(cellsInZone);
1114 
1115 
1116  if (cellsInZone[regioni] < minOverlapSize)
1117  {
1118  // Region covers too little of zone. Not good enough.
1119  regioni = -1;
1120  }
1121  else
1122  {
1123  // Check that region contains no cells that aren't in cellZone.
1124  forAll(cellRegion, celli)
1125  {
1126  if (cellRegion[celli] == regioni && existingZoneID[celli] != zoneI)
1127  {
1128  // celli in regioni but not in zoneI
1129  regioni = -1;
1130  break;
1131  }
1132  }
1133  // If one in error, all should be in error. Note that branch gets taken
1134  // on all procs.
1135  reduce(regioni, minOp<label>());
1136  }
1137 
1138  return regioni;
1139 }
1140 
1141 
1142 // Get zone per cell
1143 // - non-unique zoning
1144 // - coupled zones
1145 void getZoneID
1146 (
1147  const polyMesh& mesh,
1148  const meshCellZones& cellZones,
1149  labelList& zoneID,
1150  labelList& neiZoneID
1151 )
1152 {
1153  // Existing zoneID
1154  zoneID.setSize(mesh.nCells());
1155  zoneID = -1;
1156 
1157  forAll(cellZones, zoneI)
1158  {
1159  const cellZone& cz = cellZones[zoneI];
1160 
1161  forAll(cz, i)
1162  {
1163  label celli = cz[i];
1164  if (zoneID[celli] == -1)
1165  {
1166  zoneID[celli] = zoneI;
1167  }
1168  else
1169  {
1171  << "Cell " << celli << " with cell centre "
1172  << mesh.cellCentres()[celli]
1173  << " is multiple zones. This is not allowed." << endl
1174  << "It is in zone " << cellZones[zoneID[celli]].name()
1175  << " and in zone " << cellZones[zoneI].name()
1176  << exit(FatalError);
1177  }
1178  }
1179  }
1180 
1181  // Neighbour zoneID.
1182  neiZoneID.setSize(mesh.nFaces()-mesh.nInternalFaces());
1183 
1184  forAll(neiZoneID, i)
1185  {
1186  neiZoneID[i] = zoneID[mesh.faceOwner()[i+mesh.nInternalFaces()]];
1187  }
1188  syncTools::swapBoundaryFaceList(mesh, neiZoneID);
1189 }
1190 
1191 
1192 void matchRegions
1193 (
1194  const bool sloppyCellZones,
1195  const polyMesh& mesh,
1196 
1197  const label nCellRegions,
1198  const labelList& cellRegion,
1199 
1200  const word& defaultRegionName,
1201 
1202  labelList& regionToZone,
1204  labelList& zoneToRegion
1205 )
1206 {
1207  const meshCellZones& cellZones = mesh.cellZones();
1208 
1209  regionToZone.setSize(nCellRegions, -1);
1210  regionNames.setSize(nCellRegions);
1211  zoneToRegion.setSize(cellZones.size(), -1);
1212 
1213  // Get current per cell zoneID
1214  labelList zoneID(mesh.nCells(), -1);
1215  labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1216  getZoneID(mesh, cellZones, zoneID, neiZoneID);
1217 
1218  // Sizes per cellzone
1219  labelList zoneSizes(cellZones.size(), 0);
1220  {
1221  List<wordList> zoneNames(Pstream::nProcs());
1222  zoneNames[Pstream::myProcNo()] = cellZones.names();
1223  Pstream::gatherList(zoneNames);
1224  Pstream::scatterList(zoneNames);
1225 
1226  forAll(zoneNames, proci)
1227  {
1228  if (zoneNames[proci] != zoneNames[0])
1229  {
1231  << "cellZones not synchronised across processors." << endl
1232  << "Master has cellZones " << zoneNames[0] << endl
1233  << "Processor " << proci
1234  << " has cellZones " << zoneNames[proci]
1235  << exit(FatalError);
1236  }
1237  }
1238 
1239  forAll(cellZones, zoneI)
1240  {
1241  zoneSizes[zoneI] = returnReduce
1242  (
1243  cellZones[zoneI].size(),
1244  sumOp<label>()
1245  );
1246  }
1247  }
1248 
1249 
1250  if (sloppyCellZones)
1251  {
1252  Info<< "Trying to match regions to existing cell zones;"
1253  << " region can be subset of cell zone." << nl << endl;
1254 
1255  forAll(cellZones, zoneI)
1256  {
1257  label regioni = findCorrespondingRegion
1258  (
1259  zoneID,
1260  cellRegion,
1261  nCellRegions,
1262  zoneI,
1263  label(0.5*zoneSizes[zoneI]) // minimum overlap
1264  );
1265 
1266  if (regioni != -1)
1267  {
1268  Info<< "Sloppily matched region " << regioni
1269  //<< " size " << regionSizes[regioni]
1270  << " to zone " << zoneI << " size " << zoneSizes[zoneI]
1271  << endl;
1272  zoneToRegion[zoneI] = regioni;
1273  regionToZone[regioni] = zoneI;
1274  regionNames[regioni] = cellZones[zoneI].name();
1275  }
1276  }
1277  }
1278  else
1279  {
1280  Info<< "Trying to match regions to existing cell zones." << nl << endl;
1281 
1282  forAll(cellZones, zoneI)
1283  {
1284  label regioni = findCorrespondingRegion
1285  (
1286  zoneID,
1287  cellRegion,
1288  nCellRegions,
1289  zoneI,
1290  1 // minimum overlap
1291  );
1292 
1293  if (regioni != -1)
1294  {
1295  zoneToRegion[zoneI] = regioni;
1296  regionToZone[regioni] = zoneI;
1297  regionNames[regioni] = cellZones[zoneI].name();
1298  }
1299  }
1300  }
1301 
1302  label nUnmatchedRegions = 0;
1303 
1304  forAll(regionToZone, regioni)
1305  {
1306  if (regionToZone[regioni] == -1)
1307  {
1308  nUnmatchedRegions++;
1309  }
1310  }
1311 
1312  if (nUnmatchedRegions)
1313  {
1314  label nUnmatchedi = 1;
1315 
1316  // Allocate region names for unmatched regions.
1317  forAll(regionToZone, regioni)
1318  {
1319  if (regionToZone[regioni] == -1)
1320  {
1321  if
1322  (
1323  nUnmatchedRegions == 1
1324  && defaultRegionName != standardRegionName
1325  )
1326  {
1327  regionNames[regioni] = defaultRegionName;
1328  }
1329  else
1330  {
1331  regionNames[regioni] =
1332  defaultRegionName + Foam::name(nUnmatchedi++);
1333  }
1334  }
1335  }
1336  }
1337 }
1338 
1339 
1340 void writeCellToRegion(const fvMesh& mesh, const labelList& cellRegion)
1341 {
1342  // Write to manual decomposition option
1343  {
1344  labelIOList cellToRegion
1345  (
1346  IOobject
1347  (
1348  "cellToRegion",
1349  mesh.facesInstance(),
1350  mesh,
1353  false
1354  ),
1355  cellRegion
1356  );
1357  cellToRegion.write();
1358 
1359  Info<< "Writing region per cell file (for manual decomposition) to "
1360  << cellToRegion.relativeObjectPath() << nl << endl;
1361  }
1362  // Write for postprocessing
1363  {
1364  volScalarField cellToRegion
1365  (
1366  IOobject
1367  (
1368  "cellToRegion",
1369  mesh.time().name(),
1370  mesh,
1373  false
1374  ),
1375  mesh,
1377  zeroGradientFvPatchScalarField::typeName
1378  );
1379  forAll(cellRegion, celli)
1380  {
1381  cellToRegion[celli] = cellRegion[celli];
1382  }
1383  cellToRegion.write();
1384 
1385  Info<< "Writing region per cell as volScalarField to "
1386  << cellToRegion.relativeObjectPath() << nl << endl;
1387  }
1388 }
1389 
1390 
1391 
1392 
1393 int main(int argc, char *argv[])
1394 {
1396  (
1397  "splits mesh into multiple regions (detected by walking across faces)"
1398  );
1399  #include "addRegionOption.H"
1400  #include "addOverwriteOption.H"
1402  (
1403  "cellZones",
1404  "additionally split cellZones off into separate regions"
1405  );
1407  (
1408  "cellZonesOnly",
1409  "use cellZones only to split mesh into regions; do not use walking"
1410  );
1412  (
1413  "cellZonesFileOnly",
1414  "file",
1415  "like -cellZonesOnly, but use specified file"
1416  );
1418  (
1419  "blockedFaces",
1420  "faceSet",
1421  "specify additional region boundaries that walking does not cross"
1422  );
1424  (
1425  "makeCellZones",
1426  "place cells into cellZones instead of splitting mesh"
1427  );
1429  (
1430  "largestOnly",
1431  "only write largest region"
1432  );
1434  (
1435  "insidePoint",
1436  "point",
1437  "only write region containing point"
1438  );
1440  (
1441  "detectOnly",
1442  "do not write mesh"
1443  );
1445  (
1446  "sloppyCellZones",
1447  "try to match heuristically regions to existing cell zones"
1448  );
1450  (
1451  "useFaceZones",
1452  "use faceZones to patch inter-region faces instead of single patch"
1453  );
1455  (
1456  "prefixRegion",
1457  "prefix region name to all patches, not just coupling patches"
1458  );
1460  (
1461  "defaultRegionName",
1462  "name",
1463  "base name of the unspecified regions, defaults to \"region\""
1464  );
1466  (
1467  "noFields",
1468  "do not update fields"
1469  );
1470 
1471  #include "setRootCase.H"
1473  #include "createNamedMesh.H"
1474 
1475  const word oldInstance = mesh.pointsInstance();
1476 
1477  word blockedFacesName;
1478  if (args.optionReadIfPresent("blockedFaces", blockedFacesName))
1479  {
1480  Info<< "Reading blocked internal faces from faceSet "
1481  << blockedFacesName << nl << endl;
1482  }
1483 
1484  const bool makeCellZones = args.optionFound("makeCellZones");
1485  const bool largestOnly = args.optionFound("largestOnly");
1486  const bool insidePoint = args.optionFound("insidePoint");
1487  const bool useCellZones = args.optionFound("cellZones");
1488  const bool useCellZonesOnly = args.optionFound("cellZonesOnly");
1489  const bool useCellZonesFile = args.optionFound("cellZonesFileOnly");
1490  const bool overwrite = args.optionFound("overwrite");
1491  const bool detectOnly = args.optionFound("detectOnly");
1492  const bool sloppyCellZones = args.optionFound("sloppyCellZones");
1493  const bool useFaceZones = args.optionFound("useFaceZones");
1494  const bool prefixRegion = args.optionFound("prefixRegion");
1495  const bool fields = !args.optionFound("noFields");
1496 
1497  if
1498  (
1499  (useCellZonesOnly || useCellZonesFile)
1500  && (useCellZones || blockedFacesName.size())
1501  )
1502  {
1504  << "You cannot specify both -cellZonesOnly or -cellZonesFileOnly"
1505  << " (which specify complete split)"
1506  << " in combination with -blockedFaces or -cellZones"
1507  << " (which imply a split based on topology)"
1508  << exit(FatalError);
1509  }
1510 
1511 
1512  if (useFaceZones)
1513  {
1514  Info<< "Using current faceZones to divide inter-region interfaces"
1515  << " into multiple patches."
1516  << nl << endl;
1517  }
1518  else
1519  {
1520  Info<< "Creating single patch per inter-region interface."
1521  << nl << endl;
1522  }
1523 
1524 
1525  if (insidePoint && largestOnly)
1526  {
1528  << "You cannot specify both -largestOnly"
1529  << " (keep region with most cells)"
1530  << " and -insidePoint (keep region containing point)"
1531  << exit(FatalError);
1532  }
1533 
1534  const word defaultRegionName
1535  (
1536  args.optionLookupOrDefault("defaultRegionName", standardRegionName)
1537  );
1538 
1539  const meshCellZones& cellZones = mesh.cellZones();
1540 
1541  // Existing zoneID
1542  labelList zoneID(mesh.nCells(), -1);
1543  // Neighbour zoneID.
1544  labelList neiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1545  getZoneID(mesh, cellZones, zoneID, neiZoneID);
1546 
1547 
1548 
1549  // Determine per cell the region it belongs to
1550  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1551 
1552  // cellRegion is the labelList with the region per cell.
1553  labelList cellRegion;
1554  // Region per zone
1555  labelList regionToZone;
1556  // Name of region
1558  // Zone to region
1559  labelList zoneToRegion;
1560 
1561  label nCellRegions = 0;
1562  if (useCellZonesOnly)
1563  {
1564  Info<< "Using current cellZones to split mesh into regions."
1565  << " This requires all"
1566  << " cells to be in one and only one cellZone." << nl << endl;
1567 
1568  label unzonedCelli = findIndex(zoneID, -1);
1569  if (unzonedCelli != -1)
1570  {
1572  << "For the cellZonesOnly option all cells "
1573  << "have to be in a cellZone." << endl
1574  << "Cell " << unzonedCelli
1575  << " at" << mesh.cellCentres()[unzonedCelli]
1576  << " is not in a cellZone. There might be more unzoned cells."
1577  << exit(FatalError);
1578  }
1579  cellRegion = zoneID;
1580  nCellRegions = gMax(cellRegion)+1;
1581  regionToZone.setSize(nCellRegions);
1582  regionNames.setSize(nCellRegions);
1583  zoneToRegion.setSize(cellZones.size(), -1);
1584  for (label regioni = 0; regioni < nCellRegions; regioni++)
1585  {
1586  regionToZone[regioni] = regioni;
1587  zoneToRegion[regioni] = regioni;
1588  regionNames[regioni] = cellZones[regioni].name();
1589  }
1590  }
1591  else if (useCellZonesFile)
1592  {
1593  const word zoneFile = args.option("cellZonesFileOnly");
1594  Info<< "Reading split from cellZones file " << zoneFile << endl
1595  << "This requires all"
1596  << " cells to be in one and only one cellZone." << nl << endl;
1597 
1598  meshCellZones newCellZones
1599  (
1600  IOobject
1601  (
1602  zoneFile,
1603  mesh.facesInstance(),
1605  mesh,
1608  false
1609  ),
1610  mesh
1611  );
1612 
1613  labelList newZoneID(mesh.nCells(), -1);
1614  labelList newNeiZoneID(mesh.nFaces()-mesh.nInternalFaces());
1615  getZoneID(mesh, newCellZones, newZoneID, newNeiZoneID);
1616 
1617  label unzonedCelli = findIndex(newZoneID, -1);
1618  if (unzonedCelli != -1)
1619  {
1621  << "For the cellZonesFileOnly option all cells "
1622  << "have to be in a cellZone." << endl
1623  << "Cell " << unzonedCelli
1624  << " at" << mesh.cellCentres()[unzonedCelli]
1625  << " is not in a cellZone. There might be more unzoned cells."
1626  << exit(FatalError);
1627  }
1628  cellRegion = newZoneID;
1629  nCellRegions = gMax(cellRegion)+1;
1630  zoneToRegion.setSize(newCellZones.size(), -1);
1631  regionToZone.setSize(nCellRegions);
1632  regionNames.setSize(nCellRegions);
1633  for (label regioni = 0; regioni < nCellRegions; regioni++)
1634  {
1635  regionToZone[regioni] = regioni;
1636  zoneToRegion[regioni] = regioni;
1637  regionNames[regioni] = newCellZones[regioni].name();
1638  }
1639  }
1640  else
1641  {
1642  // Determine connected regions
1643  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~
1644 
1645  // Mark additional faces that are blocked
1646  boolList blockedFace;
1647 
1648  // Read from faceSet
1649  if (blockedFacesName.size())
1650  {
1651  faceSet blockedFaceSet(mesh, blockedFacesName);
1652  Info<< "Read "
1653  << returnReduce(blockedFaceSet.size(), sumOp<label>())
1654  << " blocked faces from set " << blockedFacesName << nl << endl;
1655 
1656  blockedFace.setSize(mesh.nFaces(), false);
1657 
1658  forAllConstIter(faceSet, blockedFaceSet, iter)
1659  {
1660  blockedFace[iter.key()] = true;
1661  }
1662  }
1663 
1664  // Imply from differing cellZones
1665  if (useCellZones)
1666  {
1667  blockedFace.setSize(mesh.nFaces(), false);
1668 
1669  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1670  {
1671  label own = mesh.faceOwner()[facei];
1672  label nei = mesh.faceNeighbour()[facei];
1673 
1674  if (zoneID[own] != zoneID[nei])
1675  {
1676  blockedFace[facei] = true;
1677  }
1678  }
1679 
1680  // Different cellZones on either side of processor patch.
1681  forAll(neiZoneID, i)
1682  {
1683  label facei = i+mesh.nInternalFaces();
1684 
1685  if (zoneID[mesh.faceOwner()[facei]] != neiZoneID[i])
1686  {
1687  blockedFace[facei] = true;
1688  }
1689  }
1690  }
1691 
1692  // Do a topological walk to determine regions
1693  regionSplit regions(mesh, blockedFace);
1694  nCellRegions = regions.nRegions();
1695  cellRegion.transfer(regions);
1696 
1697  // Make up region names. If possible match them to existing zones.
1698  matchRegions
1699  (
1700  sloppyCellZones,
1701  mesh,
1702  nCellRegions,
1703  cellRegion,
1704  defaultRegionName,
1705 
1706  regionToZone,
1707  regionNames,
1708  zoneToRegion
1709  );
1710 
1711  // Override any default region names if single region selected
1712  if (largestOnly || insidePoint)
1713  {
1714  forAll(regionToZone, regioni)
1715  {
1716  if (regionToZone[regioni] == -1)
1717  {
1718  if (overwrite)
1719  {
1721  }
1722  else if (insidePoint)
1723  {
1724  regionNames[regioni] = "insidePoint";
1725  }
1726  else if (largestOnly)
1727  {
1728  regionNames[regioni] = "largestOnly";
1729  }
1730  }
1731  }
1732  }
1733  }
1734 
1735  Info<< endl << "Number of regions:" << nCellRegions << nl << endl;
1736 
1737 
1738  // Write decomposition to file
1739  writeCellToRegion(mesh, cellRegion);
1740 
1741 
1742 
1743  // Sizes per region
1744  // ~~~~~~~~~~~~~~~~
1745 
1746  labelList regionSizes(nCellRegions, 0);
1747 
1748  forAll(cellRegion, celli)
1749  {
1750  regionSizes[cellRegion[celli]]++;
1751  }
1752  forAll(regionSizes, regioni)
1753  {
1754  reduce(regionSizes[regioni], sumOp<label>());
1755  }
1756 
1757  Info<< "Region\tCells" << nl
1758  << "------\t-----" << endl;
1759 
1760  forAll(regionSizes, regioni)
1761  {
1762  Info<< regioni << '\t' << regionSizes[regioni] << nl;
1763  }
1764  Info<< endl;
1765 
1766 
1767 
1768  // Print region to zone
1769  Info<< "Region\tZone\tName" << nl
1770  << "------\t----\t----" << endl;
1771  forAll(regionToZone, regioni)
1772  {
1773  Info<< regioni << '\t' << regionToZone[regioni] << '\t'
1774  << regionNames[regioni] << nl;
1775  }
1776  Info<< endl;
1777 
1778 
1779 
1780  // Since we're going to mess with patches and zones make sure all
1781  // is synchronised
1782  mesh.boundaryMesh().checkParallelSync(true);
1783  mesh.faceZones().checkParallelSync(true);
1784 
1785 
1786  // Interfaces
1787  // ----------
1788  // per interface:
1789  // - the two regions on either side
1790  // - the name
1791  // - the (global) size
1792  edgeList interfaces;
1793  List<Pair<word>> interfaceNames;
1794  labelList interfaceSizes;
1795  // per face the interface
1796  labelList faceToInterface;
1797 
1798  getInterfaceSizes
1799  (
1800  mesh,
1801  useFaceZones,
1802  cellRegion,
1803  regionNames,
1804 
1805  interfaces,
1806  interfaceNames,
1807  interfaceSizes,
1808  faceToInterface
1809  );
1810 
1811  Info<< "Sizes of interfaces between regions:" << nl << nl
1812  << "Interface\tRegion\tRegion\tFaces" << nl
1813  << "---------\t------\t------\t-----" << endl;
1814 
1815  forAll(interfaces, interI)
1816  {
1817  const edge& e = interfaces[interI];
1818 
1819  Info<< interI
1820  << "\t\t" << e[0] << '\t' << e[1]
1821  << '\t' << interfaceSizes[interI] << nl;
1822  }
1823  Info<< endl;
1824 
1825 
1826  if (detectOnly)
1827  {
1828  return 0;
1829  }
1830 
1831 
1832  // Read objects in time directory
1833  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1834 
1835  IOobjectList objects(mesh, runTime.name());
1836 
1837  if (fields) Info<< "Reading geometric fields" << nl << endl;
1838 
1839  #include "readVolFields.H"
1840  #include "readSurfaceFields.H"
1841  // #include "readPointFields.H"
1842 
1843  Info<< endl;
1844 
1845 
1846  // Remove any demand-driven fields ('S', 'V' etc)
1847  mesh.clearOut();
1848 
1849 
1850  if (nCellRegions == 1)
1851  {
1852  Info<< "Only one region. Doing nothing." << endl;
1853  }
1854  else if (makeCellZones)
1855  {
1856  Info<< "Putting cells into cellZones instead of splitting mesh."
1857  << endl;
1858 
1859  // Check if region overlaps with existing zone. If so keep.
1860 
1861  for (label regioni = 0; regioni < nCellRegions; regioni++)
1862  {
1863  label zoneI = regionToZone[regioni];
1864 
1865  if (zoneI != -1)
1866  {
1867  Info<< " Region " << regioni << " : corresponds to existing"
1868  << " cellZone "
1869  << zoneI << ' ' << cellZones[zoneI].name() << endl;
1870  }
1871  else
1872  {
1873  // Create new cellZone.
1874  labelList regionCells = findIndices(cellRegion, regioni);
1875 
1876  word zoneName = "region" + Foam::name(regioni);
1877 
1878  zoneI = cellZones.findZoneID(zoneName);
1879 
1880  if (zoneI == -1)
1881  {
1882  zoneI = cellZones.size();
1883  mesh.cellZones().setSize(zoneI+1);
1884  mesh.cellZones().set
1885  (
1886  zoneI,
1887  new cellZone
1888  (
1889  zoneName, // name
1890  regionCells, // addressing
1891  zoneI, // index
1892  cellZones // meshCellZones
1893  )
1894  );
1895  }
1896  else
1897  {
1898  mesh.cellZones()[zoneI].clearAddressing();
1899  mesh.cellZones()[zoneI] = regionCells;
1900  }
1901  Info<< " Region " << regioni << " : created new cellZone "
1902  << zoneI << ' ' << cellZones[zoneI].name() << endl;
1903  }
1904  }
1906 
1907  if (!overwrite)
1908  {
1909  runTime++;
1910  mesh.setInstance(runTime.name());
1911  }
1912  else
1913  {
1914  mesh.setInstance(oldInstance);
1915  }
1916 
1917  Info<< "Writing cellZones as new mesh to time " << runTime.name()
1918  << nl << endl;
1919 
1920  mesh.write();
1921 
1922 
1923  // Write cellSets for convenience
1924  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1925 
1926  Info<< "Writing cellSets corresponding to cellZones." << nl << endl;
1927 
1928  forAll(cellZones, zoneI)
1929  {
1930  const cellZone& cz = cellZones[zoneI];
1931 
1932  cellSet(mesh, cz.name(), cz).write();
1933  }
1934  }
1935  else
1936  {
1937  // Add patches for interfaces
1938  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
1939 
1940  // Add all possible patches. Empty ones get filtered later on.
1941  Info<< nl << "Adding patches" << nl << endl;
1942 
1943  labelList interfacePatches
1944  (
1945  addRegionPatches
1946  (
1947  mesh,
1948  regionNames,
1949  interfaces,
1950  interfaceNames
1951  )
1952  );
1953 
1954 
1955  if (!overwrite)
1956  {
1957  runTime++;
1958  }
1959 
1960 
1961  // Create regions
1962  // ~~~~~~~~~~~~~~
1963 
1964  if (insidePoint)
1965  {
1966  const point insidePoint = args.optionRead<point>("insidePoint");
1967 
1968  label regioni = -1;
1969 
1970  (void)mesh.tetBasePtIs();
1971 
1972  label celli = mesh.findCell(insidePoint);
1973 
1974  Info<< nl << "Found point " << insidePoint << " in cell " << celli
1975  << endl;
1976 
1977  if (celli != -1)
1978  {
1979  regioni = cellRegion[celli];
1980  }
1981 
1982  reduce(regioni, maxOp<label>());
1983 
1984  Info<< nl
1985  << "Subsetting region " << regioni
1986  << " containing point " << insidePoint << endl;
1987 
1988  if (regioni == -1)
1989  {
1991  << "Point " << insidePoint
1992  << " is not inside the mesh." << nl
1993  << "Bounding box of the mesh:" << mesh.bounds()
1994  << exit(FatalError);
1995  }
1996 
1997  createAndWriteRegion
1998  (
1999  mesh,
2000  cellRegion,
2001  regionNames,
2002  prefixRegion,
2003  faceToInterface,
2004  interfacePatches,
2005  regioni,
2006  (overwrite ? oldInstance : runTime.name())
2007  );
2008  }
2009  else if (largestOnly)
2010  {
2011  label regioni = findMax(regionSizes);
2012 
2013  Info<< nl
2014  << "Subsetting region " << regioni
2015  << " of size " << regionSizes[regioni]
2016  << " as named region " << regionNames[regioni] << endl;
2017 
2018  createAndWriteRegion
2019  (
2020  mesh,
2021  cellRegion,
2022  regionNames,
2023  prefixRegion,
2024  faceToInterface,
2025  interfacePatches,
2026  regioni,
2027  (overwrite ? oldInstance : runTime.name())
2028  );
2029  }
2030  else
2031  {
2032  // Split all
2033  for (label regioni = 0; regioni < nCellRegions; regioni++)
2034  {
2035  Info<< nl
2036  << "Region " << regioni << nl
2037  << "-------- " << endl;
2038 
2039  createAndWriteRegion
2040  (
2041  mesh,
2042  cellRegion,
2043  regionNames,
2044  prefixRegion,
2045  faceToInterface,
2046  interfacePatches,
2047  regioni,
2048  (overwrite ? oldInstance : runTime.name())
2049  );
2050  }
2051  }
2052  }
2053 
2054  Info<< "End\n" << endl;
2055 
2056  return 0;
2057 }
2058 
2059 
2060 // ************************************************************************* //
Field reading functions for post-processing utilities.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Map from edge (expressed as its endpoints) to value.
Definition: EdgeMap.H:50
Generic GeometricField class.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
An STL-conforming hash table.
Definition: HashTable.H:127
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
writeOption writeOpt() const
Definition: IOobject.H:370
const word & name() const
Return name.
Definition: IOobject.H:310
Input inter-processor communications stream.
Definition: IPstream.H:54
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: MeshZones.C:341
void clearAddressing()
Clear addressing.
Definition: MeshZones.C:387
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: MeshZones.C:221
wordList names() const
Return a list of zone names.
Definition: MeshZones.C:256
bool checkParallelSync(const bool report=false) const
Check whether all procs have all zones and in same order. Return.
Definition: MeshZones.C:428
Output inter-processor communications stream.
Definition: OPstream.H:54
virtual Ostream & write(const char)=0
Write character.
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
void append(T *)
Append an element at the end of the list.
Definition: PtrListI.H:39
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
static const word & system()
Return system name.
Definition: TimePaths.H:112
static int masterNo()
Process index of the master.
Definition: UPstream.H:417
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static int lastSlave(const label communicator=0)
Process index of last slave.
Definition: UPstream.H:452
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static int firstSlave()
Process index of first slave.
Definition: UPstream.H:446
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
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
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
T optionRead(const word &opt) const
Read a value from the named option.
Definition: argListI.H:193
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
const string & option(const word &opt) const
Return the argument string associated with the named option.
Definition: argListI.H:108
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:243
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
A collection of cell labels.
Definition: cellSet.H:51
A subset of mesh cells.
Definition: cellZone.H:64
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1169
const word & name() const
Return const reference to name.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A list of face labels.
Definition: faceSet.H:51
static tmp< VolField< Type > > interpolate(const VolField< Type > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
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 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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:402
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1736
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:286
Wall patch which holds a mapping engine to map values from another patch.
Foam::polyBoundaryMesh.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:1025
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:268
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:447
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1774
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1386
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:1077
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:1019
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1392
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:271
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:411
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:453
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Direct mesh changes based on v1.3 polyTopoChange syntax.
label nCells() const
const vectorField & faceCentres() const
const vectorField & cellCentres() const
label nInternalFaces() const
label nFaces() const
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
virtual bool write(const bool write=true) const
Write using setting from DB.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:118
Given list of cells to remove insert all the topology changes.
Definition: removeCells.H:60
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:436
A class for managing temporary objects.
Definition: tmp.H:55
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
A class for handling words, derived from string.
Definition: word.H:62
const word & name() const
Return name.
Definition: zone.H:147
Foam::word regionName
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const fvPatchList & patches
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:105
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
label findMax(const ListType &, const label start=0)
Find index of max element (and larger than given element).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurrences of given element. Linear search.
const dimensionSet dimless
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static const char nl
Definition: Ostream.H:260
Type gMax(const FieldField< Field, Type > &f)
objects
const Foam::wordList regionNames(args.optionFound("allRegions") ? runTime .controlDict().subDict("regionSolvers").toc() :wordList(1, args.optionFound("region") ? args.optionRead< word >("region") :polyMesh::defaultRegion))
Foam::argList args(argc, argv)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112