Marker.hpp 8.71 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
#pragma once
// TODO: Cleanup of copied comments
#include <string>

#include <dune/grid/common/grid.hh>

#include <amdis/AdaptInfo.hpp>
#include <amdis/Flag.hpp>
#include <amdis/Initfile.hpp>

namespace AMDiS {

  
  /**
   * \ingroup Adaption
   *
   * \brief
   * Base class for all scalar markers.
   */
20
  template <class Traits>
21
22
  class Marker
  {
23
24
25
    using GlobalBasis = typename Traits::GlobalBasis;
    using GridView = typename GlobalBasis::GridView;
    using Grid     = typename GridView::Grid;
26
27
    using IndexSet = typename GridView::IndexSet;
    using Element  = typename GridView::template Codim<0>::Entity;
28
    using EstType = std::vector<double>;
29
30
31
32
33
34

  public:
    /// Constructor.
    Marker() {}

    /// Constructor.
35
    Marker(std::string name_, std::string component_, const EstType& est_, std::shared_ptr<Grid> grid_)
36
      : name(name_),
37
        component(component_),
38
39
40
41
42
43
        grid(grid_),
        est(est_),
        maximumMarking(false),
        p(2),
        info(10),
        maxRefineLevel(-1),
44
45
46
        minRefineLevel(-1),
        refineAllowed(true),
        coarsenAllowed(false)
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
    {
      Parameters::get(name + "->p", p);
      Parameters::get(name + "->info", info);
      Parameters::get(name + "->max refinement level", maxRefineLevel);
      Parameters::get(name + "->min refinement level", minRefineLevel);
    }

    /// destructor
    virtual ~Marker() {}

    /// Marks element with newMark. If \ref maximumMarking is set, the element
    /// is marked only if the new mark is bigger than the old one. The return
    /// value specifies whether the element has been marked, or not.
    inline void mark(const Element& elem, int newMark)
    {
      int oldMark = grid->getMark(elem);

      if (!maximumMarking || (newMark > oldMark)) {
        bool (marked) = grid->mark(newMark, elem);
        if (marked) {
          if (oldMark > 0) {
            elMarkRefine--;
          } else if (oldMark < 0) {
            elMarkCoarsen--;
          }

          if (newMark > 0) {
            elMarkRefine++;
          } else if (newMark < 0) {
              elMarkCoarsen++;
          }
        } else {
          msg("Marking failed");
        }
      }
    }

    /// Can be used by sub classes. Called before traversal.
    virtual void initMarking(AdaptInfo& adaptInfo);

    /// Can be used by sub classes. Called after traversal.
    virtual void finishMarking(AdaptInfo& adaptInfo);

    /// Marks one element.
    virtual void markElement(AdaptInfo& adaptInfo, const Element& elem);

    /// Marking of the mesh.
    virtual Flag markGrid(AdaptInfo& adaptInfo);

    /// Sets \ref maximumMarking.
    inline void setMaximumMarking(bool maxMark)
    {
      maximumMarking = maxMark;
    }

    inline int getElMarkRefine()
    {
      return elMarkRefine;
    }

    inline int getElMarkCoarsen()
    {
      return elMarkCoarsen;
    }

    /// Returns \ref name of the Marker
    inline std::string getName() const
    {
      return name;
    }

    /// Creates a scalar marker depending on the strategy set in parameters.
119
    static std::shared_ptr<Marker<Traits> > createMarker(std::string name, std::string component_, const EstType& est_, std::shared_ptr<Grid> grid_);
120
121
122
123
124

  protected:
    /// Name of the scalar marker.
    std::string name;

125
126
    /// Problem component to work on
    std::string component;
127
128

    /// Pointer to the grid
129
    std::shared_ptr<Grid> grid;
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168

    /// Pointer to the local error estimates
    EstType est;

    /// estimation sum
    double estSum;

    /// estmation maximum
    double estMax;

    /// If true, elements are marked only if the new value is bigger than
    /// the current marking.
    bool maximumMarking;

    /// Lower limit for error estimation, from which an element is marked for
    /// refinement
    double markRLimit;

    /// Upper limit for error estimation, from which an element is marked for
    /// coarsening
    double markCLimit;

    /// power in estimator norm
    double p;

    /// Info level.
    int info;

    /// Counter for elements marked for refinement
    int elMarkRefine;

    /// Counter for elements marked for coarsening
    int elMarkCoarsen;

    /// Maximal level of all elements.
    int maxRefineLevel;

    /// Minimal level of all elements.
    int minRefineLevel;
169
170
171
172

    bool refineAllowed;

    bool coarsenAllowed;
173
174
175
176
177
178
179
180
181
  };


