1 #ifndef amandus_refine_strategy_h 2 #define amandus_refine_strategy_h 12 #include <amandus/amandus.h> 13 #include <amandus/integrator.h> 14 #include <deal.II/algorithms/operator.h> 15 #include <deal.II/base/smartpointer.h> 16 #include <deal.II/grid/grid_refinement.h> 17 #include <deal.II/numerics/solution_transfer.h> 32 virtual void operator()(
const dealii::Vector<double>& indicator) = 0;
56 this->
tria->set_all_refine_flags();
73 , refine_threshold(refine_threshold)
79 double threshold = this->refine_threshold * indicator.linfty_norm();
80 dealii::GridRefinement::refine(*(this->
tria), indicator, threshold);
96 MarkBulk(dealii::Triangulation<dim>&
tria,
double refine_threshold,
double coarsen_threshold = 0.)
98 , refine_threshold(refine_threshold)
99 , coarsen_threshold(coarsen_threshold)
105 dealii::Vector<double> square_indicators(indicator);
106 square_indicators.scale(indicator);
107 dealii::GridRefinement::refine_and_coarsen_fixed_fraction(
108 *(this->
tria), square_indicators, this->refine_threshold, this->coarsen_threshold);
132 dealii::GridRefinement::refine_and_coarsen_optimize(*(this->
tria), indicator, this->order);
MarkOptimize(dealii::Triangulation< dim > &tria, const unsigned int order=2)
Definition: refine_strategy.h:124
MarkMaximum(dealii::Triangulation< dim > &tria, double refine_threshold)
Definition: refine_strategy.h:71
void operator()(const dealii::Vector< double > &indicator)
Definition: refine_strategy.h:103
double refine_threshold
Definition: refine_strategy.h:84
Definition: refine_strategy.h:38
void operator()(const dealii::Vector< double > &indicator)
Definition: refine_strategy.h:130
MarkBulk(dealii::Triangulation< dim > &tria, double refine_threshold, double coarsen_threshold=0.)
Definition: refine_strategy.h:96
int order
Definition: refine_strategy.h:136
mark with maximum criterium
Definition: refine_strategy.h:68
mark with bulk criterion
Definition: refine_strategy.h:93
void operator()(const dealii::Vector< double > &indicator)
Definition: refine_strategy.h:77
virtual void operator()(const dealii::Vector< double > &indicator)=0
abstract base class AmandusRefineStrategy
Definition: refine_strategy.h:25
AmandusRefineStrategy(dealii::Triangulation< dim > &tria)
Definition: refine_strategy.h:28
mark optimize
Definition: refine_strategy.h:121
dealii::SmartPointer< dealii::Triangulation< dim >, AmandusRefineStrategy< dim > > tria
Definition: refine_strategy.h:35
double refine_threshold
Definition: refine_strategy.h:112