13#ifndef DUMUX_TWOP_ADAPTION_INDICATOR_HH
14#define DUMUX_TWOP_ADAPTION_INDICATOR_HH
17#include <dune/common/exceptions.hh>
18#include <dune/grid/common/partitionset.hh>
31template<
class TypeTag>
36 using Element =
typename GridView::template Codim<0>::Entity;
41 enum { saturationIdx = Indices::saturationIdx };
55 : gridGeometry_(gridGeometry)
56 , refineBound_(std::numeric_limits<Scalar>::max())
57 , coarsenBound_(std::numeric_limits<Scalar>::lowest())
58 , maxSaturationDelta_(gridGeometry_->gridView().size(0), 0.0)
89 void setLevels(std::size_t minLevel, std::size_t maxLevel)
105 Scalar refineTol = 0.05,
106 Scalar coarsenTol = 0.001)
109 refineBound_ = std::numeric_limits<Scalar>::max();
110 coarsenBound_ = std::numeric_limits<Scalar>::lowest();
111 maxSaturationDelta_.assign(gridGeometry_->gridView().size(0), 0.0);
114 if (minLevel_ >= maxLevel_)
118 if (coarsenTol > refineTol)
119 DUNE_THROW(Dune::InvalidStateException,
"Refine tolerance must be higher than coarsen tolerance");
122 Scalar globalMax = std::numeric_limits<Scalar>::lowest();
123 Scalar globalMin = std::numeric_limits<Scalar>::max();
126 for (
const auto& element : elements(gridGeometry_->gridView()))
129 const auto globalIdxI = gridGeometry_->elementMapper().index(element);
132 const auto geometry = element.geometry();
134 const Scalar satI =
evalSolution(element, geometry, *gridGeometry_, elemSol, geometry.center())[saturationIdx];
139 globalMin = min(satI, globalMin);
140 globalMax = max(satI, globalMax);
143 for (
const auto& intersection : intersections(gridGeometry_->gridView(), element))
146 if (intersection.neighbor())
149 const auto outside = intersection.outside();
150 const auto globalIdxJ = gridGeometry_->elementMapper().index(outside);
153 if (element.level() > outside.level() || (element.level() == outside.level() && globalIdxI < globalIdxJ))
156 const auto outsideGeometry = outside.geometry();
158 const Scalar satJ =
evalSolution(outside, outsideGeometry, *gridGeometry_, elemSolJ, outsideGeometry.center())[saturationIdx];
161 Scalar localdelta = abs(satI - satJ);
162 maxSaturationDelta_[globalIdxI] = max(maxSaturationDelta_[globalIdxI], localdelta);
163 maxSaturationDelta_[globalIdxJ] = max(maxSaturationDelta_[globalIdxJ], localdelta);
170 const auto globalDelta = globalMax - globalMin;
173 refineBound_ = refineTol*globalDelta;
174 coarsenBound_ = coarsenTol*globalDelta;
192 for (
const auto& element : elements(gridGeometry_->gridView(), Dune::Partitions::interior))
193 if (this->
operator()(element) > 0)
194 checkNeighborsRefine_(element);
208 if (element.hasFather()
209 && maxSaturationDelta_[gridGeometry_->elementMapper().index(element)] < coarsenBound_)
213 else if (element.level() < maxLevel_
214 && maxSaturationDelta_[gridGeometry_->elementMapper().index(element)] > refineBound_)
234 bool checkNeighborsRefine_(
const Element &element, std::size_t level = 1)
236 for(
const auto& intersection : intersections(gridGeometry_->gridView(), element))
238 if(!intersection.neighbor())
242 const auto outside = intersection.outside();
245 if (outside.partitionType() == Dune::GhostEntity)
248 if (outside.level() < maxLevel_ && outside.level() < element.level())
251 maxSaturationDelta_[gridGeometry_->elementMapper().index(outside)] = std::numeric_limits<Scalar>::max();
252 if(level < maxLevel_)
253 checkNeighborsRefine_(outside, ++level);
260 std::shared_ptr<const GridGeometry> gridGeometry_;
263 Scalar coarsenBound_;
264 std::vector< Scalar > maxSaturationDelta_;
265 std::size_t minLevel_;
266 std::size_t maxLevel_;
Class defining a standard, saturation dependent indicator for grid adaptation.
Definition: gridadaptindicator.hh:33
void setMinLevel(std::size_t minLevel)
Function to set the minimum allowed level.
Definition: gridadaptindicator.hh:68
int operator()(const Element &element) const
function call operator to return mark
Definition: gridadaptindicator.hh:206
void setMaxLevel(std::size_t maxLevel)
Function to set the maximum allowed level.
Definition: gridadaptindicator.hh:78
TwoPGridAdaptIndicator(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
The Constructor.
Definition: gridadaptindicator.hh:54
void calculate(const SolutionVector &sol, Scalar refineTol=0.05, Scalar coarsenTol=0.001)
Calculates the indicator used for refinement/coarsening for each grid cell.
Definition: gridadaptindicator.hh:104
void setLevels(std::size_t minLevel, std::size_t maxLevel)
Function to set the minimum/maximum allowed levels.
Definition: gridadaptindicator.hh:89
Defines all properties used in Dumux.
Element solution classes and factory functions.
free functions for the evaluation of primary variables inside elements.
PrimaryVariables evalSolution(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const CVFEElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Interpolates a given box element solution at a given global position. Uses the finite element cache o...
Definition: evalsolution.hh:152
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition: parameters.hh:149
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.