ExaDG
Toggle main menu visibility
Loading...
Searching...
No Matches
include
exadg
time_integration
time_int_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_INT_BASE_H_
23
#define EXADG_TIME_INTEGRATION_TIME_INT_BASE_H_
24
25
// C/C++
26
#include <boost/archive/binary_iarchive.hpp>
27
#include <boost/archive/binary_oarchive.hpp>
28
#include <boost/archive/text_iarchive.hpp>
29
#include <boost/archive/text_oarchive.hpp>
30
31
#include <fstream>
32
#include <sstream>
33
34
// deal.II
35
#include <deal.II/base/conditional_ostream.h>
36
#include <deal.II/base/timer.h>
37
38
// ExaDG
39
#include <exadg/time_integration/restart.h>
40
#include <exadg/time_integration/restart_data.h>
41
#include <exadg/utilities/numbers.h>
42
#include <exadg/utilities/timer_tree.h>
43
44
45
namespace
ExaDG
46
{
47
class
TimeIntBase
48
{
49
public
:
50
// Archive type used for serialization. Binary type for reduced file size.
51
typedef
boost::archive::binary_iarchive BoostInputArchiveType;
52
typedef
boost::archive::binary_oarchive BoostOutputArchiveType;
53
// Alternative human-readable type for debugging.
54
// typedef boost::archive::text_iarchive BoostInputArchiveType;
55
// typedef boost::archive::text_oarchive BoostOutputArchiveType;
56
57
TimeIntBase(
double
const
& start_time_,
58
double
const
& end_time_,
59
unsigned
int
const
max_number_of_time_steps_,
60
RestartData
const
& restart_data_,
61
MPI_Comm
const
& mpi_comm_,
62
bool
const
is_test_);
63
64
virtual
~TimeIntBase()
65
{
66
}
67
68
/*
69
* Setup of time integration scheme.
70
*/
71
virtual
void
72
setup(
bool
const
do_restart) = 0;
73
74
/*
75
* Returns true if the start time has been reached.
76
*/
77
bool
78
started()
const
;
79
80
/*
81
* returns true if the end of time loop has been reached or the maximum number of time steps
82
*/
83
bool
84
finished()
const
;
85
86
/*
87
* Performs the time loop from start_time to end_time by repeatingly calling
88
* advance_one_timestep().
89
*/
90
void
91
timeloop();
92
93
/*
94
* Perform only one time step (which is used when coupling different solvers, equations, etc.).
95
*/
96
void
97
advance_one_timestep();
98
99
/*
100
* The main sub-routines of advance_one_timestep()
101
*/
102
void
103
advance_one_timestep_pre_solve(
bool
const
print_header);
104
105
void
106
advance_one_timestep_solve();
107
108
void
109
advance_one_timestep_post_solve();
110
111
/*
112
* Reset the current time.
113
*/
114
void
115
reset_time(
double
const
& current_time);
116
117
/*
118
* In case of adaptive mesh refinement, the driver requests preparing the owned vectors
119
* for refinement and interpolation afterwards.
120
*/
121
virtual
void
122
prepare_coarsening_and_refinement();
123
124
virtual
void
125
interpolate_after_coarsening_and_refinement();
126
127
/*
128
* Get the time step size.
129
*/
130
virtual
double
131
get_time_step_size()
const
= 0;
132
133
/*
134
* Set the time step size.
135
*/
136
virtual
void
137
set_current_time_step_size(
double
const
& time_step_size) = 0;
138
139
/*
140
* Get the current time t_{n}.
141
*/
142
double
143
get_time()
const
;
144
145
/*
146
* Get time at the end of the current time step t_{n+1}.
147
*/
148
double
149
get_next_time()
const
;
150
151
/*
152
* Get number of computed time steps
153
*/
154
unsigned
int
155
get_number_of_time_steps()
const
;
156
157
std::shared_ptr<TimerTree>
158
get_timings()
const
;
159
160
protected
:
161
/*
162
* Do one time step including pre and post routines done before and after the actual solution of
163
* the current time step. Compared to the function advance_one_timestep(), do_timestep() is a raw
164
* version that does not call postprocessing routines, does not write output to pcout, and does
165
* not perform timer measurements within its sub-routines. The typical use case of do_timestep()
166
* is when using pseudo-timestepping to obtain the solution of a steady-state problem with a
167
* transient solver.
168
*/
169
void
170
do_timestep();
171
172
/*
173
* e.g., update of time integrator constants
174
*/
175
virtual
void
176
do_timestep_pre_solve(
bool
const
print_header) = 0;
177
178
/*
179
* The actual solution of the current time step
180
*/
181
virtual
void
182
do_timestep_solve() = 0;
183
184
/*
185
* e.g., update of DoF vectors, increment time, adjust time step size, etc.
186
*/
187
virtual
void
188
do_timestep_post_solve() = 0;
189
190
/*
191
* Postprocessing of solution.
192
*/
193
virtual
void
194
postprocessing()
const
= 0;
195
196
/*
197
* Get the current time step number.
198
*/
199
types::time_step
200
get_time_step_number()
const
;
201
202
/*
203
* Write solution vectors to files so that the simulation can be restart from an intermediate
204
* state. Note that the sequence of writing and reading data in `write_restart` and `read_restart`
205
* needs to be identical.
206
*/
207
void
208
write_restart()
const
;
209
210
/*
211
* Read all relevant data from restart files to start the time integrator. Note that the sequence
212
* of writing and reading data in `write_restart` and `read_restart` needs to be identical.
213
*/
214
void
215
read_restart();
216
217
/*
218
* Output solver information before solving the time step.
219
*/
220
void
221
output_solver_info_header()
const
;
222
223
224
/*
225
* Output estimated computation time until completion of the simulation.
226
*/
227
void
228
output_remaining_time()
const
;
229
230
/*
231
* Start and end times.
232
*/
233
double
start_time, end_time;
234
235
/*
236
* Physical time.
237
*/
238
double
time;
239
240
/*
241
* A small number which is much smaller than the time step size.
242
*/
243
double
const
eps;
244
245
/*
246
* Output to screen.
247
*/
248
dealii::ConditionalOStream pcout;
249
250
/*
251
* The number of the current time step starting with time_step_number = 1.
252
*/
253
types::time_step time_step_number;
254
255
/*
256
* Maximum number of time steps.
257
*/
258
unsigned
int
const
max_number_of_time_steps;
259
260
/*
261
* Restart.
262
*/
263
RestartData
const
restart_data;
264
265
/*
266
* MPI communicator.
267
*/
268
MPI_Comm
const
mpi_comm;
269
270
/*
271
* Computation time (wall clock time).
272
*/
273
dealii::Timer global_timer;
274
std::shared_ptr<TimerTree> timer_tree;
275
bool
is_test;
276
277
private
:
278
/*
279
* Write restart data.
280
*/
281
virtual
void
282
do_write_restart(std::string
const
& filename)
const
= 0;
283
284
/*
285
* Read restart data.
286
*/
287
virtual
void
288
do_read_restart(std::ifstream & in) = 0;
289
};
290
291
}
// namespace ExaDG
292
293
#endif
/* EXADG_TIME_INTEGRATION_TIME_INT_BASE_H_ */
ExaDG
Definition
driver.cpp:33
ExaDG::RestartData
Definition
restart_data.h:84
Generated by
1.17.0