8#include <dune/grid/common/grid.hh>
10#include <amdis/common/ConceptsBase.hpp>
11#include <amdis/common/TypeTraits.hpp>
13#include <amdis/gridfunctions/GridFunction.hpp>
15#include <amdis/AdaptInfo.hpp>
16#include <amdis/Flag.hpp>
17#include <amdis/Initfile.hpp>
31 using GridView =
typename Grid::LeafGridView;
32 using Element =
typename GridView::template Codim<0>::Entity;
39 Marker(std::string
const&
name, std::shared_ptr<Grid>
const& grid)
54 void mark(Element
const& elem,
int newMark);
81 std::string
const&
name()
const
144 template <
class Gr
id>
150 using Element =
typename Super::Element;
151 using Estimates = std::vector<double>;
159 std::shared_ptr<Grid>
const& grid)
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);
209 template <
class Gr
id>
214 using Element =
typename Super::Element;
215 using Estimates =
typename Super::Estimates;
219 GRMarker(std::string
const&
name, std::string
const& component, Estimates
const& est,
220 std::shared_ptr<Grid>
const& grid)
239 template <
class Gr
id>
244 using Estimates =
typename Super::Estimates;
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}
274 template <
class Gr
id>
279 using Estimates =
typename Super::Estimates;
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}
309 template <
class Gr
id>
314 using Element =
typename Super::Element;
315 using Estimates =
typename Super::Estimates;
320 std::shared_ptr<Grid>
const& grid)
321 :
Super{
name, std::move(component), est, grid}
369 template <
class Gr
id,
class Gr
idFct>
374 using Element =
typename Super::Element;
377 using IsGridFunction =
decltype(localFunction(std::declval<GF>()));
379 static_assert(Dune::Std::is_detected<IsGridFunction,GridFct>::value,
"");
405 template <
class Gr
id,
class GF>
406 GridFunctionMarker(std::string
const& name, std::shared_ptr<Grid>
const& grid, GF&& gf)
407 -> GridFunctionMarker<Grid,
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