8#include <compile_time_options.h>
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>
78 template <
typename Description,
int dim,
typename Number =
double>
105 View::n_precomputation_cycles;
118 std::map<std::string, dealii::Timer> &computing_timer,
122 const std::string &subsection =
"/HyperbolicModule");
221 template <
int stages>
224 std::array<std::reference_wrapper<const StateVector>, stages>
226 const std::array<Number, stages> stage_weights,
228 Number tau = Number(0.),
229 std::atomic<Number> tau_max = std::numeric_limits<Number>::max())
const;
238 void cfl(Number new_cfl)
const
240 Assert(cfl_ > Number(0.), dealii::ExcInternalError());
293 indicator_parameters_;
299 riemann_solver_parameters_;
310 std::map<std::
string, dealii::Timer> &computing_timer_;
312 dealii::SmartPointer<const
OfflineData<dim, Number>> offline_data_;
319 mutable
unsigned int n_restarts_;
321 mutable
unsigned int n_warnings_;
325 using ScalarVector = typename Vectors::ScalarVector<Number>;
326 mutable ScalarVector alpha_;
328 static constexpr auto n_bounds =
330 mutable Vectors::MultiComponentVector<Number, n_bounds> bounds_;
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 & alpha() const
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
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
Euler::HyperbolicSystem HyperbolicSystem