autoRefineMesh.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 | Copyright (C) 2011-2016 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 "polyTopoChanger.H"
52 #include "mapPolyMesh.H"
53 #include "labelIOList.H"
54 #include "emptyPolyPatch.H"
55 #include "removeCells.H"
56 #include "meshSearch.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().findPatchID(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<mapPolyMesh> morphMap = meshMod.changeMesh(mesh, false);
516 
517  if (morphMap().hasMotionPoints())
518  {
519  mesh.movePoints(morphMap().preMotionPoints());
520  }
521 
522  // Update topology on cellRemover
523  cellRemover.updateMesh(morphMap());
524 
525  // Update refLevel for removed cells.
526  const labelList& cellMap = morphMap().cellMap();
527 
528  labelList newRefLevel(cellMap.size());
529 
530  forAll(cellMap, i)
531  {
532  newRefLevel[i] = refLevel[cellMap[i]];
533  }
534 
535  // Transfer back to refLevel
536  refLevel.transfer(newRefLevel);
537 
538  if (writeMesh)
539  {
540  Info<< "Writing refined mesh to time " << runTime.timeName() << nl
541  << endl;
542 
544  mesh.write();
545  refLevel.write();
546  }
547 
548  // Update cutCells for removed cells.
549  cutCells.updateMesh(morphMap());
550 }
551 
552 
553 // Divide the cells according to status compared to surface. Constructs sets:
554 // - cutCells : all cells cut by surface
555 // - inside : all cells inside surface
556 // - outside : ,, outside ,,
557 // and a combined set:
558 // - selected : sum of above according to selectCut/Inside/Outside flags.
559 void classifyCells
560 (
561  const polyMesh& mesh,
562  const fileName& surfName,
563  const triSurfaceSearch& querySurf,
564  const pointField& outsidePts,
565 
566  const bool selectCut,
567  const bool selectInside,
568  const bool selectOutside,
569 
570  const label nCutLayers,
571 
572  cellSet& inside,
573  cellSet& outside,
574  cellSet& cutCells,
575  cellSet& selected
576 )
577 {
578  // Cut faces with surface and classify cells
580  (
581  mesh,
582  surfName,
583  querySurf.surface(),
584  querySurf,
585  outsidePts,
586 
587  nCutLayers,
588 
589  inside,
590  outside,
591  cutCells
592  );
593 
594  // Combine wanted parts into selected
595  if (selectCut)
596  {
597  selected.addSet(cutCells);
598  }
599  if (selectInside)
600  {
601  selected.addSet(inside);
602  }
603  if (selectOutside)
604  {
605  selected.addSet(outside);
606  }
607 
608  Info<< "Determined cell status:" << endl
609  << " inside :" << inside.size() << endl
610  << " outside :" << outside.size() << endl
611  << " cutCells:" << cutCells.size() << endl
612  << " selected:" << selected.size() << endl
613  << endl;
614 
615  writeSet(inside, "inside cells");
616  writeSet(outside, "outside cells");
617  writeSet(cutCells, "cut cells");
618  writeSet(selected, "selected cells");
619 }
620 
621 
622 
623 int main(int argc, char *argv[])
624 {
626 
627  #include "setRootCase.H"
628  #include "createTime.H"
629  #include "createPolyMesh.H"
630 
631  // If nessecary add oldInternalFaces patch
632  label newPatchi = addPatch(mesh, "oldInternalFaces");
633 
634 
635  //
636  // Read motionProperties dictionary
637  //
638 
639  Info<< "Checking for motionProperties\n" << endl;
640 
641  IOobject motionObj
642  (
643  "motionProperties",
644  runTime.constant(),
645  mesh,
648  );
649 
650  // corrector for mesh motion
651  twoDPointCorrector* correct2DPtr = NULL;
652 
653  if (motionObj.headerOk())
654  {
655  Info<< "Reading " << runTime.constant() / "motionProperties"
656  << endl << endl;
657 
658  IOdictionary motionProperties(motionObj);
659 
660  Switch twoDMotion(motionProperties.lookup("twoDMotion"));
661 
662  if (twoDMotion)
663  {
664  Info<< "Correcting for 2D motion" << endl << endl;
665  correct2DPtr = new twoDPointCorrector(mesh);
666  }
667  }
668 
669  IOdictionary refineDict
670  (
671  IOobject
672  (
673  "autoRefineMeshDict",
674  runTime.system(),
675  mesh,
678  )
679  );
680 
681  fileName surfName(refineDict.lookup("surface"));
682  surfName.expand();
683  label nCutLayers(readLabel(refineDict.lookup("nCutLayers")));
684  label cellLimit(readLabel(refineDict.lookup("cellLimit")));
685  bool selectCut(readBool(refineDict.lookup("selectCut")));
686  bool selectInside(readBool(refineDict.lookup("selectInside")));
687  bool selectOutside(readBool(refineDict.lookup("selectOutside")));
688  bool selectHanging(readBool(refineDict.lookup("selectHanging")));
689 
690  scalar minEdgeLen(readScalar(refineDict.lookup("minEdgeLen")));
691  scalar maxEdgeLen(readScalar(refineDict.lookup("maxEdgeLen")));
692  scalar curvature(readScalar(refineDict.lookup("curvature")));
693  scalar curvDist(readScalar(refineDict.lookup("curvatureDistance")));
694  pointField outsidePts(refineDict.lookup("outsidePoints"));
695  label refinementLimit(readLabel(refineDict.lookup("splitLevel")));
696  bool writeMesh(readBool(refineDict.lookup("writeMesh")));
697 
698  Info<< "Cells to be used for meshing (0=false, 1=true):" << nl
699  << " cells cut by surface : " << selectCut << nl
700  << " cells inside of surface : " << selectInside << nl
701  << " cells outside of surface : " << selectOutside << nl
702  << " hanging cells : " << selectHanging << nl
703  << endl;
704 
705 
706  if (nCutLayers > 0 && selectInside)
707  {
709  << "Illogical settings : Both nCutLayers>0 and selectInside true."
710  << endl
711  << "This would mean that inside cells get removed but should be"
712  << " included in final mesh" << endl;
713  }
714 
715  // Surface.
716  triSurface surf(surfName);
717 
718  // Search engine on surface
719  triSurfaceSearch querySurf(surf);
720 
721  // Search engine on mesh. No face decomposition since mesh unwarped.
722  meshSearch queryMesh(mesh, polyMesh::FACE_PLANES);
723 
724  // Check all 'outside' points
725  forAll(outsidePts, outsideI)
726  {
727  const point& outsidePoint = outsidePts[outsideI];
728 
729  if (queryMesh.findCell(outsidePoint, -1, false) == -1)
730  {
732  << "outsidePoint " << outsidePoint
733  << " is not inside any cell"
734  << exit(FatalError);
735  }
736  }
737 
738 
739 
740  // Current refinement level. Read if present.
741  labelIOList refLevel
742  (
743  IOobject
744  (
745  "refinementLevel",
746  runTime.timeName(),
748  mesh,
751  ),
752  labelList(mesh.nCells(), 0)
753  );
754 
755  label maxLevel = max(refLevel);
756 
757  if (maxLevel > 0)
758  {
759  Info<< "Read existing refinement level from file "
760  << refLevel.objectPath() << nl
761  << " min level : " << min(refLevel) << nl
762  << " max level : " << maxLevel << nl
763  << endl;
764  }
765  else
766  {
767  Info<< "Created zero refinement level in file "
768  << refLevel.objectPath() << nl
769  << endl;
770  }
771 
772 
773 
774 
775  // Print edge stats on original mesh. Leave out 2d-normal direction
776  direction normalDir(getNormalDir(correct2DPtr));
777  scalar meshMinEdgeLen = getEdgeStats(mesh, normalDir);
778 
779  while (meshMinEdgeLen > minEdgeLen)
780  {
781  // Get inside/outside/cutCells cellSets.
782  cellSet inside(mesh, "inside", mesh.nCells()/10);
783  cellSet outside(mesh, "outside", mesh.nCells()/10);
784  cellSet cutCells(mesh, "cutCells", mesh.nCells()/10);
785  cellSet selected(mesh, "selected", mesh.nCells()/10);
786 
787  classifyCells
788  (
789  mesh,
790  surfName,
791  querySurf,
792  outsidePts,
793 
794  selectCut, // for determination of selected
795  selectInside, // for determination of selected
796  selectOutside, // for determination of selected
797 
798  nCutLayers, // How many layers of cutCells to include
799 
800  inside,
801  outside,
802  cutCells,
803  selected // not used but determined anyway.
804  );
805 
806  Info<< " Selected " << cutCells.name() << " with "
807  << cutCells.size() << " cells" << endl;
808 
809  if ((curvDist > 0) && (meshMinEdgeLen < maxEdgeLen))
810  {
811  // Done refining enough close to surface. Now switch to more
812  // intelligent refinement where only cutCells on surface curvature
813  // are refined.
814  cellSet curveCells(mesh, "curveCells", mesh.nCells()/10);
815 
816  selectCurvatureCells
817  (
818  mesh,
819  surfName,
820  querySurf,
821  maxEdgeLen,
822  curvature,
823  curveCells
824  );
825 
826  Info<< " Selected " << curveCells.name() << " with "
827  << curveCells.size() << " cells" << endl;
828 
829  // Add neighbours to cutCells. This is if selectCut is not
830  // true and so outside and/or inside cells get exposed there is
831  // also refinement in them.
832  if (!selectCut)
833  {
834  addCutNeighbours
835  (
836  mesh,
837  selectInside,
838  selectOutside,
839  inside,
840  outside,
841  cutCells
842  );
843  }
844 
845  // Subset cutCells to only curveCells
846  cutCells.subset(curveCells);
847 
848  Info<< " Removed cells not on surface curvature. Selected "
849  << cutCells.size() << endl;
850  }
851 
852 
853  if (nCutLayers > 0)
854  {
855  // Subset mesh to remove inside cells altogether. Updates cutCells,
856  // refLevel.
857  subsetMesh(mesh, writeMesh, newPatchi, inside, cutCells, refLevel);
858  }
859 
860 
861  // Added cells from 2:1 refinement level restriction.
862  while
863  (
864  limitRefinementLevel
865  (
866  mesh,
867  refinementLimit,
868  labelHashSet(),
869  refLevel,
870  cutCells
871  )
872  )
873  {}
874 
875 
876  Info<< " Current cells : " << mesh.nCells() << nl
877  << " Selected for refinement :" << cutCells.size() << nl
878  << endl;
879 
880  if (cutCells.empty())
881  {
882  Info<< "Stopping refining since 0 cells selected to be refined ..."
883  << nl << endl;
884  break;
885  }
886 
887  if ((mesh.nCells() + 8*cutCells.size()) > cellLimit)
888  {
889  Info<< "Stopping refining since cell limit reached ..." << nl
890  << "Would refine from " << mesh.nCells()
891  << " to " << mesh.nCells() + 8*cutCells.size() << " cells."
892  << nl << endl;
893  break;
894  }
895 
896  doRefinement
897  (
898  mesh,
899  refineDict,
900  cutCells,
901  refLevel
902  );
903 
904  Info<< " After refinement:" << mesh.nCells() << nl
905  << endl;
906 
907  if (writeMesh)
908  {
909  Info<< " Writing refined mesh to time " << runTime.timeName()
910  << nl << endl;
911 
913  mesh.write();
914  refLevel.write();
915  }
916 
917  // Update mesh edge stats.
918  meshMinEdgeLen = getEdgeStats(mesh, normalDir);
919  }
920 
921 
922  if (selectHanging)
923  {
924  // Get inside/outside/cutCells cellSets.
925  cellSet inside(mesh, "inside", mesh.nCells()/10);
926  cellSet outside(mesh, "outside", mesh.nCells()/10);
927  cellSet cutCells(mesh, "cutCells", mesh.nCells()/10);
928  cellSet selected(mesh, "selected", mesh.nCells()/10);
929 
930  classifyCells
931  (
932  mesh,
933  surfName,
934  querySurf,
935  outsidePts,
936 
937  selectCut,
938  selectInside,
939  selectOutside,
940 
941  nCutLayers,
942 
943  inside,
944  outside,
945  cutCells,
946  selected
947  );
948 
949 
950  // Find any cells which have all their points on the outside of the
951  // selected set and refine them
952  labelHashSet hanging = surfaceSets::getHangingCells(mesh, selected);
953 
954  Info<< "Detected " << hanging.size() << " hanging cells"
955  << " (cells with all points on"
956  << " outside of cellSet 'selected').\nThese would probably result"
957  << " in flattened cells when snapping the mesh to the surface"
958  << endl;
959 
960  Info<< "Refining " << hanging.size() << " hanging cells" << nl
961  << endl;
962 
963  // Added cells from 2:1 refinement level restriction.
964  while
965  (
966  limitRefinementLevel
967  (
968  mesh,
969  refinementLimit,
970  labelHashSet(),
971  refLevel,
972  hanging
973  )
974  )
975  {}
976 
977  doRefinement(mesh, refineDict, hanging, refLevel);
978 
979  Info<< "Writing refined mesh to time " << runTime.timeName() << nl
980  << endl;
981 
982  // Write final mesh
984  mesh.write();
985  refLevel.write();
986 
987  }
988  else if (!writeMesh)
989  {
990  Info<< "Writing refined mesh to time " << runTime.timeName() << nl
991  << endl;
992 
993  // Write final mesh. (will have been written already if writeMesh=true)
995  mesh.write();
996  refLevel.write();
997  }
998 
999 
1000  Info<< "End\n" << endl;
1001 
1002  return 0;
1003 }
1004 
1005 
1006 // ************************************************************************* //
const Time & time() const
Return time.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search...
Definition: meshSearch.H:57
const labelListList & cellCells() const
fileName objectPath() const
Return complete path + object name.
Definition: IOobject.H:363
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1048
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
uint8_t direction
Definition: direction.H:46
A class for handling file names.
Definition: fileName.H:69
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:219
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const triSurface & surface() const
Return reference to the surface.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
const vector & planeNormal() const
Return plane normal.
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:461
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Does multiple pass refinement to refine cells in multiple directions.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Given list of cells to remove insert all the topology changes.
Definition: removeCells.H:59
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
virtual void updateMesh(const mapPolyMesh &morphMap)
Update any stored data for new labels.
Definition: cellSet.C:225
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
virtual const pointField & points() const =0
Return mesh points.
static void noParallel()
Remove the parallel options.
Definition: argList.C:146
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
bool readBool(Istream &)
Definition: boolIO.C:60
patches[0]
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
const cellList & cells() const
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
label nCells() const
static labelHashSet getHangingCells(const primitiveMesh &mesh, const labelHashSet &internalCells)
Get cells using points on &#39;outside&#39; only.
Definition: surfaceSets.C:282
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
Helper class to search on triSurface.
scalar y
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
dynamicFvMesh & mesh
virtual void addSet(const topoSet &set)
Add elements present in set.
Definition: topoSet.C:486
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Class applies a two-dimensional correction to mesh motion point field.
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
const vectorField & cellCentres() const
virtual void subset(const topoSet &set)
Subset contents. Only elements present in both sets remain.
Definition: topoSet.C:467
List< label > labelList
A List of labels.
Definition: labelList.H:56
const fileName & local() const
Definition: IOobject.H:347
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
label readLabel(Istream &is)
Definition: label.H:64
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:198
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
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::polyBoundaryMesh.
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
static const char nl
Definition: Ostream.H:262
Empty front and back plane patch. Used for 2-D geometries.
label nEdges() const
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:878
A topoSetSource to select cells based on relation to surface.
Definition: surfaceToCell.H:60
void setSize(const label)
Reset size of List.
Definition: List.C:295
label patchi
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
#define WarningInFunction
Report a warning using Foam::Warning.
string & expand(const bool allowEmpty=false)
Expand initial tildes and all occurences of environment variables.
Definition: string.C:95
A collection of cell labels.
Definition: cellSet.H:48
Direct mesh changes based on v1.3 polyTopoChange syntax.
vector vec(const pointField &) const
Return the vector (end - start)
Definition: edgeI.H:157
virtual bool write() const
Write using setting from DB.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const word & system() const
Return system name.
Definition: TimePaths.H:114
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Triangulated surface description with patch information.
Definition: triSurface.H:65
label findPatchID(const word &patchName) const
Find patch index given a name.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const fileName & instance() const
Definition: IOobject.H:337
label nInternalFaces() const
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:365
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451