ExaDG
Loading...
Searching...
No Matches
time_int_abm_base.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2023 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_ABM_BASE_H_
23#define INCLUDE_EXADG_TIME_INTEGRATION_TIME_INT_ABM_BASE_H_
24
25
26#include <exadg/time_integration/ab_constants.h>
27#include <exadg/time_integration/am_constants.h>
28#include <exadg/time_integration/push_back_vectors.h>
29#include <exadg/time_integration/restart.h>
30#include <exadg/time_integration/time_int_multistep_base.h>
31#include <exadg/utilities/print_solver_results.h>
32
33namespace ExaDG
34{
38template<typename Operator, typename VectorType>
39class TimeIntAdamsBashforthMoultonBase : public TimeIntMultistepBase
40{
41 using Number = typename VectorType::value_type;
42 using BoostInputArchiveType = TimeIntBase::BoostInputArchiveType;
43 using BoostOutputArchiveType = TimeIntBase::BoostOutputArchiveType;
44
45public:
46 TimeIntAdamsBashforthMoultonBase(std::shared_ptr<Operator> pde_operator_in,
47 double const start_time_,
48 double const end_time_,
49 unsigned int const max_number_of_time_steps_,
50 unsigned int const order_,
51 bool const start_with_low_order_,
52 bool const adaptive_time_stepping_,
53 RestartData const & restart_data_,
54 MPI_Comm const & mpi_comm_,
55 bool const is_test_)
56 : TimeIntMultistepBase(start_time_,
57 end_time_,
58 max_number_of_time_steps_,
59 order_,
60 start_with_low_order_,
61 adaptive_time_stepping_,
62 restart_data_,
63 mpi_comm_,
64 is_test_),
65 pde_operator(pde_operator_in),
66 // order of predictor can be chosen one below order of corrector
67 ab(order_ - 1, start_with_low_order_),
68 am(order_, start_with_low_order_),
69 vec_evaluated_operators(order_ - 1)
70 {
71 AssertThrow(order_ >= 1,
72 dealii::ExcMessage("Oder of ABM time integrator has to be at least 1."));
73 }
74
75 void
76 print_iterations() const
77 {
78 // explicit time integration -> no iterations
79 print_list_of_iterations(pcout, {"Adams-Bashforth-Moulton"}, {0});
80 }
81
82 void
83 ale_update()
84 {
85 AssertThrow(false, dealii::ExcMessage("not yet implemented"));
86 }
87
88 VectorType const &
89 get_solution() const
90 {
91 return solution;
92 }
93
94protected:
95 Operator const &
96 get_underlying_operator() const
97 {
98 return *pde_operator;
99 }
100
101private:
102 void
103 update_time_integrator_constants() final
104 {
105 ab.update(time_step_number, adaptive_time_stepping, time_steps);
106 am.update(time_step_number, adaptive_time_stepping, time_steps);
107 }
108
109 void
110 allocate_vectors() final
111 {
112 pde_operator->initialize_dof_vector(solution);
113 pde_operator->initialize_dof_vector(prediction);
114
115 pde_operator->initialize_dof_vector(evaluated_operator_np);
116 for(auto & evaluated_operator : vec_evaluated_operators)
117 pde_operator->initialize_dof_vector(evaluated_operator);
118 }
119
120 void
121 initialize_current_solution() final
122 {
123 pde_operator->prescribe_initial_conditions(solution, get_time());
124 }
125
126 void
127 initialize_former_multistep_dof_vectors() final
128 {
129 if(start_with_low_order)
130 {
131 if(vec_evaluated_operators.size() > 0)
132 pde_operator->evaluate(vec_evaluated_operators[0], solution, get_time());
133 }
134 else // start with high order
135 {
136 // fill evaluated operators
137 VectorType temp_sol;
138 pde_operator->initialize_dof_vector(temp_sol);
139 for(unsigned int i = 0; i < vec_evaluated_operators.size(); ++i)
140 {
141 pde_operator->prescribe_initial_conditions(temp_sol, get_previous_time(i));
142 pde_operator->evaluate(vec_evaluated_operators[i], temp_sol, get_previous_time(i));
143 }
144 }
145 }
146
147 void
148 setup_derived() final
149 {
150 }
151
152 void
153 do_timestep_predict()
154 {
155 dealii::Timer timer;
156 timer.restart();
157
158 predrict_solution(prediction, solution, vec_evaluated_operators);
159
160 timer_tree->insert({"Timeloop", "Adams-Bashforth-Moulton"}, timer.wall_time());
161 }
162
163 void
164 do_timestep_correct()
165 {
166 dealii::Timer timer;
167 timer.restart();
168
169 // evaluate operator given the predicted solution
170 pde_operator->evaluate(evaluated_operator_np, prediction, get_next_time());
171 // correct solution
172 correct_solution(solution, evaluated_operator_np, vec_evaluated_operators);
173 // correct operator by evaluating operator with correct solution
174 pde_operator->evaluate(evaluated_operator_np, solution, get_next_time());
175
176 // write output
177 if(this->print_solver_info() and not(this->is_test))
178 {
179 pcout << std::endl << "Adams-Bashforth-Moulton:";
180 print_wall_time(pcout, timer.wall_time());
181 }
182
183 timer_tree->insert({"Timeloop", "Adams-Bashforth-Moulton"}, timer.wall_time());
184 }
185
186 void
187 do_timestep_solve() final
188 {
189 do_timestep_predict();
190 do_timestep_correct();
191 }
192
193 void
194 correct_solution(VectorType & dst,
195 VectorType const & op_np,
196 std::vector<VectorType> const & ops) const
197 {
198 dst.add(static_cast<Number>(get_time_step_size() * am.get_gamma0()), op_np);
199 for(unsigned int i = 0; i < this->am.get_order() - 1; ++i)
200 dst.add(static_cast<Number>(get_time_step_size() * am.get_alpha(i)), ops[i]);
201 }
202
203 void
204 predrict_solution(VectorType & dst,
205 VectorType const & src,
206 std::vector<VectorType> const & ops) const
207 {
208 dst = src;
209 for(unsigned int i = 0; i < this->ab.get_order(); ++i)
210 dst.add(static_cast<Number>(get_time_step_size() * ab.get_alpha(i)), ops[i]);
211 }
212
213 void
214 prepare_vectors_for_next_timestep() final
215 {
216 if(vec_evaluated_operators.size() > 0)
217 {
218 push_back(vec_evaluated_operators);
219 std::swap(vec_evaluated_operators[0], evaluated_operator_np);
220 }
221 }
222
223 void
224 read_restart_vectors(BoostInputArchiveType & ia) final
225 {
226 read_write_distributed_vector(solution, ia);
227 read_write_distributed_vector(prediction, ia);
228
229 for(unsigned int i = 0; i < vec_evaluated_operators.size(); ++i)
230 {
231 read_write_distributed_vector(vec_evaluated_operators[i], ia);
232 }
233 }
234
235 void
236 write_restart_vectors(BoostOutputArchiveType & oa) const final
237 {
238 read_write_distributed_vector(solution, oa);
239 read_write_distributed_vector(prediction, oa);
240
241 for(unsigned int i = 0; i < vec_evaluated_operators.size(); ++i)
242 {
243 read_write_distributed_vector(vec_evaluated_operators[i], oa);
244 }
245 }
246
247 void
248 solve_steady_problem() final
249 {
250 AssertThrow(false, dealii::ExcMessage("Steady not implemented."));
251 }
252
253 // spatial pde operator
254 std::shared_ptr<Operator> pde_operator;
255
256 // Time integration constants
259
260 // solution vector
261 VectorType solution;
262 // temporary vector to store prediction
263 VectorType prediction;
264
265 // store evaluated operators from previous time steps
266 VectorType evaluated_operator_np;
267 std::vector<VectorType> vec_evaluated_operators;
268};
269
270} // namespace ExaDG
271
272#endif /* INCLUDE_EXADG_TIME_INTEGRATION_TIME_INT_ABM_BASE_H_*/
Definition ab_constants.h:40
Definition am_constants.h:40
Definition driver.cpp:33
void read_write_distributed_vector(VectorType &vector, BoostArchiveType &archive)
Definition restart.h:97
Definition restart_data.h:38