ExaDG
Loading...
Searching...
No Matches
time_integration_constants_base.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef EXADG_TIME_INTEGRATION_TIME_INTEGRATION_CONSTANTS_BASE_H_
23#define EXADG_TIME_INTEGRATION_TIME_INTEGRATION_CONSTANTS_BASE_H_
24
25// C/C++
26#include <vector>
27
28namespace ExaDG
29{
30class TimeIntegratorConstantsBase
31{
32public:
33 TimeIntegratorConstantsBase(unsigned int const order, bool const start_with_low_order)
34 : order(order), start_with_low_order(start_with_low_order)
35 {
36 }
37
38 virtual ~TimeIntegratorConstantsBase()
39 {
40 }
41
42 /*
43 * This function updates the time integrator constants. The argument time_steps is only used in
44 * case of adaptive time stepping.
45 */
46 void
47 update(unsigned int const current_order,
48 bool const adaptive_time_stepping,
49 std::vector<double> const & time_steps)
50 {
51 // when starting the time integrator with a low order method, ensure that
52 // the time integrator constants are set properly
53 unsigned int const update_order =
54 (current_order <= order and start_with_low_order == true) ? current_order : order;
55
56 if(adaptive_time_stepping)
57 set_adaptive_time_step(update_order, time_steps);
58 else
59 set_constant_time_step(update_order);
60 }
61
62 unsigned int
63 get_order() const
64 {
65 return order;
66 }
67
68 /*
69 * This function prints the time integrator constants
70 */
71 virtual void
72 print(dealii::ConditionalOStream & pcout) const = 0;
73
74protected:
80 void
81 disable_high_order_constants(unsigned int const current_order, std::vector<double> & constants)
82 {
83 for(unsigned int i = current_order; i < constants.size(); ++i)
84 constants[i] = 0.0;
85 }
86
87
88 // order of time integrator
89 unsigned int const order;
90
91 // use a low order time integration scheme to start the time integrator?
92 bool const start_with_low_order;
93
94private:
95 /*
96 * This function calculates the time integrator constants in case of constant time step sizes.
97 */
98 virtual void
99 set_constant_time_step(unsigned int const current_order) = 0;
100
101 /*
102 * This function calculates time integrator constants in case of varying time step sizes
103 * (adaptive time stepping).
104 */
105 virtual void
106 set_adaptive_time_step(unsigned int const current_order,
107 std::vector<double> const & time_steps) = 0;
108};
109} // namespace ExaDG
110
111#endif /* EXADG_TIME_INTEGRATION_TIME_INTEGRATION_CONSTANTS_BASE_H_ */
void disable_high_order_constants(unsigned int const current_order, std::vector< double > &constants)
Definition time_integration_constants_base.h:81
Definition driver.cpp:33