ryujin 2.1.1 revision ef0fcd4010d109b860652ece4a7b8963fb7d46b1
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
8#include <compile_time_options.h>
9
10#include <convenience_macros.h>
11#include <discretization.h>
13#include <patterns_conversion.h>
14#include <simd.h>
15#include <state_vector.h>
16
17#include <deal.II/base/parameter_acceptor.h>
18#include <deal.II/base/tensor.h>
19
20#include <array>
21
22namespace ryujin
23{
24 namespace Skeleton
25 {
26 template <int dim, typename Number>
27 class HyperbolicSystemView;
28
34 class HyperbolicSystem final : public dealii::ParameterAcceptor
35 {
36 public:
40 static inline const std::string problem_name =
41 "Skeleton Hyperbolic System";
42
46 HyperbolicSystem(const std::string &subsection = "/HyperbolicSystem")
47 : ParameterAcceptor(subsection)
48 {
49 }
50
57 template <int dim, typename Number>
58 auto view() const
59 {
61 }
62
63 private:
64 template <int dim, typename Number>
66 }; /* HyperbolicSystem */
67
68
74 template <int dim, typename Number>
76 {
77 public:
82 HyperbolicSystemView(const HyperbolicSystem &hyperbolic_system)
83 : hyperbolic_system_(hyperbolic_system)
84 {
85 }
86
90 template <int dim2, typename Number2>
91 auto view() const
92 {
93 return HyperbolicSystemView<dim2, Number2>{hyperbolic_system_};
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
212
220
222
226
230 static constexpr unsigned int n_precomputation_cycles = 0;
231
236 template <typename DISPATCH, typename SPARSITY>
237 void
238 precomputation_loop(unsigned int /*cycle*/,
239 const DISPATCH &dispatch_check,
240 const SPARSITY & /*sparsity_simd*/,
241 StateVector & /*state_vector*/,
242 unsigned int /*left*/,
243 unsigned int /*right*/,
244 const bool /*skip_constrained_dofs*/) const = delete;
245
247
251
257 bool is_admissible(const state_type & /*U*/) const
258 {
259 return true;
260 }
261
263
267
271 template <typename Lambda>
273 const dealii::types::boundary_id /*id*/,
274 const state_type &U,
275 const dealii::Tensor<1, dim, Number> & /*normal*/,
276 const Lambda & /*get_dirichlet_data*/) const
277 {
278 return U;
279 }
280
282
286
308 const InitialPrecomputedVector & /*piv*/,
309 const unsigned int /*i*/,
310 const state_type & /*U_i*/) const
311 {
312 return flux_contribution_type{};
313 }
314
317 const InitialPrecomputedVector & /*piv*/,
318 const unsigned int * /*js*/,
319 const state_type & /*U_j*/) const
320 {
321 return flux_contribution_type{};
322 }
323
330 const flux_contribution_type & /*flux_j*/,
331 const dealii::Tensor<1, dim, Number> & /*c_ij*/) const
332 {
333 return state_type{};
334 }
335
339 static constexpr bool have_high_order_flux = false;
340
344 const dealii::Tensor<1, dim, Number> &) const = delete;
345
347
351
353 static constexpr bool have_source_terms = false;
354
356 const unsigned int /*i*/,
357 const state_type & /*U_i*/,
358 const ScalarNumber /*tau*/) const = delete;
359
361 const unsigned int * /*js*/,
362 const state_type & /*U_j*/,
363 const ScalarNumber /*tau*/) const = delete;
364
366
370
381 template <typename ST>
382 state_type expand_state(const ST &state) const
383 {
384 return state;
385 }
386
391 state_type from_primitive_state(const state_type &primitive_state) const
392 {
393 return primitive_state;
394 }
395
401 {
402 return state;
403 }
404
410 template <typename Lambda>
412 const Lambda & /*lambda*/) const
413 {
414 return state;
415 }
416
417 }; /* HyperbolicSystemView */
418 } // namespace Skeleton
419} // 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 bool) 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")