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