25#ifndef DUMUX_TWOP_ADAPTION_INDICATOR_HH
26#define DUMUX_TWOP_ADAPTION_INDICATOR_HH
29#include <dune/common/exceptions.hh>
30#include <dune/grid/common/partitionset.hh>
43template<
class TypeTag>
48 using Element =
typename GridView::template Codim<0>::Entity;
53 enum { saturationIdx = Indices::saturationIdx };
67 : gridGeometry_(gridGeometry)
68 , refineBound_(std::numeric_limits<Scalar>::max())
69 , coarsenBound_(std::numeric_limits<Scalar>::lowest())
70 , maxSaturationDelta_(gridGeometry_->gridView().size(0), 0.0)
101 void setLevels(std::size_t minLevel, std::size_t maxLevel)
103 minLevel_ = minLevel;
104 maxLevel_ = maxLevel;
117 Scalar refineTol = 0.05,
118 Scalar coarsenTol = 0.001)
121 refineBound_ = std::numeric_limits<Scalar>::max();
122 coarsenBound_ = std::numeric_limits<Scalar>::lowest();
123 maxSaturationDelta_.assign(gridGeometry_->gridView().size(0), 0.0);
126 if (minLevel_ >= maxLevel_)
130 if (coarsenTol > refineTol)
131 DUNE_THROW(Dune::InvalidStateException,
"Refine tolerance must be higher than coarsen tolerance");
134 Scalar globalMax = std::numeric_limits<Scalar>::lowest();
135 Scalar globalMin = std::numeric_limits<Scalar>::max();
138 for (
const auto& element : elements(gridGeometry_->gridView()))
141 const auto globalIdxI = gridGeometry_->elementMapper().index(element);
144 const auto geometry = element.geometry();
146 const Scalar satI =
evalSolution(element, geometry, *gridGeometry_, elemSol, geometry.center())[saturationIdx];
151 globalMin = min(satI, globalMin);
152 globalMax = max(satI, globalMax);
155 for (
const auto& intersection : intersections(gridGeometry_->gridView(), element))
158 if (intersection.neighbor())
161 const auto outside = intersection.outside();
162 const auto globalIdxJ = gridGeometry_->elementMapper().index(outside);
165 if (element.level() > outside.level() || (element.level() == outside.level() && globalIdxI < globalIdxJ))
168 const auto outsideGeometry = outside.geometry();
170 const Scalar satJ =
evalSolution(outside, outsideGeometry, *gridGeometry_, elemSolJ, outsideGeometry.center())[saturationIdx];
173 Scalar localdelta = abs(satI - satJ);
174 maxSaturationDelta_[globalIdxI] = max(maxSaturationDelta_[globalIdxI], localdelta);
175 maxSaturationDelta_[globalIdxJ] = max(maxSaturationDelta_[globalIdxJ], localdelta);
182 const auto globalDelta = globalMax - globalMin;
185 refineBound_ = refineTol*globalDelta;
186 coarsenBound_ = coarsenTol*globalDelta;
204 for (
const auto& element : elements(gridGeometry_->gridView(), Dune::Partitions::interior))
205 if (this->
operator()(element) > 0)
206 checkNeighborsRefine_(element);
220 if (element.hasFather()
221 && maxSaturationDelta_[gridGeometry_->elementMapper().index(element)] < coarsenBound_)
225 else if (element.level() < maxLevel_
226 && maxSaturationDelta_[gridGeometry_->elementMapper().index(element)] > refineBound_)
246 bool checkNeighborsRefine_(
const Element &element, std::size_t level = 1)
248 for(
const auto& intersection : intersections(gridGeometry_->gridView(), element))
250 if(!intersection.neighbor())
254 const auto outside = intersection.outside();
257 if (outside.partitionType() == Dune::GhostEntity)
260 if (outside.level() < maxLevel_ && outside.level() < element.level())
263 maxSaturationDelta_[gridGeometry_->elementMapper().index(outside)] = std::numeric_limits<Scalar>::max();
264 if(level < maxLevel_)
265 checkNeighborsRefine_(outside, ++level);
272 std::shared_ptr<const GridGeometry> gridGeometry_;
275 Scalar coarsenBound_;
276 std::vector< Scalar > maxSaturationDelta_;
277 std::size_t minLevel_;
278 std::size_t maxLevel_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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 BoxElementSolution< 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:95
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition: parameters.hh:161
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Class defining a standard, saturation dependent indicator for grid adaptation.
Definition: gridadaptindicator.hh:45
void setMinLevel(std::size_t minLevel)
Function to set the minimum allowed level.
Definition: gridadaptindicator.hh:80
int operator()(const Element &element) const
function call operator to return mark
Definition: gridadaptindicator.hh:218
void setMaxLevel(std::size_t maxLevel)
Function to set the maximum allowed level.
Definition: gridadaptindicator.hh:90
TwoPGridAdaptIndicator(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
The Constructor.
Definition: gridadaptindicator.hh:66
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:116
void setLevels(std::size_t minLevel, std::size_t maxLevel)
Function to set the minumum/maximum allowed levels.
Definition: gridadaptindicator.hh:101
Declares all properties used in Dumux.