24#ifndef DUMUX_GRIDADAPTINITIALIZATIONINDICATOR_HH
25#define DUMUX_GRIDADAPTINITIALIZATIONINDICATOR_HH
29#include <dune/geometry/type.hh>
44template<
class TypeTag>
50 using Element =
typename GridView::Traits::template Codim<0>::Entity;
73 std::shared_ptr<const GridGeometry> gridGeometry,
74 std::shared_ptr<const GridVariables> gridVariables)
76 , gridGeometry_(gridGeometry)
77 , gridVariables_(gridVariables)
80 , refineAtDirichletBC_(
getParamFromGroup<bool>(problem->paramGroup(),
"Adaptive.RefineAtDirichletBC", true))
81 , refineAtFluxBC_(
getParamFromGroup<bool>(problem->paramGroup(),
"Adaptive.RefineAtFluxBC", true))
82 , refineAtSource_(
getParamFromGroup<bool>(problem->paramGroup(),
"Adaptive.RefineAtSource", true))
83 , eps_(
getParamFromGroup<Scalar>(problem->paramGroup(),
"Adaptive.BCRefinementThreshold", 1e-10))
107 minLevel_ = minLevel;
116 refineAtDirichletBC_ =
refine;
140 template<
class SolutionVector>
144 indicatorVector_.assign(gridGeometry_->gridView().size(0),
false);
146 for (
const auto& element : elements(gridGeometry_->gridView()))
148 const auto eIdx = gridGeometry_->elementMapper().index(element);
151 if (element.level() < minLevel_)
153 indicatorVector_[eIdx] =
true;
158 if (!refineAtSource_ && !refineAtFluxBC_ && !refineAtDirichletBC_)
162 if (element.level() == maxLevel_)
166 auto fvGeometry =
localView(*gridGeometry_);
167 fvGeometry.bind(element);
169 auto elemVolVars =
localView(gridVariables_->curGridVolVars());
170 elemVolVars.bind(element, fvGeometry, sol);
173 auto elemFluxVarsCache =
localView(gridVariables_->gridFluxVarsCache());
174 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
179 for (
const auto& scv : scvs(fvGeometry))
181 auto source = problem_->source(element, fvGeometry, elemVolVars, scv);
182 auto pointSource = problem_->scvPointSources(element, fvGeometry, elemVolVars, scv);
183 if (source.infinity_norm() + pointSource.infinity_norm() > eps_)
185 indicatorVector_[eIdx] =
true;
192 if (!indicatorVector_[eIdx]
193 && element.hasBoundaryIntersections()
194 && (refineAtDirichletBC_ || refineAtFluxBC_))
199 for (
const auto& scvf : scvfs(fvGeometry))
202 if (!scvf.boundary())
205 const auto bcTypes = problem_->boundaryTypes(element, scvf);
207 if(bcTypes.hasOnlyDirichlet() && refineAtDirichletBC_)
209 indicatorVector_[eIdx] =
true;
214 else if(refineAtFluxBC_)
216 const auto fluxes = problem_->neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
217 if (fluxes.infinity_norm() > eps_)
219 indicatorVector_[eIdx] =
true;
230 std::vector<BoundaryTypes> bcTypes(fvGeometry.numScv());
233 for (
const auto& scv : scvs(fvGeometry))
235 bcTypes[scv.localDofIndex()] = problem_->boundaryTypes(element, scv);
236 if (refineAtDirichletBC_ && bcTypes[scv.localDofIndex()].hasDirichlet())
238 indicatorVector_[eIdx] =
true;
244 if (!indicatorVector_[eIdx] && refineAtFluxBC_)
246 for (
const auto& scvf : scvfs(fvGeometry))
249 if (scvf.boundary() && bcTypes[scvf.insideScvIdx()].hasNeumann())
251 const auto fluxes = problem_->neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
252 if (fluxes.infinity_norm() > eps_)
254 indicatorVector_[eIdx] =
true;
278 if (indicatorVector_[gridGeometry_->elementMapper().index(element)])
284 std::shared_ptr<const Problem> problem_;
285 std::shared_ptr<const GridGeometry> gridGeometry_;
286 std::shared_ptr<const GridVariables> gridVariables_;
287 std::vector<bool> indicatorVector_;
291 bool refineAtDirichletBC_;
292 bool refineAtFluxBC_;
293 bool refineAtSource_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition: parameters.hh:374
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Class defining an initialization indicator for grid adaption. Refines at sources and boundaries....
Definition: initializationindicator.hh:46
void setRefinementAtFluxBoundaries(bool refine)
Function to set if refinement at sources is to be used.
Definition: initializationindicator.hh:130
void setLevels(std::size_t minLevel, std::size_t maxLevel)
Function to set the minumum/maximum allowed levels.
Definition: initializationindicator.hh:105
void setMaxLevel(std::size_t maxLevel)
Function to set the maximum allowed level.
Definition: initializationindicator.hh:97
void calculate(const SolutionVector &sol)
Calculates the indicator used for refinement/coarsening for each grid cell.
Definition: initializationindicator.hh:141
void setRefinementAtDirichletBC(bool refine)
Function to set if refinement at Dirichlet boundaries is to be used.
Definition: initializationindicator.hh:114
int operator()(const Element &element) const
function call operator to return mark
Definition: initializationindicator.hh:276
GridAdaptInitializationIndicator(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< const GridVariables > gridVariables)
Constructor.
Definition: initializationindicator.hh:72
bool refine(const Element &element)
Indicator function for marking of grid cells for refinement.
Definition: gridadaptinitializationindicator.hh:316
int maxLevel()
Definition: gridadaptinitializationindicator.hh:346
void setRefinementAtSources(bool refine)
Function to set if refinement at sources is to be used.
Definition: initializationindicator.hh:122
void setMinLevel(std::size_t minLevel)
Function to set the minimum allowed level.
Definition: initializationindicator.hh:89
typename Detail::template ProblemTraits< Problem, GridGeometry::discMethod >::BoundaryTypes BoundaryTypes
Definition: common/typetraits/problem.hh:46
Declares all properties used in Dumux.
Type traits for problem classes.