extrudeMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 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 or from patch read from file.
29 
30  Type of extrusion prescribed by run-time selectable model.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "argList.H"
35 #include "Time.H"
36 #include "polyTopoChange.H"
37 #include "edgeCollapser.H"
38 #include "perfectInterface.H"
39 #include "addPatchCellLayer.H"
40 #include "fvMesh.H"
41 #include "MeshedSurfaces.H"
42 #include "globalIndex.H"
43 #include "cellSet.H"
44 #include "systemDict.H"
45 
46 #include "extrudedMesh.H"
47 #include "extrudeModel.H"
48 
49 #include "wedge.H"
50 #include "wedgePolyPatch.H"
51 #include "planeExtrusion.H"
52 #include "emptyPolyPatch.H"
53 #include "processorPolyPatch.H"
54 
55 using namespace Foam;
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 enum class extrudeSurfaceType
60 {
61  mesh,
62  patch,
63  surface
64 };
65 
66 static const NamedEnum<extrudeSurfaceType, 3> extrudeSurfaceTypeNames
67 {
68  "mesh",
69  "patch",
70  "surface"
71 };
72 
73 
75 {
76  const label patchID = patches.findIndex(name);
77 
78  if (patchID == -1)
79  {
81  << "Cannot find patch " << name
82  << " in the source mesh.\n"
83  << "Valid patch names are " << patches.names()
84  << exit(FatalError);
85  }
86  return patchID;
87 }
88 
89 
90 labelList findPatchZones(polyMesh& mesh, const wordList& zoneNames)
91 {
92  labelList patchZones(zoneNames.size());
93 
94  // Add optional cellZone generation for the extrusion layer
95  forAll(zoneNames, i)
96  {
97  const label cellZonei = mesh.cellZones().findIndex(zoneNames[i]);
98  if (cellZonei >= 0)
99  {
100  patchZones[i] = cellZonei;
101  }
102  else
103  {
104  label nZones = mesh.cellZones().size();
105  mesh.cellZones().setSize(nZones + 1);
106  mesh.cellZones().set
107  (
108  nZones,
109  new cellZone
110  (
111  zoneNames[i],
112  labelList(),
113  mesh.cellZones()
114  )
115  );
116  patchZones[i] = nZones;
117  }
118  }
119 
120  return patchZones;
121 }
122 
123 
124 labelList findPatchFaces
125 (
126  polyMesh& mesh,
127  const polyBoundaryMesh& patches,
128  const wordList& patchNames,
129  const labelList& patchZones,
130  labelList& faceCellZones
131 )
132 {
133  // Count of the number of faces in the named patches
134  label n = 0;
135 
136  forAll(patchNames, i)
137  {
139  }
140 
141  labelList faceLabels(n);
142  faceCellZones.setSize(n);
143  n = 0;
144 
145  forAll(patchNames, i)
146  {
147  const polyPatch& pp = patches[findIndex(patches, patchNames[i])];
148 
149  forAll(pp, j)
150  {
151  faceLabels[n] = pp.start()+j;
152  faceCellZones[n] = patchZones[i];
153  n++;
154  }
155  }
156 
157  return faceLabels;
158 }
159 
160 
161 labelList findPatchFaces
162 (
163  const polyBoundaryMesh& patches,
164  const word& name
165 )
166 {
167  const polyPatch& pp = patches[findIndex(patches, name)];
168  labelList faceLabels(pp.size());
169 
170  forAll(pp, i)
171  {
172  faceLabels[i] = pp.start() + i;
173  }
174 
175  return faceLabels;
176 }
177 
178 
179 void updateFaceLabels(const polyTopoChangeMap& map, labelList& faceLabels)
180 {
181  const labelList& reverseMap = map.reverseFaceMap();
182 
183  labelList newFaceLabels(faceLabels.size());
184  label newI = 0;
185 
186  forAll(faceLabels, i)
187  {
188  label oldFacei = faceLabels[i];
189 
190  if (reverseMap[oldFacei] >= 0)
191  {
192  newFaceLabels[newI++] = reverseMap[oldFacei];
193  }
194  }
195  newFaceLabels.setSize(newI);
196  faceLabels.transfer(newFaceLabels);
197 }
198 
199 
200 void updateCellSet(const polyTopoChangeMap& map, labelHashSet& cellLabels)
201 {
202  const labelList& reverseMap = map.reverseCellMap();
203 
204  labelHashSet newCellLabels(2*cellLabels.size());
205 
206  forAll(cellLabels, i)
207  {
208  label oldCelli = cellLabels[i];
209 
210  if (reverseMap[oldCelli] >= 0)
211  {
212  newCellLabels.insert(reverseMap[oldCelli]);
213  }
214  }
215  cellLabels.transfer(newCellLabels);
216 }
217 
218 
219 template<class PatchType>
220 void changeFrontBackPatches
221 (
222  polyMesh& mesh,
223  const word& frontPatchName,
224  const word& backPatchName
225 )
226 {
228 
229  label frontPatchi = findIndex(patches, frontPatchName);
230  label backPatchi = findIndex(patches, backPatchName);
231 
232  DynamicList<polyPatch*> newPatches(patches.size());
233 
235  {
236  const polyPatch& pp(patches[patchi]);
237 
238  if (patchi == frontPatchi || patchi == backPatchi)
239  {
240  newPatches.append
241  (
242  new PatchType
243  (
244  pp.name(),
245  pp.size(),
246  pp.start(),
247  pp.index(),
248  mesh.boundaryMesh(),
249  PatchType::typeName
250  )
251  );
252  }
253  else
254  {
255  newPatches.append(pp.clone(mesh.boundaryMesh()).ptr());
256  }
257  }
258 
259  // Edit patches
261  mesh.addPatches(newPatches, true);
262 }
263 
264 
265 int main(int argc, char *argv[])
266 {
267  #include "addRegionOption.H"
268  #include "addMeshOption.H"
269  #include "addDictOption.H"
270 
271  #include "setRootCase.H"
272  #include "createTimeExtruded.H"
273 
274  const dictionary dict(systemDict("extrudeMeshDict", args, runTimeExtruded));
275 
276  // Point generator
278 
279  // Whether to flip normals
280  const Switch flipNormals(dict.lookup("flipNormals"));
281 
282  // What to extrude
283  const extrudeSurfaceType mode = extrudeSurfaceTypeNames.read
284  (
285  dict.lookup("constructFrom")
286  );
287 
288  // Any merging of small edges
289  const scalar mergeTol(dict.lookupOrDefault<scalar>("mergeTol", 1e-4));
290 
291  Info<< "Extruding from " << extrudeSurfaceTypeNames[mode]
292  << " using model " << model().type() << endl;
293  if (flipNormals)
294  {
295  Info<< "Flipping normals before extruding" << endl;
296  }
297  if (mergeTol > 0)
298  {
299  Info<< "Collapsing edges < " << mergeTol << " of bounding box" << endl;
300  }
301  else
302  {
303  Info<< "Not collapsing any edges after extrusion" << endl;
304  }
305  Info<< endl;
306 
307 
308  // Generated mesh (one of either)
309  autoPtr<fvMesh> meshFromMesh;
310  autoPtr<polyMesh> meshFromSurface;
311 
312  // Faces on front and back for stitching (in case of mergeFaces)
313  word frontPatchName;
314  labelList frontPatchFaces;
315  word backPatchName;
316  labelList backPatchFaces;
317 
318  // Optional added cells (get written to cellSet)
319  labelHashSet addedCellsSet;
320 
321  if (mode == extrudeSurfaceType::patch || mode == extrudeSurfaceType::mesh)
322  {
323  if (flipNormals && mode == extrudeSurfaceType::mesh)
324  {
326  << "Flipping normals not supported for extrusions from mesh."
327  << exit(FatalError);
328  }
329 
330  fileName sourceCasePath(dict.lookup("sourceCase"));
331  sourceCasePath.expand();
332  fileName sourceRootDir = sourceCasePath.path();
333  fileName sourceCaseDir = sourceCasePath.name();
334  if (Pstream::parRun())
335  {
336  sourceCaseDir =
337  sourceCaseDir
338  /"processor" + Foam::name(Pstream::myProcNo());
339  }
340  wordList sourcePatches(dict.lookup("sourcePatches"));
341 
342  if (sourcePatches.size() == 1)
343  {
344  frontPatchName = sourcePatches[0];
345  }
346 
347  wordList zoneNames(dict.lookupOrDefault("zoneNames", wordList::null()));
348  if (zoneNames.size() == 1)
349  {
350  if (zoneNames[0] == "patchNames")
351  {
352  zoneNames = sourcePatches;
353  }
354  else
355  {
356  zoneNames = wordList(sourcePatches.size(), zoneNames[0]);
357  }
358  }
359  else if (zoneNames.size() && zoneNames.size() != sourcePatches.size())
360  {
362  << "Number of zoneNames "
363  "does not equal the number of sourcePatches"
364  << exit(FatalError);
365  }
366 
367  Info<< "Extruding patches " << sourcePatches
368  << " on mesh " << sourceCasePath << nl
369  << endl;
370 
371  Time runTime
372  (
374  sourceRootDir,
375  sourceCaseDir
376  );
377 
379 
381 
382  const labelList patchZones
383  (
384  zoneNames.size()
385  ? findPatchZones(mesh, zoneNames)
386  : labelList(sourcePatches.size(), -1)
387  );
388 
389  labelList faceCellZones;
390  const labelList patchFaces
391  (
392  findPatchFaces
393  (
394  mesh,
395  patches,
396  sourcePatches,
397  patchZones,
398  faceCellZones
399  )
400  );
401 
402  // Extrusion engine. Either adding to existing mesh or
403  // creating separate mesh.
405 
406  if (mode == extrudeSurfaceType::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(patchFaces, i)
412  {
413  label patchFacei = patchFaces[i];
414 
415  label patchi = patches.whichPatch(patchFacei);
416  label own = mesh.faceOwner()[patchFacei];
417  label nei = -1;
418  if (patchi == -1)
419  {
420  nei = mesh.faceNeighbour()[patchFacei];
421  }
422 
423  meshMod.modifyFace
424  (
425  mesh.faces()[patchFacei].reverseFace(), // modified face
426  patchFacei, // label of face
427  own, // owner
428  nei, // neighbour
429  true, // face flip
430  patchi // patch for face
431  );
432  }
433 
434  // Change the mesh. without keeping old points.
435  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh);
436 
437  // Update fields
438  mesh.topoChange(map);
439  }
440 
441 
442  indirectPrimitivePatch extrudePatch
443  (
445  (
446  mesh.faces(),
447  patchFaces
448  ),
449  mesh.points()
450  );
451 
452  // Determine extrudePatch normal
453  pointField extrudePatchPointNormals
454  (
455  PatchTools::pointNormals(mesh, extrudePatch)
456  );
457 
458 
459  // Precalculate mesh edges for pp.edges.
460  const labelList meshEdges
461  (
462  extrudePatch.meshEdges
463  (
464  mesh.edges(),
465  mesh.pointEdges()
466  )
467  );
468 
469  // Global face indices engine
470  const globalIndex globalFaces(mesh.nFaces());
471 
472  // Determine extrudePatch.edgeFaces in global numbering (so across
473  // coupled patches)
474  labelListList edgeGlobalFaces
475  (
477  (
478  mesh,
479  globalFaces,
480  extrudePatch
481  )
482  );
483 
484 
485  // Determine what patches boundary edges need to get extruded into.
486  // This might actually cause edge-connected processors to become
487  // face-connected so might need to introduce new processor boundaries.
488  // Calculates:
489  // - per pp.edge the patch to extrude into
490  // - any additional processor boundaries (patchToNbrProc = map from
491  // new patchID to neighbour processor)
492  // - number of new patches (nPatches)
493 
494  labelList sidePatchID;
495  label nPatches;
496  Map<label> nbrProcToPatch;
497  Map<label> patchToNbrProc;
499  (
500  mesh,
501  globalFaces,
502  edgeGlobalFaces,
503  extrudePatch,
504 
505  sidePatchID,
506  nPatches,
507  nbrProcToPatch,
508  patchToNbrProc
509  );
510 
511 
512  // Add any patches
513 
514  label nAdded = nPatches - mesh.boundaryMesh().size();
515  reduce(nAdded, sumOp<label>());
516 
517  Info<< "Adding overall " << nAdded << " processor patches." << endl;
518 
519  if (nAdded > 0)
520  {
521  DynamicList<polyPatch*> newPatches(nPatches);
523  {
524  newPatches.append
525  (
527  (
529  ).ptr()
530  );
531  }
532  for
533  (
535  patchi < nPatches;
536  patchi++
537  )
538  {
539  label nbrProci = patchToNbrProc[patchi];
540 
541  Pout<< "Adding patch " << patchi
542  << " between " << Pstream::myProcNo()
543  << " and " << nbrProci
544  << endl;
545 
546  newPatches.append
547  (
549  (
550  0, // size
551  mesh.nFaces(), // start
552  patchi, // index
553  mesh.boundaryMesh(),// polyBoundaryMesh
554  Pstream::myProcNo(),// myProcNo
555  nbrProci // neighbProcNo
556  )
557  );
558  }
559 
560  // Add patches. Do no parallel updates.
562  mesh.addFvPatches(newPatches, true);
563  }
564 
565 
566 
567  // Only used for addPatchCellLayer into new mesh
568  labelList exposedPatchID;
569  if (mode == extrudeSurfaceType::patch)
570  {
571  dict.lookup("exposedPatchName") >> backPatchName;
572  exposedPatchID.setSize
573  (
574  extrudePatch.size(),
575  findIndex(patches, backPatchName)
576  );
577  }
578 
579  // Determine points and extrusion
580  pointField layer0Points(extrudePatch.nPoints());
581  pointField displacement(extrudePatch.nPoints());
582  forAll(displacement, pointi)
583  {
584  const vector& patchNormal = extrudePatchPointNormals[pointi];
585 
586  // layer0 point
587  layer0Points[pointi] = model()
588  (
589  extrudePatch.localPoints()[pointi],
590  patchNormal,
591  0
592  );
593  // layerN point
594  point extrudePt = model()
595  (
596  extrudePatch.localPoints()[pointi],
597  patchNormal,
598  model().nLayers()
599  );
600  displacement[pointi] = extrudePt - layer0Points[pointi];
601  }
602 
603 
604  // Check if wedge (has layer0 different from original patch points)
605  // If so move the mesh to starting position.
606  if (gMax(mag(layer0Points-extrudePatch.localPoints())) > small)
607  {
608  Info<< "Moving mesh to layer0 points since differ from original"
609  << " points - this can happen for wedge extrusions." << nl
610  << endl;
611 
612  pointField newPoints(mesh.points());
613  forAll(extrudePatch.meshPoints(), i)
614  {
615  newPoints[extrudePatch.meshPoints()[i]] = layer0Points[i];
616  }
617  mesh.setPoints(newPoints);
618  }
619 
620 
621  // Layers per face
622  labelList nFaceLayers(extrudePatch.size(), model().nLayers());
623 
624  // Layers per point
625  labelList nPointLayers(extrudePatch.nPoints(), model().nLayers());
626 
627  // Displacement for first layer
628  vectorField firstLayerDisp(displacement*model().sumThickness(1));
629 
630  // Expansion ratio not used.
631  scalarField ratio(extrudePatch.nPoints(), 1.0);
632 
633  // Topo change container. Either copy an existing mesh or start
634  // with empty storage (number of patches only needed for checking)
636  (
637  (
639  ? new polyTopoChange(mesh)
640  : new polyTopoChange(patches.size())
641  )
642  );
643 
644  layerExtrude.setRefinement
645  (
646  mesh,
647  globalFaces,
648  edgeGlobalFaces,
649 
650  ratio, // expansion ratio
651  extrudePatch, // patch faces to extrude
652  sidePatchID, // if boundary edge: patch to use
653  exposedPatchID, // if new mesh: patches for exposed faces
654  nFaceLayers,
655  nPointLayers,
656  firstLayerDisp,
657  faceCellZones,
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 actual mesh from polyTopoChange container
699  autoPtr<polyTopoChangeMap> map = meshMod().makeMesh
700  (
701  meshFromMesh,
702  IOobject
703  (
704  regionName,
705  runTimeExtruded.constant(),
706  runTimeExtruded,
709  false
710  ),
711  mesh
712  );
713 
714  layerExtrude.updateZones(meshFromMesh());
715 
716  meshFromMesh().topoChange(map);
717 
718  layerExtrude.topoChange
719  (
720  map(),
721  identityMap(extrudePatch.size()),
722  identityMap(extrudePatch.nPoints())
723  );
724 
725  // Calculate face labels for front and back.
726  frontPatchFaces = renumber
727  (
728  map().reverseFaceMap(),
729  frontPatchFaces
730  );
731  backPatchFaces = renumber
732  (
733  map().reverseFaceMap(),
734  backPatchFaces
735  );
736 
737  // Store added cells
739  {
740  const labelListList addedCells
741  (
742  layerExtrude.addedCells
743  (
744  meshFromMesh,
745  layerExtrude.layerFaces()
746  )
747  );
748  forAll(addedCells, facei)
749  {
750  const labelList& aCells = addedCells[facei];
751  forAll(aCells, i)
752  {
753  addedCellsSet.insert(aCells[i]);
754  }
755  }
756  }
757  }
758  else
759  {
760  // Read from surface
761  fileName surfName(dict.lookup("surface"));
762  surfName.expand();
763 
764  Info<< "Extruding surfaceMesh read from file " << surfName << nl
765  << endl;
766 
767  MeshedSurface<face> fMesh(surfName);
768 
769  if (flipNormals)
770  {
771  Info<< "Flipping faces." << nl << endl;
772  faceList& faces = const_cast<faceList&>(fMesh.faces());
773  forAll(faces, i)
774  {
775  faces[i] = fMesh[i].reverseFace();
776  }
777  }
778 
779  Info<< "Extruding surface with :" << nl
780  << " points : " << fMesh.points().size() << nl
781  << " faces : " << fMesh.size() << nl
782  << " normals[0] : " << fMesh.faceNormals()[0]
783  << nl
784  << endl;
785 
786  meshFromSurface.reset
787  (
788  new extrudedMesh
789  (
790  IOobject
791  (
793  runTimeExtruded.constant(),
794  runTimeExtruded
795  ),
796  fMesh,
797  model()
798  )
799  );
800 
801 
802  // Get the faces on front and back
803  frontPatchName = "originalPatch";
804  frontPatchFaces = findPatchFaces
805  (
806  meshFromSurface().boundaryMesh(),
807  frontPatchName
808  );
809  backPatchName = "otherSide";
810  backPatchFaces = findPatchFaces
811  (
812  meshFromSurface().boundaryMesh(),
813  backPatchName
814  );
815  }
816 
817 
818  polyMesh& mesh =
819  (
820  meshFromMesh.valid()
821  ? meshFromMesh()
822  : meshFromSurface()
823  );
824 
825 
826  const boundBox& bb = mesh.bounds();
827  const vector span = bb.span();
828  const scalar mergeDim = mergeTol * bb.minDim();
829 
830  Info<< "Mesh bounding box : " << bb << nl
831  << " with span : " << span << nl
832  << "Merge distance : " << mergeDim << nl
833  << endl;
834 
835 
836  // Collapse edges
837  // ~~~~~~~~~~~~~~
838 
839  if (mergeDim > 0)
840  {
841  Pout<< "Collapsing edges < " << mergeDim << " ..." << nl << endl;
842 
843  // Edge collapsing engine
844  edgeCollapser collapser(mesh);
845 
846  const edgeList& edges = mesh.edges();
847  const pointField& points = mesh.points();
848 
850  Map<point> collapsePointToLocation(mesh.nPoints());
851 
852  forAll(edges, edgeI)
853  {
854  const edge& e = edges[edgeI];
855 
856  scalar d = e.mag(points);
857 
858  if (d < mergeDim)
859  {
860  Pout<< "Merging edge " << e << " since length " << d
861  << " << " << mergeDim << nl;
862 
863  collapseEdge[edgeI] = true;
864  collapsePointToLocation.set(e[1], points[e[0]]);
865  }
866  }
867 
868  List<pointEdgeCollapse> allPointInfo;
870  labelList pointPriority(mesh.nPoints(), 0);
871 
872  collapser.consistentCollapse
873  (
874  globalPoints,
875  pointPriority,
876  collapsePointToLocation,
877  collapseEdge,
878  allPointInfo
879  );
880 
881  // Topo change container
882  polyTopoChange meshMod(mesh);
883 
884  // Put all modifications into meshMod
885  bool anyChange = collapser.setRefinement(allPointInfo, meshMod);
886  reduce(anyChange, orOp<bool>());
887 
888  if (anyChange)
889  {
890  // Construct new mesh from polyTopoChange.
891  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh);
892 
893  // Update fields
894  mesh.topoChange(map);
895 
896  // Update stored data
897  updateFaceLabels(map(), frontPatchFaces);
898  updateFaceLabels(map(), backPatchFaces);
899  updateCellSet(map(), addedCellsSet);
900  }
901  }
902 
903 
904  // Change the front and back patch types as required
905  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
906 
907  word frontBackType(word::null);
908 
909  if (isType<extrudeModels::wedge>(model()))
910  {
911  changeFrontBackPatches<wedgePolyPatch>
912  (
913  mesh,
914  frontPatchName,
915  backPatchName
916  );
917  }
918  else if (isType<extrudeModels::plane>(model()))
919  {
920  changeFrontBackPatches<emptyPolyPatch>
921  (
922  mesh,
923  frontPatchName,
924  backPatchName
925  );
926  }
927 
928 
929  // Merging front and back patch faces
930  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
931 
932  Switch mergeFaces(dict.lookup("mergeFaces"));
933  if (mergeFaces)
934  {
936  {
938  << "Cannot stitch front and back of extrusion since"
939  << " in 'mesh' mode (extrusion appended to mesh)."
940  << exit(FatalError);
941  }
942 
943  Info<< "Assuming full 360 degree axisymmetric case;"
944  << " stitching faces on patches "
945  << frontPatchName << " and "
946  << backPatchName << " together ..." << nl << endl;
947 
948  if (frontPatchFaces.size() != backPatchFaces.size())
949  {
951  << "Differing number of faces on front ("
952  << frontPatchFaces.size() << ") and back ("
953  << backPatchFaces.size() << ")"
954  << exit(FatalError);
955  }
956 
957  // Add the perfect interface mesh modifier
958  perfectInterface perfectStitcher
959  (
960  "couple",
961  mesh,
962  word::null, // dummy patch name
963  word::null // dummy patch name
964  );
965 
966  // Topo change container
967  polyTopoChange meshMod(mesh);
968 
969  perfectStitcher.setRefinement
970  (
972  (
974  (
975  mesh.faces(),
976  frontPatchFaces
977  ),
978  mesh.points()
979  ),
981  (
983  (
984  mesh.faces(),
985  backPatchFaces
986  ),
987  mesh.points()
988  ),
989  meshMod
990  );
991 
992  // Construct new mesh from polyTopoChange.
993  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh);
994 
995  // Update fields
996  mesh.topoChange(map);
997 
998  // Update local data
999  updateCellSet(map(), addedCellsSet);
1000  }
1001 
1002  mesh.setInstance(runTimeExtruded.constant());
1003  Info<< "Writing mesh to " << mesh.relativeObjectPath() << nl << endl;
1004 
1005  if (!mesh.write())
1006  {
1008  << exit(FatalError);
1009  }
1010 
1011  // Need writing cellSet
1012  label nAdded = returnReduce(addedCellsSet.size(), sumOp<label>());
1013  if (nAdded > 0)
1014  {
1015  cellSet addedCells(mesh, "addedCells", addedCellsSet);
1016  Info<< "Writing added cells to cellSet " << addedCells.name()
1017  << nl << endl;
1018  if (!addedCells.write())
1019  {
1021  << addedCells.name()
1022  << exit(FatalError);
1023  }
1024  }
1025 
1026  Info<< "End\n" << endl;
1027 
1028  return 0;
1029 }
1030 
1031 
1032 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
void transfer(HashTable< T, Key, Hash > &)
Transfer the contents of the argument table into this table.
Definition: HashTable.C:587
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
fileName relativeObjectPath() const
Return complete relativePath + object name.
Definition: IOobject.H:420
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:280
A List with indirect addressing.
Definition: IndirectList.H:105
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static const List< word > & null()
Return a null List.
Definition: ListI.H:118
void setSize(const label)
Reset size of List.
Definition: List.C:281
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
A bit-packed bool list.
static tmp< pointField > pointNormals(const polyMesh &, const PrimitivePatch< FaceList, PointField > &)
Return parallel consistent point normals for patches using mesh points.
autoPtr< T > set(const label, const word &key, T *)
Set element to pointer provided and return old element.
label findIndex(const word &key) const
Return the index of the given the key or -1 if not found.
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:208
T & first()
Return the first element of the list.
Definition: UListI.H:114
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
T & last()
Return the last element of the list.
Definition: UListI.H:128
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
Adds layers of cells to outside of polyPatch. Can optionally create stand-alone extruded mesh (addToM...
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.
static labelListList globalEdgeFaces(const polyMesh &, const globalIndex &globalFaces, const indirectPrimitivePatch &pp)
Per patch edge the pp faces (in global indices) using it. Uses.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
scalar minDim() const
Smallest length/height/width dimension.
Definition: boundBoxI.H:108
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:90
A collection of cell labels.
Definition: cellSet.H:51
Named list of cell indices representing a sub-set of the mesh.
Definition: cellZone.H:61
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
Does polyTopoChanges to remove edges. Can remove faces due to edge collapse but can not remove cells ...
Definition: edgeCollapser.H:67
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
static autoPtr< extrudeModel > New(const dictionary &)
Select null constructed.
A class for handling file names.
Definition: fileName.H:82
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:284
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:785
virtual bool write(const bool write=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1785
void removeFvBoundary()
Remove boundary patches. Warning: fvPatchFields hold ref to.
Definition: fvMesh.C:803
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1369
virtual void setPoints(const pointField &)
Reset the points.
Definition: fvMesh.C:1228
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
Calculates points shared by more than two processor patches or cyclic patches.
Definition: globalPoints.H:101
label index() const
Return the index of this patch in the boundaryMesh.
const word & name() const
Return name.
Calculates zone number for every face of patch.
Definition: patchZones.H:56
Hack of attachDetach to couple patches when they perfectly align. Does not decouple....
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:449
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:270
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1363
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:1107
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:91
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:407
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:215
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reverseCellMap() const
Reverse cell map.
const labelList & reverseFaceMap() const
Reverse face map.
Direct mesh changes based on v1.3 polyTopoChange syntax.
const labelListList & pointEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
label nPoints() const
label nEdges() const
label nFaces() const
Neighbour processor patch.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
label collapseEdge(triSurface &surf, const scalar minLen)
Keep collapsing all edges < minLen.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
const pointField & points
const fvPatchList & patches
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
const doubleScalar e
Definition: doubleScalar.H:106
mode_t mode(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file mode.
Definition: POSIX.C:461
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const word & regionName(const solver &region)
Definition: solver.H:218
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
messageStream Info
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion, const fileName &path=fileName::null)
Definition: systemDict.C:93
static const char nl
Definition: Ostream.H:267
Type gMax(const FieldField< Field, Type > &f)
wordList patchNames(nPatches)
label nPatches
Definition: readKivaGrid.H:396
dictionary dict
Foam::argList args(argc, argv)