ryujin 2.1.1 revision 0eab90fbc6e1ac9f2e0a2e6d16f9f023c13a02f7
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{
47 enum class IDViolationStrategy : std::uint8_t {
52 warn,
53
60 };
61
62
79 class Restart final
80 {
81 public:
83 };
84
85
99 class Correction final
100 {
101 };
102
103
112 template <typename Description, int dim, typename Number = double>
113 class HyperbolicModule final : public dealii::ParameterAcceptor
114 {
115 public:
120
122
123 using View =
124 typename Description::template HyperbolicSystemView<dim, Number>;
125
126 static constexpr auto problem_dimension = View::problem_dimension;
127
128 using state_type = typename View::state_type;
129
130 using precomputed_type = typename View::precomputed_type;
131
132 using initial_precomputed_type = typename View::initial_precomputed_type;
133
135
136 using InitialPrecomputedVector = typename View::InitialPrecomputedVector;
137
138 static constexpr auto n_precomputation_cycles =
139 View::n_precomputation_cycles;
140
142
146
151 const MPIEnsemble &mpi_ensemble,
152 std::map<std::string, dealii::Timer> &computing_timer,
155 const InitialValues<Description, dim, Number> &initial_values,
156 const std::string &subsection = "/HyperbolicModule");
157
163 void prepare();
164
166
170
186 void prepare_state_vector(StateVector &state_vector, Number t) const;
187
250 template <int stages>
251 Number step(
252 const StateVector &old_state_vector,
253 std::array<std::reference_wrapper<const StateVector>, stages>
254 stage_state_vectors,
255 const std::array<Number, stages> stage_weights,
256 StateVector &new_state_vector,
257 Number tau = Number(0.),
258 std::atomic<Number> tau_max = std::numeric_limits<Number>::max()) const;
259
267 void cfl(Number new_cfl) const
268 {
269 Assert(cfl_ > Number(0.), dealii::ExcInternalError());
270 cfl_ = new_cfl;
271 }
272
278
279
285 void acceptable_tau_max_ratio(Number new_acceptable_tau_max_ratio) const
286 {
287 Assert(new_acceptable_tau_max_ratio >= Number(1.),
288 dealii::ExcInternalError());
289 acceptable_tau_max_ratio_ = new_acceptable_tau_max_ratio;
290 }
291
299
300
304
309
314
321
328
335
342
343 // FIXME: refactor to function
345
346 private:
348
352 typename Description::template Indicator<dim, Number>::Parameters
353 indicator_parameters_;
354
355 typename Description::template Limiter<dim, Number>::Parameters
356 limiter_parameters_;
357
358 typename Description::template RiemannSolver<dim, Number>::Parameters
359 riemann_solver_parameters_;
360
362
364
368
369 const MPIEnsemble &mpi_ensemble_;
370 std::map<std::string, dealii::Timer> &computing_timer_;
371
372 dealii::SmartPointer<const OfflineData<dim, Number>> offline_data_;
373 dealii::SmartPointer<const HyperbolicSystem> hyperbolic_system_;
374 dealii::SmartPointer<const InitialValues<Description, dim, Number>>
375 initial_values_;
376
377 mutable Number cfl_;
378 mutable Number acceptable_tau_max_ratio_;
379
380 mutable unsigned int n_restarts_;
381 mutable unsigned int n_corrections_;
382 mutable unsigned int n_warnings_;
383
384 InitialPrecomputedVector initial_precomputed_;
385
386 using ScalarVector = typename Vectors::ScalarVector<Number>;
387 mutable ScalarVector alpha_;
388
389 static constexpr auto n_bounds =
390 Description::template Limiter<dim, Number>::n_bounds;
391 mutable Vectors::MultiComponentVector<Number, n_bounds> bounds_;
392
393 using HyperbolicVector =
394 Vectors::MultiComponentVector<Number, problem_dimension>;
395 mutable HyperbolicVector r_;
396
397 mutable SparseMatrixSIMD<Number> dij_matrix_;
398 mutable SparseMatrixSIMD<Number> lij_matrix_;
399 mutable SparseMatrixSIMD<Number> lij_matrix_next_;
400 mutable SparseMatrixSIMD<Number, problem_dimension> pij_matrix_;
401
403 };
404
405} /* 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
const auto & acceptable_tau_max_ratio() const
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