extrudeMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  extrudeMesh
29 
30 Group
31  grpMeshGenerationUtilities
32 
33 Description
34  Extrude mesh from existing patch (by default outwards facing normals;
35  optional flips faces) or from patch read from file.
36 
37  Note: Merges close points so be careful.
38 
39  Type of extrusion prescribed by run-time selectable model.
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "argList.H"
44 #include "Time.H"
45 #include "polyTopoChange.H"
46 #include "polyTopoChanger.H"
47 #include "edgeCollapser.H"
48 #include "perfectInterface.H"
49 #include "addPatchCellLayer.H"
50 #include "fvMesh.H"
51 #include "MeshedSurfaces.H"
52 #include "globalIndex.H"
53 #include "cellSet.H"
54 #include "fvMeshTools.H"
55 
56 #include "extrudedMesh.H"
57 #include "extrudeModel.H"
58 
59 #include "wedge.H"
60 #include "wedgePolyPatch.H"
61 #include "planeExtrusion.H"
62 #include "emptyPolyPatch.H"
63 #include "processorMeshes.H"
64 #include "hexRef8Data.H"
65 
66 using namespace Foam;
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 enum ExtrudeMode
71 {
72  MESH,
73  PATCH,
74  SURFACE
75 };
76 
77 static const Enum<ExtrudeMode> ExtrudeModeNames
78 {
79  { ExtrudeMode::MESH, "mesh" },
80  { ExtrudeMode::PATCH, "patch" },
81  { ExtrudeMode::SURFACE, "surface" },
82 };
83 
84 
85 label findPatchID(const polyBoundaryMesh& patches, const word& name)
86 {
87  const label patchID = patches.findPatchID(name);
88 
89  if (patchID == -1)
90  {
92  << "Cannot find patch " << name
93  << " in the source mesh.\n"
94  << "Valid patch names are " << patches.names()
95  << exit(FatalError);
96  }
97  return patchID;
98 }
99 
100 
101 labelList patchFaces(const polyBoundaryMesh& patches, const wordList& names)
102 {
103  label n = 0;
104 
105  forAll(names, i)
106  {
107  const polyPatch& pp = patches[findPatchID(patches, names[i])];
108 
109  n += pp.size();
110  }
111  labelList faceLabels(n);
112  n = 0;
113  forAll(names, i)
114  {
115  const polyPatch& pp = patches[findPatchID(patches, names[i])];
116 
117  forAll(pp, j)
118  {
119  faceLabels[n++] = pp.start()+j;
120  }
121  }
122 
123  return faceLabels;
124 }
125 
126 
127 void updateFaceLabels(const mapPolyMesh& map, labelList& faceLabels)
128 {
129  const labelList& reverseMap = map.reverseFaceMap();
130 
131  labelList newFaceLabels(faceLabels.size());
132  label newI = 0;
133 
134  forAll(faceLabels, i)
135  {
136  label oldFacei = faceLabels[i];
137 
138  if (reverseMap[oldFacei] >= 0)
139  {
140  newFaceLabels[newI++] = reverseMap[oldFacei];
141  }
142  }
143  newFaceLabels.setSize(newI);
144  faceLabels.transfer(newFaceLabels);
145 }
146 
147 
148 void updateCellSet(const mapPolyMesh& map, labelHashSet& cellLabels)
149 {
150  const labelList& reverseMap = map.reverseCellMap();
151 
152  labelHashSet newCellLabels(2*cellLabels.size());
153 
154  forAll(cellLabels, i)
155  {
156  label oldCelli = cellLabels[i];
157 
158  if (reverseMap[oldCelli] >= 0)
159  {
160  newCellLabels.insert(reverseMap[oldCelli]);
161  }
162  }
163  cellLabels.transfer(newCellLabels);
164 }
165 
166 
167 template<class PatchType>
168 void changeFrontBackPatches
169 (
170  polyMesh& mesh,
171  const word& frontPatchName,
172  const word& backPatchName
173 )
174 {
176 
177  label frontPatchi = findPatchID(patches, frontPatchName);
178  label backPatchi = findPatchID(patches, backPatchName);
179 
180  DynamicList<polyPatch*> newPatches(patches.size());
181 
182  forAll(patches, patchi)
183  {
184  const polyPatch& pp(patches[patchi]);
185 
186  if (patchi == frontPatchi || patchi == backPatchi)
187  {
188  newPatches.append
189  (
190  new PatchType
191  (
192  pp.name(),
193  pp.size(),
194  pp.start(),
195  pp.index(),
196  mesh.boundaryMesh(),
197  PatchType::typeName
198  )
199  );
200  }
201  else
202  {
203  newPatches.append(pp.clone(mesh.boundaryMesh()).ptr());
204  }
205  }
206 
207  // Edit patches
209  mesh.addPatches(newPatches, true);
210 }
211 
212 
213 int main(int argc, char *argv[])
214 {
216  (
217  "Extrude mesh from existing patch."
218  );
219 
220  #include "addRegionOption.H"
221  argList::addOption("dict", "file", "Alternative extrudeMeshDict");
222  #include "setRootCase.H"
223  #include "createTimeExtruded.H"
224 
225  // Get optional regionName
227  word regionDir;
228  if (args.readIfPresent("region", regionName))
229  {
231  Info<< "Create mesh " << regionName << " for time = "
232  << runTimeExtruded.timeName() << nl << endl;
233  }
234  else
235  {
237  Info<< "Create mesh for time = "
238  << runTimeExtruded.timeName() << nl << endl;
239  }
240 
241  const IOdictionary dict
242  (
244  (
245  IOobject
246  (
247  "extrudeMeshDict",
248  runTimeExtruded.system(),
249  runTimeExtruded,
252  ),
253  args.getOrDefault<fileName>("dict", "")
254  )
255  );
256 
257  // Point generator
259 
260  // Whether to flip normals
261  const bool flipNormals(dict.get<bool>("flipNormals"));
262 
263  // What to extrude
264  const ExtrudeMode mode = ExtrudeModeNames.get("constructFrom", dict);
265 
266  // Any merging of small edges
267  const scalar mergeTol(dict.getOrDefault<scalar>("mergeTol", 1e-4));
268 
269  Info<< "Extruding from " << ExtrudeModeNames[mode]
270  << " using model " << model().type() << endl;
271  if (flipNormals)
272  {
273  Info<< "Flipping normals before extruding" << endl;
274  }
275  if (mergeTol > 0)
276  {
277  Info<< "Collapsing edges < " << mergeTol << " of bounding box" << endl;
278  }
279  else
280  {
281  Info<< "Not collapsing any edges after extrusion" << endl;
282  }
283  Info<< endl;
284 
285 
286  // Generated mesh (one of either)
287  autoPtr<fvMesh> meshFromMesh;
288  autoPtr<polyMesh> meshFromSurface;
289 
290  // Faces on front and back for stitching (in case of mergeFaces)
291  word frontPatchName;
292  labelList frontPatchFaces;
293  word backPatchName;
294  labelList backPatchFaces;
295 
296  // Optional added cells (get written to cellSet)
297  labelHashSet addedCellsSet;
298 
299  // Optional refinement data
300  autoPtr<hexRef8Data> refDataPtr;
301 
302  if (mode == PATCH || mode == MESH)
303  {
304  if (flipNormals && mode == MESH)
305  {
307  << "Flipping normals not supported for extrusions from mesh."
308  << exit(FatalError);
309  }
310 
311  fileName sourceCasePath(dict.get<fileName>("sourceCase").expand());
312  fileName sourceRootDir = sourceCasePath.path();
313  fileName sourceCaseDir = sourceCasePath.name();
314  if (Pstream::parRun())
315  {
316  sourceCaseDir =
317  sourceCaseDir/("processor" + Foam::name(Pstream::myProcNo()));
318  }
319  wordList sourcePatches;
320  dict.readEntry("sourcePatches", sourcePatches);
321 
322  if (sourcePatches.size() == 1)
323  {
324  frontPatchName = sourcePatches[0];
325  }
326 
327  Info<< "Extruding patches " << sourcePatches
328  << " on mesh " << sourceCasePath << nl
329  << endl;
330 
331  Time runTime
332  (
334  sourceRootDir,
335  sourceCaseDir
336  );
337 
338  #include "createNamedMesh.H"
339 
341 
342 
343  // Extrusion engine. Either adding to existing mesh or
344  // creating separate mesh.
345  addPatchCellLayer layerExtrude(mesh, (mode == MESH));
346 
347  const labelList meshFaces(patchFaces(patches, sourcePatches));
348 
349  if (mode == PATCH && flipNormals)
350  {
351  // Cheat. Flip patch faces in mesh. This invalidates the
352  // mesh (open cells) but does produce the correct extrusion.
353  polyTopoChange meshMod(mesh);
354  forAll(meshFaces, i)
355  {
356  label meshFacei = meshFaces[i];
357 
358  label patchi = patches.whichPatch(meshFacei);
359  label own = mesh.faceOwner()[meshFacei];
360  label nei = -1;
361  if (patchi == -1)
362  {
363  nei = mesh.faceNeighbour()[meshFacei];
364  }
365 
366  label zoneI = mesh.faceZones().whichZone(meshFacei);
367  bool zoneFlip = false;
368  if (zoneI != -1)
369  {
370  label index = mesh.faceZones()[zoneI].whichFace(meshFacei);
371  zoneFlip = mesh.faceZones()[zoneI].flipMap()[index];
372  }
373 
374  meshMod.modifyFace
375  (
376  mesh.faces()[meshFacei].reverseFace(), // modified face
377  meshFacei, // label of face
378  own, // owner
379  nei, // neighbour
380  true, // face flip
381  patchi, // patch for face
382  zoneI, // zone for face
383  zoneFlip // face flip in zone
384  );
385  }
386 
387  // Change the mesh. No inflation.
388  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
389 
390  // Update fields
391  mesh.updateMesh(map());
392 
393  // Move mesh (since morphing does not do this)
394  if (map().hasMotionPoints())
395  {
396  mesh.movePoints(map().preMotionPoints());
397  }
398  }
399 
400 
401  indirectPrimitivePatch extrudePatch
402  (
404  (
405  mesh.faces(),
406  meshFaces
407  ),
408  mesh.points()
409  );
410 
411  // Determine extrudePatch normal
412  pointField extrudePatchPointNormals
413  (
414  PatchTools::pointNormals(mesh, extrudePatch)
415  );
416 
417 
418  // Precalculate mesh edges for pp.edges.
419  const labelList meshEdges
420  (
421  extrudePatch.meshEdges
422  (
423  mesh.edges(),
424  mesh.pointEdges()
425  )
426  );
427 
428  // Global face indices engine
429  const globalIndex globalFaces(mesh.nFaces());
430 
431  // Determine extrudePatch.edgeFaces in global numbering (so across
432  // coupled patches)
433  labelListList edgeGlobalFaces
434  (
436  (
437  mesh,
438  globalFaces,
439  extrudePatch
440  )
441  );
442 
443 
444  // Determine what patches boundary edges need to get extruded into.
445  // This might actually cause edge-connected processors to become
446  // face-connected so might need to introduce new processor boundaries.
447  // Calculates:
448  // - per pp.edge the patch to extrude into
449  // - any additional processor boundaries (patchToNbrProc = map from
450  // new patchID to neighbour processor)
451  // - number of new patches (nPatches)
452 
453  labelList edgePatchID;
454  labelList edgeZoneID;
455  boolList edgeFlip;
456  labelList inflateFaceID;
457  label nPatches;
458  Map<label> nbrProcToPatch;
459  Map<label> patchToNbrProc;
461  (
462  true, // for internal edges get zone info from any face
463 
464  mesh,
465  globalFaces,
466  edgeGlobalFaces,
467  extrudePatch,
468 
469  edgePatchID,
470  nPatches,
471  nbrProcToPatch,
472  patchToNbrProc,
473 
474  edgeZoneID,
475  edgeFlip,
476  inflateFaceID
477  );
478 
479 
480  // Add any patches.
481 
482  label nAdded = nPatches - mesh.boundaryMesh().size();
483  reduce(nAdded, sumOp<label>());
484 
485  Info<< "Adding overall " << nAdded << " processor patches." << endl;
486 
487  if (nAdded > 0)
488  {
489  DynamicList<polyPatch*> newPatches(nPatches);
490  forAll(mesh.boundaryMesh(), patchi)
491  {
492  newPatches.append
493  (
494  mesh.boundaryMesh()[patchi].clone
495  (
497  ).ptr()
498  );
499  }
500  for
501  (
502  label patchi = mesh.boundaryMesh().size();
503  patchi < nPatches;
504  patchi++
505  )
506  {
507  label nbrProci = patchToNbrProc[patchi];
508 
509  Pout<< "Adding patch " << patchi
510  << " between " << Pstream::myProcNo()
511  << " and " << nbrProci
512  << endl;
513 
514  newPatches.append
515  (
517  (
518  0, // size
519  mesh.nFaces(), // start
520  patchi, // index
521  mesh.boundaryMesh(),// polyBoundaryMesh
522  Pstream::myProcNo(),// myProcNo
523  nbrProci // neighbProcNo
524  )
525  );
526  }
527 
528  // Add patches. Do no parallel updates.
530  mesh.addFvPatches(newPatches, true);
531  }
532 
533 
534 
535  // Only used for addPatchCellLayer into new mesh
536  labelList exposedPatchID;
537  if (mode == PATCH)
538  {
539  dict.readEntry("exposedPatchName", backPatchName);
540  exposedPatchID.setSize
541  (
542  extrudePatch.size(),
543  findPatchID(patches, backPatchName)
544  );
545  }
546 
547  // Determine points and extrusion
548  pointField layer0Points(extrudePatch.nPoints());
549  pointField displacement(extrudePatch.nPoints());
550  forAll(displacement, pointi)
551  {
552  const vector& patchNormal = extrudePatchPointNormals[pointi];
553 
554  // layer0 point
555  layer0Points[pointi] = model()
556  (
557  extrudePatch.localPoints()[pointi],
558  patchNormal,
559  0
560  );
561  // layerN point
562  point extrudePt = model()
563  (
564  extrudePatch.localPoints()[pointi],
565  patchNormal,
566  model().nLayers()
567  );
568  displacement[pointi] = extrudePt - layer0Points[pointi];
569  }
570 
571 
572  // Check if wedge (has layer0 different from original patch points)
573  // If so move the mesh to starting position.
574  if (gMax(mag(layer0Points-extrudePatch.localPoints())) > SMALL)
575  {
576  Info<< "Moving mesh to layer0 points since differ from original"
577  << " points - this can happen for wedge extrusions." << nl
578  << endl;
579 
580  pointField newPoints(mesh.points());
581  forAll(extrudePatch.meshPoints(), i)
582  {
583  newPoints[extrudePatch.meshPoints()[i]] = layer0Points[i];
584  }
585  mesh.movePoints(newPoints);
586  }
587 
588 
589  // Layers per face
590  labelList nFaceLayers(extrudePatch.size(), model().nLayers());
591 
592  // Layers per point
593  labelList nPointLayers(extrudePatch.nPoints(), model().nLayers());
594 
595  // Displacement for first layer
596  vectorField firstLayerDisp(displacement*model().sumThickness(1));
597 
598  // Expansion ratio not used.
599  scalarField ratio(extrudePatch.nPoints(), 1.0);
600 
601 
602  // Load any refinement data
603  if (mode == MESH)
604  {
605  // Tricky: register hexRef8 onto the database
606  // since the mesh does not yet exist. It should not be registered
607  // onto the input mesh.
608  refDataPtr.reset
609  (
610  new hexRef8Data
611  (
612  IOobject
613  (
614  "dummy",
617  runTimeExtruded, //mesh,
620  false
621  )
622  )
623  );
624  }
625 
626 
627  // Topo change container. Either copy an existing mesh or start
628  // with empty storage (number of patches only needed for checking)
630  (
631  (
632  mode == MESH
633  ? new polyTopoChange(mesh)
634  : new polyTopoChange(patches.size())
635  )
636  );
637 
638  layerExtrude.setRefinement
639  (
640  globalFaces,
641  edgeGlobalFaces,
642 
643  ratio, // expansion ratio
644  extrudePatch, // patch faces to extrude
645 
646  edgePatchID, // if boundary edge: patch for extruded face
647  edgeZoneID, // optional zone for extruded face
648  edgeFlip,
649  inflateFaceID, // mesh face that zone/patch info is from
650 
651  exposedPatchID, // if new mesh: patches for exposed faces
652  nFaceLayers,
653  nPointLayers,
654  firstLayerDisp,
655  meshMod()
656  );
657 
658  // Reset points according to extrusion model
659  forAll(layerExtrude.addedPoints(), pointi)
660  {
661  const labelList& pPoints = layerExtrude.addedPoints()[pointi];
662  forAll(pPoints, pPointi)
663  {
664  label meshPointi = pPoints[pPointi];
665 
666  point modelPt
667  (
668  model()
669  (
670  extrudePatch.localPoints()[pointi],
671  extrudePatchPointNormals[pointi],
672  pPointi+1 // layer
673  )
674  );
675 
676  const_cast<DynamicList<point>&>
677  (
678  meshMod().points()
679  )[meshPointi] = modelPt;
680  }
681  }
682 
683  // Store faces on front and exposed patch (if mode=patch there are
684  // only added faces so cannot used map to old faces)
685  const labelListList& layerFaces = layerExtrude.layerFaces();
686  backPatchFaces.setSize(layerFaces.size());
687  frontPatchFaces.setSize(layerFaces.size());
688  forAll(backPatchFaces, patchFacei)
689  {
690  backPatchFaces[patchFacei] = layerFaces[patchFacei].first();
691  frontPatchFaces[patchFacei] = layerFaces[patchFacei].last();
692  }
693 
694 
695  // Create dummy fvSchemes, fvSolution
697 
698  // Create actual mesh from polyTopoChange container
699  autoPtr<mapPolyMesh> map = meshMod().makeMesh
700  (
701  meshFromMesh,
702  IOobject
703  (
704  regionName,
705  runTimeExtruded.constant(),
706  runTimeExtruded,
707  IOobject::READ_IF_PRESENT, // Read fv* if present
709  false
710  ),
711  mesh
712  );
713 
714  layerExtrude.updateMesh
715  (
716  map(),
717  identity(extrudePatch.size()),
718  identity(extrudePatch.nPoints())
719  );
720 
721  // Calculate face labels for front and back.
722  frontPatchFaces = renumber
723  (
724  map().reverseFaceMap(),
725  frontPatchFaces
726  );
727  backPatchFaces = renumber
728  (
729  map().reverseFaceMap(),
730  backPatchFaces
731  );
732 
733  // Update
734  if (refDataPtr)
735  {
736  refDataPtr->updateMesh(map());
737  }
738 
739  // Store added cells
740  if (mode == MESH)
741  {
742  const labelListList addedCells
743  (
744  layerExtrude.addedCells
745  (
746  *meshFromMesh,
747  layerExtrude.layerFaces()
748  )
749  );
750  forAll(addedCells, facei)
751  {
752  const labelList& aCells = addedCells[facei];
753  addedCellsSet.insert(aCells);
754  }
755  }
756  }
757  else
758  {
759  // Read from surface
760  fileName surfName(dict.get<fileName>("surface").expand());
761 
762  Info<< "Extruding surfaceMesh read from file " << surfName << nl
763  << endl;
764 
765  MeshedSurface<face> fMesh(surfName);
766 
767  if (flipNormals)
768  {
769  Info<< "Flipping faces." << nl << endl;
770  faceList& faces = const_cast<faceList&>(fMesh.surfFaces());
771  forAll(faces, i)
772  {
773  faces[i] = fMesh[i].reverseFace();
774  }
775  }
776 
777  Info<< "Extruding surface with :" << nl
778  << " points : " << fMesh.points().size() << nl
779  << " faces : " << fMesh.size() << nl
780  << " normals[0] : " << fMesh.faceNormals()[0]
781  << nl
782  << endl;
783 
784  meshFromSurface.reset
785  (
786  new extrudedMesh
787  (
788  IOobject
789  (
791  runTimeExtruded.constant(),
792  runTimeExtruded
793  ),
794  fMesh,
795  model()
796  )
797  );
798 
799 
800  // Get the faces on front and back
801  frontPatchName = "originalPatch";
802  frontPatchFaces = patchFaces
803  (
804  meshFromSurface().boundaryMesh(),
805  wordList(1, frontPatchName)
806  );
807  backPatchName = "otherSide";
808  backPatchFaces = patchFaces
809  (
810  meshFromSurface().boundaryMesh(),
811  wordList(1, backPatchName)
812  );
813  }
814 
815 
816  polyMesh& mesh =
817  (
818  meshFromMesh
819  ? *meshFromMesh
820  : *meshFromSurface
821  );
822 
823 
824  const boundBox& bb = mesh.bounds();
825  const vector span = bb.span();
826  const scalar mergeDim = mergeTol * bb.minDim();
827 
828  Info<< "Mesh bounding box : " << bb << nl
829  << " with span : " << span << nl
830  << "Merge distance : " << mergeDim << nl
831  << endl;
832 
833 
834  // Collapse edges
835  // ~~~~~~~~~~~~~~
836 
837  if (mergeDim > 0)
838  {
839  Info<< "Collapsing edges < " << mergeDim << " ..." << nl << endl;
840 
841  // Edge collapsing engine
842  edgeCollapser collapser(mesh);
843 
844  const edgeList& edges = mesh.edges();
845  const pointField& points = mesh.points();
846 
848  Map<point> collapsePointToLocation(mesh.nPoints());
849 
850  forAll(edges, edgeI)
851  {
852  const edge& e = edges[edgeI];
853 
854  scalar d = e.mag(points);
855 
856  if (d < mergeDim)
857  {
858  Info<< "Merging edge " << e << " since length " << d
859  << " << " << mergeDim << nl;
860 
861  collapseEdge.set(edgeI);
862  collapsePointToLocation.set(e[1], points[e[0]]);
863  }
864  }
865 
866  List<pointEdgeCollapse> allPointInfo;
868  labelList pointPriority(mesh.nPoints(), Zero);
869 
870  collapser.consistentCollapse
871  (
872  globalPoints,
873  pointPriority,
874  collapsePointToLocation,
875  collapseEdge,
876  allPointInfo
877  );
878 
879  // Topo change container
880  polyTopoChange meshMod(mesh);
881 
882  // Put all modifications into meshMod
883  bool anyChange = collapser.setRefinement(allPointInfo, meshMod);
884  reduce(anyChange, orOp<bool>());
885 
886  if (anyChange)
887  {
888  // Construct new mesh from polyTopoChange.
889  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
890 
891  // Update fields
892  mesh.updateMesh(map());
893 
894  // Update stored data
895  updateFaceLabels(map(), frontPatchFaces);
896  updateFaceLabels(map(), backPatchFaces);
897  updateCellSet(map(), addedCellsSet);
898 
899  if (refDataPtr)
900  {
901  refDataPtr->updateMesh(map());
902  }
903 
904  // Move mesh (if inflation used)
905  if (map().hasMotionPoints())
906  {
907  mesh.movePoints(map().preMotionPoints());
908  }
909  }
910  }
911 
912 
913  // Change the front and back patch types as required
914  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
915 
916  word frontBackType(word::null);
917 
918  if (isType<extrudeModels::wedge>(model()))
919  {
920  changeFrontBackPatches<wedgePolyPatch>
921  (
922  mesh,
923  frontPatchName,
924  backPatchName
925  );
926  }
927  else if (isType<extrudeModels::plane>(model()))
928  {
929  changeFrontBackPatches<emptyPolyPatch>
930  (
931  mesh,
932  frontPatchName,
933  backPatchName
934  );
935  }
936 
937 
938  // Merging front and back patch faces
939  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
940 
941  if (dict.get<bool>("mergeFaces"))
942  {
943  if (mode == MESH)
944  {
946  << "Cannot stitch front and back of extrusion since"
947  << " in 'mesh' mode (extrusion appended to mesh)."
948  << exit(FatalError);
949  }
950 
951  Info<< "Assuming full 360 degree axisymmetric case;"
952  << " stitching faces on patches "
953  << frontPatchName << " and "
954  << backPatchName << " together ..." << nl << endl;
955 
956  if (frontPatchFaces.size() != backPatchFaces.size())
957  {
959  << "Differing number of faces on front ("
960  << frontPatchFaces.size() << ") and back ("
961  << backPatchFaces.size() << ")"
962  << exit(FatalError);
963  }
964 
965 
966  polyTopoChanger stitcher(mesh);
967  stitcher.setSize(1);
968 
969  const word cutZoneName("originalCutFaceZone");
970 
971  List<faceZone*> fz
972  (
973  1,
974  new faceZone
975  (
976  cutZoneName,
977  frontPatchFaces,
978  false, // none are flipped
979  0,
980  mesh.faceZones()
981  )
982  );
983 
985 
986  // Add the perfect interface mesh modifier
987  perfectInterface perfectStitcher
988  (
989  "couple",
990  0,
991  stitcher,
992  cutZoneName,
993  word::null, // dummy patch name
994  word::null // dummy patch name
995  );
996 
997  // Topo change container
998  polyTopoChange meshMod(mesh);
999 
1000  perfectStitcher.setRefinement
1001  (
1003  (
1005  (
1006  mesh.faces(),
1007  frontPatchFaces
1008  ),
1009  mesh.points()
1010  ),
1012  (
1014  (
1015  mesh.faces(),
1016  backPatchFaces
1017  ),
1018  mesh.points()
1019  ),
1020  meshMod
1021  );
1022 
1023  // Construct new mesh from polyTopoChange.
1024  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
1025 
1026  // Update fields
1027  mesh.updateMesh(map());
1028 
1029  // Update local data
1030  updateCellSet(map(), addedCellsSet);
1031 
1032  if (refDataPtr)
1033  {
1034  refDataPtr->updateMesh(map());
1035  }
1036 
1037  // Move mesh (if inflation used)
1038  if (map().hasMotionPoints())
1039  {
1040  mesh.movePoints(map().preMotionPoints());
1041  }
1042  }
1043 
1044  mesh.setInstance(runTimeExtruded.constant());
1045  Info<< "Writing mesh to " << mesh.objectPath() << nl << endl;
1046 
1047  if (!mesh.write())
1048  {
1050  << exit(FatalError);
1051  }
1052  // Remove any left-over files
1054 
1055  // Need writing cellSet
1056  label nAdded = returnReduce(addedCellsSet.size(), sumOp<label>());
1057  if (nAdded > 0)
1058  {
1059  cellSet addedCells(mesh, "addedCells", addedCellsSet);
1060  Info<< "Writing added cells to cellSet " << addedCells.name()
1061  << nl << endl;
1062  if (!addedCells.write())
1063  {
1065  << addedCells.name()
1066  << exit(FatalError);
1067  }
1068  }
1069 
1070  if (refDataPtr)
1071  {
1072  refDataPtr->write();
1073  }
1074 
1075 
1076  Info<< "End\n" << endl;
1077 
1078  return 0;
1079 }
1080 
1081 
1082 // ************************************************************************* //
planeExtrusion.H
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::polyMesh::addPatches
void addPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Definition: polyMesh.C:954
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
runTime
engineTime & runTime
Definition: createEngineTime.H:13
addPatchCellLayer.H
Foam::polyMesh::points
virtual const pointField & points() const
Definition: polyMesh.C:1062
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
Foam::PatchTools::pointNormals
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &)
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:53
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:190
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Foam::globalPoints
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:98
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::fileName
A class for handling file names.
Definition: fileName.H:71
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Definition: fvMesh.C:1034
Foam::fvMeshTools::createDummyFvMeshFiles
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Definition: fvMeshTools.C:753
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:88
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:59
Foam::polyMesh::defaultRegion
static word defaultRegion
Definition: polyMesh.H:314
Foam::argList::getOrDefault
T getOrDefault(const word &optName, const T &deflt) const
Definition: argListI.H:300
Foam::fileName::path
static std::string path(const std::string &str)
Definition: fileNameI.H:169
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::addPatchCellLayer::globalEdgeFaces
static labelListList globalEdgeFaces(const polyMesh &, const globalIndex &globalFaces, const indirectPrimitivePatch &pp)
Definition: addPatchCellLayer.C:640
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:51
wedgePolyPatch.H
globalIndex.H
Foam::polyTopoChanger
List of mesh modifiers defining the mesh dynamics.
Definition: polyTopoChanger.H:58
Foam::fileName::name
static std::string name(const std::string &str)
Definition: fileNameI.H:192
polyTopoChanger.H
Foam::polyMesh::meshSubDir
static word meshSubDir
Definition: polyMesh.H:317
polyTopoChange.H
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:95
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:59
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
Foam::primitiveMesh::nEdges
label nEdges() const
Definition: primitiveMeshI.H:60
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Definition: polyMesh.C:845
Foam::Map< label >
Foam::fvMesh::removeFvBoundary
void removeFvBoundary()
Definition: fvMesh.C:628
Foam::fvMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Definition: fvMesh.C:858
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Definition: polyMesh.H:440
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:100
Foam::primitiveMesh::edges
const edgeList & edges() const
Definition: primitiveMeshEdges.C:498
Foam::Pout
prefixOSstream Pout
Foam::primitiveMesh::pointEdges
const labelListList & pointEdges() const
Definition: primitiveMeshEdges.C:509
Foam::HashSet
A HashTable with keys but without contents that is similar to std::unordered_set.
Definition: HashSet.H:73
processorMeshes.H
Foam::Time::controlDictName
static word controlDictName
Definition: Time.H:222
Foam::argList::readIfPresent
bool readIfPresent(const word &optName, T &val) const
Definition: argListI.H:316
wedge.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
Foam::polyBoundaryMesh::names
wordList names() const
Definition: polyBoundaryMesh.C:548
Foam::sumOp
Definition: ops.H:207
Foam::primitiveMesh::nPoints
label nPoints() const noexcept
Definition: primitiveMeshI.H:30
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:58
Foam::mode
mode_t mode(const fileName &name, const bool followLink=true)
Definition: POSIX.C:692
Foam::edgeCollapser
Does polyTopoChanges to remove edges. Can remove faces due to edge collapse but can not remove cells ...
Definition: edgeCollapser.H:64
Foam::boundBox::span
vector span() const
Definition: boundBoxI.H:120
n
label n
Definition: TABSMDCalcMethod2.H:31
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:45
nPatches
const label nPatches
Definition: printMeshSummary.H:24
createTimeExtruded.H
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::extrudedMesh
Definition: extrudedMesh.H:48
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
edgeCollapser.H
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:64
Foam::List::setSize
void setSize(const label n)
Definition: List.H:218
extrudeModel.H
argList.H
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Definition: polyMesh.C:1100
Foam::List::transfer
void transfer(List< T > &list)
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:183
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:52
addRegionOption.H
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Definition: fvMesh.C:946
Foam::boundBox::minDim
scalar minDim() const
Definition: boundBoxI.H:138
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:295
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Definition: polyMesh.H:482
Foam::extrudeModel::New
static autoPtr< extrudeModel > New(const dictionary &dict)
Definition: extrudeModelNew.C:27
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Definition: polyBoundaryMesh.C:805
Foam::processorPolyPatch
Neighbour processor patch.
Definition: processorPolyPatch.H:54
Foam::IOobject::selectIO
static IOobject selectIO(const IOobject &io, const fileName &altFile, const word &ioName="")
Definition: IOobject.C:231
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Definition: polyBoundaryMesh.C:758
Foam::hexRef8Data
Various for reading/decomposing/reconstructing/distributing refinement data.
Definition: hexRef8Data.H:56
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Definition: ZoneMesh.C:282
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::fvMesh::addFvPatches
void addFvPatches(PtrList< polyPatch > &plist, const bool validBoundary=true)
Definition: fvMesh.C:598
Foam::PtrList::clone
PtrList< T > clone(Args &&... args) const
Foam::perfectInterface
Hack of attachDetach to couple patches when they perfectly align. Does not decouple....
Definition: perfectInterface.H:54
createNamedMesh.H
Required Variables.
Foam::addPatchCellLayer
Adds layers of cells to outside of polyPatch. Can optionally create stand-alone extruded mesh (addToM...
Definition: addPatchCellLayer.H:124
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
hexRef8Data.H
fvMesh.H
Foam
Definition: atmBoundaryLayer.C:26
Foam::patchIdentifier::index
label index() const noexcept
Definition: patchIdentifier.H:143
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:47
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
Foam::polyPatch::start
label start() const
Definition: polyPatch.H:357
emptyPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::polyMesh::removeBoundary
void removeBoundary()
Definition: polyMeshClear.C:32
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Definition: mapPolyMesh.H:497
Foam::polyMesh::bounds
const boundBox & bounds() const
Definition: polyMesh.H:446
Foam::autoPtr::reset
void reset(autoPtr< T > &&other) noexcept
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:49
setRootCase.H
Foam::renumber
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Definition: ListOpsTemplates.C:30
Foam::mapPolyMesh::reverseCellMap
const labelList & reverseCellMap() const
Definition: mapPolyMesh.H:528
Foam::polyMesh::faces
virtual const faceList & faces() const
Definition: polyMesh.C:1087
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
extrudedMesh.H
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Definition: UPstream.H:459
Foam::nl
constexpr char nl
Definition: Ostream.H:424
fvMeshTools.H
Foam::Vector< scalar >
Foam::UPstream::parRun
static bool & parRun() noexcept
Definition: UPstream.H:429
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::addPatchCellLayer::calcExtrudeInfo
static void calcExtrudeInfo(const bool zoneFromAnyFace, const polyMesh &, const globalIndex &globalFaces, const labelListList &globalEdgeFaces, const indirectPrimitivePatch &pp, labelList &edgePatchID, label &nPatches, Map< label > &nbrProcToPatch, Map< label > &patchToNbrProc, labelList &edgeZoneID, boolList &edgeFlip, labelList &inflateFaceID)
Definition: addPatchCellLayer.C:1083
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::string::expand
string & expand(const bool allowEmpty=false)
Definition: string.C:166
Foam::processorMeshes::removeFiles
static void removeFiles(const polyMesh &mesh)
Definition: processorMeshes.C:267
MeshedSurfaces.H
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::PDRpatchDef
Bookkeeping for patch definitions.
Definition: PDRpatchDef.H:49
Foam::identity
labelList identity(const label len, label start=0)
Definition: labelList.C:31
Foam::HashSet::insert
bool insert(const Key &key)
Definition: HashSet.H:191
Foam::word::null
static const word null
Definition: word.H:78
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:182
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:59
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:57
Foam::polyPatch::clone
virtual autoPtr< polyPatch > clone(const labelList &faceCells) const
Definition: polyPatch.H:231
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
perfectInterface.H
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:58
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Definition: primitiveMeshI.H:83
Foam::patchIdentifier::name
const word & name() const noexcept
Definition: patchIdentifier.H:131
cellSet.H
Foam::PATCH
PDRpatchDef PATCH
Definition: PDRpatchDef.H:108
collapseEdge
label collapseEdge(triSurface &surf, const scalar minLen)
regionDir
const word & regionDir
Definition: findMeshDefinitionDict.H:28
Foam::IOobject::objectPath
fileName objectPath() const
Definition: IOobjectI.H:207
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:141
Foam::MeshedSurface< face >
Foam::argList::addOption
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Definition: argList.C:328
args
Foam::argList args(argc, argv)
Foam::orOp
Definition: ops.H:228
Foam::PtrListOps::names
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Foam::polyMesh::addZones
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Definition: polyMesh.C:992
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:585
Foam::polyMesh::setInstance
void setInstance(const fileName &instance, const IOobject::writeOption wOpt=IOobject::AUTO_WRITE)
Definition: polyMeshIO.C:29
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Definition: polyMesh.C:1106
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:75