ryujin 2.1.1 revision 1c453cc82f1d29edf537280cd96267402ac73e60
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 MPI_Comm &mpi_communicator,
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_communicator,
26 computing_timer,
27 hyperbolic_system,
28 parabolic_system,
29 offline_data,
30 initial_values,
31 subsection)
32 , n_restarts_(0)
33 , n_warnings_(0)
34 {
35 }
36
37
38 template <typename Description, int dim, typename Number>
40 {
41#ifdef DEBUG_OUTPUT
42 std::cout << "ParabolicModule<Description, dim, Number>::prepare()"
43 << std::endl;
44#endif
45 if constexpr (!ParabolicSystem::is_identity)
46 parabolic_solver_.prepare();
47
48 cycle_ = 0;
49 }
50
51
52 template <typename Description, int dim, typename Number>
53 template <int stages>
55 const StateVector &old_state_vector,
56 const Number old_t,
57 std::array<std::reference_wrapper<const StateVector>,
58 stages> /*stage_state_vectors*/,
59 const std::array<Number, stages> /*stage_weights*/,
60 StateVector &new_state_vector,
61 Number tau) const
62 {
63 if constexpr (ParabolicSystem::is_identity) {
64 AssertThrow(
65 false,
66 dealii::ExcMessage("The parabolic system is the identity. This "
67 "function should have never been called."));
68 __builtin_trap();
69
70 } else {
71
72 static_assert(stages == 0, "high order fluxes are not implemented");
73
74 const bool reinit_gmg = cycle_++ % 4 == 0;
75 parabolic_solver_.backward_euler_step(old_state_vector,
76 old_t,
77 new_state_vector,
78 tau,
79 id_violation_strategy_,
80 reinit_gmg);
81 n_restarts_ = parabolic_solver_.n_restarts();
82 n_warnings_ = parabolic_solver_.n_warnings();
83 }
84 }
85
86
87 template <typename Description, int dim, typename Number>
89 std::ostream &output) const
90 {
91 if constexpr (!ParabolicSystem::is_identity) {
92 parabolic_solver_.print_solver_statistics(output);
93 }
94 }
95
96} /* namespace ryujin */
typename View::StateVector StateVector
ParabolicModule(const MPI_Comm &mpi_communicator, 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::ParabolicSystem ParabolicSystem
typename Description::HyperbolicSystem HyperbolicSystem
void 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
void print_solver_statistics(std::ostream &output) const