extrudeMesh.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  extrudeMesh
26 
27 Description
28  Extrude mesh from existing patch (by default outwards facing normals;
29  optional flips faces) or from patch read from file.
30 
31  Note: Merges close points so be careful.
32 
33  Type of extrusion prescribed by run-time selectable model.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "argList.H"
38 #include "Time.H"
39 #include "polyTopoChange.H"
40 #include "polyTopoChanger.H"
41 #include "edgeCollapser.H"
42 #include "perfectInterface.H"
43 #include "addPatchCellLayer.H"
44 #include "fvMesh.H"
45 #include "MeshedSurfaces.H"
46 #include "globalIndex.H"
47 #include "cellSet.H"
48 #include "systemDict.H"
49 
50 #include "extrudedMesh.H"
51 #include "extrudeModel.H"
52 
53 #include "wedge.H"
54 #include "wedgePolyPatch.H"
55 #include "planeExtrusion.H"
56 #include "emptyPolyPatch.H"
57 
58 using namespace Foam;
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 enum ExtrudeMode
63 {
64  MESH,
65  PATCH,
66  SURFACE
67 };
68 
69 namespace Foam
70 {
71  template<>
72  const char* NamedEnum<ExtrudeMode, 3>::names[] =
73  {
74  "mesh",
75  "patch",
76  "surface"
77  };
78 }
79 
80 static const NamedEnum<ExtrudeMode, 3> ExtrudeModeNames;
81 
82 
83 void createDummyFvMeshFiles(const polyMesh& mesh, const word& regionName)
84 {
85  // Create dummy system/fv*
86  {
87  IOobject io
88  (
89  "fvSchemes",
90  mesh.time().system(),
91  regionName,
92  mesh,
95  false
96  );
97 
98  Info<< "Testing:" << io.objectPath() << endl;
99 
100  if (!io.typeHeaderOk<IOdictionary>(false))
101  {
102  Info<< "Writing dummy " << regionName/io.name() << endl;
103  dictionary dummyDict;
104  dictionary divDict;
105  dummyDict.add("divSchemes", divDict);
106  dictionary gradDict;
107  dummyDict.add("gradSchemes", gradDict);
108  dictionary laplDict;
109  dummyDict.add("laplacianSchemes", laplDict);
110 
111  IOdictionary(io, dummyDict).regIOobject::write();
112  }
113  }
114  {
115  IOobject io
116  (
117  "fvSolution",
118  mesh.time().system(),
119  regionName,
120  mesh,
123  false
124  );
125 
126  if (!io.typeHeaderOk<IOdictionary>(false))
127  {
128  Info<< "Writing dummy " << regionName/io.name() << endl;
129  dictionary dummyDict;
130  IOdictionary(io, dummyDict).regIOobject::write();
131  }
132  }
133 }
134 
135 
136 label findPatchID(const polyBoundaryMesh& patches, const word& name)
137 {
138  const label patchID = patches.findPatchID(name);
139 
140  if (patchID == -1)
141  {
143  << "Cannot find patch " << name
144  << " in the source mesh.\n"
145  << "Valid patch names are " << patches.names()
146  << exit(FatalError);
147  }
148  return patchID;
149 }
150 
151 
152 labelList patchFaces(const polyBoundaryMesh& patches, const wordList& names)
153 {
154  label n = 0;
155 
156  forAll(names, i)
157  {
158  const polyPatch& pp = patches[findPatchID(patches, names[i])];
159 
160  n += pp.size();
161  }
162  labelList faceLabels(n);
163  n = 0;
164  forAll(names, i)
165  {
166  const polyPatch& pp = patches[findPatchID(patches, names[i])];
167 
168  forAll(pp, j)
169  {
170  faceLabels[n++] = pp.start()+j;
171  }
172  }
173 
174  return faceLabels;
175 }
176 
177 
178 void updateFaceLabels(const mapPolyMesh& map, labelList& faceLabels)
179 {
180  const labelList& reverseMap = map.reverseFaceMap();
181 
182  labelList newFaceLabels(faceLabels.size());
183  label newI = 0;
184 
185  forAll(faceLabels, i)
186  {
187  label oldFacei = faceLabels[i];
188 
189  if (reverseMap[oldFacei] >= 0)
190  {
191  newFaceLabels[newI++] = reverseMap[oldFacei];
192  }
193  }
194  newFaceLabels.setSize(newI);
195  faceLabels.transfer(newFaceLabels);
196 }
197 
198 
199 void updateCellSet(const mapPolyMesh& map, labelHashSet& cellLabels)
200 {
201  const labelList& reverseMap = map.reverseCellMap();
202 
203  labelHashSet newCellLabels(2*cellLabels.size());
204 
205  forAll(cellLabels, i)
206  {
207  label oldCelli = cellLabels[i];
208 
209  if (reverseMap[oldCelli] >= 0)
210  {
211  newCellLabels.insert(reverseMap[oldCelli]);
212  }
213  }
214  cellLabels.transfer(newCellLabels);
215 }
216 
217 
218 template<class PatchType>
219 void changeFrontBackPatches
220 (
221  polyMesh& mesh,
222  const word& frontPatchName,
223  const word& backPatchName
224 )
225 {
226  const polyBoundaryMesh& patches = mesh.boundaryMesh();
227 
228  label frontPatchi = findPatchID(patches, frontPatchName);
229  label backPatchi = findPatchID(patches, backPatchName);
230 
231  DynamicList<polyPatch*> newPatches(patches.size());
232 
233  forAll(patches, patchi)
234  {
235  const polyPatch& pp(patches[patchi]);
236 
237  if (patchi == frontPatchi || patchi == backPatchi)
238  {
239  newPatches.append
240  (
241  new PatchType
242  (
243  pp.name(),
244  pp.size(),
245  pp.start(),
246  pp.index(),
247  mesh.boundaryMesh(),
248  PatchType::typeName
249  )
250  );
251  }
252  else
253  {
254  newPatches.append(pp.clone(mesh.boundaryMesh()).ptr());
255  }
256  }
257 
258  // Edit patches
259  mesh.removeBoundary();
260  mesh.addPatches(newPatches, true);
261 }
262 
263 
264 int main(int argc, char *argv[])
265 {
266  #include "addRegionOption.H"
267  #include "addDictOption.H"
268 
269  #include "setRootCase.H"
270  #include "createTimeExtruded.H"
271 
272  // Get optional regionName
274  word regionDir;
275  if (args.optionReadIfPresent("region", regionName))
276  {
277  regionDir = regionName;
278  Info<< "Create mesh " << regionName << " for time = "
279  << runTimeExtruded.timeName() << nl << endl;
280  }
281  else
282  {
283  regionName = fvMesh::defaultRegion;
284  Info<< "Create mesh for time = "
285  << runTimeExtruded.timeName() << nl << endl;
286  }
287 
288  const dictionary dict(systemDict("extrudeMeshDict", args, runTimeExtruded));
289 
290  // Point generator
292 
293  // Whether to flip normals
294  const Switch flipNormals(dict.lookup("flipNormals"));
295 
296  // What to extrude
297  const ExtrudeMode mode = ExtrudeModeNames.read
298  (
299  dict.lookup("constructFrom")
300  );
301 
302  // Any merging of small edges
303  const scalar mergeTol(dict.lookupOrDefault<scalar>("mergeTol", 1e-4));
304 
305  Info<< "Extruding from " << ExtrudeModeNames[mode]
306  << " using model " << model().type() << endl;
307  if (flipNormals)
308  {
309  Info<< "Flipping normals before extruding" << endl;
310  }
311  if (mergeTol > 0)
312  {
313  Info<< "Collapsing edges < " << mergeTol << " of bounding box" << endl;
314  }
315  else
316  {
317  Info<< "Not collapsing any edges after extrusion" << endl;
318  }
319  Info<< endl;
320 
321 
322  // Generated mesh (one of either)
323  autoPtr<fvMesh> meshFromMesh;
324  autoPtr<polyMesh> meshFromSurface;
325 
326  // Faces on front and back for stitching (in case of mergeFaces)
327  word frontPatchName;
328  labelList frontPatchFaces;
329  word backPatchName;
330  labelList backPatchFaces;
331 
332  // Optional added cells (get written to cellSet)
333  labelHashSet addedCellsSet;
334 
335  if (mode == PATCH || mode == MESH)
336  {
337  if (flipNormals && mode == MESH)
338  {
340  << "Flipping normals not supported for extrusions from mesh."
341  << exit(FatalError);
342  }
343 
344  fileName sourceCasePath(dict.lookup("sourceCase"));
345  sourceCasePath.expand();
346  fileName sourceRootDir = sourceCasePath.path();
347  fileName sourceCaseDir = sourceCasePath.name();
348  if (Pstream::parRun())
349  {
350  sourceCaseDir =
351  sourceCaseDir
352  /"processor" + Foam::name(Pstream::myProcNo());
353  }
354  wordList sourcePatches;
355  dict.lookup("sourcePatches") >> sourcePatches;
356 
357  if (sourcePatches.size() == 1)
358  {
359  frontPatchName = sourcePatches[0];
360  }
361 
362  Info<< "Extruding patches " << sourcePatches
363  << " on mesh " << sourceCasePath << nl
364  << endl;
365 
366  Time runTime
367  (
369  sourceRootDir,
370  sourceCaseDir
371  );
372 
373  #include "createMesh.H"
374 
375  const polyBoundaryMesh& patches = mesh.boundaryMesh();
376 
377 
378  // Extrusion engine. Either adding to existing mesh or
379  // creating separate mesh.
380  addPatchCellLayer layerExtrude(mesh, (mode == MESH));
381 
382  const labelList meshFaces(patchFaces(patches, sourcePatches));
383 
384  if (mode == PATCH && flipNormals)
385  {
386  // Cheat. Flip patch faces in mesh. This invalidates the
387  // mesh (open cells) but does produce the correct extrusion.
388  polyTopoChange meshMod(mesh);
389  forAll(meshFaces, i)
390  {
391  label meshFacei = meshFaces[i];
392 
393  label patchi = patches.whichPatch(meshFacei);
394  label own = mesh.faceOwner()[meshFacei];
395  label nei = -1;
396  if (patchi == -1)
397  {
398  nei = mesh.faceNeighbour()[meshFacei];
399  }
400 
401  label zoneI = mesh.faceZones().whichZone(meshFacei);
402  bool zoneFlip = false;
403  if (zoneI != -1)
404  {
405  label index = mesh.faceZones()[zoneI].whichFace(meshFacei);
406  zoneFlip = mesh.faceZones()[zoneI].flipMap()[index];
407  }
408 
409  meshMod.modifyFace
410  (
411  mesh.faces()[meshFacei].reverseFace(), // modified face
412  meshFacei, // label of face
413  own, // owner
414  nei, // neighbour
415  true, // face flip
416  patchi, // patch for face
417  zoneI, // zone for face
418  zoneFlip // face flip in zone
419  );
420  }
421 
422  // Change the mesh. No inflation.
423  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
424 
425  // Update fields
426  mesh.updateMesh(map);
427 
428  // Move mesh (since morphing does not do this)
429  if (map().hasMotionPoints())
430  {
431  mesh.movePoints(map().preMotionPoints());
432  }
433  }
434 
435 
436 
437  indirectPrimitivePatch extrudePatch
438  (
440  (
441  mesh.faces(),
442  meshFaces
443  ),
444  mesh.points()
445  );
446 
447  // Determine extrudePatch normal
448  pointField extrudePatchPointNormals
449  (
450  PatchTools::pointNormals(mesh, extrudePatch)
451  );
452 
453 
454  // Precalculate mesh edges for pp.edges.
455  const labelList meshEdges
456  (
457  extrudePatch.meshEdges
458  (
459  mesh.edges(),
460  mesh.pointEdges()
461  )
462  );
463 
464  // Global face indices engine
465  const globalIndex globalFaces(mesh.nFaces());
466 
467  // Determine extrudePatch.edgeFaces in global numbering (so across
468  // coupled patches)
469  labelListList edgeGlobalFaces
470  (
472  (
473  mesh,
474  globalFaces,
475  extrudePatch
476  )
477  );
478 
479 
480  // Determine what patches boundary edges need to get extruded into.
481  // This might actually cause edge-connected processors to become
482  // face-connected so might need to introduce new processor boundaries.
483  // Calculates:
484  // - per pp.edge the patch to extrude into
485  // - any additional processor boundaries (patchToNbrProc = map from
486  // new patchID to neighbour processor)
487  // - number of new patches (nPatches)
488 
489  labelList sidePatchID;
490  label nPatches;
491  Map<label> nbrProcToPatch;
492  Map<label> patchToNbrProc;
494  (
495  mesh,
496  globalFaces,
497  edgeGlobalFaces,
498  extrudePatch,
499 
500  sidePatchID,
501  nPatches,
502  nbrProcToPatch,
503  patchToNbrProc
504  );
505 
506 
507  // Add any patches.
508 
509  label nAdded = nPatches - mesh.boundaryMesh().size();
510  reduce(nAdded, sumOp<label>());
511 
512  Info<< "Adding overall " << nAdded << " processor patches." << endl;
513 
514  if (nAdded > 0)
515  {
516  DynamicList<polyPatch*> newPatches(nPatches);
517  forAll(mesh.boundaryMesh(), patchi)
518  {
519  newPatches.append
520  (
521  mesh.boundaryMesh()[patchi].clone
522  (
523  mesh.boundaryMesh()
524  ).ptr()
525  );
526  }
527  for
528  (
529  label patchi = mesh.boundaryMesh().size();
530  patchi < nPatches;
531  patchi++
532  )
533  {
534  label nbrProci = patchToNbrProc[patchi];
535 
536  Pout<< "Adding patch " << patchi
537  << " between " << Pstream::myProcNo()
538  << " and " << nbrProci
539  << endl;
540 
541  newPatches.append
542  (
544  (
545  0, // size
546  mesh.nFaces(), // start
547  patchi, // index
548  mesh.boundaryMesh(),// polyBoundaryMesh
549  Pstream::myProcNo(),// myProcNo
550  nbrProci // neighbProcNo
551  )
552  );
553  }
554 
555  // Add patches. Do no parallel updates.
556  mesh.removeFvBoundary();
557  mesh.addFvPatches(newPatches, true);
558  }
559 
560 
561 
562  // Only used for addPatchCellLayer into new mesh
563  labelList exposedPatchID;
564  if (mode == PATCH)
565  {
566  dict.lookup("exposedPatchName") >> backPatchName;
567  exposedPatchID.setSize
568  (
569  extrudePatch.size(),
570  findPatchID(patches, backPatchName)
571  );
572  }
573 
574  // Determine points and extrusion
575  pointField layer0Points(extrudePatch.nPoints());
576  pointField displacement(extrudePatch.nPoints());
577  forAll(displacement, pointi)
578  {
579  const vector& patchNormal = extrudePatchPointNormals[pointi];
580 
581  // layer0 point
582  layer0Points[pointi] = model()
583  (
584  extrudePatch.localPoints()[pointi],
585  patchNormal,
586  0
587  );
588  // layerN point
589  point extrudePt = model()
590  (
591  extrudePatch.localPoints()[pointi],
592  patchNormal,
593  model().nLayers()
594  );
595  displacement[pointi] = extrudePt - layer0Points[pointi];
596  }
597 
598 
599  // Check if wedge (has layer0 different from original patch points)
600  // If so move the mesh to starting position.
601  if (gMax(mag(layer0Points-extrudePatch.localPoints())) > small)
602  {
603  Info<< "Moving mesh to layer0 points since differ from original"
604  << " points - this can happen for wedge extrusions." << nl
605  << endl;
606 
607  pointField newPoints(mesh.points());
608  forAll(extrudePatch.meshPoints(), i)
609  {
610  newPoints[extrudePatch.meshPoints()[i]] = layer0Points[i];
611  }
612  mesh.movePoints(newPoints);
613  }
614 
615 
616  // Layers per face
617  labelList nFaceLayers(extrudePatch.size(), model().nLayers());
618 
619  // Layers per point
620  labelList nPointLayers(extrudePatch.nPoints(), model().nLayers());
621 
622  // Displacement for first layer
623  vectorField firstLayerDisp(displacement*model().sumThickness(1));
624 
625  // Expansion ratio not used.
626  scalarField ratio(extrudePatch.nPoints(), 1.0);
627 
628  // Topo change container. Either copy an existing mesh or start
629  // with empty storage (number of patches only needed for checking)
631  (
632  (
633  mode == MESH
634  ? new polyTopoChange(mesh)
635  : new polyTopoChange(patches.size())
636  )
637  );
638 
639  layerExtrude.setRefinement
640  (
641  globalFaces,
642  edgeGlobalFaces,
643 
644  ratio, // expansion ratio
645  extrudePatch, // patch faces to extrude
646  sidePatchID, // if boundary edge: patch to use
647  exposedPatchID, // if new mesh: patches for exposed faces
648  nFaceLayers,
649  nPointLayers,
650  firstLayerDisp,
651  meshMod()
652  );
653 
654  // Reset points according to extrusion model
655  forAll(layerExtrude.addedPoints(), pointi)
656  {
657  const labelList& pPoints = layerExtrude.addedPoints()[pointi];
658  forAll(pPoints, pPointi)
659  {
660  label meshPointi = pPoints[pPointi];
661 
662  point modelPt
663  (
664  model()
665  (
666  extrudePatch.localPoints()[pointi],
667  extrudePatchPointNormals[pointi],
668  pPointi+1 // layer
669  )
670  );
671 
672  const_cast<DynamicList<point>&>
673  (
674  meshMod().points()
675  )[meshPointi] = modelPt;
676  }
677  }
678 
679  // Store faces on front and exposed patch (if mode=patch there are
680  // only added faces so cannot used map to old faces)
681  const labelListList& layerFaces = layerExtrude.layerFaces();
682  backPatchFaces.setSize(layerFaces.size());
683  frontPatchFaces.setSize(layerFaces.size());
684  forAll(backPatchFaces, patchFacei)
685  {
686  backPatchFaces[patchFacei] = layerFaces[patchFacei].first();
687  frontPatchFaces[patchFacei] = layerFaces[patchFacei].last();
688  }
689 
690 
691  // Create dummy fvSchemes, fvSolution
692  createDummyFvMeshFiles(mesh, regionDir);
693 
694  // Create actual mesh from polyTopoChange container
695  autoPtr<mapPolyMesh> map = meshMod().makeMesh
696  (
697  meshFromMesh,
698  IOobject
699  (
700  regionName,
701  runTimeExtruded.constant(),
702  runTimeExtruded,
705  false
706  ),
707  mesh
708  );
709 
710  layerExtrude.updateMesh
711  (
712  map(),
713  identity(extrudePatch.size()),
714  identity(extrudePatch.nPoints())
715  );
716 
717  // Calculate face labels for front and back.
718  frontPatchFaces = renumber
719  (
720  map().reverseFaceMap(),
721  frontPatchFaces
722  );
723  backPatchFaces = renumber
724  (
725  map().reverseFaceMap(),
726  backPatchFaces
727  );
728 
729  // Store added cells
730  if (mode == MESH)
731  {
732  const labelListList addedCells
733  (
734  layerExtrude.addedCells
735  (
736  meshFromMesh,
737  layerExtrude.layerFaces()
738  )
739  );
740  forAll(addedCells, facei)
741  {
742  const labelList& aCells = addedCells[facei];
743  forAll(aCells, i)
744  {
745  addedCellsSet.insert(aCells[i]);
746  }
747  }
748  }
749  }
750  else
751  {
752  // Read from surface
753  fileName surfName(dict.lookup("surface"));
754  surfName.expand();
755 
756  Info<< "Extruding surfaceMesh read from file " << surfName << nl
757  << endl;
758 
759  MeshedSurface<face> fMesh(surfName);
760 
761  if (flipNormals)
762  {
763  Info<< "Flipping faces." << nl << endl;
764  faceList& faces = const_cast<faceList&>(fMesh.faces());
765  forAll(faces, i)
766  {
767  faces[i] = fMesh[i].reverseFace();
768  }
769  }
770 
771  Info<< "Extruding surface with :" << nl
772  << " points : " << fMesh.points().size() << nl
773  << " faces : " << fMesh.size() << nl
774  << " normals[0] : " << fMesh.faceNormals()[0]
775  << nl
776  << endl;
777 
778  meshFromSurface.reset
779  (
780  new extrudedMesh
781  (
782  IOobject
783  (
785  runTimeExtruded.constant(),
786  runTimeExtruded
787  ),
788  fMesh,
789  model()
790  )
791  );
792 
793 
794  // Get the faces on front and back
795  frontPatchName = "originalPatch";
796  frontPatchFaces = patchFaces
797  (
798  meshFromSurface().boundaryMesh(),
799  wordList(1, frontPatchName)
800  );
801  backPatchName = "otherSide";
802  backPatchFaces = patchFaces
803  (
804  meshFromSurface().boundaryMesh(),
805  wordList(1, backPatchName)
806  );
807  }
808 
809 
810  polyMesh& mesh =
811  (
812  meshFromMesh.valid()
813  ? meshFromMesh()
814  : meshFromSurface()
815  );
816 
817 
818  const boundBox& bb = mesh.bounds();
819  const vector span = bb.span();
820  const scalar mergeDim = mergeTol * bb.minDim();
821 
822  Info<< "Mesh bounding box : " << bb << nl
823  << " with span : " << span << nl
824  << "Merge distance : " << mergeDim << nl
825  << endl;
826 
827 
828  // Collapse edges
829  // ~~~~~~~~~~~~~~
830 
831  if (mergeDim > 0)
832  {
833  Pout<< "Collapsing edges < " << mergeDim << " ..." << nl << endl;
834 
835  // Edge collapsing engine
836  edgeCollapser collapser(mesh);
837 
838  const edgeList& edges = mesh.edges();
839  const pointField& points = mesh.points();
840 
842  Map<point> collapsePointToLocation(mesh.nPoints());
843 
844  forAll(edges, edgeI)
845  {
846  const edge& e = edges[edgeI];
847 
848  scalar d = e.mag(points);
849 
850  if (d < mergeDim)
851  {
852  Pout<< "Merging edge " << e << " since length " << d
853  << " << " << mergeDim << nl;
854 
855  collapseEdge[edgeI] = true;
856  collapsePointToLocation.set(e[1], points[e[0]]);
857  }
858  }
859 
860  List<pointEdgeCollapse> allPointInfo;
861  const globalIndex globalPoints(mesh.nPoints());
862  labelList pointPriority(mesh.nPoints(), 0);
863 
864  collapser.consistentCollapse
865  (
866  globalPoints,
867  pointPriority,
868  collapsePointToLocation,
869  collapseEdge,
870  allPointInfo
871  );
872 
873  // Topo change container
874  polyTopoChange meshMod(mesh);
875 
876  // Put all modifications into meshMod
877  bool anyChange = collapser.setRefinement(allPointInfo, meshMod);
878  reduce(anyChange, orOp<bool>());
879 
880  if (anyChange)
881  {
882  // Construct new mesh from polyTopoChange.
883  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
884 
885  // Update fields
886  mesh.updateMesh(map);
887 
888  // Update stored data
889  updateFaceLabels(map(), frontPatchFaces);
890  updateFaceLabels(map(), backPatchFaces);
891  updateCellSet(map(), addedCellsSet);
892 
893  // Move mesh (if inflation used)
894  if (map().hasMotionPoints())
895  {
896  mesh.movePoints(map().preMotionPoints());
897  }
898  }
899  }
900 
901 
902  // Change the front and back patch types as required
903  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
904 
905  word frontBackType(word::null);
906 
907  if (isType<extrudeModels::wedge>(model()))
908  {
909  changeFrontBackPatches<wedgePolyPatch>
910  (
911  mesh,
912  frontPatchName,
913  backPatchName
914  );
915  }
916  else if (isType<extrudeModels::plane>(model()))
917  {
918  changeFrontBackPatches<emptyPolyPatch>
919  (
920  mesh,
921  frontPatchName,
922  backPatchName
923  );
924  }
925 
926 
927  // Merging front and back patch faces
928  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
929 
930  Switch mergeFaces(dict.lookup("mergeFaces"));
931  if (mergeFaces)
932  {
933  if (mode == MESH)
934  {
936  << "Cannot stitch front and back of extrusion since"
937  << " in 'mesh' mode (extrusion appended to mesh)."
938  << exit(FatalError);
939  }
940 
941  Info<< "Assuming full 360 degree axisymmetric case;"
942  << " stitching faces on patches "
943  << frontPatchName << " and "
944  << backPatchName << " together ..." << nl << endl;
945 
946  if (frontPatchFaces.size() != backPatchFaces.size())
947  {
949  << "Differing number of faces on front ("
950  << frontPatchFaces.size() << ") and back ("
951  << backPatchFaces.size() << ")"
952  << exit(FatalError);
953  }
954 
955 
956 
957  polyTopoChanger stitcher(mesh);
958  stitcher.setSize(1);
959 
960  const word cutZoneName("originalCutFaceZone");
961 
962  List<faceZone*> fz
963  (
964  1,
965  new faceZone
966  (
967  cutZoneName,
968  frontPatchFaces,
969  boolList(frontPatchFaces.size(), false),
970  0,
971  mesh.faceZones()
972  )
973  );
974 
975  mesh.addZones(List<pointZone*>(0), fz, List<cellZone*>(0));
976 
977  // Add the perfect interface mesh modifier
978  perfectInterface perfectStitcher
979  (
980  "couple",
981  0,
982  stitcher,
983  cutZoneName,
984  word::null, // dummy patch name
985  word::null // dummy patch name
986  );
987 
988  // Topo change container
989  polyTopoChange meshMod(mesh);
990 
991  perfectStitcher.setRefinement
992  (
994  (
996  (
997  mesh.faces(),
998  frontPatchFaces
999  ),
1000  mesh.points()
1001  ),
1003  (
1005  (
1006  mesh.faces(),
1007  backPatchFaces
1008  ),
1009  mesh.points()
1010  ),
1011  meshMod
1012  );
1013 
1014  // Construct new mesh from polyTopoChange.
1015  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
1016 
1017  // Update fields
1018  mesh.updateMesh(map);
1019 
1020  // Update local data
1021  updateCellSet(map(), addedCellsSet);
1022 
1023  // Move mesh (if inflation used)
1024  if (map().hasMotionPoints())
1025  {
1026  mesh.movePoints(map().preMotionPoints());
1027  }
1028  }
1029 
1030  mesh.setInstance(runTimeExtruded.constant());
1031  Info<< "Writing mesh to " << mesh.localObjectPath() << nl << endl;
1032 
1033  if (!mesh.write())
1034  {
1036  << exit(FatalError);
1037  }
1038 
1039  // Need writing cellSet
1040  label nAdded = returnReduce(addedCellsSet.size(), sumOp<label>());
1041  if (nAdded > 0)
1042  {
1043  cellSet addedCells(mesh, "addedCells", addedCellsSet);
1044  Info<< "Writing added cells to cellSet " << addedCells.name()
1045  << nl << endl;
1046  if (!addedCells.write())
1047  {
1049  << addedCells.name()
1050  << exit(FatalError);
1051  }
1052  }
1053 
1054  Info<< "End\n" << endl;
1055 
1056  return 0;
1057 }
1058 
1059 
1060 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:402
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
dictionary dict
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:276
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1287
#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
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 class for handling file names.
Definition: fileName.H:79
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#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
const labelListList & pointEdges() const
label nFaces() const
Foam::word regionName
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
engineTime & runTime
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Does polyTopoChanges to remove edges. Can remove faces due to edge collapse but can not remove cells ...
Definition: edgeCollapser.H:66
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
label findPatchID(const word &patchName) const
Find patch index given a name.
label collapseEdge(triSurface &surf, const scalar minLen)
Keep collapsing all edges < minLen.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
void transfer(HashTable< T, Key, Hash > &)
Transfer the contents of the argument table into this table.
Definition: HashTable.C:513
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
T & first()
Return the first element of the list.
Definition: UListI.H:114
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1133
Neighbour processor patch.
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
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:228
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
const word & regionDir(const word &regionName)
dynamicFvMesh & mesh
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
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
List of mesh modifiers defining the mesh dynamics.
Adds layers of cells to outside of polyPatch. Can optionally create stand-alone extruded mesh (addToM...
A class for handling words, derived from string.
Definition: word.H:59
wordList names() const
Return a list of patch names.
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
static const word null
An empty word.
Definition: word.H:77
fileName localObjectPath() const
Return complete localPath + object name.
Definition: IOobject.H:425
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: edgeI.H:181
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
void addZones(const List< pointZone *> &pz, const List< faceZone *> &fz, const List< cellZone *> &cz)
Add mesh zones.
Definition: polyMesh.C:953
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:199
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
void addPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:909
const word & system() const
Return system name.
Definition: TimePaths.H:113
Foam::polyBoundaryMesh.
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
static const char nl
Definition: Ostream.H:260
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:521
Type gMax(const FieldField< Field, Type > &f)
const Time & time() const
Return time.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
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 nEdges() const
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &)
Return parallel consistent point normals for patches using mesh points.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
List< word > wordList
A List of words.
Definition: fileName.H:54
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:100
void setSize(const label)
Reset size of List.
Definition: List.C:281
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: MeshZones.C:221
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:440
A bit-packed bool list.
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
string & expand(const bool allowEmpty=false)
Expand initial tildes and all occurrences of environment variables.
Definition: string.C:103
A collection of cell labels.
Definition: cellSet.H:48
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:476
label index() const
Return the index of this patch in the boundaryMesh.
Direct mesh changes based on v1.3 polyTopoChange syntax.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:309
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
static labelListList globalEdgeFaces(const polyMesh &, const globalIndex &globalFaces, const indirectPrimitivePatch &pp)
Per patch edge the pp faces (in global indices) using it. Uses.
label nPoints() const
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:61
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:265
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
Hack of attachDetach to couple patches when they perfectly align. Does not decouple. Used by stitchMesh app. Does geometric matching.
virtual bool write(const bool write=true) const
Write using setting from DB.
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
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static void calcSidePatch(const polyMesh &, const globalIndex &globalFaces, const labelListList &globalEdgeFaces, const indirectPrimitivePatch &pp, labelList &sidePatchID, label &nPatches, Map< label > &nbrProcToPatch, Map< label > &patchToNbrProc)
Boundary edges get extruded into boundary faces. Determine patch.
A List with indirect addressing.
Definition: IndirectList.H:101
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
scalar minDim() const
Smallest length/height/width dimension.
Definition: boundBoxI.H:102
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:490
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