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