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