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)
160 :
Super{std::move(name), grid}
161 , component_(std::move(component))
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);
185 double estSum_ = 0.0;
188 double estMax_ = 0.0;
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}
260 double msGamma_ = 0.5;
263 double msGammaC_ = 0.1;
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}
295 double esTheta_ = 0.9;
298 double esThetaC_ = 0.2;
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}
332 void markElementForRefinement(
AdaptInfo& adaptInfo, Element
const& elem);
335 void markElementForCoarsening(
AdaptInfo& adaptInfo, Element
const& elem);
339 double gersSum_ = 0.0;
342 double oldErrSum_ = 0.0;
345 double gersThetaStar_ = 0.6;
348 double gersNu_ = 0.1;
351 double gersThetaC_ = 0.1;
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>
412 #include "Marker.inc.hpp" GERSMarker(std::string const &name, std::string component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:319
The Flag class encapsulates flags which represents simple information. Used e.g. while mesh traversal...
Definition: Flag.hpp:13
Base class for all estimator-based markers.
Definition: Marker.hpp:145
int elMarkCoarsen() const
Returns elMarkCoarsen_ of the Marker.
Definition: Marker.hpp:75
Estimates est_
Pointer to the local error estimates.
Definition: Marker.hpp:182
virtual void finishMarking(AdaptInfo &adaptInfo)
Called after traversal.
Definition: Marker.inc.hpp:40
void markElement(AdaptInfo &, Element const &elem) override
Marks one element.
Definition: Marker.hpp:224
EstimatorMarker(std::string name, std::string component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:158
bool maximumMarking_
Definition: Marker.hpp:113
Maximum strategy.
Definition: Marker.hpp:240
void allowCoarsening(bool allow)
Sets coarsenAllowed_.
Definition: Marker.hpp:99
decltype(auto) makeGridFunction(PreGridFct const &preGridFct, GridView const &gridView)
Generator for Gridfunctions from Expressions (PreGridfunctions)
Definition: GridFunction.hpp:154
Marker based on an indicator given as grid-function.
Definition: Marker.hpp:370
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
MSMarker(std::string const &name, std::string component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:248
Base class for all markers.
Definition: Marker.hpp:28
Global refinement.
Definition: Marker.hpp:210
int maxRefineLevel_
Maximal level of all elements.
Definition: Marker.hpp:125
Marker()=default
Constructor.
static std::optional< T > get(std::string const &key)
Get parameter-values from parameter-tree.
Definition: Initfile.hpp:25
ESMarker(std::string const &name, std::string component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:283
double markCLimit_
Definition: Marker.hpp:196
Marker(std::string const &name, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:39
std::string name_
Name of the marker.
Definition: Marker.hpp:106
void setMaximumMarking(bool maxMark)
Sets maximumMarking_.
Definition: Marker.hpp:87
GridFunctionMarker(std::string const &name, std::shared_ptr< Grid > const &grid, GF &&gf)
Constructor.
Definition: Marker.hpp:384
void mark(Element const &elem, int newMark)
Definition: Marker.inc.hpp:6
virtual Flag markGrid(AdaptInfo &adaptInfo)
Marking of the whole grid.
Definition: Marker.inc.hpp:50
virtual void initMarking(AdaptInfo &adaptInfo)
Called before traversal.
Definition: Marker.inc.hpp:32
bool refineAllowed_
Allow elements to be marked for refinement.
Definition: Marker.hpp:131
Holds adapt parameters and infos about the problem.
Definition: AdaptInfo.hpp:25
double markRLimit_
Definition: Marker.hpp:192
std::string component_
Problem component to work on.
Definition: Marker.hpp:179
Guaranteed error reduction strategy.
Definition: Marker.hpp:310
int minRefineLevel_
Minimal level of all elements.
Definition: Marker.hpp:128
int elMarkCoarsen_
Counter for elements marked for coarsening.
Definition: Marker.hpp:122
GRMarker(std::string const &name, std::string const &component, Estimates const &est, std::shared_ptr< Grid > const &grid)
Constructor.
Definition: Marker.hpp:219
int elMarkRefine() const
Returns elMarkRefine_ of the Marker.
Definition: Marker.hpp:69
std::shared_ptr< Grid > grid_
Pointer to the grid.
Definition: Marker.hpp:109
bool coarsenAllowed_
Allow elements to be marked for coarsening.
Definition: Marker.hpp:134
virtual void markElement(AdaptInfo &adaptInfo, Element const &elem)=0
Marks one element.
void allowRefinement(bool allow)
Sets refineAllowed_.
Definition: Marker.hpp:93
std::string const & name() const
Returns name_ of the Marker.
Definition: Marker.hpp:81
int elMarkRefine_
Counter for elements marked for refinement.
Definition: Marker.hpp:119
Equidistribution strategy.
Definition: Marker.hpp:275
int info_
Info level.
Definition: Marker.hpp:116
virtual ~Marker()=default
Destructor.
void markElement(AdaptInfo &, Element const &) final
Implementation of Marker::markElement. Does nothing since marking is done in markGrid().
Definition: Marker.hpp:391