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-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  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 "polyModifyFace.H"
58 #include "polyAddFace.H"
59 #include "regionSplit.H"
60 #include "Tuple2.H"
61 #include "cyclicFvPatch.H"
62 
63 using namespace Foam;
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 void modifyOrAddFace
68 (
69  polyTopoChange& meshMod,
70  const face& f,
71  const label facei,
72  const label own,
73  const bool flipFaceFlux,
74  const label newPatchi,
75  const label zoneID,
76  const bool zoneFlip,
77 
78  PackedBoolList& modifiedFace
79 )
80 {
81  if (!modifiedFace[facei])
82  {
83  // First usage of face. Modify.
84  meshMod.setAction
85  (
87  (
88  f, // modified face
89  facei, // label of face
90  own, // owner
91  -1, // neighbour
92  flipFaceFlux, // face flip
93  newPatchi, // patch for face
94  false, // remove from zone
95  zoneID, // zone for face
96  zoneFlip // face flip in zone
97  )
98  );
99  modifiedFace[facei] = 1;
100  }
101  else
102  {
103  // Second or more usage of face. Add.
104  meshMod.setAction
105  (
107  (
108  f, // modified face
109  own, // owner
110  -1, // neighbour
111  -1, // master point
112  -1, // master edge
113  facei, // master face
114  flipFaceFlux, // face flip
115  newPatchi, // patch for face
116  zoneID, // zone for face
117  zoneFlip // face flip in zone
118  )
119  );
120  }
121 }
122 
123 
124 template<class Type>
125 void subsetVolFields
126 (
127  const fvMeshSubset& subsetter,
128  const IOobjectList& objectsList,
129  const label patchi,
130  const Type& exposedValue,
131  const word GeomVolType,
132  PtrList<VolField<Type>>& subFields
133 )
134 {
135  const fvMesh& baseMesh = subsetter.baseMesh();
136 
137  label i = 0;
138 
139  forAllConstIter(IOobjectList , objectsList, iter)
140  {
141  if (iter()->headerClassName() == GeomVolType)
142  {
143  const word fieldName = iter()->name();
144 
145  Info<< "Subsetting field " << fieldName << endl;
146 
148  (
149  *iter(),
150  baseMesh
151  );
152 
153  subFields.set(i, subsetter.interpolate(volField));
154 
155  // Explicitly set exposed faces (in patchi) to exposedValue.
156  if (patchi >= 0)
157  {
159  subFields[i++].boundaryFieldRef()[patchi];
160 
161  label newStart = fld.patch().patch().start();
162 
163  label oldPatchi = subsetter.patchMap()[patchi];
164 
165  if (oldPatchi == -1)
166  {
167  // New patch. Reset whole value.
168  fld = exposedValue;
169  }
170  else
171  {
172  // Reset those faces that originate from different patch
173  // or internal faces.
174  label oldSize = volField.boundaryField()[oldPatchi].size();
175  label oldStart = volField.boundaryField()
176  [
177  oldPatchi
178  ].patch().patch().start();
179 
180  forAll(fld, j)
181  {
182  label oldFacei = subsetter.faceMap()[newStart+j];
183 
184  if (oldFacei < oldStart || oldFacei >= oldStart+oldSize)
185  {
186  fld[j] = exposedValue;
187  }
188  }
189  }
190  }
191  }
192  }
193 }
194 
195 
196 template<class Type>
197 void subsetSurfaceFields
198 (
199  const fvMeshSubset& subsetter,
200  const IOobjectList& objectsList,
201  const label patchi,
202  const Type& exposedValue,
203  const word GeomSurfType,
204  PtrList<SurfaceField<Type>>& subFields
205 )
206 {
207  const fvMesh& baseMesh = subsetter.baseMesh();
208 
209  label i(0);
210 
211  forAllConstIter(IOobjectList , objectsList, iter)
212  {
213  if (iter()->headerClassName() == GeomSurfType)
214  {
215  const word& fieldName = iter.key();
216 
217  Info<< "Subsetting field " << fieldName << endl;
218 
220  (
221  *iter(),
222  baseMesh
223  );
224 
225  subFields.set(i, subsetter.interpolate(volField));
226 
227 
228  // Explicitly set exposed faces (in patchi) to exposedValue.
229  if (patchi >= 0)
230  {
232  subFields[i++].boundaryFieldRef()[patchi];
233 
234  label newStart = fld.patch().patch().start();
235 
236  label oldPatchi = subsetter.patchMap()[patchi];
237 
238  if (oldPatchi == -1)
239  {
240  // New patch. Reset whole value.
241  fld = exposedValue;
242  }
243  else
244  {
245  // Reset those faces that originate from different patch
246  // or internal faces.
247  label oldSize = volField.boundaryField()[oldPatchi].size();
248  label oldStart = volField.boundaryField()
249  [
250  oldPatchi
251  ].patch().patch().start();
252 
253  forAll(fld, j)
254  {
255  label oldFacei = subsetter.faceMap()[newStart+j];
256 
257  if (oldFacei < oldStart || oldFacei >= oldStart+oldSize)
258  {
259  fld[j] = exposedValue;
260  }
261  }
262  }
263  }
264  }
265  }
266 }
267 
268 
269 // Faces introduced into zero-sized patches don't get a value at all.
270 // This is hack to set them to an initial value.
271 template<class GeoField>
272 void initCreatedPatches
273 (
274  const fvMesh& mesh,
275  const polyTopoChangeMap& map,
276  const typename GeoField::value_type initValue
277 )
278 {
280  (
281  mesh.objectRegistry::lookupClass<GeoField>()
282  );
283 
284  for
285  (
286  typename HashTable<const GeoField*>::
287  iterator fieldIter = fields.begin();
288  fieldIter != fields.end();
289  ++fieldIter
290  )
291  {
292  GeoField& field = const_cast<GeoField&>(*fieldIter());
293 
294  typename GeoField::Boundary& fieldBf =
295  field.boundaryFieldRef();
296 
297  forAll(fieldBf, patchi)
298  {
299  if (map.oldPatchSizes()[patchi] == 0)
300  {
301  // Not mapped.
302  fieldBf[patchi] = initValue;
303 
304  if (fieldBf[patchi].fixesValue())
305  {
306  fieldBf[patchi] == initValue;
307  }
308  }
309  }
310  }
311 }
312 
313 
314 void createCoupledBaffles
315 (
316  fvMesh& mesh,
317  const labelList& coupledWantedPatch,
318  polyTopoChange& meshMod,
319  PackedBoolList& modifiedFace
320 )
321 {
322  const meshFaceZones& faceZones = mesh.faceZones();
323 
324  forAll(coupledWantedPatch, facei)
325  {
326  if (coupledWantedPatch[facei] != -1)
327  {
328  const face& f = mesh.faces()[facei];
329  label zoneID = faceZones.whichZone(facei);
330  bool zoneFlip = false;
331 
332  if (zoneID >= 0)
333  {
334  const faceZone& fZone = faceZones[zoneID];
335  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
336  }
337 
338  // Use owner side of face
339  modifyOrAddFace
340  (
341  meshMod,
342  f, // modified face
343  facei, // label of face
344  mesh.faceOwner()[facei], // owner
345  false, // face flip
346  coupledWantedPatch[facei], // patch for face
347  zoneID, // zone for face
348  zoneFlip, // face flip in zone
349  modifiedFace // modify or add status
350  );
351 
352  if (mesh.isInternalFace(facei))
353  {
354  label zoneID = faceZones.whichZone(facei);
355  bool zoneFlip = false;
356 
357  if (zoneID >= 0)
358  {
359  const faceZone& fZone = faceZones[zoneID];
360  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
361  }
362  // Use neighbour side of face
363  modifyOrAddFace
364  (
365  meshMod,
366  f.reverseFace(), // modified face
367  facei, // label of face
368  mesh.faceNeighbour()[facei],// owner
369  false, // face flip
370  coupledWantedPatch[facei], // patch for face
371  zoneID, // zone for face
372  zoneFlip, // face flip in zone
373  modifiedFace // modify or add status
374  );
375  }
376  }
377  }
378 }
379 
380 
381 void createCyclicCoupledBaffles
382 (
383  fvMesh& mesh,
384  const labelList& cyclicMasterPatch,
385  const labelList& cyclicSlavePatch,
386  polyTopoChange& meshMod,
387  PackedBoolList& modifiedFace
388 )
389 {
390  const meshFaceZones& faceZones = mesh.faceZones();
391 
392  forAll(cyclicMasterPatch, facei)
393  {
394  if (cyclicMasterPatch[facei] != -1)
395  {
396  const face& f = mesh.faces()[facei];
397 
398  label zoneID = faceZones.whichZone(facei);
399  bool zoneFlip = false;
400 
401  if (zoneID >= 0)
402  {
403  const faceZone& fZone = faceZones[zoneID];
404  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
405  }
406 
407  modifyOrAddFace
408  (
409  meshMod,
410  f.reverseFace(), // modified face
411  facei, // label of face
412  mesh.faceNeighbour()[facei], // owner
413  false, // face flip
414  cyclicMasterPatch[facei], // patch for face
415  zoneID, // zone for face
416  zoneFlip, // face flip in zone
417  modifiedFace // modify or add
418  );
419  }
420  }
421 
422  forAll(cyclicSlavePatch, facei)
423  {
424  if (cyclicSlavePatch[facei] != -1)
425  {
426  const face& f = mesh.faces()[facei];
427  if (mesh.isInternalFace(facei))
428  {
429  label zoneID = faceZones.whichZone(facei);
430  bool zoneFlip = false;
431 
432  if (zoneID >= 0)
433  {
434  const faceZone& fZone = faceZones[zoneID];
435  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
436  }
437  // Use owner side of face
438  modifyOrAddFace
439  (
440  meshMod,
441  f, // modified face
442  facei, // label of face
443  mesh.faceOwner()[facei], // owner
444  false, // face flip
445  cyclicSlavePatch[facei], // patch for face
446  zoneID, // zone for face
447  zoneFlip, // face flip in zone
448  modifiedFace // modify or add status
449  );
450  }
451  }
452  }
453 }
454 
455 
456 void createBaffles
457 (
458  fvMesh& mesh,
459  const labelList& wantedPatch,
460  polyTopoChange& meshMod
461 )
462 {
463  const meshFaceZones& faceZones = mesh.faceZones();
464  Info << "faceZone:createBaffle " << faceZones << endl;
465  forAll(wantedPatch, facei)
466  {
467  if (wantedPatch[facei] != -1)
468  {
469  const face& f = mesh.faces()[facei];
470 
471  label zoneID = faceZones.whichZone(facei);
472  bool zoneFlip = false;
473 
474  if (zoneID >= 0)
475  {
476  const faceZone& fZone = faceZones[zoneID];
477  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
478  }
479 
480  meshMod.setAction
481  (
483  (
484  f, // modified face
485  facei, // label of face
486  mesh.faceOwner()[facei], // owner
487  -1, // neighbour
488  false, // face flip
489  wantedPatch[facei], // patch for face
490  false, // remove from zone
491  zoneID, // zone for face
492  zoneFlip // face flip in zone
493  )
494  );
495 
496  if (mesh.isInternalFace(facei))
497  {
498  label zoneID = faceZones.whichZone(facei);
499  bool zoneFlip = false;
500 
501  if (zoneID >= 0)
502  {
503  const faceZone& fZone = faceZones[zoneID];
504  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
505  }
506 
507  meshMod.setAction
508  (
510  (
511  f.reverseFace(), // modified face
512  mesh.faceNeighbour()[facei],// owner
513  -1, // neighbour
514  -1, // masterPointID
515  -1, // masterEdgeID
516  facei, // masterFaceID,
517  false, // face flip
518  wantedPatch[facei], // patch for face
519  zoneID, // zone for face
520  zoneFlip // face flip in zone
521  )
522  );
523  }
524  }
525  }
526 }
527 
528 
529 // Wrapper around find patch. Also makes sure same patch in parallel.
530 label findPatch(const polyBoundaryMesh& patches, const word& patchName)
531 {
532  label patchi = patches.findPatchID(patchName);
533 
534  if (patchi == -1)
535  {
537  << "Illegal patch " << patchName
538  << nl << "Valid patches are " << patches.names()
539  << exit(FatalError);
540  }
541 
542  // Check same patch for all procs
543  {
544  label newPatch = patchi;
545  reduce(newPatch, minOp<label>());
546 
547  if (newPatch != patchi)
548  {
550  << "Patch " << patchName
551  << " should have the same patch index on all processors." << nl
552  << "On my processor it has index " << patchi
553  << " ; on some other processor it has index " << newPatch
554  << exit(FatalError);
555  }
556  }
557  return patchi;
558 }
559 
560 
561 
562 int main(int argc, char *argv[])
563 {
564  #include "addOverwriteOption.H"
565  #include "setRootCase.H"
567  #include "createMeshNoChangers.H"
568 
569  // Read control dictionary
570  // ~~~~~~~~~~~~~~~~~~~~~~~
571 
573  (
574  IOobject
575  (
576  "PDRMeshDict",
577  runTime.system(),
578  mesh,
581  )
582  );
583 
584  // Per faceSet the patch to put the baffles into
585  const List<Pair<word>> setsAndPatches(dict.lookup("blockedFaces"));
586 
587  // Per faceSet the patch to put the coupled baffles into
588  DynamicList<FixedList<word, 3>> coupledAndPatches(10);
589  const dictionary& functionDicts = dict.subDict("coupledFaces");
590  forAllConstIter(dictionary, functionDicts, iter)
591  {
592  // safety:
593  if (!iter().isDict())
594  {
595  continue;
596  }
597  const word& key = iter().keyword();
598 
599  const dictionary& dict = iter().dict();
600  const word cyclicName = dict.lookup("cyclicMasterPatchName");
601  const word wallName = dict.lookup("wallPatchName");
602  FixedList<word, 3> nameAndType;
603  nameAndType[0] = key;
604  nameAndType[1] = wallName;
605  nameAndType[2] = cyclicName;
606  coupledAndPatches.append(nameAndType);
607  }
608 
609  forAll(setsAndPatches, setI)
610  {
611  Info<< "Faces in faceSet " << setsAndPatches[setI][0]
612  << " become baffles in patch " << setsAndPatches[setI][1]
613  << endl;
614  }
615 
616  forAll(coupledAndPatches, setI)
617  {
618  Info<< "Faces in faceSet " << coupledAndPatches[setI][0]
619  << " become coupled baffles in patch " << coupledAndPatches[setI][1]
620  << endl;
621  }
622 
623  // All exposed faces that are not explicitly marked to be put into a patch
624  const word defaultPatch(dict.lookup("defaultPatch"));
625 
626  Info<< "Faces that get exposed become boundary faces in patch "
627  << defaultPatch << endl;
628 
629  const word blockedSetName(dict.lookup("blockedCells"));
630 
631  Info<< "Reading blocked cells from cellSet " << blockedSetName
632  << endl;
633 
634  const bool overwrite = args.optionFound("overwrite");
635 
636 
637  // Read faceSets, lookup patches
638  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
639 
640  // Check that face sets don't have coincident faces
641  labelList wantedPatch(mesh.nFaces(), -1);
642  forAll(setsAndPatches, setI)
643  {
644  faceSet fSet(mesh, setsAndPatches[setI][0]);
645 
646  label patchi = findPatch
647  (
648  mesh.boundaryMesh(),
649  setsAndPatches[setI][1]
650  );
651 
652  forAllConstIter(faceSet, fSet, iter)
653  {
654  if (wantedPatch[iter.key()] != -1)
655  {
657  << "Face " << iter.key()
658  << " is in faceSet " << setsAndPatches[setI][0]
659  << " destined for patch " << setsAndPatches[setI][1]
660  << " but also in patch " << wantedPatch[iter.key()]
661  << exit(FatalError);
662  }
663  wantedPatch[iter.key()] = patchi;
664  }
665  }
666 
667  // Per face the patch for coupled baffle or -1.
668  labelList coupledWantedPatch(mesh.nFaces(), -1);
669  labelList cyclicWantedPatch_half0(mesh.nFaces(), -1);
670  labelList cyclicWantedPatch_half1(mesh.nFaces(), -1);
671 
672  forAll(coupledAndPatches, setI)
673  {
674  const polyBoundaryMesh& patches = mesh.boundaryMesh();
675  const label cyclicId =
676  findPatch(patches, coupledAndPatches[setI][2]);
677 
678  const label cyclicSlaveId = findPatch
679  (
680  patches,
681  refCast<const cyclicFvPatch>
682  (
683  mesh.boundary()[cyclicId]
684  ).neighbFvPatch().name()
685  );
686 
687  faceSet fSet(mesh, coupledAndPatches[setI][0]);
688  label patchi = findPatch(patches, coupledAndPatches[setI][1]);
689 
690  forAllConstIter(faceSet, fSet, iter)
691  {
692  if (coupledWantedPatch[iter.key()] != -1)
693  {
695  << "Face " << iter.key()
696  << " is in faceSet " << coupledAndPatches[setI][0]
697  << " destined for patch " << coupledAndPatches[setI][1]
698  << " but also in patch " << coupledWantedPatch[iter.key()]
699  << exit(FatalError);
700  }
701  coupledWantedPatch[iter.key()] = patchi;
702  cyclicWantedPatch_half0[iter.key()] = cyclicId;
703  cyclicWantedPatch_half1[iter.key()] = cyclicSlaveId;
704  }
705  }
706 
707  // Exposed faces patch
708  label defaultPatchi = findPatch(mesh.boundaryMesh(), defaultPatch);
709 
710 
711  //
712  // Removing blockedCells (like subsetMesh)
713  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
714  //
715 
716  // Create mesh subsetting engine
717  fvMeshSubset subsetter(mesh);
718 
719  {
720 
721  cellSet blockedCells(mesh, blockedSetName);
722 
723  // invert
724  blockedCells.invert(mesh.nCells());
725 
726  // Create subsetted mesh.
727  subsetter.setLargeCellSubset(blockedCells, defaultPatchi, true);
728  }
729 
730 
731  // Subset wantedPatch. Note that might also include boundary faces
732  // that have been exposed by subsetting.
733  wantedPatch = IndirectList<label>(wantedPatch, subsetter.faceMap())();
734 
735  coupledWantedPatch = IndirectList<label>
736  (
737  coupledWantedPatch,
738  subsetter.faceMap()
739  )();
740 
741  cyclicWantedPatch_half0 = IndirectList<label>
742  (
743  cyclicWantedPatch_half0,
744  subsetter.faceMap()
745  )();
746 
747  cyclicWantedPatch_half1 = IndirectList<label>
748  (
749  cyclicWantedPatch_half1,
750  subsetter.faceMap()
751  )();
752 
753  // Read all fields in time and constant directories
754  IOobjectList objects(mesh, runTime.name());
755  IOobjectList timeObjects(IOobjectList(mesh, mesh.facesInstance()));
756  forAllConstIter(IOobjectList, timeObjects, iter)
757  {
758  if
759  (
760  iter()->headerClassName() == volScalarField::typeName
761  || iter()->headerClassName() == volVectorField::typeName
762  || iter()->headerClassName() == volSphericalTensorField::typeName
763  || iter()->headerClassName() == volTensorField::typeName
764  || iter()->headerClassName() == volSymmTensorField::typeName
765  || iter()->headerClassName() == surfaceScalarField::typeName
766  || iter()->headerClassName() == surfaceVectorField::typeName
767  || iter()->headerClassName()
769  || iter()->headerClassName() == surfaceSymmTensorField::typeName
770  || iter()->headerClassName() == surfaceTensorField::typeName
771  )
772  {
773  objects.add(*iter());
774  }
775  }
776 
777  // Read vol fields and subset.
778 
779  wordList scalarNames(objects.names(volScalarField::typeName));
780  PtrList<volScalarField> scalarFlds(scalarNames.size());
781  subsetVolFields
782  (
783  subsetter,
784  objects,
785  defaultPatchi,
786  scalar(Zero),
788  scalarFlds
789  );
790 
791  wordList vectorNames(objects.names(volVectorField::typeName));
792  PtrList<volVectorField> vectorFlds(vectorNames.size());
793  subsetVolFields
794  (
795  subsetter,
796  objects,
797  defaultPatchi,
798  vector(Zero),
800  vectorFlds
801  );
802 
803  wordList sphericalTensorNames
804  (
806  );
807  PtrList<volSphericalTensorField> sphericalTensorFlds
808  (
809  sphericalTensorNames.size()
810  );
811  subsetVolFields
812  (
813  subsetter,
814  objects,
815  defaultPatchi,
818  sphericalTensorFlds
819  );
820 
821  wordList symmTensorNames(objects.names(volSymmTensorField::typeName));
822  PtrList<volSymmTensorField> symmTensorFlds(symmTensorNames.size());
823  subsetVolFields
824  (
825  subsetter,
826  objects,
827  defaultPatchi,
828  symmTensor(Zero),
830  symmTensorFlds
831  );
832 
833  wordList tensorNames(objects.names(volTensorField::typeName));
834  PtrList<volTensorField> tensorFlds(tensorNames.size());
835  subsetVolFields
836  (
837  subsetter,
838  objects,
839  defaultPatchi,
840  tensor(Zero),
842  tensorFlds
843  );
844 
845  // Read surface fields and subset.
846 
847  wordList surfScalarNames(objects.names(surfaceScalarField::typeName));
848  PtrList<surfaceScalarField> surfScalarFlds(surfScalarNames.size());
849  subsetSurfaceFields
850  (
851  subsetter,
852  objects,
853  defaultPatchi,
854  scalar(Zero),
856  surfScalarFlds
857  );
858 
859  wordList surfVectorNames(objects.names(surfaceVectorField::typeName));
860  PtrList<surfaceVectorField> surfVectorFlds(surfVectorNames.size());
861  subsetSurfaceFields
862  (
863  subsetter,
864  objects,
865  defaultPatchi,
866  vector(Zero),
868  surfVectorFlds
869  );
870 
871  wordList surfSphericalTensorNames
872  (
874  );
875  PtrList<surfaceSphericalTensorField> surfSphericalTensorFlds
876  (
877  surfSphericalTensorNames.size()
878  );
879  subsetSurfaceFields
880  (
881  subsetter,
882  objects,
883  defaultPatchi,
886  surfSphericalTensorFlds
887  );
888 
889  wordList surfSymmTensorNames
890  (
892  );
893 
894  PtrList<surfaceSymmTensorField> surfSymmTensorFlds
895  (
896  surfSymmTensorNames.size()
897  );
898 
899  subsetSurfaceFields
900  (
901  subsetter,
902  objects,
903  defaultPatchi,
904  symmTensor(Zero),
906  surfSymmTensorFlds
907  );
908 
909  wordList surfTensorNames(objects.names(surfaceTensorField::typeName));
910  PtrList<surfaceTensorField> surfTensorFlds(surfTensorNames.size());
911  subsetSurfaceFields
912  (
913  subsetter,
914  objects,
915  defaultPatchi,
916  tensor(Zero),
918  surfTensorFlds
919  );
920 
921  if (!overwrite)
922  {
923  runTime++;
924  }
925 
926  Info<< "Writing mesh without blockedCells to time "
927  << runTime.userTimeName() << endl;
928 
929  // Subsetting adds 'subset' prefix. Rename field to be like original.
930  forAll(scalarFlds, i)
931  {
932  scalarFlds[i].rename(scalarNames[i]);
933  scalarFlds[i].writeOpt() = IOobject::AUTO_WRITE;
934  scalarFlds[i].checkIn();
935  }
936  forAll(vectorFlds, i)
937  {
938  vectorFlds[i].rename(vectorNames[i]);
939  vectorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
940  vectorFlds[i].checkIn();
941  }
942  forAll(sphericalTensorFlds, i)
943  {
944  sphericalTensorFlds[i].rename(sphericalTensorNames[i]);
945  sphericalTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
946  sphericalTensorFlds[i].checkIn();
947  }
948  forAll(symmTensorFlds, i)
949  {
950  symmTensorFlds[i].rename(symmTensorNames[i]);
951  symmTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
952  symmTensorFlds[i].checkIn();
953  }
954  forAll(tensorFlds, i)
955  {
956  tensorFlds[i].rename(tensorNames[i]);
957  tensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
958  tensorFlds[i].checkIn();
959  }
960 
961  // Surface ones.
962  forAll(surfScalarFlds, i)
963  {
964  surfScalarFlds[i].rename(surfScalarNames[i]);
965  surfScalarFlds[i].writeOpt() = IOobject::AUTO_WRITE;
966  surfScalarFlds[i].checkIn();
967  }
968  forAll(surfVectorFlds, i)
969  {
970  surfVectorFlds[i].rename(surfVectorNames[i]);
971  surfVectorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
972  surfVectorFlds[i].checkIn();
973  }
974  forAll(surfSphericalTensorFlds, i)
975  {
976  surfSphericalTensorFlds[i].rename(surfSphericalTensorNames[i]);
977  surfSphericalTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
978  surfSphericalTensorFlds[i].checkIn();
979  }
980  forAll(surfSymmTensorFlds, i)
981  {
982  surfSymmTensorFlds[i].rename(surfSymmTensorNames[i]);
983  surfSymmTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
984  surfSymmTensorFlds[i].checkIn();
985  }
986  forAll(surfTensorNames, i)
987  {
988  surfTensorFlds[i].rename(surfTensorNames[i]);
989  surfTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
990  surfTensorFlds[i].checkIn();
991  }
992 
993  subsetter.subMesh().write();
994 
995 
996  //
997  // Splitting blockedFaces
998  // ~~~~~~~~~~~~~~~~~~~~~~
999  //
1000 
1001  // Synchronise wantedPatch across coupled patches.
1003  (
1004  subsetter.subMesh(),
1005  wantedPatch,
1006  maxEqOp<label>()
1007  );
1008 
1009  // Synchronise coupledWantedPatch across coupled patches.
1011  (
1012  subsetter.subMesh(),
1013  coupledWantedPatch,
1014  maxEqOp<label>()
1015  );
1016 
1017  // Synchronise cyclicWantedPatch across coupled patches.
1019  (
1020  subsetter.subMesh(),
1021  cyclicWantedPatch_half0,
1022  maxEqOp<label>()
1023  );
1024 
1025  // Synchronise cyclicWantedPatch across coupled patches.
1027  (
1028  subsetter.subMesh(),
1029  cyclicWantedPatch_half1,
1030  maxEqOp<label>()
1031  );
1032 
1033  // Topochange container
1034  polyTopoChange meshMod(subsetter.subMesh());
1035 
1036 
1037  // Whether first use of face (modify) or consecutive (add)
1038  PackedBoolList modifiedFace(mesh.nFaces());
1039 
1040  // Create coupled wall-side baffles
1041  createCoupledBaffles
1042  (
1043  subsetter.subMesh(),
1044  coupledWantedPatch,
1045  meshMod,
1046  modifiedFace
1047  );
1048 
1049  // Create coupled master/slave cyclic baffles
1050  createCyclicCoupledBaffles
1051  (
1052  subsetter.subMesh(),
1053  cyclicWantedPatch_half0,
1054  cyclicWantedPatch_half1,
1055  meshMod,
1056  modifiedFace
1057  );
1058 
1059  // Create wall baffles
1060  createBaffles
1061  (
1062  subsetter.subMesh(),
1063  wantedPatch,
1064  meshMod
1065  );
1066 
1067  if (!overwrite)
1068  {
1069  runTime++;
1070  }
1071 
1072  // Change the mesh. Change points directly (no inflation).
1074  meshMod.changeMesh(subsetter.subMesh(), false);
1075 
1076  // Update fields
1077  subsetter.subMesh().topoChange(map);
1078 
1079  // Fix faces that get mapped to zero-sized patches (these don't get any
1080  // value)
1081  initCreatedPatches<volScalarField>
1082  (
1083  subsetter.subMesh(),
1084  map,
1085  0.0
1086  );
1087  initCreatedPatches<volVectorField>
1088  (
1089  subsetter.subMesh(),
1090  map,
1091  Zero
1092  );
1093  initCreatedPatches<volSphericalTensorField>
1094  (
1095  subsetter.subMesh(),
1096  map,
1097  Zero
1098  );
1099  initCreatedPatches<volSymmTensorField>
1100  (
1101  subsetter.subMesh(),
1102  map,
1103  Zero
1104  );
1105  initCreatedPatches<volTensorField>
1106  (
1107  subsetter.subMesh(),
1108  map,
1109  Zero
1110  );
1111 
1112  initCreatedPatches<surfaceScalarField>
1113  (
1114  subsetter.subMesh(),
1115  map,
1116  0.0
1117  );
1118  initCreatedPatches<surfaceVectorField>
1119  (
1120  subsetter.subMesh(),
1121  map,
1122  Zero
1123  );
1124  initCreatedPatches<surfaceSphericalTensorField>
1125  (
1126  subsetter.subMesh(),
1127  map,
1128  Zero
1129  );
1130  initCreatedPatches<surfaceSymmTensorField>
1131  (
1132  subsetter.subMesh(),
1133  map,
1134  Zero
1135  );
1136  initCreatedPatches<surfaceTensorField>
1137  (
1138  subsetter.subMesh(),
1139  map,
1140  Zero
1141  );
1142 
1143 
1144  // Move mesh (since morphing might not do this)
1145  if (map().hasMotionPoints())
1146  {
1147  subsetter.subMesh().setPoints(map().preMotionPoints());
1148  }
1149 
1150  Info<< "Writing mesh with split blockedFaces to time "
1151  << runTime.userTimeName() << endl;
1152 
1153  subsetter.subMesh().write();
1154 
1155 
1156  //
1157  // Removing inaccessible regions
1158  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1159  //
1160 
1161  // Determine connected regions. regionSplit is the labelList with the
1162  // region per cell.
1163  regionSplit cellRegion(subsetter.subMesh());
1164 
1165  if (cellRegion.nRegions() > 1)
1166  {
1168  << "Removing blocked faces and cells created "
1169  << cellRegion.nRegions()
1170  << " regions that are not connected via a face." << nl
1171  << " This is not supported in solvers." << nl
1172  << " Use" << nl << nl
1173  << " splitMeshRegions <root> <case> -largestOnly" << nl << nl
1174  << " to extract a single region of the mesh." << nl
1175  << " This mesh will be written to a new timedirectory"
1176  << " so might have to be moved back to constant/" << nl
1177  << endl;
1178 
1179  word startFrom(runTime.controlDict().lookup("startFrom"));
1180 
1181  if (startFrom != "latestTime")
1182  {
1184  << "To run splitMeshRegions please set your"
1185  << " startFrom entry to latestTime" << endl;
1186  }
1187  }
1188 
1189  Info << nl << "End" << endl;
1190 
1191  return 0;
1192 }
1193 
1194 
1195 // ************************************************************************* //
#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:105
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
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: MeshZones.C:221
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:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
A list of face labels.
Definition: faceSet.H:51
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:68
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:307
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:896
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:101
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:895
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1736
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1236
virtual void setPoints(const pointField &)
Reset the points.
Definition: fvMesh.C:1133
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:54
Foam::polyBoundaryMesh.
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:1025
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:447
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1373
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1386
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1392
Class describing modification of a face.
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 inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
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:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const fvPatchList & patches
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
#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:251
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:260
labelList f(nPoints)
objects
dictionary dict
Foam::argList args(argc, argv)