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