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