extrudeToRegionMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  extrudeToRegionMesh
26 
27 Description
28  Extrude faceZones (internal or boundary faces) or faceSets (boundary faces
29  only) into a separate mesh (as a different region).
30 
31  - used to e.g. extrude baffles (extrude internal faces) or create
32  liquid film regions.
33  - if extruding internal faces:
34  - create baffles in original mesh with mappedWall patches
35  - if extruding boundary faces:
36  - convert boundary faces to mappedWall patches
37  - extrude edges of faceZone as a <zone>_sidePatch
38  - extrudes into master direction (i.e. away from the owner cell
39  if flipMap is false)
40 
41 \verbatim
42 
43 Internal face extrusion
44 -----------------------
45 
46  +-------------+
47  | |
48  | |
49  +---AAAAAAA---+
50  | |
51  | |
52  +-------------+
53 
54  AAA=faceZone to extrude.
55 
56 
57 For the case of no flipMap the extrusion starts at owner and extrudes
58 into the space of the neighbour:
59 
60  +CCCCCCC+
61  | | <= extruded mesh
62  +BBBBBBB+
63 
64  +-------------+
65  | |
66  | (neighbour) |
67  |___CCCCCCC___| <= original mesh (with 'baffles' added)
68  | BBBBBBB |
69  |(owner side) |
70  | |
71  +-------------+
72 
73  BBB=mapped between owner on original mesh and new extrusion.
74  (zero offset)
75  CCC=mapped between neighbour on original mesh and new extrusion
76  (offset due to the thickness of the extruded mesh)
77 
78 For the case of flipMap the extrusion is the other way around: from the
79 neighbour side into the owner side.
80 
81 
82 Boundary face extrusion
83 -----------------------
84 
85  +--AAAAAAA--+
86  | |
87  | |
88  +-----------+
89 
90  AAA=faceZone to extrude. E.g. slave side is owner side (no flipmap)
91 
92 becomes
93 
94  +CCCCCCC+
95  | | <= extruded mesh
96  +BBBBBBB+
97 
98  +--BBBBBBB--+
99  | | <= original mesh
100  | |
101  +-----------+
102 
103  BBB=mapped between original mesh and new extrusion
104  CCC=polypatch
105 
106 
107  Notes:
108  - when extruding cyclics with only one cell in between it does not
109  detect this as a cyclic since the face is the same face. It will
110  only work if the coupled edge extrudes a different face so if there
111  are more than 1 cell in between.
112 
113 \endverbatim
114 
115 \*---------------------------------------------------------------------------*/
116 
117 #include "argList.H"
118 #include "polyTopoChange.H"
119 #include "mappedWallPolyPatch.H"
121 #include "createShellMesh.H"
122 #include "syncTools.H"
123 #include "cyclicPolyPatch.H"
124 #include "wedgePolyPatch.H"
125 #include "extrudeModel.H"
126 #include "faceSet.H"
127 #include "fvMeshTools.H"
128 #include "OBJstream.H"
129 #include "PatchTools.H"
130 #include "systemDict.H"
131 
132 using namespace Foam;
133 
134 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135 
136 label findIndex(const List<polyPatch*>& newPatches, const word& name)
137 {
138  forAll(newPatches, i)
139  {
140  if (newPatches[i]->name() == name)
141  {
142  return i;
143  }
144  }
145  return -1;
146 }
147 
148 
149 label addPatch
150 (
151  const polyBoundaryMesh& patches,
152  const word& patchName,
153  const word& patchType,
154  const dictionary& dict,
155  DynamicList<polyPatch*>& newPatches
156 )
157 {
158  label patchi = findIndex(newPatches, patchName);
159 
160  if (patchi != -1)
161  {
162  if (newPatches[patchi]->type() == patchType)
163  {
164  return patchi;
165  }
166  else
167  {
169  << "Already have patch " << patchName
170  << " but of type " << newPatches[patchi]->type()
171  << exit(FatalError);
172  }
173  }
174 
175  patchi = newPatches.size();
176 
177  label startFacei = 0;
178  if (patchi > 0)
179  {
180  const polyPatch& pp = *newPatches.last();
181  startFacei = pp.start()+pp.size();
182  }
183 
184  dictionary patchDict(dict);
185  patchDict.set("type", patchType);
186  patchDict.set("nFaces", 0);
187  patchDict.set("startFace", startFacei);
188 
189  newPatches.append
190  (
192  (
193  patchName,
194  patchDict,
195  patchi,
196  patches
197  ).ptr()
198  );
199 
200  return patchi;
201 }
202 
203 
204 // Remove zero-sized patches
205 void deleteEmptyPatches(fvMesh& mesh)
206 {
208 
209  wordList masterNames;
210  if (Pstream::master())
211  {
212  masterNames = patches.names();
213  }
214  Pstream::scatter(masterNames);
215 
216 
217  labelList oldToNew(patches.size(), -1);
218  label usedI = 0;
219  label notUsedI = patches.size();
220 
221  // Add all the non-empty, non-processor patches
222  forAll(masterNames, masterI)
223  {
224  label patchi = patches.findIndex(masterNames[masterI]);
225 
226  if (patchi != -1)
227  {
228  if (isA<processorPolyPatch>(patches[patchi]))
229  {
230  // Similar named processor patch? Not 'possible'.
231  if (patches[patchi].size() == 0)
232  {
233  Pout<< "Deleting processor patch " << patchi
234  << " name:" << patches[patchi].name()
235  << endl;
236  oldToNew[patchi] = --notUsedI;
237  }
238  else
239  {
240  oldToNew[patchi] = usedI++;
241  }
242  }
243  else
244  {
245  // Common patch.
246  if (returnReduce(patches[patchi].size(), sumOp<label>()) == 0)
247  {
248  Pout<< "Deleting patch " << patchi
249  << " name:" << patches[patchi].name()
250  << endl;
251  oldToNew[patchi] = --notUsedI;
252  }
253  else
254  {
255  oldToNew[patchi] = usedI++;
256  }
257  }
258  }
259  }
260 
261  // Add remaining patches at the end
263  {
264  if (oldToNew[patchi] == -1)
265  {
266  // Unique to this processor. Note: could check that these are
267  // only processor patches.
268  if (patches[patchi].size() == 0)
269  {
270  Pout<< "Deleting processor patch " << patchi
271  << " name:" << patches[patchi].name()
272  << endl;
273  oldToNew[patchi] = --notUsedI;
274  }
275  else
276  {
277  oldToNew[patchi] = usedI++;
278  }
279  }
280  }
281 
282  fvMeshTools::reorderPatches(mesh, oldToNew, usedI, true);
283 }
284 
285 
286 enum class connectivity
287 {
288  unset,
289  boundary,
290  internal,
291  error
292 };
293 
294 
295 struct isInternalEqOp
296 {
297  void operator()(connectivity& x, const connectivity& y) const
298  {
299  if (x == connectivity::unset)
300  {
301  // Set x if not yet set
302  x = y;
303  }
304  else if (y == connectivity::unset)
305  {
306  // Do not set x if y is not yet set
307  }
308  else if (x == y)
309  {
310  // Do nothing if x and y are equal
311  }
312  else
313  {
314  // Set an error value if x and y are both set and unequal
315  x = connectivity::error;
316  }
317  }
318 };
319 
320 
321 Ostream& operator<<(Ostream& os, const connectivity& p)
322 {
323  return os << static_cast<label>(p);
324 }
325 
326 
327 Istream& operator>>(Istream& is, connectivity& p)
328 {
329  label l;
330  is >> l;
331  p = static_cast<connectivity>(l);
332  return is;
333 }
334 
335 
336 // Return a bool list indicating whether or not the zone contains internal or
337 // boundary faces
338 boolList calcZoneIsInternal
339 (
340  const polyMesh& mesh,
341  const wordList& zoneNames,
342  const labelList& extrudeFaces,
343  const labelList& extrudeFaceZoneIDs
344 )
345 {
346  List<connectivity> zoneConnectivity(zoneNames.size(), connectivity::unset);
347 
348  forAll(extrudeFaces, i)
349  {
350  const label facei = extrudeFaces[i];
351  const label zonei = extrudeFaceZoneIDs[i];
352  isInternalEqOp()
353  (
354  zoneConnectivity[zonei],
355  mesh.isInternalFace(facei)
356  ? connectivity::internal
358  );
359  }
360 
361  Pstream::listCombineGather(zoneConnectivity, isInternalEqOp());
362  Pstream::listCombineScatter(zoneConnectivity);
363 
364  boolList zoneIsInternal(zoneNames.size());
365  forAll(zoneConnectivity, zonei)
366  {
367  switch (zoneConnectivity[zonei])
368  {
369  case connectivity::unset:
370  // This zone doesn't contain any faces, so it doesn't matter
371  // what we do here
372  zoneIsInternal[zonei] = false;
373  break;
375  zoneIsInternal[zonei] = false;
376  break;
377  case connectivity::internal:
378  zoneIsInternal[zonei] = true;
379  break;
380  case connectivity::error:
382  << "Zone " << zoneNames[zonei]
383  << " contains both internal and boundary faces."
384  << " This is not allowed." << exit(FatalError);
385  break;
386  }
387  }
388 
389  return zoneIsInternal;
390 }
391 
392 
393 // To combineReduce a labelList. Filters out duplicates.
394 struct uniqueEqOp
395 {
396  void operator()(labelList& x, const labelList& y) const
397  {
398  if (x.empty())
399  {
400  if (y.size())
401  {
402  x = y;
403  }
404  }
405  else
406  {
407  forAll(y, yi)
408  {
409  if (findIndex(x, y[yi]) == -1)
410  {
411  label sz = x.size();
412  x.setSize(sz+1);
413  x[sz] = y[yi];
414  }
415  }
416  }
417  }
418 };
419 
420 
421 // Calculate global pp faces per pp edge.
422 labelListList globalEdgeFaces
423 (
424  const polyMesh& mesh,
425  const globalIndex& globalFaces,
426  const primitiveFacePatch& pp,
427  const labelList& ppMeshEdges
428 )
429 {
430  // From mesh edge to global pp face labels.
431  labelListList globalEdgeFaces(ppMeshEdges.size());
432 
433  const labelListList& edgeFaces = pp.edgeFaces();
434 
435  forAll(edgeFaces, edgeI)
436  {
437  const labelList& eFaces = edgeFaces[edgeI];
438 
439  // Store pp face and processor as unique tag.
440  labelList& globalEFaces = globalEdgeFaces[edgeI];
441  globalEFaces.setSize(eFaces.size());
442  forAll(eFaces, i)
443  {
444  globalEFaces[i] = globalFaces.toGlobal(eFaces[i]);
445  }
446  }
447 
448  // Synchronise across coupled edges.
450  (
451  mesh,
452  ppMeshEdges,
453  globalEdgeFaces,
454  uniqueEqOp(),
455  labelList() // null value
456  );
457 
458  return globalEdgeFaces;
459 }
460 
461 
462 // Find a patch face that is not extruded. Return -1 if not found.
463 label findUncoveredPatchFace
464 (
465  const fvMesh& mesh,
466  const UIndirectList<label>& extrudeFaces, // mesh faces that are extruded
467  const label meshEdgeI // mesh edge
468 )
469 {
470  // Make set of extruded faces.
471  labelHashSet extrudeFaceSet(extrudeFaces.size());
472  forAll(extrudeFaces, i)
473  {
474  extrudeFaceSet.insert(extrudeFaces[i]);
475  }
476 
477  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
478  const labelList& eFaces = mesh.edgeFaces()[meshEdgeI];
479  forAll(eFaces, i)
480  {
481  label facei = eFaces[i];
482  label patchi = pbm.whichPatch(facei);
483 
484  if
485  (
486  patchi != -1
487  && !pbm[patchi].coupled()
488  && !extrudeFaceSet.found(facei)
489  )
490  {
491  return facei;
492  }
493  }
494  return -1;
495 }
496 
497 
498 // Same as findUncoveredPatchFace, except explicitly checks for cyclic faces
499 label findUncoveredCyclicPatchFace
500 (
501  const fvMesh& mesh,
502  const UIndirectList<label>& extrudeFaces, // mesh faces that are extruded
503  const label meshEdgeI // mesh edge
504 )
505 {
506  // Make set of extruded faces.
507  labelHashSet extrudeFaceSet(extrudeFaces.size());
508  forAll(extrudeFaces, i)
509  {
510  extrudeFaceSet.insert(extrudeFaces[i]);
511  }
512 
513  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
514  const labelList& eFaces = mesh.edgeFaces()[meshEdgeI];
515  forAll(eFaces, i)
516  {
517  label facei = eFaces[i];
518  label patchi = pbm.whichPatch(facei);
519 
520  if
521  (
522  patchi != -1
523  && isA<cyclicPolyPatch>(pbm[patchi])
524  && !extrudeFaceSet.found(facei)
525  )
526  {
527  return facei;
528  }
529  }
530  return -1;
531 }
532 
533 
534 // Add coupled patches into the mesh
535 void addCouplingPatches
536 (
537  const fvMesh& mesh,
538  const bool isShellMesh,
539  const bool intrude,
540  const word& regionName,
541  const word& nbrRegionName,
542  const wordList& zoneNames,
543  const wordList& oppositeZoneNames,
544  const boolList& zoneIsInternal,
545  const dictionary& dict,
546  DynamicList<polyPatch*>& newPatches,
547  labelList& zoneTopPatch,
548  labelList& zoneBottomPatch
549 )
550 {
551  Pout<< "Adding coupling patches:" << nl << nl
552  << "patchID\tpatch\ttype" << nl
553  << "-------\t-----\t----"
554  << endl;
555 
557  (
558  dict.lookupOrDefault("patchNames", wordList())
559  );
560  wordList regionPatchNames
561  (
562  dict.lookupOrDefault("regionPatchNames", wordList())
563  );
564  if (isShellMesh)
565  {
566  patchNames.swap(regionPatchNames);
567  }
568 
569  const wordList patchTypes
570  (
571  isShellMesh
572  ? dict.lookupOrDefault("regionPatchTypes", wordList())
573  : dict.lookupOrDefault("patchTypes", wordList())
574  );
575 
576  wordList oppositePatchNames
577  (
578  dict.lookupOrDefault("oppositePatchNames", wordList())
579  );
580  wordList regionOppositePatchNames
581  (
582  dict.lookupOrDefault("regionOppositePatchNames", wordList())
583  );
584  if (isShellMesh)
585  {
586  oppositePatchNames.swap(regionOppositePatchNames);
587  }
588 
589  const wordList oppositePatchTypes
590  (
591  isShellMesh
592  ? dict.lookupOrDefault("regionOppositePatchTypes", wordList())
593  : dict.lookupOrDefault("oppositePatchType", wordList())
594  );
595 
596  zoneTopPatch.setSize(zoneNames.size(), -1);
597  zoneBottomPatch.setSize(zoneNames.size(), -1);
598 
599  dictionary patchDict;
600  patchDict.add("neighbourRegion", nbrRegionName);
601 
602  label nOldPatches = newPatches.size();
603  forAll(zoneNames, zonei)
604  {
605  const word patchNamePrefix =
606  regionName + "_to_" + nbrRegionName + '_';
607 
608  const word nbrPatchNamePrefix =
609  nbrRegionName + "_to_" + regionName + '_';
610 
611  word bottomPatchName;
612  word bottomNbrPatchName;
613  word topPatchName;
614  word topNbrPatchName;
615 
616  if (zoneIsInternal[zonei])
617  {
618  bottomPatchName =
619  patchNamePrefix + zoneNames[zonei] + "_bottom";
620  bottomNbrPatchName =
621  nbrPatchNamePrefix + zoneNames[zonei] + "_bottom";
622  topPatchName =
623  patchNamePrefix + zoneNames[zonei] + "_top";
624  topNbrPatchName =
625  nbrPatchNamePrefix + zoneNames[zonei] + "_top";
626  }
627  else if (!oppositeZoneNames[zonei].empty())
628  {
629  bottomPatchName = patchNamePrefix + zoneNames[zonei];
630  bottomNbrPatchName = nbrPatchNamePrefix + zoneNames[zonei];
631  topPatchName = patchNamePrefix + oppositeZoneNames[zonei];
632  topNbrPatchName = nbrPatchNamePrefix + oppositeZoneNames[zonei];
633  }
634  else
635  {
636  bottomPatchName = patchNamePrefix + zoneNames[zonei];
637  bottomNbrPatchName = nbrPatchNamePrefix + zoneNames[zonei];
638  topPatchName = zoneNames[zonei] + "_top";
639  topNbrPatchName = word::null;
640  }
641 
642  if (patchNames.size())
643  {
644  bottomPatchName = patchNames[zonei];
645  }
646  if (regionPatchNames.size())
647  {
648  bottomNbrPatchName = regionPatchNames[zonei];
649  }
650  if (oppositePatchNames.size())
651  {
652  topPatchName = oppositePatchNames[zonei];
653  }
654  if (regionOppositePatchNames.size())
655  {
656  topNbrPatchName = regionOppositePatchNames[zonei];
657  }
658 
659  const bool bothMapped =
660  zoneIsInternal[zonei] || !oppositeZoneNames[zonei].empty();
661 
662  const bool bottomMapped = bothMapped || !(intrude && isShellMesh);
663  const bool topMapped = bothMapped || intrude;
664 
665  const bool bottomExtruded = !isShellMesh && intrude;
666  const bool topExtruded = !bottomExtruded;
667 
668  dictionary bottomPatchDict;
669  if (bottomMapped)
670  {
671  bottomPatchDict = patchDict;
672  bottomPatchDict.add
673  (
674  "neighbourPatch",
675  !intrude ? bottomNbrPatchName : topNbrPatchName
676  );
677  if (bottomExtruded)
678  {
679  bottomPatchDict.add("isExtrudedRegion", isShellMesh);
680  }
681  }
682 
683  zoneBottomPatch[zonei] = addPatch
684  (
685  mesh.boundaryMesh(),
686  bottomPatchName,
687  patchTypes.size() ? patchTypes[zonei]
688  : !bottomMapped ? polyPatch::typeName
689  : !bottomExtruded
690  ? mappedWallPolyPatch::typeName
691  : mappedExtrudedWallPolyPatch::typeName,
692  bottomPatchDict,
693  newPatches
694  );
695 
696  Pout<< zoneBottomPatch[zonei]
697  << '\t' << newPatches[zoneBottomPatch[zonei]]->name()
698  << '\t' << newPatches[zoneBottomPatch[zonei]]->type()
699  << nl;
700 
701  if (!isShellMesh && !bothMapped) continue;
702 
703  dictionary topPatchDict;
704  if (topMapped)
705  {
706  topPatchDict = patchDict;
707  topPatchDict.add
708  (
709  "neighbourPatch",
710  !intrude ? topNbrPatchName : bottomNbrPatchName
711  );
712  if (topExtruded)
713  {
714  topPatchDict.add("isExtrudedRegion", isShellMesh);
715  }
716  }
717 
718  zoneTopPatch[zonei] = addPatch
719  (
720  mesh.boundaryMesh(),
721  topPatchName,
722  oppositePatchTypes.size() ? oppositePatchTypes[zonei]
723  : !topMapped ? polyPatch::typeName
724  : !topExtruded
725  ? mappedWallPolyPatch::typeName
726  : mappedExtrudedWallPolyPatch::typeName,
727  topPatchDict,
728  newPatches
729  );
730 
731  Pout<< zoneTopPatch[zonei]
732  << '\t' << newPatches[zoneTopPatch[zonei]]->name()
733  << '\t' << newPatches[zoneTopPatch[zonei]]->type()
734  << nl;
735  }
736 
737  Pout<< "Added " << newPatches.size()-nOldPatches
738  << " inter-region patches." << nl
739  << endl;
740 }
741 
742 
743 // Count the number of faces in patches that need to be created
744 labelList countExtrudePatches
745 (
746  const fvMesh& mesh,
747  const label nZones,
748  const primitiveFacePatch& extrudePatch,
749  const labelList& extrudeFaces,
750  const labelList& extrudeFaceZoneIDs,
751  const labelList& extrudeEdgeMeshEdges,
752  const labelListList& extrudeEdgeGlobalFaces
753 )
754 {
755  labelList zoneSideNFaces(nZones, 0);
756 
757  forAll(extrudePatch.edgeFaces(), edgeI)
758  {
759  const labelList& eFaces = extrudePatch.edgeFaces()[edgeI];
760 
761  if (eFaces.size() == 2)
762  {
763  // Internal edge
764  }
765  else if
766  (
767  eFaces.size() == 1
768  && extrudeEdgeGlobalFaces[edgeI].size() == 2
769  )
770  {
771  // Coupled edge
772  }
773  else
774  {
775  // Perimeter edge. Check whether we are on a mesh edge with
776  // external patches. If so choose any uncovered one. If none found
777  // put face in undetermined zone 'side' patch.
778  const label facei = findUncoveredPatchFace
779  (
780  mesh,
781  UIndirectList<label>(extrudeFaces, eFaces),
782  extrudeEdgeMeshEdges[edgeI]
783  );
784 
785  if (facei == -1)
786  {
787  forAll(extrudePatch.edgeFaces()[edgeI], i)
788  {
789  const label facei = extrudePatch.edgeFaces()[edgeI][i];
790  zoneSideNFaces[extrudeFaceZoneIDs[facei]] ++;
791  }
792  }
793  }
794  }
795 
796  Pstream::listCombineGather(zoneSideNFaces, plusEqOp<label>());
797  Pstream::listCombineScatter(zoneSideNFaces);
798 
799  return zoneSideNFaces;
800 }
801 
802 
803 // Add side patches into the mesh
804 labelList addZoneSidePatches
805 (
806  const fvMesh& mesh,
807  const wordList& zoneNames,
808  const labelList& zoneSideNFaces,
809  const word& oneDPolyPatchType,
810  DynamicList<polyPatch*>& newPatches
811 )
812 {
813  Pout<< "Adding patches for sides on zones:" << nl << nl
814  << "patchID\tpatch" << nl << "-------\t-----" << endl;
815 
816  labelList zoneSidePatches(zoneNames.size(), -labelMax);
817 
818  const label nOldPatches = newPatches.size();
819 
820  if (oneDPolyPatchType != word::null)
821  {
822  forAll(zoneSideNFaces, zoneI)
823  {
824  word patchName;
825 
826  if (oneDPolyPatchType == "empty")
827  {
828  patchName = "oneDEmptyPatch";
829  zoneSidePatches[zoneI] = addPatch
830  (
831  mesh.boundaryMesh(),
832  patchName,
833  emptyPolyPatch::typeName,
834  dictionary(),
835  newPatches
836  );
837  }
838  else if (oneDPolyPatchType == "wedge")
839  {
840  patchName = "oneDWedgePatch";
841  zoneSidePatches[zoneI] = addPatch
842  (
843  mesh.boundaryMesh(),
844  patchName,
845  wedgePolyPatch::typeName,
846  dictionary(),
847  newPatches
848  );
849  }
850  else
851  {
853  << "Type " << oneDPolyPatchType << " does not exist "
854  << exit(FatalError);
855  }
856 
857  Pout<< zoneSidePatches[zoneI] << '\t' << patchName << nl;
858  }
859  }
860  else
861  {
862  forAll(zoneSideNFaces, zoneI)
863  {
864  if (zoneSideNFaces[zoneI] > 0)
865  {
866  word patchName = zoneNames[zoneI] + "_" + "side";
867 
868  zoneSidePatches[zoneI] = addPatch
869  (
870  mesh.boundaryMesh(),
871  patchName,
872  polyPatch::typeName,
873  dictionary(),
874  newPatches
875  );
876 
877  Pout<< zoneSidePatches[zoneI] << '\t' << patchName << nl;
878  }
879  }
880  }
881 
882  Pout<< "Added " << newPatches.size() - nOldPatches << " zone-side patches."
883  << nl << endl;
884 
885  return zoneSidePatches;
886 }
887 
888 
889 // Add any additional coupled side patches that might be necessary. Return
890 // correspondence between extruded edges and their side patches.
891 labelList addExtrudeEdgeSidePatches
892 (
893  const fvMesh& mesh,
894  const primitiveFacePatch& extrudePatch,
895  const labelList& extrudeFaces,
896  const labelList& extrudeEdgeMeshEdges,
897  const distributionMap& extrudeEdgeFacesMap,
898  const labelListList& extrudeEdgeGlobalFaces,
899  DynamicList<polyPatch*>& newPatches
900 )
901 {
902  // Get procID in extrudeEdgeGlobalFaces order
903  labelList procID(extrudeEdgeGlobalFaces.size(), Pstream::myProcNo());
904  extrudeEdgeFacesMap.distribute(procID);
905 
906  // Calculate opposite processor for coupled edges (only if shared by
907  // two procs). Note: Could have saved original globalEdgeFaces structure.
908  labelList minProcID(extrudeEdgeGlobalFaces.size(), labelMax);
909  labelList maxProcID(extrudeEdgeGlobalFaces.size(), labelMin);
910  forAll(extrudeEdgeGlobalFaces, edgeI)
911  {
912  const labelList& eFaces = extrudeEdgeGlobalFaces[edgeI];
913  if (eFaces.size())
914  {
915  forAll(eFaces, i)
916  {
917  label proci = procID[eFaces[i]];
918  minProcID[edgeI] = min(minProcID[edgeI], proci);
919  maxProcID[edgeI] = max(maxProcID[edgeI], proci);
920  }
921  }
922  }
924  (
925  mesh,
926  extrudeEdgeMeshEdges,
927  minProcID,
928  minEqOp<label>(),
929  labelMax // null value
930  );
932  (
933  mesh,
934  extrudeEdgeMeshEdges,
935  maxProcID,
936  maxEqOp<label>(),
937  labelMin // null value
938  );
939 
940  Pout<< "Adding processor or cyclic patches:" << nl << nl
941  << "patchID\tpatch" << nl << "-------\t-----" << endl;
942 
943  const label nOldPatches = newPatches.size();
944 
945  labelList extrudeEdgeSidePatches(extrudePatch.edgeFaces().size(), -1);
946  forAll(extrudePatch.edgeFaces(), edgeI)
947  {
948  const labelList& eFaces = extrudePatch.edgeFaces()[edgeI];
949 
950  if
951  (
952  eFaces.size() == 1
953  && extrudeEdgeGlobalFaces[edgeI].size() == 2
954  )
955  {
956  // Coupled boundary edge. Find matching patch.
957  label nbrProci = minProcID[edgeI];
958  if (nbrProci == Pstream::myProcNo())
959  {
960  nbrProci = maxProcID[edgeI];
961  }
962 
963  if (nbrProci == Pstream::myProcNo())
964  {
965  // Cyclic patch since both procs the same. This cyclic should
966  // already exist in newPatches so no adding necessary.
967 
968  const label facei = findUncoveredCyclicPatchFace
969  (
970  mesh,
971  UIndirectList<label>(extrudeFaces, eFaces),
972  extrudeEdgeMeshEdges[edgeI]
973  );
974 
975  if (facei != -1)
976  {
978 
979  const label newPatchi = findIndex
980  (
981  newPatches,
982  patches[patches.whichPatch(facei)].name()
983  );
984 
985  extrudeEdgeSidePatches[edgeI] = newPatchi;
986  }
987  else
988  {
990  << "Unable to determine coupled patch addressing"
991  << abort(FatalError);
992  }
993  }
994  else
995  {
996  // Processor patch
997  const word name =
999 
1000  extrudeEdgeSidePatches[edgeI] = findIndex(newPatches, name);
1001 
1002  if (extrudeEdgeSidePatches[edgeI] == -1)
1003  {
1004  dictionary patchDict;
1005  patchDict.add("myProcNo", Pstream::myProcNo());
1006  patchDict.add("neighbProcNo", nbrProci);
1007 
1008  extrudeEdgeSidePatches[edgeI] = addPatch
1009  (
1010  mesh.boundaryMesh(),
1011  name,
1012  processorPolyPatch::typeName,
1013  patchDict,
1014  newPatches
1015  );
1016 
1017  Pout<< extrudeEdgeSidePatches[edgeI] << '\t' << name << nl;
1018  }
1019  }
1020  }
1021  }
1022 
1023  Pout<< "Added " << newPatches.size() - nOldPatches
1024  << " coupled patches." << nl << endl;
1025 
1026  return extrudeEdgeSidePatches;
1027 }
1028 
1029 
1030 int main(int argc, char *argv[])
1031 {
1032  argList::addNote("Create region mesh by extruding a faceZone or faceSet");
1033 
1034  #include "addMeshOption.H"
1035  #include "addRegionOption.H"
1036  #include "addNoOverwriteOption.H"
1037  #include "addDictOption.H"
1038  #include "setRootCase.H"
1039  #include "createTime.H"
1041 
1042  if (mesh.boundaryMesh().checkParallelSync(true))
1043  {
1044  List<wordList> allNames(Pstream::nProcs());
1045  allNames[Pstream::myProcNo()] = mesh.boundaryMesh().names();
1046  Pstream::gatherList(allNames);
1047  Pstream::scatterList(allNames);
1048 
1050  << "Patches are not synchronised on all processors."
1051  << " Per processor patches " << allNames
1052  << exit(FatalError);
1053  }
1054 
1055  #include "setNoOverwrite.H"
1056 
1057  const dictionary dict(systemDict("extrudeToRegionMeshDict", args, mesh));
1058 
1059  // Region to create by extrusion
1060  const word shellRegionName(dict.lookup("region"));
1061 
1062  // Should the extruded region overlap the existing region, i.e. "intrude"?
1063  const Switch intrude(dict.lookupOrDefault("intrude", false));
1064 
1065  if (shellRegionName == regionName)
1066  {
1067  FatalIOErrorIn(args.executable().c_str(), dict)
1068  << "Cannot extrude into same region as mesh." << endl
1069  << "Mesh region : " << regionName << endl
1070  << "Shell region : " << shellRegionName
1071  << exit(FatalIOError);
1072  }
1073 
1074  // Select faces to extrude
1075  enum class zoneSourceType
1076  {
1077  zone,
1078  set,
1079  patch
1080  };
1081 
1082  static const wordList zoneSourceTypeNames
1083  {
1084  "faceZone",
1085  "faceSet",
1086  "patch"
1087  };
1088 
1089  static const wordList zoneSourcesTypeNames
1090  {
1091  "faceZones",
1092  "faceSets",
1093  "patches"
1094  };
1095 
1096  wordList zoneNames;
1097  wordList oppositeZoneNames;
1098  List<zoneSourceType> zoneSourceTypes;
1099 
1100  auto lookupZones = [&](const zoneSourceType& type)
1101  {
1102  const word& keyword = zoneSourcesTypeNames[unsigned(type)];
1103 
1104  if (dict.found(keyword))
1105  {
1106  zoneNames.append(dict.lookup<wordList>(keyword));
1107 
1108  oppositeZoneNames.append
1109  (
1111  (
1112  {
1113  "opposite" + keyword.capitalise(),
1114  keyword + "Shadow"
1115  },
1116  wordList
1117  (
1118  zoneNames.size() - oppositeZoneNames.size(),
1119  word::null
1120  )
1121  )
1122  );
1123 
1124  zoneSourceTypes.setSize(zoneNames.size(), type);
1125  }
1126  };
1127  lookupZones(zoneSourceType::zone);
1128  lookupZones(zoneSourceType::set);
1129  lookupZones(zoneSourceType::patch);
1130 
1131  Info<< nl << "Extruding:" << nl << incrIndent;
1132  forAll(zoneNames, zonei)
1133  {
1134  const unsigned typei = unsigned(zoneSourceTypes[zonei]);
1135 
1136  if (oppositeZoneNames[zonei].empty())
1137  {
1138  Info<< indent << "From " << zoneSourceTypeNames[typei] << " \""
1139  << zoneNames[zonei] << "\"" << nl;
1140  }
1141  else
1142  {
1143  Info<< indent << "Between " << zoneSourcesTypeNames[typei] << " \""
1144  << zoneNames[zonei] << "\" and \"" << oppositeZoneNames[zonei]
1145  << "\"" << nl;
1146  }
1147  }
1148  Info<< endl << decrIndent;
1149 
1150  // One-dimensional extrusion settings
1151  const Switch oneD(dict.lookupOrDefault("oneD", false));
1152  Switch oneDNonManifoldEdges(false);
1153  word oneDPatchType(emptyPolyPatch::typeName);
1154  if (oneD)
1155  {
1156  oneDNonManifoldEdges = dict.lookupOrDefault("nonManifold", false);
1157  dict.lookup("oneDPolyPatchType") >> oneDPatchType;
1158  }
1159 
1160  if (oneD)
1161  {
1162  if (oneDNonManifoldEdges)
1163  {
1164  Info<< "Extruding as 1D columns with sides in patch type "
1165  << oneDPatchType
1166  << " and connected points (except on non-manifold areas)."
1167  << endl;
1168  }
1169  else
1170  {
1171  Info<< "Extruding as 1D columns with sides in patch type "
1172  << oneDPatchType
1173  << " and duplicated points (overlapping volumes)."
1174  << endl;
1175  }
1176  }
1177 
1178  // Construct the point generator
1180 
1181  // Change the primary mesh?
1182  const Switch adaptMesh(dict.lookup("adaptMesh"));
1183 
1184 
1185  // Determine output instance
1186  word meshInstance;
1187  if (!overwrite)
1188  {
1189  runTime++;
1190  meshInstance = runTime.name();
1191  }
1192  else
1193  {
1194  meshInstance = mesh.pointsInstance();
1195  }
1196  Info<< "Writing meshes to " << meshInstance << nl << endl;
1197 
1198 
1199  // Map from extrude zone to mesh zone, or -1 if not a mesh zone
1200  labelList zoneMeshZoneID(zoneNames.size(), -1);
1201  labelList oppositeZoneMeshZoneID(zoneNames.size(), -1);
1202  forAll(zoneNames, zonei)
1203  {
1204  if (zoneSourceTypes[zonei] != zoneSourceType::zone) continue;
1205 
1206  zoneMeshZoneID[zonei] =
1207  mesh.faceZones().findIndex(zoneNames[zonei]);
1208 
1209  if (zoneMeshZoneID[zonei] == -1)
1210  {
1211  FatalIOErrorIn(args.executable().c_str(), dict)
1212  << "Cannot find zone " << zoneNames[zonei]
1213  << endl << "Valid zones are " << mesh.faceZones().toc()
1214  << exit(FatalIOError);
1215  }
1216 
1217  if (!oppositeZoneNames[zonei].empty())
1218  {
1219  oppositeZoneMeshZoneID[zonei] =
1220  mesh.faceZones().findIndex(oppositeZoneNames[zonei]);
1221 
1222  if (oppositeZoneMeshZoneID[zonei] == -1)
1223  {
1224  FatalIOErrorIn(args.executable().c_str(), dict)
1225  << "Cannot find opposite zone " << oppositeZoneNames[zonei]
1226  << endl << "Valid zones are " << mesh.faceZones().toc()
1227  << exit(FatalIOError);
1228  }
1229  }
1230  }
1231 
1232 
1233  // Extract faces to extrude
1234  labelList extrudeFaces, oppositeExtrudeFaces;
1235  labelList extrudeFaceZoneIDs, oppositeExtrudeFaceZoneIDs;
1236  boolList extrudeFaceFlips, oppositeExtrudeFaceFlips;
1237  {
1238  // Load any faceSets that we need
1239  PtrList<faceSet> zoneSets(zoneNames.size());
1240  PtrList<faceSet> oppositeZoneSets(zoneNames.size());
1241  forAll(zoneNames, zonei)
1242  {
1243  if (zoneSourceTypes[zonei] == zoneSourceType::set)
1244  {
1245  zoneSets.set(zonei, new faceSet(mesh, zoneNames[zonei]));
1246  if (!oppositeZoneNames.empty())
1247  {
1248  oppositeZoneSets.set
1249  (
1250  zonei,
1251  new faceSet(mesh, oppositeZoneNames[zonei])
1252  );
1253  }
1254  }
1255  }
1256 
1257  // Create dynamic face lists
1258  DynamicList<label> facesDyn, oppositeFacesDyn;
1259  DynamicList<label> zoneIDsDyn, oppositeZoneIDsDyn;
1260  DynamicList<bool> flipsDyn, oppositeFlipsDyn;
1261  forAll(zoneNames, zonei)
1262  {
1263  switch (zoneSourceTypes[zonei])
1264  {
1265  case zoneSourceType::zone:
1266  {
1267  const faceZone& fz =
1268  mesh.faceZones()[zoneMeshZoneID[zonei]];
1269  facesDyn.append(fz);
1270  zoneIDsDyn.append(labelList(fz.size(), zonei));
1271  flipsDyn.append(fz.flipMap());
1272 
1273  if (!oppositeZoneNames[zonei].empty())
1274  {
1275  const faceZone& sfz =
1276  mesh.faceZones()[oppositeZoneMeshZoneID[zonei]];
1277  if (sfz.size() != fz.size())
1278  {
1279  FatalIOErrorIn(args.executable().c_str(), dict)
1280  << "Opposite zone " << oppositeZoneNames[zonei]
1281  << "is a different size from its "
1282  << "corresponding zone " << zoneNames[zonei]
1283  << exit(FatalIOError);
1284  }
1285  oppositeFacesDyn.append(sfz);
1286  oppositeZoneIDsDyn.append(labelList(sfz.size(), zonei));
1287  oppositeFlipsDyn.append(sfz.flipMap());
1288  }
1289  else
1290  {
1291  oppositeFacesDyn.append(labelList(fz.size(), -1));
1292  oppositeZoneIDsDyn.append(labelList(fz.size(), -1));
1293  oppositeFlipsDyn.append(boolList(fz.size(), false));
1294  }
1295  break;
1296  }
1297  case zoneSourceType::set:
1298  {
1299  const faceSet& fs = zoneSets[zonei];
1300  facesDyn.append(fs.toc());
1301  zoneIDsDyn.append(labelList(fs.size(), zonei));
1302  flipsDyn.append(boolList(fs.size(), false));
1303 
1304  if (!oppositeZoneNames[zonei].empty())
1305  {
1306  const faceSet& sfs = oppositeZoneSets[zonei];
1307  if (sfs.size() != fs.size())
1308  {
1309  FatalIOErrorIn(args.executable().c_str(), dict)
1310  << "Opposite set " << oppositeZoneNames[zonei]
1311  << "is a different size from its "
1312  << "corresponding zone " << zoneNames[zonei]
1313  << exit(FatalIOError);
1314  }
1315  oppositeFacesDyn.append(sfs.toc());
1316  oppositeZoneIDsDyn.append(labelList(sfs.size(), zonei));
1317  oppositeFlipsDyn.append(boolList(sfs.size(), false));
1318  }
1319  else
1320  {
1321  oppositeFacesDyn.append(labelList(fs.size(), -1));
1322  oppositeZoneIDsDyn.append(labelList(fs.size(), -1));
1323  oppositeFlipsDyn.append(boolList(fs.size(), false));
1324  }
1325  break;
1326  }
1327  case zoneSourceType::patch:
1328  {
1329  const polyPatch& pp =
1330  mesh.boundaryMesh()[zoneNames[zonei]];
1331  facesDyn.append(pp.start() + identityMap(pp.size()));
1332  zoneIDsDyn.append(labelList(pp.size(), zonei));
1333  flipsDyn.append(boolList(pp.size(), false));
1334 
1335  if (!oppositeZoneNames[zonei].empty())
1336  {
1337  const polyPatch& spp =
1338  mesh.boundaryMesh()[oppositeZoneNames[zonei]];
1339  if (spp.size() != pp.size())
1340  {
1341  FatalIOErrorIn(args.executable().c_str(), dict)
1342  << "Opposite patch " << oppositeZoneNames[zonei]
1343  << "is a different size from its "
1344  << "corresponding zone " << zoneNames[zonei]
1345  << exit(FatalIOError);
1346  }
1347  oppositeFacesDyn.append
1348  (
1349  spp.start() + identityMap(spp.size())
1350  );
1351  oppositeZoneIDsDyn.append(labelList(spp.size(), zonei));
1352  oppositeFlipsDyn.append(boolList(spp.size(), false));
1353  }
1354  else
1355  {
1356  oppositeFacesDyn.append(labelList(pp.size(), -1));
1357  oppositeZoneIDsDyn.append(labelList(pp.size(), -1));
1358  oppositeFlipsDyn.append(boolList(pp.size(), false));
1359  }
1360  break;
1361  }
1362  }
1363  }
1364 
1365  // Transfer to non-dynamic storage
1366  extrudeFaces.transfer(facesDyn);
1367  extrudeFaceZoneIDs.transfer(zoneIDsDyn);
1368  extrudeFaceFlips.transfer(flipsDyn);
1369  oppositeExtrudeFaces.transfer(oppositeFacesDyn);
1370  oppositeExtrudeFaceZoneIDs.transfer(oppositeZoneIDsDyn);
1371  oppositeExtrudeFaceFlips.transfer(oppositeFlipsDyn);
1372  }
1373 
1374  faceList extrudeFaceList(UIndirectList<face>(mesh.faces(), extrudeFaces));
1375 
1376  forAll(extrudeFaceList, facei)
1377  {
1378  if (intrude != extrudeFaceFlips[facei])
1379  {
1380  extrudeFaceList[facei].flip();
1381  }
1382  }
1383 
1384  // Create a primitive patch of the extruded faces
1385  const primitiveFacePatch extrudePatch
1386  (
1387  extrudeFaceList,
1388  mesh.points()
1389  );
1390 
1391  // Check zone either all internal or all external faces
1392  const boolList zoneIsInternal
1393  (
1394  calcZoneIsInternal(mesh, zoneNames, extrudeFaces, extrudeFaceZoneIDs)
1395  );
1396 
1397 
1398  Info<< "extrudePatch :"
1399  << " faces:" << returnReduce(extrudePatch.size(), sumOp<label>())
1400  << " points:" << returnReduce(extrudePatch.nPoints(), sumOp<label>())
1401  << " edges:" << returnReduce(extrudePatch.nEdges(), sumOp<label>())
1402  << nl << endl;
1403 
1404 
1405  // Determine per-extrude-edge info
1406  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1407 
1408  // Corresponding mesh edges
1409  const labelList extrudeEdgeMeshEdges
1410  (
1411  extrudePatch.meshEdges
1412  (
1413  mesh.edges(),
1414  mesh.pointEdges()
1415  )
1416  );
1417 
1418  const globalIndex globalExtrudeFaces(extrudePatch.size());
1419 
1420  // Global pp faces per pp edge.
1421  labelListList extrudeEdgeGlobalFaces
1422  (
1423  globalEdgeFaces
1424  (
1425  mesh,
1426  globalExtrudeFaces,
1427  extrudePatch,
1428  extrudeEdgeMeshEdges
1429  )
1430  );
1431  List<Map<label>> compactMap;
1432  const distributionMap extrudeEdgeFacesMap
1433  (
1434  globalExtrudeFaces,
1435  extrudeEdgeGlobalFaces,
1436  compactMap
1437  );
1438 
1439 
1440  // Copy all non-local patches since these are used on boundary edges of
1441  // the extrusion
1442  DynamicList<polyPatch*> regionPatches(mesh.boundaryMesh().size());
1444  {
1445  if (!isA<processorPolyPatch>(mesh.boundaryMesh()[patchi]))
1446  {
1447  regionPatches.append
1448  (
1450  (
1451  mesh.boundaryMesh(),
1452  regionPatches.size(),
1453  0, // size
1454  0 // start
1455  ).ptr()
1456  );
1457  }
1458  }
1459 
1460 
1461  // Add interface patches
1462  // ~~~~~~~~~~~~~~~~~~~~~
1463 
1464  // From zone to interface patch (region side)
1465  labelList zoneTopPatch, zoneBottomPatch;
1466  addCouplingPatches
1467  (
1468  mesh,
1469  true,
1470  intrude,
1471  shellRegionName,
1472  regionName,
1473  zoneNames,
1474  oppositeZoneNames,
1475  zoneIsInternal,
1476  dict,
1477  regionPatches,
1478  zoneTopPatch,
1479  zoneBottomPatch
1480  );
1481 
1482  // From zone to interface patch (mesh side)
1483  labelList interMeshTopPatch;
1484  labelList interMeshBottomPatch;
1485  if (adaptMesh)
1486  {
1488 
1489  // Clone existing non-processor patches
1490  DynamicList<polyPatch*> newPatches(patches.size());
1492  {
1493  if (!isA<processorPolyPatch>(patches[patchi]))
1494  {
1495  newPatches.append(patches[patchi].clone(patches).ptr());
1496  }
1497  }
1498 
1499  // Add new patches
1500  addCouplingPatches
1501  (
1502  mesh,
1503  false,
1504  intrude,
1505  regionName,
1506  shellRegionName,
1507  zoneNames,
1508  oppositeZoneNames,
1509  zoneIsInternal,
1510  dict,
1511  newPatches,
1512  interMeshTopPatch,
1513  interMeshBottomPatch
1514  );
1515 
1516  // Clone existing processor patches
1518  {
1519  if (isA<processorPolyPatch>(patches[patchi]))
1520  {
1521  newPatches.append
1522  (
1524  (
1525  patches,
1526  newPatches.size(),
1527  patches[patchi].size(),
1528  patches[patchi].start()
1529  ).ptr()
1530  );
1531  }
1532  }
1533 
1534  // Add to mesh
1535  mesh.clearOut();
1537  mesh.addFvPatches(newPatches, true);
1538 
1539  // Note: from this point on mesh patches differs from regionPatches
1540  }
1541 
1542 
1543  // Patch per extruded face
1544  labelList extrudeFaceTopPatchID(extrudePatch.size());
1545  labelList extrudeFaceBottomPatchID(extrudePatch.size());
1546  forAll(extrudeFaceZoneIDs, facei)
1547  {
1548  extrudeFaceTopPatchID[facei] =
1549  zoneTopPatch[extrudeFaceZoneIDs[facei]];
1550  extrudeFaceBottomPatchID[facei] =
1551  zoneBottomPatch[extrudeFaceZoneIDs[facei]];
1552  }
1553 
1554 
1555  // Count how many patches on special edges of extrudePatch are necessary
1556  const labelList zoneSideNFaces
1557  (
1558  countExtrudePatches
1559  (
1560  mesh,
1561  zoneNames.size(),
1562  extrudePatch, // patch
1563  extrudeFaces, // mesh face per patch face
1564  extrudeFaceZoneIDs, // ...
1565  extrudeEdgeMeshEdges, // mesh edge per patch edge
1566  extrudeEdgeGlobalFaces // global indexing per patch edge
1567  )
1568  );
1569 
1570 
1571  // Add the zone-side patches.
1572  const labelList zoneSidePatches
1573  (
1574  addZoneSidePatches
1575  (
1576  mesh,
1577  zoneNames,
1578  zoneSideNFaces,
1579  (oneD ? oneDPatchType : word::null),
1580  regionPatches
1581  )
1582  );
1583 
1584 
1585  // Sets extrudeEdgeSidePatches[edgei] to interprocessor patch. Adds any
1586  // interprocessor or cyclic patches if necessary.
1587  const labelList extrudeEdgeSidePatches
1588  (
1589  addExtrudeEdgeSidePatches
1590  (
1591  mesh,
1592  extrudePatch,
1593  extrudeFaces,
1594  extrudeEdgeMeshEdges,
1595  extrudeEdgeFacesMap,
1596  extrudeEdgeGlobalFaces,
1597  regionPatches
1598  )
1599  );
1600 
1601 
1602  // Set patches to use for edges to be extruded into boundary faces
1603  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1604  // In order of edgeFaces: per edge, per originating face the
1605  // patch to use for the side face (from the extruded edge).
1606  // If empty size create an internal face.
1607  labelListList extrudeEdgePatches(extrudePatch.nEdges());
1608 
1609  // Is edge a non-manifold edge
1610  PackedBoolList nonManifoldEdge(extrudePatch.nEdges(), false);
1611 
1612  // Note: logic has to be same as in countExtrudePatches.
1613  forAll(extrudePatch.edgeFaces(), edgeI)
1614  {
1615  const labelList& eFaces = extrudePatch.edgeFaces()[edgeI];
1616  labelList& ePatches = extrudeEdgePatches[edgeI];
1617 
1618  if (oneD)
1619  {
1620  ePatches.setSize(eFaces.size());
1621  forAll(eFaces, i)
1622  {
1623  ePatches[i] =
1624  zoneSidePatches[extrudeFaceZoneIDs[eFaces[i]]];
1625  }
1626 
1627  if (oneDNonManifoldEdges)
1628  {
1629  if (eFaces.size() != 2)
1630  {
1631  nonManifoldEdge[edgeI] = true;
1632  }
1633  }
1634  else
1635  {
1636  nonManifoldEdge[edgeI] = true;
1637  }
1638  }
1639  else if (eFaces.size() == 2)
1640  {
1641  // Internal edge
1642  }
1643  else if (extrudeEdgeSidePatches[edgeI] != -1)
1644  {
1645  // Coupled edge
1646  ePatches.setSize(eFaces.size());
1647  forAll(eFaces, i)
1648  {
1649  ePatches[i] = extrudeEdgeSidePatches[edgeI];
1650  }
1651  }
1652  else
1653  {
1654  // Perimeter edge
1655  const label facei = findUncoveredPatchFace
1656  (
1657  mesh,
1658  UIndirectList<label>(extrudeFaces, eFaces),
1659  extrudeEdgeMeshEdges[edgeI]
1660  );
1661 
1662  if (facei != -1)
1663  {
1664  const label newPatchi =
1665  findIndex
1666  (
1667  regionPatches,
1668  mesh.boundaryMesh()
1669  [
1670  mesh.boundaryMesh().whichPatch(facei)
1671  ].name()
1672  );
1673  ePatches.setSize(eFaces.size(), newPatchi);
1674  }
1675  else
1676  {
1677  ePatches.setSize(eFaces.size());
1678  forAll(eFaces, i)
1679  {
1680  ePatches[i] =
1681  zoneSidePatches[extrudeFaceZoneIDs[eFaces[i]]];
1682  }
1683  }
1684 
1685  nonManifoldEdge[edgeI] = true;
1686  }
1687  }
1688 
1689 
1690 
1691  // Assign point regions
1692  // ~~~~~~~~~~~~~~~~~~~~
1693 
1694  // Per face, per point the region number.
1695  faceList pointGlobalRegions;
1696  faceList pointLocalRegions;
1697  labelList localToGlobalRegion;
1698 
1700  (
1701  mesh.globalData(),
1702  extrudePatch,
1703  nonManifoldEdge,
1704  false, // keep cyclic separated regions apart
1705  pointGlobalRegions,
1706  pointLocalRegions,
1707  localToGlobalRegion
1708  );
1709 
1710  // Per local region an originating point
1711  labelList localRegionPoints(localToGlobalRegion.size());
1712  forAll(pointLocalRegions, facei)
1713  {
1714  const face& f = extrudePatch.localFaces()[facei];
1715  const face& pRegions = pointLocalRegions[facei];
1716  forAll(pRegions, fp)
1717  {
1718  localRegionPoints[pRegions[fp]] = f[fp];
1719  }
1720  }
1721 
1722  // Calculate region normals by reducing local region normals
1723  pointField localRegionNormals(localToGlobalRegion.size());
1724  {
1725  pointField localSum(localToGlobalRegion.size(), Zero);
1726 
1727  forAll(pointLocalRegions, facei)
1728  {
1729  const face& pRegions = pointLocalRegions[facei];
1730  forAll(pRegions, fp)
1731  {
1732  const label localRegionI = pRegions[fp];
1733 
1734  // Add a small amount of the face-centre-to-point vector in
1735  // order to stabilise the computation of normals on the edges
1736  // of baffles
1737  localSum[localRegionI] +=
1738  rootSmall
1739  *(
1740  extrudePatch.points()[extrudePatch[facei][fp]]
1741  - extrudePatch.faceCentres()[facei]
1742  )
1743  + (1 - rootSmall)
1744  *(
1745  extrudePatch.faceNormals()[facei]
1746  );
1747  }
1748  }
1749 
1750  Map<point> globalSum(2*localToGlobalRegion.size());
1751 
1752  forAll(localSum, localRegionI)
1753  {
1754  label globalRegionI = localToGlobalRegion[localRegionI];
1755  globalSum.insert(globalRegionI, localSum[localRegionI]);
1756  }
1757 
1758  // Reduce
1760  Pstream::mapCombineScatter(globalSum);
1761 
1762  forAll(localToGlobalRegion, localRegionI)
1763  {
1764  label globalRegionI = localToGlobalRegion[localRegionI];
1765  localRegionNormals[localRegionI] = globalSum[globalRegionI];
1766  }
1767  localRegionNormals /= mag(localRegionNormals) ;
1768  }
1769 
1770 
1771  // Use model to create displacements of first layer
1772  vectorField firstDisp(localRegionNormals.size());
1773  forAll(firstDisp, regionI)
1774  {
1775  const point& regionPt = extrudePatch.points()
1776  [
1777  extrudePatch.meshPoints()
1778  [
1779  localRegionPoints[regionI]
1780  ]
1781  ];
1782  const vector& n = localRegionNormals[regionI];
1783  firstDisp[regionI] = model()(regionPt, n, 1) - regionPt;
1784  }
1785 
1786 
1787  // Create a new mesh
1788  // ~~~~~~~~~~~~~~~~~
1789 
1790  fvMesh regionMesh
1791  (
1792  IOobject
1793  (
1794  shellRegionName,
1795  meshInstance,
1796  meshPath,
1797  runTime,
1800  false
1801  ),
1802  pointField(),
1803  faceList(),
1804  labelList(),
1805  labelList(),
1806  false
1807  );
1808 
1809  // Add the new patches
1810  forAll(regionPatches, patchi)
1811  {
1812  polyPatch* ppPtr = regionPatches[patchi];
1813  regionPatches[patchi] = ppPtr->clone(regionMesh.boundaryMesh()).ptr();
1814  delete ppPtr;
1815  }
1816  regionMesh.clearOut();
1817  regionMesh.removeFvBoundary();
1818  regionMesh.addFvPatches(regionPatches, true);
1819 
1820  // Extrude
1821  {
1822  polyTopoChange meshMod(regionPatches.size());
1823 
1824  createShellMesh extruder
1825  (
1826  extrudePatch,
1827  pointLocalRegions,
1828  localRegionPoints
1829  );
1830 
1831  extruder.setRefinement
1832  (
1833  firstDisp, // first displacement
1834  model().expansionRatio(),
1835  model().nLayers(), // nLayers
1836  extrudeFaceTopPatchID,
1837  extrudeFaceBottomPatchID,
1838  extrudeEdgePatches,
1839  meshMod
1840  );
1841 
1842  // Enforce actual point positions according to extrudeModel (model), as
1843  // the extruder only does a fixed expansionRatio. The regionPoints and
1844  // nLayers are looped in the same way as in createShellMesh.
1845  DynamicList<point>& newPoints =
1846  const_cast<DynamicList<point>&>(meshMod.points());
1847  label meshPointi = extrudePatch.localPoints().size();
1848  forAll(localRegionPoints, regioni)
1849  {
1850  const label pointi = localRegionPoints[regioni];
1851  const point& pt = extrudePatch.localPoints()[pointi];
1852  const vector& n = localRegionNormals[regioni];
1853 
1854  for (label layeri = 1; layeri <= model().nLayers(); layeri++)
1855  {
1856  newPoints[meshPointi++] = model()(pt, n, layeri);
1857  }
1858  }
1859 
1860  meshMod.changeMesh(regionMesh);
1861  }
1862 
1863  // Remove any unused patches
1864  deleteEmptyPatches(regionMesh);
1865 
1866  // Set region mesh instance and write options
1867  regionMesh.setInstance(meshInstance);
1868 
1869  Info<< "Writing mesh " << regionMesh.name() << " to "
1870  << regionMesh.facesInstance() << nl << endl;
1871 
1872  regionMesh.write();
1873 
1874 
1875  // Insert baffles into original mesh
1876  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1877 
1878  autoPtr<polyTopoChangeMap> addBafflesMap;
1879 
1880  if (adaptMesh)
1881  {
1882  polyTopoChange meshMod(mesh);
1883 
1884  // Modify faces to be in bottom (= always coupled) patch
1885  forAll(extrudeFaces, facei)
1886  {
1887  const label meshFacei = extrudeFaces[facei];
1888  const label zonei = extrudeFaceZoneIDs[facei];
1889  const bool flip = extrudeFaceFlips[facei];
1890  const face& f = mesh.faces()[meshFacei];
1891 
1892  if (!flip)
1893  {
1894  meshMod.modifyFace
1895  (
1896  f, // modified face
1897  meshFacei, // label of face being modified
1898  mesh.faceOwner()[meshFacei],// owner
1899  -1, // neighbour
1900  false, // face flip
1901  interMeshBottomPatch[zonei] // patch for face
1902  );
1903  }
1904  else if (mesh.isInternalFace(meshFacei))
1905  {
1906  meshMod.modifyFace
1907  (
1908  f.reverseFace(), // modified face
1909  meshFacei, // label of modified face
1910  mesh.faceNeighbour()[meshFacei],// owner
1911  -1, // neighbour
1912  true, // face flip
1913  interMeshBottomPatch[zonei] // patch for face
1914  );
1915  }
1916  }
1917 
1918  forAll(extrudeFaces, facei)
1919  {
1920  if (oppositeExtrudeFaces[facei] != -1)
1921  {
1922  const label meshFacei = oppositeExtrudeFaces[facei];
1923  const label zonei = oppositeExtrudeFaceZoneIDs[facei];
1924  const bool flip = oppositeExtrudeFaceFlips[facei];
1925  const face& f = mesh.faces()[meshFacei];
1926 
1927  if (!flip)
1928  {
1929  meshMod.modifyFace
1930  (
1931  f, // modified face
1932  meshFacei, // face being modified
1933  mesh.faceOwner()[meshFacei],// owner
1934  -1, // neighbour
1935  false, // face flip
1936  interMeshTopPatch[zonei] // patch for face
1937  );
1938  }
1939  else if (mesh.isInternalFace(meshFacei))
1940  {
1941  meshMod.modifyFace
1942  (
1943  f.reverseFace(), // modified face
1944  meshFacei, // label modified face
1945  mesh.faceNeighbour()[meshFacei],// owner
1946  -1, // neighbour
1947  true, // face flip
1948  interMeshTopPatch[zonei] // patch for face
1949  );
1950  }
1951  }
1952  else
1953  {
1954  const label meshFacei = extrudeFaces[facei];
1955  const label zonei = extrudeFaceZoneIDs[facei];
1956  const bool flip = extrudeFaceFlips[facei];
1957  const face& f = mesh.faces()[meshFacei];
1958 
1959  if (!flip)
1960  {
1961  if (mesh.isInternalFace(meshFacei))
1962  {
1963  meshMod.addFace
1964  (
1965  f.reverseFace(), // modified face
1966  mesh.faceNeighbour()[meshFacei],// owner
1967  -1, // neighbour
1968  meshFacei, // master face
1969  true, // flip flux
1970  interMeshTopPatch[zonei] // patch for face
1971  );
1972  }
1973  }
1974  else
1975  {
1976  meshMod.addFace
1977  (
1978  f, // face
1979  mesh.faceOwner()[meshFacei], // owner
1980  -1, // neighbour
1981  meshFacei, // master face
1982  false, // flip flux
1983  interMeshTopPatch[zonei] // patch for face
1984  );
1985  }
1986  }
1987  }
1988 
1989  // Change the mesh. Change points directly (without keeping old points).
1990  addBafflesMap = meshMod.changeMesh(mesh);
1991 
1992  // Update fields
1993  mesh.topoChange(addBafflesMap);
1994 
1995  // Remove any unused patches
1996  deleteEmptyPatches(mesh);
1997 
1998  mesh.setInstance(meshInstance);
1999 
2000  Info<< "Writing mesh " << mesh.name() << " to "
2001  << mesh.facesInstance() << nl << endl;
2002 
2003  mesh.write();
2004  }
2005 
2006  Info << "End\n" << endl;
2007 
2008  return 0;
2009 }
2010 
2011 
2012 // ************************************************************************* //
scalar y
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
wordList toc() const
Return the table of contents.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:227
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:280
const word & name() const
Return name.
Definition: IOobject.H:307
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
A HashTable to objects of type <T> with a label key.
Definition: Map.H:52
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
A bit-packed bool list.
A list of faces which address into the list of points.
label nEdges() const
Return number of edges in patch.
label nPoints() const
Return number of points supporting patch faces.
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
const Field< PointType > & points() const
Return reference to global points.
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
const labelListList & edgeFaces() const
Return edge-face addressing.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
const Field< PointType > & faceNormals() const
Return face normals for patch.
const Field< PointType > & faceCentres() const
Return face centres for patch.
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
label findIndex(const word &key) const
Return the index of the given the key or -1 if not found.
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
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
A List with indirect addressing.
Definition: UIndirectList.H:60
label size() const
Return the number of elements in the list.
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
void swap(UList< T > &)
Swap two ULists of the same type in constant time.
Definition: UList.C:90
T & last()
Return the last element of the list.
Definition: UListI.H:128
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
const word & executable() const
Name of executable without the path.
Definition: argListI.H:36
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Creates mesh by extruding a patch.
static void calcPointRegions(const globalMeshData &globalData, const primitiveFacePatch &patch, const PackedBoolList &nonManifoldEdge, const bool syncNonCollocated, faceList &pointGlobalRegions, faceList &pointLocalRegions, labelList &localToGlobalRegion)
Helper: calculate point regions. The point region is the.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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:740
T lookupOrDefaultBackwardsCompatible(const wordList &, const T &) const
Find and return a T, trying a list of keywords in sequence,.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1020
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
Class containing processor-to-processor mapping information.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:71
static autoPtr< extrudeModel > New(const dictionary &)
Select null constructed.
A list of face labels.
Definition: faceSet.H:51
Named list of face indices representing a sub-set of the mesh faces.
Definition: faceZone.H:66
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:262
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches. If validBoundary call is.
Definition: fvMeshTools.C:104
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:785
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1785
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:803
const word & name() const
Return reference to name.
Definition: fvMesh.H:434
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1369
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:313
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
Foam::polyBoundaryMesh.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
wordList names() const
Return the list of patch names.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:994
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:988
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:443
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1363
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:91
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:32
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:215
Direct mesh changes based on v1.3 polyTopoChange syntax.
const labelListList & pointEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const labelListList & edgeFaces() const
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
static word newName(const label myProcNo, const label neighbProcNo)
Return the name of a processorPolyPatch.
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh edges.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:340
#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
const fvPatchList & patches
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
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const word & regionName(const solver &region)
Definition: solver.H:218
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Istream & operator>>(Istream &, pistonPointEdgeData &)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
T clone(const T &t)
Definition: List.H:55
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
static const label labelMax
Definition: label.H:62
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion, const fileName &path=fileName::null)
Definition: systemDict.C:93
List< face > faceList
Definition: faceListFwd.H:41
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
static const char nl
Definition: Ostream.H:267
static const label labelMin
Definition: label.H:61
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
wordList patchTypes(nPatches)
wordList patchNames(nPatches)
faceListList boundary(nPatches)
labelList f(nPoints)
dictionary dict
word meshPath
Definition: setMeshPath.H:1
const bool overwrite
Definition: setNoOverwrite.H:1
Foam::argList args(argc, argv)
volScalarField & p