ryujin 2.1.1 revision 0eab90fbc6e1ac9f2e0a2e6d16f9f023c13a02f7
parabolic_module.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2023 - 2024 by the ryujin authors
4//
5
6#pragma once
7
8#include "hyperbolic_module.h"
9
10#include <compile_time_options.h>
11
12#include "convenience_macros.h"
13#include "initial_values.h"
14#include "mpi_ensemble.h"
15#include "offline_data.h"
16
17#include <deal.II/base/parameter_acceptor.h>
18#include <deal.II/base/timer.h>
19#include <deal.II/lac/sparse_matrix.templates.h>
20#include <deal.II/lac/vector.h>
21
22#include <functional>
23
24namespace ryujin
25{
32 template <typename Description, int dim, typename Number = double>
33 class ParabolicModule final : public dealii::ParameterAcceptor
34 {
35 public:
40
42
43 using View =
44 typename Description::template HyperbolicSystemView<dim, Number>;
45
47
49 typename Description::template ParabolicSolver<dim, Number>;
50
52
54
58
63 const MPIEnsemble &mpi_ensemble,
64 std::map<std::string, dealii::Timer> &computing_timer,
65 const OfflineData<dim, Number> &offline_data,
66 const HyperbolicSystem &hyperbolic_system,
67 const ParabolicSystem &parabolic_system,
68 const InitialValues<Description, dim, Number> &initial_values,
69 const std::string &subsection = "/ParabolicModule");
70
76 void prepare();
77
79
83
91 void prepare_state_vector(StateVector &state_vector, Number t) const;
92
102 template <int stages>
103 void
104 backward_euler_step(const StateVector &old_state_vector,
105 const Number old_t,
106 std::array<std::reference_wrapper<const StateVector>,
107 stages> stage_state_vectors,
108 const std::array<Number, stages> stage_weights,
109 StateVector &new_state_vector,
110 Number tau) const;
111
120 void crank_nicolson_step(const StateVector &old_state_vector,
121 const Number old_t,
122 StateVector &new_state_vector,
123 Number tau) const;
124
130 void print_solver_statistics(std::ostream &output) const;
131
133
141
142
143
149
155
156 // FIXME: refactor to function
158
159 private:
161
165
166 // FIXME: refactor contents into this class.
167 ParabolicSolver parabolic_solver_;
168 mutable unsigned int cycle_;
169
170 mutable unsigned int n_restarts_;
171 mutable unsigned int n_corrections_;
172 mutable unsigned int n_warnings_;
173
175 };
176
177} /* namespace ryujin */
void prepare_state_vector(StateVector &state_vector, Number t) const
typename View::StateVector StateVector
void crank_nicolson_step(const StateVector &old_state_vector, const Number old_t, StateVector &new_state_vector, Number tau) const
const auto & n_warnings() const
const auto & n_corrections() const
IDViolationStrategy id_violation_strategy_
typename Description::ParabolicSystem ParabolicSystem
const auto & n_restarts() const
typename Description::template ParabolicSolver< dim, Number > ParabolicSolver
typename Description::template HyperbolicSystemView< dim, Number > View
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
#define ACCESSOR_READ_ONLY(member)
std::tuple< MultiComponentVector< Number, problem_dim >, MultiComponentVector< Number, prec_dim >, BlockVector< Number > > StateVector
Definition: state_vector.h:51
Euler::HyperbolicSystem HyperbolicSystem
Definition: description.h:32
ryujin::StubParabolicSystem ParabolicSystem
Definition: description.h:37