ryujin 2.1.1 revision 0b194b984a74af675d09b5e928529ca8c7b634f2
hyperbolic_module.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2020 - 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 "initial_values.h"
12#include "offline_data.h"
13#include "sparse_matrix_simd.h"
14#include "state_vector.h"
15
16#include <deal.II/base/parameter_acceptor.h>
17#include <deal.II/base/timer.h>
18#include <deal.II/lac/sparse_matrix.templates.h>
19#include <deal.II/lac/vector.h>
20
21#include <functional>
22
23namespace ryujin
24{
41 warn,
42
49 };
50
51
58 class Restart final
59 {
60 };
61
62
71 template <typename Description, int dim, typename Number = double>
72 class HyperbolicModule final : public dealii::ParameterAcceptor
73 {
74 public:
79
81
82 using View =
83 typename Description::template HyperbolicSystemView<dim, Number>;
84
85 static constexpr auto problem_dimension = View::problem_dimension;
86
87 using state_type = typename View::state_type;
88
89 using precomputed_type = typename View::precomputed_type;
90
91 using initial_precomputed_type = typename View::initial_precomputed_type;
92
94
95 using InitialPrecomputedVector = typename View::InitialPrecomputedVector;
96
97 static constexpr auto n_precomputation_cycles =
98 View::n_precomputation_cycles;
99
101
105
110 const MPI_Comm &mpi_communicator,
111 std::map<std::string, dealii::Timer> &computing_timer,
114 const InitialValues<Description, dim, Number> &initial_values,
115 const std::string &subsection = "/HyperbolicModule");
116
122 void prepare();
123
125
129
150 void prepare_state_vector(StateVector &state_vector, Number t) const;
151
214 template <int stages>
215 Number step(
216 const StateVector &old_state_vector,
217 std::array<std::reference_wrapper<const StateVector>, stages>
218 stage_state_vectors,
219 const std::array<Number, stages> stage_weights,
220 StateVector &new_state_vector,
221 Number tau = Number(0.),
222 std::atomic<Number> tau_max = std::numeric_limits<Number>::max()) const;
223
231 void cfl(Number new_cfl) const
232 {
233 Assert(cfl_ > Number(0.), dealii::ExcInternalError());
234 cfl_ = new_cfl;
235 }
236
242
243
247
252
257
264
269
275
276 // FIXME: refactor to function
278
279 private:
281
285 typename Description::template Indicator<dim, Number>::Parameters
286 indicator_parameters_;
287
288 typename Description::template Limiter<dim, Number>::Parameters
289 limiter_parameters_;
290
291 typename Description::template RiemannSolver<dim, Number>::Parameters
292 riemann_solver_parameters_;
293
295
297
301
302 const MPI_Comm &mpi_communicator_;
303 std::map<std::string, dealii::Timer> &computing_timer_;
304
305 dealii::SmartPointer<const OfflineData<dim, Number>> offline_data_;
306 dealii::SmartPointer<const HyperbolicSystem> hyperbolic_system_;
307 dealii::SmartPointer<const InitialValues<Description, dim, Number>>
308 initial_values_;
309
310 mutable Number cfl_;
311
312 mutable unsigned int n_restarts_;
313
314 mutable unsigned int n_warnings_;
315
316 InitialPrecomputedVector initial_precomputed_;
317
318 using ScalarVector = typename Vectors::ScalarVector<Number>;
319 mutable ScalarVector alpha_;
320
321 static constexpr auto n_bounds =
322 Description::template Limiter<dim, Number>::n_bounds;
323 mutable Vectors::MultiComponentVector<Number, n_bounds> bounds_;
324
325 using HyperbolicVector =
326 Vectors::MultiComponentVector<Number, problem_dimension>;
327 mutable HyperbolicVector r_;
328
329 mutable SparseMatrixSIMD<Number> dij_matrix_;
330 mutable SparseMatrixSIMD<Number> lij_matrix_;
331 mutable SparseMatrixSIMD<Number> lij_matrix_next_;
332 mutable SparseMatrixSIMD<Number, problem_dimension> pij_matrix_;
333
335 };
336
337} /* namespace ryujin */
auto & hyperbolic_system() const
static constexpr auto problem_dimension
typename View::initial_precomputed_type initial_precomputed_type
static constexpr auto n_precomputation_cycles
typename Description::HyperbolicSystem HyperbolicSystem
typename View::StateVector StateVector
HyperbolicModule(const MPI_Comm &mpi_communicator, std::map< std::string, dealii::Timer > &computing_timer, const OfflineData< dim, Number > &offline_data, const HyperbolicSystem &hyperbolic_system, const InitialValues< Description, dim, Number > &initial_values, const std::string &subsection="/HyperbolicModule")
typename View::precomputed_type precomputed_type
typename View::state_type state_type
typename Description::template HyperbolicSystemView< dim, Number > View
typename View::InitialPrecomputedVector InitialPrecomputedVector
Number step(const StateVector &old_state_vector, std::array< std::reference_wrapper< const StateVector >, stages > stage_state_vectors, const std::array< Number, stages > stage_weights, StateVector &new_state_vector, Number tau=Number(0.), std::atomic< Number > tau_max=std::numeric_limits< Number >::max()) const
void cfl(Number new_cfl) const
IDViolationStrategy id_violation_strategy_
auto & initial_precomputed() const
void prepare_state_vector(StateVector &state_vector, Number t) 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