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