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 | Copyright (C) 2011-2015 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 mixed-coupled
38  duplicate baffles
39 
40  Subsets out the blocked cells and splits the blockedFaces and updates
41  fields.
42 
43  The face splitting is done by duplicating the faces. No duplication of
44  points for now (so checkMesh will show a lot of 'duplicate face' messages)
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #include "fvMeshSubset.H"
49 #include "argList.H"
50 #include "cellSet.H"
51 #include "IOobjectList.H"
52 #include "volFields.H"
53 #include "mapPolyMesh.H"
54 #include "faceSet.H"
55 #include "cellSet.H"
56 #include "syncTools.H"
57 #include "polyTopoChange.H"
58 #include "polyModifyFace.H"
59 #include "polyAddFace.H"
60 #include "regionSplit.H"
61 #include "Tuple2.H"
62 #include "cyclicFvPatch.H"
63 
64 using namespace Foam;
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 void modifyOrAddFace
69 (
70  polyTopoChange& meshMod,
71  const face& f,
72  const label faceI,
73  const label own,
74  const bool flipFaceFlux,
75  const label newPatchI,
76  const label zoneID,
77  const bool zoneFlip,
78 
79  PackedBoolList& modifiedFace
80 )
81 {
82  if (!modifiedFace[faceI])
83  {
84  // First usage of face. Modify.
85  meshMod.setAction
86  (
88  (
89  f, // modified face
90  faceI, // label of face
91  own, // owner
92  -1, // neighbour
93  flipFaceFlux, // face flip
94  newPatchI, // patch for face
95  false, // remove from zone
96  zoneID, // zone for face
97  zoneFlip // face flip in zone
98  )
99  );
100  modifiedFace[faceI] = 1;
101  }
102  else
103  {
104  // Second or more usage of face. Add.
105  meshMod.setAction
106  (
108  (
109  f, // modified face
110  own, // owner
111  -1, // neighbour
112  -1, // master point
113  -1, // master edge
114  faceI, // master face
115  flipFaceFlux, // face flip
116  newPatchI, // patch for face
117  zoneID, // zone for face
118  zoneFlip // face flip in zone
119  )
120  );
121  }
122 }
123 
124 
125 template<class Type>
126 void subsetVolFields
127 (
128  const fvMeshSubset& subsetter,
129  const IOobjectList& objectsList,
130  const label patchI,
131  const Type& exposedValue,
132  const word GeomVolType,
134 )
135 {
136  const fvMesh& baseMesh = subsetter.baseMesh();
137 
138  label i = 0;
139 
140  forAllConstIter(IOobjectList , objectsList, iter)
141  {
142  if (iter()->headerClassName() == GeomVolType)
143  {
144  const word fieldName = iter()->name();
145 
146  Info<< "Subsetting field " << fieldName << endl;
147 
149  (
150  *iter(),
151  baseMesh
152  );
153 
154  subFields.set(i, subsetter.interpolate(volField));
155 
156  // Explicitly set exposed faces (in patchI) to exposedValue.
157  if (patchI >= 0)
158  {
160  subFields[i++].boundaryField()[patchI];
161 
162  label newStart = fld.patch().patch().start();
163 
164  label oldPatchI = subsetter.patchMap()[patchI];
165 
166  if (oldPatchI == -1)
167  {
168  // New patch. Reset whole value.
169  fld = exposedValue;
170  }
171  else
172  {
173  // Reset those faces that originate from different patch
174  // or internal faces.
175  label oldSize = volField.boundaryField()[oldPatchI].size();
176  label oldStart = volField.boundaryField()
177  [
178  oldPatchI
179  ].patch().patch().start();
180 
181  forAll(fld, j)
182  {
183  label oldFaceI = subsetter.faceMap()[newStart+j];
184 
185  if (oldFaceI < oldStart || oldFaceI >= oldStart+oldSize)
186  {
187  fld[j] = exposedValue;
188  }
189  }
190  }
191  }
192  }
193  }
194 }
195 
196 
197 template<class Type>
199 (
200  const fvMeshSubset& subsetter,
201  const IOobjectList& objectsList,
202  const label patchI,
203  const Type& exposedValue,
204  const word GeomSurfType,
206 )
207 {
208  const fvMesh& baseMesh = subsetter.baseMesh();
209 
210  label i(0);
211 
212  forAllConstIter(IOobjectList , objectsList, iter)
213  {
214  if (iter()->headerClassName() == GeomSurfType)
215  {
216  const word& fieldName = iter.key();
217 
218  Info<< "Subsetting field " << fieldName << endl;
219 
221  (
222  *iter(),
223  baseMesh
224  );
225 
226  subFields.set(i, subsetter.interpolate(volField));
227 
228 
229  // Explicitly set exposed faces (in patchI) to exposedValue.
230  if (patchI >= 0)
231  {
233  subFields[i++].boundaryField()[patchI];
234 
235  label newStart = fld.patch().patch().start();
236 
237  label oldPatchI = subsetter.patchMap()[patchI];
238 
239  if (oldPatchI == -1)
240  {
241  // New patch. Reset whole value.
242  fld = exposedValue;
243  }
244  else
245  {
246  // Reset those faces that originate from different patch
247  // or internal faces.
248  label oldSize = volField.boundaryField()[oldPatchI].size();
249  label oldStart = volField.boundaryField()
250  [
251  oldPatchI
252  ].patch().patch().start();
253 
254  forAll(fld, j)
255  {
256  label oldFaceI = subsetter.faceMap()[newStart+j];
257 
258  if (oldFaceI < oldStart || oldFaceI >= oldStart+oldSize)
259  {
260  fld[j] = exposedValue;
261  }
262  }
263  }
264  }
265  }
266  }
267 }
268 
269 
270 // Faces introduced into zero-sized patches don't get a value at all.
271 // This is hack to set them to an initial value.
272 template<class GeoField>
273 void initCreatedPatches
274 (
275  const fvMesh& mesh,
276  const mapPolyMesh& map,
277  const typename GeoField::value_type initValue
278 )
279 {
281  (
282  mesh.objectRegistry::lookupClass<GeoField>()
283  );
284 
285  for
286  (
287  typename HashTable<const GeoField*>::
288  iterator fieldIter = fields.begin();
289  fieldIter != fields.end();
290  ++fieldIter
291  )
292  {
293  GeoField& field = const_cast<GeoField&>(*fieldIter());
294 
295  forAll(field.boundaryField(), patchi)
296  {
297  if (map.oldPatchSizes()[patchi] == 0)
298  {
299  // Not mapped.
300  field.boundaryField()[patchi] = initValue;
301 
302  if (field.boundaryField()[patchi].fixesValue())
303  {
304  field.boundaryField()[patchi] == initValue;
305  }
306  }
307  }
308  }
309 }
310 
311 
312 void createCoupledBaffles
313 (
314  fvMesh& mesh,
315  const labelList& coupledWantedPatch,
316  polyTopoChange& meshMod,
317  PackedBoolList& modifiedFace
318 )
319 {
320  const faceZoneMesh& faceZones = mesh.faceZones();
321 
322  forAll(coupledWantedPatch, faceI)
323  {
324  if (coupledWantedPatch[faceI] != -1)
325  {
326  const face& f = mesh.faces()[faceI];
327  label zoneID = faceZones.whichZone(faceI);
328  bool zoneFlip = false;
329 
330  if (zoneID >= 0)
331  {
332  const faceZone& fZone = faceZones[zoneID];
333  zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
334  }
335 
336  // Use owner side of face
337  modifyOrAddFace
338  (
339  meshMod,
340  f, // modified face
341  faceI, // label of face
342  mesh.faceOwner()[faceI], // owner
343  false, // face flip
344  coupledWantedPatch[faceI], // patch for face
345  zoneID, // zone for face
346  zoneFlip, // face flip in zone
347  modifiedFace // modify or add status
348  );
349 
350  if (mesh.isInternalFace(faceI))
351  {
352  label zoneID = faceZones.whichZone(faceI);
353  bool zoneFlip = false;
354 
355  if (zoneID >= 0)
356  {
357  const faceZone& fZone = faceZones[zoneID];
358  zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
359  }
360  // Use neighbour side of face
361  modifyOrAddFace
362  (
363  meshMod,
364  f.reverseFace(), // modified face
365  faceI, // label of face
366  mesh.faceNeighbour()[faceI],// owner
367  false, // face flip
368  coupledWantedPatch[faceI], // patch for face
369  zoneID, // zone for face
370  zoneFlip, // face flip in zone
371  modifiedFace // modify or add status
372  );
373  }
374  }
375  }
376 }
377 
378 
379 void createCyclicCoupledBaffles
380 (
381  fvMesh& mesh,
382  const labelList& cyclicMasterPatch,
383  const labelList& cyclicSlavePatch,
384  polyTopoChange& meshMod,
385  PackedBoolList& modifiedFace
386 )
387 {
388  const faceZoneMesh& faceZones = mesh.faceZones();
389 
390  forAll(cyclicMasterPatch, faceI)
391  {
392  if (cyclicMasterPatch[faceI] != -1)
393  {
394  const face& f = mesh.faces()[faceI];
395 
396  label zoneID = faceZones.whichZone(faceI);
397  bool zoneFlip = false;
398 
399  if (zoneID >= 0)
400  {
401  const faceZone& fZone = faceZones[zoneID];
402  zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
403  }
404 
405  modifyOrAddFace
406  (
407  meshMod,
408  f.reverseFace(), // modified face
409  faceI, // label of face
410  mesh.faceNeighbour()[faceI], // owner
411  false, // face flip
412  cyclicMasterPatch[faceI], // patch for face
413  zoneID, // zone for face
414  zoneFlip, // face flip in zone
415  modifiedFace // modify or add
416  );
417  }
418  }
419 
420  forAll(cyclicSlavePatch, faceI)
421  {
422  if (cyclicSlavePatch[faceI] != -1)
423  {
424  const face& f = mesh.faces()[faceI];
425  if (mesh.isInternalFace(faceI))
426  {
427  label zoneID = faceZones.whichZone(faceI);
428  bool zoneFlip = false;
429 
430  if (zoneID >= 0)
431  {
432  const faceZone& fZone = faceZones[zoneID];
433  zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
434  }
435  // Use owner side of face
436  modifyOrAddFace
437  (
438  meshMod,
439  f, // modified face
440  faceI, // label of face
441  mesh.faceOwner()[faceI], // owner
442  false, // face flip
443  cyclicSlavePatch[faceI], // patch for face
444  zoneID, // zone for face
445  zoneFlip, // face flip in zone
446  modifiedFace // modify or add status
447  );
448  }
449  }
450  }
451 }
452 
453 
454 void createBaffles
455 (
456  fvMesh& mesh,
457  const labelList& wantedPatch,
458  polyTopoChange& meshMod
459 )
460 {
461  const faceZoneMesh& faceZones = mesh.faceZones();
462  Info << "faceZone:createBaffle " << faceZones << endl;
463  forAll(wantedPatch, faceI)
464  {
465  if (wantedPatch[faceI] != -1)
466  {
467  const face& f = mesh.faces()[faceI];
468 
469  label zoneID = faceZones.whichZone(faceI);
470  bool zoneFlip = false;
471 
472  if (zoneID >= 0)
473  {
474  const faceZone& fZone = faceZones[zoneID];
475  zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
476  }
477 
478  meshMod.setAction
479  (
481  (
482  f, // modified face
483  faceI, // label of face
484  mesh.faceOwner()[faceI], // owner
485  -1, // neighbour
486  false, // face flip
487  wantedPatch[faceI], // patch for face
488  false, // remove from zone
489  zoneID, // zone for face
490  zoneFlip // face flip in zone
491  )
492  );
493 
494  if (mesh.isInternalFace(faceI))
495  {
496  label zoneID = faceZones.whichZone(faceI);
497  bool zoneFlip = false;
498 
499  if (zoneID >= 0)
500  {
501  const faceZone& fZone = faceZones[zoneID];
502  zoneFlip = fZone.flipMap()[fZone.whichFace(faceI)];
503  }
504 
505  meshMod.setAction
506  (
508  (
509  f.reverseFace(), // modified face
510  mesh.faceNeighbour()[faceI],// owner
511  -1, // neighbour
512  -1, // masterPointID
513  -1, // masterEdgeID
514  faceI, // masterFaceID,
515  false, // face flip
516  wantedPatch[faceI], // patch for face
517  zoneID, // zone for face
518  zoneFlip // face flip in zone
519  )
520  );
521  }
522  }
523  }
524 }
525 
526 
527 // Wrapper around find patch. Also makes sure same patch in parallel.
528 label findPatch(const polyBoundaryMesh& patches, const word& patchName)
529 {
530  label patchI = patches.findPatchID(patchName);
531 
532  if (patchI == -1)
533  {
535  << "Illegal patch " << patchName
536  << nl << "Valid patches are " << patches.names()
537  << exit(FatalError);
538  }
539 
540  // Check same patch for all procs
541  {
542  label newPatch = patchI;
543  reduce(newPatch, minOp<label>());
544 
545  if (newPatch != patchI)
546  {
548  << "Patch " << patchName
549  << " should have the same patch index on all processors." << nl
550  << "On my processor it has index " << patchI
551  << " ; on some other processor it has index " << newPatch
552  << exit(FatalError);
553  }
554  }
555  return patchI;
556 }
557 
558 
559 
560 int main(int argc, char *argv[])
561 {
562  #include "addOverwriteOption.H"
563  #include "setRootCase.H"
564  #include "createTime.H"
565  runTime.functionObjects().off();
566  #include "createMesh.H"
567 
568  // Read control dictionary
569  // ~~~~~~~~~~~~~~~~~~~~~~~
570 
572  (
573  IOobject
574  (
575  "PDRMeshDict",
576  runTime.system(),
577  mesh,
580  )
581  );
582 
583  // Per faceSet the patch to put the baffles into
584  const List<Pair<word> > setsAndPatches(dict.lookup("blockedFaces"));
585 
586  // Per faceSet the patch to put the coupled baffles into
587  DynamicList<FixedList<word, 3> > coupledAndPatches(10);
588  const dictionary& functionDicts = dict.subDict("coupledFaces");
589  forAllConstIter(dictionary, functionDicts, iter)
590  {
591  // safety:
592  if (!iter().isDict())
593  {
594  continue;
595  }
596  const word& key = iter().keyword();
597 
598  const dictionary& dict = iter().dict();
599  const word cyclicName = dict.lookup("cyclicMasterPatchName");
600  const word wallName = dict.lookup("wallPatchName");
601  FixedList<word, 3> nameAndType;
602  nameAndType[0] = key;
603  nameAndType[1] = wallName;
604  nameAndType[2] = cyclicName;
605  coupledAndPatches.append(nameAndType);
606  }
607 
608  forAll(setsAndPatches, setI)
609  {
610  Info<< "Faces in faceSet " << setsAndPatches[setI][0]
611  << " become baffles in patch " << setsAndPatches[setI][1]
612  << endl;
613  }
614 
615  forAll(coupledAndPatches, setI)
616  {
617  Info<< "Faces in faceSet " << coupledAndPatches[setI][0]
618  << " become coupled baffles in patch " << coupledAndPatches[setI][1]
619  << endl;
620  }
621 
622  // All exposed faces that are not explicitly marked to be put into a patch
623  const word defaultPatch(dict.lookup("defaultPatch"));
624 
625  Info<< "Faces that get exposed become boundary faces in patch "
626  << defaultPatch << endl;
627 
628  const word blockedSetName(dict.lookup("blockedCells"));
629 
630  Info<< "Reading blocked cells from cellSet " << blockedSetName
631  << endl;
632 
633  const bool overwrite = args.optionFound("overwrite");
634 
635 
636  // Read faceSets, lookup patches
637  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
638 
639  // Check that face sets don't have coincident faces
640  labelList wantedPatch(mesh.nFaces(), -1);
641  forAll(setsAndPatches, setI)
642  {
643  faceSet fSet(mesh, setsAndPatches[setI][0]);
644 
645  label patchI = findPatch
646  (
647  mesh.boundaryMesh(),
648  setsAndPatches[setI][1]
649  );
650 
651  forAllConstIter(faceSet, fSet, iter)
652  {
653  if (wantedPatch[iter.key()] != -1)
654  {
656  << "Face " << iter.key()
657  << " is in faceSet " << setsAndPatches[setI][0]
658  << " destined for patch " << setsAndPatches[setI][1]
659  << " but also in patch " << wantedPatch[iter.key()]
660  << exit(FatalError);
661  }
662  wantedPatch[iter.key()] = patchI;
663  }
664  }
665 
666  // Per face the patch for coupled baffle or -1.
667  labelList coupledWantedPatch(mesh.nFaces(), -1);
668  labelList cyclicWantedPatch_half0(mesh.nFaces(), -1);
669  labelList cyclicWantedPatch_half1(mesh.nFaces(), -1);
670 
671  forAll(coupledAndPatches, setI)
672  {
674  const label cyclicId =
675  findPatch(patches, coupledAndPatches[setI][2]);
676 
677  const label cyclicSlaveId = findPatch
678  (
679  patches,
680  refCast<const cyclicFvPatch>
681  (
682  mesh.boundary()[cyclicId]
683  ).neighbFvPatch().name()
684  );
685 
686  faceSet fSet(mesh, coupledAndPatches[setI][0]);
687  label patchI = findPatch(patches, coupledAndPatches[setI][1]);
688 
689  forAllConstIter(faceSet, fSet, iter)
690  {
691  if (coupledWantedPatch[iter.key()] != -1)
692  {
694  << "Face " << iter.key()
695  << " is in faceSet " << coupledAndPatches[setI][0]
696  << " destined for patch " << coupledAndPatches[setI][1]
697  << " but also in patch " << coupledWantedPatch[iter.key()]
698  << exit(FatalError);
699  }
700  coupledWantedPatch[iter.key()] = patchI;
701  cyclicWantedPatch_half0[iter.key()] = cyclicId;
702  cyclicWantedPatch_half1[iter.key()] = cyclicSlaveId;
703  }
704  }
705 
706  // Exposed faces patch
707  label defaultPatchI = findPatch(mesh.boundaryMesh(), defaultPatch);
708 
709 
710  //
711  // Removing blockedCells (like subsetMesh)
712  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
713  //
714 
715  // Create mesh subsetting engine
716  fvMeshSubset subsetter(mesh);
717 
718  {
719 
720  cellSet blockedCells(mesh, blockedSetName);
721 
722  // invert
723  blockedCells.invert(mesh.nCells());
724 
725  // Create subsetted mesh.
726  subsetter.setLargeCellSubset(blockedCells, defaultPatchI, true);
727  }
728 
729 
730  // Subset wantedPatch. Note that might also include boundary faces
731  // that have been exposed by subsetting.
732  wantedPatch = IndirectList<label>(wantedPatch, subsetter.faceMap())();
733 
734  coupledWantedPatch = IndirectList<label>
735  (
736  coupledWantedPatch,
737  subsetter.faceMap()
738  )();
739 
740  cyclicWantedPatch_half0 = IndirectList<label>
741  (
742  cyclicWantedPatch_half0,
743  subsetter.faceMap()
744  )();
745 
746  cyclicWantedPatch_half1 = IndirectList<label>
747  (
748  cyclicWantedPatch_half1,
749  subsetter.faceMap()
750  )();
751 
752  // Read all fields in time and constant directories
753  IOobjectList objects(mesh, runTime.timeName());
755  forAllConstIter(IOobjectList, timeObjects, iter)
756  {
757  if
758  (
759  iter()->headerClassName() == volScalarField::typeName
760  || iter()->headerClassName() == volVectorField::typeName
761  || iter()->headerClassName() == volSphericalTensorField::typeName
762  || iter()->headerClassName() == volTensorField::typeName
763  || iter()->headerClassName() == volSymmTensorField::typeName
764  || iter()->headerClassName() == surfaceScalarField::typeName
765  || iter()->headerClassName() == surfaceVectorField::typeName
766  || iter()->headerClassName()
767  == surfaceSphericalTensorField::typeName
768  || iter()->headerClassName() == surfaceSymmTensorField::typeName
769  || iter()->headerClassName() == surfaceTensorField::typeName
770  )
771  {
772  objects.add(*iter());
773  }
774  }
775 
776  // Read vol fields and subset.
777 
778  wordList scalarNames(objects.names(volScalarField::typeName));
779  PtrList<volScalarField> scalarFlds(scalarNames.size());
781  (
782  subsetter,
783  objects,
784  defaultPatchI,
786  volScalarField::typeName,
787  scalarFlds
788  );
789 
790  wordList vectorNames(objects.names(volVectorField::typeName));
791  PtrList<volVectorField> vectorFlds(vectorNames.size());
793  (
794  subsetter,
795  objects,
796  defaultPatchI,
798  volVectorField::typeName,
799  vectorFlds
800  );
801 
802  wordList sphericalTensorNames
803  (
804  objects.names(volSphericalTensorField::typeName)
805  );
806  PtrList<volSphericalTensorField> sphericalTensorFlds
807  (
808  sphericalTensorNames.size()
809  );
811  (
812  subsetter,
813  objects,
814  defaultPatchI,
816  volSphericalTensorField::typeName,
817  sphericalTensorFlds
818  );
819 
820  wordList symmTensorNames(objects.names(volSymmTensorField::typeName));
821  PtrList<volSymmTensorField> symmTensorFlds(symmTensorNames.size());
823  (
824  subsetter,
825  objects,
826  defaultPatchI,
828  volSymmTensorField::typeName,
829  symmTensorFlds
830  );
831 
832  wordList tensorNames(objects.names(volTensorField::typeName));
833  PtrList<volTensorField> tensorFlds(tensorNames.size());
835  (
836  subsetter,
837  objects,
838  defaultPatchI,
840  volTensorField::typeName,
841  tensorFlds
842  );
843 
844  // Read surface fields and subset.
845 
846  wordList surfScalarNames(objects.names(surfaceScalarField::typeName));
847  PtrList<surfaceScalarField> surfScalarFlds(surfScalarNames.size());
849  (
850  subsetter,
851  objects,
852  defaultPatchI,
854  surfaceScalarField::typeName,
855  surfScalarFlds
856  );
857 
858  wordList surfVectorNames(objects.names(surfaceVectorField::typeName));
859  PtrList<surfaceVectorField> surfVectorFlds(surfVectorNames.size());
861  (
862  subsetter,
863  objects,
864  defaultPatchI,
866  surfaceVectorField::typeName,
867  surfVectorFlds
868  );
869 
870  wordList surfSphericalTensorNames
871  (
872  objects.names(surfaceSphericalTensorField::typeName)
873  );
874  PtrList<surfaceSphericalTensorField> surfSphericalTensorFlds
875  (
876  surfSphericalTensorNames.size()
877  );
879  (
880  subsetter,
881  objects,
882  defaultPatchI,
884  surfaceSphericalTensorField::typeName,
885  surfSphericalTensorFlds
886  );
887 
888  wordList surfSymmTensorNames
889  (
890  objects.names(surfaceSymmTensorField::typeName)
891  );
892 
893  PtrList<surfaceSymmTensorField> surfSymmTensorFlds
894  (
895  surfSymmTensorNames.size()
896  );
897 
899  (
900  subsetter,
901  objects,
902  defaultPatchI,
904  surfaceSymmTensorField::typeName,
905  surfSymmTensorFlds
906  );
907 
908  wordList surfTensorNames(objects.names(surfaceTensorField::typeName));
909  PtrList<surfaceTensorField> surfTensorFlds(surfTensorNames.size());
911  (
912  subsetter,
913  objects,
914  defaultPatchI,
916  surfaceTensorField::typeName,
917  surfTensorFlds
918  );
919 
920  if (!overwrite)
921  {
922  runTime++;
923  }
924 
925  Info<< "Writing mesh without blockedCells to time " << runTime.value()
926  << endl;
927 
928  // Subsetting adds 'subset' prefix. Rename field to be like original.
929  forAll(scalarFlds, i)
930  {
931  scalarFlds[i].rename(scalarNames[i]);
932  scalarFlds[i].writeOpt() = IOobject::AUTO_WRITE;
933  }
934  forAll(vectorFlds, i)
935  {
936  vectorFlds[i].rename(vectorNames[i]);
937  vectorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
938  }
939  forAll(sphericalTensorFlds, i)
940  {
941  sphericalTensorFlds[i].rename(sphericalTensorNames[i]);
942  sphericalTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
943  }
944  forAll(symmTensorFlds, i)
945  {
946  symmTensorFlds[i].rename(symmTensorNames[i]);
947  symmTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
948  }
949  forAll(tensorFlds, i)
950  {
951  tensorFlds[i].rename(tensorNames[i]);
952  tensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
953  }
954 
955  // Surface ones.
956  forAll(surfScalarFlds, i)
957  {
958  surfScalarFlds[i].rename(surfScalarNames[i]);
959  surfScalarFlds[i].writeOpt() = IOobject::AUTO_WRITE;
960  }
961  forAll(surfVectorFlds, i)
962  {
963  surfVectorFlds[i].rename(surfVectorNames[i]);
964  surfVectorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
965  }
966  forAll(surfSphericalTensorFlds, i)
967  {
968  surfSphericalTensorFlds[i].rename(surfSphericalTensorNames[i]);
969  surfSphericalTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
970  }
971  forAll(surfSymmTensorFlds, i)
972  {
973  surfSymmTensorFlds[i].rename(surfSymmTensorNames[i]);
974  surfSymmTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
975  }
976  forAll(surfTensorNames, i)
977  {
978  surfTensorFlds[i].rename(surfTensorNames[i]);
979  surfTensorFlds[i].writeOpt() = IOobject::AUTO_WRITE;
980  }
981 
982  subsetter.subMesh().write();
983 
984 
985  //
986  // Splitting blockedFaces
987  // ~~~~~~~~~~~~~~~~~~~~~~
988  //
989 
990  // Synchronize wantedPatch across coupled patches.
992  (
993  subsetter.subMesh(),
994  wantedPatch,
996  );
997 
998  // Synchronize coupledWantedPatch across coupled patches.
1000  (
1001  subsetter.subMesh(),
1002  coupledWantedPatch,
1003  maxEqOp<label>()
1004  );
1005 
1006  // Synchronize cyclicWantedPatch across coupled patches.
1008  (
1009  subsetter.subMesh(),
1010  cyclicWantedPatch_half0,
1011  maxEqOp<label>()
1012  );
1013 
1014  // Synchronize cyclicWantedPatch across coupled patches.
1016  (
1017  subsetter.subMesh(),
1018  cyclicWantedPatch_half1,
1019  maxEqOp<label>()
1020  );
1021 
1022  // Topochange container
1023  polyTopoChange meshMod(subsetter.subMesh());
1024 
1025 
1026  // Whether first use of face (modify) or consecutive (add)
1027  PackedBoolList modifiedFace(mesh.nFaces());
1028 
1029  // Create coupled wall-side baffles
1030  createCoupledBaffles
1031  (
1032  subsetter.subMesh(),
1033  coupledWantedPatch,
1034  meshMod,
1035  modifiedFace
1036  );
1037 
1038  // Create coupled master/slave cyclic baffles
1039  createCyclicCoupledBaffles
1040  (
1041  subsetter.subMesh(),
1042  cyclicWantedPatch_half0,
1043  cyclicWantedPatch_half1,
1044  meshMod,
1045  modifiedFace
1046  );
1047 
1048  // Create wall baffles
1049  createBaffles
1050  (
1051  subsetter.subMesh(),
1052  wantedPatch,
1053  meshMod
1054  );
1055 
1056  if (!overwrite)
1057  {
1058  runTime++;
1059  }
1060 
1061  // Change the mesh. Change points directly (no inflation).
1062  autoPtr<mapPolyMesh> map = meshMod.changeMesh(subsetter.subMesh(), false);
1063 
1064  // Update fields
1065  subsetter.subMesh().updateMesh(map);
1066 
1067  // Fix faces that get mapped to zero-sized patches (these don't get any
1068  // value)
1069  initCreatedPatches<volScalarField>
1070  (
1071  subsetter.subMesh(),
1072  map,
1073  0.0
1074  );
1075  initCreatedPatches<volVectorField>
1076  (
1077  subsetter.subMesh(),
1078  map,
1079  vector::zero
1080  );
1081  initCreatedPatches<volSphericalTensorField>
1082  (
1083  subsetter.subMesh(),
1084  map,
1086  );
1087  initCreatedPatches<volSymmTensorField>
1088  (
1089  subsetter.subMesh(),
1090  map,
1092  );
1093  initCreatedPatches<volTensorField>
1094  (
1095  subsetter.subMesh(),
1096  map,
1097  tensor::zero
1098  );
1099 
1100  initCreatedPatches<surfaceScalarField>
1101  (
1102  subsetter.subMesh(),
1103  map,
1104  0.0
1105  );
1106  initCreatedPatches<surfaceVectorField>
1107  (
1108  subsetter.subMesh(),
1109  map,
1110  vector::zero
1111  );
1112  initCreatedPatches<surfaceSphericalTensorField>
1113  (
1114  subsetter.subMesh(),
1115  map,
1117  );
1118  initCreatedPatches<surfaceSymmTensorField>
1119  (
1120  subsetter.subMesh(),
1121  map,
1123  );
1124  initCreatedPatches<surfaceTensorField>
1125  (
1126  subsetter.subMesh(),
1127  map,
1128  tensor::zero
1129  );
1130 
1131 
1132  // Move mesh (since morphing might not do this)
1133  if (map().hasMotionPoints())
1134  {
1135  subsetter.subMesh().movePoints(map().preMotionPoints());
1136  }
1137 
1138  Info<< "Writing mesh with split blockedFaces to time " << runTime.value()
1139  << endl;
1140 
1141  subsetter.subMesh().write();
1142 
1143 
1144  //
1145  // Removing inaccessible regions
1146  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1147  //
1148 
1149  // Determine connected regions. regionSplit is the labelList with the
1150  // region per cell.
1151  regionSplit cellRegion(subsetter.subMesh());
1152 
1153  if (cellRegion.nRegions() > 1)
1154  {
1156  << "Removing blocked faces and cells created "
1157  << cellRegion.nRegions()
1158  << " regions that are not connected via a face." << nl
1159  << " This is not supported in solvers." << nl
1160  << " Use" << nl << nl
1161  << " splitMeshRegions <root> <case> -largestOnly" << nl << nl
1162  << " to extract a single region of the mesh." << nl
1163  << " This mesh will be written to a new timedirectory"
1164  << " so might have to be moved back to constant/" << nl
1165  << endl;
1166 
1167  word startFrom(runTime.controlDict().lookup("startFrom"));
1168 
1169  if (startFrom != "latestTime")
1170  {
1172  << "To run splitMeshRegions please set your"
1173  << " startFrom entry to latestTime" << endl;
1174  }
1175  }
1176 
1177  Info << nl << "End" << endl;
1178 
1179  return 0;
1180 }
1181 
1182 
1183 // ************************************************************************* //
Foam::fvMeshSubset::setLargeCellSubset
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:795
Foam::fvPatchField< Type >
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::fvMeshSubset::interpolate
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
Definition: fvMeshSubsetInterpolate.C:43
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::PackedBoolList
A bit-packed bool list.
Definition: PackedBoolList.H:63
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
cyclicFvPatch.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Tuple2.H
Foam::fvMeshSubset
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:73
Foam::polyTopoChange::changeMesh
autoPtr< mapPolyMesh > 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.
Definition: polyTopoChange.C:3039
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
fields
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
polyAddFace.H
fvMeshSubset.H
addOverwriteOption.H
mapPolyMesh.H
subsetVolFields
void subsetVolFields(const fvMeshSubset &subsetter, const wordList &fieldNames, PtrList< GeometricField< Type, fvPatchField, volMesh > > &subFields)
Definition: subsetMesh.C:124
Foam::minOp
Definition: ops.H:173
polyTopoChange.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::fvsPatchField< Type >
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:97
Foam::maxEqOp
Definition: ops.H:77
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
Foam::faceSet
A list of face labels.
Definition: faceSet.H:48
Foam::polyModifyFace
Class describing modification of a face.
Definition: polyModifyFace.H:48
IOobjectList.H
Foam::fvMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:726
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mapPolyMesh::oldPatchSizes
const labelList & oldPatchSizes() const
Return list of the old patch sizes.
Definition: mapPolyMesh.H:622
syncTools.H
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:109
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::label
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
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
argList.H
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
faceSet.H
regionSplit.H
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:102
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:802
Foam::ZoneMesh< faceZone, polyMesh >
Foam::regionSplit
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:119
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:226
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::Tensor::zero
static const Tensor zero
Definition: Tensor.H:80
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::SymmTensor< scalar >::zero
static const SymmTensor zero
Definition: SymmTensor.H:77
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:314
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:48
Foam::fvMeshSubset::patchMap
const labelList & patchMap() const
Return patch map.
Definition: fvMeshSubset.C:1550
subsetSurfaceFields
void subsetSurfaceFields(const fvMeshSubset &subsetter, const wordList &fieldNames, PtrList< GeometricField< Type, fvsPatchField, surfaceMesh > > &subFields)
Definition: subsetMesh.C:158
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
polyModifyFace.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:550
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
setRootCase.H
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::polyTopoChange::setAction
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
Definition: polyTopoChange.C:2530
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
f
labelList f(nPoints)
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:70
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
createMesh.H
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::fvMeshSubset::subMesh
const fvMesh & subMesh() const
Return reference to subset mesh.
Definition: fvMeshSubset.C:1472
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::polyAddFace
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:51
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
cellSet.H
Foam::syncTools::syncFaceList
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:381
Foam::fvMeshSubset::baseMesh
const fvMesh & baseMesh() const
Original mesh.
Definition: fvMeshSubset.H:216
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
args
Foam::argList args(argc, argv)
volField
Foam::tmp< Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > > volField(const Foam::fvMeshSubset &, const Foam::GeometricField< Type, Foam::fvPatchField, Foam::volMesh > &vf)
Wrapper to get hold of the field or the subsetted field.
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::SphericalTensor::zero
static const SphericalTensor zero
Definition: SphericalTensor.H:74
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
Foam::fvMeshSubset::faceMap
const labelList & faceMap() const
Return face map.
Definition: fvMeshSubset.C:1496