PDRMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  PDRMesh
26 
27 Description
28  Mesh and field preparation utility for PDR type simulations.
29 
30  Reads
31  - cellSet giving blockedCells
32  - faceSets giving blockedFaces and the patch they should go into
33 
34  NOTE: To avoid exposing wrong fields values faceSets should include
35  faces contained in the blockedCells cellset.
36 
37  - coupledFaces reads coupledFacesSet to introduces mixe-coupled baffles
38 
39  Subsets out the blocked cells and splits the blockedFaces and updates
40  fields.
41 
42  The face splitting is done by duplicating the faces. No duplication of
43  points for now (so checkMesh will show a lot of 'duplicate face' messages)
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #include "fvMeshSubset.H"
48 #include "argList.H"
49 #include "cellSet.H"
50 #include "IOobjectList.H"
51 #include "volFields.H"
52 #include "polyTopoChangeMap.H"
53 #include "faceSet.H"
54 #include "cellSet.H"
55 #include "syncTools.H"
56 #include "polyTopoChange.H"
57 #include "regionSplit.H"
58 #include "Tuple2.H"
59 #include "cyclicFvPatch.H"
60 
61 using namespace Foam;
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 void modifyOrAddFace
66 (
67  polyTopoChange& meshMod,
68  const face& f,
69  const label facei,
70  const label own,
71  const bool flipFaceFlux,
72  const label newPatchi,
73 
74  PackedBoolList& modifiedFace
75 )
76 {
77  if (!modifiedFace[facei])
78  {
79  // First usage of face. Modify.
80  meshMod.modifyFace
81  (
82  f, // modified face
83  facei, // label of face
84  own, // owner
85  -1, // neighbour
86  flipFaceFlux, // face flip
87  newPatchi // patch for face
88  );
89  modifiedFace[facei] = 1;
90  }
91  else
92  {
93  // Second or more usage of face. Add.
94  meshMod.addFace
95  (
96  f, // modified face
97  own, // owner
98  -1, // neighbour
99  facei, // master face
100  flipFaceFlux, // face flip
101  newPatchi // patch for face
102  );
103  }
104 }
105 
106 
107 template<class Type>
108 void subsetVolFields
109 (
110  const fvMeshSubset& subsetter,
111  const IOobjectList& objectsList,
112  const label patchi,
113  const Type& exposedValue,
114  const word GeomVolType,
115  PtrList<VolField<Type>>& subFields
116 )
117 {
118  const fvMesh& baseMesh = subsetter.baseMesh();
119 
120  label i = 0;
121 
122  forAllConstIter(IOobjectList , objectsList, iter)
123  {
124  if (iter()->headerClassName() == GeomVolType)
125  {
126  const word fieldName = iter()->name();
127 
128  Info<< "Subsetting field " << fieldName << endl;
129 
131  (
132  *iter(),
133  baseMesh
134  );
135 
136  subFields.set(i, subsetter.interpolate(volField));
137 
138  // Explicitly set exposed faces (in patchi) to exposedValue.
139  if (patchi >= 0)
140  {
142  subFields[i++].boundaryFieldRef()[patchi];
143 
144  label newStart = fld.patch().patch().start();
145 
146  label oldPatchi = subsetter.patchMap()[patchi];
147 
148  if (oldPatchi == -1)
149  {
150  // New patch. Reset whole value.
151  fld = exposedValue;
152  }
153  else
154  {
155  // Reset those faces that originate from different patch
156  // or internal faces.
157  label oldSize = volField.boundaryField()[oldPatchi].size();
158  label oldStart = volField.boundaryField()
159  [
160  oldPatchi
161  ].patch().patch().start();
162 
163  forAll(fld, j)
164  {
165  label oldFacei = subsetter.faceMap()[newStart+j];
166 
167  if (oldFacei < oldStart || oldFacei >= oldStart+oldSize)
168  {
169  fld[j] = exposedValue;
170  }
171  }
172  }
173  }
174  }
175  }
176 }
177 
178 
179 template<class Type>
180 void subsetSurfaceFields
181 (
182  const fvMeshSubset& subsetter,
183  const IOobjectList& objectsList,
184  const label patchi,
185  const Type& exposedValue,
186  const word GeomSurfType,
187  PtrList<SurfaceField<Type>>& subFields
188 )
189 {
190  const fvMesh& baseMesh = subsetter.baseMesh();
191 
192  label i(0);
193 
194  forAllConstIter(IOobjectList , objectsList, iter)
195  {
196  if (iter()->headerClassName() == GeomSurfType)
197  {
198  const word& fieldName = iter.key();
199 
200  Info<< "Subsetting field " << fieldName << endl;
201 
203  (
204  *iter(),
205  baseMesh
206  );
207 
208  subFields.set(i, subsetter.interpolate(volField));
209 
210 
211  // Explicitly set exposed faces (in patchi) to exposedValue.
212  if (patchi >= 0)
213  {
215  subFields[i++].boundaryFieldRef()[patchi];
216 
217  label newStart = fld.patch().patch().start();
218 
219  label oldPatchi = subsetter.patchMap()[patchi];
220 
221  if (oldPatchi == -1)
222  {
223  // New patch. Reset whole value.
224  fld = exposedValue;
225  }
226  else
227  {
228  // Reset those faces that originate from different patch
229  // or internal faces.
230  label oldSize = volField.boundaryField()[oldPatchi].size();
231  label oldStart = volField.boundaryField()
232  [
233  oldPatchi
234  ].patch().patch().start();
235 
236  forAll(fld, j)
237  {
238  label oldFacei = subsetter.faceMap()[newStart+j];
239 
240  if (oldFacei < oldStart || oldFacei >= oldStart+oldSize)
241  {
242  fld[j] = exposedValue;
243  }
244  }
245  }
246  }
247  }
248  }
249 }
250 
251 
252 // Faces introduced into zero-sized patches don't get a value at all.
253 // This is hack to set them to an initial value.
254 template<class GeoField>
255 void initCreatedPatches
256 (
257  const fvMesh& mesh,
258  const polyTopoChangeMap& map,
259  const typename GeoField::value_type initValue
260 )
261 {
263  (
264  mesh.objectRegistry::lookupClass<GeoField>()
265  );
266 
267  for
268  (
269  typename HashTable<const GeoField*>::
270  iterator fieldIter = fields.begin();
271  fieldIter != fields.end();
272  ++fieldIter
273  )
274  {
275  GeoField& field = const_cast<GeoField&>(*fieldIter());
276 
277  typename GeoField::Boundary& fieldBf =
278  field.boundaryFieldRef();
279 
280  forAll(fieldBf, patchi)
281  {
282  if (map.oldPatchSizes()[patchi] == 0)
283  {
284  // Not mapped.
285  fieldBf[patchi] = initValue;
286 
287  if (fieldBf[patchi].fixesValue())
288  {
289  fieldBf[patchi] == initValue;
290  }
291  }
292  }
293  }
294 }
295 
296 
297 void createCoupledBaffles
298 (
299  fvMesh& mesh,
300  const labelList& coupledWantedPatch,
301  polyTopoChange& meshMod,
302  PackedBoolList& modifiedFace
303 )
304 {
305  forAll(coupledWantedPatch, facei)
306  {
307  if (coupledWantedPatch[facei] != -1)
308  {
309  const face& f = mesh.faces()[facei];
310 
311  // Use owner side of face
312  modifyOrAddFace
313  (
314  meshMod,
315  f, // modified face
316  facei, // label of face
317  mesh.faceOwner()[facei], // owner
318  false, // face flip
319  coupledWantedPatch[facei], // patch for face
320  modifiedFace // modify or add status
321  );
322 
323  if (mesh.isInternalFace(facei))
324  {
325  // Use neighbour side of face
326  modifyOrAddFace
327  (
328  meshMod,
329  f.reverseFace(), // modified face
330  facei, // label of face
331  mesh.faceNeighbour()[facei],// owner
332  false, // face flip
333  coupledWantedPatch[facei], // patch for face
334  modifiedFace // modify or add status
335  );
336  }
337  }
338  }
339 }
340 
341 
342 void createCyclicCoupledBaffles
343 (
344  fvMesh& mesh,
345  const labelList& cyclicMasterPatch,
346  const labelList& cyclicSlavePatch,
347  polyTopoChange& meshMod,
348  PackedBoolList& modifiedFace
349 )
350 {
351  forAll(cyclicMasterPatch, facei)
352  {
353  if (cyclicMasterPatch[facei] != -1)
354  {
355  const face& f = mesh.faces()[facei];
356 
357  modifyOrAddFace
358  (
359  meshMod,
360  f.reverseFace(), // modified face
361  facei, // label of face
362  mesh.faceNeighbour()[facei], // owner
363  false, // face flip
364  cyclicMasterPatch[facei], // patch for face
365  modifiedFace // modify or add
366  );
367  }
368  }
369 
370  forAll(cyclicSlavePatch, facei)
371  {
372  if (cyclicSlavePatch[facei] != -1)
373  {
374  const face& f = mesh.faces()[facei];
375  if (mesh.isInternalFace(facei))
376  {
377  // Use owner side of face
378  modifyOrAddFace
379  (
380  meshMod,
381  f, // modified face
382  facei, // label of face
383  mesh.faceOwner()[facei], // owner
384  false, // face flip
385  cyclicSlavePatch[facei], // patch for face
386  modifiedFace // modify or add status
387  );
388  }
389  }
390  }
391 }
392 
393 
394 void createBaffles
395 (
396  fvMesh& mesh,
397  const labelList& wantedPatch,
398  polyTopoChange& meshMod
399 )
400 {
401  forAll(wantedPatch, facei)
402  {
403  if (wantedPatch[facei] != -1)
404  {
405  const face& f = mesh.faces()[facei];
406 
407  meshMod.modifyFace
408  (
409  f, // modified face
410  facei, // label of face
411  mesh.faceOwner()[facei], // owner
412  -1, // neighbour
413  false, // face flip
414  wantedPatch[facei] // patch for face
415  );
416 
417  if (mesh.isInternalFace(facei))
418  {
419  meshMod.addFace
420  (
421  f.reverseFace(), // modified face
422  mesh.faceNeighbour()[facei],// owner
423  -1, // neighbour
424  facei, // masterFaceID,
425  false, // face flip
426  wantedPatch[facei] // patch for face
427  );
428  }
429  }
430  }
431 }
432 
433 
434 // Wrapper around find patch. Also makes sure same patch in parallel.
435 label findPatch(const polyBoundaryMesh& patches, const word& patchName)
436 {
437  label patchi = patches.findIndex(patchName);
438 
439  if (patchi == -1)
440  {
442  << "Illegal patch " << patchName
443  << nl << "Valid patches are " << patches.names()
444  << exit(FatalError);
445  }
446 
447  // Check same patch for all procs
448  {
449  label newPatch = patchi;
450  reduce(newPatch, minOp<label>());
451 
452  if (newPatch != patchi)
453  {
455  << "Patch " << patchName
456  << " should have the same patch index on all processors." << nl
457  << "On my processor it has index " << patchi
458  << " ; on some other processor it has index " << newPatch
459  << exit(FatalError);
460  }
461  }
462  return patchi;
463 }
464 
465 
466 
467 int main(int argc, char *argv[])
468 {
469  #include "addOverwriteOption.H"
470  #include "setRootCase.H"
472  #include "createMeshNoChangers.H"
473 
474  // Read control dictionary
475  // ~~~~~~~~~~~~~~~~~~~~~~~
476 
478  (
479  IOobject
480  (
481  "PDRMeshDict",
482  runTime.system(),
483  mesh,
486  )
487  );
488 
489  // Per faceSet the patch to put the baffles into
490  const List<Pair<word>> setsAndPatches(dict.lookup("blockedFaces"));
491 
492  // Per faceSet the patch to put the coupled baffles into
493  DynamicList<FixedList<word, 3>> coupledAndPatches(10);
494  const dictionary& functionDicts = dict.subDict("coupledFaces");
495  forAllConstIter(dictionary, functionDicts, iter)
496  {
497  // safety:
498  if (!iter().isDict())
499  {
500  continue;
501  }
502  const word& key = iter().keyword();
503 
504  const dictionary& dict = iter().dict();
505  const word cyclicName = dict.lookup("cyclicMasterPatchName");
506  const word wallName = dict.lookup("wallPatchName");
507  FixedList<word, 3> nameAndType;
508  nameAndType[0] = key;
509  nameAndType[1] = wallName;
510  nameAndType[2] = cyclicName;
511  coupledAndPatches.append(nameAndType);
512  }
513 
514  forAll(setsAndPatches, setI)
515  {
516  Info<< "Faces in faceSet " << setsAndPatches[setI][0]
517  << " become baffles in patch " << setsAndPatches[setI][1]
518  << endl;
519  }
520 
521  forAll(coupledAndPatches, setI)
522  {
523  Info<< "Faces in faceSet " << coupledAndPatches[setI][0]
524  << " become coupled baffles in patch " << coupledAndPatches[setI][1]
525  << endl;
526  }
527 
528  // All exposed faces that are not explicitly marked to be put into a patch
529  const word defaultPatch(dict.lookup("defaultPatch"));
530 
531  Info<< "Faces that get exposed become boundary faces in patch "
532  << defaultPatch << endl;
533 
534  const word blockedSetName(dict.lookup("blockedCells"));
535 
536  Info<< "Reading blocked cells from cellSet " << blockedSetName
537  << endl;
538 
539  const bool overwrite = args.optionFound("overwrite");
540 
541 
542  // Read faceSets, lookup patches
543  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
544 
545  // Check that face sets don't have coincident faces
546  labelList wantedPatch(mesh.nFaces(), -1);
547  forAll(setsAndPatches, setI)
548  {
549  faceSet fSet(mesh, setsAndPatches[setI][0]);
550 
551  label patchi = findPatch
552  (
553  mesh.boundaryMesh(),
554  setsAndPatches[setI][1]
555  );
556 
557  forAllConstIter(faceSet, fSet, iter)
558  {
559  if (wantedPatch[iter.key()] != -1)
560  {
562  << "Face " << iter.key()
563  << " is in faceSet " << setsAndPatches[setI][0]
564  << " destined for patch " << setsAndPatches[setI][1]
565  << " but also in patch " << wantedPatch[iter.key()]
566  << exit(FatalError);
567  }
568  wantedPatch[iter.key()] = patchi;
569  }
570  }
571 
572  // Per face the patch for coupled baffle or -1.
573  labelList coupledWantedPatch(mesh.nFaces(), -1);
574  labelList cyclicWantedPatch_half0(mesh.nFaces(), -1);
575  labelList cyclicWantedPatch_half1(mesh.nFaces(), -1);
576 
577  forAll(coupledAndPatches, setI)
578  {
579  const polyBoundaryMesh& patches = mesh.boundaryMesh();
580  const label cyclicId =
581  findPatch(patches, coupledAndPatches[setI][2]);
582 
583  const label cyclicSlaveId = findPatch
584  (
585  patches,
586  refCast<const cyclicFvPatch>
587  (
588  mesh.boundary()[cyclicId]
589  ).neighbFvPatch().name()
590  );
591 
592  faceSet fSet(mesh, coupledAndPatches[setI][0]);
593  label patchi = findPatch(patches, coupledAndPatches[setI][1]);
594 
595  forAllConstIter(faceSet, fSet, iter)
596  {
597  if (coupledWantedPatch[iter.key()] != -1)
598  {
600  << "Face " << iter.key()
601  << " is in faceSet " << coupledAndPatches[setI][0]
602  << " destined for patch " << coupledAndPatches[setI][1]
603  << " but also in patch " << coupledWantedPatch[iter.key()]
604  << exit(FatalError);
605  }
606  coupledWantedPatch[iter.key()] = patchi;
607  cyclicWantedPatch_half0[iter.key()] = cyclicId;
608  cyclicWantedPatch_half1[iter.key()] = cyclicSlaveId;
609  }
610  }
611 
612  // Exposed faces patch
613  label defaultPatchi = findPatch(mesh.boundaryMesh(), defaultPatch);
614 
615 
616  //
617  // Removing blockedCells (like subsetMesh)
618  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
619  //
620 
621  // Create mesh subsetting engine
622  fvMeshSubset subsetter(mesh);
623 
624  {
625 
626  cellSet blockedCells(mesh, blockedSetName);
627 
628  // invert
629  blockedCells.invert(mesh.nCells());
630 
631  // Create subsetted mesh.
632  subsetter.setLargeCellSubset(blockedCells, defaultPatchi, true);
633  }
634 
635 
636  // Subset wantedPatch. Note that might also include boundary faces
637  // that have been exposed by subsetting.
638  wantedPatch = IndirectList<label>(wantedPatch, subsetter.faceMap())();
639 
640  coupledWantedPatch = IndirectList<label>
641  (
642  coupledWantedPatch,
643  subsetter.faceMap()
644  )();
645 
646  cyclicWantedPatch_half0 = IndirectList<label>
647  (
648  cyclicWantedPatch_half0,
649  subsetter.faceMap()
650  )();
651 
652  cyclicWantedPatch_half1 = IndirectList<label>
653  (
654  cyclicWantedPatch_half1,
655  subsetter.faceMap()
656  )();
657 
658  // Read all fields in time and constant directories
659  IOobjectList objects(mesh, runTime.name());
660  IOobjectList timeObjects(IOobjectList(mesh, mesh.facesInstance()));
661  forAllConstIter(IOobjectList, timeObjects, iter)
662  {
663  if
664  (
665  iter()->headerClassName() == volScalarField::typeName
666  || iter()->headerClassName() == volVectorField::typeName
667  || iter()->headerClassName() == volSphericalTensorField::typeName
668  || iter()->headerClassName() == volTensorField::typeName
669  || iter()->headerClassName() == volSymmTensorField::typeName
670  || iter()->headerClassName() == surfaceScalarField::typeName
671  || iter()->headerClassName() == surfaceVectorField::typeName
672  || iter()->headerClassName()
674  || iter()->headerClassName() == surfaceSymmTensorField::typeName
675  || iter()->headerClassName() == surfaceTensorField::typeName
676  )
677  {
678  objects.add(*iter());
679  }
680  }
681 
682  // Read vol fields and subset.
683 
684  wordList scalarNames(objects.names(volScalarField::typeName));
685  PtrList<volScalarField> scalarFlds(scalarNames.size());
686  subsetVolFields
687  (
688  subsetter,
689  objects,
690  defaultPatchi,
691  scalar(Zero),
693  scalarFlds
694  );
695 
696  wordList vectorNames(objects.names(volVectorField::typeName));
697  PtrList<volVectorField> vectorFlds(vectorNames.size());
698  subsetVolFields
699  (
700  subsetter,
701  objects,
702  defaultPatchi,
703  vector(Zero),
705  vectorFlds
706  );
707 
708  wordList sphericalTensorNames
709  (
711  );
712  PtrList<volSphericalTensorField> sphericalTensorFlds
713  (
714  sphericalTensorNames.size()
715  );
716  subsetVolFields
717  (
718  subsetter,
719  objects,
720  defaultPatchi,
723  sphericalTensorFlds
724  );
725 
726  wordList symmTensorNames(objects.names(volSymmTensorField::typeName));
727  PtrList<volSymmTensorField> symmTensorFlds(symmTensorNames.size());
728  subsetVolFields
729  (
730  subsetter,
731  objects,
732  defaultPatchi,
733  symmTensor(Zero),
735  symmTensorFlds
736  );
737 
738  wordList tensorNames(objects.names(volTensorField::typeName));
739  PtrList<volTensorField> tensorFlds(tensorNames.size());
740  subsetVolFields
741  (
742  subsetter,
743  objects,
744  defaultPatchi,
745  tensor(Zero),
747  tensorFlds
748  );
749 
750  // Read surface fields and subset.
751 
752  wordList surfScalarNames(objects.names(surfaceScalarField::typeName));
753  PtrList<surfaceScalarField> surfScalarFlds(surfScalarNames.size());
754  subsetSurfaceFields
755  (
756  subsetter,
757  objects,
758  defaultPatchi,
759  scalar(Zero),
761  surfScalarFlds
762  );
763 
764  wordList surfVectorNames(objects.names(surfaceVectorField::typeName));
765  PtrList<surfaceVectorField> surfVectorFlds(surfVectorNames.size());
766  subsetSurfaceFields
767  (
768  subsetter,
769  objects,
770  defaultPatchi,
771  vector(Zero),
773  surfVectorFlds
774  );
775 
776  wordList surfSphericalTensorNames
777  (
779  );
780  PtrList<surfaceSphericalTensorField> surfSphericalTensorFlds
781  (
782  surfSphericalTensorNames.size()
783  );
784  subsetSurfaceFields
785  (
786  subsetter,
787  objects,
788  defaultPatchi,
791  surfSphericalTensorFlds
792  );
793 
794  wordList surfSymmTensorNames
795  (
797  );
798 
799  PtrList<surfaceSymmTensorField> surfSymmTensorFlds
800  (
801  surfSymmTensorNames.size()
802  );
803 
804  subsetSurfaceFields
805  (
806  subsetter,
807  objects,
808  defaultPatchi,
809  symmTensor(Zero),
811  surfSymmTensorFlds
812  );
813 
814  wordList surfTensorNames(objects.names(surfaceTensorField::typeName));
815  PtrList<surfaceTensorField> surfTensorFlds(surfTensorNames.size());
816  subsetSurfaceFields
817  (
818  subsetter,
819  objects,
820  defaultPatchi,
821  tensor(Zero),
823  surfTensorFlds
824  );
825 
826  if (!overwrite)
827  {
828  runTime++;
829  }
830 
831  Info<< "Writing mesh without blockedCells to time "
832  << runTime.userTimeName() << endl;
833 
834  // Subsetting adds 'subset' prefix. Rename field to be like original.
835  forAll(scalarFlds, i)
836  {
837  scalarFlds[i].rename(scalarNames[i]);
838  scalarFlds[i].writeOpt() = IOobject::AUTO_WRITE;
839  scalarFlds[i].checkIn();
840  }
841  forAll(vectorFlds, i)
842  {
843  vectorFlds[i].rename(vectorNames[i]);
844  vectorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
845  vectorFlds[i].checkIn();
846  }
847  forAll(sphericalTensorFlds, i)
848  {
849  sphericalTensorFlds[i].rename(sphericalTensorNames[i]);
850  sphericalTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
851  sphericalTensorFlds[i].checkIn();
852  }
853  forAll(symmTensorFlds, i)
854  {
855  symmTensorFlds[i].rename(symmTensorNames[i]);
856  symmTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
857  symmTensorFlds[i].checkIn();
858  }
859  forAll(tensorFlds, i)
860  {
861  tensorFlds[i].rename(tensorNames[i]);
862  tensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
863  tensorFlds[i].checkIn();
864  }
865 
866  // Surface ones.
867  forAll(surfScalarFlds, i)
868  {
869  surfScalarFlds[i].rename(surfScalarNames[i]);
870  surfScalarFlds[i].writeOpt() = IOobject::AUTO_WRITE;
871  surfScalarFlds[i].checkIn();
872  }
873  forAll(surfVectorFlds, i)
874  {
875  surfVectorFlds[i].rename(surfVectorNames[i]);
876  surfVectorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
877  surfVectorFlds[i].checkIn();
878  }
879  forAll(surfSphericalTensorFlds, i)
880  {
881  surfSphericalTensorFlds[i].rename(surfSphericalTensorNames[i]);
882  surfSphericalTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
883  surfSphericalTensorFlds[i].checkIn();
884  }
885  forAll(surfSymmTensorFlds, i)
886  {
887  surfSymmTensorFlds[i].rename(surfSymmTensorNames[i]);
888  surfSymmTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
889  surfSymmTensorFlds[i].checkIn();
890  }
891  forAll(surfTensorNames, i)
892  {
893  surfTensorFlds[i].rename(surfTensorNames[i]);
894  surfTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
895  surfTensorFlds[i].checkIn();
896  }
897 
898  subsetter.subMesh().write();
899 
900 
901  //
902  // Splitting blockedFaces
903  // ~~~~~~~~~~~~~~~~~~~~~~
904  //
905 
906  // Synchronise wantedPatch across coupled patches.
908  (
909  subsetter.subMesh(),
910  wantedPatch,
912  );
913 
914  // Synchronise coupledWantedPatch across coupled patches.
916  (
917  subsetter.subMesh(),
918  coupledWantedPatch,
920  );
921 
922  // Synchronise cyclicWantedPatch across coupled patches.
924  (
925  subsetter.subMesh(),
926  cyclicWantedPatch_half0,
928  );
929 
930  // Synchronise cyclicWantedPatch across coupled patches.
932  (
933  subsetter.subMesh(),
934  cyclicWantedPatch_half1,
936  );
937 
938  // Topochange container
939  polyTopoChange meshMod(subsetter.subMesh());
940 
941 
942  // Whether first use of face (modify) or consecutive (add)
943  PackedBoolList modifiedFace(mesh.nFaces());
944 
945  // Create coupled wall-side baffles
946  createCoupledBaffles
947  (
948  subsetter.subMesh(),
949  coupledWantedPatch,
950  meshMod,
951  modifiedFace
952  );
953 
954  // Create coupled master/slave cyclic baffles
955  createCyclicCoupledBaffles
956  (
957  subsetter.subMesh(),
958  cyclicWantedPatch_half0,
959  cyclicWantedPatch_half1,
960  meshMod,
961  modifiedFace
962  );
963 
964  // Create wall baffles
965  createBaffles
966  (
967  subsetter.subMesh(),
968  wantedPatch,
969  meshMod
970  );
971 
972  if (!overwrite)
973  {
974  runTime++;
975  }
976 
977  // Change the mesh. Change points directly
978  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(subsetter.subMesh());
979 
980  // Update fields
981  subsetter.subMesh().topoChange(map);
982 
983  // Fix faces that get mapped to zero-sized patches (these don't get any
984  // value)
985  initCreatedPatches<volScalarField>
986  (
987  subsetter.subMesh(),
988  map,
989  0.0
990  );
991  initCreatedPatches<volVectorField>
992  (
993  subsetter.subMesh(),
994  map,
995  Zero
996  );
997  initCreatedPatches<volSphericalTensorField>
998  (
999  subsetter.subMesh(),
1000  map,
1001  Zero
1002  );
1003  initCreatedPatches<volSymmTensorField>
1004  (
1005  subsetter.subMesh(),
1006  map,
1007  Zero
1008  );
1009  initCreatedPatches<volTensorField>
1010  (
1011  subsetter.subMesh(),
1012  map,
1013  Zero
1014  );
1015 
1016  initCreatedPatches<surfaceScalarField>
1017  (
1018  subsetter.subMesh(),
1019  map,
1020  0.0
1021  );
1022  initCreatedPatches<surfaceVectorField>
1023  (
1024  subsetter.subMesh(),
1025  map,
1026  Zero
1027  );
1028  initCreatedPatches<surfaceSphericalTensorField>
1029  (
1030  subsetter.subMesh(),
1031  map,
1032  Zero
1033  );
1034  initCreatedPatches<surfaceSymmTensorField>
1035  (
1036  subsetter.subMesh(),
1037  map,
1038  Zero
1039  );
1040  initCreatedPatches<surfaceTensorField>
1041  (
1042  subsetter.subMesh(),
1043  map,
1044  Zero
1045  );
1046 
1047  Info<< "Writing mesh with split blockedFaces to time "
1048  << runTime.userTimeName() << endl;
1049 
1050  subsetter.subMesh().write();
1051 
1052 
1053  //
1054  // Removing inaccessible regions
1055  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1056  //
1057 
1058  // Determine connected regions. regionSplit is the labelList with the
1059  // region per cell.
1060  regionSplit cellRegion(subsetter.subMesh());
1061 
1062  if (cellRegion.nRegions() > 1)
1063  {
1065  << "Removing blocked faces and cells created "
1066  << cellRegion.nRegions()
1067  << " regions that are not connected via a face." << nl
1068  << " This is not supported in solvers." << nl
1069  << " Use" << nl << nl
1070  << " splitMeshRegions <root> <case> -largestOnly" << nl << nl
1071  << " to extract a single region of the mesh." << nl
1072  << " This mesh will be written to a new timedirectory"
1073  << " so might have to be moved back to constant/" << nl
1074  << endl;
1075 
1076  word startFrom(runTime.controlDict().lookup("startFrom"));
1077 
1078  if (startFrom != "latestTime")
1079  {
1081  << "To run splitMeshRegions please set your"
1082  << " startFrom entry to latestTime" << endl;
1083  }
1084  }
1085 
1086  Info << nl << "End" << endl;
1087 
1088  return 0;
1089 }
1090 
1091 
1092 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
static const char *const typeName
Definition: Field.H:106
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:78
Generic GeometricField class.
An STL-conforming hash table.
Definition: HashTable.H:127
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
A List with indirect addressing.
Definition: IndirectList.H:105
A bit-packed bool list.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
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
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
A list of face labels.
Definition: faceSet.H:51
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:74
const labelList & faceMap() const
Return face map.
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:222
void setLargeCellSubset(const labelList &region, const label currentRegion, const label patchID=-1, const bool syncCouples=true)
Set the subset from all cells with region == currentRegion.
Definition: fvMeshSubset.C:893
const fvMesh & subMesh() const
Return reference to subset mesh.
const labelList & patchMap() const
Return patch map.
static tmp< VolField< Type > > interpolate(const VolField< Type > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1737
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1264
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
Foam::polyBoundaryMesh.
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:971
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1339
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1345
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & oldPatchSizes() const
Return list of the old patch sizes.
Direct mesh changes based on v1.3 polyTopoChange syntax.
autoPtr< polyTopoChangeMap > changeMesh(polyMesh &mesh, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
void modifyFace(const face &f, const label facei, const label own, const label nei, const bool flipFaceFlux, const label patchID)
Modify vertices or cell of face.
label addFace(const face &f, const label own, const label nei, const label masterFaceID, const bool flipFaceFlux, const label patchID)
Add face to cells and return new face index.
label nCells() const
label nFaces() const
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:118
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronise values on all mesh faces.
Definition: syncTools.H:387
A class for handling words, derived from string.
Definition: word.H:62
Foam::tmp< Foam::VolField< Type > > volField(const Foam::fvMeshSubset &, const Foam::VolField< Type > &vf)
Wrapper to get hold of the field or the subsetted field.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const fvPatchList & patches
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:228
#define WarningInFunction
Report a warning using Foam::Warning.
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
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars.
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
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:48
messageStream Info
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
error FatalError
static const char nl
Definition: Ostream.H:266
labelList f(nPoints)
objects
dictionary dict
Foam::argList args(argc, argv)