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