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