AdaptInfo.hpp 14.8 KB
Newer Older
1
2
3
#pragma once

// std c++ headers
4
#include <algorithm>
5
#include <cmath>
6
#include <forward_list>
7
#include <limits>
8
#include <map>
9
10
11
12
#include <string>
#include <vector>

// AMDiS includes
13
14
#include <amdis/Output.hpp>
#include <amdis/common/Math.hpp>
15
16
17
18
19
20
21
22

namespace AMDiS
{

  /**
   * \ingroup Adaption
   *
   * \brief
23
   * Holds adapt parameters and infos about the problem.
24
25
26
   */
  class AdaptInfo
  {
27
28
    using List = typename std::vector<std::string>;

29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
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
  protected:
    /** \brief
     * Stores adapt infos for a scalar problem or for one component of a
     * vector valued problem.
     */
    class ScalContent
    {
    public:
      /// Constructor.
      explicit ScalContent(std::string prefix);

      /// Sum of all error estimates
      double est_sum = 0.0;

      /// Sum of all time error estimates
      double est_t_sum = 0.0;

      /// maximal local error estimate.
      double est_max = 0.0;

      /// Maximum of all time error estimates
      double est_t_max = 0.0;

      /// factors to combine max and integral time estimate
      double fac_max = 0.0, fac_sum = 1.0;

      /// Tolerance for the (absolute or relative) error
      double spaceTolerance = 0.0;

      /// Time tolerance.
      double timeTolerance = 0.0;

      /// Relative time tolerance
      double timeRelativeTolerance = 0.0;

      /// Lower bound for the time error.
      double timeErrLow = 0.0;

      /// true if coarsening is allowed, false otherwise.
      int coarsenAllowed = 0;

      /// true if refinement is allowed, false otherwise.
      int refinementAllowed = 1;
    };

  public:
    /// Constructor.
76
    explicit AdaptInfo(std::string name, List compNames = {"0"});
77
78

    /// Destructor.
79
    virtual ~AdaptInfo() {}
80
81
82
83
84
85
86

    /// Resets all variables to zero (or something equivalent)
    void reset();

    /// Returns whether space tolerance is reached.
    virtual bool spaceToleranceReached() const
    {
87
88
      for (auto it = scalContents.begin(); it != scalContents.end(); it++)      {
        if (!(it->second->est_sum < it->second->spaceTolerance))
89
90
91
92
93
94
          return false;
      }

      return true;
    }

95
96
    /// Returns whether space tolerance of component associated with key is reached.
    virtual bool spaceToleranceReached(std::string key) const
97
    {
98
      if (!(scalContents.at(key)->est_sum < scalContents.at(key)->spaceTolerance))
99
100
101
102
103
104
105
106
        return false;
      else
        return true;
    }

    /// Returns whether time tolerance is reached.
    virtual bool timeToleranceReached() const
    {
107
108
      for (auto it = scalContents.begin(); it != scalContents.end(); it++)
        if (!(getTimeEstCombined(it->first) < it->second->timeTolerance))
109
110
111
112
113
          return false;

      return true;
    }

114
115
    /// Returns whether time tolerance of component associated with key is reached.
    virtual bool timeToleranceReached(std::string key) const
116
    {
117
      if (!(getTimeEstCombined(key) < scalContents.at(key)->timeTolerance))
118
119
120
121
122
123
124
125
        return false;
      else
        return true;
    }

    /// Returns whether time error is under its lower bound.
    virtual bool timeErrorLow() const
    {
126
127
      for (auto it = scalContents.begin(); it != scalContents.end(); it++)
        if (!(getTimeEstCombined(it->first) < it->second->timeErrLow))
128
129
130
131
132
133
          return false;

      return true;
    }
    /// Returns the time estimation as a combination
    /// of maximal and integral time error
134
    double getTimeEstCombined(std::string key) const
135
136
    {
      return
137
138
        scalContents.at(key)->est_t_max * scalContents.at(key)->fac_max +
        scalContents.at(key)->est_t_sum * scalContents.at(key)->fac_sum;
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
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
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
254
255
256
257
258
259
260
261
262
263
264
265
    }


    /// Print debug information about time error and its bound.
    void printTimeErrorLowInfo() const;

