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-2015 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  {
244  FatalErrorIn("dynamicRefineFvMesh::refine(const labelList&)")
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  (
299  "dynamicRefineFvMesh::refine(const labelList&)"
300  ) << "Problem: should not have removed faces"
301  << " when refining."
302  << nl << "face:" << faceI << abort(FatalError);
303  }
304  else if (masterFaceI != faceI)
305  {
306  masterFaces.insert(masterFaceI);
307  }
308  }
309  }
310  if (debug)
311  {
312  Pout<< "Found " << masterFaces.size() << " split faces " << endl;
313  }
314 
316  (
317  lookupClass<surfaceScalarField>()
318  );
320  {
321  if (!correctFluxes_.found(iter.key()))
322  {
323  WarningIn("dynamicRefineFvMesh::refine(const labelList&)")
324  << "Cannot find surfaceScalarField " << iter.key()
325  << " in user-provided flux mapping table "
326  << correctFluxes_ << endl
327  << " The flux mapping table is used to recreate the"
328  << " flux on newly created faces." << endl
329  << " Either add the entry if it is a flux or use ("
330  << iter.key() << " none) to suppress this warning."
331  << endl;
332  continue;
333  }
334 
335  const word& UName = correctFluxes_[iter.key()];
336 
337  if (UName == "none")
338  {
339  continue;
340  }
341 
342  if (UName == "NaN")
343  {
344  Pout<< "Setting surfaceScalarField " << iter.key()
345  << " to NaN" << endl;
346 
347  surfaceScalarField& phi = *iter();
348 
350 
351  continue;
352  }
353 
354  if (debug)
355  {
356  Pout<< "Mapping flux " << iter.key()
357  << " using interpolated flux " << UName
358  << endl;
359  }
360 
361  surfaceScalarField& phi = *iter();
362  const surfaceScalarField phiU
363  (
365  (
366  lookupObject<volVectorField>(UName)
367  )
368  & Sf()
369  );
370 
371  // Recalculate new internal faces.
372  for (label faceI = 0; faceI < nInternalFaces(); faceI++)
373  {
374  label oldFaceI = faceMap[faceI];
375 
376  if (oldFaceI == -1)
377  {
378  // Inflated/appended
379  phi[faceI] = phiU[faceI];
380  }
381  else if (reverseFaceMap[oldFaceI] != faceI)
382  {
383  // face-from-masterface
384  phi[faceI] = phiU[faceI];
385  }
386  }
387 
388  // Recalculate new boundary faces.
389  surfaceScalarField::GeometricBoundaryField& bphi =
390  phi.boundaryField();
391  forAll(bphi, patchI)
392  {
393  fvsPatchScalarField& patchPhi = bphi[patchI];
394  const fvsPatchScalarField& patchPhiU =
395  phiU.boundaryField()[patchI];
396 
397  label faceI = patchPhi.patch().start();
398 
399  forAll(patchPhi, i)
400  {
401  label oldFaceI = faceMap[faceI];
402 
403  if (oldFaceI == -1)
404  {
405  // Inflated/appended
406  patchPhi[i] = patchPhiU[i];
407  }
408  else if (reverseFaceMap[oldFaceI] != faceI)
409  {
410  // face-from-masterface
411  patchPhi[i] = patchPhiU[i];
412  }
413 
414  faceI++;
415  }
416  }
417 
418  // Update master faces
419  forAllConstIter(labelHashSet, masterFaces, iter)
420  {
421  label faceI = iter.key();
422 
423  if (isInternalFace(faceI))
424  {
425  phi[faceI] = phiU[faceI];
426  }
427  else
428  {
429  label patchI = boundaryMesh().whichPatch(faceI);
430  label i = faceI - boundaryMesh()[patchI].start();
431 
432  const fvsPatchScalarField& patchPhiU =
433  phiU.boundaryField()[patchI];
434 
435  fvsPatchScalarField& patchPhi = bphi[patchI];
436 
437  patchPhi[i] = patchPhiU[i];
438  }
439  }
440  }
441  }
442 
443 
444 
445  // Update numbering of cells/vertices.
446  meshCutter_.updateMesh(map);
447 
448  // Update numbering of protectedCell_
449  if (protectedCell_.size())
450  {
451  PackedBoolList newProtectedCell(nCells());
452 
453  forAll(newProtectedCell, cellI)
454  {
455  label oldCellI = map().cellMap()[cellI];
456  newProtectedCell.set(cellI, protectedCell_.get(oldCellI));
457  }
458  protectedCell_.transfer(newProtectedCell);
459  }
460 
461  // Debug: Check refinement levels (across faces only)
462  meshCutter_.checkRefinementLevels(-1, labelList(0));
463 
464  return map;
465 }
466 
467 
468 // Combines previously split cells, maps fields and recalculates
469 // (an approximate) flux
472 (
473  const labelList& splitPoints
474 )
475 {
476  polyTopoChange meshMod(*this);
477 
478  // Play refinement commands into mesh changer.
479  meshCutter_.setUnrefinement(splitPoints, meshMod);
480 
481 
482  // Save information on faces that will be combined
483  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
484 
485  // Find the faceMidPoints on cells to be combined.
486  // for each face resulting of split of face into four store the
487  // midpoint
488  Map<label> faceToSplitPoint(3*splitPoints.size());
489 
490  {
491  forAll(splitPoints, i)
492  {
493  label pointI = splitPoints[i];
494 
495  const labelList& pEdges = pointEdges()[pointI];
496 
497  forAll(pEdges, j)
498  {
499  label otherPointI = edges()[pEdges[j]].otherVertex(pointI);
500 
501  const labelList& pFaces = pointFaces()[otherPointI];
502 
503  forAll(pFaces, pFaceI)
504  {
505  faceToSplitPoint.insert(pFaces[pFaceI], otherPointI);
506  }
507  }
508  }
509  }
510 
511 
512  // Change mesh and generate map.
513  //autoPtr<mapPolyMesh> map = meshMod.changeMesh(*this, true);
514  autoPtr<mapPolyMesh> map = meshMod.changeMesh(*this, false);
515 
516  Info<< "Unrefined from "
517  << returnReduce(map().nOldCells(), sumOp<label>())
518  << " to " << globalData().nTotalCells() << " cells."
519  << endl;
520 
521  // Update fields
522  updateMesh(map);
523 
524 
525  // Move mesh
526  /*
527  pointField newPoints;
528  if (map().hasMotionPoints())
529  {
530  newPoints = map().preMotionPoints();
531  }
532  else
533  {
534  newPoints = points();
535  }
536  movePoints(newPoints);
537  */
538 
539  // Correct the flux for modified faces.
540  {
541  const labelList& reversePointMap = map().reversePointMap();
542  const labelList& reverseFaceMap = map().reverseFaceMap();
543 
545  (
546  lookupClass<surfaceScalarField>()
547  );
549  {
550  if (!correctFluxes_.found(iter.key()))
551  {
552  WarningIn("dynamicRefineFvMesh::refine(const labelList&)")
553  << "Cannot find surfaceScalarField " << iter.key()
554  << " in user-provided flux mapping table "
555  << correctFluxes_ << endl
556  << " The flux mapping table is used to recreate the"
557  << " flux on newly created faces." << endl
558  << " Either add the entry if it is a flux or use ("
559  << iter.key() << " none) to suppress this warning."
560  << endl;
561  continue;
562  }
563 
564  const word& UName = correctFluxes_[iter.key()];
565 
566  if (UName == "none")
567  {
568  continue;
569  }
570 
571  if (debug)
572  {
573  Info<< "Mapping flux " << iter.key()
574  << " using interpolated flux " << UName
575  << endl;
576  }
577 
578  surfaceScalarField& phi = *iter();
579  surfaceScalarField::GeometricBoundaryField& bphi =
580  phi.boundaryField();
581 
582  const surfaceScalarField phiU
583  (
585  (
586  lookupObject<volVectorField>(UName)
587  )
588  & Sf()
589  );
590 
591 
592  forAllConstIter(Map<label>, faceToSplitPoint, iter)
593  {
594  label oldFaceI = iter.key();
595  label oldPointI = iter();
596 
597  if (reversePointMap[oldPointI] < 0)
598  {
599  // midpoint was removed. See if face still exists.
600  label faceI = reverseFaceMap[oldFaceI];
601 
602  if (faceI >= 0)
603  {
604  if (isInternalFace(faceI))
605  {
606  phi[faceI] = phiU[faceI];
607  }
608  else
609  {
610  label patchI = boundaryMesh().whichPatch(faceI);
611  label i = faceI - boundaryMesh()[patchI].start();
612 
613  const fvsPatchScalarField& patchPhiU =
614  phiU.boundaryField()[patchI];
615  fvsPatchScalarField& patchPhi = bphi[patchI];
616  patchPhi[i] = patchPhiU[i];
617  }
618  }
619  }
620  }
621  }
622  }
623 
624 
625  // Update numbering of cells/vertices.
626  meshCutter_.updateMesh(map);
627 
628  // Update numbering of protectedCell_
629  if (protectedCell_.size())
630  {
631  PackedBoolList newProtectedCell(nCells());
632 
633  forAll(newProtectedCell, cellI)
634  {
635  label oldCellI = map().cellMap()[cellI];
636  if (oldCellI >= 0)
637  {
638  newProtectedCell.set(cellI, protectedCell_.get(oldCellI));
639  }
640  }
641  protectedCell_.transfer(newProtectedCell);
642  }
643 
644  // Debug: Check refinement levels (across faces only)
645  meshCutter_.checkRefinementLevels(-1, labelList(0));
646 
647  return map;
648 }
649 
650 
651 // Get max of connected point
654 {
655  scalarField vFld(nCells(), -GREAT);
656 
657  forAll(pointCells(), pointI)
658  {
659  const labelList& pCells = pointCells()[pointI];
660 
661  forAll(pCells, i)
662  {
663  vFld[pCells[i]] = max(vFld[pCells[i]], pFld[pointI]);
664  }
665  }
666  return vFld;
667 }
668 
669 
670 // Get max of connected cell
673 {
674  scalarField pFld(nPoints(), -GREAT);
675 
676  forAll(pointCells(), pointI)
677  {
678  const labelList& pCells = pointCells()[pointI];
679 
680  forAll(pCells, i)
681  {
682  pFld[pointI] = max(pFld[pointI], vFld[pCells[i]]);
683  }
684  }
685  return pFld;
686 }
687 
688 
689 // Simple (non-parallel) interpolation by averaging.
692 {
693  scalarField pFld(nPoints());
694 
695  forAll(pointCells(), pointI)
696  {
697  const labelList& pCells = pointCells()[pointI];
698 
699  scalar sum = 0.0;
700  forAll(pCells, i)
701  {
702  sum += vFld[pCells[i]];
703  }
704  pFld[pointI] = sum/pCells.size();
705  }
706  return pFld;
707 }
708 
709 
710 // Calculate error. Is < 0 or distance to minLevel, maxLevel
712 (
713  const scalarField& fld,
714  const scalar minLevel,
715  const scalar maxLevel
716 ) const
717 {
718  scalarField c(fld.size(), -1);
719 
720  forAll(fld, i)
721  {
722  scalar err = min(fld[i]-minLevel, maxLevel-fld[i]);
723 
724  if (err >= 0)
725  {
726  c[i] = err;
727  }
728  }
729  return c;
730 }
731 
732 
734 (
735  const scalar lowerRefineLevel,
736  const scalar upperRefineLevel,
737  const scalarField& vFld,
738  PackedBoolList& candidateCell
739 ) const
740 {
741  // Get error per cell. Is -1 (not to be refined) to >0 (to be refined,
742  // higher more desirable to be refined).
743  scalarField cellError
744  (
745  maxPointField
746  (
747  error
748  (
749  cellToPoint(vFld),
750  lowerRefineLevel,
751  upperRefineLevel
752  )
753  )
754  );
755 
756  // Mark cells that are candidates for refinement.
757  forAll(cellError, cellI)
758  {
759  if (cellError[cellI] > 0)
760  {
761  candidateCell.set(cellI, 1);
762  }
763  }
764 }
765 
766 
768 (
769  const label maxCells,
770  const label maxRefinement,
771  const PackedBoolList& candidateCell
772 ) const
773 {
774  // Every refined cell causes 7 extra cells
775  label nTotToRefine = (maxCells - globalData().nTotalCells()) / 7;
776 
777  const labelList& cellLevel = meshCutter_.cellLevel();
778 
779  // Mark cells that cannot be refined since they would trigger refinement
780  // of protected cells (since 2:1 cascade)
781  PackedBoolList unrefineableCell;
782  calculateProtectedCells(unrefineableCell);
783 
784  // Count current selection
785  label nLocalCandidates = count(candidateCell, 1);
786  label nCandidates = returnReduce(nLocalCandidates, sumOp<label>());
787 
788  // Collect all cells
789  DynamicList<label> candidates(nLocalCandidates);
790 
791  if (nCandidates < nTotToRefine)
792  {
793  forAll(candidateCell, cellI)
794  {
795  if
796  (
797  cellLevel[cellI] < maxRefinement
798  && candidateCell.get(cellI)
799  && (
800  unrefineableCell.empty()
801  || !unrefineableCell.get(cellI)
802  )
803  )
804  {
805  candidates.append(cellI);
806  }
807  }
808  }
809  else
810  {
811  // Sort by error? For now just truncate.
812  for (label level = 0; level < maxRefinement; level++)
813  {
814  forAll(candidateCell, cellI)
815  {
816  if
817  (
818  cellLevel[cellI] == level
819  && candidateCell.get(cellI)
820  && (
821  unrefineableCell.empty()
822  || !unrefineableCell.get(cellI)
823  )
824  )
825  {
826  candidates.append(cellI);
827  }
828  }
829 
830  if (returnReduce(candidates.size(), sumOp<label>()) > nTotToRefine)
831  {
832  break;
833  }
834  }
835  }
836 
837  // Guarantee 2:1 refinement after refinement
838  labelList consistentSet
839  (
840  meshCutter_.consistentRefinement
841  (
842  candidates.shrink(),
843  true // Add to set to guarantee 2:1
844  )
845  );
846 
847  Info<< "Selected " << returnReduce(consistentSet.size(), sumOp<label>())
848  << " cells for refinement out of " << globalData().nTotalCells()
849  << "." << endl;
850 
851  return consistentSet;
852 }
853 
854 
856 (
857  const scalar unrefineLevel,
858  const PackedBoolList& markedCell,
859  const scalarField& pFld
860 ) const
861 {
862  // All points that can be unrefined
863  const labelList splitPoints(meshCutter_.getSplitPoints());
864 
865  DynamicList<label> newSplitPoints(splitPoints.size());
866 
867  forAll(splitPoints, i)
868  {
869  label pointI = splitPoints[i];
870 
871  if (pFld[pointI] < unrefineLevel)
872  {
873  // Check that all cells are not marked
874  const labelList& pCells = pointCells()[pointI];
875 
876  bool hasMarked = false;
877 
878  forAll(pCells, pCellI)
879  {
880  if (markedCell.get(pCells[pCellI]))
881  {
882  hasMarked = true;
883  break;
884  }
885  }
886 
887  if (!hasMarked)
888  {
889  newSplitPoints.append(pointI);
890  }
891  }
892  }
893 
894 
895  newSplitPoints.shrink();
896 
897  // Guarantee 2:1 refinement after unrefinement
898  labelList consistentSet
899  (
900  meshCutter_.consistentUnrefinement
901  (
902  newSplitPoints,
903  false
904  )
905  );
906  Info<< "Selected " << returnReduce(consistentSet.size(), sumOp<label>())
907  << " split points out of a possible "
908  << returnReduce(splitPoints.size(), sumOp<label>())
909  << "." << endl;
910 
911  return consistentSet;
912 }
913 
914 
916 (
917  PackedBoolList& markedCell
918 ) const
919 {
920  // Mark faces using any marked cell
921  boolList markedFace(nFaces(), false);
922 
923  forAll(markedCell, cellI)
924  {
925  if (markedCell.get(cellI))
926  {
927  const cell& cFaces = cells()[cellI];
928 
929  forAll(cFaces, i)
930  {
931  markedFace[cFaces[i]] = true;
932  }
933  }
934  }
935 
936  syncTools::syncFaceList(*this, markedFace, orEqOp<bool>());
937 
938  // Update cells using any markedFace
939  for (label faceI = 0; faceI < nInternalFaces(); faceI++)
940  {
941  if (markedFace[faceI])
942  {
943  markedCell.set(faceOwner()[faceI], 1);
944  markedCell.set(faceNeighbour()[faceI], 1);
945  }
946  }
947  for (label faceI = nInternalFaces(); faceI < nFaces(); faceI++)
948  {
949  if (markedFace[faceI])
950  {
951  markedCell.set(faceOwner()[faceI], 1);
952  }
953  }
954 }
955 
956 
958 (
959  PackedBoolList& protectedCell,
960  label& nProtected
961 ) const
962 {
963  const labelList& cellLevel = meshCutter_.cellLevel();
964  const labelList& pointLevel = meshCutter_.pointLevel();
965 
966  labelList nAnchorPoints(nCells(), 0);
967 
968  forAll(pointLevel, pointI)
969  {
970  const labelList& pCells = pointCells(pointI);
971 
972  forAll(pCells, pCellI)
973  {
974  label cellI = pCells[pCellI];
975 
976  if (pointLevel[pointI] <= cellLevel[cellI])
977  {
978  // Check if cell has already 8 anchor points -> protect cell
979  if (nAnchorPoints[cellI] == 8)
980  {
981  if (protectedCell.set(cellI, true))
982  {
983  nProtected++;
984  }
985  }
986 
987  if (!protectedCell[cellI])
988  {
989  nAnchorPoints[cellI]++;
990  }
991  }
992  }
993  }
994 
995 
996  forAll(protectedCell, cellI)
997  {
998  if (!protectedCell[cellI] && nAnchorPoints[cellI] != 8)
999  {
1000  protectedCell.set(cellI, true);
1001  nProtected++;
1002  }
1003  }
1004 }
1005 
1006 
1007 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1008 
1009 Foam::dynamicRefineFvMesh::dynamicRefineFvMesh(const IOobject& io)
1010 :
1011  dynamicFvMesh(io),
1012  meshCutter_(*this),
1013  dumpLevel_(false),
1014  nRefinementIterations_(0),
1015  protectedCell_(nCells(), 0)
1016 {
1017  // Read static part of dictionary
1018  readDict();
1019 
1020 
1021  const labelList& cellLevel = meshCutter_.cellLevel();
1022  const labelList& pointLevel = meshCutter_.pointLevel();
1023 
1024  // Set cells that should not be refined.
1025  // This is currently any cell which does not have 8 anchor points or
1026  // uses any face which does not have 4 anchor points.
1027  // Note: do not use cellPoint addressing
1028 
1029  // Count number of points <= cellLevel
1030  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1031 
1032  labelList nAnchors(nCells(), 0);
1033 
1034  label nProtected = 0;
1035 
1036  forAll(pointCells(), pointI)
1037  {
1038  const labelList& pCells = pointCells()[pointI];
1039 
1040  forAll(pCells, i)
1041  {
1042  label cellI = pCells[i];
1043 
1044  if (!protectedCell_.get(cellI))
1045  {
1046  if (pointLevel[pointI] <= cellLevel[cellI])
1047  {
1048  nAnchors[cellI]++;
1049 
1050  if (nAnchors[cellI] > 8)
1051  {
1052  protectedCell_.set(cellI, 1);
1053  nProtected++;
1054  }
1055  }
1056  }
1057  }
1058  }
1059 
1060 
1061  // Count number of points <= faceLevel
1062  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1063  // Bit tricky since proc face might be one more refined than the owner since
1064  // the coupled one is refined.
1065 
1066  {
1067  labelList neiLevel(nFaces());
1068 
1069  for (label faceI = 0; faceI < nInternalFaces(); faceI++)
1070  {
1071  neiLevel[faceI] = cellLevel[faceNeighbour()[faceI]];
1072  }
1073  for (label faceI = nInternalFaces(); faceI < nFaces(); faceI++)
1074  {
1075  neiLevel[faceI] = cellLevel[faceOwner()[faceI]];
1076  }
1077  syncTools::swapFaceList(*this, neiLevel);
1078 
1079 
1080  boolList protectedFace(nFaces(), false);
1081 
1082  forAll(faceOwner(), faceI)
1083  {
1084  label faceLevel = max
1085  (
1086  cellLevel[faceOwner()[faceI]],
1087  neiLevel[faceI]
1088  );
1089 
1090  const face& f = faces()[faceI];
1091 
1092  label nAnchors = 0;
1093 
1094  forAll(f, fp)
1095  {
1096  if (pointLevel[f[fp]] <= faceLevel)
1097  {
1098  nAnchors++;
1099 
1100  if (nAnchors > 4)
1101  {
1102  protectedFace[faceI] = true;
1103  break;
1104  }
1105  }
1106  }
1107  }
1108 
1109  syncTools::syncFaceList(*this, protectedFace, orEqOp<bool>());
1110 
1111  for (label faceI = 0; faceI < nInternalFaces(); faceI++)
1112  {
1113  if (protectedFace[faceI])
1114  {
1115  protectedCell_.set(faceOwner()[faceI], 1);
1116  nProtected++;
1117  protectedCell_.set(faceNeighbour()[faceI], 1);
1118  nProtected++;
1119  }
1120  }
1121  for (label faceI = nInternalFaces(); faceI < nFaces(); faceI++)
1122  {
1123  if (protectedFace[faceI])
1124  {
1125  protectedCell_.set(faceOwner()[faceI], 1);
1126  nProtected++;
1127  }
1128  }
1129 
1130  // Also protect any cells that are less than hex
1131  forAll(cells(), cellI)
1132  {
1133  const cell& cFaces = cells()[cellI];
1134 
1135  if (cFaces.size() < 6)
1136  {
1137  if (protectedCell_.set(cellI, 1))
1138  {
1139  nProtected++;
1140  }
1141  }
1142  else
1143  {
1144  forAll(cFaces, cFaceI)
1145  {
1146  if (faces()[cFaces[cFaceI]].size() < 4)
1147  {
1148  if (protectedCell_.set(cellI, 1))
1149  {
1150  nProtected++;
1151  }
1152  break;
1153  }
1154  }
1155  }
1156  }
1157 
1158  // Check cells for 8 corner points
1160  }
1161 
1162  if (returnReduce(nProtected, sumOp<label>()) == 0)
1163  {
1165  }
1166  else
1167  {
1168 
1169  cellSet protectedCells(*this, "protectedCells", nProtected);
1170  forAll(protectedCell_, cellI)
1171  {
1172  if (protectedCell_[cellI])
1173  {
1174  protectedCells.insert(cellI);
1175  }
1176  }
1177 
1178  Info<< "Detected " << returnReduce(nProtected, sumOp<label>())
1179  << " cells that are protected from refinement."
1180  << " Writing these to cellSet "
1181  << protectedCells.name()
1182  << "." << endl;
1183 
1184  protectedCells.write();
1185  }
1186 }
1187 
1188 
1189 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1190 
1192 {}
1193 
1194 
1195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1196 
1198 {
1199  // Re-read dictionary. Choosen since usually -small so trivial amount
1200  // of time compared to actual refinement. Also very useful to be able
1201  // to modify on-the-fly.
1202  dictionary refineDict
1203  (
1204  IOdictionary
1205  (
1206  IOobject
1207  (
1208  "dynamicMeshDict",
1209  time().constant(),
1210  *this,
1213  false
1214  )
1215  ).subDict(typeName + "Coeffs")
1216  );
1217 
1218  label refineInterval = readLabel(refineDict.lookup("refineInterval"));
1219 
1220  bool hasChanged = false;
1221 
1222  if (refineInterval == 0)
1223  {
1224  topoChanging(hasChanged);
1225 
1226  return false;
1227  }
1228  else if (refineInterval < 0)
1229  {
1230  FatalErrorIn("dynamicRefineFvMesh::update()")
1231  << "Illegal refineInterval " << refineInterval << nl
1232  << "The refineInterval setting in the dynamicMeshDict should"
1233  << " be >= 1." << nl
1234  << exit(FatalError);
1235  }
1236 
1237 
1238 
1239 
1240  // Note: cannot refine at time 0 since no V0 present since mesh not
1241  // moved yet.
1242 
1243  if (time().timeIndex() > 0 && time().timeIndex() % refineInterval == 0)
1244  {
1245  label maxCells = readLabel(refineDict.lookup("maxCells"));
1246 
1247  if (maxCells <= 0)
1248  {
1249  FatalErrorIn("dynamicRefineFvMesh::update()")
1250  << "Illegal maximum number of cells " << maxCells << nl
1251  << "The maxCells setting in the dynamicMeshDict should"
1252  << " be > 0." << nl
1253  << exit(FatalError);
1254  }
1255 
1256  label maxRefinement = readLabel(refineDict.lookup("maxRefinement"));
1257 
1258  if (maxRefinement <= 0)
1259  {
1260  FatalErrorIn("dynamicRefineFvMesh::update()")
1261  << "Illegal maximum refinement level " << maxRefinement << nl
1262  << "The maxCells setting in the dynamicMeshDict should"
1263  << " be > 0." << nl
1264  << exit(FatalError);
1265  }
1266 
1267  const word fieldName(refineDict.lookup("field"));
1268 
1269  const volScalarField& vFld = lookupObject<volScalarField>(fieldName);
1270 
1271  const scalar lowerRefineLevel =
1272  readScalar(refineDict.lookup("lowerRefineLevel"));
1273  const scalar upperRefineLevel =
1274  readScalar(refineDict.lookup("upperRefineLevel"));
1275  const scalar unrefineLevel = refineDict.lookupOrDefault<scalar>
1276  (
1277  "unrefineLevel",
1278  GREAT
1279  );
1280  const label nBufferLayers =
1281  readLabel(refineDict.lookup("nBufferLayers"));
1282 
1283  // Cells marked for refinement or otherwise protected from unrefinement.
1285 
1286  // Determine candidates for refinement (looking at field only)
1288  (
1289  lowerRefineLevel,
1290  upperRefineLevel,
1291  vFld,
1292  refineCell
1293  );
1294 
1295  if (globalData().nTotalCells() < maxCells)
1296  {
1297  // Select subset of candidates. Take into account max allowable
1298  // cells, refinement level, protected cells.
1299  labelList cellsToRefine
1300  (
1302  (
1303  maxCells,
1304  maxRefinement,
1305  refineCell
1306  )
1307  );
1308 
1309  label nCellsToRefine = returnReduce
1310  (
1311  cellsToRefine.size(), sumOp<label>()
1312  );
1313 
1314  if (nCellsToRefine > 0)
1315  {
1316  // Refine/update mesh and map fields
1317  autoPtr<mapPolyMesh> map = refine(cellsToRefine);
1318 
1319  // Update refineCell. Note that some of the marked ones have
1320  // not been refined due to constraints.
1321  {
1322  const labelList& cellMap = map().cellMap();
1323  const labelList& reverseCellMap = map().reverseCellMap();
1324 
1325  PackedBoolList newRefineCell(cellMap.size());
1326 
1327  forAll(cellMap, cellI)
1328  {
1329  label oldCellI = cellMap[cellI];
1330 
1331  if (oldCellI < 0)
1332  {
1333  newRefineCell.set(cellI, 1);
1334  }
1335  else if (reverseCellMap[oldCellI] != cellI)
1336  {
1337  newRefineCell.set(cellI, 1);
1338  }
1339  else
1340  {
1341  newRefineCell.set(cellI, refineCell.get(oldCellI));
1342  }
1343  }
1344  refineCell.transfer(newRefineCell);
1345  }
1346 
1347  // Extend with a buffer layer to prevent neighbouring points
1348  // being unrefined.
1349  for (label i = 0; i < nBufferLayers; i++)
1350  {
1351  extendMarkedCells(refineCell);
1352  }
1353 
1354  hasChanged = true;
1355  }
1356  }
1357 
1358 
1359  {
1360  // Select unrefineable points that are not marked in refineCell
1361  labelList pointsToUnrefine
1362  (
1364  (
1365  unrefineLevel,
1366  refineCell,
1367  maxCellField(vFld)
1368  )
1369  );
1370 
1371  label nSplitPoints = returnReduce
1372  (
1373  pointsToUnrefine.size(),
1374  sumOp<label>()
1375  );
1376 
1377  if (nSplitPoints > 0)
1378  {
1379  // Refine/update mesh
1380  unrefine(pointsToUnrefine);
1381 
1382  hasChanged = true;
1383  }
1384  }
1385 
1386 
1387  if ((nRefinementIterations_ % 10) == 0)
1388  {
1389  // Compact refinement history occassionally (how often?).
1390  // Unrefinement causes holes in the refinementHistory.
1391  const_cast<refinementHistory&>(meshCutter().history()).compact();
1392  }
1394  }
1395 
1396  topoChanging(hasChanged);
1397  if (hasChanged)
1398  {
1399  // Reset moving flag (if any). If not using inflation we'll not move,
1400  // if are using inflation any follow on movePoints will set it.
1401  moving(false);
1402  }
1403 
1404  return hasChanged;
1405 }
1406 
1407 
1413 ) const
1414 {
1415  // Force refinement data to go to the current time directory.
1416  const_cast<hexRef8&>(meshCutter_).setInstance(time().timeName());
1417 
1418  bool writeOk =
1419  (
1420  dynamicFvMesh::writeObjects(fmt, ver, cmp)
1421  && meshCutter_.write()
1422  );
1423 
1424  if (dumpLevel_)
1425  {
1426  volScalarField scalarCellLevel
1427  (
1428  IOobject
1429  (
1430  "cellLevel",
1431  time().timeName(),
1432  *this,
1435  false
1436  ),
1437  *this,
1438  dimensionedScalar("level", dimless, 0)
1439  );
1440 
1441  const labelList& cellLevel = meshCutter_.cellLevel();
1442 
1443  forAll(cellLevel, cellI)
1444  {
1445  scalarCellLevel[cellI] = cellLevel[cellI];
1446  }
1447 
1448  writeOk = writeOk && scalarCellLevel.write();
1449  }
1450 
1451  return writeOk;
1452 }
1453 
1454 
1455 // ************************************************************************* //
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:52
Foam::surfaceFields.
const refinementHistory & history() const
Definition: hexRef8.H:399
label nFaces() const
PackedBoolList protectedCell_
Protected cells (usually since not hexes)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
virtual bool writeObjects(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:864
autoPtr< mapPolyMesh > refine(const labelList &)
Refine cells. Update mesh and fields.
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:390
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
labelList f(nPoints)
scalarField maxCellField(const volScalarField &) const
Get point max of connected cell.
const fvPatch & patch() const
Return patch.
void extendMarkedCells(PackedBoolList &markedCell) const
Extend markedCell with cell-face-cell.
An STL-conforming hash table.
Definition: HashTable.H:61
void set(const PackedList< 1 > &)
Set specified bits.
#define forAllIter(Container, container, iter)
Definition: UList.H:440
A collection of cell labels.
Definition: cellSet.H:48
label timeIndex
Definition: getTimeIndex.H:4
A bit-packed bool list.
A class for handling words, derived from string.
Definition: word.H:59
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:59
label nRefinementIterations_
Number of refinement/unrefinement steps done so far.
bool topoChanging() const
Is mesh topology changing.
Definition: polyMesh.H:507
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:310
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
InternalField & internalField()
Return internal field.
label nCells() const
const cellList & cells() const
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:64
Surface Interpolation.
const labelListList & pointCells() const
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
messageStream Info
All refinement history. Used in unrefinement.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
virtual void selectRefineCandidates(const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, PackedBoolList &candidateCell) const
Select candidate cells for refinement.
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:380
Namespace for OpenFOAM.
label readLabel(Istream &is)
Definition: label.H:64
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
void calculateProtectedCells(PackedBoolList &unrefineableCell) const
Calculate cells that cannot be refined since would trigger.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
virtual bool write() const
Write using setting from DB.
const hexRef8 & meshCutter() const
Direct access to the refinement engine.
label n
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
scalarField cellToPoint(const scalarField &vFld) const
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp) const
Write using given format, version and compression.
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const labelIOList & cellLevel() const
Definition: hexRef8.H:389
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
#define WarningIn(functionName)
Report a warning using Foam::Warning.
const labelIOList & pointLevel() const
Definition: hexRef8.H:394
bool write() const
Force writing refinement+history to polyMesh directory.
Definition: hexRef8.C:5832
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
Switch dumpLevel_
Dump cellLevel for postprocessing.
static label count(const PackedBoolList &, const unsigned int)
Count set/unset elements in packedlist.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
static void fillNan(UList< scalar > &)
Fill block of data with NaN.
Definition: sigFpe.C:48
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:258
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:429
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
label size() const
Return number of elements in table.
void readDict()
Read the projection parameters from dictionary.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const word & name() const
Return name.
Definition: IOobject.H:260
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1073
virtual labelList selectRefineCells(const label maxCells, const label maxRefinement, const PackedBoolList &candidateCell) const
Subset candidate cells for refinement.
surfaceScalarField & phi
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
const cellShapeList & cells
label nInternalFaces() const
virtual ~dynamicRefineFvMesh()
Destructor.
A topoSetSource to select points based on usage in cells.
Definition: cellToPoint.H:49
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:50
hexRef8 meshCutter_
Mesh cutting engine.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:66
void clear()
Clear the list, i.e. set addressable size to zero.
Definition: PackedListI.H:900
error FatalError
Constant dispersed-phase particle diameter model.
scalarField maxPointField(const scalarField &) const
Get per cell max of connected point.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
List< label > labelList
A List of labels.
Definition: labelList.H:56
const dimensionedScalar c
Speed of light in a vacuum.
Direct mesh changes based on v1.3 polyTopoChange syntax.
virtual bool update()
Update the mesh for both mesh motion and topology change.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1060
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
Abstract base class for geometry and/or topology changing fvMesh.
Definition: dynamicFvMesh.H:51
scalarField error(const scalarField &fld, const scalar minLevel, const scalar maxLevel) const
void checkEightAnchorPoints(PackedBoolList &protectedCell, label &nProtected) const
Check all cells have 8 anchor points.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Version number type.
Definition: IOstream.H:96
label nPoints
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1079
static void swapFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled face values.
Definition: syncTools.H:462
autoPtr< mapPolyMesh > unrefine(const labelList &)
Unrefine cells. Gets passed in centre points of cells to combine.
bool empty() const
Return true if the list is empty (ie, size() is zero).
Definition: PackedListI.H:729
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic 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:117
virtual labelList selectUnrefinePoints(const scalar unrefineLevel, const PackedBoolList &markedCell, const scalarField &pFld) const
Select points that can be unrefined.
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.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1200
defineTypeNameAndDebug(combustionModel, 0)
word timeName
Definition: getTimeIndex.H:3
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
unsigned int get(const label) const
Get value at index I.
Definition: PackedListI.H:972
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116