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  }
112  labelList faceLabels(n);
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,
132  labelList& faceLabels,
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  }
143  faceLabels.setSize(n);
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  // Get optional regionName
260  if (args.readIfPresent("region", regionName))
261  {
262  Info<< "Create mesh " << regionName << " for time = "
263  << runTimeExtruded.timeName() << nl << endl;
264  }
265  else
266  {
267  Info<< "Create mesh for time = "
268  << runTimeExtruded.timeName() << nl << endl;
269  }
270 
271 
272  const IOdictionary dict
273  (
275  (
276  IOobject
277  (
278  "extrudeMeshDict",
279  runTimeExtruded.system(),
280  runTimeExtruded,
283  ),
284  args.getOrDefault<fileName>("dict", "")
285  )
286  );
287 
288  // Point generator
290 
291  // Whether to flip normals
292  const bool flipNormals(dict.get<bool>("flipNormals"));
293 
294  // What to extrude
295  const ExtrudeMode mode = ExtrudeModeNames.get("constructFrom", dict);
296 
297  // Any merging of small edges
298  const scalar mergeTol(dict.getOrDefault<scalar>("mergeTol", 1e-4));
299 
300  Info<< "Extruding from " << ExtrudeModeNames[mode]
301  << " using model " << model().type() << endl;
302  if (flipNormals)
303  {
304  Info<< "Flipping normals before extruding" << endl;
305  }
306  if (mergeTol > 0)
307  {
308  Info<< "Collapsing edges < " << mergeTol << " of bounding box" << endl;
309  }
310  else
311  {
312  Info<< "Not collapsing any edges after extrusion" << endl;
313  }
314  Info<< endl;
315 
316 
317  // Generated mesh (one of either)
318  autoPtr<fvMesh> meshFromMesh;
319  autoPtr<polyMesh> meshFromSurface;
320 
321  // Faces on front and back for stitching (in case of mergeFaces)
322  word frontPatchName;
323  labelList frontPatchFaces;
324  word backPatchName;
325  labelList backPatchFaces;
326 
327  // Optional added cells (get written to cellSet)
328  labelHashSet addedCellsSet;
329 
330  // Optional refinement data
331  autoPtr<hexRef8Data> refDataPtr;
332 
333  if (mode == PATCH || mode == MESH)
334  {
335  if (flipNormals && mode == MESH)
336  {
338  << "Flipping normals not supported for extrusions from mesh."
339  << exit(FatalError);
340  }
341 
342  fileName sourceCasePath(dict.get<fileName>("sourceCase").expand());
343  fileName sourceRootDir = sourceCasePath.path();
344  fileName sourceCaseDir = sourceCasePath.name();
345  if (Pstream::parRun())
346  {
347  sourceCaseDir =
348  sourceCaseDir/("processor" + Foam::name(Pstream::myProcNo()));
349  }
350  wordList sourcePatches;
351  wordList sourceFaceZones;
352  if
353  (
355  (
356  "sourcePatches",
357  sourcePatches,
359  )
360  )
361  {
362  if (sourcePatches.size() == 1)
363  {
364  frontPatchName = sourcePatches[0];
365  }
366 
367  Info<< "Extruding patches " << sourcePatches
368  << " on mesh " << sourceCasePath << nl
369  << endl;
370  }
371  else
372  {
373  dict.readEntry("sourceFaceZones", sourceFaceZones);
374  Info<< "Extruding faceZones " << sourceFaceZones
375  << " on mesh " << sourceCasePath << nl
376  << endl;
377  }
378 
379 
380  Time runTime
381  (
383  sourceRootDir,
384  sourceCaseDir
385  );
386 
387  #include "createNamedMesh.H"
388 
390 
391 
392  // Extrusion engine. Either adding to existing mesh or
393  // creating separate mesh.
394  addPatchCellLayer layerExtrude(mesh, (mode == MESH));
395 
396  labelList meshFaces;
397  bitSet faceFlips;
398  if (sourceFaceZones.size())
399  {
400  zoneFaces(mesh.faceZones(), sourceFaceZones, meshFaces, faceFlips);
401  }
402  else
403  {
404  meshFaces = patchFaces(patches, sourcePatches);
405  faceFlips.setSize(meshFaces.size());
406  }
407 
408  if (mode == PATCH && flipNormals)
409  {
410  // Cheat. Flip patch faces in mesh. This invalidates the
411  // mesh (open cells) but does produce the correct extrusion.
412  polyTopoChange meshMod(mesh);
413  forAll(meshFaces, i)
414  {
415  label meshFacei = meshFaces[i];
416 
417  label patchi = patches.whichPatch(meshFacei);
418  label own = mesh.faceOwner()[meshFacei];
419  label nei = -1;
420  if (patchi == -1)
421  {
422  nei = mesh.faceNeighbour()[meshFacei];
423  }
424 
425  label zoneI = mesh.faceZones().whichZone(meshFacei);
426  bool zoneFlip = false;
427  if (zoneI != -1)
428  {
429  label index = mesh.faceZones()[zoneI].whichFace(meshFacei);
430  zoneFlip = mesh.faceZones()[zoneI].flipMap()[index];
431  }
432 
433  meshMod.modifyFace
434  (
435  mesh.faces()[meshFacei].reverseFace(), // modified face
436  meshFacei, // label of face
437  own, // owner
438  nei, // neighbour
439  true, // face flip
440  patchi, // patch for face
441  zoneI, // zone for face
442  zoneFlip // face flip in zone
443  );
444  }
445 
446  // Change the mesh. No inflation.
447  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
448 
449  // Update fields
450  mesh.updateMesh(map());
451 
452  // Move mesh (since morphing does not do this)
453  if (map().hasMotionPoints())
454  {
455  mesh.movePoints(map().preMotionPoints());
456  }
457  }
458 
459 
460  indirectPrimitivePatch extrudePatch
461  (
463  (
464  mesh.faces(),
465  meshFaces
466  ),
467  mesh.points()
468  );
469 
470  // Determine extrudePatch normal
471  pointField extrudePatchPointNormals
472  (
473  PatchTools::pointNormals(mesh, extrudePatch, faceFlips)
474  );
475 
476 
477  // Precalculate mesh edges for pp.edges.
478  const labelList meshEdges
479  (
480  extrudePatch.meshEdges
481  (
482  mesh.edges(),
483  mesh.pointEdges()
484  )
485  );
486 
487  // Global face indices engine
488  const globalIndex globalFaces(mesh.nFaces());
489 
490  // Determine extrudePatch.edgeFaces in global numbering (so across
491  // coupled patches)
492  labelListList edgeGlobalFaces
493  (
495  (
496  mesh,
497  globalFaces,
498  extrudePatch
499  )
500  );
501 
502 
503  // Determine what patches boundary edges need to get extruded into.
504  // This might actually cause edge-connected processors to become
505  // face-connected so might need to introduce new processor boundaries.
506  // Calculates:
507  // - per pp.edge the patch to extrude into
508  // - any additional processor boundaries (patchToNbrProc = map from
509  // new patchID to neighbour processor)
510  // - number of new patches (nPatches)
511 
512  labelList edgePatchID;
513  labelList edgeZoneID;
514  boolList edgeFlip;
515  labelList inflateFaceID;
516  label nPatches;
517  Map<label> nbrProcToPatch;
518  Map<label> patchToNbrProc;
520  (
521  true, // for internal edges get zone info from any face
522 
523  mesh,
524  globalFaces,
525  edgeGlobalFaces,
526  extrudePatch,
527 
528  edgePatchID,
529  nPatches,
530  nbrProcToPatch,
531  patchToNbrProc,
532 
533  edgeZoneID,
534  edgeFlip,
535  inflateFaceID
536  );
537 
538 
539  // Add any patches.
540 
541  const label nAdded = returnReduce
542  (
544  sumOp<label>()
545  );
546 
547  Info<< "Adding overall " << nAdded << " processor patches." << endl;
548 
549  if (nAdded > 0)
550  {
551  DynamicList<polyPatch*> newPatches(nPatches);
552  forAll(mesh.boundaryMesh(), patchi)
553  {
554  newPatches.append
555  (
556  mesh.boundaryMesh()[patchi].clone
557  (
559  ).ptr()
560  );
561  }
562  for
563  (
564  label patchi = mesh.boundaryMesh().size();
565  patchi < nPatches;
566  patchi++
567  )
568  {
569  label nbrProci = patchToNbrProc[patchi];
570 
571  Pout<< "Adding patch " << patchi
572  << " between " << Pstream::myProcNo()
573  << " and " << nbrProci
574  << endl;
575 
576  newPatches.append
577  (
579  (
580  0, // size
581  mesh.nFaces(), // start
582  patchi, // index
583  mesh.boundaryMesh(),// polyBoundaryMesh
584  Pstream::myProcNo(),// myProcNo
585  nbrProci // neighbProcNo
586  )
587  );
588  }
589 
590  // Add patches. Do no parallel updates.
592  mesh.addFvPatches(newPatches, true);
593  }
594 
595 
596 
597  // Only used for addPatchCellLayer into new mesh
598  labelList exposedPatchID;
599  if (mode == PATCH)
600  {
601  dict.readEntry("exposedPatchName", backPatchName);
602  exposedPatchID.setSize
603  (
604  extrudePatch.size(),
605  findPatchID(patches, backPatchName)
606  );
607  }
608 
609  // Determine points and extrusion
610  pointField layer0Points(extrudePatch.nPoints());
611  pointField displacement(extrudePatch.nPoints());
612  forAll(displacement, pointi)
613  {
614  const vector& patchNormal = extrudePatchPointNormals[pointi];
615 
616  // layer0 point
617  layer0Points[pointi] = model()
618  (
619  extrudePatch.localPoints()[pointi],
620  patchNormal,
621  0
622  );
623  // layerN point
624  point extrudePt = model()
625  (
626  extrudePatch.localPoints()[pointi],
627  patchNormal,
628  model().nLayers()
629  );
630  displacement[pointi] = extrudePt - layer0Points[pointi];
631  }
632 
633 
634  // Check if wedge (has layer0 different from original patch points)
635  // If so move the mesh to starting position.
636  if (gMax(mag(layer0Points-extrudePatch.localPoints())) > SMALL)
637  {
638  Info<< "Moving mesh to layer0 points since differ from original"
639  << " points - this can happen for wedge extrusions." << nl
640  << endl;
641 
642  pointField newPoints(mesh.points());
643  forAll(extrudePatch.meshPoints(), i)
644  {
645  newPoints[extrudePatch.meshPoints()[i]] = layer0Points[i];
646  }
647  mesh.movePoints(newPoints);
648  }
649 
650 
651  // Layers per face
652  labelList nFaceLayers(extrudePatch.size(), model().nLayers());
653 
654  // Layers per point
655  labelList nPointLayers(extrudePatch.nPoints(), model().nLayers());
656 
657  // Displacement for first layer
658  vectorField firstLayerDisp(displacement*model().sumThickness(1));
659 
660  // Expansion ratio not used.
661  scalarField ratio(extrudePatch.nPoints(), 1.0);
662 
663 
664  // Load any refinement data
665  if (mode == MESH)
666  {
667  // Tricky: register hexRef8 onto the database
668  // since the mesh does not yet exist. It should not be registered
669  // onto the input mesh.
670  refDataPtr.reset
671  (
672  new hexRef8Data
673  (
674  IOobject
675  (
676  "dummy",
679  runTimeExtruded, //mesh,
682  false
683  )
684  )
685  );
686  }
687 
688 
689  // Topo change container. Either copy an existing mesh or start
690  // with empty storage (number of patches only needed for checking)
692  (
693  (
694  mode == MESH
695  ? new polyTopoChange(mesh)
696  : new polyTopoChange(patches.size())
697  )
698  );
699 
700  layerExtrude.setRefinement
701  (
702  globalFaces,
703  edgeGlobalFaces,
704 
705  ratio, // expansion ratio
706  extrudePatch, // patch faces to extrude
707  faceFlips, // side to extrude (for internal/coupled faces)
708 
709  edgePatchID, // if boundary edge: patch for extruded face
710  edgeZoneID, // optional zone for extruded face
711  edgeFlip,
712  inflateFaceID, // mesh face that zone/patch info is from
713 
714  exposedPatchID, // if new mesh: patches for exposed faces
715  nFaceLayers,
716  nPointLayers,
717  firstLayerDisp,
718  meshMod()
719  );
720 
721  // Reset points according to extrusion model
722  forAll(layerExtrude.addedPoints(), pointi)
723  {
724  const labelList& pPoints = layerExtrude.addedPoints()[pointi];
725  forAll(pPoints, pPointi)
726  {
727  label meshPointi = pPoints[pPointi];
728 
729  point modelPt
730  (
731  model()
732  (
733  extrudePatch.localPoints()[pointi],
734  extrudePatchPointNormals[pointi],
735  pPointi+1 // layer
736  )
737  );
738 
739  const_cast<DynamicList<point>&>
740  (
741  meshMod().points()
742  )[meshPointi] = modelPt;
743  }
744  }
745 
746  // Store faces on front and exposed patch (if mode=patch there are
747  // only added faces so cannot used map to old faces)
748  const labelListList& layerFaces = layerExtrude.layerFaces();
749  backPatchFaces.setSize(layerFaces.size());
750  frontPatchFaces.setSize(layerFaces.size());
751  forAll(backPatchFaces, patchFacei)
752  {
753  backPatchFaces[patchFacei] = layerFaces[patchFacei].first();
754  frontPatchFaces[patchFacei] = layerFaces[patchFacei].last();
755  }
756 
757 
758  // Create dummy fvSchemes, fvSolution
760  (
761  mesh,
763  true
764  );
765 
766  // Create actual mesh from polyTopoChange container
767  autoPtr<mapPolyMesh> map = meshMod().makeMesh
768  (
769  meshFromMesh,
770  IOobject
771  (
772  regionName,
773  runTimeExtruded.constant(),
774  runTimeExtruded,
775  IOobject::READ_IF_PRESENT, // Read fv* if present
777  false
778  ),
779  mesh
780  );
781 
782  layerExtrude.updateMesh
783  (
784  map(),
785  identity(extrudePatch.size()),
786  identity(extrudePatch.nPoints())
787  );
788 
789  // Calculate face labels for front and back.
790  frontPatchFaces = renumber
791  (
792  map().reverseFaceMap(),
793  frontPatchFaces
794  );
795  backPatchFaces = renumber
796  (
797  map().reverseFaceMap(),
798  backPatchFaces
799  );
800 
801  // Update
802  if (refDataPtr)
803  {
804  refDataPtr->updateMesh(map());
805  }
806 
807  // Store added cells
808  if (mode == MESH)
809  {
810  const labelListList addedCells
811  (
812  layerExtrude.addedCells
813  (
814  *meshFromMesh,
815  layerExtrude.layerFaces()
816  )
817  );
818  forAll(addedCells, facei)
819  {
820  const labelList& aCells = addedCells[facei];
821  addedCellsSet.insert(aCells);
822  }
823  }
824  }
825  else
826  {
827  // Read from surface
828  fileName surfName(dict.get<fileName>("surface").expand());
829 
830  Info<< "Extruding surfaceMesh read from file " << surfName << nl
831  << endl;
832 
833  MeshedSurface<face> fMesh(surfName);
834 
835  if (flipNormals)
836  {
837  Info<< "Flipping faces." << nl << endl;
838  faceList& faces = const_cast<faceList&>(fMesh.surfFaces());
839  forAll(faces, i)
840  {
841  faces[i] = fMesh[i].reverseFace();
842  }
843  }
844 
845  Info<< "Extruding surface with :" << nl
846  << " points : " << fMesh.points().size() << nl
847  << " faces : " << fMesh.size() << nl
848  << " normals[0] : " << fMesh.faceNormals()[0]
849  << nl
850  << endl;
851 
852  meshFromSurface.reset
853  (
854  new extrudedMesh
855  (
856  IOobject
857  (
859  runTimeExtruded.constant(),
860  runTimeExtruded
861  ),
862  fMesh,
863  model()
864  )
865  );
866 
867 
868  // Get the faces on front and back
869  frontPatchName = "originalPatch";
870  frontPatchFaces = patchFaces
871  (
872  meshFromSurface().boundaryMesh(),
873  wordList(1, frontPatchName)
874  );
875  backPatchName = "otherSide";
876  backPatchFaces = patchFaces
877  (
878  meshFromSurface().boundaryMesh(),
879  wordList(1, backPatchName)
880  );
881  }
882 
883 
884  polyMesh& mesh =
885  (
886  meshFromMesh
887  ? *meshFromMesh
888  : *meshFromSurface
889  );
890 
891 
892  const boundBox& bb = mesh.bounds();
893  const scalar mergeDim = mergeTol * bb.minDim();
894 
895  Info<< "Mesh bounding box : " << bb << nl
896  << " with span : " << bb.span() << nl
897  << "Merge distance : " << mergeDim << nl
898  << endl;
899 
900 
901  // Collapse edges
902  // ~~~~~~~~~~~~~~
903 
904  if (mergeDim > 0)
905  {
906  Info<< "Collapsing edges < " << mergeDim << " ..." << nl << endl;
907 
908  // Edge collapsing engine
909  edgeCollapser collapser(mesh);
910 
911  const edgeList& edges = mesh.edges();
912  const pointField& points = mesh.points();
913 
915  Map<point> collapsePointToLocation(mesh.nPoints());
916 
917  forAll(edges, edgeI)
918  {
919  const edge& e = edges[edgeI];
920 
921  scalar d = e.mag(points);
922 
923  if (d < mergeDim)
924  {
925  Info<< "Merging edge " << e << " since length " << d
926  << " << " << mergeDim << nl;
927 
928  collapseEdge.set(edgeI);
929  collapsePointToLocation.set(e[1], points[e[0]]);
930  }
931  }
932 
933  List<pointEdgeCollapse> allPointInfo;
935  labelList pointPriority(mesh.nPoints(), Zero);
936 
937  collapser.consistentCollapse
938  (
939  globalPoints,
940  pointPriority,
941  collapsePointToLocation,
942  collapseEdge,
943  allPointInfo
944  );
945 
946  // Topo change container
947  polyTopoChange meshMod(mesh);
948 
949  // Put all modifications into meshMod
950  bool anyChange = collapser.setRefinement(allPointInfo, meshMod);
951 
952  if (returnReduceOr(anyChange))
953  {
954  // Construct new mesh from polyTopoChange.
955  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
956 
957  // Update fields
958  mesh.updateMesh(map());
959 
960  // Update stored data
961  updateFaceLabels(map(), frontPatchFaces);
962  updateFaceLabels(map(), backPatchFaces);
963  updateCellSet(map(), addedCellsSet);
964 
965  if (refDataPtr)
966  {
967  refDataPtr->updateMesh(map());
968  }
969 
970  // Move mesh (if inflation used)
971  if (map().hasMotionPoints())
972  {
973  mesh.movePoints(map().preMotionPoints());
974  }
975  }
976  }
977 
978 
979  // Change the front and back patch types as required
980  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
981 
982  word frontBackType(word::null);
983 
984  if (isType<extrudeModels::wedge>(model()))
985  {
986  changeFrontBackPatches<wedgePolyPatch>
987  (
988  mesh,
989  frontPatchName,
990  backPatchName
991  );
992  }
993  else if (isType<extrudeModels::plane>(model()))
994  {
995  changeFrontBackPatches<emptyPolyPatch>
996  (
997  mesh,
998  frontPatchName,
999  backPatchName
1000  );
1001  }
1002 
1003 
1004  // Merging front and back patch faces
1005  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1006 
1007  if (dict.get<bool>("mergeFaces"))
1008  {
1009  if (mode == MESH)
1010  {
1012  << "Cannot stitch front and back of extrusion since"
1013  << " in 'mesh' mode (extrusion appended to mesh)."
1014  << exit(FatalError);
1015  }
1016 
1017  Info<< "Assuming full 360 degree axisymmetric case;"
1018  << " stitching faces on patches "
1019  << frontPatchName << " and "
1020  << backPatchName << " together ..." << nl << endl;
1021 
1022  if (frontPatchFaces.size() != backPatchFaces.size())
1023  {
1025  << "Differing number of faces on front ("
1026  << frontPatchFaces.size() << ") and back ("
1027  << backPatchFaces.size() << ")"
1028  << exit(FatalError);
1029  }
1030 
1031 
1032  polyTopoChanger stitcher(mesh);
1033  stitcher.setSize(1);
1034 
1035  const word cutZoneName("originalCutFaceZone");
1036 
1037  List<faceZone*> fz
1038  (
1039  1,
1040  new faceZone
1041  (
1042  cutZoneName,
1043  frontPatchFaces,
1044  false, // none are flipped
1045  0,
1046  mesh.faceZones()
1047  )
1048  );
1049 
1051 
1052  // Add the perfect interface mesh modifier
1053  perfectInterface perfectStitcher
1054  (
1055  "couple",
1056  0,
1057  stitcher,
1058  cutZoneName,
1059  word::null, // dummy patch name
1060  word::null // dummy patch name
1061  );
1062 
1063  // Topo change container
1064  polyTopoChange meshMod(mesh);
1065 
1066  perfectStitcher.setRefinement
1067  (
1069  (
1071  (
1072  mesh.faces(),
1073  frontPatchFaces
1074  ),
1075  mesh.points()
1076  ),
1078  (
1080  (
1081  mesh.faces(),
1082  backPatchFaces
1083  ),
1084  mesh.points()
1085  ),
1086  meshMod
1087  );
1088 
1089  // Construct new mesh from polyTopoChange.
1090  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
1091 
1092  // Update fields
1093  mesh.updateMesh(map());
1094 
1095  // Update local data
1096  updateCellSet(map(), addedCellsSet);
1097 
1098  if (refDataPtr)
1099  {
1100  refDataPtr->updateMesh(map());
1101  }
1102 
1103  // Move mesh (if inflation used)
1104  if (map().hasMotionPoints())
1105  {
1106  mesh.movePoints(map().preMotionPoints());
1107  }
1108  }
1109 
1110  mesh.setInstance(runTimeExtruded.constant());
1111  Info<< "Writing mesh to " << mesh.objectPath() << nl << endl;
1112 
1113  if (!mesh.write())
1114  {
1116  << exit(FatalError);
1117  }
1118  // Remove any left-over files
1120 
1121  // Need writing cellSet
1122  if (returnReduceOr(addedCellsSet.size()))
1123  {
1124  cellSet addedCells(mesh, "addedCells", addedCellsSet);
1125  Info<< "Writing added cells to cellSet " << addedCells.name()
1126  << nl << endl;
1127  if (!addedCells.write())
1128  {
1130  << addedCells.name()
1131  << exit(FatalError);
1132  }
1133  }
1134 
1135  if (refDataPtr)
1136  {
1137  refDataPtr->write();
1138  }
1139 
1140 
1141  Info<< "End\n" << endl;
1142 
1143  return 0;
1144 }
1145 
1146 
1147 // ************************************************************************* //
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.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
dictionary dict
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:453
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
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:655
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:71
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:853
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:439
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:578
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1110
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:841
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1072
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:49
void setSize(const label n, unsigned int val=0u)
Alias for resize()
Definition: PackedList.H:765
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:402
T & first()
Access first element of the list, position [0].
Definition: UList.H:798
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:487
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:888
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:639
Does polyTopoChanges to remove edges. Can remove faces due to edge collapse but can not remove cells ...
Definition: edgeCollapser.H:64
Required Variables.
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
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:45
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)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:688
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:227
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:239
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:1066
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
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...
Foam::word regionName(Foam::polyMesh::defaultRegion)
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
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:289
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:964
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
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:513
An edge is a list of two point labels. The functionality it provides supports the discretisation 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)
Definition: labelList.C:31
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 Foam::string.
Definition: word.H:63
Reading required, file watched for runTime modification.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:397
wordList names() const
Return a list of patch names.
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:99
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
Definition: polyMesh.C:997
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1104
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:376
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1091
String literal.
Definition: keyType.H:82
virtual autoPtr< polyPatch > clone(const labelList &faceCells) const
Construct and return a clone, setting faceCells.
Definition: polyPatch.H:284
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:275
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Create additional fvSchemes/fvSolution files.
Definition: fvMeshTools.C:1057
void updateMesh(const mapPolyMesh &)
Update local numbering for changed mesh.
Definition: hexRef8Data.C:286
label nEdges() const
Number of mesh edges.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:32
const word & name() const noexcept
The patch name.
bool write() const
Write.
Definition: hexRef8Data.C:377
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:646
void addFvPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:625
static autoPtr< extrudeModel > New(const dictionary &dict)
Select null constructed.
T & last()
Access last element of the list, position [size()-1].
Definition: UList.H:812
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:276
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
List< word > wordList
A List of words.
Definition: fileName.H:58
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:98
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:592
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...
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
A collection of cell labels.
Definition: cellSet.H:47
Direct mesh changes based on v1.3 polyTopoChange syntax.
label whichPatch(const label faceIndex) const
Return patch index for a given mesh face index.
void transfer(HashTable< T, Key, Hash > &rhs)
Transfer contents into this table.
Definition: HashTable.C:712
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:441
const polyBoundaryMesh & patches
Automatically write from objectRegistry::writeObject()
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:185
messageStream Info
Information stream (stdout output on master, null elsewhere)
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:73
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
bool readIfPresent(const word &optName, T &val) const
Read a value from the named option if present.
Definition: argListI.H:316
label index() const noexcept
The index of this patch in the boundaryMesh.
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:959
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:166
mode_t mode(const fileName &name, const bool followLink=true)
Return the file mode, normally following symbolic links.
Definition: POSIX.C:726
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:209
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:617
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:231
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157