ryujin 2.1.1 revision 6dc06e5864abd5d99e5d7ab641dbe621936411d9
parabolic_module.template.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2023 by the ryujin authors
4//
5
6#pragma once
7
8#include "parabolic_module.h"
9
10namespace ryujin
11{
12 using namespace dealii;
13
14 template <typename Description, int dim, typename Number>
16 const MPIEnsemble &mpi_ensemble,
17 std::map<std::string, dealii::Timer> &computing_timer,
18 const OfflineData<dim, Number> &offline_data,
19 const HyperbolicSystem &hyperbolic_system,
20 const ParabolicSystem &parabolic_system,
21 const InitialValues<Description, dim, Number> &initial_values,
22 const std::string &subsection /*= "ParabolicModule"*/)
23 : ParameterAcceptor(subsection)
24 , id_violation_strategy_(IDViolationStrategy::warn)
25 , parabolic_solver_(mpi_ensemble,
26 computing_timer,
27 hyperbolic_system,
28 parabolic_system,
29 offline_data,
30 initial_values,
31 subsection)
32 , n_restarts_(0)
33 , n_corrections_(0)
34 , n_warnings_(0)
35 {
36 }
37
38
39 template <typename Description, int dim, typename Number>
41 {
42#ifdef DEBUG_OUTPUT
43 std::cout << "ParabolicModule<Description, dim, Number>::prepare()"
44 << std::endl;
45#endif
46 if constexpr (!ParabolicSystem::is_identity)
47 parabolic_solver_.prepare();
48
49 cycle_ = 0;
50 }
51
52
53 template <typename Description, int dim, typename Number>
54 template <int stages>
56 const StateVector &old_state_vector,
57 const Number old_t,
58 std::array<std::reference_wrapper<const StateVector>,
59 stages> /*stage_state_vectors*/,
60 const std::array<Number, stages> /*stage_weights*/,
61 StateVector &new_state_vector,
62 Number tau) const
63 {
64 if constexpr (ParabolicSystem::is_identity) {
65 AssertThrow(
66 false,
67 dealii::ExcMessage("The parabolic system is the identity. This "
68 "function should have never been called."));
69 __builtin_trap();
70
71 } else {
72
73 AssertThrow(stages == 0,
74 dealii::ExcMessage("Although IMEX schemes are implemented, "
75 "the high order fluxes are not. "));
77 const bool reinit_gmg = cycle_++ % 4 == 0;
78 parabolic_solver_.backward_euler_step(old_state_vector,
79 old_t,
80 new_state_vector,
81 tau,
82 id_violation_strategy_,
83 reinit_gmg);
84
85 /* Update statistics: */
86 n_restarts_ = parabolic_solver_.n_restarts();
87 n_corrections_ = parabolic_solver_.n_corrections();
88 n_warnings_ = parabolic_solver_.n_warnings();
89 }
90 }
91
92
93 template <typename Description, int dim, typename Number>
95 const StateVector &old_state_vector,
96 const Number old_t,
97 StateVector &new_state_vector,
98 Number tau) const
99 {
100 if constexpr (ParabolicSystem::is_identity) {
101 AssertThrow(
102 false,
103 dealii::ExcMessage("The parabolic system is the identity. This "
104 "function should have never been called."));
105 __builtin_trap();
106
107 } else {
108
109 const bool reinit_gmg = cycle_++ % 4 == 0;
110 parabolic_solver_.crank_nicolson_step(old_state_vector,
111 old_t,
112 new_state_vector,
113 tau,
114 id_violation_strategy_,
115 reinit_gmg);
116
117 /* Update statistics: */
118 n_restarts_ = parabolic_solver_.n_restarts();
119 n_corrections_ = parabolic_solver_.n_corrections();
120 n_warnings_ = parabolic_solver_.n_warnings();
122 }
123
124
125 template <typename Description, int dim, typename Number>
127 std::ostream &output) const
128 {
129 if constexpr (!ParabolicSystem::is_identity) {
130 parabolic_solver_.print_solver_statistics(output);
131 }
132 }
133
134} /* namespace ryujin */
typename View::StateVector StateVector
void crank_nicolson_step(const StateVector &old_state_vector, const Number old_t, StateVector &new_state_vector, Number tau) const
typename Description::ParabolicSystem ParabolicSystem
ParabolicModule(const MPIEnsemble &mpi_ensemble, std::map< std::string, dealii::Timer > &computing_timer, const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system, const ParabolicSystem &parabolic_system, const InitialValues< Description, dim, Number > &initial_values, const std::string &subsection="/ParabolicModule")
typename Description::HyperbolicSystem HyperbolicSystem
void print_solver_statistics(std::ostream &output) const
void backward_euler_step(const StateVector &old_state_vector, const Number old_t, std::array< std::reference_wrapper< const StateVector >, stages > stage_state_vectors, const std::array< Number, stages > stage_weights, StateVector &new_state_vector, Number tau) const