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