ExaDG
Loading...
Searching...
No Matches
time_int_multistep_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 INCLUDE_EXADG_TIME_INTEGRATION_TIME_INT_MULTISTEP_BASE_H_
23#define INCLUDE_EXADG_TIME_INTEGRATION_TIME_INT_MULTISTEP_BASE_H_
24
25// ExaDG
26#include <exadg/time_integration/time_int_base.h>
27
28namespace ExaDG
29{
30class TimeIntMultistepBase : public TimeIntBase
31{
32public:
33 using BoostInputArchiveType = TimeIntBase::BoostInputArchiveType;
34 using BoostOutputArchiveType = TimeIntBase::BoostOutputArchiveType;
35
36 /*
37 * Constructor.
38 */
39 TimeIntMultistepBase(double const start_time_,
40 double const end_time_,
41 unsigned int const max_number_of_time_steps_,
42 unsigned const order_,
43 bool const start_with_low_order_,
44 bool const adaptive_time_stepping_,
45 RestartData const & restart_data_,
46 MPI_Comm const & mpi_comm_,
47 bool const is_test_);
48
49 /*
50 * Destructor.
51 */
52 virtual ~TimeIntMultistepBase()
53 {
54 }
55
56 /*
57 * Setup function where allocations/initializations are done. Calls another function
58 * setup_derived() in which the setup of derived classes can be performed.
59 */
60 void
61 setup(bool const do_restart) final;
62
63 /*
64 * Pseudo-time-stepping for steady-state problems.
65 */
66 void
67 timeloop_steady_problem();
68
69 /*
70 * Get the time step size.
71 */
72 double
73 get_time_step_size() const final;
74
75 double
76 get_time_step_size(int const index) const;
77
78 /*
79 * Set the time step size. Note that the push-back of time step sizes in case of adaptive time
80 * stepping may not be done here because calling this public function several times would falsify
81 * the results. Hence, set_time_step_size() is only allowed to overwrite the current time step
82 * size.
83 */
84 void
85 set_current_time_step_size(double const & time_step_size) final;
86
87 /*
88 * Get time at the end of the current time step.
89 */
90 double
91 get_previous_time(int const i /* t_{n-i} */) const;
92
93protected:
94 /*
95 * Get current order of time integrator
96 */
97 unsigned int
98 get_current_order() const;
99
100 /*
101 * Do one time step including different updates before and after the actual solution of the
102 * current time step.
103 */
104 void
105 do_timestep_pre_solve(bool const print_header) final;
106
107 void
108 do_timestep_post_solve() final;
109
110 /*
111 * Update the time integrator constants.
112 */
113 virtual void
114 update_time_integrator_constants() = 0;
115
116 /*
117 * Get reference to vector with time step sizes
118 */
119 std::vector<double>
120 get_time_step_vector() const;
121
122 /*
123 * Update of time step sizes in case of variable time steps.
124 */
125 void
126 push_back_time_step_sizes();
127
128 /*
129 * Calculate time step size.
130 */
131 virtual double
132 calculate_time_step_size() = 0;
133
134 /*
135 * returns whether solver info has to be written in the current time step.
136 */
137 virtual bool
138 print_solver_info() const = 0;
139
140 /*
141 * Order of time integration scheme.
142 */
143 unsigned int const order;
144
145 /*
146 * Start with low order (1st order) time integration scheme in first time step.
147 */
148 bool const start_with_low_order;
149
150 /*
151 * Use adaptive time stepping?
152 */
153 bool const adaptive_time_stepping;
154
155 /*
156 * Vector with time step sizes.
157 */
158 std::vector<double> time_steps;
159
160private:
161 /*
162 * Allocate solution vectors (has to be implemented by derived classes).
163 */
164 virtual void
165 allocate_vectors() = 0;
166
167 /*
168 * Initializes the solution vectors by prescribing initial conditions or reading data from
169 * restart files and initializes the time step size.
170 */
171 virtual void
172 initialize_solution_and_time_step_size(bool do_restart);
173
174 /*
175 * Initializes the solution vectors at time t
176 */
177 virtual void
178 initialize_current_solution() = 0;
179
180 /*
181 * Initializes solutions (BDF) or evaluated operators (Adams) at time t - dt[1], t - dt[1] -
182 * dt[2], etc.
183 */
184 virtual void
185 initialize_former_multistep_dof_vectors() = 0;
186
187 /*
188 * Setup of derived classes.
189 */
190 virtual void
191 setup_derived() = 0;
192
193 /*
194 * This function prepares the solution vectors for the next time step, e.g., by switching pointers
195 * to the solution vectors (called push back here).
196 */
197 virtual void
198 prepare_vectors_for_next_timestep() = 0;
199
200 /*
201 * Solve for a steady-state solution using pseudo-time-stepping.
202 */
203 virtual void
204 solve_steady_problem();
205
206 /*
207 * Postprocessing of solution.
208 */
209 virtual void
210 postprocessing_steady_problem() const;
211
212 /*
213 * Restart: read solution vectors (has to be implemented in derived classes).
214 */
215 void
216 do_read_restart(std::ifstream & in) final;
217
218 void
219 read_restart_preamble(BoostInputArchiveType & ia);
220
221 virtual void
222 read_restart_vectors(BoostInputArchiveType & ia) = 0;
223
224 /*
225 * Write solution vectors to files so that the simulation can be restart from an intermediate
226 * state.
227 */
228 void
229 do_write_restart(std::string const & filename) const final;
230
231 void
232 write_restart_preamble(BoostOutputArchiveType & oa) const;
233
234 virtual void
235 write_restart_vectors(BoostOutputArchiveType & oa) const = 0;
236
237 /*
238 * Recalculate the time step size after each time step in case of adaptive time stepping.
239 */
240 virtual double
241 recalculate_time_step_size() const = 0;
242};
243
244} // namespace ExaDG
245
246#endif /* INCLUDE_EXADG_TIME_INTEGRATION_TIME_INT_MULTISTEP_BASE_H_ */
Definition driver.cpp:33
Definition restart_data.h:38