version 3.10-dev
initializationindicator.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_GRIDADAPTINITIALIZATIONINDICATOR_HH
13#define DUMUX_GRIDADAPTINITIALIZATIONINDICATOR_HH
14
15#include <memory>
16
17#include <dune/geometry/type.hh>
22
23namespace Dumux {
24
32template<class TypeTag>
34{
38 using Element = typename GridView::Traits::template Codim<0>::Entity;
39
42
44
45public:
46
60 GridAdaptInitializationIndicator(std::shared_ptr<const Problem> problem,
61 std::shared_ptr<const GridGeometry> gridGeometry,
62 std::shared_ptr<const GridVariables> gridVariables)
63 : problem_(problem)
64 , gridGeometry_(gridGeometry)
65 , gridVariables_(gridVariables)
66 , minLevel_(getParamFromGroup<std::size_t>(problem->paramGroup(), "Adaptive.MinLevel"))
67 , maxLevel_(getParamFromGroup<std::size_t>(problem->paramGroup(), "Adaptive.MaxLevel"))
68 , refineAtDirichletBC_(getParamFromGroup<bool>(problem->paramGroup(), "Adaptive.RefineAtDirichletBC", true))
69 , refineAtFluxBC_(getParamFromGroup<bool>(problem->paramGroup(), "Adaptive.RefineAtFluxBC", true))
70 , refineAtSource_(getParamFromGroup<bool>(problem->paramGroup(), "Adaptive.RefineAtSource", true))
71 , eps_(getParamFromGroup<Scalar>(problem->paramGroup(), "Adaptive.BCRefinementThreshold", 1e-10))
72 {}
73
77 void setMinLevel(std::size_t minLevel)
78 {
79 minLevel_ = minLevel;
80 }
81
85 void setMaxLevel(std::size_t maxLevel)
86 {
87 maxLevel_ = maxLevel;
88 }
89
93 void setLevels(std::size_t minLevel, std::size_t maxLevel)
94 {
95 minLevel_ = minLevel;
96 maxLevel_ = maxLevel;
97 }
98
103 {
104 refineAtDirichletBC_ = refine;
105 }
106
110 void setRefinementAtSources(bool refine)
111 {
112 refineAtSource_ = refine;
113 }
114
119 {
120 refineAtFluxBC_ = refine;
121 }
122
128 template<class SolutionVector>
129 void calculate(const SolutionVector& sol)
130 {
132 indicatorVector_.assign(gridGeometry_->gridView().size(0), false);
133
134 // get the fvGeometry and elementVolVars needed for the bc and source interfaces
135 auto fvGeometry = localView(*gridGeometry_);
136 auto elemVolVars = localView(gridVariables_->curGridVolVars());
137 // elemFluxVarsCache for neumann interface
138 auto elemFluxVarsCache = localView(gridVariables_->gridFluxVarsCache());
139
140 for (const auto& element : elements(gridGeometry_->gridView()))
141 {
142 const auto eIdx = gridGeometry_->elementMapper().index(element);
143
145 if (element.level() < minLevel_)
146 {
147 indicatorVector_[eIdx] = true;
148 continue; // proceed to the next element
149 }
150
151 // If refinement at sources/BCs etc is deactivated, skip the rest
152 if (!refineAtSource_ && !refineAtFluxBC_ && !refineAtDirichletBC_)
153 continue;
154
155 // if the element is already on the maximum permissive level, skip rest
156 if (element.level() == maxLevel_)
157 continue;
158
159 // Bind all of the local views
160 fvGeometry.bind(element);
161 elemVolVars.bind(element, fvGeometry, sol);
162 elemFluxVarsCache.bind(element, fvGeometry, elemVolVars);
163
164
166 if (refineAtSource_)
167 {
168 for (const auto& scv : scvs(fvGeometry))
169 {
170 auto source = problem_->source(element, fvGeometry, elemVolVars, scv);
171 auto pointSource = problem_->scvPointSources(element, fvGeometry, elemVolVars, scv);
172 if (source.infinity_norm() + pointSource.infinity_norm() > eps_)
173 {
174 indicatorVector_[eIdx] = true;
175 break; // element is marked, escape scv loop
176 }
177 }
178 }
179
181 if (!indicatorVector_[eIdx] // proceed if element is not already marked
182 && element.hasBoundaryIntersections() // proceed if element is on boundary
183 && (refineAtDirichletBC_ || refineAtFluxBC_)) // proceed if boundary refinement is active
184 {
185 // cell-centered schemes
186 if (!isBox)
187 {
188 for (const auto& scvf : scvfs(fvGeometry))
189 {
190 // skip non-boundary scvfs
191 if (!scvf.boundary())
192 continue;
193
194 const auto bcTypes = problem_->boundaryTypes(element, scvf);
195 // We assume pure BCs, mixed boundary types are not allowed anyway!
196 if(bcTypes.hasOnlyDirichlet() && refineAtDirichletBC_)
197 {
198 indicatorVector_[eIdx] = true;
199 break; // element is marked, escape scvf loop
200 }
201
202 // we are on a pure Neumann boundary
203 else if(refineAtFluxBC_)
204 {
205 const auto fluxes = problem_->neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
206 if (fluxes.infinity_norm() > eps_)
207 {
208 indicatorVector_[eIdx] = true;
209 break; // element is marked, escape scvf loop
210 }
211 }
212 }
213 }
214 // box-scheme
215 else
216 {
217 // container to store bcTypes
219 std::vector<BoundaryTypes> bcTypes(fvGeometry.numScv());
220
221 // Get bcTypes and maybe mark for refinement on Dirichlet boundaries
222 for (const auto& scv : scvs(fvGeometry))
223 {
224 bcTypes[scv.localDofIndex()] = problem_->boundaryTypes(element, scv);
225 if (refineAtDirichletBC_ && bcTypes[scv.localDofIndex()].hasDirichlet())
226 {
227 indicatorVector_[eIdx] = true;
228 break; // element is marked, escape scv loop
229 }
230 }
231
232 // If element hasn't been marked because of Dirichlet BCS, check Neumann BCs
233 if (!indicatorVector_[eIdx] && refineAtFluxBC_)
234 {
235 for (const auto& scvf : scvfs(fvGeometry))
236 {
238 if (scvf.boundary() && bcTypes[scvf.insideScvIdx()].hasNeumann())
239 {
240 const auto fluxes = problem_->neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
241 if (fluxes.infinity_norm() > eps_)
242 {
243 indicatorVector_[eIdx] = true;
244 break; // element is marked, escape scvf loop
245 }
246 }
247 }
248 }
249 }
250 }
251 }
252 }
253
265 int operator() (const Element& element) const
266 {
267 if (indicatorVector_[gridGeometry_->elementMapper().index(element)])
268 return 1;
269 return 0;
270 }
271
272private:
273 std::shared_ptr<const Problem> problem_;
274 std::shared_ptr<const GridGeometry> gridGeometry_;
275 std::shared_ptr<const GridVariables> gridVariables_;
276 std::vector<bool> indicatorVector_;
277
278 int minLevel_;
279 int maxLevel_;
280 bool refineAtDirichletBC_;
281 bool refineAtFluxBC_;
282 bool refineAtSource_;
283 Scalar eps_;
284};
285
286} // end namespace Dumux
287
288#endif
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:26
Class defining an initialization indicator for grid adaption. Refines at sources and boundaries....
Definition: initializationindicator.hh:34
void setRefinementAtFluxBoundaries(bool refine)
Function to set if refinement at sources is to be used.
Definition: initializationindicator.hh:118
void setLevels(std::size_t minLevel, std::size_t maxLevel)
Function to set the minimum/maximum allowed levels.
Definition: initializationindicator.hh:93
void setMaxLevel(std::size_t maxLevel)
Function to set the maximum allowed level.
Definition: initializationindicator.hh:85
void calculate(const SolutionVector &sol)
Calculates the indicator used for refinement/coarsening for each grid cell.
Definition: initializationindicator.hh:129
void setRefinementAtDirichletBC(bool refine)
Function to set if refinement at Dirichlet boundaries is to be used.
Definition: initializationindicator.hh:102
int operator()(const Element &element) const
function call operator to return mark
Definition: initializationindicator.hh:265
GridAdaptInitializationIndicator(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< const GridVariables > gridVariables)
Constructor.
Definition: initializationindicator.hh:60
void setRefinementAtSources(bool refine)
Function to set if refinement at sources is to be used.
Definition: initializationindicator.hh:110
void setMinLevel(std::size_t minLevel)
Function to set the minimum allowed level.
Definition: initializationindicator.hh:77
Defines all properties used in Dumux.
Type traits for problem classes.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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 available discretization methods in Dumux.
constexpr Box box
Definition: method.hh:147
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
typename Detail::template ProblemTraits< Problem, typename GridGeometry::DiscretizationMethod >::BoundaryTypes BoundaryTypes
Definition: common/typetraits/problem.hh:34