    /// Returns \ref spaceIteration.
    int getSpaceIteration() const
    {
      return spaceIteration;
    }

    /// Sets \ref spaceIteration.
    void setSpaceIteration(int it)
    {
      spaceIteration = it;
    }

    /// Returns \ref maxSpaceIteration.
    int getMaxSpaceIteration() const
    {
      return maxSpaceIteration;
    }

    /// Sets \ref maxSpaceIteration.
    void setMaxSpaceIteration(int it)
    {
      maxSpaceIteration = it;
    }

    /// Increments \ref spaceIteration by 1;
    void incSpaceIteration()
    {
      spaceIteration++;
    }

    /// Sets \ref timestepIteration.
    void setTimestepIteration(int it)
    {
      timestepIteration = it;
    }

    /// Returns \ref timestepIteration.
    int getTimestepIteration() const
    {
      return timestepIteration;
    }

    /// Increments \ref timestepIteration by 1;
    void incTimestepIteration()
    {
      timestepIteration++;
    }

    /// Returns \ref maxTimestepIteration.
    int getMaxTimestepIteration() const
    {
      return maxTimestepIteration;
    }

    /// Sets \ref maxTimestepIteration.
    void setMaxTimestepIteration(int it)
    {
      maxTimestepIteration = it;
    }

    /// Sets \ref timeIteration.
    void setTimeIteration(int it)
    {
      timeIteration = it;
    }

    /// Returns \ref timeIteration.
    int getTimeIteration() const
    {
      return timeIteration;
    }

    /// Increments \ref timesIteration by 1;
    void incTimeIteration()
    {
      timeIteration++;
    }

    /// Returns \ref maxTimeIteration.
    int getMaxTimeIteration() const
    {
      return maxTimeIteration;
    }

    /// Sets \ref maxTimeIteration.
    void setMaxTimeIteration(int it)
    {
      maxTimeIteration = it;
    }

    /// Returns \ref timestepNumber.
    int getTimestepNumber() const
    {
      return timestepNumber;
    }

    /// Sets \ref timestepNumber.
    void setTimestepNumber(int num)
    {
      timestepNumber = std::min(nTimesteps, num);
    }

    /// Returns \ref nTimesteps.
    int getNumberOfTimesteps() const
    {
      return nTimesteps;
    }

    /// Sets \ref nTimesteps.
    void setNumberOfTimesteps(int num)
    {
      nTimesteps = std::max(0, num);
    }

    /// Increments \ref timestepNumber by 1;
    void incTimestepNumber()
    {
      timestepNumber++;
    }

    /// Sets \ref est_sum.
266
    void setEstSum(double e, std::string key)
267
    {
268
      scalContents.at(key)->est_sum = e;
269
270
271
    }

    /// Sets \ref est_max.
272
    void setEstMax(double e, std::string key)
273
    {
274
      scalContents.at(key)->est_max = e;
275
276
277
    }

    /// Sets \ref est_max.
278
    void setTimeEstMax(double e, std::string key)
279
    {
280
      scalContents.at(key)->est_t_max = e;
281
282
283
    }

    /// Sets \ref est_t_sum.
284
    void setTimeEstSum(double e, std::string key)
285
    {
286
      scalContents.at(key)->est_t_sum = e;
287
288
289
    }

    /// Returns \ref est_sum.
290
    double getEstSum(std::string key) const
291
    {
292
      AMDIS_FUNCNAME_DBG("AdaptInfo::getEstSum()");
293
      test_exit_dbg(scalContents.count(key) == 1, "Wrong key for adaptInfo!\n");
294

295
      return scalContents.at(key)->est_sum;
296
297
298
    }

    /// Returns \ref est_t_sum.
299
    double getEstTSum(std::string key) const
300
    {
301
      return scalContents.at(key)->est_t_sum;
302
303
304
    }

    /// Returns \ref est_max.
305
    double getEstMax(std::string key) const
306
    {
307
      AMDIS_FUNCNAME_DBG("AdaptInfo::getEstSum()");
308
      test_exit_dbg(scalContents.count(key) == 1, "Wrong key for adaptInfo!\n");
309

310
      return scalContents.at(key)->est_max;
311
312
313
    }

