Amandus: Simulations based on multilevel Schwarz methods
refine_strategy.h
Go to the documentation of this file.
1 #ifndef amandus_refine_strategy_h
2 #define amandus_refine_strategy_h
3 
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>
18 
24 template <int dim>
26 {
27 public:
28  AmandusRefineStrategy(dealii::Triangulation<dim>& tria)
29  : tria(&tria)
30  {
31  }
32  virtual void operator()(const dealii::Vector<double>& indicator) = 0;
33 
34 protected:
35  dealii::SmartPointer<dealii::Triangulation<dim>, AmandusRefineStrategy<dim>> tria;
36 };
37 
38 namespace RefineStrategy
39 {
45 template <int dim>
47 {
48 public:
49  MarkUniform(dealii::Triangulation<dim>& tria)
50  : AmandusRefineStrategy<dim>(tria)
51  {
52  }
53 
54  void operator()(const dealii::Vector<double>& /*indicator*/)
55  {
56  this->tria->set_all_refine_flags();
57  }
58 
59  void operator()() { this->tria->set_all_refine_flags(); }
60 };
61 
67 template <int dim>
69 {
70 public:
71  MarkMaximum(dealii::Triangulation<dim>& tria, double refine_threshold)
72  : AmandusRefineStrategy<dim>(tria)
73  , refine_threshold(refine_threshold)
74  {
75  }
76 
77  void operator()(const dealii::Vector<double>& indicator)
78  {
79  double threshold = this->refine_threshold * indicator.linfty_norm();
80  dealii::GridRefinement::refine(*(this->tria), indicator, threshold);
81  }
82 
83 protected:
85 };
86 
92 template <int dim>
93 class MarkBulk : public AmandusRefineStrategy<dim>
94 {
95 public:
96  MarkBulk(dealii::Triangulation<dim>& tria, double refine_threshold, double coarsen_threshold = 0.)
97  : AmandusRefineStrategy<dim>(tria)
98  , refine_threshold(refine_threshold)
99  , coarsen_threshold(coarsen_threshold)
100  {
101  }
102 
103  void operator()(const dealii::Vector<double>& indicator)
104  {
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);
109  }
110 
111 protected:
112  double refine_threshold, coarsen_threshold;
113 };
114 
120 template <int dim>
122 {
123 public:
124  MarkOptimize(dealii::Triangulation<dim>& tria, const unsigned int order = 2)
125  : AmandusRefineStrategy<dim>(tria)
126  , order(order)
127  {
128  }
129 
130  void operator()(const dealii::Vector<double>& indicator)
131  {
132  dealii::GridRefinement::refine_and_coarsen_optimize(*(this->tria), indicator, this->order);
133  }
134 
135 protected:
136  int order;
137 };
138 }
139 
140 #endif
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()()
Definition: refine_strategy.h:59
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
MarkUniform(dealii::Triangulation< dim > &tria)
Definition: refine_strategy.h:49
AmandusRefineStrategy(dealii::Triangulation< dim > &tria)
Definition: refine_strategy.h:28
mark all, i.e. uniform refinement
Definition: refine_strategy.h:46
mark optimize
Definition: refine_strategy.h:121
dealii::SmartPointer< dealii::Triangulation< dim >, AmandusRefineStrategy< dim > > tria
Definition: refine_strategy.h:35
void operator()(const dealii::Vector< double > &)
Definition: refine_strategy.h:54
double refine_threshold
Definition: refine_strategy.h:112