24#ifndef DUMUX_GRIDADAPTIONINDICATOR2PLOCAL_HH
25#define DUMUX_GRIDADAPTIONINDICATOR2PLOCAL_HH
38template<
class TypeTag>
45 using Element =
typename GridView::Traits::template Codim<0>::Entity;
47 using SolutionTypes =
typename GET_PROP(TypeTag, SolutionTypes);
48 using ScalarSolutionType =
typename SolutionTypes::ScalarSolution;
50 using PrimaryVariables =
typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
51 using BoundaryTypes =
typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
53 using Indices =
typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
57 sw = Indices::saturationW,
58 sn = Indices::saturationNw
62 wPhaseIdx = Indices::wPhaseIdx,
63 nPhaseIdx = Indices::nPhaseIdx
75 if(indicatorVector_.size() != problem_.variables().cellDataGlobal().size())
77 indicatorVector_.resize(problem_.variables().cellDataGlobal().size());
79 indicatorVector_ = -1e100;
81 Scalar globalMax = -1e100;
82 Scalar globalMin = 1e100;
84 Scalar maxLocalDelta = 0;
88 for (
const auto& element : elements(problem_.gridView()))
92 int globalIdxI = problem_.variables().index(element);
96 PrimaryVariables source(0.0);
97 problem_.source(source, element);
98 for (
int i = 0; i < 2; i++)
101 if (abs(source[i]) > 1e-10)
103 indicatorVector_[globalIdxI] = 10;
109 if (indicatorVector_[globalIdxI] == 10)
115 switch (saturationType_)
118 satI = problem_.variables().cellData(globalIdxI).saturation(wPhaseIdx);
121 satI = problem_.variables().cellData(globalIdxI).saturation(nPhaseIdx);
126 globalMin = min(satI, globalMin);
128 globalMax = max(satI, globalMax);
131 for (
const auto& intersection :
intersections(problem_.gridView(), element))
133 if (indicatorVector_[globalIdxI] == 10)
139 if (intersection.boundary())
141 BoundaryTypes bcTypes;
142 problem_.boundaryTypes(bcTypes, intersection);
144 for (
int i = 0; i < 2; i++)
146 if (bcTypes.isNeumann(i))
148 PrimaryVariables flux(0.0);
149 problem_.neumann(flux, intersection);
151 bool fluxBound =
false;
152 for (
int j = 0; j < 2; j++)
155 if (abs(flux[j]) > 1e-10)
159 indicatorVector_[globalIdxI] = 10;
168 else if (bcTypes.isDirichlet(i))
170 if (refineAtDirichletBC_)
172 indicatorVector_[globalIdxI] = 10;
181 auto outside = intersection.outside();
182 int globalIdxJ = problem_.variables().index(outside);
185 if (element.level() > outside.level()
186 || (element.level() == outside.level() && globalIdxI < globalIdxJ))
189 switch (saturationType_)
192 satJ = problem_.variables().cellData(globalIdxJ).saturation(wPhaseIdx);
195 satJ = problem_.variables().cellData(globalIdxJ).saturation(nPhaseIdx);
200 Scalar localdelta = abs(satI - satJ);
202 indicatorVector_[globalIdxI][0] = max(indicatorVector_[globalIdxI][0], localdelta);
203 indicatorVector_[globalIdxJ][0] = max(indicatorVector_[globalIdxJ][0], localdelta);
205 maxLocalDelta = max(maxLocalDelta, localdelta);
211 Scalar globaldelta = globalMax - globalMin;
213 if (maxLocalDelta > 0.)
215 indicatorVector_ /= maxLocalDelta;
217 if (globaldelta > 0.)
219 maxLocalDelta /= globaldelta;
222 refineBound_ = refinetol_*maxLocalDelta;
223 coarsenBound_ = coarsentol_*maxLocalDelta;
235 return (indicatorVector_[problem_.elementMapper().index(element)] > refineBound_);
247 return (indicatorVector_[problem_.elementMapper().index(element)] < coarsenBound_);
266 indicatorVector_[idx] = coarsenBound_+1;
276 indicatorVector_[idx] = refineBound_ - 1;
292 refinetol_ = getParam<Scalar>(
"GridAdapt.RefineTolerance");
293 coarsentol_ = getParam<Scalar>(
"GridAdapt.CoarsenTolerance");
294 refineAtDirichletBC_ = getParam<bool>(
"GridAdapt.RefineAtDirichletBC");
295 refineAtFluxBC_ = getParam<bool>(
"GridAdapt.RefineAtFluxBC");
296 refineAtSource_ = getParam<bool>(
"GridAdapt.RefineAtSource");
304 Scalar coarsenBound_;
305 ScalarSolutionType indicatorVector_;
307 bool refineAtDirichletBC_;
308 bool refineAtFluxBC_;
309 bool refineAtSource_;
#define GET_PROP_VALUE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:282
#define GET_PROP(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:281
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
Base file for properties related to sequential IMPET algorithms.
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag SaturationFormulation
The formulation of the saturation model.
Definition: porousmediumflow/2p/sequential/properties.hh:53
Class defining a standard, saturation dependent indicator for grid adaption.
Definition: gridadaptionindicatorlocal.hh:40
void setIndicatorToCoarse(int idx)
Function for changing the indicatorVector values for coarsening.
Definition: gridadaptionindicatorlocal.hh:274
bool refine(const Element &element)
Indicator function for marking of grid cells for refinement.
Definition: gridadaptionindicatorlocal.hh:233
void calculateIndicator()
Calculates the indicator used for refinement/coarsening for each grid cell.
Definition: gridadaptionindicatorlocal.hh:72
void init()
Initializes the adaption indicator class.
Definition: gridadaptionindicatorlocal.hh:253
GridAdaptionIndicator2PLocal(Problem &problem)
Constructs a GridAdaptionIndicator instance.
Definition: gridadaptionindicatorlocal.hh:289
bool coarsen(const Element &element)
Indicator function for marking of grid cells for coarsening.
Definition: gridadaptionindicatorlocal.hh:245
void setIndicatorToRefine(int idx)
Function for changing the indicatorVector values for refinement.
Definition: gridadaptionindicatorlocal.hh:264
Defines the properties required for (immiscible) two-phase sequential models.