    /// Returns \ref est_max.
314
    double getTimeEstMax(std::string key) const
315
    {
316
      return scalContents.at(key)->est_t_max;
317
318
319
    }

    /// Returns \ref est_t_sum.
320
    double getTimeEstSum(std::string key) const
321
    {
322
      return scalContents.at(key)->est_t_sum;
323
324
325
326
327
328
329
330
331
332
333
334
335
336
    }

    /// Returns \ref est_t the estimated overall time error
    double getTimeEst() const
    {
      return est_t;
    }

    void setTimeEst(double value)
    {
      est_t = value;
    }

    /// Returns \ref spaceTolerance.
337
    double getSpaceTolerance(std::string key) const
338
    {
339
      return scalContents.at(key)->spaceTolerance;
340
341
342
    }

    /// Sets \ref spaceTolerance.
343
    void setSpaceTolerance(std::string key, double tol)
344
    {
345
      scalContents.at(key)->spaceTolerance = tol;
346
347
348
    }

    /// Returns \ref timeTolerance.
349
    double getTimeTolerance(std::string key) const
350
    {
351
      return scalContents.at(key)->timeTolerance;
352
353
354
    }

    /// Returns \ref timeRelativeTolerance.
355
    double getTimeRelativeTolerance(std::string key) const
356
    {
357
      return scalContents.at(key)->timeRelativeTolerance;
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
    }

    /// Sets \ref time
    double setTime(double t)
    {
      time = t;
      if (time > endTime)
        time = endTime;
      if (time < startTime)
        time = startTime;

      return time;
    }

    /// Gets \ref time
    double getTime() const
    {
      return time;
    }

    /// Gets \ref &time
    double* getTimePtr()
    {
      return &time;
    }

    /// Sets \ref timestep
    double setTimestep(double t)
    {
      timestep = t;
      if (timestep > maxTimestep)
        timestep = maxTimestep;
      if (timestep < minTimestep)
        timestep = minTimestep;
      if (time + timestep > endTime)
        timestep = endTime - time;

      return timestep;
    }
    /// Gets \ref timestep
    double getTimestep() const
    {
      return timestep;
    }

    void setLastProcessedTimestep(double t)
    {
      lastProcessedTimestep = t;
    }

    double getLastProcessedTimestep() const
    {
      return lastProcessedTimestep;
    }

    /// Returns true, if the end time is reached and no more timestep
    /// computations must be done.
    bool reachedEndTime() const
    {
      if (nTimesteps > 0)
        return !(timestepNumber < nTimesteps);

420
      return !(std::abs(time - endTime) > std::numeric_limits<double>::epsilon());
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
    }


    /// Sets \ref minTimestep
    void setMinTimestep(double t)
    {
      minTimestep = t;
    }

    /// Gets \ref minTimestep
    double getMinTimestep() const
    {
      return minTimestep;
    }

    /// Sets \ref maxTimestep
    void setMaxTimestep(double t)
    {
      maxTimestep = t;
    }

    /// Gets \ref maxTimestep
    double getMaxTimestep() const
    {
      return maxTimestep;
    }

    /// Gets \ref &timestep
    double* getTimestepPtr()
    {
      return &timestep;
    }

    /// Sets \ref startTime = time
    void setStartTime(double time)
    {
      startTime = time;
    }

    /// Sets \ref endTime = time
    void setEndTime(double time)
    {
      endTime = time;
    }

    /// Returns \ref startTime
    double getStartTime() const
    {
      return startTime;
    }

    /// Returns \ref endTime
    double getEndTime() const
    {
      return endTime;
    }

    /// Returns \ref timeErrLow.
479
    double getTimeErrLow(std::string key) const
480
    {
481
      return scalContents.at(key)->timeErrLow;
482
483
484
    }

    /// Returns whether coarsening is allowed or not.
485
    bool isCoarseningAllowed(std::string key) const
486
    {
487
      return (scalContents.at(key)->coarsenAllowed == 1);
488
489
490
    }

    /// Returns whether coarsening is allowed or not.
491
    bool isRefinementAllowed(std::string key) const
492
    {
493
      return (scalContents.at(key)->refinementAllowed == 1);
494
495
496
    }

