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-2019 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  - extrude edges in between different faceZones as a
39  (nonuniformTransform)cyclic <zoneA>_<zoneB>
40  - extrudes into master direction (i.e. away from the owner cell
41  if flipMap is false)
42 
43 \verbatim
44 
45 Internal face extrusion
46 -----------------------
47 
48  +-------------+
49  | |
50  | |
51  +---AAAAAAA---+
52  | |
53  | |
54  +-------------+
55 
56  AAA=faceZone to extrude.
57 
58 
59 For the case of no flipMap the extrusion starts at owner and extrudes
60 into the space of the neighbour:
61 
62  +CCCCCCC+
63  | | <= extruded mesh
64  +BBBBBBB+
65 
66  +-------------+
67  | |
68  | (neighbour) |
69  |___CCCCCCC___| <= original mesh (with 'baffles' added)
70  | BBBBBBB |
71  |(owner side) |
72  | |
73  +-------------+
74 
75  BBB=mapped between owner on original mesh and new extrusion.
76  (zero offset)
77  CCC=mapped between neighbour on original mesh and new extrusion
78  (offset due to the thickness of the extruded mesh)
79 
80 For the case of flipMap the extrusion is the other way around: from the
81 neighbour side into the owner side.
82 
83 
84 Boundary face extrusion
85 -----------------------
86 
87  +--AAAAAAA--+
88  | |
89  | |
90  +-----------+
91 
92  AAA=faceZone to extrude. E.g. slave side is owner side (no flipmap)
93 
94 becomes
95 
96  +CCCCCCC+
97  | | <= extruded mesh
98  +BBBBBBB+
99 
100  +--BBBBBBB--+
101  | | <= original mesh
102  | |
103  +-----------+
104 
105  BBB=mapped between original mesh and new extrusion
106  CCC=polypatch
107 
108 
109 Notes:
110  - when extruding cyclics with only one cell in between it does not
111  detect this as a cyclic since the face is the same face. It will
112  only work if the coupled edge extrudes a different face so if there
113  are more than 1 cell in between.
114 
115 \endverbatim
116 
117 \*---------------------------------------------------------------------------*/
118 
119 #include "argList.H"
120 #include "polyTopoChange.H"
121 #include "mappedWallPolyPatch.H"
122 #include "createShellMesh.H"
123 #include "syncTools.H"
124 #include "cyclicPolyPatch.H"
125 #include "wedgePolyPatch.H"
127 #include "extrudeModel.H"
128 #include "faceSet.H"
129 #include "fvMeshTools.H"
130 #include "OBJstream.H"
131 #include "PatchTools.H"
132 
133 using namespace Foam;
134 
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
137 label findPatchID(const List<polyPatch*>& newPatches, const word& name)
138 {
139  forAll(newPatches, i)
140  {
141  if (newPatches[i]->name() == name)
142  {
143  return i;
144  }
145  }
146  return -1;
147 }
148 
149 
150 template<class PatchType>
151 label addPatch
152 (
153  const polyBoundaryMesh& patches,
154  const word& patchName,
155  DynamicList<polyPatch*>& newPatches
156 )
157 {
158  label patchi = findPatchID(newPatches, patchName);
159 
160  if (patchi != -1)
161  {
162  if (isA<PatchType>(*newPatches[patchi]))
163  {
164  // Already there
165  return patchi;
166  }
167  else
168  {
170  << "Already have patch " << patchName
171  << " but of type " << newPatches[patchi]->type()
172  << exit(FatalError);
173  }
174  }
175 
176 
177  patchi = newPatches.size();
178 
179  label startFacei = 0;
180  if (patchi > 0)
181  {
182  const polyPatch& pp = *newPatches.last();
183  startFacei = pp.start()+pp.size();
184  }
185 
186 
187  newPatches.append
188  (
190  (
191  PatchType::typeName,
192  patchName,
193  0, // size
194  startFacei, // nFaces
195  patchi,
196  patches
197  ).ptr()
198  );
199 
200  return patchi;
201 }
202 
203 
204 template<class PatchType>
205 label addPatch
206 (
207  const polyBoundaryMesh& patches,
208  const word& patchName,
209  const dictionary& dict,
210  DynamicList<polyPatch*>& newPatches
211 )
212 {
213  label patchi = findPatchID(newPatches, patchName);
214 
215  if (patchi != -1)
216  {
217  if (isA<PatchType>(*newPatches[patchi]))
218  {
219  // Already there
220  return patchi;
221  }
222  else
223  {
225  << "Already have patch " << patchName
226  << " but of type " << newPatches[patchi]->type()
227  << exit(FatalError);
228  }
229  }
230 
231 
232  patchi = newPatches.size();
233 
234  label startFacei = 0;
235  if (patchi > 0)
236  {
237  const polyPatch& pp = *newPatches.last();
238  startFacei = pp.start()+pp.size();
239  }
240 
241  dictionary patchDict(dict);
242  patchDict.set("type", PatchType::typeName);
243  patchDict.set("nFaces", 0);
244  patchDict.set("startFace", startFacei);
245 
246  newPatches.append
247  (
249  (
250  patchName,
251  patchDict,
252  patchi,
253  patches
254  ).ptr()
255  );
256 
257  return patchi;
258 }
259 
260 
261 // Remove zero-sized patches
262 void deleteEmptyPatches(fvMesh& mesh)
263 {
264  const polyBoundaryMesh& patches = mesh.boundaryMesh();
265 
266  wordList masterNames;
267  if (Pstream::master())
268  {
269  masterNames = patches.names();
270  }
271  Pstream::scatter(masterNames);
272 
273 
274  labelList oldToNew(patches.size(), -1);
275  label usedI = 0;
276  label notUsedI = patches.size();
277 
278  // Add all the non-empty, non-processor patches
279  forAll(masterNames, masterI)
280  {
281  label patchi = patches.findPatchID(masterNames[masterI]);
282 
283  if (patchi != -1)
284  {
285  if (isA<processorPolyPatch>(patches[patchi]))
286  {
287  // Similar named processor patch? Not 'possible'.
288  if (patches[patchi].size() == 0)
289  {
290  Pout<< "Deleting processor patch " << patchi
291  << " name:" << patches[patchi].name()
292  << endl;
293  oldToNew[patchi] = --notUsedI;
294  }
295  else
296  {
297  oldToNew[patchi] = usedI++;
298  }
299  }
300  else
301  {
302  // Common patch.
303  if (returnReduce(patches[patchi].size(), sumOp<label>()) == 0)
304  {
305  Pout<< "Deleting patch " << patchi
306  << " name:" << patches[patchi].name()
307  << endl;
308  oldToNew[patchi] = --notUsedI;
309  }
310  else
311  {
312  oldToNew[patchi] = usedI++;
313  }
314  }
315  }
316  }
317 
318  // Add remaining patches at the end
319  forAll(patches, patchi)
320  {
321  if (oldToNew[patchi] == -1)
322  {
323  // Unique to this processor. Note: could check that these are
324  // only processor patches.
325  if (patches[patchi].size() == 0)
326  {
327  Pout<< "Deleting processor patch " << patchi
328  << " name:" << patches[patchi].name()
329  << endl;
330  oldToNew[patchi] = --notUsedI;
331  }
332  else
333  {
334  oldToNew[patchi] = usedI++;
335  }
336  }
337  }
338 
339  fvMeshTools::reorderPatches(mesh, oldToNew, usedI, true);
340 }
341 
342 
343 void createDummyFvMeshFiles(const polyMesh& mesh, const word& regionName)
344 {
345  // Create dummy system/fv*
346  {
347  IOobject io
348  (
349  "fvSchemes",
350  mesh.time().system(),
351  regionName,
352  mesh,
355  false
356  );
357 
358  Info<< "Testing:" << io.objectPath() << endl;
359 
360  if (!io.typeHeaderOk<IOdictionary>(true))
361  {
362  Info<< "Writing dummy " << regionName/io.name() << endl;
363  dictionary dummyDict;
364  dictionary divDict;
365  dummyDict.add("divSchemes", divDict);
366  dictionary gradDict;
367  dummyDict.add("gradSchemes", gradDict);
368  dictionary laplDict;
369  dummyDict.add("laplacianSchemes", laplDict);
370 
371  IOdictionary(io, dummyDict).regIOobject::write();
372  }
373  }
374  {
375  IOobject io
376  (
377  "fvSolution",
378  mesh.time().system(),
379  regionName,
380  mesh,
383  false
384  );
385 
386  if (!io.typeHeaderOk<IOdictionary>(true))
387  {
388  Info<< "Writing dummy " << regionName/io.name() << endl;
389  dictionary dummyDict;
390  IOdictionary(io, dummyDict).regIOobject::write();
391  }
392  }
393 }
394 
395 
396 // Check zone either all internal or all external faces
397 void checkZoneInside
398 (
399  const polyMesh& mesh,
400  const wordList& zoneNames,
401  const labelList& zoneID,
402  const labelList& extrudeMeshFaces,
403  const boolList& isInternal
404 )
405 {
406  forAll(zoneNames, i)
407  {
408  if (isInternal[i])
409  {
410  Info<< "Zone " << zoneNames[i] << " has internal faces" << endl;
411  }
412  else
413  {
414  Info<< "Zone " << zoneNames[i] << " has boundary faces" << endl;
415  }
416  }
417 
418  forAll(extrudeMeshFaces, i)
419  {
420  label facei = extrudeMeshFaces[i];
421  label zoneI = zoneID[i];
422  if (isInternal[zoneI] != mesh.isInternalFace(facei))
423  {
425  << "Zone " << zoneNames[zoneI]
426  << " is not consistently all internal or all boundary faces."
427  << " Face " << facei << " at " << mesh.faceCentres()[facei]
428  << " is the first occurrence."
429  << exit(FatalError);
430  }
431  }
432 }
433 
434 
435 // To combineReduce a labelList. Filters out duplicates.
436 class uniqueEqOp
437 {
438 
439 public:
440 
441  void operator()(labelList& x, const labelList& y) const
442  {
443  if (x.empty())
444  {
445  if (y.size())
446  {
447  x = y;
448  }
449  }
450  else
451  {
452  forAll(y, yi)
453  {
454  if (findIndex(x, y[yi]) == -1)
455  {
456  label sz = x.size();
457  x.setSize(sz+1);
458  x[sz] = y[yi];
459  }
460  }
461  }
462  }
463 };
464 
465 
466 // Calculate global pp faces per pp edge.
467 labelListList globalEdgeFaces
468 (
469  const polyMesh& mesh,
470  const globalIndex& globalFaces,
471  const primitiveFacePatch& pp,
472  const labelList& ppMeshEdges
473 )
474 {
475  // From mesh edge to global pp face labels.
476  labelListList globalEdgeFaces(ppMeshEdges.size());
477 
478  const labelListList& edgeFaces = pp.edgeFaces();
479 
480  forAll(edgeFaces, edgeI)
481  {
482  const labelList& eFaces = edgeFaces[edgeI];
483 
484  // Store pp face and processor as unique tag.
485  labelList& globalEFaces = globalEdgeFaces[edgeI];
486  globalEFaces.setSize(eFaces.size());
487  forAll(eFaces, i)
488  {
489  globalEFaces[i] = globalFaces.toGlobal(eFaces[i]);
490  }
491  }
492 
493  // Synchronise across coupled edges.
495  (
496  mesh,
497  ppMeshEdges,
498  globalEdgeFaces,
499  uniqueEqOp(),
500  labelList() // null value
501  );
502 
503  return globalEdgeFaces;
504 }
505 
506 
507 // Find a patch face that is not extruded. Return -1 if not found.
508 label findUncoveredPatchFace
509 (
510  const fvMesh& mesh,
511  const UIndirectList<label>& extrudeMeshFaces,// mesh faces that are extruded
512  const label meshEdgeI // mesh edge
513 )
514 {
515  // Make set of extruded faces.
516  labelHashSet extrudeFaceSet(extrudeMeshFaces.size());
517  forAll(extrudeMeshFaces, i)
518  {
519  extrudeFaceSet.insert(extrudeMeshFaces[i]);
520  }
521 
522  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
523  const labelList& eFaces = mesh.edgeFaces()[meshEdgeI];
524  forAll(eFaces, i)
525  {
526  label facei = eFaces[i];
527  label patchi = pbm.whichPatch(facei);
528 
529  if
530  (
531  patchi != -1
532  && !pbm[patchi].coupled()
533  && !extrudeFaceSet.found(facei)
534  )
535  {
536  return facei;
537  }
538  }
539  return -1;
540 }
541 
542 
543 // Same as findUncoveredPatchFace, except explicitly checks for cyclic faces
544 label findUncoveredCyclicPatchFace
545 (
546  const fvMesh& mesh,
547  const UIndirectList<label>& extrudeMeshFaces,// mesh faces that are extruded
548  const label meshEdgeI // mesh edge
549 )
550 {
551  // Make set of extruded faces.
552  labelHashSet extrudeFaceSet(extrudeMeshFaces.size());
553  forAll(extrudeMeshFaces, i)
554  {
555  extrudeFaceSet.insert(extrudeMeshFaces[i]);
556  }
557 
558  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
559  const labelList& eFaces = mesh.edgeFaces()[meshEdgeI];
560  forAll(eFaces, i)
561  {
562  label facei = eFaces[i];
563  label patchi = pbm.whichPatch(facei);
564 
565  if
566  (
567  patchi != -1
568  && isA<cyclicPolyPatch>(pbm[patchi])
569  && !extrudeFaceSet.found(facei)
570  )
571  {
572  return facei;
573  }
574  }
575  return -1;
576 }
577 
578 
579 // Calculate per edge min and max zone
580 void calcEdgeMinMaxZone
581 (
582  const fvMesh& mesh,
583  const primitiveFacePatch& extrudePatch,
584  const labelList& extrudeMeshEdges,
585  const labelList& zoneID,
586  const mapDistribute& extrudeEdgeFacesMap,
587  const labelListList& extrudeEdgeGlobalFaces,
588 
589  labelList& minZoneID,
590  labelList& maxZoneID
591 )
592 {
593  // Get zoneIDs in extrudeEdgeGlobalFaces order
594  labelList mappedZoneID(zoneID);
595  extrudeEdgeFacesMap.distribute(mappedZoneID);
596 
597  // Get min and max zone per edge
598  minZoneID.setSize(extrudeEdgeGlobalFaces.size(), labelMax);
599  maxZoneID.setSize(extrudeEdgeGlobalFaces.size(), labelMin);
600 
601  forAll(extrudeEdgeGlobalFaces, edgeI)
602  {
603  const labelList& eFaces = extrudeEdgeGlobalFaces[edgeI];
604  if (eFaces.size())
605  {
606  forAll(eFaces, i)
607  {
608  label zoneI = mappedZoneID[eFaces[i]];
609  minZoneID[edgeI] = min(minZoneID[edgeI], zoneI);
610  maxZoneID[edgeI] = max(maxZoneID[edgeI], zoneI);
611  }
612  }
613  }
615  (
616  mesh,
617  extrudeMeshEdges,
618  minZoneID,
619  minEqOp<label>(),
620  labelMax // null value
621  );
623  (
624  mesh,
625  extrudeMeshEdges,
626  maxZoneID,
627  maxEqOp<label>(),
628  labelMin // null value
629  );
630 }
631 
632 
633 // Count the number of faces in patches that need to be created. Calculates:
634 // zoneSidePatch[zoneI] : the number of side faces to be created
635 // zoneZonePatch[zoneA,zoneB] : the number of faces in between zoneA and B
636 // Since this only counts we're not taking the processor patches into
637 // account.
638 void countExtrudePatches
639 (
640  const fvMesh& mesh,
641  const label nZones,
642  const primitiveFacePatch& extrudePatch,
643  const labelList& extrudeMeshFaces,
644  const labelList& extrudeMeshEdges,
645 
646  const labelListList& extrudeEdgeGlobalFaces,
647  const labelList& minZoneID,
648  const labelList& maxZoneID,
649 
650  labelList& zoneSidePatch,
651  labelList& zoneZonePatch
652 )
653 {
654  // Check on master edge for use of zones. Since we only want to know
655  // whether they are being used at all no need to accurately count on slave
656  // edge as well. Just add all together at the end of this routine so it
657  // gets detected at least.
658 
659  forAll(extrudePatch.edgeFaces(), edgeI)
660  {
661  const labelList& eFaces = extrudePatch.edgeFaces()[edgeI];
662 
663  if (eFaces.size() == 2)
664  {
665  // Internal edge - check if in between different zones.
666  if (minZoneID[edgeI] != maxZoneID[edgeI])
667  {
668  zoneZonePatch[minZoneID[edgeI]*nZones+maxZoneID[edgeI]]++;
669  }
670  }
671  else if
672  (
673  eFaces.size() == 1
674  && extrudeEdgeGlobalFaces[edgeI].size() == 2
675  )
676  {
677  // Coupled edge - check if in between different zones.
678  if (minZoneID[edgeI] != maxZoneID[edgeI])
679  {
680  const edge& e = extrudePatch.edges()[edgeI];
681  const pointField& pts = extrudePatch.localPoints();
683  << "Edge " << edgeI
684  << "at " << pts[e[0]] << pts[e[1]]
685  << " is a coupled edge and in between two different zones "
686  << minZoneID[edgeI] << " and " << maxZoneID[edgeI] << endl
687  << " This is currently not supported." << endl;
688 
689  zoneZonePatch[minZoneID[edgeI]*nZones+maxZoneID[edgeI]]++;
690  }
691  }
692  else
693  {
694  // One or more than two edge-faces.
695  // Check whether we are on a mesh edge with external patches. If
696  // so choose any uncovered one. If none found put face in
697  // undetermined zone 'side' patch
698 
699  label facei = findUncoveredPatchFace
700  (
701  mesh,
702  UIndirectList<label>(extrudeMeshFaces, eFaces),
703  extrudeMeshEdges[edgeI]
704  );
705 
706  if (facei == -1)
707  {
708  zoneSidePatch[minZoneID[edgeI]]++;
709  }
710  }
711  }
712  // Synchronise decistion. Actual numbers are not important, just make
713  // sure that they're > 0 on all processors.
715  Pstream::listCombineScatter(zoneSidePatch);
717  Pstream::listCombineScatter(zoneZonePatch);
718 }
719 
720 
721 void addCouplingPatches
722 (
723  const fvMesh& mesh,
724  const word& regionName,
725  const word& shellRegionName,
726  const wordList& zoneNames,
727  const wordList& zoneShadowNames,
728  const boolList& isInternal,
729  const labelList& zoneIDs,
730 
731  DynamicList<polyPatch*>& newPatches,
732  labelList& interRegionTopPatch,
733  labelList& interRegionBottomPatch
734 )
735 {
736  Pout<< "Adding coupling patches:" << nl << nl
737  << "patchID\tpatch\ttype" << nl
738  << "-------\t-----\t----"
739  << endl;
740 
741  interRegionTopPatch.setSize(zoneNames.size(), -1);
742  interRegionBottomPatch.setSize(zoneNames.size(), -1);
743 
744  label nOldPatches = newPatches.size();
745  forAll(zoneNames, zoneI)
746  {
747  word interName
748  (
749  regionName
750  +"_to_"
751  +shellRegionName
752  +'_'
753  +zoneNames[zoneI]
754  );
755 
756  if (isInternal[zoneI])
757  {
758  interRegionTopPatch[zoneI] = addPatch<mappedWallPolyPatch>
759  (
760  mesh.boundaryMesh(),
761  interName + "_top",
762  newPatches
763  );
764  Pout<< interRegionTopPatch[zoneI]
765  << '\t' << newPatches[interRegionTopPatch[zoneI]]->name()
766  << '\t' << newPatches[interRegionTopPatch[zoneI]]->type()
767  << nl;
768 
769  interRegionBottomPatch[zoneI] = addPatch<mappedWallPolyPatch>
770  (
771  mesh.boundaryMesh(),
772  interName + "_bottom",
773  newPatches
774  );
775  Pout<< interRegionBottomPatch[zoneI]
776  << '\t' << newPatches[interRegionBottomPatch[zoneI]]->name()
777  << '\t' << newPatches[interRegionBottomPatch[zoneI]]->type()
778  << nl;
779  }
780  else if (zoneShadowNames.size() == 0)
781  {
782  interRegionTopPatch[zoneI] = addPatch<polyPatch>
783  (
784  mesh.boundaryMesh(),
785  zoneNames[zoneI] + "_top",
786  newPatches
787  );
788  Pout<< interRegionTopPatch[zoneI]
789  << '\t' << newPatches[interRegionTopPatch[zoneI]]->name()
790  << '\t' << newPatches[interRegionTopPatch[zoneI]]->type()
791  << nl;
792 
793  interRegionBottomPatch[zoneI] = addPatch<mappedWallPolyPatch>
794  (
795  mesh.boundaryMesh(),
796  interName,
797  newPatches
798  );
799  Pout<< interRegionBottomPatch[zoneI]
800  << '\t' << newPatches[interRegionBottomPatch[zoneI]]->name()
801  << '\t' << newPatches[interRegionBottomPatch[zoneI]]->type()
802  << nl;
803  }
804  else // patch using shadow face zones.
805  {
806  interRegionTopPatch[zoneI] = addPatch<mappedWallPolyPatch>
807  (
808  mesh.boundaryMesh(),
809  zoneShadowNames[zoneI] + "_top",
810  newPatches
811  );
812  Pout<< interRegionTopPatch[zoneI]
813  << '\t' << newPatches[interRegionTopPatch[zoneI]]->name()
814  << '\t' << newPatches[interRegionTopPatch[zoneI]]->type()
815  << nl;
816 
817  interRegionBottomPatch[zoneI] = addPatch<mappedWallPolyPatch>
818  (
819  mesh.boundaryMesh(),
820  interName,
821  newPatches
822  );
823  Pout<< interRegionBottomPatch[zoneI]
824  << '\t' << newPatches[interRegionBottomPatch[zoneI]]->name()
825  << '\t' << newPatches[interRegionBottomPatch[zoneI]]->type()
826  << nl;
827  }
828  }
829  Pout<< "Added " << newPatches.size()-nOldPatches
830  << " inter-region patches." << nl
831  << endl;
832 }
833 
834 
835 // Sets sidePatch[edgeI] to interprocessor or cyclic patch. Adds any
836 // coupled patches if necessary.
837 void addCoupledPatches
838 (
839  const fvMesh& mesh,
840  const primitiveFacePatch& extrudePatch,
841  const labelList& extrudeMeshFaces,
842  const labelList& extrudeMeshEdges,
843  const mapDistribute& extrudeEdgeFacesMap,
844  const labelListList& extrudeEdgeGlobalFaces,
845 
846  labelList& sidePatchID,
847  DynamicList<polyPatch*>& newPatches
848 )
849 {
850  // Calculate opposite processor for coupled edges (only if shared by
851  // two procs). Note: could have saved original globalEdgeFaces structure.
852 
853  // Get procID in extrudeEdgeGlobalFaces order
854  labelList procID(extrudeEdgeGlobalFaces.size(), Pstream::myProcNo());
855  extrudeEdgeFacesMap.distribute(procID);
856 
857  labelList minProcID(extrudeEdgeGlobalFaces.size(), labelMax);
858  labelList maxProcID(extrudeEdgeGlobalFaces.size(), labelMin);
859 
860  forAll(extrudeEdgeGlobalFaces, edgeI)
861  {
862  const labelList& eFaces = extrudeEdgeGlobalFaces[edgeI];
863  if (eFaces.size())
864  {
865  forAll(eFaces, i)
866  {
867  label proci = procID[eFaces[i]];
868  minProcID[edgeI] = min(minProcID[edgeI], proci);
869  maxProcID[edgeI] = max(maxProcID[edgeI], proci);
870  }
871  }
872  }
874  (
875  mesh,
876  extrudeMeshEdges,
877  minProcID,
878  minEqOp<label>(),
879  labelMax // null value
880  );
882  (
883  mesh,
884  extrudeMeshEdges,
885  maxProcID,
886  maxEqOp<label>(),
887  labelMin // null value
888  );
889 
890  Pout<< "Adding processor or cyclic patches:" << nl << nl
891  << "patchID\tpatch" << nl
892  << "-------\t-----"
893  << endl;
894 
895  label nOldPatches = newPatches.size();
896 
897  sidePatchID.setSize(extrudePatch.edgeFaces().size(), -1);
898  forAll(extrudePatch.edgeFaces(), edgeI)
899  {
900  const labelList& eFaces = extrudePatch.edgeFaces()[edgeI];
901 
902  if
903  (
904  eFaces.size() == 1
905  && extrudeEdgeGlobalFaces[edgeI].size() == 2
906  )
907  {
908  // coupled boundary edge. Find matching patch.
909  label nbrProci = minProcID[edgeI];
910  if (nbrProci == Pstream::myProcNo())
911  {
912  nbrProci = maxProcID[edgeI];
913  }
914 
915 
916  if (nbrProci == Pstream::myProcNo())
917  {
918  // Cyclic patch since both procs the same. This cyclic should
919  // already exist in newPatches so no adding necessary.
920 
921  label facei = findUncoveredCyclicPatchFace
922  (
923  mesh,
924  UIndirectList<label>(extrudeMeshFaces, eFaces),
925  extrudeMeshEdges[edgeI]
926  );
927 
928  if (facei != -1)
929  {
930  const polyBoundaryMesh& patches = mesh.boundaryMesh();
931 
932  label newPatchi = findPatchID
933  (
934  newPatches,
935  patches[patches.whichPatch(facei)].name()
936  );
937 
938  sidePatchID[edgeI] = newPatchi;
939  }
940  else
941  {
943  << "Unable to determine coupled patch addressing"
944  << abort(FatalError);
945  }
946  }
947  else
948  {
949  // Processor patch
950  word name
951  (
953  );
954 
955  sidePatchID[edgeI] = findPatchID(newPatches, name);
956 
957  if (sidePatchID[edgeI] == -1)
958  {
959  dictionary patchDict;
960  patchDict.add("myProcNo", Pstream::myProcNo());
961  patchDict.add("neighbProcNo", nbrProci);
962 
963  sidePatchID[edgeI] = addPatch<processorPolyPatch>
964  (
965  mesh.boundaryMesh(),
966  name,
967  patchDict,
968  newPatches
969  );
970 
971  Pout<< sidePatchID[edgeI] << '\t' << name
972  << nl;
973  }
974  }
975  }
976  }
977  Pout<< "Added " << newPatches.size()-nOldPatches
978  << " coupled patches." << nl
979  << endl;
980 }
981 
982 
983 void addZoneSidePatches
984 (
985  const fvMesh& mesh,
986  const wordList& zoneNames,
987  const word& oneDPolyPatchType,
988 
989  DynamicList<polyPatch*>& newPatches,
990  labelList& zoneSidePatch
991 )
992 {
993  Pout<< "Adding patches for sides on zones:" << nl << nl
994  << "patchID\tpatch" << nl
995  << "-------\t-----"
996  << endl;
997 
998  label nOldPatches = newPatches.size();
999 
1000  forAll(zoneSidePatch, zoneI)
1001  {
1002  if (oneDPolyPatchType != word::null)
1003  {
1004  // Reuse single empty patch.
1005  word patchName;
1006  if (oneDPolyPatchType == "empty")
1007  {
1008  patchName = "oneDEmptyPatch";
1009  zoneSidePatch[zoneI] = addPatch<emptyPolyPatch>
1010  (
1011  mesh.boundaryMesh(),
1012  patchName,
1013  newPatches
1014  );
1015  }
1016  else if (oneDPolyPatchType == "wedge")
1017  {
1018  patchName = "oneDWedgePatch";
1019  zoneSidePatch[zoneI] = addPatch<wedgePolyPatch>
1020  (
1021  mesh.boundaryMesh(),
1022  patchName,
1023  newPatches
1024  );
1025  }
1026  else
1027  {
1029  << "Type " << oneDPolyPatchType << " does not exist "
1030  << exit(FatalError);
1031  }
1032 
1033  Pout<< zoneSidePatch[zoneI] << '\t' << patchName << nl;
1034  }
1035  else if (zoneSidePatch[zoneI] > 0)
1036  {
1037  word patchName = zoneNames[zoneI] + "_" + "side";
1038 
1039  zoneSidePatch[zoneI] = addPatch<polyPatch>
1040  (
1041  mesh.boundaryMesh(),
1042  patchName,
1043  newPatches
1044  );
1045 
1046  Pout<< zoneSidePatch[zoneI] << '\t' << patchName << nl;
1047  }
1048  }
1049  Pout<< "Added " << newPatches.size()-nOldPatches << " zone-side patches."
1050  << nl << endl;
1051 }
1052 
1053 
1054 void addInterZonePatches
1055 (
1056  const fvMesh& mesh,
1057  const wordList& zoneNames,
1058  const bool oneD,
1059 
1060  labelList& zoneZonePatch_min,
1061  labelList& zoneZonePatch_max,
1062  DynamicList<polyPatch*>& newPatches
1063 )
1064 {
1065  Pout<< "Adding inter-zone patches:" << nl << nl
1066  << "patchID\tpatch" << nl
1067  << "-------\t-----"
1068  << endl;
1069 
1070  dictionary transformDict;
1071  transformDict.add
1072  (
1073  "transform",
1074  cyclicPolyPatch::transformTypeNames[cyclicPolyPatch::NOORDERING]
1075  );
1076 
1077  // Allow arbitrary matching for nonuniformTransformCyclicPolyPatch
1078  transformDict.add("matchTolerance", 1e6);
1079 
1080  label nOldPatches = newPatches.size();
1081 
1082  if (!oneD)
1083  {
1084  forAll(zoneZonePatch_min, minZone)
1085  {
1086  for (label maxZone = minZone; maxZone < zoneNames.size(); maxZone++)
1087  {
1088  label index = minZone*zoneNames.size()+maxZone;
1089 
1090  if (zoneZonePatch_min[index] > 0)
1091  {
1092  word minToMax =
1093  zoneNames[minZone]
1094  + "_to_"
1095  + zoneNames[maxZone];
1096  word maxToMin =
1097  zoneNames[maxZone]
1098  + "_to_"
1099  + zoneNames[minZone];
1100 
1101  {
1102  transformDict.set("neighbourPatch", maxToMin);
1103  zoneZonePatch_min[index] =
1104  addPatch<nonuniformTransformCyclicPolyPatch>
1105  (
1106  mesh.boundaryMesh(),
1107  minToMax,
1108  transformDict,
1109  newPatches
1110  );
1111  Pout<< zoneZonePatch_min[index] << '\t' << minToMax
1112  << nl;
1113  }
1114  {
1115  transformDict.set("neighbourPatch", minToMax);
1116  zoneZonePatch_max[index] =
1117  addPatch<nonuniformTransformCyclicPolyPatch>
1118  (
1119  mesh.boundaryMesh(),
1120  maxToMin,
1121  transformDict,
1122  newPatches
1123  );
1124  Pout<< zoneZonePatch_max[index] << '\t' << maxToMin
1125  << nl;
1126  }
1127 
1128  }
1129  }
1130  }
1131  }
1132  Pout<< "Added " << newPatches.size()-nOldPatches << " inter-zone patches."
1133  << nl << endl;
1134 }
1135 
1136 
1137 tmp<pointField> calcOffset
1138 (
1139  const primitiveFacePatch& extrudePatch,
1140  const createShellMesh& extruder,
1141  const polyPatch& pp
1142 )
1143 {
1145 
1146  tmp<pointField> toffsets(new pointField(fc.size()));
1147  pointField& offsets = toffsets.ref();
1148 
1149  forAll(fc, i)
1150  {
1151  label meshFacei = pp.start()+i;
1152  label patchFacei = mag(extruder.faceToFaceMap()[meshFacei])-1;
1153  point patchFc = extrudePatch[patchFacei].centre
1154  (
1155  extrudePatch.points()
1156  );
1157  offsets[i] = patchFc - fc[i];
1158  }
1159  return toffsets;
1160 }
1161 
1162 
1163 void setCouplingInfo
1164 (
1165  fvMesh& mesh,
1166  const labelList& zoneToPatch,
1167  const word& sampleRegion,
1169  const List<pointField>& offsets
1170 )
1171 {
1172  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1173 
1174  List<polyPatch*> newPatches
1175  (
1176  patches.size(),
1177  static_cast<polyPatch*>(nullptr)
1178  );
1179 
1180  forAll(zoneToPatch, zoneI)
1181  {
1182  label patchi = zoneToPatch[zoneI];
1183 
1184  if (patchi != -1)
1185  {
1186  const polyPatch& pp = patches[patchi];
1187 
1188  if (isA<mappedWallPolyPatch>(pp))
1189  {
1190  newPatches[patchi] = new mappedWallPolyPatch
1191  (
1192  pp.name(),
1193  pp.size(),
1194  pp.start(),
1195  patchi,
1196  sampleRegion, // sampleRegion
1197  mode, // sampleMode
1198  pp.name(), // samplePatch
1199  offsets[zoneI], // offset
1200  patches
1201  );
1202  }
1203  }
1204  }
1205 
1206  forAll(newPatches, patchi)
1207  {
1208  if (!newPatches[patchi])
1209  {
1210  newPatches[patchi] = patches[patchi].clone(patches).ptr();
1211  }
1212  }
1213 
1214  mesh.removeFvBoundary();
1215  mesh.addFvPatches(newPatches, true);
1216 }
1217 
1218 
1219 // Extrude and write geometric properties
1220 void extrudeGeometricProperties
1221 (
1222  const polyMesh& mesh,
1223  const primitiveFacePatch& extrudePatch,
1224  const createShellMesh& extruder,
1225  const polyMesh& regionMesh,
1226  const extrudeModel& model
1227 )
1228 {
1229  const pointIOField patchFaceCentres
1230  (
1231  IOobject
1232  (
1233  "patchFaceCentres",
1234  mesh.pointsInstance(),
1235  mesh.meshSubDir,
1236  mesh,
1238  )
1239  );
1240 
1241  const pointIOField patchEdgeCentres
1242  (
1243  IOobject
1244  (
1245  "patchEdgeCentres",
1246  mesh.pointsInstance(),
1247  mesh.meshSubDir,
1248  mesh,
1250  )
1251  );
1252 
1253  // forAll(extrudePatch.edges(), edgeI)
1254  //{
1255  // const edge& e = extrudePatch.edges()[edgeI];
1256  // Pout<< "Edge:" << e.centre(extrudePatch.localPoints()) << nl
1257  // << "read:" << patchEdgeCentres[edgeI]
1258  // << endl;
1259  //}
1260 
1261 
1262  // Determine edge normals on original patch
1263  labelList patchEdges;
1264  labelList coupledEdges;
1265  PackedBoolList sameEdgeOrientation;
1267  (
1268  extrudePatch,
1269  mesh.globalData().coupledPatch(),
1270  patchEdges,
1271  coupledEdges,
1272  sameEdgeOrientation
1273  );
1274 
1275  pointField patchEdgeNormals
1276  (
1278  (
1279  mesh,
1280  extrudePatch,
1281  patchEdges,
1282  coupledEdges
1283  )
1284  );
1285 
1286 
1287  pointIOField faceCentres
1288  (
1289  IOobject
1290  (
1291  "faceCentres",
1292  regionMesh.pointsInstance(),
1293  regionMesh.meshSubDir,
1294  regionMesh,
1297  false
1298  ),
1299  regionMesh.nFaces()
1300  );
1301 
1302 
1303  // Work out layers. Guaranteed in columns so no fancy parallel bits.
1304 
1305 
1306  forAll(extruder.faceToFaceMap(), facei)
1307  {
1308  if (extruder.faceToFaceMap()[facei] != 0)
1309  {
1310  // 'horizontal' face
1311  label patchFacei = mag(extruder.faceToFaceMap()[facei])-1;
1312 
1313  label celli = regionMesh.faceOwner()[facei];
1314  if (regionMesh.isInternalFace(facei))
1315  {
1316  celli = max(celli, regionMesh.faceNeighbour()[facei]);
1317  }
1318 
1319  // Calculate layer from cell numbering (see createShellMesh)
1320  label layerI = (celli % model.nLayers());
1321 
1322  if
1323  (
1324  !regionMesh.isInternalFace(facei)
1325  && extruder.faceToFaceMap()[facei] > 0
1326  )
1327  {
1328  // Top face
1329  layerI++;
1330  }
1331 
1332 
1333  // Recalculate based on extrusion model
1334  faceCentres[facei] = model
1335  (
1336  patchFaceCentres[patchFacei],
1337  extrudePatch.faceNormals()[patchFacei],
1338  layerI
1339  );
1340  }
1341  else
1342  {
1343  // 'vertical face
1344  label patchEdgeI = extruder.faceToEdgeMap()[facei];
1345  label layerI =
1346  (
1347  regionMesh.faceOwner()[facei]
1348  % model.nLayers()
1349  );
1350 
1351  // Extrude patch edge centre to this layer
1352  point pt0 = model
1353  (
1354  patchEdgeCentres[patchEdgeI],
1355  patchEdgeNormals[patchEdgeI],
1356  layerI
1357  );
1358  // Extrude patch edge centre to next layer
1359  point pt1 = model
1360  (
1361  patchEdgeCentres[patchEdgeI],
1362  patchEdgeNormals[patchEdgeI],
1363  layerI+1
1364  );
1365 
1366  // Interpolate
1367  faceCentres[facei] = 0.5*(pt0+pt1);
1368  }
1369  }
1370 
1371  pointIOField cellCentres
1372  (
1373  IOobject
1374  (
1375  "cellCentres",
1376  regionMesh.pointsInstance(),
1377  regionMesh.meshSubDir,
1378  regionMesh,
1381  false
1382  ),
1383  regionMesh.nCells()
1384  );
1385 
1386  forAll(extruder.cellToFaceMap(), celli)
1387  {
1388  label patchFacei = extruder.cellToFaceMap()[celli];
1389 
1390  // Calculate layer from cell numbering (see createShellMesh)
1391  label layerI = (celli % model.nLayers());
1392 
1393  // Recalculate based on extrusion model
1394  point pt0 = model
1395  (
1396  patchFaceCentres[patchFacei],
1397  extrudePatch.faceNormals()[patchFacei],
1398  layerI
1399  );
1400  point pt1 = model
1401  (
1402  patchFaceCentres[patchFacei],
1403  extrudePatch.faceNormals()[patchFacei],
1404  layerI+1
1405  );
1406 
1407  // Interpolate
1408  cellCentres[celli] = 0.5*(pt0+pt1);
1409  }
1410 
1411 
1412  // Bit of checking
1413  if (false)
1414  {
1415  OBJstream faceStr(regionMesh.time().path()/"faceCentres.obj");
1416  OBJstream cellStr(regionMesh.time().path()/"cellCentres.obj");
1417 
1418  forAll(faceCentres, facei)
1419  {
1420  Pout<< "Model :" << faceCentres[facei] << endl
1421  << "regionMesh:" << regionMesh.faceCentres()[facei] << endl;
1422  faceStr.write
1423  (
1424  linePointRef
1425  (
1426  faceCentres[facei],
1427  regionMesh.faceCentres()[facei]
1428  )
1429  );
1430  }
1431  forAll(cellCentres, celli)
1432  {
1433  Pout<< "Model :" << cellCentres[celli] << endl
1434  << "regionMesh:" << regionMesh.cellCentres()[celli] << endl;
1435  cellStr.write
1436  (
1437  linePointRef
1438  (
1439  cellCentres[celli],
1440  regionMesh.cellCentres()[celli]
1441  )
1442  );
1443  }
1444  }
1445 
1446 
1447 
1448  Info<< "Writing geometric properties for mesh " << regionMesh.name()
1449  << " to " << regionMesh.pointsInstance() << nl
1450  << endl;
1451 
1452  bool ok = faceCentres.write() && cellCentres.write();
1453 
1454  if (!ok)
1455  {
1457  << "Failed writing " << faceCentres.objectPath()
1458  << " and " << cellCentres.objectPath()
1459  << exit(FatalError);
1460  }
1461 }
1462 
1463 
1464 
1465 
1466 int main(int argc, char *argv[])
1467 {
1468  argList::addNote("Create region mesh by extruding a faceZone or faceSet");
1469 
1470  #include "addRegionOption.H"
1471  #include "addOverwriteOption.H"
1472  #include "addDictOption.H"
1473  #include "setRootCase.H"
1474  #include "createTime.H"
1475  #include "createNamedMesh.H"
1476 
1477  if (mesh.boundaryMesh().checkParallelSync(true))
1478  {
1479  List<wordList> allNames(Pstream::nProcs());
1480  allNames[Pstream::myProcNo()] = mesh.boundaryMesh().names();
1481  Pstream::gatherList(allNames);
1482  Pstream::scatterList(allNames);
1483 
1485  << "Patches are not synchronised on all processors."
1486  << " Per processor patches " << allNames
1487  << exit(FatalError);
1488  }
1489 
1490 
1491  const word oldInstance = mesh.pointsInstance();
1492  bool overwrite = args.optionFound("overwrite");
1493 
1494 
1495  const word dictName("extrudeToRegionMeshDict");
1496 
1497  #include "setSystemMeshDictionaryIO.H"
1498 
1500 
1501 
1502  // Point generator
1504 
1505 
1506  // Region
1507  const word shellRegionName(dict.lookup("region"));
1508 
1509  // Faces to extrude - either faceZones or faceSets (boundary faces only)
1510  wordList zoneNames;
1511  wordList zoneShadowNames;
1512 
1513  bool hasZones = dict.found("faceZones");
1514  if (hasZones)
1515  {
1516  dict.lookup("faceZones") >> zoneNames;
1517  dict.readIfPresent("faceZonesShadow", zoneShadowNames);
1518 
1519  // Check
1520  if (dict.found("faceSets"))
1521  {
1522  FatalIOErrorIn(args.executable().c_str(), dict)
1523  << "Please supply faces to extrude either through 'faceZones'"
1524  << " or 'faceSets' entry. Found both."
1525  << exit(FatalIOError);
1526  }
1527  }
1528  else
1529  {
1530  dict.lookup("faceSets") >> zoneNames;
1531  dict.readIfPresent("faceSetsShadow", zoneShadowNames);
1532  }
1533 
1534 
1535  mappedPatchBase::sampleMode sampleMode =
1536  mappedPatchBase::sampleModeNames_[dict.lookup("sampleMode")];
1537 
1538  const Switch oneD(dict.lookup("oneD"));
1539  Switch oneDNonManifoldEdges(false);
1540  word oneDPatchType(emptyPolyPatch::typeName);
1541  if (oneD)
1542  {
1543  oneDNonManifoldEdges = dict.lookupOrDefault("nonManifold", false);
1544  dict.lookup("oneDPolyPatchType") >> oneDPatchType;
1545  }
1546 
1547  const Switch adaptMesh(dict.lookup("adaptMesh"));
1548 
1549  if (hasZones)
1550  {
1551  Info<< "Extruding zones " << zoneNames
1552  << " on mesh " << regionName
1553  << " into shell mesh " << shellRegionName
1554  << endl;
1555  }
1556  else
1557  {
1558  Info<< "Extruding faceSets " << zoneNames
1559  << " on mesh " << regionName
1560  << " into shell mesh " << shellRegionName
1561  << endl;
1562  }
1563 
1564  if (shellRegionName == regionName)
1565  {
1566  FatalIOErrorIn(args.executable().c_str(), dict)
1567  << "Cannot extrude into same region as mesh." << endl
1568  << "Mesh region : " << regionName << endl
1569  << "Shell region : " << shellRegionName
1570  << exit(FatalIOError);
1571  }
1572 
1573 
1574  if (oneD)
1575  {
1576  if (oneDNonManifoldEdges)
1577  {
1578  Info<< "Extruding as 1D columns with sides in patch type "
1579  << oneDPatchType
1580  << " and connected points (except on non-manifold areas)."
1581  << endl;
1582  }
1583  else
1584  {
1585  Info<< "Extruding as 1D columns with sides in patch type "
1586  << oneDPatchType
1587  << " and duplicated points (overlapping volumes)."
1588  << endl;
1589  }
1590  }
1591 
1592 
1593  // Create dummy fv* files
1594  createDummyFvMeshFiles(mesh, shellRegionName);
1595 
1596 
1597  word meshInstance;
1598  if (!overwrite)
1599  {
1600  runTime++;
1601  meshInstance = runTime.timeName();
1602  }
1603  else
1604  {
1605  meshInstance = oldInstance;
1606  }
1607  Info<< "Writing meshes to " << meshInstance << nl << endl;
1608 
1609 
1610  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1611 
1612 
1613  // Extract faces to extrude
1614  // ~~~~~~~~~~~~~~~~~~~~~~~~
1615  // Note: zoneID are regions of extrusion. They are not mesh.faceZones
1616  // indices.
1617 
1618  // From extrude zone to mesh zone (or -1 if extruding faceSets)
1619  labelList meshZoneID;
1620  // Per extrude zone whether contains internal or external faces
1621  boolList isInternal(zoneNames.size(), false);
1622 
1623  labelList extrudeMeshFaces;
1624  faceList zoneFaces;
1625  labelList zoneID;
1626  boolList zoneFlipMap;
1627  // Shadow
1628  labelList zoneShadowIDs; // from extrude shadow zone to mesh zone
1629  labelList extrudeMeshShadowFaces;
1630  boolList zoneShadowFlipMap;
1631  labelList zoneShadowID;
1632 
1633  if (hasZones)
1634  {
1635  const faceZoneMesh& faceZones = mesh.faceZones();
1636 
1637  meshZoneID.setSize(zoneNames.size());
1638  forAll(zoneNames, i)
1639  {
1640  meshZoneID[i] = faceZones.findZoneID(zoneNames[i]);
1641  if (meshZoneID[i] == -1)
1642  {
1643  FatalIOErrorIn(args.executable().c_str(), dict)
1644  << "Cannot find zone " << zoneNames[i] << endl
1645  << "Valid zones are " << faceZones.names()
1646  << exit(FatalIOError);
1647  }
1648  }
1649  // Collect per face information
1650  label nExtrudeFaces = 0;
1651  forAll(meshZoneID, i)
1652  {
1653  nExtrudeFaces += faceZones[meshZoneID[i]].size();
1654  }
1655  extrudeMeshFaces.setSize(nExtrudeFaces);
1656  zoneFaces.setSize(nExtrudeFaces);
1657  zoneID.setSize(nExtrudeFaces);
1658  zoneFlipMap.setSize(nExtrudeFaces);
1659  nExtrudeFaces = 0;
1660  forAll(meshZoneID, i)
1661  {
1662  const faceZone& fz = faceZones[meshZoneID[i]];
1663  const primitiveFacePatch& fzp = fz();
1664  forAll(fz, j)
1665  {
1666  extrudeMeshFaces[nExtrudeFaces] = fz[j];
1667  zoneFaces[nExtrudeFaces] = fzp[j];
1668  zoneID[nExtrudeFaces] = i;
1669  zoneFlipMap[nExtrudeFaces] = fz.flipMap()[j];
1670  nExtrudeFaces++;
1671 
1672  if (mesh.isInternalFace(fz[j]))
1673  {
1674  isInternal[i] = true;
1675  }
1676  }
1677  }
1678 
1679  // Shadow zone
1680  // ~~~~~~~~~~~
1681 
1682  if (zoneShadowNames.size())
1683  {
1684  zoneShadowIDs.setSize(zoneShadowNames.size());
1685  forAll(zoneShadowNames, i)
1686  {
1687  zoneShadowIDs[i] = faceZones.findZoneID(zoneShadowNames[i]);
1688  if (zoneShadowIDs[i] == -1)
1689  {
1690  FatalIOErrorIn(args.executable().c_str(), dict)
1691  << "Cannot find zone " << zoneShadowNames[i] << endl
1692  << "Valid zones are " << faceZones.names()
1693  << exit(FatalIOError);
1694  }
1695  }
1696 
1697  label nShadowFaces = 0;
1698  forAll(zoneShadowIDs, i)
1699  {
1700  nShadowFaces += faceZones[zoneShadowIDs[i]].size();
1701  }
1702 
1703  extrudeMeshShadowFaces.setSize(nShadowFaces);
1704  zoneShadowFlipMap.setSize(nShadowFaces);
1705  zoneShadowID.setSize(nShadowFaces);
1706 
1707  nShadowFaces = 0;
1708  forAll(zoneShadowIDs, i)
1709  {
1710  const faceZone& fz = faceZones[zoneShadowIDs[i]];
1711  forAll(fz, j)
1712  {
1713  extrudeMeshShadowFaces[nShadowFaces] = fz[j];
1714  zoneShadowFlipMap[nShadowFaces] = fz.flipMap()[j];
1715  zoneShadowID[nShadowFaces] = i;
1716  nShadowFaces++;
1717  }
1718  }
1719  }
1720  }
1721  else
1722  {
1723  meshZoneID.setSize(zoneNames.size(), -1);
1724  // Load faceSets
1725  PtrList<faceSet> zones(zoneNames.size());
1726  forAll(zoneNames, i)
1727  {
1728  Info<< "Loading faceSet " << zoneNames[i] << endl;
1729  zones.set(i, new faceSet(mesh, zoneNames[i]));
1730  }
1731 
1732 
1733  // Collect per face information
1734  label nExtrudeFaces = 0;
1735  forAll(zones, i)
1736  {
1737  nExtrudeFaces += zones[i].size();
1738  }
1739  extrudeMeshFaces.setSize(nExtrudeFaces);
1740  zoneFaces.setSize(nExtrudeFaces);
1741  zoneID.setSize(nExtrudeFaces);
1742  zoneFlipMap.setSize(nExtrudeFaces);
1743 
1744  nExtrudeFaces = 0;
1745  forAll(zones, i)
1746  {
1747  const faceSet& fz = zones[i];
1748  forAllConstIter(faceSet, fz, iter)
1749  {
1750  label facei = iter.key();
1751  if (mesh.isInternalFace(facei))
1752  {
1753  FatalIOErrorIn(args.executable().c_str(), dict)
1754  << "faceSet " << fz.name()
1755  << "contains internal faces."
1756  << " This is not permitted."
1757  << exit(FatalIOError);
1758  }
1759  extrudeMeshFaces[nExtrudeFaces] = facei;
1760  zoneFaces[nExtrudeFaces] = mesh.faces()[facei];
1761  zoneID[nExtrudeFaces] = i;
1762  zoneFlipMap[nExtrudeFaces] = false;
1763  nExtrudeFaces++;
1764 
1765  if (mesh.isInternalFace(facei))
1766  {
1767  isInternal[i] = true;
1768  }
1769  }
1770  }
1771 
1772 
1773  // Shadow zone
1774  // ~~~~~~~~~~~
1775 
1776  PtrList<faceSet> shadowZones(zoneShadowNames.size());
1777  if (zoneShadowNames.size())
1778  {
1779  zoneShadowIDs.setSize(zoneShadowNames.size(), -1);
1780  forAll(zoneShadowNames, i)
1781  {
1782  shadowZones.set(i, new faceSet(mesh, zoneShadowNames[i]));
1783  }
1784 
1785  label nShadowFaces = 0;
1786  forAll(shadowZones, i)
1787  {
1788  nShadowFaces += shadowZones[i].size();
1789  }
1790 
1791  if (nExtrudeFaces != nShadowFaces)
1792  {
1793  FatalIOErrorIn(args.executable().c_str(), dict)
1794  << "Extruded faces " << nExtrudeFaces << endl
1795  << "is different from shadow faces. " << nShadowFaces
1796  << "This is not permitted " << endl
1797  << exit(FatalIOError);
1798  }
1799 
1800  extrudeMeshShadowFaces.setSize(nShadowFaces);
1801  zoneShadowFlipMap.setSize(nShadowFaces);
1802  zoneShadowID.setSize(nShadowFaces);
1803 
1804  nShadowFaces = 0;
1805  forAll(shadowZones, i)
1806  {
1807  const faceSet& fz = shadowZones[i];
1808  forAllConstIter(faceSet, fz, iter)
1809  {
1810  label facei = iter.key();
1811  if (mesh.isInternalFace(facei))
1812  {
1813  FatalIOErrorIn(args.executable().c_str(), dict)
1814  << "faceSet " << fz.name()
1815  << "contains internal faces."
1816  << " This is not permitted."
1817  << exit(FatalIOError);
1818  }
1819  extrudeMeshShadowFaces[nShadowFaces] = facei;
1820  zoneShadowFlipMap[nShadowFaces] = false;
1821  zoneShadowID[nShadowFaces] = i;
1822  nShadowFaces++;
1823  }
1824  }
1825  }
1826  }
1827  const primitiveFacePatch extrudePatch(move(zoneFaces), mesh.points());
1828 
1829 
1831  Pstream::listCombineScatter(isInternal);
1832 
1833  // Check zone either all internal or all external faces
1834  checkZoneInside(mesh, zoneNames, zoneID, extrudeMeshFaces, isInternal);
1835 
1836 
1837 
1838  const pointField& extrudePoints = extrudePatch.localPoints();
1839  const faceList& extrudeFaces = extrudePatch.localFaces();
1840  const labelListList& edgeFaces = extrudePatch.edgeFaces();
1841 
1842 
1843  Info<< "extrudePatch :"
1844  << " faces:" << extrudePatch.size()
1845  << " points:" << extrudePatch.nPoints()
1846  << " edges:" << extrudePatch.nEdges()
1847  << nl
1848  << endl;
1849 
1850 
1851  // Determine per-extrude-edge info
1852  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1853 
1854  // Corresponding mesh edges
1855  const labelList extrudeMeshEdges
1856  (
1857  extrudePatch.meshEdges
1858  (
1859  mesh.edges(),
1860  mesh.pointEdges()
1861  )
1862  );
1863 
1864  const globalIndex globalExtrudeFaces(extrudePatch.size());
1865 
1866  // Global pp faces per pp edge.
1867  labelListList extrudeEdgeGlobalFaces
1868  (
1869  globalEdgeFaces
1870  (
1871  mesh,
1872  globalExtrudeFaces,
1873  extrudePatch,
1874  extrudeMeshEdges
1875  )
1876  );
1877  List<Map<label>> compactMap;
1878  const mapDistribute extrudeEdgeFacesMap
1879  (
1880  globalExtrudeFaces,
1881  extrudeEdgeGlobalFaces,
1882  compactMap
1883  );
1884 
1885 
1886  // Determine min and max zone per edge
1887  labelList edgeMinZoneID;
1888  labelList edgeMaxZoneID;
1889  calcEdgeMinMaxZone
1890  (
1891  mesh,
1892  extrudePatch,
1893  extrudeMeshEdges,
1894  zoneID,
1895  extrudeEdgeFacesMap,
1896  extrudeEdgeGlobalFaces,
1897 
1898  edgeMinZoneID,
1899  edgeMaxZoneID
1900  );
1901 
1902 
1903 
1904 
1905  DynamicList<polyPatch*> regionPatches(patches.size());
1906  // Copy all non-local patches since these are used on boundary edges of
1907  // the extrusion
1908  forAll(patches, patchi)
1909  {
1910  if (!isA<processorPolyPatch>(patches[patchi]))
1911  {
1912  label newPatchi = regionPatches.size();
1913  regionPatches.append
1914  (
1915  patches[patchi].clone
1916  (
1917  patches,
1918  newPatchi,
1919  0, // size
1920  0 // start
1921  ).ptr()
1922  );
1923  }
1924  }
1925 
1926 
1927  // Add interface patches
1928  // ~~~~~~~~~~~~~~~~~~~~~
1929 
1930  // From zone to interface patch (region side)
1931  labelList interRegionTopPatch;
1932  labelList interRegionBottomPatch;
1933 
1934  addCouplingPatches
1935  (
1936  mesh,
1937  regionName,
1938  shellRegionName,
1939  zoneNames,
1940  zoneShadowNames,
1941  isInternal,
1942  meshZoneID,
1943 
1944  regionPatches,
1945  interRegionTopPatch,
1946  interRegionBottomPatch
1947  );
1948 
1949 
1950  // From zone to interface patch (mesh side)
1951  labelList interMeshTopPatch;
1952  labelList interMeshBottomPatch;
1953 
1954  if (adaptMesh)
1955  {
1956  // Add coupling patches to mesh
1957 
1958  // Clone existing patches
1959  DynamicList<polyPatch*> newPatches(patches.size());
1960  forAll(patches, patchi)
1961  {
1962  newPatches.append(patches[patchi].clone(patches).ptr());
1963  }
1964 
1965  // Add new patches
1966  addCouplingPatches
1967  (
1968  mesh,
1969  regionName,
1970  shellRegionName,
1971  zoneNames,
1972  zoneShadowNames,
1973  isInternal,
1974  meshZoneID,
1975 
1976  newPatches,
1977  interMeshTopPatch,
1978  interMeshBottomPatch
1979  );
1980 
1981  // Add to mesh
1982  mesh.clearOut();
1983  mesh.removeFvBoundary();
1984  mesh.addFvPatches(newPatches, true);
1985 
1987  }
1988 
1989 
1990  // Patch per extruded face
1991  labelList extrudeTopPatchID(extrudePatch.size());
1992  labelList extrudeBottomPatchID(extrudePatch.size());
1993 
1994  forAll(zoneID, facei)
1995  {
1996  extrudeTopPatchID[facei] = interRegionTopPatch[zoneID[facei]];
1997  extrudeBottomPatchID[facei] = interRegionBottomPatch[zoneID[facei]];
1998  }
1999 
2000 
2001 
2002  // Count how many patches on special edges of extrudePatch are necessary
2003  // - zoneXXX_sides
2004  // - zoneXXX_zoneYYY
2005  labelList zoneSidePatch(zoneNames.size(), 0);
2006  // Patch to use for minZone
2007  labelList zoneZonePatch_min(zoneNames.size()*zoneNames.size(), 0);
2008  // Patch to use for maxZone
2009  labelList zoneZonePatch_max(zoneNames.size()*zoneNames.size(), 0);
2010 
2011  countExtrudePatches
2012  (
2013  mesh,
2014  zoneNames.size(),
2015 
2016  extrudePatch, // patch
2017  extrudeMeshFaces, // mesh face per patch face
2018  extrudeMeshEdges, // mesh edge per patch edge
2019 
2020  extrudeEdgeGlobalFaces, // global indexing per patch edge
2021  edgeMinZoneID, // minZone per patch edge
2022  edgeMaxZoneID, // maxZone per patch edge
2023 
2024  zoneSidePatch, // per zone-side num edges that extrude into it
2025  zoneZonePatch_min // per zone-zone num edges that extrude into it
2026  );
2027 
2028  // Now we'll have:
2029  // zoneSidePatch[zoneA] : number of faces needed on the side of zoneA
2030  // zoneZonePatch_min[zoneA,zoneB] : number of faces needed in between A,B
2031 
2032 
2033  // Add the zone-side patches.
2034  addZoneSidePatches
2035  (
2036  mesh,
2037  zoneNames,
2038  (oneD ? oneDPatchType : word::null),
2039 
2040  regionPatches,
2041  zoneSidePatch
2042  );
2043 
2044 
2045  // Add the patches in between zones
2046  addInterZonePatches
2047  (
2048  mesh,
2049  zoneNames,
2050  oneD,
2051 
2052  zoneZonePatch_min,
2053  zoneZonePatch_max,
2054  regionPatches
2055  );
2056 
2057 
2058  // Sets sidePatchID[edgeI] to interprocessor patch. Adds any
2059  // interprocessor or cyclic patches if necessary.
2060  labelList sidePatchID;
2061  addCoupledPatches
2062  (
2063  mesh,
2064  extrudePatch,
2065  extrudeMeshFaces,
2066  extrudeMeshEdges,
2067  extrudeEdgeFacesMap,
2068  extrudeEdgeGlobalFaces,
2069 
2070  sidePatchID,
2071  regionPatches
2072  );
2073 
2074 
2075 // // Add all the newPatches to the mesh and fields
2076 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2077 // {
2078 // forAll(newPatches, patchi)
2079 // {
2080 // Pout<< "Adding patch " << patchi
2081 // << " name:" << newPatches[patchi]->name()
2082 // << endl;
2083 // }
2084 // // label nOldPatches = mesh.boundary().size();
2085 // mesh.clearOut();
2086 // mesh.removeFvBoundary();
2087 // mesh.addFvPatches(newPatches, true);
2088 // //// Add calculated fvPatchFields for the added patches
2089 // // for
2090 // //(
2091 // // label patchi = nOldPatches;
2092 // // patchi < mesh.boundary().size();
2093 // // patchi++
2094 // //)
2095 // //{
2096 // // Pout<< "ADDing calculated to patch " << patchi
2097 // // << endl;
2098 // // addCalculatedPatchFields(mesh);
2099 // //}
2100 // // Pout<< "** Added " << mesh.boundary().size()-nOldPatches
2101 // // << " patches." << endl;
2102 // }
2103 
2104 
2105  // Set patches to use for edges to be extruded into boundary faces
2106  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2107  // In order of edgeFaces: per edge, per originating face the
2108  // patch to use for the side face (from the extruded edge).
2109  // If empty size create an internal face.
2110  labelListList extrudeEdgePatches(extrudePatch.nEdges());
2111 
2112  // Is edge a non-manifold edge
2113  PackedBoolList nonManifoldEdge(extrudePatch.nEdges());
2114 
2115  // Note: logic has to be same as in countExtrudePatches.
2116  forAll(edgeFaces, edgeI)
2117  {
2118  const labelList& eFaces = edgeFaces[edgeI];
2119 
2120  labelList& ePatches = extrudeEdgePatches[edgeI];
2121 
2122  if (oneD)
2123  {
2124  ePatches.setSize(eFaces.size());
2125  forAll(eFaces, i)
2126  {
2127  ePatches[i] = zoneSidePatch[zoneID[eFaces[i]]];
2128  }
2129 
2130  if (oneDNonManifoldEdges)
2131  {
2132  //- Set nonManifoldEdge[edgeI] for non-manifold edges only
2133  // The other option is to have non-manifold edges everywhere
2134  // and generate space overlapping columns of cells.
2135  if (eFaces.size() != 2)
2136  {
2137  nonManifoldEdge[edgeI] = 1;
2138  }
2139  }
2140  else
2141  {
2142  nonManifoldEdge[edgeI] = 1;
2143  }
2144  }
2145  else if (eFaces.size() == 2)
2146  {
2147  label zone0 = zoneID[eFaces[0]];
2148  label zone1 = zoneID[eFaces[1]];
2149 
2150  if (zone0 != zone1) // || (cos(angle) > blabla))
2151  {
2152  label minZone = min(zone0,zone1);
2153  label maxZone = max(zone0,zone1);
2154  label index = minZone*zoneNames.size()+maxZone;
2155 
2156  ePatches.setSize(eFaces.size());
2157 
2158  if (zone0 == minZone)
2159  {
2160  ePatches[0] = zoneZonePatch_min[index];
2161  ePatches[1] = zoneZonePatch_max[index];
2162  }
2163  else
2164  {
2165  ePatches[0] = zoneZonePatch_max[index];
2166  ePatches[1] = zoneZonePatch_min[index];
2167  }
2168 
2169  nonManifoldEdge[edgeI] = 1;
2170  }
2171  }
2172  else if (sidePatchID[edgeI] != -1)
2173  {
2174  // Coupled extrusion
2175  ePatches.setSize(eFaces.size());
2176  forAll(eFaces, i)
2177  {
2178  ePatches[i] = sidePatchID[edgeI];
2179  }
2180  }
2181  else
2182  {
2183  label facei = findUncoveredPatchFace
2184  (
2185  mesh,
2186  UIndirectList<label>(extrudeMeshFaces, eFaces),
2187  extrudeMeshEdges[edgeI]
2188  );
2189 
2190  if (facei != -1)
2191  {
2192  label newPatchi = findPatchID
2193  (
2194  regionPatches,
2195  patches[patches.whichPatch(facei)].name()
2196  );
2197  ePatches.setSize(eFaces.size(), newPatchi);
2198  }
2199  else
2200  {
2201  ePatches.setSize(eFaces.size());
2202  forAll(eFaces, i)
2203  {
2204  ePatches[i] = zoneSidePatch[zoneID[eFaces[i]]];
2205  }
2206  }
2207  nonManifoldEdge[edgeI] = 1;
2208  }
2209  }
2210 
2211 
2212 
2213  // Assign point regions
2214  // ~~~~~~~~~~~~~~~~~~~~
2215 
2216  // Per face, per point the region number.
2217  faceList pointGlobalRegions;
2218  faceList pointLocalRegions;
2219  labelList localToGlobalRegion;
2220 
2222  (
2223  mesh.globalData(),
2224  extrudePatch,
2225  nonManifoldEdge,
2226  false, // keep cyclic separated regions apart
2227 
2228  pointGlobalRegions,
2229  pointLocalRegions,
2230  localToGlobalRegion
2231  );
2232 
2233  // Per local region an originating point
2234  labelList localRegionPoints(localToGlobalRegion.size());
2235  forAll(pointLocalRegions, facei)
2236  {
2237  const face& f = extrudePatch.localFaces()[facei];
2238  const face& pRegions = pointLocalRegions[facei];
2239  forAll(pRegions, fp)
2240  {
2241  localRegionPoints[pRegions[fp]] = f[fp];
2242  }
2243  }
2244 
2245  // Calculate region normals by reducing local region normals
2246  pointField localRegionNormals(localToGlobalRegion.size());
2247  {
2248  pointField localSum(localToGlobalRegion.size(), Zero);
2249 
2250  forAll(pointLocalRegions, facei)
2251  {
2252  const face& pRegions = pointLocalRegions[facei];
2253  forAll(pRegions, fp)
2254  {
2255  label localRegionI = pRegions[fp];
2256  localSum[localRegionI] += extrudePatch.faceNormals()[facei];
2257  }
2258  }
2259 
2260  Map<point> globalSum(2*localToGlobalRegion.size());
2261 
2262  forAll(localSum, localRegionI)
2263  {
2264  label globalRegionI = localToGlobalRegion[localRegionI];
2265  globalSum.insert(globalRegionI, localSum[localRegionI]);
2266  }
2267 
2268  // Reduce
2270  Pstream::mapCombineScatter(globalSum);
2271 
2272  forAll(localToGlobalRegion, localRegionI)
2273  {
2274  label globalRegionI = localToGlobalRegion[localRegionI];
2275  localRegionNormals[localRegionI] = globalSum[globalRegionI];
2276  }
2277  localRegionNormals /= mag(localRegionNormals);
2278  }
2279 
2280 
2281  // For debugging: dump hedgehog plot of normals
2282  if (false)
2283  {
2284  OFstream str(runTime.path()/"localRegionNormals.obj");
2285  label vertI = 0;
2286 
2287  scalar thickness = model().sumThickness(1);
2288 
2289  forAll(pointLocalRegions, facei)
2290  {
2291  const face& f = extrudeFaces[facei];
2292 
2293  forAll(f, fp)
2294  {
2295  label region = pointLocalRegions[facei][fp];
2296  const point& pt = extrudePoints[f[fp]];
2297 
2298  meshTools::writeOBJ(str, pt);
2299  vertI++;
2301  (
2302  str,
2303  pt+thickness*localRegionNormals[region]
2304  );
2305  vertI++;
2306  str << "l " << vertI-1 << ' ' << vertI << nl;
2307  }
2308  }
2309  }
2310 
2311 
2312  // Use model to create displacements of first layer
2313  vectorField firstDisp(localRegionNormals.size());
2314  forAll(firstDisp, regionI)
2315  {
2316  // const point& regionPt = regionCentres[regionI];
2317  const point& regionPt = extrudePatch.points()
2318  [
2319  extrudePatch.meshPoints()
2320  [
2321  localRegionPoints[regionI]
2322  ]
2323  ];
2324  const vector& n = localRegionNormals[regionI];
2325  firstDisp[regionI] = model()(regionPt, n, 1) - regionPt;
2326  }
2327 
2328 
2329  // Create a new mesh
2330  // ~~~~~~~~~~~~~~~~~
2331 
2332  createShellMesh extruder
2333  (
2334  extrudePatch,
2335  pointLocalRegions,
2336  localRegionPoints
2337  );
2338 
2339 
2340  autoPtr<mapPolyMesh> shellMap;
2341  fvMesh regionMesh
2342  (
2343  IOobject
2344  (
2345  shellRegionName,
2346  meshInstance,
2347  runTime,
2348  IOobject::NO_READ,
2350  false
2351  ),
2352  pointField(),
2353  faceList(),
2354  labelList(),
2355  labelList(),
2356  false
2357  );
2358 
2359  // Add the new patches
2360  forAll(regionPatches, patchi)
2361  {
2362  polyPatch* ppPtr = regionPatches[patchi];
2363  regionPatches[patchi] = ppPtr->clone(regionMesh.boundaryMesh()).ptr();
2364  delete ppPtr;
2365  }
2366  regionMesh.clearOut();
2367  regionMesh.removeFvBoundary();
2368  regionMesh.addFvPatches(regionPatches, true);
2369 
2370  {
2371  polyTopoChange meshMod(regionPatches.size());
2372 
2373  extruder.setRefinement
2374  (
2375  firstDisp, // first displacement
2376  model().expansionRatio(),
2377  model().nLayers(), // nLayers
2378  extrudeTopPatchID,
2379  extrudeBottomPatchID,
2380  extrudeEdgePatches,
2381  meshMod
2382  );
2383 
2384  // Enforce actual point posititions according to extrudeModel (model)
2385  // (extruder.setRefinement only does fixed expansionRatio)
2386  // The regionPoints and nLayers are looped in the same way as in
2387  // createShellMesh
2388  DynamicList<point>& newPoints = const_cast<DynamicList<point>&>
2389  (
2390  meshMod.points()
2391  );
2392  label meshPointi = extrudePatch.localPoints().size();
2393  forAll(localRegionPoints, regionI)
2394  {
2395  label pointi = localRegionPoints[regionI];
2396  point pt = extrudePatch.localPoints()[pointi];
2397  const vector& n = localRegionNormals[regionI];
2398 
2399  for (label layerI = 1; layerI <= model().nLayers(); layerI++)
2400  {
2401  newPoints[meshPointi++] = model()(pt, n, layerI);
2402  }
2403  }
2404 
2405  shellMap = meshMod.changeMesh
2406  (
2407  regionMesh, // mesh to change
2408  false // inflate
2409  );
2410  }
2411 
2412  // Necessary?
2413  regionMesh.setInstance(meshInstance);
2414 
2415 
2416  // Update numbering on extruder.
2417  extruder.updateMesh(shellMap);
2418 
2419 
2420  // Calculate offsets from shell mesh back to original mesh
2421  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2422 
2423  List<pointField> topOffsets(zoneNames.size());
2424  List<pointField> bottomOffsets(zoneNames.size());
2425 
2426  forAll(regionMesh.boundaryMesh(), patchi)
2427  {
2428  const polyPatch& pp = regionMesh.boundaryMesh()[patchi];
2429 
2430  if (isA<mappedWallPolyPatch>(pp))
2431  {
2432  if (findIndex(interRegionTopPatch, patchi) != -1)
2433  {
2434  label zoneI = findIndex(interRegionTopPatch, patchi);
2435  topOffsets[zoneI] = calcOffset(extrudePatch, extruder, pp);
2436  }
2437  else if (findIndex(interRegionBottomPatch, patchi) != -1)
2438  {
2439  label zoneI = findIndex(interRegionBottomPatch, patchi);
2440  bottomOffsets[zoneI] = calcOffset(extrudePatch, extruder, pp);
2441  }
2442  }
2443  }
2444 
2445 
2446  // Change top and bottom boundary conditions on regionMesh
2447  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2448 
2449  {
2450  // Correct top patches for offset
2451  setCouplingInfo
2452  (
2453  regionMesh,
2454  interRegionTopPatch,
2455  regionName, // name of main mesh
2456  sampleMode, // sampleMode
2457  topOffsets
2458  );
2459 
2460  // Correct bottom patches for offset
2461  setCouplingInfo
2462  (
2463  regionMesh,
2464  interRegionBottomPatch,
2465  regionName,
2466  sampleMode, // sampleMode
2467  bottomOffsets
2468  );
2469 
2470  // Remove any unused patches
2471  deleteEmptyPatches(regionMesh);
2472  }
2473 
2474  // Change top and bottom boundary conditions on main mesh
2475  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2476 
2477  if (adaptMesh)
2478  {
2479  // Correct top patches for offset
2480  setCouplingInfo
2481  (
2482  mesh,
2483  interMeshTopPatch,
2484  shellRegionName, // name of shell mesh
2485  sampleMode, // sampleMode
2486  -topOffsets
2487  );
2488 
2489  // Correct bottom patches for offset
2490  setCouplingInfo
2491  (
2492  mesh,
2493  interMeshBottomPatch,
2494  shellRegionName,
2495  sampleMode,
2496  -bottomOffsets
2497  );
2498  }
2499 
2500 
2501 
2502  // Write addressing from region mesh back to originating patch
2503  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2504 
2505  labelIOList cellToPatchFaceAddressing
2506  (
2507  IOobject
2508  (
2509  "cellToPatchFaceAddressing",
2510  regionMesh.facesInstance(),
2511  regionMesh.meshSubDir,
2512  regionMesh,
2515  false
2516  ),
2517  extruder.cellToFaceMap()
2518  );
2519  cellToPatchFaceAddressing.note() = "cell to patch face addressing";
2520 
2521  labelIOList faceToPatchFaceAddressing
2522  (
2523  IOobject
2524  (
2525  "faceToPatchFaceAddressing",
2526  regionMesh.facesInstance(),
2527  regionMesh.meshSubDir,
2528  regionMesh,
2531  false
2532  ),
2533  extruder.faceToFaceMap()
2534  );
2535  faceToPatchFaceAddressing.note() =
2536  "front/back face + turning index to patch face addressing";
2537 
2538  labelIOList faceToPatchEdgeAddressing
2539  (
2540  IOobject
2541  (
2542  "faceToPatchEdgeAddressing",
2543  regionMesh.facesInstance(),
2544  regionMesh.meshSubDir,
2545  regionMesh,
2548  false
2549  ),
2550  extruder.faceToEdgeMap()
2551  );
2552  faceToPatchEdgeAddressing.note() =
2553  "side face to patch edge addressing";
2554 
2555  labelIOList pointToPatchPointAddressing
2556  (
2557  IOobject
2558  (
2559  "pointToPatchPointAddressing",
2560  regionMesh.facesInstance(),
2561  regionMesh.meshSubDir,
2562  regionMesh,
2565  false
2566  ),
2567  extruder.pointToPointMap()
2568  );
2569  pointToPatchPointAddressing.note() =
2570  "point to patch point addressing";
2571 
2572 
2573  Info<< "Writing mesh " << regionMesh.name()
2574  << " to " << regionMesh.facesInstance() << nl
2575  << endl;
2576 
2577  bool ok =
2578  regionMesh.write()
2579  && cellToPatchFaceAddressing.write()
2580  && faceToPatchFaceAddressing.write()
2581  && faceToPatchEdgeAddressing.write()
2582  && pointToPatchPointAddressing.write();
2583 
2584  if (!ok)
2585  {
2587  << "Failed writing mesh " << regionMesh.name()
2588  << " at location " << regionMesh.facesInstance()
2589  << exit(FatalError);
2590  }
2591 
2592 
2593  // See if we need to extrude coordinates as well
2594  {
2595  autoPtr<pointIOField> patchFaceCentresPtr;
2596 
2597  IOobject io
2598  (
2599  "patchFaceCentres",
2600  mesh.pointsInstance(),
2601  mesh.meshSubDir,
2602  mesh,
2604  );
2605  if (io.typeHeaderOk<pointIOField>(true))
2606  {
2607  // Read patchFaceCentres and patchEdgeCentres
2608  Info<< "Reading patch face,edge centres : "
2609  << io.name() << " and patchEdgeCentres" << endl;
2610 
2611  extrudeGeometricProperties
2612  (
2613  mesh,
2614  extrudePatch,
2615  extruder,
2616  regionMesh,
2617  model()
2618  );
2619  }
2620  }
2621 
2622 
2623 
2624 
2625  // Insert baffles into original mesh
2626  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2627 
2628  autoPtr<mapPolyMesh> addBafflesMap;
2629 
2630  if (adaptMesh)
2631  {
2632  polyTopoChange meshMod(mesh);
2633 
2634  // Modify faces to be in bottom (= always coupled) patch
2635  forAll(extrudeMeshFaces, zoneFacei)
2636  {
2637  label meshFacei = extrudeMeshFaces[zoneFacei];
2638  label zoneI = zoneID[zoneFacei];
2639  bool flip = zoneFlipMap[zoneFacei];
2640  const face& f = mesh.faces()[meshFacei];
2641 
2642  if (!flip)
2643  {
2644  meshMod.modifyFace
2645  (
2646  f, // modified face
2647  meshFacei, // label of face being modified
2648  mesh.faceOwner()[meshFacei],// owner
2649  -1, // neighbour
2650  false, // face flip
2651  interMeshBottomPatch[zoneI],// patch for face
2652  meshZoneID[zoneI], // zone for face
2653  flip // face flip in zone
2654  );
2655  }
2656  else if (mesh.isInternalFace(meshFacei))
2657  {
2658  meshMod.modifyFace
2659  (
2660  f.reverseFace(), // modified face
2661  meshFacei, // label of modified face
2662  mesh.faceNeighbour()[meshFacei],// owner
2663  -1, // neighbour
2664  true, // face flip
2665  interMeshBottomPatch[zoneI], // patch for face
2666  meshZoneID[zoneI], // zone for face
2667  !flip // face flip in zone
2668  );
2669  }
2670  }
2671 
2672  if (zoneShadowNames.size() > 0) // if there is a top faceZone specified
2673  {
2674  forAll(extrudeMeshFaces, zoneFacei)
2675  {
2676  label meshFacei = extrudeMeshShadowFaces[zoneFacei];
2677  label zoneI = zoneShadowID[zoneFacei];
2678  bool flip = zoneShadowFlipMap[zoneFacei];
2679  const face& f = mesh.faces()[meshFacei];
2680 
2681  if (!flip)
2682  {
2683  meshMod.modifyFace
2684  (
2685  f, // modified face
2686  meshFacei, // face being modified
2687  mesh.faceOwner()[meshFacei],// owner
2688  -1, // neighbour
2689  false, // face flip
2690  interMeshTopPatch[zoneI], // patch for face
2691  meshZoneID[zoneI], // zone for face
2692  flip // face flip in zone
2693  );
2694  }
2695  else if (mesh.isInternalFace(meshFacei))
2696  {
2697  meshMod.modifyFace
2698  (
2699  f.reverseFace(), // modified face
2700  meshFacei, // label modified face
2701  mesh.faceNeighbour()[meshFacei],// owner
2702  -1, // neighbour
2703  true, // face flip
2704  interMeshTopPatch[zoneI], // patch for face
2705  meshZoneID[zoneI], // zone for face
2706  !flip // face flip in zone
2707  );
2708  }
2709  }
2710  }
2711  else
2712  {
2713  // Add faces (using same points) to be in top patch
2714  forAll(extrudeMeshFaces, zoneFacei)
2715  {
2716  label meshFacei = extrudeMeshFaces[zoneFacei];
2717  label zoneI = zoneID[zoneFacei];
2718  bool flip = zoneFlipMap[zoneFacei];
2719  const face& f = mesh.faces()[meshFacei];
2720 
2721  if (!flip)
2722  {
2723  if (mesh.isInternalFace(meshFacei))
2724  {
2725  meshMod.addFace
2726  (
2727  f.reverseFace(), // modified face
2728  mesh.faceNeighbour()[meshFacei],// owner
2729  -1, // neighbour
2730  -1, // master point
2731  -1, // master edge
2732  meshFacei, // master face
2733  true, // flip flux
2734  interMeshTopPatch[zoneI], // patch for face
2735  -1, // zone for face
2736  false // face flip in zone
2737  );
2738  }
2739  }
2740  else
2741  {
2742  meshMod.addFace
2743  (
2744  f, // face
2745  mesh.faceOwner()[meshFacei], // owner
2746  -1, // neighbour
2747  -1, // master point
2748  -1, // master edge
2749  meshFacei, // master face
2750  false, // flip flux
2751  interMeshTopPatch[zoneI], // patch for face
2752  -1, // zone for face
2753  false // zone flip
2754  );
2755  }
2756  }
2757  }
2758 
2759  // Change the mesh. Change points directly (no inflation).
2760  addBafflesMap = meshMod.changeMesh(mesh, false);
2761 
2762  // Update fields
2763  mesh.updateMesh(addBafflesMap);
2764 
2765 
2766 //XXXXXX
2767 // Update maps! e.g. faceToPatchFaceAddressing
2768 //XXXXXX
2769 
2770  // Move mesh (since morphing might not do this)
2771  if (addBafflesMap().hasMotionPoints())
2772  {
2773  mesh.movePoints(addBafflesMap().preMotionPoints());
2774  }
2775 
2776  mesh.setInstance(meshInstance);
2777 
2778  // Remove any unused patches
2779  deleteEmptyPatches(mesh);
2780 
2781  Info<< "Writing mesh " << mesh.name()
2782  << " to " << mesh.facesInstance() << nl
2783  << endl;
2784 
2785  if (!mesh.write())
2786  {
2788  << "Failed writing mesh " << mesh.name()
2789  << " at location " << mesh.facesInstance()
2790  << exit(FatalError);
2791  }
2792  }
2793 
2794  Info << "End\n" << endl;
2795 
2796  return 0;
2797 }
2798 
2799 
2800 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
label nPoints() const
Return number of points supporting patch faces.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
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.
dictionary dict
static void reorderPatches(fvMesh &, const labelList &oldToNew, const label nPatches, const bool validBoundary)
Reorder and remove trailing patches. If validBoundary call is parallel.
Definition: fvMeshTools.C:140
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:268
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
const word & executable() const
Name of executable without the path.
Definition: argListI.H:30
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1287
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
mode_t mode(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file mode.
Definition: POSIX.C:461
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:313
fileName path() const
Return path.
Definition: Time.H:266
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:473
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
const word & name() const
Return name.
A line primitive.
Definition: line.H:56
const word & name() const
Return name.
Definition: IOobject.H:295
A list of face labels.
Definition: faceSet.H:48
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:476
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:808
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void addFvPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:455
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1175
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
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
const labelListList & pointEdges() const
label nFaces() const
Foam::word regionName
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
Output to file stream.
Definition: OFstream.H:82
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
static void matchEdges(const PrimitivePatch< FaceList1, PointField1 > &p1, const PrimitivePatch< FaceList2, PointField2 > &p2, labelList &p1EdgeLabels, labelList &p2EdgeLabels, PackedBoolList &sameOrientation)
Find corresponding edges on patches sharing the same points.
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:248
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
List< face > faceList
Definition: faceListFwd.H:43
label nCells() const
engineTime & runTime
static word newName(const label myProcNo, const label neighbProcNo)
Return the name of a processorPolyPatch.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
label findPatchID(const word &patchName) const
Find patch index given a name.
const labelList & faceToEdgeMap() const
From region side-face to patch edge. -1 for non-edge faces.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
patches[0]
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Pre-declare related SubField type.
Definition: Field.H:60
static tmp< pointField > edgeNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &, const labelList &patchEdges, const labelList &coupledEdges)
Return parallel consistent edge normals for patches using mesh points.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:821
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
Top level extrusion model class.
Definition: extrudeModel.H:51
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:222
string & note()
Return non-constant access to the optional note.
Definition: IOobject.H:313
scalar y
A list of faces which address into the list of points.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
face reverseFace() const
Return face with reverse direction.
Definition: face.C:619
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
T clone(const T &t)
Definition: List.H:54
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:802
dynamicFvMesh & mesh
static const NamedEnum< transformType, 5 > transformTypeNames
const labelList & cellToFaceMap() const
From region cell to patch face. Consecutively added so.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:341
static autoPtr< extrudeModel > New(const dictionary &)
Select null constructed.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh edges.
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
A class for handling words, derived from string.
Definition: word.H:59
label size() const
Return the number of elements in the list.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
wordList names() const
Return a list of patch names.
label nLayers() const
Definition: extrudeModel.C:59
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
const labelList & faceToFaceMap() const
From region face to patch face. Contains turning index:
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
Creates mesh by extruding a patch.
const Field< PointType > & points() const
Return reference to global points.
static const word null
An empty word.
Definition: word.H:77
const labelListList & edgeFaces() const
Return edge-face addressing.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
static const label labelMax
Definition: label.H:62
List< label > labelList
A List of labels.
Definition: labelList.H:56
const word dictName("particleTrackDict")
static const zero Zero
Definition: zero.H:97
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
const Field< PointType > & faceNormals() const
Return face normals for patch.
const vectorField & cellCentres() const
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const word & system() const
Return system name.
Definition: TimePaths.H:114
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
Foam::polyBoundaryMesh.
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
bool checkParallelSync(const bool report=false) const
Check whether all procs have all patches and in same order. Return.
label nEdges() const
Return number of edges in patch.
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
static const char nl
Definition: Ostream.H:265
void updateMesh(const mapPolyMesh &)
Update any locally stored mesh information.
const Time & time() const
Return time.
scalar sumThickness(const label layer) const
Helper: calculate cumulative relative thickness for layer.
Definition: extrudeModel.C:71
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
static const NamedEnum< sampleMode, 6 > sampleModeNames_
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
const vectorField & faceCentres() const
void setSize(const label)
Reset size of List.
Definition: List.C:281
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
A bit-packed bool list.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
label patchi
Class containing processor-to-processor mapping information.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
const labelList & pointToPointMap() const
From region point to patch point.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing.
Definition: fvMatrix.H:106
Direct mesh changes based on v1.3 polyTopoChange syntax.
sampleMode
Mesh items to sample.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:303
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:948
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:256
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
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
void setRefinement(const pointField &firstLayerThickness, const scalar expansionRatio, const label nLayers, const labelList &topPatchID, const labelList &bottomPatchID, const labelListList &extrudeEdgePatches, polyTopoChange &meshMod)
Play commands into polyTopoChange to create layer mesh.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual Ostream & write(const token &)=0
Write next token to stream.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:158
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: PtrList.H:53
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
T & last()
Return the last element of the list.
Definition: UListI.H:128
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
const labelListList & edgeFaces() const
static const label labelMin
Definition: label.H:61
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49
IOerror FatalIOError
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:284