ryujin 2.1.1 revision 6dc06e5864abd5d99e5d7ab641dbe621936411d9
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 "mpi_ensemble.h"
13#include "offline_data.h"
14#include "sparse_matrix_simd.h"
15#include "state_vector.h"
16
17#include <deal.II/base/parameter_acceptor.h>
18#include <deal.II/base/smartpointer.h>
19#include <deal.II/base/timer.h>
20#include <deal.II/lac/sparse_matrix.templates.h>
21#include <deal.II/lac/vector.h>
22
23#include <functional>
24
25namespace ryujin
26{
43 enum class IDViolationStrategy : std::uint8_t {
48 warn,
49
56 };
57
58
72 class Restart final
73 {
74 };
75
76
90 class Correction final
91 {
92 };
93
94
103 template <typename Description, int dim, typename Number = double>
104 class HyperbolicModule final : public dealii::ParameterAcceptor
105 {
106 public:
111
113
114 using View =
115 typename Description::template HyperbolicSystemView<dim, Number>;
116
117 static constexpr auto problem_dimension = View::problem_dimension;
118
119 using state_type = typename View::state_type;
120
121 using precomputed_type = typename View::precomputed_type;
122
123 using initial_precomputed_type = typename View::initial_precomputed_type;
124
126
127 using InitialPrecomputedVector = typename View::InitialPrecomputedVector;
128
129 static constexpr auto n_precomputation_cycles =
130 View::n_precomputation_cycles;
131
133
137
142 const MPIEnsemble &mpi_ensemble,
143 std::map<std::string, dealii::Timer> &computing_timer,
146 const InitialValues<Description, dim, Number> &initial_values,
147 const std::string &subsection = "/HyperbolicModule");
148
154 void prepare();
155
157
161
182 void prepare_state_vector(StateVector &state_vector, Number t) const;
183
246 template <int stages>
247 Number step(
248 const StateVector &old_state_vector,
249 std::array<std::reference_wrapper<const StateVector>, stages>
250 stage_state_vectors,
251 const std::array<Number, stages> stage_weights,
252 StateVector &new_state_vector,
253 Number tau = Number(0.),
254 std::atomic<Number> tau_max = std::numeric_limits<Number>::max()) const;
255
263 void cfl(Number new_cfl) const
264 {
265 Assert(cfl_ > Number(0.), dealii::ExcInternalError());
266 cfl_ = new_cfl;
267 }
268
274
275
279
284
289
296
303
310
317
318 // FIXME: refactor to function
320
321 private:
323
327 typename Description::template Indicator<dim, Number>::Parameters
328 indicator_parameters_;
329
330 typename Description::template Limiter<dim, Number>::Parameters
331 limiter_parameters_;
332
333 typename Description::template RiemannSolver<dim, Number>::Parameters
334 riemann_solver_parameters_;
335
337
339
343
344 const MPIEnsemble &mpi_ensemble_;
345 std::map<std::string, dealii::Timer> &computing_timer_;
346
347 dealii::SmartPointer<const OfflineData<dim, Number>> offline_data_;
348 dealii::SmartPointer<const HyperbolicSystem> hyperbolic_system_;
349 dealii::SmartPointer<const InitialValues<Description, dim, Number>>
350 initial_values_;
351
352 mutable Number cfl_;
353
354 mutable unsigned int n_restarts_;
355 mutable unsigned int n_corrections_;
356 mutable unsigned int n_warnings_;
357
358 InitialPrecomputedVector initial_precomputed_;
359
360 using ScalarVector = typename Vectors::ScalarVector<Number>;
361 mutable ScalarVector alpha_;
362
363 static constexpr auto n_bounds =
364 Description::template Limiter<dim, Number>::n_bounds;
365 mutable Vectors::MultiComponentVector<Number, n_bounds> bounds_;
366
367 using HyperbolicVector =
368 Vectors::MultiComponentVector<Number, problem_dimension>;
369 mutable HyperbolicVector r_;
370
371 mutable SparseMatrixSIMD<Number> dij_matrix_;
372 mutable SparseMatrixSIMD<Number> lij_matrix_;
373 mutable SparseMatrixSIMD<Number> lij_matrix_next_;
374 mutable SparseMatrixSIMD<Number, problem_dimension> pij_matrix_;
375
377 };
378
379} /* namespace ryujin */
const auto & n_warnings() const
HyperbolicModule(const MPIEnsemble &mpi_ensemble, 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")
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
typename View::precomputed_type precomputed_type
typename View::state_type state_type
const auto & cfl() const
const auto & alpha() const
typename Description::template HyperbolicSystemView< dim, Number > View
const auto & n_corrections() const
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
const auto & hyperbolic_system() const
void cfl(Number new_cfl) const
IDViolationStrategy id_violation_strategy_
const auto & offline_data() const
void prepare_state_vector(StateVector &state_vector, Number t) const
const auto & n_restarts() const
const auto & initial_precomputed() 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