ryujin 2.1.1 revision 336b16a72e829721302c626ec7071b92032b8248
hyperbolic_system.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
9#include <discretization.h>
11#include <patterns_conversion.h>
12#include <simd.h>
13#include <state_vector.h>
14
15#include <deal.II/base/parameter_acceptor.h>
16#include <deal.II/base/tensor.h>
17
18#include <array>
19#include <functional>
20
21namespace ryujin
22{
23 namespace Skeleton
24 {
25 template <int dim, typename Number>
26 class HyperbolicSystemView;
27
33 class HyperbolicSystem final : public dealii::ParameterAcceptor
34 {
35 public:
39 static inline const std::string problem_name =
40 "Skeleton Hyperbolic System";
41
45 HyperbolicSystem(const std::string &subsection = "/HyperbolicSystem")
46 : ParameterAcceptor(subsection)
47 {
48 }
49
56 template <int dim, typename Number>
57 auto view() const
58 {
60 }
61
62 private:
63 template <int dim, typename Number>
65 }; /* HyperbolicSystem */
66
67
73 template <int dim, typename Number>
75 {
76 public:
81 HyperbolicSystemView(const HyperbolicSystem &hyperbolic_system)
82 : hyperbolic_system_(hyperbolic_system)
83 {
84 }
85
89 template <int dim2, typename Number2>
90 auto view() const
91 {
92 return HyperbolicSystemView<dim2, Number2>{hyperbolic_system_};
93 }
94
95
96 private:
97 const HyperbolicSystem &hyperbolic_system_;
98
99
100 public:
105
110
114 static constexpr unsigned int problem_dimension = 1;
115
119 using state_type = dealii::Tensor<1, problem_dimension, Number>;
120
124 using flux_type =
125 dealii::Tensor<1, problem_dimension, dealii::Tensor<1, dim, Number>>;
126
131
136 static inline const auto component_names =
137 []() -> std::array<std::string, problem_dimension> {
138 if constexpr (dim == 1)
139 return {"u"};
140 else if constexpr (dim == 2)
141 return {"u"};
142 else if constexpr (dim == 3)
143 return {"u"};
144 __builtin_trap();
145 }();
146
151 static inline const auto primitive_component_names =
152 []() -> std::array<std::string, problem_dimension> {
153 if constexpr (dim == 1)
154 return {"u"};
155 else if constexpr (dim == 2)
156 return {"u"};
157 else if constexpr (dim == 3)
158 return {"u"};
159 __builtin_trap();
160 }();
161
165 static constexpr unsigned int n_precomputed_values = 0;
166
170 using precomputed_type = std::array<Number, n_precomputed_values>;
171
175 static inline const auto precomputed_names =
176 std::array<std::string, n_precomputed_values>{};
177
181 static constexpr unsigned int n_initial_precomputed_values = 0;
182
187 std::array<Number, n_initial_precomputed_values>;
188
192 static inline const auto initial_precomputed_names =
193 std::array<std::string, n_initial_precomputed_values>{};
194
198 using StateVector = Vectors::
199 StateVector<ScalarNumber, problem_dimension, n_precomputed_values>;
200
206
214
216
220
224 static constexpr unsigned int n_precomputation_cycles = 0;
225
230 template <typename DISPATCH, typename SPARSITY>
231 void precomputation_loop(unsigned int /*cycle*/,
232 const DISPATCH &dispatch_check,
233 const SPARSITY & /*sparsity_simd*/,
234 StateVector & /*state_vector*/,
235 unsigned int /*left*/,
236 unsigned int /*right*/) const = delete;
237
239
243
249 bool is_admissible(const state_type & /*U*/) const
250 {
251 return true;
252 }
253
255
259
263 template <typename Lambda>
265 const dealii::types::boundary_id /*id*/,
266 const state_type &U,
267 const dealii::Tensor<1, dim, Number> & /*normal*/,
268 const Lambda & /*get_dirichlet_data*/) const
269 {
270 return U;
271 }
272
274
278
300 const InitialPrecomputedVector & /*piv*/,
301 const unsigned int /*i*/,
302 const state_type & /*U_i*/) const
303 {
304 return flux_contribution_type{};
305 }
306
309 const InitialPrecomputedVector & /*piv*/,
310 const unsigned int * /*js*/,
311 const state_type & /*U_j*/) const
312 {
313 return flux_contribution_type{};
314 }
315
322 const flux_contribution_type & /*flux_j*/,
323 const dealii::Tensor<1, dim, Number> & /*c_ij*/) const
324 {
325 return state_type{};
326 }
327
331 static constexpr bool have_high_order_flux = false;
332
336 const dealii::Tensor<1, dim, Number> &) const = delete;
337
339
343
345 static constexpr bool have_source_terms = false;
346
348 const unsigned int /*i*/,
349 const state_type & /*U_i*/,
350 const ScalarNumber /*tau*/) const = delete;
351
353 const unsigned int * /*js*/,
354 const state_type & /*U_j*/,
355 const ScalarNumber /*tau*/) const = delete;
356
358
362
373 template <typename ST>
374 state_type expand_state(const ST &state) const
375 {
376 return state;
377 }
378
383 state_type from_primitive_state(const state_type &primitive_state) const
384 {
385 return primitive_state;
386 }
387
393 {
394 return state;
395 }
396
402 template <typename Lambda>
404 const Lambda & /*lambda*/) const
405 {
406 return state;
407 }
408
409 }; /* HyperbolicSystemView */
410 } // namespace Skeleton
411} // namespace ryujin
std::array< Number, n_precomputed_values > precomputed_type
state_type apply_galilei_transform(const state_type &state, const Lambda &) const
bool is_admissible(const state_type &) const
static constexpr unsigned int n_precomputed_values
state_type nodal_source(const PrecomputedVector &, const unsigned int *, const state_type &, const ScalarNumber) const =delete
state_type nodal_source(const PrecomputedVector &, const unsigned int, const state_type &, const ScalarNumber) const =delete
dealii::Tensor< 1, problem_dimension, Number > state_type
HyperbolicSystemView(const HyperbolicSystem &hyperbolic_system)
state_type expand_state(const ST &state) const
typename get_value_type< Number >::type ScalarNumber
dealii::Tensor< 1, problem_dimension, dealii::Tensor< 1, dim, Number > > flux_type
state_type apply_boundary_conditions(const dealii::types::boundary_id, const state_type &U, const dealii::Tensor< 1, dim, Number > &, const Lambda &) const
Vectors::StateVector< ScalarNumber, problem_dimension, n_precomputed_values > StateVector
void precomputation_loop(unsigned int, const DISPATCH &dispatch_check, const SPARSITY &, StateVector &, unsigned int, unsigned int) const =delete
state_type high_order_flux_divergence(const flux_contribution_type &, const flux_contribution_type &, const dealii::Tensor< 1, dim, Number > &) const =delete
flux_contribution_type flux_contribution(const PrecomputedVector &, const InitialPrecomputedVector &, const unsigned int *, const state_type &) const
state_type flux_divergence(const flux_contribution_type &, const flux_contribution_type &, const dealii::Tensor< 1, dim, Number > &) const
static constexpr unsigned int n_precomputation_cycles
flux_contribution_type flux_contribution(const PrecomputedVector &, const InitialPrecomputedVector &, const unsigned int, const state_type &) const
static constexpr unsigned int problem_dimension
static constexpr unsigned int n_initial_precomputed_values
state_type to_primitive_state(const state_type &state) const
state_type from_primitive_state(const state_type &primitive_state) const
std::array< Number, n_initial_precomputed_values > initial_precomputed_type
static const std::string problem_name
HyperbolicSystem(const std::string &subsection="/HyperbolicSystem")