  /**
   * \ingroup Adaption
   *
   * \brief
   * Global refinement.
   */
182
183
  template <class Traits>
  class GRMarker : public Marker<Traits>
184
  {
185
186
187
    using GlobalBasis = typename Traits::GlobalBasis;
    using GridView = typename GlobalBasis::GridView;
    using Grid     = typename GridView::Grid;
188
189
    using IndexSet = typename GridView::IndexSet;
    using Element  = typename GridView::template Codim<0>::Entity;
190
    using EstType = std::vector<double>;
191
192
193

  public:
    /// Constructor.
194
    GRMarker(std::string name_, std::string component_, const EstType& est_, std::shared_ptr<Grid> grid_)
195
      : Marker<Traits>(name_, component_, est_, grid_)
196
197
198
199
200
    {}

    /// Implementation of Marker::markElement().
    virtual void markElement(AdaptInfo& adaptInfo, const Element& elem)
    {
201
      if (this->refineAllowed)
202
203
204
205
206
207
208
209
210
211
212
213
        this->mark(elem, 1);
    }
  };


  /**
   * \ingroup Adaption
   *
   * \brief
   * Maximum strategy.
   */
  
214
215
  template <class Traits>
  class MSMarker : public Marker<Traits>
216
  {
217
218
219
    using GlobalBasis = typename Traits::GlobalBasis;
    using GridView = typename GlobalBasis::GridView;
    using Grid     = typename GridView::Grid;
220
221
    using IndexSet = typename GridView::IndexSet;
    using Element  = typename GridView::template Codim<0>::Entity;
222
    using EstType = std::vector<double>;
223
224
225

  public:
    /// Constructor.
226
    MSMarker(std::string name_, std::string component_, const EstType& est_, std::shared_ptr<Grid> grid_)
227
      : Marker<Traits>(name_, component_, est_, grid_),
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
        MSGamma(0.5),
        MSGammaC(0.1)
    {
      Parameters::get(this->name + "->MSGamma", MSGamma);
      Parameters::get(this->name + "->MSGammaC", MSGammaC);
    }

    /// Implementation of Marker::initMarking().
    void initMarking(AdaptInfo& adaptInfo);

  protected:
    /// Marking parameter.
    double MSGamma;

    /// Marking parameter.
    double MSGammaC;
  };


  /**
   * \ingroup Adaption
   *
   * \brief
   * Equidistribution strategy
   */
  
254
255
  template <class Traits>
  class ESMarker : public Marker<Traits>
256
  {
257
258
259
    using GlobalBasis = typename Traits::GlobalBasis;
    using GridView = typename GlobalBasis::GridView;
    using Grid     = typename GridView::Grid;
260
261
    using IndexSet = typename GridView::IndexSet;
    using Element  = typename GridView::template Codim<0>::Entity;
262
    using EstType = std::vector<double>;
263
264
265

  public:
    /// Constructor.
266
    ESMarker(std::string name_, std::string component_, const EstType& est_, std::shared_ptr<Grid> grid_)
267
      : Marker<Traits>(name_, component_, est_, grid_),
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
        ESTheta(0.9),
        ESThetaC(0.2)
    {
      Parameters::get(this->name + "->ESTheta", ESTheta);
      Parameters::get(this->name + "->ESThetaC", ESThetaC);
    }

    /// Implementation of Marker::initMarking().
    virtual void initMarking(AdaptInfo& adaptInfo);

  protected:
    /// Marking parameter.
    double ESTheta;

    /// Marking parameter.
    double ESThetaC;
  };


  /**
   * \ingroup Adaption
   *
   * \brief
   * Guaranteed error reduction strategy
   */
  
294
295
  template <class Traits>
  class GERSMarker : public Marker<Traits>
296
  {
297
298
299
    using GlobalBasis = typename Traits::GlobalBasis;
    using GridView = typename GlobalBasis::GridView;
    using Grid     = typename GridView::Grid;
300
301
    using IndexSet = typename GridView::IndexSet;
    using Element  = typename GridView::template Codim<0>::Entity;
302
    using EstType = std::vector<double>;
303
304
305

  public:
    /// Constructor.
306
    GERSMarker(std::string name_, std::string component_, const EstType& est_, std::shared_ptr<Grid> grid_)
307
      : Marker<Traits>(name_, component_, est_, grid_),
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
        oldErrSum(0.0),
        GERSThetaStar(0.6),
        GERSNu(0.1),
        GERSThetaC(0.1)
    {
      Parameters::get(this->name + "->GERSThetaStar", GERSThetaStar);
      Parameters::get(this->name + "->GERSNu", GERSNu);
      Parameters::get(this->name + "->GERSThetaC", GERSThetaC);
    }

    /// Implementation of Marker::markGrid().
    virtual Flag markGrid(AdaptInfo& adaptInfo);

  protected:
    /// Refinement marking function.
    void markElementForRefinement(AdaptInfo& adaptInfo, const Element& elem);

    /// Coarsening marking function.
    void markElementForCoarsening(AdaptInfo& adaptInfo, const Element& elem);

  protected:
    /// Marking parameter.
    double GERSSum;

    /// Marking parameter.
    double oldErrSum;

    /// Marking parameter.
    double GERSThetaStar;

    /// Marking parameter.
    double GERSNu;

    /// Marking parameter.
    double GERSThetaC;
  };

}
346
347

#include "Marker.inc.hpp"