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-2022 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 "processorPolyPatch.H"
64 #include "processorMeshes.H"
65 #include "hexRef8Data.H"
66 
67 using namespace Foam;
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 enum ExtrudeMode
72 {
73  MESH,
74  PATCH,
75  SURFACE
76 };
77 
78 static const Enum<ExtrudeMode> ExtrudeModeNames
79 {
80  { ExtrudeMode::MESH, "mesh" },
81  { ExtrudeMode::PATCH, "patch" },
82  { ExtrudeMode::SURFACE, "surface" },
83 };
84 
85 
86 label findPatchID(const polyBoundaryMesh& patches, const word& name)
87 {
88  const label patchID = patches.findPatchID(name);
89 
90  if (patchID == -1)
91  {
93  << "Cannot find patch " << name
94  << " in the source mesh.\n"
95  << "Valid patch names are " << patches.names()
96  << exit(FatalError);
97  }
98  return patchID;
99 }
100 
101 
102 labelList patchFaces(const polyBoundaryMesh& patches, const wordList& names)
103 {
104  label n = 0;
105 
106  forAll(names, i)
107  {
108  const polyPatch& pp = patches[findPatchID(patches, names[i])];
109 
110  n += pp.size();
111  }
113  n = 0;
114  forAll(names, i)
115  {
116  const polyPatch& pp = patches[findPatchID(patches, names[i])];
117 
118  forAll(pp, j)
119  {
120  faceLabels[n++] = pp.start()+j;
121  }
122  }
123 
124  return faceLabels;
125 }
126 
127 
128 void zoneFaces
129 (
130  const faceZoneMesh& fzs,
131  const wordList& names,
133  bitSet& faceFlip
134 )
135 {
136  label n = 0;
137 
138  forAll(names, i)
139  {
140  const auto& pp = fzs[fzs.findZoneID(names[i])];
141  n += pp.size();
142  }
144  faceFlip.setSize(n);
145  n = 0;
146  forAll(names, i)
147  {
148  const auto& pp = fzs[fzs.findZoneID(names[i])];
149  const boolList& ppFlip = pp.flipMap();
150  forAll(pp, i)
151  {
152  faceLabels[n] = pp[i];
153  faceFlip[n] = ppFlip[i];
154  n++;
155  }
156  }
157 }
158 
159 
160 void updateFaceLabels(const mapPolyMesh& map, labelList& faceLabels)
161 {
162  const labelList& reverseMap = map.reverseFaceMap();
163 
164  labelList newFaceLabels(faceLabels.size());
165  label newI = 0;
166 
167  forAll(faceLabels, i)
168  {
169  label oldFacei = faceLabels[i];
170 
171  if (reverseMap[oldFacei] >= 0)
172  {
173  newFaceLabels[newI++] = reverseMap[oldFacei];
174  }
175  }
176  newFaceLabels.setSize(newI);
177  faceLabels.transfer(newFaceLabels);
178 }
179 
180 
181 void updateCellSet(const mapPolyMesh& map, labelHashSet& cellLabels)
182 {
183  const labelList& reverseMap = map.reverseCellMap();
184 
185  labelHashSet newCellLabels(2*cellLabels.size());
186 
187  forAll(cellLabels, i)
188  {
189  label oldCelli = cellLabels[i];
190 
191  if (reverseMap[oldCelli] >= 0)
192  {
193  newCellLabels.insert(reverseMap[oldCelli]);
194  }
195  }
196  cellLabels.transfer(newCellLabels);
197 }
198 
199 
200 template<class PatchType>
201 void changeFrontBackPatches
202 (
203  polyMesh& mesh,
204  const word& frontPatchName,
205  const word& backPatchName
206 )
207 {
209 
210  label frontPatchi = findPatchID(patches, frontPatchName);
211  label backPatchi = findPatchID(patches, backPatchName);
212 
213  DynamicList<polyPatch*> newPatches(patches.size());
214 
215  forAll(patches, patchi)
216  {
217  const polyPatch& pp(patches[patchi]);
218 
219  if (patchi == frontPatchi || patchi == backPatchi)
220  {
221  newPatches.append
222  (
223  new PatchType
224  (
225  pp.name(),
226  pp.size(),
227  pp.start(),
228  pp.index(),
229  mesh.boundaryMesh(),
230  PatchType::typeName
231  )
232  );
233  }
234  else
235  {
236  newPatches.append(pp.clone(mesh.boundaryMesh()).ptr());
237  }
238  }
239 
240  // Edit patches
242  mesh.addPatches(newPatches, true);
243 }
244 
245 
246 int main(int argc, char *argv[])
247 {
249  (
250  "Extrude mesh from existing patch."
251  );
252 
253  #include "addRegionOption.H"
254  argList::addOption("dict", "file", "Alternative extrudeMeshDict");
255  #include "setRootCase.H"
256  #include "createTimeExtruded.H"
257 
258  // Specified region or default region
259  #include "getRegionOption.H"
260 
261  {
262  Info<< "Create mesh";
263  if (!polyMesh::regionName(regionName).empty())
264  {
265  Info<< ' ' << regionName;
266  }
267  Info<< " for time = " << runTimeExtruded.timeName() << nl << nl;
268  }
269 
270  const IOdictionary dict
271  (
273  (
274  IOobject
275  (
276  "extrudeMeshDict",
277  runTimeExtruded.system(),
278  runTimeExtruded,
281  ),
282  args.getOrDefault<fileName>("dict", "")
283  )
284  );
285 
286  // Point generator
288 
289  // Whether to flip normals
290  const bool flipNormals(dict.get<bool>("flipNormals"));
291 
292  // What to extrude
293  const ExtrudeMode mode = ExtrudeModeNames.get("constructFrom", dict);
294 
295  // Any merging of small edges
296  const scalar mergeTol(dict.getOrDefault<scalar>("mergeTol", 1e-4));
297 
298  Info<< "Extruding from " << ExtrudeModeNames[mode]
299  << " using model " << model().type() << endl;
300  if (flipNormals)
301  {
302  Info<< "Flipping normals before extruding" << endl;
303  }
304  if (mergeTol > 0)
305  {
306  Info<< "Collapsing edges < " << mergeTol << " of bounding box" << endl;
307  }
308  else
309  {
310  Info<< "Not collapsing any edges after extrusion" << endl;
311  }
312  Info<< endl;
313 
314 
315  // Generated mesh (one of either)
316  autoPtr<fvMesh> meshFromMesh;
317  autoPtr<polyMesh> meshFromSurface;
318 
319  // Faces on front and back for stitching (in case of mergeFaces)
320  word frontPatchName;
321  labelList frontPatchFaces;
322  word backPatchName;
323  labelList backPatchFaces;
324 
325  // Optional added cells (get written to cellSet)
326  labelHashSet addedCellsSet;
327 
328  // Optional refinement data
329  autoPtr<hexRef8Data> refDataPtr;
330 
331  if (mode == PATCH || mode == MESH)
332  {
333  if (flipNormals && mode == MESH)
334  {
336  << "Flipping normals not supported for extrusions from mesh."
337  << exit(FatalError);
338  }
339 
340  fileName sourceCasePath(dict.get<fileName>("sourceCase").expand());
341  fileName sourceRootDir = sourceCasePath.path();
342  fileName sourceCaseDir = sourceCasePath.name();
343  if (Pstream::parRun())
344  {
345  sourceCaseDir =
346  sourceCaseDir/("processor" + Foam::name(Pstream::myProcNo()));
347  }
348  wordList sourcePatches;
349  wordList sourceFaceZones;
350  if
351  (
353  (
354  "sourcePatches",
355  sourcePatches,
357  )
358  )
359  {
360  if (sourcePatches.size() == 1)
361  {
362  frontPatchName = sourcePatches[0];
363  }
364 
365  Info<< "Extruding patches " << sourcePatches
366  << " on mesh " << sourceCasePath << nl
367  << endl;
368  }
369  else
370  {
371  dict.readEntry("sourceFaceZones", sourceFaceZones);
372  Info<< "Extruding faceZones " << sourceFaceZones
373  << " on mesh " << sourceCasePath << nl
374  << endl;
375  }
376 
377 
378  Time runTime
379  (
381  sourceRootDir,
382  sourceCaseDir
383  );
384 
385  #include "createNamedMesh.H"
386 
388 
389 
390  // Extrusion engine. Either adding to existing mesh or
391  // creating separate mesh.
392  addPatchCellLayer layerExtrude(mesh, (mode == MESH));
393 
394  labelList meshFaces;
395  bitSet faceFlips;
396  if (sourceFaceZones.size())
397  {
398  zoneFaces(mesh.faceZones(), sourceFaceZones, meshFaces, faceFlips);
399  }
400  else
401  {
402  meshFaces = patchFaces(patches, sourcePatches);
403  faceFlips.setSize(meshFaces.size());
404  }
405 
406  if (mode == PATCH && flipNormals)
407  {
408  // Cheat. Flip patch faces in mesh. This invalidates the
409  // mesh (open cells) but does produce the correct extrusion.
410  polyTopoChange meshMod(mesh);
411  forAll(meshFaces, i)
412  {
413  label meshFacei = meshFaces[i];
414 
415  label patchi = patches.whichPatch(meshFacei);
416  label own = mesh.faceOwner()[meshFacei];
417  label nei = -1;
418  if (patchi == -1)
419  {
420  nei = mesh.faceNeighbour()[meshFacei];
421  }
422 
423  label zoneI = mesh.faceZones().whichZone(meshFacei);
424  bool zoneFlip = false;
425  if (zoneI != -1)
426  {
427  label index = mesh.faceZones()[zoneI].whichFace(meshFacei);
428  zoneFlip = mesh.faceZones()[zoneI].flipMap()[index];
429  }
430 
431  meshMod.modifyFace
432  (
433  mesh.faces()[meshFacei].reverseFace(), // modified face
434  meshFacei, // label of face
435  own, // owner
436  nei, // neighbour
437  true, // face flip
438  patchi, // patch for face
439  zoneI, // zone for face
440  zoneFlip // face flip in zone
441  );
442  }
443 
444  // Change the mesh. No inflation.
445  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
446 
447  // Update fields
448  mesh.updateMesh(map());
449 
450  // Move mesh (since morphing does not do this)
451  if (map().hasMotionPoints())
452  {
453  mesh.movePoints(map().preMotionPoints());
454  }
455  }
456 
457 
458  indirectPrimitivePatch extrudePatch
459  (
461  (
462  mesh.faces(),
463  meshFaces
464  ),
465  mesh.points()
466  );
467 
468  // Determine extrudePatch normal
469  pointField extrudePatchPointNormals
470  (
471  PatchTools::pointNormals(mesh, extrudePatch, faceFlips)
472  );
473 
474 
475  // Precalculate mesh edges for pp.edges.
476  const labelList meshEdges
477  (
478  extrudePatch.meshEdges
479  (
480  mesh.edges(),
481  mesh.pointEdges()
482  )
483  );
484 
485  // Global face indices engine
486  const globalIndex globalFaces(mesh.nFaces());
487 
488  // Determine extrudePatch.edgeFaces in global numbering (so across
489  // coupled patches)
490  labelListList edgeGlobalFaces
491  (
493  (
494  mesh,
495  globalFaces,
496  extrudePatch
497  )
498  );
499 
500 
501  // Determine what patches boundary edges need to get extruded into.
502  // This might actually cause edge-connected processors to become
503  // face-connected so might need to introduce new processor boundaries.
504  // Calculates:
505  // - per pp.edge the patch to extrude into
506  // - any additional processor boundaries (patchToNbrProc = map from
507  // new patchID to neighbour processor)
508  // - number of new patches (nPatches)
509 
510  labelList edgePatchID;
511  labelList edgeZoneID;
512  boolList edgeFlip;
513  labelList inflateFaceID;
514  label nPatches;
515  Map<label> nbrProcToPatch;
516  Map<label> patchToNbrProc;
518  (
519  true, // for internal edges get zone info from any face
520 
521  mesh,
522  globalFaces,
523  edgeGlobalFaces,
524  extrudePatch,
525 
526  edgePatchID,
527  nPatches,
528  nbrProcToPatch,
529  patchToNbrProc,
530 
531  edgeZoneID,
532  edgeFlip,
533  inflateFaceID
534  );
535 
536 
537  // Add any patches.
538 
539  const label nAdded = returnReduce
540  (
542  sumOp<label>()
543  );
544 
545  Info<< "Adding overall " << nAdded << " processor patches." << endl;
546 
547  if (nAdded > 0)
548  {
549  DynamicList<polyPatch*> newPatches(nPatches);
550  forAll(mesh.boundaryMesh(), patchi)
551  {
552  newPatches.append
553  (
554  mesh.boundaryMesh()[patchi].clone
555  (
557  ).ptr()
558  );
559  }
560  for
561  (
562  label patchi = mesh.boundaryMesh().size();
563  patchi < nPatches;
564  patchi++
565  )
566  {
567  label nbrProci = patchToNbrProc[patchi];
568 
569  Pout<< "Adding patch " << patchi
570  << " between " << Pstream::myProcNo()
571  << " and " << nbrProci
572  << endl;
573 
574  newPatches.append
575  (
577  (
578  0, // size
579  mesh.nFaces(), // start
580  patchi, // index
581  mesh.boundaryMesh(),// polyBoundaryMesh
582  Pstream::myProcNo(),// myProcNo
583  nbrProci // neighbProcNo
584  )
585  );
586  }
587 
588  // Add patches. Do no parallel updates.
590  mesh.addFvPatches(newPatches, true);
591  }
592 
593 
594 
595  // Only used for addPatchCellLayer into new mesh
596  labelList exposedPatchID;
597  if (mode == PATCH)
598  {
599  dict.readEntry("exposedPatchName", backPatchName);
600  exposedPatchID.setSize
601  (
602  extrudePatch.size(),
603  findPatchID(patches, backPatchName)
604  );
605  }
606 
607  // Determine points and extrusion
608  pointField layer0Points(extrudePatch.nPoints());
609  pointField displacement(extrudePatch.nPoints());
610  forAll(displacement, pointi)
611  {
612  const vector& patchNormal = extrudePatchPointNormals[pointi];
613 
614  // layer0 point
615  layer0Points[pointi] = model()
616  (
617  extrudePatch.localPoints()[pointi],
618  patchNormal,
619  0
620  );
621  // layerN point
622  point extrudePt = model()
623  (
624  extrudePatch.localPoints()[pointi],
625  patchNormal,
626  model().nLayers()
627  );
628  displacement[pointi] = extrudePt - layer0Points[pointi];
629  }
630 
631 
632  // Check if wedge (has layer0 different from original patch points)
633  // If so move the mesh to starting position.
634  if (gMax(mag(layer0Points-extrudePatch.localPoints())) > SMALL)
635  {
636  Info<< "Moving mesh to layer0 points since differ from original"
637  << " points - this can happen for wedge extrusions." << nl
638  << endl;
639 
640  pointField newPoints(mesh.points());
641  forAll(extrudePatch.meshPoints(), i)
642  {
643  newPoints[extrudePatch.meshPoints()[i]] = layer0Points[i];
644  }
645  mesh.movePoints(newPoints);
646  }
647 
648 
649  // Layers per face
650  labelList nFaceLayers(extrudePatch.size(), model().nLayers());
651 
652  // Layers per point
653  labelList nPointLayers(extrudePatch.nPoints(), model().nLayers());
654 
655  // Displacement for first layer
656  vectorField firstLayerDisp(displacement*model().sumThickness(1));
657 
658  // Expansion ratio not used.
659  scalarField ratio(extrudePatch.nPoints(), 1.0);
660 
661 
662  // Load any refinement data
663  if (mode == MESH)
664  {
665  // Tricky: register hexRef8 onto the database
666  // since the mesh does not yet exist. It should not be registered
667  // onto the input mesh.
668  refDataPtr.reset
669  (
670  new hexRef8Data
671  (
672  IOobject
673  (
674  "dummy",
677  runTimeExtruded, //mesh,
681  )
682  )
683  );
684  }
685 
686 
687  // Topo change container. Either copy an existing mesh or start
688  // with empty storage (number of patches only needed for checking)
690  (
691  (
692  mode == MESH
693  ? new polyTopoChange(mesh)
694  : new polyTopoChange(patches.size())
695  )
696  );
697 
698  layerExtrude.setRefinement
699  (
700  globalFaces,
701  edgeGlobalFaces,
702 
703  ratio, // expansion ratio
704  extrudePatch, // patch faces to extrude
705  faceFlips, // side to extrude (for internal/coupled faces)
706 
707  edgePatchID, // if boundary edge: patch for extruded face
708  edgeZoneID, // optional zone for extruded face
709  edgeFlip,
710  inflateFaceID, // mesh face that zone/patch info is from
711 
712  exposedPatchID, // if new mesh: patches for exposed faces
713  nFaceLayers,
714  nPointLayers,
715  firstLayerDisp,
716  meshMod()
717  );
718 
719  // Reset points according to extrusion model
720  forAll(layerExtrude.addedPoints(), pointi)
721  {
722  const labelList& pPoints = layerExtrude.addedPoints()[pointi];
723  forAll(pPoints, pPointi)
724  {
725  label meshPointi = pPoints[pPointi];
726 
727  point modelPt
728  (
729  model()
730  (
731  extrudePatch.localPoints()[pointi],
732  extrudePatchPointNormals[pointi],
733  pPointi+1 // layer
734  )
735  );
736 
737  const_cast<DynamicList<point>&>
738  (
739  meshMod().points()
740  )[meshPointi] = modelPt;
741  }
742  }
743 
744  // Store faces on front and exposed patch (if mode=patch there are
745  // only added faces so cannot used map to old faces)
746  const labelListList& layerFaces = layerExtrude.layerFaces();
747  backPatchFaces.setSize(layerFaces.size());
748  frontPatchFaces.setSize(layerFaces.size());
749  forAll(backPatchFaces, patchFacei)
750  {
751  backPatchFaces[patchFacei] = layerFaces[patchFacei].first();
752  frontPatchFaces[patchFacei] = layerFaces[patchFacei].last();
753  }
754 
755 
756  // Create dummy fvSchemes, fvSolution
758  (
759  mesh,
761  true
762  );
763 
764  // Create actual mesh from polyTopoChange container
765  autoPtr<mapPolyMesh> map = meshMod().makeMesh
766  (
767  meshFromMesh,
768  IOobject
769  (
770  regionName,
771  runTimeExtruded.constant(),
772  runTimeExtruded,
773  IOobject::READ_IF_PRESENT, // Read fv* if present
776  ),
777  mesh
778  );
779 
780  layerExtrude.updateMesh
781  (
782  map(),
783  identity(extrudePatch.size()),
784  identity(extrudePatch.nPoints())
785  );
786 
787  // Calculate face labels for front and back.
788  frontPatchFaces = renumber
789  (
790  map().reverseFaceMap(),
791  frontPatchFaces
792  );
793  backPatchFaces = renumber
794  (
795  map().reverseFaceMap(),
796  backPatchFaces
797  );
798 
799  // Update
800  if (refDataPtr)
801  {
802  refDataPtr->updateMesh(map());
803  }
804 
805  // Store added cells
806  if (mode == MESH)
807  {
808  const labelListList addedCells
809  (
810  layerExtrude.addedCells
811  (
812  *meshFromMesh,
813  layerExtrude.layerFaces()
814  )
815  );
816  forAll(addedCells, facei)
817  {
818  const labelList& aCells = addedCells[facei];
819  addedCellsSet.insert(aCells);
820  }
821  }
822  }
823  else
824  {
825  // Read from surface
826  fileName surfName(dict.get<fileName>("surface").expand());
827 
828  Info<< "Extruding surfaceMesh read from file " << surfName << nl
829  << endl;
830 
831  MeshedSurface<face> fMesh(surfName);
832 
833  if (flipNormals)
834  {
835  Info<< "Flipping faces." << nl << endl;
836  faceList& faces = const_cast<faceList&>(fMesh.surfFaces());
837  forAll(faces, i)
838  {
839  faces[i] = fMesh[i].reverseFace();
840  }
841  }
842 
843  Info<< "Extruding surface with :" << nl
844  << " points : " << fMesh.points().size() << nl
845  << " faces : " << fMesh.size() << nl
846  << " normals[0] : " << fMesh.faceNormals()[0]
847  << nl
848  << endl;
849 
850  meshFromSurface.reset
851  (
852  new extrudedMesh
853  (
854  IOobject
855  (
857  runTimeExtruded.constant(),
858  runTimeExtruded
859  ),
860  fMesh,
861  model()
862  )
863  );
864 
865 
866  // Get the faces on front and back
867  frontPatchName = "originalPatch";
868  frontPatchFaces = patchFaces
869  (
870  meshFromSurface().boundaryMesh(),
871  wordList(1, frontPatchName)
872  );
873  backPatchName = "otherSide";
874  backPatchFaces = patchFaces
875  (
876  meshFromSurface().boundaryMesh(),
877  wordList(1, backPatchName)
878  );
879  }
880 
881 
882  polyMesh& mesh =
883  (
884  meshFromMesh
885  ? *meshFromMesh
886  : *meshFromSurface
887  );
888 
889 
890  const boundBox& bb = mesh.bounds();
891  const scalar mergeDim = mergeTol * bb.minDim();
892 
893  Info<< "Mesh bounding box : " << bb << nl
894  << " with span : " << bb.span() << nl
895  << "Merge distance : " << mergeDim << nl
896  << endl;
897 
898 
899  // Collapse edges
900  // ~~~~~~~~~~~~~~
901 
902  if (mergeDim > 0)
903  {
904  Info<< "Collapsing edges < " << mergeDim << " ..." << nl << endl;
905 
906  // Edge collapsing engine
907  edgeCollapser collapser(mesh);
908 
909  const edgeList& edges = mesh.edges();
910  const pointField& points = mesh.points();
911 
913  Map<point> collapsePointToLocation(mesh.nPoints());
914 
915  forAll(edges, edgeI)
916  {
917  const edge& e = edges[edgeI];
918 
919  scalar d = e.mag(points);
920 
921  if (d < mergeDim)
922  {
923  Info<< "Merging edge " << e << " since length " << d
924  << " << " << mergeDim << nl;
925 
926  collapseEdge.set(edgeI);
927  collapsePointToLocation.set(e[1], points[e[0]]);
928  }
929  }
930 
931  List<pointEdgeCollapse> allPointInfo;
933  labelList pointPriority(mesh.nPoints(), Zero);
934 
935  collapser.consistentCollapse
936  (
937  globalPoints,
938  pointPriority,
939  collapsePointToLocation,
940  collapseEdge,
941  allPointInfo
942  );
943 
944  // Topo change container
945  polyTopoChange meshMod(mesh);
946 
947  // Put all modifications into meshMod
948  bool anyChange = collapser.setRefinement(allPointInfo, meshMod);
949 
950  if (returnReduceOr(anyChange))
951  {
952  // Construct new mesh from polyTopoChange.
953  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
954 
955  // Update fields
956  mesh.updateMesh(map());
957 
958  // Update stored data
959  updateFaceLabels(map(), frontPatchFaces);
960  updateFaceLabels(map(), backPatchFaces);
961  updateCellSet(map(), addedCellsSet);
962 
963  if (refDataPtr)
964  {
965  refDataPtr->updateMesh(map());
966  }
967 
968  // Move mesh (if inflation used)
969  if (map().hasMotionPoints())
970  {
971  mesh.movePoints(map().preMotionPoints());
972  }
973  }
974  }
975 
976 
977  // Change the front and back patch types as required
978  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
979 
980  word frontBackType;
981 
982  if (isType<extrudeModels::wedge>(model()))
983  {
984  changeFrontBackPatches<wedgePolyPatch>
985  (
986  mesh,
987  frontPatchName,
988  backPatchName
989  );
990  }
991  else if (isType<extrudeModels::plane>(model()))
992  {
993  changeFrontBackPatches<emptyPolyPatch>
994  (
995  mesh,
996  frontPatchName,
997  backPatchName
998  );
999  }
1000 
1001 
1002  // Merging front and back patch faces
1003  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1004 
1005  if (dict.get<bool>("mergeFaces"))
1006  {
1007  if (mode == MESH)
1008  {
1010  << "Cannot stitch front and back of extrusion since"
1011  << " in 'mesh' mode (extrusion appended to mesh)."
1012  << exit(FatalError);
1013  }
1014 
1015  Info<< "Assuming full 360 degree axisymmetric case;"
1016  << " stitching faces on patches "
1017  << frontPatchName << " and "
1018  << backPatchName << " together ..." << nl << endl;
1019 
1020  if (frontPatchFaces.size() != backPatchFaces.size())
1021  {
1023  << "Differing number of faces on front ("
1024  << frontPatchFaces.size() << ") and back ("
1025  << backPatchFaces.size() << ")"
1026  << exit(FatalError);
1027  }
1028 
1029 
1030  polyTopoChanger stitcher(mesh);
1031  stitcher.setSize(1);
1032 
1033  const word cutZoneName("originalCutFaceZone");
1034 
1035  List<faceZone*> fz
1036  (
1037  1,
1038  new faceZone
1039  (
1040  cutZoneName,
1041  frontPatchFaces,
1042  false, // none are flipped
1043  0,
1044  mesh.faceZones()
1045  )
1046  );
1047 
1049 
1050  // Add the perfect interface mesh modifier
1051  perfectInterface perfectStitcher
1052  (
1053  "couple",
1054  0,
1055  stitcher,
1056  cutZoneName,
1057  word::null, // dummy patch name
1058  word::null // dummy patch name
1059  );
1060 
1061  // Topo change container
1062  polyTopoChange meshMod(mesh);
1063 
1064  perfectStitcher.setRefinement
1065  (
1067  (
1069  (
1070  mesh.faces(),
1071  frontPatchFaces
1072  ),
1073  mesh.points()
1074  ),
1076  (
1078  (
1079  mesh.faces(),
1080  backPatchFaces
1081  ),
1082  mesh.points()
1083  ),
1084  meshMod
1085  );
1086 
1087  // Construct new mesh from polyTopoChange.
1088  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
1089 
1090  // Update fields
1091  mesh.updateMesh(map());
1092 
1093  // Update local data
1094  updateCellSet(map(), addedCellsSet);
1095 
1096  if (refDataPtr)
1097  {
1098  refDataPtr->updateMesh(map());
1099  }
1100 
1101  // Move mesh (if inflation used)
1102  if (map().hasMotionPoints())
1103  {
1104  mesh.movePoints(map().preMotionPoints());
1105  }
1106  }
1107 
1108  mesh.setInstance(runTimeExtruded.constant());
1109  Info<< "Writing mesh to " << mesh.objectPath() << nl << endl;
1110 
1111  if (!mesh.write())
1112  {
1114  << exit(FatalError);
1115  }
1116  // Remove any left-over files
1118 
1119  // Need writing cellSet
1120  if (returnReduceOr(addedCellsSet.size()))
1121  {
1122  cellSet addedCells(mesh, "addedCells", addedCellsSet);
1123  Info<< "Writing added cells to cellSet " << addedCells.name()
1124  << nl << endl;
1125  if (!addedCells.write())
1126  {
1128  << addedCells.name()
1129  << exit(FatalError);
1130  }
1131  }
1132 
1133  if (refDataPtr)
1134  {
1135  refDataPtr->write();
1136  }
1137 
1138 
1139  Info<< "End\n" << endl;
1140 
1141  return 0;
1142 }
1143 
1144 
1145 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
dictionary dict
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to these fvPatches.
Definition: fvMesh.C:692
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
A class for handling file names.
Definition: fileName.H:72
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:859
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:847
const labelListList & pointEdges() const
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values (not the indices) of a list.
PtrList< T > clone(Args &&... args) const
Make a copy by cloning each of the list elements.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void setSize(const label n, unsigned int val=0u)
Alias for resize()
Definition: PackedList.H:785
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:410
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:918
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)
Determine extrude information per patch edge:
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
Does polyTopoChanges to remove edges. Can remove faces due to edge collapse but can not remove cells ...
Definition: edgeCollapser.H:64
Required Classes.
label collapseEdge(triSurface &surf, const scalar minLen)
Keep collapsing all edges < minLen.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Required Classes.
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:169
Ignore writing from objectRegistry::writeObject()
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:300
label nFaces() const noexcept
Number of mesh faces.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
labelList faceLabels(nFaceLabels)
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:284
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
Neighbour processor patch.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:37
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:342
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:61
Required Classes.
PDRpatchDef PATCH
Definition: PDRpatchDef.H:120
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:994
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:629
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
List of mesh modifiers defining the mesh dynamics.
Adds layers of cells to outside of polyPatch. Can optionally create stand-alone extruded mesh (addToM...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:405
wordList names() const
Return a list of patch names.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
Definition: polyMesh.C:1009
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
Reading is optional [identical to LAZY_READ].
static const word null
An empty word.
Definition: word.H:84
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:385
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
String literal.
Definition: keyType.H:82
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:268
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Create additional fvSchemes/fvSolution files.
Definition: fvMeshTools.C:1212
void updateMesh(const mapPolyMesh &)
Update local numbering for changed mesh.
Definition: hexRef8Data.C:308
label nEdges() const
Number of mesh edges.
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1102
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:32
bool write() const
Write.
Definition: hexRef8Data.C:399
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:192
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:653
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Type gMax(const FieldField< Field, Type > &f)
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:662
static autoPtr< extrudeModel > New(const dictionary &dict)
Select null constructed.
T & last()
Access last element of the list, position [size()-1].
Definition: UList.H:867
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:346
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
List< word > wordList
List of word.
Definition: fileName.H:59
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:98
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
Definition: polyMeshIO.C:29
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
string & expand(const bool allowEmpty=false)
Inplace expand initial tags, tildes, and all occurrences of environment variables as per stringOps::e...
Definition: string.C:166
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
A collection of cell labels.
Definition: cellSet.H:47
Direct mesh changes based on v1.3 polyTopoChange syntax.
void transfer(HashTable< T, Key, Hash > &rhs)
Transfer contents into this table.
Definition: HashTable.C:777
const polyBoundaryMesh & patches
Automatically write from objectRegistry::writeObject()
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:192
messageStream Info
Information stream (stdout output on master, null elsewhere)
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition: polyMesh.H:616
static labelListList globalEdgeFaces(const polyMesh &, const globalIndex &globalFaces, const indirectPrimitivePatch &pp)
Per patch edge the pp faces (in global indices) using it.
label n
Various for reading/decomposing/reconstructing/distributing refinement data.
Definition: hexRef8Data.H:56
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Hack of attachDetach to couple patches when they perfectly align. Does not decouple. Used by stitchMesh app. Does geometric matching.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &, const bitSet &flipMap=bitSet::null())
Return parallel consistent point normals for patches using mesh points.
void addPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:971
Foam::argList args(argc, argv)
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: POSIX.C:773
Bookkeeping for patch definitions.
Definition: PDRpatchDef.H:49
A List with indirect addressing.
Definition: IndirectList.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
scalar minDim() const
Smallest length/height/width dimension.
Definition: boundBoxI.H:216
Do not request registration (bool: false)
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:617
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static IOobject selectIO(const IOobject &io, const fileName &altFile, const word &ioName="")
Return the IOobject, but also consider an alternative file name.
Definition: IOobject.C:256
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127