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)
78 , minLevel_(
getParamFromGroup<std::size_t>(problem->paramGroup(),
"Adaptive.MinLevel"))
79 , maxLevel_(
getParamFromGroup<std::size_t>(problem->paramGroup(),
"Adaptive.MaxLevel"))
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))
105 void setLevels(std::size_t minLevel, std::size_t maxLevel)
107 minLevel_ = minLevel;
108 maxLevel_ = maxLevel;
116 refineAtDirichletBC_ = refine;
124 refineAtSource_ = refine;
132 refineAtFluxBC_ = refine;
140 template<
class SolutionVector>
144 indicatorVector_.assign(gridGeometry_->gridView().size(0),
false);
147 auto fvGeometry =
localView(*gridGeometry_);
148 auto elemVolVars =
localView(gridVariables_->curGridVolVars());
150 auto elemFluxVarsCache =
localView(gridVariables_->gridFluxVarsCache());
152 for (
const auto& element : elements(gridGeometry_->gridView()))
154 const auto eIdx = gridGeometry_->elementMapper().index(element);
157 if (element.level() < minLevel_)
159 indicatorVector_[eIdx] =
true;
164 if (!refineAtSource_ && !refineAtFluxBC_ && !refineAtDirichletBC_)
168 if (element.level() == maxLevel_)
172 fvGeometry.bind(element);
173 elemVolVars.bind(element, fvGeometry, sol);
174 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
180 for (
const auto& scv : scvs(fvGeometry))
182 auto source = problem_->source(element, fvGeometry, elemVolVars, scv);
183 auto pointSource = problem_->scvPointSources(element, fvGeometry, elemVolVars, scv);
184 if (source.infinity_norm() + pointSource.infinity_norm() > eps_)
186 indicatorVector_[eIdx] =
true;
193 if (!indicatorVector_[eIdx]
194 && element.hasBoundaryIntersections()
195 && (refineAtDirichletBC_ || refineAtFluxBC_))
200 for (
const auto& scvf : scvfs(fvGeometry))
203 if (!scvf.boundary())
206 const auto bcTypes = problem_->boundaryTypes(element, scvf);
208 if(bcTypes.hasOnlyDirichlet() && refineAtDirichletBC_)
210 indicatorVector_[eIdx] =
true;
215 else if(refineAtFluxBC_)
217 const auto fluxes = problem_->neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
218 if (fluxes.infinity_norm() > eps_)
220 indicatorVector_[eIdx] =
true;
231 std::vector<BoundaryTypes> bcTypes(fvGeometry.numScv());
234 for (
const auto& scv : scvs(fvGeometry))
236 bcTypes[scv.localDofIndex()] = problem_->boundaryTypes(element, scv);
237 if (refineAtDirichletBC_ && bcTypes[scv.localDofIndex()].hasDirichlet())
239 indicatorVector_[eIdx] =
true;
245 if (!indicatorVector_[eIdx] && refineAtFluxBC_)
247 for (
const auto& scvf : scvfs(fvGeometry))
250 if (scvf.boundary() && bcTypes[scvf.insideScvIdx()].hasNeumann())
252 const auto fluxes = problem_->neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
253 if (fluxes.infinity_norm() > eps_)
255 indicatorVector_[eIdx] =
true;
279 if (indicatorVector_[gridGeometry_->elementMapper().index(element)])
285 std::shared_ptr<const Problem> problem_;
286 std::shared_ptr<const GridGeometry> gridGeometry_;
287 std::shared_ptr<const GridVariables> gridVariables_;
288 std::vector<bool> indicatorVector_;
292 bool refineAtDirichletBC_;
293 bool refineAtFluxBC_;
294 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:161
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr Box box
Definition: method.hh:136
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 minimum/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:277
GridAdaptInitializationIndicator(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< const GridVariables > gridVariables)
Constructor.
Definition: initializationindicator.hh:72
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
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:38
typename Detail::template ProblemTraits< Problem, typename GridGeometry::DiscretizationMethod >::BoundaryTypes BoundaryTypes
Definition: common/typetraits/problem.hh:46
Declares all properties used in Dumux.
Type traits for problem classes.