    ///
497
    void allowRefinement(bool allow, std::string key)
498
    {
499
      scalContents.at(key)->refinementAllowed = allow;
500
501
502
    }

    ///
503
    void allowCoarsening(bool allow, std::string key)
504
    {
505
      scalContents.at(key)->coarsenAllowed = allow;
506
507
508
509
    }

    int getSize() const
    {
510
      return int(scalContents.size());
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
    }

    // TODO: remove from AdaptInfo
    bool getRosenbrockMode() const
    {
      return rosenbrockMode;
    }

    void setSolverIterations(int it)
    {
      solverIterations = it;
    }

    int getSolverIterations() const
    {
      return solverIterations;
    }

    void setMaxSolverIterations(int it)
    {
      maxSolverIterations = it;
    }

    int getMaxSolverIterations() const
    {
      return maxSolverIterations;
    }

    void setSolverTolerance(double tol)
    {
      solverTolerance = tol;
    }

    double getSolverTolerance() const
    {
      return solverTolerance;
    }

    void setSolverResidual(double res)
    {
      solverResidual = res;
    }

    double getSolverResidual() const
    {
      return solverResidual;
    }

    void setGlobalTimeTolerance(double tol)
    {
      globalTimeTolerance = tol;
    }

    double getGlobalTimeTolerance() const
    {
      return globalTimeTolerance;
    }

    // TODO: remove from AdaptInfo
    void setRosenbrockMode(bool b)
    {
      rosenbrockMode = b;
    }

    /// Creates new scalContents with the given size.
576
    void setScalContents(List names);
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626

    /** \brief
     * Resets timestep, current time and time boundaries without
     * any check. Is used by the parareal algorithm.
     */
    void resetTimeValues(double newTimeStep,
                         double newStartTime,
                         double newEndTime)
    {
      time = newStartTime;
      startTime = newStartTime;
      endTime = newEndTime;
      timestep = newTimeStep;
      timestepNumber = 0;
    }

  protected:
    /// Name.
    std::string name;

    /// Current space iteration
    int spaceIteration = -1;

    /** \brief
     * maximal allowed number of iterations of the adaptive procedure; if
     * maxIteration <= 0, no iteration bound is used
     */
    int maxSpaceIteration = -1;

    /// Current timestep iteration
    int timestepIteration = 0;

    /// Maximal number of iterations for choosing a timestep
    int maxTimestepIteration = 30;

    /// Current time iteration
    int timeIteration = 0;

    /// Maximal number of time iterations
    int maxTimeIteration = 30;

    /// Actual time, end of time interval for current time step
    double time = 0.0;

    /// Initial time
    double startTime = 0.0;

    /// Final time
    double endTime = 1.0;

627
    /// Time step size to be used
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
    double timestep = 0.0;

    /// Last processed time step size of finished iteration
    double lastProcessedTimestep = 0.0;

    /// Minimal step size
    double minTimestep = 0.0;

    /// Maximal step size
    double maxTimestep = 1.0;

    /// Number of current time step
    int timestepNumber = 0;

    /** \brief
     * Per default this value is 0 and not used. If it is set to a non-zero value,
     * the computation of the stationary problem is done nTimesteps times with a
     * fixed timestep.
     */
    int nTimesteps = 0;

    /// number of iterations needed of linear or nonlinear solver
    int solverIterations = 0;

    /// maximal number of iterations needed of linear or nonlinear solver
    int maxSolverIterations = 0;

    ///
    double solverTolerance = 1.e-8;

    ///
    double solverResidual = 0.0;

    /// tolerance for the overall time error
    double globalTimeTolerance = 1.0;

    /// Scalar adapt infos.
665
    std::map<std::string, std::unique_ptr<ScalContent>> scalContents;
666
667
668
669
670
671
672
673
674
675
676
677

    /// Is true, if the adaptive procedure was deserialized from a file. TODO: remove deserialization
    bool deserialized = false;

    /// Is true, if the time adaption is controlled by a Rosenbrock Method.
    bool rosenbrockMode = false;

    /// overall time error estimate
    double est_t = 0.0;
  };

} // end namespace AMDiS