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