Loading...
Searching...
No Matches
snappyRefineMesh.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) 2018-2024 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 snappyRefineMesh
29
30Group
31 grpMeshAdvancedUtilities
32
33Description
34 Refine cells near to a surface.
35
36\*---------------------------------------------------------------------------*/
37
38#include "argList.H"
39#include "Time.H"
40#include "polyMesh.H"
41#include "twoDPointCorrector.H"
42#include "OFstream.H"
43#include "multiDirRefinement.H"
44
45#include "triSurface.H"
46#include "triSurfaceSearch.H"
47
48#include "cellSet.H"
49#include "pointSet.H"
50#include "surfaceToCell.H"
51#include "surfaceToPoint.H"
52#include "cellToPoint.H"
53#include "pointToCell.H"
54#include "cellToCell.H"
55#include "surfaceSets.H"
56#include "polyTopoChange.H"
57#include "polyTopoChanger.H"
58#include "mapPolyMesh.H"
59#include "labelIOList.H"
60#include "emptyPolyPatch.H"
61#include "removeCells.H"
62#include "meshSearch.H"
63
64using namespace Foam;
65
66// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67
68
69// Max cos angle for edges to be considered aligned with axis.
70static const scalar edgeTol = 1e-3;
71
72
73void writeSet(const cellSet& cells, const string& msg)
74{
75 Info<< "Writing " << msg << " (" << cells.size() << ") to cellSet "
76 << cells.instance()/cells.local()/cells.name()
77 << endl;
78 cells.write();
79}
80
81
82direction getNormalDir(const twoDPointCorrector* correct2DPtr)
83{
84 if (correct2DPtr)
85 {
86 const vector& normal = correct2DPtr->planeNormal();
87
88 if (mag(normal.x()) > 1-edgeTol)
89 {
90 return vector::X;
91 }
92 else if (mag(normal.y()) > 1-edgeTol)
93 {
94 return vector::Y;
95 }
96 else if (mag(normal.z()) > 1-edgeTol)
97 {
98 return vector::Z;
99 }
100 }
101
102 return direction(3);
103}
104
105
106
107// Calculate some edge statistics on mesh. Return min. edge length over all
108// directions but exclude component (0=x, 1=y, 2=z, other=none)
109scalar getEdgeStats(const primitiveMesh& mesh, const direction excludeCmpt)
110{
111 label nAny(0);
112 label nX(0);
113 label nY(0);
114 label nZ(0);
115
116 scalarMinMax limitsAny(GREAT, -GREAT);
117 scalarMinMax limitsX(limitsAny);
118 scalarMinMax limitsY(limitsAny);
119 scalarMinMax limitsZ(limitsAny);
120
121 const edgeList& edges = mesh.edges();
122
123 for (const edge& e : edges)
124 {
125 vector eVec(e.vec(mesh.points()));
126
127 scalar eMag = mag(eVec);
128
129 eVec /= eMag;
130
131 if (mag(eVec.x()) > 1-edgeTol)
132 {
133 limitsX.add(eMag);
134 nX++;
135 }
136 else if (mag(eVec.y()) > 1-edgeTol)
137 {
138 limitsY.add(eMag);
139 nY++;
140 }
141 else if (mag(eVec.z()) > 1-edgeTol)
142 {
143 limitsZ.add(eMag);
144 nZ++;
145 }
146 else
147 {
148 limitsAny.add(eMag);
149 nAny++;
150 }
151 }
152
153 Info<< "Mesh bounding box:" << boundBox(mesh.points()) << nl << nl
154 << "Mesh edge statistics:" << nl
155 << " x aligned : number:" << nX
156 << "\tminLen:" << limitsX.min() << "\tmaxLen:" << limitsX.max() << nl
157 << " y aligned : number:" << nY
158 << "\tminLen:" << limitsY.min() << "\tmaxLen:" << limitsY.max() << nl
159 << " z aligned : number:" << nZ
160 << "\tminLen:" << limitsZ.min() << "\tmaxLen:" << limitsZ.max() << nl
161 << " other : number:" << nAny
162 << "\tminLen:" << limitsAny.min()
163 << "\tmaxLen:" << limitsAny.max() << nl << endl;
164
165 if (excludeCmpt == vector::X)
166 {
167 return Foam::min
168 (
169 limitsAny.min(),
170 Foam::min(limitsY.min(), limitsZ.min())
171 );
172 }
173 else if (excludeCmpt == vector::Y)
174 {
175 return Foam::min
176 (
177 limitsAny.min(),
178 Foam::min(limitsX.min(), limitsZ.min())
179 );
180 }
181 else if (excludeCmpt == vector::Z)
182 {
183 return Foam::min
184 (
185 limitsAny.min(),
186 Foam::min(limitsX.min(), limitsY.min())
187 );
188 }
189 else
190 {
191 return Foam::min
192 (
193 limitsAny.min(),
195 (
196 limitsX.min(),
197 Foam::min(limitsY.min(), limitsZ.min())
198 )
199 );
200 }
201}
202
203
204// Adds empty patch if not yet there. Returns patchID.
205label addPatch(polyMesh& mesh, const word& patchName)
206{
207 label patchi = mesh.boundaryMesh().findPatchID(patchName);
208
209 if (patchi == -1)
210 {
211 const polyBoundaryMesh& patches = mesh.boundaryMesh();
212
213 polyPatchList newPatches(patches.size() + 1);
214
215 label nPatches = 0;
216
217 // Add empty patch as 0th entry (Note: only since subsetMesh wants this)
218 patchi = 0;
219 {
220 newPatches.set
221 (
222 nPatches,
224 (
225 patchName,
226 0,
227 mesh.nInternalFaces(),
228 nPatches,
229 patches,
230 emptyPolyPatch::typeName
231 )
232 );
233 ++nPatches;
234 }
235
236 for (const polyPatch& pp : patches)
237 {
238 newPatches.set
239 (
240 nPatches,
241 pp.clone
242 (
243 patches,
244 nPatches,
245 pp.size(),
246 pp.start()
247 )
248 );
249 ++nPatches;
250 }
251
252 mesh.removeBoundary();
253 mesh.addPatches(newPatches);
254
255 Info<< "Created patch oldInternalFaces at " << patchi << endl;
256 }
257 else
258 {
259 Info<< "Reusing patch oldInternalFaces at " << patchi << endl;
260 }
261
262
263 return patchi;
264}
265
266
267// Take surface and select cells based on surface curvature.
268void selectCurvatureCells
269(
270 const polyMesh& mesh,
271 const fileName& surfName,
272 const triSurfaceSearch& querySurf,
273 const scalar nearDist,
274 const scalar curvature,
276)
277{
278 // Use surfaceToCell to do actual calculation.
279
280 // Since we're adding make sure set is on disk.
281 cells.write();
282
283 // Take centre of cell 0 as outside point since info not used.
284
285 surfaceToCell cutSource
286 (
287 mesh,
288 surfName,
289 querySurf.surface(),
290 querySurf,
291 pointField(1, mesh.cellCentres()[0]),
292 false, // includeCut
293 false, // includeInside
294 false, // includeOutside
295 false, // geometricOnly
296 nearDist,
297 curvature
298 );
299
300 cutSource.applyToSet(topoSetSource::ADD, cells);
301}
302
303
304// cutCells contains currently selected cells to be refined. Add neighbours
305// on the inside or outside to them.
306void addCutNeighbours
307(
308 const polyMesh& mesh,
309 const bool selectInside,
310 const bool selectOutside,
311 const labelHashSet& inside,
312 const labelHashSet& outside,
313 labelHashSet& cutCells
314)
315{
316 // Pick up face neighbours of cutCells
317
318 labelHashSet addCutFaces(cutCells.size());
319
320 for (const label celli : cutCells)
321 {
322 const labelList& cFaces = mesh.cells()[celli];
323
324 forAll(cFaces, i)
325 {
326 const label facei = cFaces[i];
327
328 if (mesh.isInternalFace(facei))
329 {
330 label nbr = mesh.faceOwner()[facei];
331
332 if (nbr == celli)
333 {
334 nbr = mesh.faceNeighbour()[facei];
335 }
336
337 if (selectInside && inside.found(nbr))
338 {
339 addCutFaces.insert(nbr);
340 }
341 else if (selectOutside && outside.found(nbr))
342 {
343 addCutFaces.insert(nbr);
344 }
345 }
346 }
347 }
348
349 Info<< " Selected an additional " << addCutFaces.size()
350 << " neighbours of cutCells to refine" << endl;
351
352 for (const label facei : addCutFaces)
353 {
354 cutCells.insert(facei);
355 }
356}
357
358
359// Return true if any cells had to be split to keep a difference between
360// neighbouring refinement levels < limitDiff.
361// Gets cells which will be refined (so increase the refinement level) and
362// updates it.
363bool limitRefinementLevel
364(
365 const primitiveMesh& mesh,
366 const label limitDiff,
367 const labelHashSet& excludeCells,
368 const labelList& refLevel,
369 labelHashSet& cutCells
370)
371{
372 // Do simple check on validity of refinement level.
373 forAll(refLevel, celli)
374 {
375 if (!excludeCells.found(celli))
376 {
377 const labelList& cCells = mesh.cellCells()[celli];
378
379 forAll(cCells, i)
380 {
381 label nbr = cCells[i];
382
383 if (!excludeCells.found(nbr))
384 {
385 if (refLevel[celli] - refLevel[nbr] >= limitDiff)
386 {
388 << "Level difference between neighbouring cells "
389 << celli << " and " << nbr
390 << " greater than or equal to " << limitDiff << endl
391 << "refLevels:" << refLevel[celli] << ' '
392 << refLevel[nbr] << abort(FatalError);
393 }
394 }
395 }
396 }
397 }
398
399
400 labelHashSet addCutCells(cutCells.size());
401
402 for (const label celli : cutCells)
403 {
404 // celli will be refined.
405 const labelList& cCells = mesh.cellCells()[celli];
406
407 forAll(cCells, i)
408 {
409 const label nbr = cCells[i];
410
411 if (!excludeCells.found(nbr) && !cutCells.found(nbr))
412 {
413 if (refLevel[celli] + 1 - refLevel[nbr] >= limitDiff)
414 {
415 addCutCells.insert(nbr);
416 }
417 }
418 }
419 }
420
421 if (addCutCells.size())
422 {
423 // Add cells to cutCells.
424
425 Info<< "Added an additional " << addCutCells.size() << " cells"
426 << " to satisfy 1:" << limitDiff << " refinement level"
427 << endl;
428
429 for (const label celli : addCutCells)
430 {
431 cutCells.insert(celli);
432 }
433 return true;
434 }
435
436 Info<< "Added no additional cells"
437 << " to satisfy 1:" << limitDiff << " refinement level"
438 << endl;
439
440 return false;
441}
442
443
444// Do all refinement (i.e. refCells) according to refineDict and update
445// refLevel afterwards for added cells
446void doRefinement
447(
448 polyMesh& mesh,
449 const dictionary& refineDict,
450 const labelHashSet& refCells,
451 labelList& refLevel
452)
453{
454 label oldCells = mesh.nCells();
455
456 // Multi-iteration, multi-direction topology modifier.
457 multiDirRefinement multiRef
458 (
459 mesh,
460 refCells.toc(),
461 refineDict
462 );
463
464 //
465 // Update refLevel for split cells
466 //
467
468 refLevel.setSize(mesh.nCells());
469
470 for (label celli = oldCells; celli < mesh.nCells(); celli++)
471 {
472 refLevel[celli] = 0;
473 }
474
475 const labelListList& addedCells = multiRef.addedCells();
476
477 forAll(addedCells, oldCelli)
478 {
479 const labelList& added = addedCells[oldCelli];
480
481 if (added.size())
482 {
483 // Give all cells resulting from split the refinement level
484 // of the master.
485 label masterLevel = ++refLevel[oldCelli];
486
487 forAll(added, i)
488 {
489 refLevel[added[i]] = masterLevel;
490 }
491 }
492 }
493}
494
495
496// Subset mesh and update refLevel and cutCells
497void subsetMesh
498(
499 polyMesh& mesh,
500 const label writeMesh,
501 const label patchi, // patchID for exposed faces
502 const labelHashSet& cellsToRemove,
503 cellSet& cutCells,
504 labelIOList& refLevel
505)
506{
507 removeCells cellRemover(mesh);
508
509 labelList cellLabels(cellsToRemove.toc());
510
511 Info<< "Mesh has:" << mesh.nCells() << " cells."
512 << " Removing:" << cellLabels.size() << " cells" << endl;
513
514 // exposed faces
515 labelList exposedFaces(cellRemover.getExposedFaces(cellLabels));
516
517 polyTopoChange meshMod(mesh);
518 cellRemover.setRefinement
519 (
520 cellLabels,
521 exposedFaces,
522 labelList(exposedFaces.size(), patchi),
523 meshMod
524 );
525
526 // Do all changes
527 Info<< "Morphing ..." << endl;
528
529 const Time& runTime = mesh.time();
530
531 autoPtr<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
532
533 if (morphMap().hasMotionPoints())
534 {
535 mesh.movePoints(morphMap().preMotionPoints());
536 }
537
538 // Update topology on cellRemover
539 cellRemover.updateMesh(morphMap());
540
541 // Update refLevel for removed cells.
542 const labelList& cellMap = morphMap().cellMap();
543
544 labelList newRefLevel(cellMap.size());
545
546 forAll(cellMap, i)
547 {
548 newRefLevel[i] = refLevel[cellMap[i]];
549 }
550
551 // Transfer back to refLevel
552 refLevel.transfer(newRefLevel);
553
554 if (writeMesh)
555 {
556 Info<< "Writing refined mesh to time " << runTime.timeName() << nl
557 << endl;
558
559 // More precision (for points data)
560 const auto oldPrec = IOstream::minPrecision(10);
561
562 mesh.write();
563 refLevel.write();
564
566 }
567
568 // Update cutCells for removed cells.
569 cutCells.updateMesh(morphMap());
570}
571
572
573// Divide the cells according to status compared to surface. Constructs sets:
574// - cutCells : all cells cut by surface
575// - inside : all cells inside surface
576// - outside : ,, outside ,,
577// and a combined set:
578// - selected : sum of above according to selectCut/Inside/Outside flags.
579void classifyCells
580(
581 const polyMesh& mesh,
582 const fileName& surfName,
583 const triSurfaceSearch& querySurf,
584 const pointField& outsidePts,
585
586 const bool selectCut,
587 const bool selectInside,
588 const bool selectOutside,
589
590 const label nCutLayers,
591
592 cellSet& inside,
593 cellSet& outside,
594 cellSet& cutCells,
595 cellSet& selected
596)
597{
598 // Cut faces with surface and classify cells
600 (
601 mesh,
602 surfName,
603 querySurf.surface(),
604 querySurf,
605 outsidePts,
606
607 nCutLayers,
608
609 inside,
610 outside,
611 cutCells
612 );
613
614 // Combine wanted parts into selected
615 if (selectCut)
616 {
617 selected.addSet(cutCells);
618 }
619 if (selectInside)
620 {
621 selected.addSet(inside);
622 }
623 if (selectOutside)
624 {
625 selected.addSet(outside);
626 }
627
628 Info<< "Determined cell status:" << endl
629 << " inside :" << inside.size() << endl
630 << " outside :" << outside.size() << endl
631 << " cutCells:" << cutCells.size() << endl
632 << " selected:" << selected.size() << endl
633 << endl;
634
635 writeSet(inside, "inside cells");
636 writeSet(outside, "outside cells");
637 writeSet(cutCells, "cut cells");
638 writeSet(selected, "selected cells");
639}
640
641
642
643int main(int argc, char *argv[])
644{
646 (
647 "Refine cells near to a surface"
648 );
650
651 #include "setRootCase.H"
652 #include "createTime.H"
653 #include "createPolyMesh.H"
654
655 // If necessary add oldInternalFaces patch
656 label newPatchi = addPatch(mesh, "oldInternalFaces");
657
658
659 //
660 // Read motionProperties dictionary
661 //
662
663 Info<< "Checking for motionProperties\n" << endl;
664
665 IOobject motionObj
666 (
667 "motionProperties",
668 runTime.constant(),
669 mesh,
672 );
673
674 // corrector for mesh motion
675 twoDPointCorrector* correct2DPtr = nullptr;
676
677 if (motionObj.typeHeaderOk<IOdictionary>(true))
678 {
679 Info<< "Reading " << runTime.constant() / "motionProperties"
680 << endl << endl;
681
682 IOdictionary motionProperties(motionObj);
683
684 if (motionProperties.get<bool>("twoDMotion"))
685 {
686 Info<< "Correcting for 2D motion" << endl << endl;
687 correct2DPtr = new twoDPointCorrector(mesh);
688 }
689 }
690
691 IOdictionary refineDict
692 (
694 (
695 "snappyRefineMeshDict",
696 runTime.system(),
697 mesh,
700 )
701 );
702
703 fileName surfName(refineDict.get<fileName>("surface"));
704 surfName.expand();
705 const label nCutLayers(refineDict.get<label>("nCutLayers"));
706 const label cellLimit(refineDict.get<label>("cellLimit"));
707 const bool selectCut(refineDict.get<bool>("selectCut"));
708 const bool selectInside(refineDict.get<bool>("selectInside"));
709 const bool selectOutside(refineDict.get<bool>("selectOutside"));
710 const bool selectHanging(refineDict.get<bool>("selectHanging"));
711
712 const scalar minEdgeLen(refineDict.get<scalar>("minEdgeLen"));
713 const scalar maxEdgeLen(refineDict.get<scalar>("maxEdgeLen"));
714 const scalar curvature(refineDict.get<scalar>("curvature"));
715 const scalar curvDist(refineDict.get<scalar>("curvatureDistance"));
716 pointField outsidePts(refineDict.lookup("outsidePoints"));
717 const label refinementLimit(refineDict.get<label>("splitLevel"));
718 const bool writeMesh(refineDict.get<bool>("writeMesh"));
719
720 Info<< "Cells to be used for meshing (0=false, 1=true):" << nl
721 << " cells cut by surface : " << selectCut << nl
722 << " cells inside of surface : " << selectInside << nl
723 << " cells outside of surface : " << selectOutside << nl
724 << " hanging cells : " << selectHanging << nl
725 << endl;
726
727
728 if (nCutLayers > 0 && selectInside)
729 {
731 << "Illogical settings : Both nCutLayers>0 and selectInside true."
732 << endl
733 << "This would mean that inside cells get removed but should be"
734 << " included in final mesh" << endl;
735 }
736
737 // Surface.
738 triSurface surf(surfName);
739
740 // Search engine on surface
741 triSurfaceSearch querySurf(surf);
742
743 // Search engine on mesh. No face decomposition since mesh unwarped.
745
746 // Check all 'outside' points
747 forAll(outsidePts, outsidei)
748 {
749 const point& outsidePoint = outsidePts[outsidei];
750
751 if (queryMesh.findCell(outsidePoint, -1, false) == -1)
752 {
754 << "outsidePoint " << outsidePoint
755 << " is not inside any cell"
756 << exit(FatalError);
757 }
758 }
759
760
761
762 // Current refinement level. Read if present.
763 labelIOList refLevel
764 (
766 (
767 "refinementLevel",
768 runTime.timeName(),
770 mesh,
773 ),
774 labelList(mesh.nCells(), Zero)
775 );
776
777 label maxLevel = max(refLevel);
778
779 if (maxLevel > 0)
780 {
781 Info<< "Read existing refinement level from file "
782 << refLevel.objectPath() << nl
783 << " min level : " << min(refLevel) << nl
784 << " max level : " << maxLevel << nl
785 << endl;
786 }
787 else
788 {
789 Info<< "Created zero refinement level in file "
790 << refLevel.objectPath() << nl
791 << endl;
792 }
793
794
795
796
797 // Print edge stats on original mesh. Leave out 2d-normal direction
798 direction normalDir(getNormalDir(correct2DPtr));
799 scalar meshMinEdgeLen = getEdgeStats(mesh, normalDir);
800
801 while (meshMinEdgeLen > minEdgeLen)
802 {
803 // Get inside/outside/cutCells cellSets.
804 cellSet inside(mesh, "inside", mesh.nCells()/10);
805 cellSet outside(mesh, "outside", mesh.nCells()/10);
806 cellSet cutCells(mesh, "cutCells", mesh.nCells()/10);
807 cellSet selected(mesh, "selected", mesh.nCells()/10);
808
809 classifyCells
810 (
811 mesh,
812 surfName,
813 querySurf,
814 outsidePts,
815
816 selectCut, // for determination of selected
817 selectInside, // for determination of selected
818 selectOutside, // for determination of selected
819
820 nCutLayers, // How many layers of cutCells to include
821
822 inside,
823 outside,
824 cutCells,
825 selected // not used but determined anyway.
826 );
827
828 Info<< " Selected " << cutCells.name() << " with "
829 << cutCells.size() << " cells" << endl;
830
831 if ((curvDist > 0) && (meshMinEdgeLen < maxEdgeLen))
832 {
833 // Done refining enough close to surface. Now switch to more
834 // intelligent refinement where only cutCells on surface curvature
835 // are refined.
836 cellSet curveCells(mesh, "curveCells", mesh.nCells()/10);
837
838 selectCurvatureCells
839 (
840 mesh,
841 surfName,
842 querySurf,
843 maxEdgeLen,
844 curvature,
845 curveCells
846 );
847
848 Info<< " Selected " << curveCells.name() << " with "
849 << curveCells.size() << " cells" << endl;
850
851 // Add neighbours to cutCells. This is if selectCut is not
852 // true and so outside and/or inside cells get exposed there is
853 // also refinement in them.
854 if (!selectCut)
855 {
856 addCutNeighbours
857 (
858 mesh,
859 selectInside,
860 selectOutside,
861 inside,
862 outside,
863 cutCells
864 );
865 }
866
867 // Subset cutCells to only curveCells
868 cutCells.subset(curveCells);
869
870 Info<< " Removed cells not on surface curvature. Selected "
871 << cutCells.size() << endl;
872 }
873
874
875 if (nCutLayers > 0)
876 {
877 // Subset mesh to remove inside cells altogether. Updates cutCells,
878 // refLevel.
879 subsetMesh(mesh, writeMesh, newPatchi, inside, cutCells, refLevel);
880 }
881
882
883 // Added cells from 2:1 refinement level restriction.
884 while
885 (
886 limitRefinementLevel
887 (
888 mesh,
889 refinementLimit,
890 labelHashSet(),
891 refLevel,
892 cutCells
893 )
894 )
895 {}
896
897
898 Info<< " Current cells : " << mesh.nCells() << nl
899 << " Selected for refinement :" << cutCells.size() << nl
900 << endl;
901
902 if (cutCells.empty())
903 {
904 Info<< "Stopping refining since 0 cells selected to be refined ..."
905 << nl << endl;
906 break;
907 }
908
909 if ((mesh.nCells() + 8*cutCells.size()) > cellLimit)
910 {
911 Info<< "Stopping refining since cell limit reached ..." << nl
912 << "Would refine from " << mesh.nCells()
913 << " to " << mesh.nCells() + 8*cutCells.size() << " cells."
914 << nl << endl;
915 break;
916 }
917
918 doRefinement
919 (
920 mesh,
921 refineDict,
922 cutCells,
923 refLevel
924 );
925
926 Info<< " After refinement:" << mesh.nCells() << nl
927 << endl;
928
929 if (writeMesh)
930 {
931 Info<< " Writing refined mesh to time " << runTime.timeName()
932 << nl << endl;
933
934 // More precision (for points data)
935 const auto oldPrec = IOstream::minPrecision(10);
936
937 mesh.write();
938 refLevel.write();
939
941 }
942
943 // Update mesh edge stats.
944 meshMinEdgeLen = getEdgeStats(mesh, normalDir);
945 }
946
947
948 if (selectHanging)
949 {
950 // Get inside/outside/cutCells cellSets.
951 cellSet inside(mesh, "inside", mesh.nCells()/10);
952 cellSet outside(mesh, "outside", mesh.nCells()/10);
953 cellSet cutCells(mesh, "cutCells", mesh.nCells()/10);
954 cellSet selected(mesh, "selected", mesh.nCells()/10);
955
956 classifyCells
957 (
958 mesh,
959 surfName,
960 querySurf,
961 outsidePts,
962
963 selectCut,
964 selectInside,
965 selectOutside,
966
967 nCutLayers,
968
969 inside,
970 outside,
971 cutCells,
972 selected
973 );
974
975
976 // Find any cells which have all their points on the outside of the
977 // selected set and refine them
979
980 Info<< "Detected " << hanging.size() << " hanging cells"
981 << " (cells with all points on"
982 << " outside of cellSet 'selected').\nThese would probably result"
983 << " in flattened cells when snapping the mesh to the surface"
984 << endl;
985
986 Info<< "Refining " << hanging.size() << " hanging cells" << nl
987 << endl;
988
989 // Added cells from 2:1 refinement level restriction.
990 while
991 (
992 limitRefinementLevel
993 (
994 mesh,
995 refinementLimit,
996 labelHashSet(),
997 refLevel,
998 hanging
999 )
1000 )
1001 {}
1002
1003 doRefinement(mesh, refineDict, hanging, refLevel);
1004
1005 Info<< "Writing refined mesh to time " << runTime.timeName() << nl
1006 << endl;
1007
1008 // Write final mesh
1009
1010 // More precision (for points data)
1011 const auto oldPrec = IOstream::minPrecision(10);
1012
1013 mesh.write();
1014 refLevel.write();
1015
1017 }
1018 else if (!writeMesh)
1019 {
1020 // Write final mesh. (will have been written already if writeMesh=true)
1021
1022 Info<< "Writing refined mesh to time " << runTime.timeName() << nl
1023 << endl;
1024
1025 // More precision (for points data)
1026 const auto oldPrec = IOstream::minPrecision(10);
1027
1028 mesh.write();
1029 refLevel.write();
1030
1032 }
1033
1034
1035 Info<< "End\n" << endl;
1036
1037 return 0;
1038}
1039
1040
1041// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition HashTable.C:141
bool empty() const noexcept
True if the hash table is empty.
Definition HashTable.H:353
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
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...
@ 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
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
fileName objectPath() const
The complete path + object name.
Definition IOobjectI.H:313
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
Definition IOstream.H:459
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition IOstream.H:437
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void setSize(label n)
Alias for resize().
Definition List.H:536
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static void noParallel()
Remove the parallel options.
Definition argList.C:599
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
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
A collection of cell labels.
Definition cellSet.H:50
virtual void updateMesh(const mapPolyMesh &morphMap)
Update any stored data for new labels.
Definition cellSet.C:208
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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.
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream.
Definition dictionary.C:367
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
Empty front and back plane patch. Used for 2-D geometries.
A class for handling file names.
Definition fileName.H:75
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition meshSearch.H:57
Does multiple pass refinement to refine cells in multiple directions.
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 word defaultRegion
Return the default region name.
Definition polyMesh.H:406
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.
Cell-face mesh analysis engine.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Given list of cells to remove, insert all the topology changes.
Definition removeCells.H:60
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
static labelHashSet getHangingCells(const primitiveMesh &mesh, const labelHashSet &internalCells)
Get cells using points on 'outside' only.
static void getSurfaceSets(const polyMesh &mesh, const fileName &surfName, const triSurface &surf, const triSurfaceSearch &querySurf, const pointField &outsidePts, const label nCutLayers, labelHashSet &inside, labelHashSet &outside, labelHashSet &cut)
Divide cells into cut,inside and outside.
A topoSetCellSource to select cells based on relation to a surface given by an external file.
@ ADD
Add elements to current set.
virtual void subset(const labelUList &elems)
Subset contents. Only elements present in both sets remain.
Definition topoSet.C:597
virtual void addSet(const labelUList &elems)
Add given elements to the set.
Definition topoSet.C:625
Helper class to search on triSurface.
const triSurface & surface() const
Return reference to the surface.
Triangulated surface description with patch information.
Definition triSurface.H:74
Class applies a two-dimensional correction to mesh motion point field.
const vector & planeNormal() const
Return plane normal.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const cellShapeList & cells
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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
IOList< label > labelIOList
IO for a List of label.
Definition labelIOList.H:32
messageStream Info
Information stream (stdout output on master, null elsewhere).
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label nPatches
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299