dynamicRefineFvMesh.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-2019 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 \*---------------------------------------------------------------------------*/
25 
26 #include "dynamicRefineFvMesh.H"
28 #include "surfaceInterpolate.H"
29 #include "volFields.H"
30 #include "polyTopoChange.H"
31 #include "surfaceFields.H"
32 #include "syncTools.H"
33 #include "pointFields.H"
34 #include "sigFpe.H"
35 #include "cellSet.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(dynamicRefineFvMesh, 0);
42  addToRunTimeSelectionTable(dynamicFvMesh, dynamicRefineFvMesh, IOobject);
43 }
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 (
49  const PackedBoolList& l,
50  const unsigned int val
51 )
52 {
53  label n = 0;
54  forAll(l, i)
55  {
56  if (l.get(i) == val)
57  {
58  n++;
59  }
60 
61  // debug also serves to get-around Clang compiler trying to optimise
62  // out this forAll loop under O3 optimisation
63  if (debug)
64  {
65  Info<< "n=" << n << endl;
66  }
67  }
68 
69  return n;
70 }
71 
72 
74 (
75  PackedBoolList& unrefineableCell
76 ) const
77 {
78  if (protectedCell_.empty())
79  {
80  unrefineableCell.clear();
81  return;
82  }
83 
84  const labelList& cellLevel = meshCutter_.cellLevel();
85 
86  unrefineableCell = protectedCell_;
87 
88  // Get neighbouring cell level
89  labelList neiLevel(nFaces()-nInternalFaces());
90 
91  for (label facei = nInternalFaces(); facei < nFaces(); facei++)
92  {
93  neiLevel[facei-nInternalFaces()] = cellLevel[faceOwner()[facei]];
94  }
95  syncTools::swapBoundaryFaceList(*this, neiLevel);
96 
97 
98  while (true)
99  {
100  // Pick up faces on border of protected cells
101  boolList seedFace(nFaces(), false);
102 
103  forAll(faceNeighbour(), facei)
104  {
105  label own = faceOwner()[facei];
106  bool ownProtected = unrefineableCell.get(own);
107  label nei = faceNeighbour()[facei];
108  bool neiProtected = unrefineableCell.get(nei);
109 
110  if (ownProtected && (cellLevel[nei] > cellLevel[own]))
111  {
112  seedFace[facei] = true;
113  }
114  else if (neiProtected && (cellLevel[own] > cellLevel[nei]))
115  {
116  seedFace[facei] = true;
117  }
118  }
119  for (label facei = nInternalFaces(); facei < nFaces(); facei++)
120  {
121  label own = faceOwner()[facei];
122  bool ownProtected = unrefineableCell.get(own);
123  if
124  (
125  ownProtected
126  && (neiLevel[facei-nInternalFaces()] > cellLevel[own])
127  )
128  {
129  seedFace[facei] = true;
130  }
131  }
132 
133  syncTools::syncFaceList(*this, seedFace, orEqOp<bool>());
134 
135 
136  // Extend unrefineableCell
137  bool hasExtended = false;
138 
139  for (label facei = 0; facei < nInternalFaces(); facei++)
140  {
141  if (seedFace[facei])
142  {
143  label own = faceOwner()[facei];
144  if (unrefineableCell.get(own) == 0)
145  {
146  unrefineableCell.set(own, 1);
147  hasExtended = true;
148  }
149 
150  label nei = faceNeighbour()[facei];
151  if (unrefineableCell.get(nei) == 0)
152  {
153  unrefineableCell.set(nei, 1);
154  hasExtended = true;
155  }
156  }
157  }
158  for (label facei = nInternalFaces(); facei < nFaces(); facei++)
159  {
160  if (seedFace[facei])
161  {
162  label own = faceOwner()[facei];
163  if (unrefineableCell.get(own) == 0)
164  {
165  unrefineableCell.set(own, 1);
166  hasExtended = true;
167  }
168  }
169  }
170 
171  if (!returnReduce(hasExtended, orOp<bool>()))
172  {
173  break;
174  }
175  }
176 }
177 
178 
180 {
181  const dictionary refineDict
182  (
183  dynamicMeshDict().optionalSubDict(typeName + "Coeffs")
184  );
185 
186  List<Pair<word>> fluxVelocities = List<Pair<word>>
187  (
188  refineDict.lookup("correctFluxes")
189  );
190  // Rework into hashtable.
191  correctFluxes_.resize(fluxVelocities.size());
192  forAll(fluxVelocities, i)
193  {
194  correctFluxes_.insert(fluxVelocities[i][0], fluxVelocities[i][1]);
195  }
196 
197  dumpLevel_ = Switch(refineDict.lookup("dumpLevel"));
198 }
199 
200 
201 // Refines cells, maps fields and recalculates (an approximate) flux
204 (
205  const labelList& cellsToRefine
206 )
207 {
208  // Mesh changing engine.
209  polyTopoChange meshMod(*this);
210 
211  // Play refinement commands into mesh changer.
212  meshCutter_.setRefinement(cellsToRefine, meshMod);
213 
214  // Create mesh (with inflation), return map from old to new mesh.
215  // autoPtr<mapPolyMesh> map = meshMod.changeMesh(*this, true);
216  autoPtr<mapPolyMesh> map = meshMod.changeMesh(*this, false);
217 
218  Info<< "Refined from "
219  << returnReduce(map().nOldCells(), sumOp<label>())
220  << " to " << globalData().nTotalCells() << " cells." << endl;
221 
222  if (debug)
223  {
224  // Check map.
225  for (label facei = 0; facei < nInternalFaces(); facei++)
226  {
227  label oldFacei = map().faceMap()[facei];
228 
229  if (oldFacei >= nInternalFaces())
230  {
232  << "New internal face:" << facei
233  << " fc:" << faceCentres()[facei]
234  << " originates from boundary oldFace:" << oldFacei
235  << abort(FatalError);
236  }
237  }
238  }
239 
240  // // Remove the stored tet base points
241  // tetBasePtIsPtr_.clear();
242  // // Remove the cell tree
243  // cellTreePtr_.clear();
244 
245  // Update fields
246  updateMesh(map);
247 
248 
249  // Move mesh
250  /*
251  pointField newPoints;
252  if (map().hasMotionPoints())
253  {
254  newPoints = map().preMotionPoints();
255  }
256  else
257  {
258  newPoints = points();
259  }
260  movePoints(newPoints);
261  */
262 
263  // Correct the flux for modified/added faces. All the faces which only
264  // have been renumbered will already have been handled by the mapping.
265  {
266  const labelList& faceMap = map().faceMap();
267  const labelList& reverseFaceMap = map().reverseFaceMap();
268 
269  // Storage for any master faces. These will be the original faces
270  // on the coarse cell that get split into four (or rather the
271  // master face gets modified and three faces get added from the master)
272  labelHashSet masterFaces(4*cellsToRefine.size());
273 
274  forAll(faceMap, facei)
275  {
276  label oldFacei = faceMap[facei];
277 
278  if (oldFacei >= 0)
279  {
280  label masterFacei = reverseFaceMap[oldFacei];
281 
282  if (masterFacei < 0)
283  {
285  << "Problem: should not have removed faces"
286  << " when refining."
287  << nl << "face:" << facei << abort(FatalError);
288  }
289  else if (masterFacei != facei)
290  {
291  masterFaces.insert(masterFacei);
292  }
293  }
294  }
295  if (debug)
296  {
297  Pout<< "Found " << masterFaces.size() << " split faces " << endl;
298  }
299 
301  (
302  lookupClass<surfaceScalarField>()
303  );
305  {
306  if (!correctFluxes_.found(iter.key()))
307  {
309  << "Cannot find surfaceScalarField " << iter.key()
310  << " in user-provided flux mapping table "
311  << correctFluxes_ << endl
312  << " The flux mapping table is used to recreate the"
313  << " flux on newly created faces." << endl
314  << " Either add the entry if it is a flux or use ("
315  << iter.key() << " none) to suppress this warning."
316  << endl;
317  continue;
318  }
319 
320  const word& UName = correctFluxes_[iter.key()];
321 
322  if (UName == "none")
323  {
324  continue;
325  }
326 
327  if (UName == "NaN")
328  {
329  Pout<< "Setting surfaceScalarField " << iter.key()
330  << " to NaN" << endl;
331 
332  surfaceScalarField& phi = *iter();
333 
335 
336  continue;
337  }
338 
339  if (debug)
340  {
341  Pout<< "Mapping flux " << iter.key()
342  << " using interpolated flux " << UName
343  << endl;
344  }
345 
346  surfaceScalarField& phi = *iter();
347  const surfaceScalarField phiU
348  (
350  (
351  lookupObject<volVectorField>(UName)
352  )
353  & Sf()
354  );
355 
356  // Recalculate new internal faces.
357  for (label facei = 0; facei < nInternalFaces(); facei++)
358  {
359  label oldFacei = faceMap[facei];
360 
361  if (oldFacei == -1)
362  {
363  // Inflated/appended
364  phi[facei] = phiU[facei];
365  }
366  else if (reverseFaceMap[oldFacei] != facei)
367  {
368  // face-from-masterface
369  phi[facei] = phiU[facei];
370  }
371  }
372 
373  // Recalculate new boundary faces.
374  surfaceScalarField::Boundary& phiBf =
375  phi.boundaryFieldRef();
376  forAll(phiBf, patchi)
377  {
378  fvsPatchScalarField& patchPhi = phiBf[patchi];
379  const fvsPatchScalarField& patchPhiU =
380  phiU.boundaryField()[patchi];
381 
382  label facei = patchPhi.patch().start();
383 
384  forAll(patchPhi, i)
385  {
386  label oldFacei = faceMap[facei];
387 
388  if (oldFacei == -1)
389  {
390  // Inflated/appended
391  patchPhi[i] = patchPhiU[i];
392  }
393  else if (reverseFaceMap[oldFacei] != facei)
394  {
395  // face-from-masterface
396  patchPhi[i] = patchPhiU[i];
397  }
398 
399  facei++;
400  }
401  }
402 
403  // Update master faces
404  forAllConstIter(labelHashSet, masterFaces, iter)
405  {
406  label facei = iter.key();
407 
408  if (isInternalFace(facei))
409  {
410  phi[facei] = phiU[facei];
411  }
412  else
413  {
414  label patchi = boundaryMesh().whichPatch(facei);
415  label i = facei - boundaryMesh()[patchi].start();
416 
417  const fvsPatchScalarField& patchPhiU =
418  phiU.boundaryField()[patchi];
419 
420  fvsPatchScalarField& patchPhi = phiBf[patchi];
421 
422  patchPhi[i] = patchPhiU[i];
423  }
424  }
425  }
426  }
427 
428 
429 
430  // Update numbering of cells/vertices.
431  meshCutter_.updateMesh(map);
432 
433  // Update numbering of protectedCell_
434  if (protectedCell_.size())
435  {
436  PackedBoolList newProtectedCell(nCells());
437 
438  forAll(newProtectedCell, celli)
439  {
440  label oldCelli = map().cellMap()[celli];
441  newProtectedCell.set(celli, protectedCell_.get(oldCelli));
442  }
443  protectedCell_.transfer(newProtectedCell);
444  }
445 
446  // Debug: Check refinement levels (across faces only)
447  meshCutter_.checkRefinementLevels(-1, labelList(0));
448 
449  return map;
450 }
451 
452 
455 (
456  const labelList& splitPoints
457 )
458 {
459  polyTopoChange meshMod(*this);
460 
461  // Play refinement commands into mesh changer.
462  meshCutter_.setUnrefinement(splitPoints, meshMod);
463 
464 
465  // Save information on faces that will be combined
466  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
467 
468  // Find the faceMidPoints on cells to be combined.
469  // for each face resulting of split of face into four store the
470  // midpoint
471  Map<label> faceToSplitPoint(3*splitPoints.size());
472 
473  {
474  forAll(splitPoints, i)
475  {
476  label pointi = splitPoints[i];
477 
478  const labelList& pEdges = pointEdges()[pointi];
479 
480  forAll(pEdges, j)
481  {
482  label otherPointi = edges()[pEdges[j]].otherVertex(pointi);
483 
484  const labelList& pFaces = pointFaces()[otherPointi];
485 
486  forAll(pFaces, pFacei)
487  {
488  faceToSplitPoint.insert(pFaces[pFacei], otherPointi);
489  }
490  }
491  }
492  }
493 
494 
495  // Change mesh and generate map.
496  // autoPtr<mapPolyMesh> map = meshMod.changeMesh(*this, true);
497  autoPtr<mapPolyMesh> map = meshMod.changeMesh(*this, false);
498 
499  Info<< "Unrefined from "
500  << returnReduce(map().nOldCells(), sumOp<label>())
501  << " to " << globalData().nTotalCells() << " cells."
502  << endl;
503 
504  // Update fields
505  updateMesh(map);
506 
507 
508  // Move mesh
509  /*
510  pointField newPoints;
511  if (map().hasMotionPoints())
512  {
513  newPoints = map().preMotionPoints();
514  }
515  else
516  {
517  newPoints = points();
518  }
519  movePoints(newPoints);
520  */
521 
522  // Correct the flux for modified faces.
523  {
524  const labelList& reversePointMap = map().reversePointMap();
525  const labelList& reverseFaceMap = map().reverseFaceMap();
526 
528  (
529  lookupClass<surfaceScalarField>()
530  );
532  {
533  if (!correctFluxes_.found(iter.key()))
534  {
536  << "Cannot find surfaceScalarField " << iter.key()
537  << " in user-provided flux mapping table "
538  << correctFluxes_ << endl
539  << " The flux mapping table is used to recreate the"
540  << " flux on newly created faces." << endl
541  << " Either add the entry if it is a flux or use ("
542  << iter.key() << " none) to suppress this warning."
543  << endl;
544  continue;
545  }
546 
547  const word& UName = correctFluxes_[iter.key()];
548 
549  if (UName == "none")
550  {
551  continue;
552  }
553 
554  if (debug)
555  {
556  Info<< "Mapping flux " << iter.key()
557  << " using interpolated flux " << UName
558  << endl;
559  }
560 
561  surfaceScalarField& phi = *iter();
562  surfaceScalarField::Boundary& phiBf =
563  phi.boundaryFieldRef();
564 
565  const surfaceScalarField phiU
566  (
568  (
569  lookupObject<volVectorField>(UName)
570  )
571  & Sf()
572  );
573 
574 
575  forAllConstIter(Map<label>, faceToSplitPoint, iter)
576  {
577  label oldFacei = iter.key();
578  label oldPointi = iter();
579 
580  if (reversePointMap[oldPointi] < 0)
581  {
582  // midpoint was removed. See if face still exists.
583  label facei = reverseFaceMap[oldFacei];
584 
585  if (facei >= 0)
586  {
587  if (isInternalFace(facei))
588  {
589  phi[facei] = phiU[facei];
590  }
591  else
592  {
593  label patchi = boundaryMesh().whichPatch(facei);
594  label i = facei - boundaryMesh()[patchi].start();
595 
596  const fvsPatchScalarField& patchPhiU =
597  phiU.boundaryField()[patchi];
598  fvsPatchScalarField& patchPhi = phiBf[patchi];
599  patchPhi[i] = patchPhiU[i];
600  }
601  }
602  }
603  }
604  }
605  }
606 
607 
608  // Update numbering of cells/vertices.
609  meshCutter_.updateMesh(map);
610 
611  // Update numbering of protectedCell_
612  if (protectedCell_.size())
613  {
614  PackedBoolList newProtectedCell(nCells());
615 
616  forAll(newProtectedCell, celli)
617  {
618  label oldCelli = map().cellMap()[celli];
619  if (oldCelli >= 0)
620  {
621  newProtectedCell.set(celli, protectedCell_.get(oldCelli));
622  }
623  }
624  protectedCell_.transfer(newProtectedCell);
625  }
626 
627  // Debug: Check refinement levels (across faces only)
628  meshCutter_.checkRefinementLevels(-1, labelList(0));
629 
630  return map;
631 }
632 
633 
636 {
637  scalarField vFld(nCells(), -great);
638 
639  forAll(pointCells(), pointi)
640  {
641  const labelList& pCells = pointCells()[pointi];
642 
643  forAll(pCells, i)
644  {
645  vFld[pCells[i]] = max(vFld[pCells[i]], pFld[pointi]);
646  }
647  }
648  return vFld;
649 }
650 
651 
654 {
655  scalarField pFld(nPoints(), -great);
656 
657  forAll(pointCells(), pointi)
658  {
659  const labelList& pCells = pointCells()[pointi];
660 
661  forAll(pCells, i)
662  {
663  pFld[pointi] = max(pFld[pointi], vFld[pCells[i]]);
664  }
665  }
666  return pFld;
667 }
668 
669 
672 {
673  scalarField pFld(nPoints());
674 
675  forAll(pointCells(), pointi)
676  {
677  const labelList& pCells = pointCells()[pointi];
678 
679  scalar sum = 0.0;
680  forAll(pCells, i)
681  {
682  sum += vFld[pCells[i]];
683  }
684  pFld[pointi] = sum/pCells.size();
685  }
686  return pFld;
687 }
688 
689 
691 (
692  const scalarField& fld,
693  const scalar minLevel,
694  const scalar maxLevel
695 ) const
696 {
697  scalarField c(fld.size(), -1);
698 
699  forAll(fld, i)
700  {
701  scalar err = min(fld[i]-minLevel, maxLevel-fld[i]);
702 
703  if (err >= 0)
704  {
705  c[i] = err;
706  }
707  }
708  return c;
709 }
710 
711 
713 (
714  const scalar lowerRefineLevel,
715  const scalar upperRefineLevel,
716  const scalarField& vFld,
717  PackedBoolList& candidateCell
718 ) const
719 {
720  // Get error per cell. Is -1 (not to be refined) to >0 (to be refined,
721  // higher more desirable to be refined).
722  scalarField cellError
723  (
724  maxPointField
725  (
726  error
727  (
728  cellToPoint(vFld),
729  lowerRefineLevel,
730  upperRefineLevel
731  )
732  )
733  );
734 
735  // Mark cells that are candidates for refinement.
736  forAll(cellError, celli)
737  {
738  if (cellError[celli] > 0)
739  {
740  candidateCell.set(celli, 1);
741  }
742  }
743 }
744 
745 
747 (
748  const label maxCells,
749  const label maxRefinement,
750  const PackedBoolList& candidateCell
751 ) const
752 {
753  // Every refined cell causes 7 extra cells
754  label nTotToRefine = (maxCells - globalData().nTotalCells()) / 7;
755 
756  const labelList& cellLevel = meshCutter_.cellLevel();
757 
758  // Mark cells that cannot be refined since they would trigger refinement
759  // of protected cells (since 2:1 cascade)
760  PackedBoolList unrefineableCell;
761  calculateProtectedCells(unrefineableCell);
762 
763  // Count current selection
764  label nLocalCandidates = count(candidateCell, 1);
765  label nCandidates = returnReduce(nLocalCandidates, sumOp<label>());
766 
767  // Collect all cells
768  DynamicList<label> candidates(nLocalCandidates);
769 
770  if (nCandidates < nTotToRefine)
771  {
772  forAll(candidateCell, celli)
773  {
774  if
775  (
776  cellLevel[celli] < maxRefinement
777  && candidateCell.get(celli)
778  && (
779  unrefineableCell.empty()
780  || !unrefineableCell.get(celli)
781  )
782  )
783  {
784  candidates.append(celli);
785  }
786  }
787  }
788  else
789  {
790  // Sort by error? For now just truncate.
791  for (label level = 0; level < maxRefinement; level++)
792  {
793  forAll(candidateCell, celli)
794  {
795  if
796  (
797  cellLevel[celli] == level
798  && candidateCell.get(celli)
799  && (
800  unrefineableCell.empty()
801  || !unrefineableCell.get(celli)
802  )
803  )
804  {
805  candidates.append(celli);
806  }
807  }
808 
809  if (returnReduce(candidates.size(), sumOp<label>()) > nTotToRefine)
810  {
811  break;
812  }
813  }
814  }
815 
816  // Guarantee 2:1 refinement after refinement
817  labelList consistentSet
818  (
819  meshCutter_.consistentRefinement
820  (
821  candidates.shrink(),
822  true // Add to set to guarantee 2:1
823  )
824  );
825 
826  Info<< "Selected " << returnReduce(consistentSet.size(), sumOp<label>())
827  << " cells for refinement out of " << globalData().nTotalCells()
828  << "." << endl;
829 
830  return consistentSet;
831 }
832 
833 
835 (
836  const scalar unrefineLevel,
837  const PackedBoolList& markedCell,
838  const scalarField& pFld
839 ) const
840 {
841  // All points that can be unrefined
842  const labelList splitPoints(meshCutter_.getSplitPoints());
843 
844  DynamicList<label> newSplitPoints(splitPoints.size());
845 
846  forAll(splitPoints, i)
847  {
848  label pointi = splitPoints[i];
849 
850  if (pFld[pointi] < unrefineLevel)
851  {
852  // Check that all cells are not marked
853  const labelList& pCells = pointCells()[pointi];
854 
855  bool hasMarked = false;
856 
857  forAll(pCells, pCelli)
858  {
859  if (markedCell.get(pCells[pCelli]))
860  {
861  hasMarked = true;
862  break;
863  }
864  }
865 
866  if (!hasMarked)
867  {
868  newSplitPoints.append(pointi);
869  }
870  }
871  }
872 
873 
874  newSplitPoints.shrink();
875 
876  // Guarantee 2:1 refinement after unrefinement
877  labelList consistentSet
878  (
879  meshCutter_.consistentUnrefinement
880  (
881  newSplitPoints,
882  false
883  )
884  );
885  Info<< "Selected " << returnReduce(consistentSet.size(), sumOp<label>())
886  << " split points out of a possible "
887  << returnReduce(splitPoints.size(), sumOp<label>())
888  << "." << endl;
889 
890  return consistentSet;
891 }
892 
893 
895 (
896  PackedBoolList& markedCell
897 ) const
898 {
899  // Mark faces using any marked cell
900  boolList markedFace(nFaces(), false);
901 
902  forAll(markedCell, celli)
903  {
904  if (markedCell.get(celli))
905  {
906  const cell& cFaces = cells()[celli];
907 
908  forAll(cFaces, i)
909  {
910  markedFace[cFaces[i]] = true;
911  }
912  }
913  }
914 
915  syncTools::syncFaceList(*this, markedFace, orEqOp<bool>());
916 
917  // Update cells using any markedFace
918  for (label facei = 0; facei < nInternalFaces(); facei++)
919  {
920  if (markedFace[facei])
921  {
922  markedCell.set(faceOwner()[facei], 1);
923  markedCell.set(faceNeighbour()[facei], 1);
924  }
925  }
926  for (label facei = nInternalFaces(); facei < nFaces(); facei++)
927  {
928  if (markedFace[facei])
929  {
930  markedCell.set(faceOwner()[facei], 1);
931  }
932  }
933 }
934 
935 
937 (
938  PackedBoolList& protectedCell,
939  label& nProtected
940 ) const
941 {
942  const labelList& cellLevel = meshCutter_.cellLevel();
943  const labelList& pointLevel = meshCutter_.pointLevel();
944 
945  labelList nAnchorPoints(nCells(), 0);
946 
947  forAll(pointLevel, pointi)
948  {
949  const labelList& pCells = pointCells(pointi);
950 
951  forAll(pCells, pCelli)
952  {
953  label celli = pCells[pCelli];
954 
955  if (pointLevel[pointi] <= cellLevel[celli])
956  {
957  // Check if cell has already 8 anchor points -> protect cell
958  if (nAnchorPoints[celli] == 8)
959  {
960  if (protectedCell.set(celli, true))
961  {
962  nProtected++;
963  }
964  }
965 
966  if (!protectedCell[celli])
967  {
968  nAnchorPoints[celli]++;
969  }
970  }
971  }
972  }
973 
974 
975  forAll(protectedCell, celli)
976  {
977  if (!protectedCell[celli] && nAnchorPoints[celli] != 8)
978  {
979  protectedCell.set(celli, true);
980  nProtected++;
981  }
982  }
983 }
984 
985 
986 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
987 
989 :
990  dynamicFvMesh(io),
991  meshCutter_(*this),
992  dumpLevel_(false),
993  nRefinementIterations_(0),
994  protectedCell_(nCells(), 0)
995 {
996  // Read static part of dictionary
997  readDict();
998 
999 
1000  const labelList& cellLevel = meshCutter_.cellLevel();
1001  const labelList& pointLevel = meshCutter_.pointLevel();
1002 
1003  // Set cells that should not be refined.
1004  // This is currently any cell which does not have 8 anchor points or
1005  // uses any face which does not have 4 anchor points.
1006  // Note: do not use cellPoint addressing
1007 
1008  // Count number of points <= cellLevel
1009  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1010 
1011  labelList nAnchors(nCells(), 0);
1012 
1013  label nProtected = 0;
1014 
1015  forAll(pointCells(), pointi)
1016  {
1017  const labelList& pCells = pointCells()[pointi];
1018 
1019  forAll(pCells, i)
1020  {
1021  label celli = pCells[i];
1022 
1023  if (!protectedCell_.get(celli))
1024  {
1025  if (pointLevel[pointi] <= cellLevel[celli])
1026  {
1027  nAnchors[celli]++;
1028 
1029  if (nAnchors[celli] > 8)
1030  {
1031  protectedCell_.set(celli, 1);
1032  nProtected++;
1033  }
1034  }
1035  }
1036  }
1037  }
1038 
1039 
1040  // Count number of points <= faceLevel
1041  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1042  // Bit tricky since proc face might be one more refined than the owner since
1043  // the coupled one is refined.
1044 
1045  {
1046  labelList neiLevel(nFaces());
1047 
1048  for (label facei = 0; facei < nInternalFaces(); facei++)
1049  {
1050  neiLevel[facei] = cellLevel[faceNeighbour()[facei]];
1051  }
1052  for (label facei = nInternalFaces(); facei < nFaces(); facei++)
1053  {
1054  neiLevel[facei] = cellLevel[faceOwner()[facei]];
1055  }
1056  syncTools::swapFaceList(*this, neiLevel);
1057 
1058 
1059  boolList protectedFace(nFaces(), false);
1060 
1061  forAll(faceOwner(), facei)
1062  {
1063  label faceLevel = max
1064  (
1065  cellLevel[faceOwner()[facei]],
1066  neiLevel[facei]
1067  );
1068 
1069  const face& f = faces()[facei];
1070 
1071  label nAnchors = 0;
1072 
1073  forAll(f, fp)
1074  {
1075  if (pointLevel[f[fp]] <= faceLevel)
1076  {
1077  nAnchors++;
1078 
1079  if (nAnchors > 4)
1080  {
1081  protectedFace[facei] = true;
1082  break;
1083  }
1084  }
1085  }
1086  }
1087 
1088  syncTools::syncFaceList(*this, protectedFace, orEqOp<bool>());
1089 
1090  for (label facei = 0; facei < nInternalFaces(); facei++)
1091  {
1092  if (protectedFace[facei])
1093  {
1094  protectedCell_.set(faceOwner()[facei], 1);
1095  nProtected++;
1096  protectedCell_.set(faceNeighbour()[facei], 1);
1097  nProtected++;
1098  }
1099  }
1100  for (label facei = nInternalFaces(); facei < nFaces(); facei++)
1101  {
1102  if (protectedFace[facei])
1103  {
1104  protectedCell_.set(faceOwner()[facei], 1);
1105  nProtected++;
1106  }
1107  }
1108 
1109  // Also protect any cells that are less than hex
1110  forAll(cells(), celli)
1111  {
1112  const cell& cFaces = cells()[celli];
1113 
1114  if (cFaces.size() < 6)
1115  {
1116  if (protectedCell_.set(celli, 1))
1117  {
1118  nProtected++;
1119  }
1120  }
1121  else
1122  {
1123  forAll(cFaces, cFacei)
1124  {
1125  if (faces()[cFaces[cFacei]].size() < 4)
1126  {
1127  if (protectedCell_.set(celli, 1))
1128  {
1129  nProtected++;
1130  }
1131  break;
1132  }
1133  }
1134  }
1135  }
1136 
1137  // Check cells for 8 corner points
1139  }
1140 
1141  if (returnReduce(nProtected, sumOp<label>()) == 0)
1142  {
1144  }
1145  else
1146  {
1147 
1148  cellSet protectedCells(*this, "protectedCells", nProtected);
1149  forAll(protectedCell_, celli)
1150  {
1151  if (protectedCell_[celli])
1152  {
1153  protectedCells.insert(celli);
1154  }
1155  }
1156 
1157  Info<< "Detected " << returnReduce(nProtected, sumOp<label>())
1158  << " cells that are protected from refinement."
1159  << " Writing these to cellSet "
1160  << protectedCells.name()
1161  << "." << endl;
1162 
1163  protectedCells.write();
1164  }
1165 }
1166 
1167 
1168 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1169 
1171 {}
1172 
1173 
1174 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1175 
1177 {
1178  // Re-read dictionary. Chosen since usually -small so trivial amount
1179  // of time compared to actual refinement. Also very useful to be able
1180  // to modify on-the-fly.
1181  const dictionary refineDict
1182  (
1183  dynamicMeshDict().optionalSubDict(typeName + "Coeffs")
1184  );
1185 
1186  label refineInterval = refineDict.lookup<label>("refineInterval");
1187 
1188  bool hasChanged = false;
1189 
1190  if (refineInterval == 0)
1191  {
1192  topoChanging(hasChanged);
1193 
1194  return false;
1195  }
1196  else if (refineInterval < 0)
1197  {
1199  << "Illegal refineInterval " << refineInterval << nl
1200  << "The refineInterval setting in the dynamicMeshDict should"
1201  << " be >= 1." << nl
1202  << exit(FatalError);
1203  }
1204 
1205  // Note: cannot refine at time 0 since no V0 present since mesh not
1206  // moved yet.
1207 
1208  if (time().timeIndex() > 0 && time().timeIndex() % refineInterval == 0)
1209  {
1210  label maxCells = refineDict.lookup<label>("maxCells");
1211 
1212  if (maxCells <= 0)
1213  {
1215  << "Illegal maximum number of cells " << maxCells << nl
1216  << "The maxCells setting in the dynamicMeshDict should"
1217  << " be > 0." << nl
1218  << exit(FatalError);
1219  }
1220 
1221  label maxRefinement = refineDict.lookup<label>("maxRefinement");
1222 
1223  if (maxRefinement <= 0)
1224  {
1226  << "Illegal maximum refinement level " << maxRefinement << nl
1227  << "The maxCells setting in the dynamicMeshDict should"
1228  << " be > 0." << nl
1229  << exit(FatalError);
1230  }
1231 
1232  const word fieldName(refineDict.lookup("field"));
1233 
1234  const volScalarField& vFld = lookupObject<volScalarField>(fieldName);
1235 
1236  const scalar lowerRefineLevel =
1237  refineDict.lookup<scalar>("lowerRefineLevel");
1238  const scalar upperRefineLevel =
1239  refineDict.lookup<scalar>("upperRefineLevel");
1240  const scalar unrefineLevel = refineDict.lookupOrDefault<scalar>
1241  (
1242  "unrefineLevel",
1243  great
1244  );
1245  const label nBufferLayers =
1246  refineDict.lookup<label>("nBufferLayers");
1247 
1248  // Cells marked for refinement or otherwise protected from unrefinement.
1250 
1251  // Determine candidates for refinement (looking at field only)
1253  (
1254  lowerRefineLevel,
1255  upperRefineLevel,
1256  vFld,
1257  refineCell
1258  );
1259 
1260  if (globalData().nTotalCells() < maxCells)
1261  {
1262  // Select subset of candidates. Take into account max allowable
1263  // cells, refinement level, protected cells.
1264  labelList cellsToRefine
1265  (
1267  (
1268  maxCells,
1269  maxRefinement,
1270  refineCell
1271  )
1272  );
1273 
1274  label nCellsToRefine = returnReduce
1275  (
1276  cellsToRefine.size(), sumOp<label>()
1277  );
1278 
1279  if (nCellsToRefine > 0)
1280  {
1281  // Refine/update mesh and map fields
1282  autoPtr<mapPolyMesh> map = refine(cellsToRefine);
1283 
1284  // Update refineCell. Note that some of the marked ones have
1285  // not been refined due to constraints.
1286  {
1287  const labelList& cellMap = map().cellMap();
1288  const labelList& reverseCellMap = map().reverseCellMap();
1289 
1290  PackedBoolList newRefineCell(cellMap.size());
1291 
1292  forAll(cellMap, celli)
1293  {
1294  label oldCelli = cellMap[celli];
1295 
1296  if (oldCelli < 0)
1297  {
1298  newRefineCell.set(celli, 1);
1299  }
1300  else if (reverseCellMap[oldCelli] != celli)
1301  {
1302  newRefineCell.set(celli, 1);
1303  }
1304  else
1305  {
1306  newRefineCell.set(celli, refineCell.get(oldCelli));
1307  }
1308  }
1309  refineCell.transfer(newRefineCell);
1310  }
1311 
1312  // Extend with a buffer layer to prevent neighbouring points
1313  // being unrefined.
1314  for (label i = 0; i < nBufferLayers; i++)
1315  {
1316  extendMarkedCells(refineCell);
1317  }
1318 
1319  hasChanged = true;
1320  }
1321  }
1322 
1323 
1324  {
1325  // Select unrefineable points that are not marked in refineCell
1326  labelList pointsToUnrefine
1327  (
1329  (
1330  unrefineLevel,
1331  refineCell,
1332  maxCellField(vFld)
1333  )
1334  );
1335 
1336  label nSplitPoints = returnReduce
1337  (
1338  pointsToUnrefine.size(),
1339  sumOp<label>()
1340  );
1341 
1342  if (nSplitPoints > 0)
1343  {
1344  // Refine/update mesh
1345  unrefine(pointsToUnrefine);
1346 
1347  hasChanged = true;
1348  }
1349  }
1350 
1351 
1352  if ((nRefinementIterations_ % 10) == 0)
1353  {
1354  // Compact refinement history occasionally (how often?).
1355  // Unrefinement causes holes in the refinementHistory.
1356  const_cast<refinementHistory&>(meshCutter().history()).compact();
1357  }
1359  }
1360 
1361  topoChanging(hasChanged);
1362  if (hasChanged)
1363  {
1364  // Reset moving flag (if any). If not using inflation we'll not move,
1365  // if are using inflation any follow on movePoints will set it.
1366  moving(false);
1367  }
1368 
1369  return hasChanged;
1370 }
1371 
1372 
1378  const bool write
1379 ) const
1380 {
1381  // Force refinement data to go to the current time directory.
1382  const_cast<hexRef8&>(meshCutter_).setInstance(time().timeName());
1383 
1384  bool writeOk =
1385  (
1386  dynamicFvMesh::writeObject(fmt, ver, cmp, write)
1387  && meshCutter_.write(write)
1388  );
1389 
1390  if (dumpLevel_)
1391  {
1392  volScalarField scalarCellLevel
1393  (
1394  IOobject
1395  (
1396  "cellLevel",
1397  time().timeName(),
1398  *this,
1401  false
1402  ),
1403  *this,
1405  );
1406 
1407  const labelList& cellLevel = meshCutter_.cellLevel();
1408 
1409  forAll(cellLevel, celli)
1410  {
1411  scalarCellLevel[celli] = cellLevel[celli];
1412  }
1413 
1414  writeOk = writeOk && scalarCellLevel.write();
1415  }
1416 
1417  return writeOk;
1418 }
1419 
1420 
1421 // ************************************************************************* //
Foam::surfaceFields.
const refinementHistory & history() const
Definition: hexRef8.H:395
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual labelList selectRefineCells(const label maxCells, const label maxRefinement, const PackedBoolList &candidateCell) const
Subset candidate cells for refinement.
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
const word & name() const
Return name.
Definition: IOobject.H:303
scalarField maxPointField(const scalarField &) const
Get per cell max of connected point.
static label count(const PackedBoolList &, const unsigned int)
Count set/unset elements in packedlist.
bool moving() const
Is mesh moving.
Definition: polyMesh.H:512
const hexRef8 & meshCutter() const
Direct access to the refinement engine.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
PackedBoolList protectedCell_
Protected cells (usually since not hexes)
unsigned int get(const label) const
Get value at index I.
Definition: PackedListI.H:954
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const labelIOList & pointLevel() const
Definition: hexRef8.H:390
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:149
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1175
const Boundary & boundaryField() const
Return const-reference to the boundary field.
label nFaces() const
scalarField cellToPoint(const scalarField &vFld) const
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label nCells() const
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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/any.
Definition: Switch.H:60
void checkEightAnchorPoints(PackedBoolList &protectedCell, label &nProtected) const
Check all cells have 8 anchor points.
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write using given format, version and compression.
const cellList & cells() const
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
const dimensionedScalar & c
Speed of light in a vacuum.
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
void calculateProtectedCells(PackedBoolList &unrefineableCell) const
Calculate cells that cannot be refined since would trigger.
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Macros for easy insertion into run-time selection tables.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool write=true) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:1012
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
phi
Definition: pEqn.H:104
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:66
virtual bool update()
Update the mesh for both mesh motion and topology change.
const cellShapeList & cells
void set(const PackedList< 1 > &)
Set specified bits.
autoPtr< mapPolyMesh > unrefine(const labelList &)
Unrefine cells. Gets passed in centre points of cells to combine.
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1001
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:64
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
hexRef8 meshCutter_
Mesh cutting engine.
void extendMarkedCells(PackedBoolList &markedCell) const
Extend markedCell with cell-face-cell.
bool write(const bool write=true) const
Force writing refinement+history to polyMesh directory.
Definition: hexRef8.C:5746
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
autoPtr< mapPolyMesh > refine(const labelList &)
Refine cells. Update mesh and fields.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
static void swapFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled face values.
Definition: syncTools.H:463
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
static void fillNan(UList< scalar > &)
Fill block of data with NaN.
Definition: sigFpe.C:46
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:56
An STL-conforming hash table.
Definition: HashTable.H:61
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void clear()
Clear the list, i.e. set addressable size to zero.
Definition: PackedListI.H:889
const labelListList & pointCells() const
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< small) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &mergedCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:252
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
A topoSetSource to select points based on usage in cells.
Definition: cellToPoint.H:49
labelList f(nPoints)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
label nRefinementIterations_
Number of refinement/unrefinement steps done so far.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
virtual void selectRefineCandidates(const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, PackedBoolList &candidateCell) const
Select candidate cells for refinement.
const fvPatch & patch() const
Return patch.
void readDict()
Read the projection parameters from dictionary.
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:381
A bit-packed bool list.
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
virtual labelList selectUnrefinePoints(const scalar unrefineLevel, const PackedBoolList &markedCell, const scalarField &pFld) const
Select points that can be unrefined.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dynamicRefineFvMesh(const IOobject &io)
Construct from IOobject.
Abstract base class for geometry and/or topology changing fvMesh.
Definition: dynamicFvMesh.H:49
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
A collection of cell labels.
Definition: cellSet.H:48
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Direct mesh changes based on v1.3 polyTopoChange syntax.
bool topoChanging() const
Is mesh topology changing.
Definition: polyMesh.H:526
label timeIndex
Definition: getTimeIndex.H:4
Version number type.
Definition: IOstream.H:96
virtual ~dynamicRefineFvMesh()
Destructor.
messageStream Info
label n
bool empty() const
Return true if the list is empty (ie, size() is zero).
Definition: PackedListI.H:718
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
virtual bool write(const bool write=true) const
Write using setting from DB.
All refinement history. Used in unrefinement.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
scalarField error(const scalarField &fld, const scalar minLevel, const scalar maxLevel) const
const labelIOList & cellLevel() const
Definition: hexRef8.H:385
Switch dumpLevel_
Dump cellLevel for postprocessing.
Namespace for OpenFOAM.
scalarField maxCellField(const volScalarField &) const
Get point max of connected cell.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
const dictionary & dynamicMeshDict() const
Return the dynamicMeshDict.