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