AMDiS 2.10
The Adaptive Multi-Dimensional Simulation Toolbox
Marker.hpp
1#pragma once
2
3#include <limits>
4#include <optional>
5#include <string>
6#include <utility>
7
8#include <dune/grid/common/grid.hh>
9
10#include <amdis/common/ConceptsBase.hpp>
11#include <amdis/common/TypeTraits.hpp>
12
13#include <amdis/gridfunctions/GridFunction.hpp>
14
15#include <amdis/AdaptInfo.hpp>
16#include <amdis/Flag.hpp>
17#include <amdis/Initfile.hpp>
18
19namespace AMDiS
20{
27 template <class Grid>
28 class Marker
29 {
30 protected:
31 using GridView = typename Grid::LeafGridView;
32 using Element = typename GridView::template Codim<0>::Entity;
33
34 public:
36 Marker() = default;
37
39 Marker(std::string const& name, std::shared_ptr<Grid> const& grid)
40 : name_(name)
41 , grid_(grid)
42 {
43 Parameters::get(name + "->info", info_);
44 Parameters::get(name + "->max refinement level", maxRefineLevel_);
45 Parameters::get(name + "->min refinement level", minRefineLevel_);
46 }
47
49 virtual ~Marker() = default;
50
54 void mark(Element const& elem, int newMark);
55
57 virtual void initMarking(AdaptInfo& adaptInfo);
58
60 virtual void finishMarking(AdaptInfo& adaptInfo);
61
63 virtual void markElement(AdaptInfo& adaptInfo, Element const& elem) = 0;
64
66 virtual Flag markGrid(AdaptInfo& adaptInfo);
67
69 int elMarkRefine() const
70 {
71 return elMarkRefine_;
72 }
73
75 int elMarkCoarsen() const
76 {
77 return elMarkCoarsen_;
78 }
79
81 std::string const& name() const
82 {
83 return name_;
84 }
85
87 void setMaximumMarking(bool maxMark)
88 {
89 maximumMarking_ = maxMark;
90 }
91
93 void allowRefinement(bool allow)
94 {
95 refineAllowed_ = allow;
96 }
97
99 void allowCoarsening(bool allow)
100 {
101 coarsenAllowed_ = allow;
102 }
103
104 protected:
106 std::string name_;
107
109 std::shared_ptr<Grid> grid_;
110
113 bool maximumMarking_ = false;
114
116 int info_ = 0;
117
120
123
125 int maxRefineLevel_ = std::numeric_limits<int>::max();
126
129
131 bool refineAllowed_ = true;
132
134 bool coarsenAllowed_ = true;
135 };
136
137
144 template <class Grid>
146 : public Marker<Grid>
147 {
148 protected:
149 using Super = Marker<Grid>;
150 using Element = typename Super::Element;
151 using Estimates = std::vector<double>;
152
153 public:
155 EstimatorMarker() = default;
156
158 EstimatorMarker(std::string name, std::string component, Estimates const& est,
159 std::shared_ptr<Grid> const& grid)
160 : Super{std::move(name), grid}
161 , component_(std::move(component))
162 , est_(est)
163 {
164 Parameters::get(this->name_ + "->p", p_);
165 }
166
168 void initMarking(AdaptInfo& adaptInfo) override;
169
171 void markElement(AdaptInfo& adaptInfo, Element const& elem) override;
172
174 static std::unique_ptr<EstimatorMarker<Grid> > createMarker(std::string const& name,
175 std::string const& component, Estimates const& est, std::shared_ptr<Grid> const& grid);
176
177 protected:
179 std::string component_;
180
182 Estimates est_;
183
185 double estSum_ = 0.0;
186
188 double estMax_ = 0.0;
189
193
197
199 double p_ = 2.0;
200 };
201
202
209 template <class Grid>
211 : public EstimatorMarker<Grid>
212 {
214 using Element = typename Super::Element;
215 using Estimates = typename Super::Estimates;
216
217 public:
219 GRMarker(std::string const& name, std::string const& component, Estimates const& est,
220 std::shared_ptr<Grid> const& grid)
221 : Super{name, component, est, grid}
222 {}
223
224 void markElement(AdaptInfo& /*adaptInfo*/, Element const& elem) override
225 {
226 if (this->refineAllowed_)
227 this->mark(elem, 1);
228 }
229 };
230
231
239 template <class Grid>
241 : public EstimatorMarker<Grid>
242 {
244 using Estimates = typename Super::Estimates;
245
246 public:
248 MSMarker(std::string const& name, std::string component, Estimates const& est,
249 std::shared_ptr<Grid> const& grid)
250 : Super{name, std::move(component), est, grid}
251 {
252 Parameters::get(name + "->MSGamma", msGamma_);
253 Parameters::get(name + "->MSGammaC", msGammaC_);
254 }
255
256 void initMarking(AdaptInfo& adaptInfo) override;
257
258 protected:
260 double msGamma_ = 0.5;
261
263 double msGammaC_ = 0.1;
264 };
265
266
274 template <class Grid>
276 : public EstimatorMarker<Grid>
277 {
279 using Estimates = typename Super::Estimates;
280
281 public:
283 ESMarker(std::string const& name, std::string component, Estimates const& est,
284 std::shared_ptr<Grid> const& grid)
285 : Super{name, std::move(component), est, grid}
286 {
287 Parameters::get(name + "->ESTheta", esTheta_);
288 Parameters::get(name + "->ESThetaC", esThetaC_);
289 }
290
291 void initMarking(AdaptInfo& adaptInfo) override;
292
293 protected:
295 double esTheta_ = 0.9;
296
298 double esThetaC_ = 0.2;
299 };
300
301
309 template <class Grid>
311 : public EstimatorMarker<Grid>
312 {
314 using Element = typename Super::Element;
315 using Estimates = typename Super::Estimates;
316
317 public:
319 GERSMarker(std::string const& name, std::string component, Estimates const& est,
320 std::shared_ptr<Grid> const& grid)
321 : Super{name, std::move(component), est, grid}
322 {
323 Parameters::get(name + "->GERSThetaStar", gersThetaStar_);
324 Parameters::get(name + "->GERSNu", gersNu_);
325 Parameters::get(name + "->GERSThetaC", gersThetaC_);
326 }
327
328 Flag markGrid(AdaptInfo& adaptInfo) override;
329
330 protected:
332 void markElementForRefinement(AdaptInfo& adaptInfo, Element const& elem);
333
335 void markElementForCoarsening(AdaptInfo& adaptInfo, Element const& elem);
336
337 protected:
339 double gersSum_ = 0.0;
340
342 double oldErrSum_ = 0.0;
343
345 double gersThetaStar_ = 0.6;
346
348 double gersNu_ = 0.1;
349
351 double gersThetaC_ = 0.1;
352 };
353
354
369 template <class Grid, class GridFct>
371 : public Marker<Grid>
372 {
373 using Super = Marker<Grid>;
374 using Element = typename Super::Element;
375
376 template <class GF>
377 using IsGridFunction = decltype(localFunction(std::declval<GF>()));
378
379 static_assert(Dune::Std::is_detected<IsGridFunction,GridFct>::value, "");
380
381 public:
383 template <class GF>
384 GridFunctionMarker(std::string const& name, std::shared_ptr<Grid> const& grid, GF&& gf)
385 : Super{name, grid}
386 , gridFct_{makeGridFunction(FWD(gf), grid->leafGridView())}
387 {}
388
391 void markElement(AdaptInfo& /*adaptInfo*/, Element const& /*elem*/) final { /* do nothing */}
392
396 Flag markGrid(AdaptInfo& adaptInfo) override;
397
398 private:
400 GridFct gridFct_;
401 };
402
403
404 // Deduction guide for GridFunctionMarker class
405 template <class Grid, class GF>
406 GridFunctionMarker(std::string const& name, std::shared_ptr<Grid> const& grid, GF&& gf)
407 -> GridFunctionMarker<Grid,
408 TYPEOF( makeGridFunction(FWD(gf), grid->leafGridView()) )>;
409
410} // end namespace AMDiS
411
412#include "Marker.inc.hpp"
Holds adapt parameters and infos about the problem.
Definition: AdaptInfo.hpp:26
Equidistribution strategy.
Definition: Marker.hpp:277
double esTheta_
Marking parameter.
Definition: Marker.hpp:295
void initMarking(AdaptInfo &adaptInfo) override
Can be used by sub classes. Called before traversal.
Definition: Marker.inc.hpp:148
double esThetaC_
Marking parameter.
Definition: Marker.hpp:298
ESMarker(std::string const &name, std::string component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:283
Base class for all estimator-based markers.
Definition: Marker.hpp:147
double markRLimit_
Definition: Marker.hpp:192
EstimatorMarker()=default
Constructor.
Estimates est_
Pointer to the local error estimates.
Definition: Marker.hpp:182
double markCLimit_
Definition: Marker.hpp:196
EstimatorMarker(std::string name, std::string component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:158
void markElement(AdaptInfo &adaptInfo, Element const &elem) override
Marks one element.
Definition: Marker.inc.hpp:85
static std::unique_ptr< EstimatorMarker< Grid > > createMarker(std::string const &name, std::string const &component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Creates a scalar marker depending on the strategy set in parameters.
Definition: Marker.inc.hpp:102
void initMarking(AdaptInfo &adaptInfo) override
Can be used by sub classes. Called before traversal.
Definition: Marker.inc.hpp:74
double p_
Power in estimator norm.
Definition: Marker.hpp:199
double estMax_
Estmation maximum.
Definition: Marker.hpp:188
std::string component_
Problem component to work on.
Definition: Marker.hpp:179
double estSum_
Estimation sum.
Definition: Marker.hpp:185
The Flag class encapsulates flags which represents simple information. Used e.g. while mesh traversal...
Definition: Flag.hpp:14
Guaranteed error reduction strategy.
Definition: Marker.hpp:312
Flag markGrid(AdaptInfo &adaptInfo) override
Marking of the whole grid.
Definition: Marker.inc.hpp:169
double gersNu_
Marking parameter.
Definition: Marker.hpp:348
double gersSum_
Marking parameter.
Definition: Marker.hpp:339
void markElementForRefinement(AdaptInfo &adaptInfo, Element const &elem)
Refinement marking function.
Definition: Marker.inc.hpp:243
GERSMarker(std::string const &name, std::string component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:319
double oldErrSum_
Marking parameter.
Definition: Marker.hpp:342
double gersThetaStar_
Marking parameter.
Definition: Marker.hpp:345
double gersThetaC_
Marking parameter.
Definition: Marker.hpp:351
void markElementForCoarsening(AdaptInfo &adaptInfo, Element const &elem)
Coarsening marking function.
Definition: Marker.inc.hpp:255
Global refinement.
Definition: Marker.hpp:212
void markElement(AdaptInfo &, Element const &elem) override
Marks one element.
Definition: Marker.hpp:224
GRMarker(std::string const &name, std::string const &component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:219
Marker based on an indicator given as grid-function.
Definition: Marker.hpp:372
Flag markGrid(AdaptInfo &adaptInfo) override
Definition: Marker.inc.hpp:271
GridFunctionMarker(std::string const &name, std::shared_ptr< Grid > const &grid, GF &&gf)
Constructor.
Definition: Marker.hpp:384
void markElement(AdaptInfo &, Element const &) final
Implementation of Marker::markElement. Does nothing since marking is done in markGrid().
Definition: Marker.hpp:391
static std::optional< T > get(std::string const &key)
Get parameter-values from parameter-tree.
Definition: Initfile.hpp:25
Maximum strategy.
Definition: Marker.hpp:242
double msGammaC_
Marking parameter.
Definition: Marker.hpp:263
MSMarker(std::string const &name, std::string component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:248
void initMarking(AdaptInfo &adaptInfo) override
Can be used by sub classes. Called before traversal.
Definition: Marker.inc.hpp:132
double msGamma_
Marking parameter.
Definition: Marker.hpp:260
Base class for all markers.
Definition: Marker.hpp:29
bool maximumMarking_
Definition: Marker.hpp:113
std::shared_ptr< Grid > grid_
Pointer to the grid.
Definition: Marker.hpp:109
std::string const & name() const
Returns name_ of the Marker.
Definition: Marker.hpp:81
virtual Flag markGrid(AdaptInfo &adaptInfo)
Marking of the whole grid.
Definition: Marker.inc.hpp:50
virtual void finishMarking(AdaptInfo &adaptInfo)
Called after traversal.
Definition: Marker.inc.hpp:40
int maxRefineLevel_
Maximal level of all elements.
Definition: Marker.hpp:125
int info_
Info level.
Definition: Marker.hpp:116
std::string name_
Name of the marker.
Definition: Marker.hpp:106
virtual ~Marker()=default
Destructor.
bool refineAllowed_
Allow elements to be marked for refinement.
Definition: Marker.hpp:131
void mark(Element const &elem, int newMark)
Definition: Marker.inc.hpp:6
virtual void initMarking(AdaptInfo &adaptInfo)
Called before traversal.
Definition: Marker.inc.hpp:32
virtual void markElement(AdaptInfo &adaptInfo, Element const &elem)=0
Marks one element.
int elMarkRefine_
Counter for elements marked for refinement.
Definition: Marker.hpp:119
void allowCoarsening(bool allow)
Sets coarsenAllowed_.
Definition: Marker.hpp:99
void setMaximumMarking(bool maxMark)
Sets maximumMarking_.
Definition: Marker.hpp:87
Marker()=default
Constructor.
Marker(std::string const &name, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:39
void allowRefinement(bool allow)
Sets refineAllowed_.
Definition: Marker.hpp:93
int minRefineLevel_
Minimal level of all elements.
Definition: Marker.hpp:128
int elMarkCoarsen_
Counter for elements marked for coarsening.
Definition: Marker.hpp:122
int elMarkRefine() const
Returns elMarkRefine_ of the Marker.
Definition: Marker.hpp:69
bool coarsenAllowed_
Allow elements to be marked for coarsening.
Definition: Marker.hpp:134
int elMarkCoarsen() const
Returns elMarkCoarsen_ of the Marker.
Definition: Marker.hpp:75
decltype(auto) makeGridFunction(PreGridFct const &preGridFct, GridView const &gridView)
Generator for Gridfunctions from Expressions (PreGridfunctions)
Definition: GridFunction.hpp:168