PDRMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  PDRMesh
29 
30 Group
31  grpMeshAdvancedUtilities
32 
33 Description
34  Mesh and field preparation utility for PDR type simulations.
35 
36  Reads
37  - cellSet giving blockedCells
38  - faceSets giving blockedFaces and the patch they should go into
39 
40  NOTE: To avoid exposing wrong fields values faceSets should include
41  faces contained in the blockedCells cellset.
42 
43  - coupledFaces reads coupledFacesSet to introduces mixed-coupled
44  duplicate baffles
45 
46  Subsets out the blocked cells and splits the blockedFaces and updates
47  fields.
48 
49  The face splitting is done by duplicating the faces. No duplication of
50  points for now (so checkMesh will show a lot of 'duplicate face' messages)
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #include "fvMeshSubset.H"
55 #include "argList.H"
56 #include "cellSet.H"
57 #include "BitOps.H"
58 #include "IOobjectList.H"
59 #include "volFields.H"
60 #include "mapPolyMesh.H"
61 #include "faceSet.H"
62 #include "cellSet.H"
63 #include "pointSet.H"
64 #include "syncTools.H"
65 #include "ReadFields.H"
66 #include "polyTopoChange.H"
67 #include "polyModifyFace.H"
68 #include "polyAddFace.H"
69 #include "regionSplit.H"
70 #include "Tuple2.H"
71 #include "cyclicFvPatch.H"
72 
73 using namespace Foam;
74 
75 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76 
77 void modifyOrAddFace
78 (
79  polyTopoChange& meshMod,
80  const face& f,
81  const label facei,
82  const label own,
83  const bool flipFaceFlux,
84  const label newPatchi,
85  const label zoneID,
86  const bool zoneFlip,
87 
88  bitSet& modifiedFace
89 )
90 {
91  if (modifiedFace.set(facei))
92  {
93  // First usage of face. Modify.
94  meshMod.setAction
95  (
97  (
98  f, // modified face
99  facei, // label of face
100  own, // owner
101  -1, // neighbour
102  flipFaceFlux, // face flip
103  newPatchi, // patch for face
104  false, // remove from zone
105  zoneID, // zone for face
106  zoneFlip // face flip in zone
107  )
108  );
109  }
110  else
111  {
112  // Second or more usage of face. Add.
113  meshMod.setAction
114  (
116  (
117  f, // modified face
118  own, // owner
119  -1, // neighbour
120  -1, // master point
121  -1, // master edge
122  facei, // master face
123  flipFaceFlux, // face flip
124  newPatchi, // patch for face
125  zoneID, // zone for face
126  zoneFlip // face flip in zone
127  )
128  );
129  }
130 }
131 
132 
133 template<class Type>
135 (
136  const fvMeshSubset& subsetter,
137  const IOobjectList& objects,
138  const label patchi,
139  const Type& exposedValue
140 )
141 {
143 
144  const fvMesh& baseMesh = subsetter.baseMesh();
145 
146  const UPtrList<const IOobject> fieldObjects
147  (
148  objects.csorted<GeoField>()
149  );
150 
151  PtrList<GeoField> subFields(fieldObjects.size());
152 
153  label nFields = 0;
154  for (const IOobject& io : fieldObjects)
155  {
156  if (!nFields)
157  {
158  Info<< "Subsetting " << GeoField::typeName << " (";
159  }
160  else
161  {
162  Info<< ' ';
163  }
164  Info<< " " << io.name() << endl;
165 
166  // Read unregistered
168  GeoField origField(rio, baseMesh);
169 
170  subFields.set(nFields, subsetter.interpolate(origField));
171  auto& subField = subFields[nFields];
172  ++nFields;
173 
174 
175  // Subsetting adds 'subset' prefix. Rename field to be like original.
176  subField.rename(io.name());
177  subField.writeOpt(IOobjectOption::AUTO_WRITE);
178 
179 
180  // Explicitly set exposed faces (in patchi) to exposedValue.
181  if (patchi >= 0)
182  {
183  fvPatchField<Type>& fld = subField.boundaryFieldRef()[patchi];
184 
185  const label newStart = fld.patch().patch().start();
186  const label oldPatchi = subsetter.patchMap()[patchi];
187 
188  if (oldPatchi == -1)
189  {
190  // New patch. Reset whole value.
191  fld = exposedValue;
192  }
193  else
194  {
195  // Reset faces that originate from different patch
196  // or internal faces.
197 
198  const fvPatchField<Type>& origPfld =
199  origField.boundaryField()[oldPatchi];
200 
201  const label oldSize = origPfld.size();
202  const label oldStart = origPfld.patch().patch().start();
203 
204  forAll(fld, j)
205  {
206  const label oldFacei = subsetter.faceMap()[newStart+j];
207 
208  if (oldFacei < oldStart || oldFacei >= oldStart+oldSize)
209  {
210  fld[j] = exposedValue;
211  }
212  }
213  }
214  }
215  }
216 
217  if (nFields)
218  {
219  Info<< ')' << nl;
220  }
221 
222  return subFields;
223 }
224 
225 
226 template<class Type>
228 (
229  const fvMeshSubset& subsetter,
230  const IOobjectList& objects,
231  const label patchi,
232  const Type& exposedValue
233 )
234 {
236 
237  const fvMesh& baseMesh = subsetter.baseMesh();
238 
239  const UPtrList<const IOobject> fieldObjects
240  (
241  objects.csorted<GeoField>()
242  );
243 
244  PtrList<GeoField> subFields(fieldObjects.size());
245 
246  label nFields = 0;
247  for (const IOobject& io : fieldObjects)
248  {
249  if (!nFields)
250  {
251  Info<< "Subsetting " << GeoField::typeName << " (";
252  }
253  else
254  {
255  Info<< ' ';
256  }
257  Info<< io.name();
258 
259  // Read unregistered
261  GeoField origField(rio, baseMesh);
262 
263  subFields.set(nFields, subsetter.interpolate(origField));
264  auto& subField = subFields[nFields];
265  ++nFields;
266 
267  // Subsetting adds 'subset' prefix. Rename field to be like original.
268  subField.rename(io.name());
269  subField.writeOpt(IOobjectOption::AUTO_WRITE);
270 
271 
272  // Explicitly set exposed faces (in patchi) to exposedValue.
273  if (patchi >= 0)
274  {
275  fvsPatchField<Type>& fld = subField.boundaryFieldRef()[patchi];
276 
277  const label newStart = fld.patch().patch().start();
278  const label oldPatchi = subsetter.patchMap()[patchi];
279 
280  if (oldPatchi == -1)
281  {
282  // New patch. Reset whole value.
283  fld = exposedValue;
284  }
285  else
286  {
287  // Reset faces that originate from different patch
288  // or internal faces.
289 
290  const fvsPatchField<Type>& origPfld =
291  origField.boundaryField()[oldPatchi];
292 
293  const label oldSize = origPfld.size();
294  const label oldStart = origPfld.patch().patch().start();
295 
296  forAll(fld, j)
297  {
298  const label oldFacei = subsetter.faceMap()[newStart+j];
299 
300  if (oldFacei < oldStart || oldFacei >= oldStart+oldSize)
301  {
302  fld[j] = exposedValue;
303  }
304  }
305  }
306  }
307  }
308 
309  if (nFields)
310  {
311  Info<< ')' << nl;
312  }
313 
314  return subFields;
315 }
316 
317 
318 // Faces introduced into zero-sized patches don't get a value at all.
319 // This is hack to set them to an initial value.
320 template<class GeoField>
321 void initCreatedPatches
322 (
323  const fvMesh& mesh,
324  const mapPolyMesh& map,
325  const typename GeoField::value_type initValue
326 )
327 {
328  for (const GeoField& field : mesh.objectRegistry::csorted<GeoField>())
329  {
330  auto& fieldBf = const_cast<GeoField&>(field).boundaryFieldRef();
331 
332  forAll(fieldBf, patchi)
333  {
334  if (map.oldPatchSizes()[patchi] == 0)
335  {
336  // Not mapped.
337  fieldBf[patchi] = initValue;
338 
339  if (fieldBf[patchi].fixesValue())
340  {
341  fieldBf[patchi] == initValue;
342  }
343  }
344  }
345  }
346 }
347 
348 
349 template<class TopoSet>
350 void subsetTopoSets
351 (
352  const fvMesh& mesh,
353  const IOobjectList& objects,
354  const labelList& map,
355  const fvMesh& subMesh,
356  PtrList<TopoSet>& subSets
357 )
358 {
359  // Read original sets
360  PtrList<TopoSet> sets;
361  ReadFields<TopoSet>(objects, sets);
362 
363  subSets.resize_null(sets.size());
364 
365  forAll(sets, seti)
366  {
367  const TopoSet& set = sets[seti];
368 
369  Info<< "Subsetting " << set.type() << ' ' << set.name() << endl;
370 
371  labelHashSet subset(2*min(set.size(), map.size()));
372 
373  // Map the data
374  forAll(map, i)
375  {
376  if (set.contains(map[i]))
377  {
378  subset.insert(i);
379  }
380  }
381 
382  subSets.set
383  (
384  seti,
385  new TopoSet
386  (
387  subMesh,
388  set.name(),
389  std::move(subset),
391  )
392  );
393  }
394 }
395 
396 
397 void createCoupledBaffles
398 (
399  fvMesh& mesh,
400  const labelList& coupledWantedPatch,
401  polyTopoChange& meshMod,
402  bitSet& modifiedFace
403 )
404 {
405  const faceZoneMesh& faceZones = mesh.faceZones();
406 
407  forAll(coupledWantedPatch, facei)
408  {
409  if (coupledWantedPatch[facei] != -1)
410  {
411  const face& f = mesh.faces()[facei];
412  label zoneID = faceZones.whichZone(facei);
413  bool zoneFlip = false;
414 
415  if (zoneID >= 0)
416  {
417  const faceZone& fZone = faceZones[zoneID];
418  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
419  }
420 
421  // Use owner side of face
422  modifyOrAddFace
423  (
424  meshMod,
425  f, // modified face
426  facei, // label of face
427  mesh.faceOwner()[facei], // owner
428  false, // face flip
429  coupledWantedPatch[facei], // patch for face
430  zoneID, // zone for face
431  zoneFlip, // face flip in zone
432  modifiedFace // modify or add status
433  );
434 
435  if (mesh.isInternalFace(facei))
436  {
437  label zoneID = faceZones.whichZone(facei);
438  bool zoneFlip = false;
439 
440  if (zoneID >= 0)
441  {
442  const faceZone& fZone = faceZones[zoneID];
443  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
444  }
445  // Use neighbour side of face
446  modifyOrAddFace
447  (
448  meshMod,
449  f.reverseFace(), // modified face
450  facei, // label of face
451  mesh.faceNeighbour()[facei],// owner
452  false, // face flip
453  coupledWantedPatch[facei], // patch for face
454  zoneID, // zone for face
455  zoneFlip, // face flip in zone
456  modifiedFace // modify or add status
457  );
458  }
459  }
460  }
461 }
462 
463 
464 void createCyclicCoupledBaffles
465 (
466  fvMesh& mesh,
467  const labelList& cyclicMasterPatch,
468  const labelList& cyclicSlavePatch,
469  polyTopoChange& meshMod,
470  bitSet& modifiedFace
471 )
472 {
473  const faceZoneMesh& faceZones = mesh.faceZones();
474 
475  forAll(cyclicMasterPatch, facei)
476  {
477  if (cyclicMasterPatch[facei] != -1)
478  {
479  const face& f = mesh.faces()[facei];
480 
481  label zoneID = faceZones.whichZone(facei);
482  bool zoneFlip = false;
483 
484  if (zoneID >= 0)
485  {
486  const faceZone& fZone = faceZones[zoneID];
487  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
488  }
489 
490  modifyOrAddFace
491  (
492  meshMod,
493  f.reverseFace(), // modified face
494  facei, // label of face
495  mesh.faceNeighbour()[facei], // owner
496  false, // face flip
497  cyclicMasterPatch[facei], // patch for face
498  zoneID, // zone for face
499  zoneFlip, // face flip in zone
500  modifiedFace // modify or add
501  );
502  }
503  }
504 
505  forAll(cyclicSlavePatch, facei)
506  {
507  if (cyclicSlavePatch[facei] != -1)
508  {
509  const face& f = mesh.faces()[facei];
510  if (mesh.isInternalFace(facei))
511  {
512  label zoneID = faceZones.whichZone(facei);
513  bool zoneFlip = false;
514 
515  if (zoneID >= 0)
516  {
517  const faceZone& fZone = faceZones[zoneID];
518  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
519  }
520  // Use owner side of face
521  modifyOrAddFace
522  (
523  meshMod,
524  f, // modified face
525  facei, // label of face
526  mesh.faceOwner()[facei], // owner
527  false, // face flip
528  cyclicSlavePatch[facei], // patch for face
529  zoneID, // zone for face
530  zoneFlip, // face flip in zone
531  modifiedFace // modify or add status
532  );
533  }
534  }
535  }
536 }
537 
538 
539 void createBaffles
540 (
541  fvMesh& mesh,
542  const labelList& wantedPatch,
543  polyTopoChange& meshMod
544 )
545 {
546  const faceZoneMesh& faceZones = mesh.faceZones();
547  Info << "faceZone:createBaffle " << faceZones << endl;
548  forAll(wantedPatch, facei)
549  {
550  if (wantedPatch[facei] != -1)
551  {
552  const face& f = mesh.faces()[facei];
553 
554  label zoneID = faceZones.whichZone(facei);
555  bool zoneFlip = false;
556 
557  if (zoneID >= 0)
558  {
559  const faceZone& fZone = faceZones[zoneID];
560  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
561  }
562 
563  meshMod.setAction
564  (
566  (
567  f, // modified face
568  facei, // label of face
569  mesh.faceOwner()[facei], // owner
570  -1, // neighbour
571  false, // face flip
572  wantedPatch[facei], // patch for face
573  false, // remove from zone
574  zoneID, // zone for face
575  zoneFlip // face flip in zone
576  )
577  );
578 
579  if (mesh.isInternalFace(facei))
580  {
581  label zoneID = faceZones.whichZone(facei);
582  bool zoneFlip = false;
583 
584  if (zoneID >= 0)
585  {
586  const faceZone& fZone = faceZones[zoneID];
587  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
588  }
589 
590  meshMod.setAction
591  (
593  (
594  f.reverseFace(), // modified face
595  mesh.faceNeighbour()[facei],// owner
596  -1, // neighbour
597  -1, // masterPointID
598  -1, // masterEdgeID
599  facei, // masterFaceID,
600  false, // face flip
601  wantedPatch[facei], // patch for face
602  zoneID, // zone for face
603  zoneFlip // face flip in zone
604  )
605  );
606  }
607  }
608  }
609 }
610 
611 
612 // Wrapper around find patch. Also makes sure same patch in parallel.
613 label findPatch(const polyBoundaryMesh& patches, const word& patchName)
614 {
615  const label patchi = patches.findPatchID(patchName);
616 
617  if (patchi == -1)
618  {
620  << "Illegal patch " << patchName
621  << nl << "Valid patches are " << patches.names()
622  << exit(FatalError);
623  }
624 
625  // Check same patch for all procs
626  {
627  const label newPatchi = returnReduce(patchi, minOp<label>());
628 
629  if (newPatchi != patchi)
630  {
632  << "Patch " << patchName
633  << " should have the same patch index on all processors." << nl
634  << "On my processor it has index " << patchi
635  << " ; on some other processor it has index " << newPatchi
636  << exit(FatalError);
637  }
638  }
639  return patchi;
640 }
641 
642 
643 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
644 
645 int main(int argc, char *argv[])
646 {
648  (
649  "Mesh and field preparation utility for PDR type simulations."
650  );
651  #include "addOverwriteOption.H"
652 
653  argList::noFunctionObjects(); // Never use function objects
654 
655  #include "setRootCase.H"
656  #include "createTime.H"
657  #include "createNamedMesh.H"
658 
659  // Read control dictionary
660  // ~~~~~~~~~~~~~~~~~~~~~~~
661 
663  (
664  IOobject
665  (
666  "PDRMeshDict",
667  runTime.system(),
668  mesh,
671  )
672  );
673 
674  // Per faceSet the patch to put the baffles into
675  const List<Pair<word>> setsAndPatches(dict.lookup("blockedFaces"));
676 
677  // Per faceSet the patch to put the coupled baffles into
678  DynamicList<FixedList<word, 3>> coupledAndPatches(10);
679 
680  const dictionary& functionDicts = dict.subDict("coupledFaces");
681 
682  for (const entry& dEntry : functionDicts)
683  {
684  if (!dEntry.isDict()) // Safety
685  {
686  continue;
687  }
688 
689  const word& key = dEntry.keyword();
690  const dictionary& dict = dEntry.dict();
691 
692  const word cyclicName = dict.get<word>("cyclicMasterPatch");
693  const word wallName = dict.get<word>("wallPatch");
694  FixedList<word, 3> nameAndType;
695  nameAndType[0] = key;
696  nameAndType[1] = wallName;
697  nameAndType[2] = cyclicName;
698  coupledAndPatches.append(nameAndType);
699  }
700 
701  forAll(setsAndPatches, setI)
702  {
703  Info<< "Faces in faceSet " << setsAndPatches[setI][0]
704  << " become baffles in patch " << setsAndPatches[setI][1]
705  << endl;
706  }
707 
708  forAll(coupledAndPatches, setI)
709  {
710  Info<< "Faces in faceSet " << coupledAndPatches[setI][0]
711  << " become coupled baffles in patch " << coupledAndPatches[setI][1]
712  << endl;
713  }
714 
715  // All exposed faces that are not explicitly marked to be put into a patch
716  const word defaultPatch(dict.get<word>("defaultPatch"));
717 
718  Info<< "Faces that get exposed become boundary faces in patch "
719  << defaultPatch << endl;
720 
721  const word blockedSetName(dict.get<word>("blockedCells"));
722 
723  Info<< "Reading blocked cells from cellSet " << blockedSetName
724  << endl;
725 
726  const bool overwrite = args.found("overwrite");
727 
728 
729  // Read faceSets, lookup patches
730  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
731 
732  // Check that face sets don't have coincident faces
733  labelList wantedPatch(mesh.nFaces(), -1);
734  forAll(setsAndPatches, setI)
735  {
736  faceSet fSet(mesh, setsAndPatches[setI][0]);
737 
738  label patchi = findPatch
739  (
740  mesh.boundaryMesh(),
741  setsAndPatches[setI][1]
742  );
743 
744  for (const label facei : fSet)
745  {
746  if (wantedPatch[facei] != -1)
747  {
749  << "Face " << facei
750  << " is in faceSet " << setsAndPatches[setI][0]
751  << " destined for patch " << setsAndPatches[setI][1]
752  << " but also in patch " << wantedPatch[facei]
753  << exit(FatalError);
754  }
755  wantedPatch[facei] = patchi;
756  }
757  }
758 
759  // Per face the patch for coupled baffle or -1.
760  labelList coupledWantedPatch(mesh.nFaces(), -1);
761  labelList cyclicWantedPatch_half0(mesh.nFaces(), -1);
762  labelList cyclicWantedPatch_half1(mesh.nFaces(), -1);
763 
764  forAll(coupledAndPatches, setI)
765  {
767  const label cyclicId =
768  findPatch(patches, coupledAndPatches[setI][2]);
769 
770  const label cyclicSlaveId = findPatch
771  (
772  patches,
773  refCast<const cyclicFvPatch>
774  (
775  mesh.boundary()[cyclicId]
776  ).neighbFvPatch().name()
777  );
778 
779  faceSet fSet(mesh, coupledAndPatches[setI][0]);
780  label patchi = findPatch(patches, coupledAndPatches[setI][1]);
781 
782  for (const label facei : fSet)
783  {
784  if (coupledWantedPatch[facei] != -1)
785  {
787  << "Face " << facei
788  << " is in faceSet " << coupledAndPatches[setI][0]
789  << " destined for patch " << coupledAndPatches[setI][1]
790  << " but also in patch " << coupledWantedPatch[facei]
791  << exit(FatalError);
792  }
793 
794  coupledWantedPatch[facei] = patchi;
795  cyclicWantedPatch_half0[facei] = cyclicId;
796  cyclicWantedPatch_half1[facei] = cyclicSlaveId;
797  }
798  }
799 
800  // Exposed faces patch
801  label defaultPatchi = findPatch(mesh.boundaryMesh(), defaultPatch);
802 
803 
804  //
805  // Removing blockedCells (like subsetMesh)
806  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
807  //
808 
809  // Mesh subsetting engine
810  fvMeshSubset subsetter
811  (
812  mesh,
814  (
815  mesh.nCells(),
816  cellSet(mesh, blockedSetName), // Blocked cells as labelHashSet
817  false // on=false: invert logic => retain the unblocked cells
818  ),
819  defaultPatchi,
820  true
821  );
822 
823 
824  // Subset wantedPatch. Note that might also include boundary faces
825  // that have been exposed by subsetting.
826  wantedPatch = IndirectList<label>(wantedPatch, subsetter.faceMap())();
827 
828  coupledWantedPatch = IndirectList<label>
829  (
830  coupledWantedPatch,
831  subsetter.faceMap()
832  )();
833 
834  cyclicWantedPatch_half0 = IndirectList<label>
835  (
836  cyclicWantedPatch_half0,
837  subsetter.faceMap()
838  )();
839 
840  cyclicWantedPatch_half1 = IndirectList<label>
841  (
842  cyclicWantedPatch_half1,
843  subsetter.faceMap()
844  )();
845 
846  // Read all fields in time and constant directories
847  IOobjectList objects(mesh, runTime.timeName());
848  {
849  IOobjectList timeObjects(mesh, mesh.facesInstance());
850 
851  // Transfer specific types
852  forAllIters(timeObjects, iter)
853  {
854  autoPtr<IOobject> objPtr(timeObjects.remove(iter));
855  const auto& obj = *objPtr;
856 
857  if
858  (
859  obj.isHeaderClass<volScalarField>()
860  || obj.isHeaderClass<volVectorField>()
861  || obj.isHeaderClass<volSphericalTensorField>()
862  || obj.isHeaderClass<volTensorField>()
863  || obj.isHeaderClass<volSymmTensorField>()
864  || obj.isHeaderClass<surfaceScalarField>()
865  || obj.isHeaderClass<surfaceVectorField>()
866  || obj.isHeaderClass<surfaceSphericalTensorField>()
867  || obj.isHeaderClass<surfaceSymmTensorField>()
868  || obj.isHeaderClass<surfaceTensorField>()
869  )
870  {
871  objects.add(objPtr);
872  }
873  }
874  }
875 
876  // Read vol fields and subset.
877  PtrList<volScalarField> scalarFlds
878  (
879  subsetVolFields<scalar>
880  (
881  subsetter,
882  objects,
883  defaultPatchi,
884  scalar(Zero)
885  )
886  );
887 
888  PtrList<volVectorField> vectorFlds
889  (
890  subsetVolFields<vector>
891  (
892  subsetter,
893  objects,
894  defaultPatchi,
895  vector(Zero)
896  )
897  );
898 
900  (
901  subsetVolFields<sphericalTensor>
902  (
903  subsetter,
904  objects,
905  defaultPatchi,
907  )
908  );
909 
910  PtrList<volSymmTensorField> symmTensorFlds
911  (
912  subsetVolFields<symmTensor>
913  (
914  subsetter,
915  objects,
916  defaultPatchi,
918  )
919  );
920 
921  PtrList<volTensorField> tensorFlds
922  (
923  subsetVolFields<tensor>
924  (
925  subsetter,
926  objects,
927  defaultPatchi,
928  tensor(Zero)
929  )
930  );
931 
932  // Read surface fields and subset.
933  PtrList<surfaceScalarField> surfScalarFlds
934  (
935  subsetSurfaceFields<scalar>
936  (
937  subsetter,
938  objects,
939  defaultPatchi,
940  scalar(Zero)
941  )
942  );
943 
944  PtrList<surfaceVectorField> surfVectorFlds
945  (
946  subsetSurfaceFields<vector>
947  (
948  subsetter,
949  objects,
950  defaultPatchi,
951  vector(Zero)
952  )
953  );
954 
955  PtrList<surfaceSphericalTensorField> surfSphericalTensorFlds
956  (
957  subsetSurfaceFields<sphericalTensor>
958  (
959  subsetter,
960  objects,
961  defaultPatchi,
963  )
964  );
965 
966  PtrList<surfaceSymmTensorField> surfSymmTensorFlds
967  (
968  subsetSurfaceFields<symmTensor>
969  (
970  subsetter,
971  objects,
972  defaultPatchi,
974  )
975  );
976 
977  PtrList<surfaceTensorField> surfTensorFlds
978  (
979  subsetSurfaceFields<tensor>
980  (
981  subsetter,
982  objects,
983  defaultPatchi,
984  tensor(Zero)
985  )
986  );
987 
988 
989  // Set handling
990  PtrList<cellSet> cellSets;
991  PtrList<faceSet> faceSets;
992  PtrList<pointSet> pointSets;
993  {
994  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
995  subsetTopoSets
996  (
997  mesh,
998  objects,
999  subsetter.cellMap(),
1000  subsetter.subMesh(),
1001  cellSets
1002  );
1003  subsetTopoSets
1004  (
1005  mesh,
1006  objects,
1007  subsetter.faceMap(),
1008  subsetter.subMesh(),
1009  faceSets
1010  );
1011  subsetTopoSets
1012  (
1013  mesh,
1014  objects,
1015  subsetter.pointMap(),
1016  subsetter.subMesh(),
1017  pointSets
1018  );
1019  }
1020 
1021 
1022  if (!overwrite)
1023  {
1024  ++runTime;
1025  }
1026 
1027  Info<< "Writing mesh without blockedCells to time "
1028  << runTime.value() << endl;
1029 
1030  subsetter.subMesh().write();
1031 
1032 
1033  //
1034  // Splitting blockedFaces
1035  // ~~~~~~~~~~~~~~~~~~~~~~
1036  //
1037 
1038  // Synchronize wantedPatch across coupled patches.
1040  (
1041  subsetter.subMesh(),
1042  wantedPatch,
1043  maxEqOp<label>()
1044  );
1045 
1046  // Synchronize coupledWantedPatch across coupled patches.
1048  (
1049  subsetter.subMesh(),
1050  coupledWantedPatch,
1051  maxEqOp<label>()
1052  );
1053 
1054  // Synchronize cyclicWantedPatch across coupled patches.
1056  (
1057  subsetter.subMesh(),
1058  cyclicWantedPatch_half0,
1059  maxEqOp<label>()
1060  );
1061 
1062  // Synchronize cyclicWantedPatch across coupled patches.
1064  (
1065  subsetter.subMesh(),
1066  cyclicWantedPatch_half1,
1067  maxEqOp<label>()
1068  );
1069 
1070  // Topochange container
1071  polyTopoChange meshMod(subsetter.subMesh());
1072 
1073 
1074  // Whether first use of face (modify) or consecutive (add)
1075  bitSet modifiedFace(mesh.nFaces());
1076 
1077  // Create coupled wall-side baffles
1078  createCoupledBaffles
1079  (
1080  subsetter.subMesh(),
1081  coupledWantedPatch,
1082  meshMod,
1083  modifiedFace
1084  );
1085 
1086  // Create coupled master/slave cyclic baffles
1087  createCyclicCoupledBaffles
1088  (
1089  subsetter.subMesh(),
1090  cyclicWantedPatch_half0,
1091  cyclicWantedPatch_half1,
1092  meshMod,
1093  modifiedFace
1094  );
1095 
1096  // Create wall baffles
1097  createBaffles
1098  (
1099  subsetter.subMesh(),
1100  wantedPatch,
1101  meshMod
1102  );
1103 
1104  if (!overwrite)
1105  {
1106  ++runTime;
1107  }
1108 
1109  // Change the mesh. Change points directly (no inflation).
1110  autoPtr<mapPolyMesh> mapPtr =
1111  meshMod.changeMesh(subsetter.subMesh(), false);
1112  mapPolyMesh& map = *mapPtr;
1113 
1114  // Update fields
1115  subsetter.subMesh().updateMesh(map);
1116 
1117  // Fix faces that get mapped to zero-sized patches (these don't get any
1118  // value)
1119  initCreatedPatches<volScalarField>
1120  (
1121  subsetter.subMesh(),
1122  map,
1123  Zero
1124  );
1125  initCreatedPatches<volVectorField>
1126  (
1127  subsetter.subMesh(),
1128  map,
1129  Zero
1130  );
1131  initCreatedPatches<volSphericalTensorField>
1132  (
1133  subsetter.subMesh(),
1134  map,
1135  Zero
1136  );
1137  initCreatedPatches<volSymmTensorField>
1138  (
1139  subsetter.subMesh(),
1140  map,
1141  Zero
1142  );
1143  initCreatedPatches<volTensorField>
1144  (
1145  subsetter.subMesh(),
1146  map,
1147  Zero
1148  );
1149 
1150  initCreatedPatches<surfaceScalarField>
1151  (
1152  subsetter.subMesh(),
1153  map,
1154  Zero
1155  );
1156  initCreatedPatches<surfaceVectorField>
1157  (
1158  subsetter.subMesh(),
1159  map,
1160  Zero
1161  );
1162  initCreatedPatches<surfaceSphericalTensorField>
1163  (
1164  subsetter.subMesh(),
1165  map,
1166  Zero
1167  );
1168  initCreatedPatches<surfaceSymmTensorField>
1169  (
1170  subsetter.subMesh(),
1171  map,
1172  Zero
1173  );
1174  initCreatedPatches<surfaceTensorField>
1175  (
1176  subsetter.subMesh(),
1177  map,
1178  Zero
1179  );
1180 
1181  // Update numbering of topoSets
1182  topoSet::updateMesh(subsetter.subMesh().facesInstance(), map, cellSets);
1183  topoSet::updateMesh(subsetter.subMesh().facesInstance(), map, faceSets);
1184  topoSet::updateMesh(subsetter.subMesh().facesInstance(), map, pointSets);
1185 
1186 
1187  // Move mesh (since morphing might not do this)
1188  if (map.hasMotionPoints())
1189  {
1190  subsetter.subMesh().movePoints(map.preMotionPoints());
1191  }
1192 
1193  Info<< "Writing mesh with split blockedFaces to time " << runTime.value()
1194  << endl;
1195 
1196  subsetter.subMesh().write();
1197 
1198 
1199  //
1200  // Removing inaccessible regions
1201  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1202  //
1203 
1204  // Determine connected regions. regionSplit is the labelList with the
1205  // region per cell.
1206  regionSplit cellRegion(subsetter.subMesh());
1207 
1208  if (cellRegion.nRegions() > 1)
1209  {
1211  << "Removing blocked faces and cells created "
1212  << cellRegion.nRegions()
1213  << " regions that are not connected via a face." << nl
1214  << " This is not supported in solvers." << nl
1215  << " Use" << nl << nl
1216  << " splitMeshRegions <root> <case> -largestOnly" << nl << nl
1217  << " to extract a single region of the mesh." << nl
1218  << " This mesh will be written to a new timedirectory"
1219  << " so might have to be moved back to constant/" << nl
1220  << endl;
1221 
1222  const word startFrom(runTime.controlDict().get<word>("startFrom"));
1223 
1224  if (startFrom != "latestTime")
1225  {
1227  << "To run splitMeshRegions please set your"
1228  << " startFrom entry to latestTime" << endl;
1229  }
1230  }
1231 
1232  Info<< "\nEnd\n" << endl;
1233 
1234  return 0;
1235 }
1236 
1237 
1238 // ************************************************************************* //
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:136
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream. ...
Definition: dictionary.C:367
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition: syncTools.H:465
const Type & value() const noexcept
Return const reference to value.
static tmp< DimensionedField< Type, volMesh > > interpolate(const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelUList &cellMap)
Map volume internal (dimensioned) field.
dictionary dict
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:547
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
rDeltaTY field()
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:498
A list of face labels.
Definition: faceSet.H:47
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Class describing modification of a face.
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:859
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:446
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:929
Required Classes.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const fvPatch & patch() const noexcept
Return the patch.
Generic GeometricField class.
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:397
Field reading functions for post-processing utilities.
const fvMesh & baseMesh() const noexcept
Original mesh.
Definition: fvMeshSubsetI.H:23
Ignore writing from objectRegistry::writeObject()
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
label nFaces() const noexcept
Number of mesh faces.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:47
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
const labelList & faceMap() const
Return face map.
Definition: fvMeshSubsetI.H:65
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const labelList & oldPatchSizes() const noexcept
Return list of the old patch sizes.
Definition: mapPolyMesh.H:776
const fvMesh & subMesh() const
Return reference to subset mesh.
Definition: fvMeshSubsetI.H:41
const labelList & patchMap() const
Return patch map.
Definition: fvMeshSubsetI.H:92
void resize_null(const label newLen)
Set the addressed list to the given size, deleting all existing entries. Afterwards the list contains...
Definition: PtrListI.H:96
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1005
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:55
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
wordList names() const
Return a list of patch names.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const word & system() const noexcept
Return system name.
Definition: TimePathsI.H:118
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:336
const labelList & pointMap() const
Return point map.
Definition: fvMeshSubsetI.H:57
const labelList & cellMap() const
Return cell map.
Definition: fvMeshSubsetI.H:84
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
const dictionary & controlDict() const noexcept
Return read access to the controlDict dictionary.
Definition: Time.H:539
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1113
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:671
labelList f(nPoints)
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: fvMeshSubset.H:75
bool hasMotionPoints() const noexcept
Has valid preMotionPoints?
Definition: mapPolyMesh.H:767
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:406
List< T > subset(const BoolListType &select, const UList< T > &input, const bool invert=false)
Extract elements of the input list when select is true.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
label nCells() const noexcept
Number of mesh cells.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
A collection of cell labels.
Definition: cellSet.H:47
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Direct mesh changes based on v1.3 polyTopoChange syntax.
const polyBoundaryMesh & patches
Automatically write from objectRegistry::writeObject()
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:367
messageStream Info
Information stream (stdout output on master, null elsewhere)
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
bitSet create(const label n, const labelHashSet &locations, const bool on=true)
Create a bitSet with length n with the specified on locations.
Definition: BitOps.C:228
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition: fvPatch.H:202
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
virtual void updateMesh(const mapPolyMesh &morphMap)
Update any stored data for new labels. Not implemented.
Definition: topoSet.C:687
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
A List with indirect addressing.
Definition: IndirectList.H:60
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
const pointField & preMotionPoints() const noexcept
Pre-motion point positions.
Definition: mapPolyMesh.H:759
Do not request registration (bool: false)
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
Namespace for OpenFOAM.
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:63
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127