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