ryujin 2.1.1 revision 284c1ee75f38e11ba18808dd0d991d905cf7c040
time_integrator.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2022 - 2024 by the ryujin authors
4//
5
6#pragma once
7
8#include <compile_time_options.h>
9
10#include "convenience_macros.h"
11#include "hyperbolic_module.h"
12#include "offline_data.h"
13#include "parabolic_module.h"
14#include "patterns_conversion.h"
15
16namespace ryujin
17{
26 none,
27
35 };
36
37
41 enum class TimeSteppingScheme {
69
74 erk_11,
75
87 erk_22,
88
101 erk_33,
102
115 erk_43,
116
121 erk_54,
122
128
134
140 };
141} // namespace ryujin
142
143#ifndef DOXYGEN
144DECLARE_ENUM(ryujin::CFLRecoveryStrategy,
147 "bang bang control"}));
148
149DECLARE_ENUM(
151 LIST({ryujin::TimeSteppingScheme::ssprk_22, "ssprk 22"},
160 {ryujin::TimeSteppingScheme::strang_erk_43_cn, "strang erk 43 cn"}, ));
161#endif
162
163namespace ryujin
164{
171 template <typename Description, int dim, typename Number = double>
172 class TimeIntegrator final : public dealii::ParameterAcceptor
173 {
174 public:
179
184
188 using View =
189 typename Description::template HyperbolicSystemView<dim, Number>;
190
194 static constexpr unsigned int problem_dimension = View::problem_dimension;
195
199 static constexpr unsigned int n_precomputed_values =
200 View::n_precomputed_values;
201
206
211
216 const MPI_Comm &mpi_communicator,
217 std::map<std::string, dealii::Timer> &computing_timer,
218 const OfflineData<dim, Number> &offline_data,
219 const HyperbolicModule<Description, dim, Number> &hyperbolic_module,
220 const ParabolicModule<Description, dim, Number> &parabolic_module,
221 const std::string &subsection = "/TimeIntegrator");
222
228 void prepare();
229
234
247 Number step(vector_type &U, Number t);
248
253
261
262 protected:
273 Number step_ssprk_22(vector_type &U, Number t);
274
285 Number step_ssprk_33(vector_type &U, Number t);
286
292 Number step_erk_11(vector_type &U, Number t);
293
299 Number step_erk_22(vector_type &U, Number t);
300
306 Number step_erk_33(vector_type &U, Number t);
307
313 Number step_erk_43(vector_type &U, Number t);
314
321 Number step_erk_54(vector_type &U, Number t);
322
330 Number step_strang_ssprk_33_cn(vector_type &U, Number t);
331
339 Number step_strang_erk_33_cn(vector_type &U, Number t);
340
348 Number step_strang_erk_43_cn(vector_type &U, Number t);
349
350 private:
352
356
357 Number cfl_min_;
358 Number cfl_max_;
359
360 CFLRecoveryStrategy cfl_recovery_strategy_;
361
362 TimeSteppingScheme time_stepping_scheme_;
363 double efficiency_;
364
366
368
372
373 const MPI_Comm &mpi_communicator_;
374 std::map<std::string, dealii::Timer> &computing_timer_;
375
376 dealii::SmartPointer<const OfflineData<dim, Number>> offline_data_;
377 dealii::SmartPointer<const HyperbolicModule<Description, dim, Number>>
378 hyperbolic_module_;
379 dealii::SmartPointer<const ParabolicModule<Description, dim, Number>>
380 parabolic_module_;
381
382 std::vector<vector_type> U_;
383 std::vector<precomputed_type> precomputed_;
384
386 };
387
388} /* namespace ryujin */
static constexpr unsigned int problem_dimension
Number step_erk_22(vector_type &U, Number t)
typename Description::HyperbolicSystem HyperbolicSystem
TimeIntegrator(const MPI_Comm &mpi_communicator, std::map< std::string, dealii::Timer > &computing_timer, const OfflineData< dim, Number > &offline_data, const HyperbolicModule< Description, dim, Number > &hyperbolic_module, const ParabolicModule< Description, dim, Number > &parabolic_module, const std::string &subsection="/TimeIntegrator")
Number step_erk_11(vector_type &U, Number t)
Number step_erk_33(vector_type &U, Number t)
static constexpr unsigned int n_precomputed_values
Number step(vector_type &U, Number t)
typename Description::ParabolicSystem ParabolicSystem
auto & efficiency() const
Number step_strang_erk_43_cn(vector_type &U, Number t)
Number step_strang_ssprk_33_cn(vector_type &U, Number t)
auto & time_stepping_scheme() const
Number step_erk_54(vector_type &U, Number t)
typename Description::template HyperbolicSystemView< dim, Number > View
Number step_ssprk_22(vector_type &U, Number t)
Number step_erk_43(vector_type &U, Number t)
Number step_ssprk_33(vector_type &U, Number t)
Number step_strang_erk_33_cn(vector_type &U, Number t)
#define ACCESSOR_READ_ONLY(member)
Euler::ParabolicSystem ParabolicSystem
Definition: description.h:37
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:32