ryujin 2.1.1 revision 350e54cc11f3d21282bddcf3e6517944dbc508bf
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 private:
96 const HyperbolicSystem &hyperbolic_system_;
97
98
99 public:
104
109
113 static constexpr unsigned int problem_dimension = 1;
114
118 using state_type = dealii::Tensor<1, problem_dimension, Number>;
119
123 using flux_type =
124 dealii::Tensor<1, problem_dimension, dealii::Tensor<1, dim, Number>>;
125
130
135 static inline const auto component_names =
136 []() -> std::array<std::string, problem_dimension> {
137 if constexpr (dim == 1)
138 return {"u"};
139 else if constexpr (dim == 2)
140 return {"u"};
141 else if constexpr (dim == 3)
142 return {"u"};
143 __builtin_trap();
144 }();
145
150 static inline const auto primitive_component_names =
151 []() -> std::array<std::string, problem_dimension> {
152 if constexpr (dim == 1)
153 return {"u"};
154 else if constexpr (dim == 2)
155 return {"u"};
156 else if constexpr (dim == 3)
157 return {"u"};
158 __builtin_trap();
159 }();
160
164 static constexpr unsigned int n_precomputed_values = 0;
165
169 using precomputed_type = std::array<Number, n_precomputed_values>;
170
174 static inline const auto precomputed_names =
175 std::array<std::string, n_precomputed_values>{};
176
180 static constexpr unsigned int n_initial_precomputed_values = 0;
181
186 std::array<Number, n_initial_precomputed_values>;
187
191 static inline const auto initial_precomputed_names =
192 std::array<std::string, n_initial_precomputed_values>{};
193
197 using StateVector = Vectors::
198 StateVector<ScalarNumber, problem_dimension, n_precomputed_values>;
199
205
211
219
221
225
229 static constexpr unsigned int n_precomputation_cycles = 0;
230
235 template <typename DISPATCH, typename SPARSITY>
236 void
237 precomputation_loop(unsigned int /*cycle*/,
238 const DISPATCH &dispatch_check,
239 const SPARSITY & /*sparsity_simd*/,
240 StateVector & /*state_vector*/,
241 unsigned int /*left*/,
242 unsigned int /*right*/,
243 const bool /*skip_constrained_dofs*/) const = delete;
244
246
250
256 bool is_admissible(const state_type & /*U*/) const
257 {
258 return true;
259 }
260
262
266
270 template <typename Lambda>
272 const dealii::types::boundary_id /*id*/,
273 const state_type &U,
274 const dealii::Tensor<1, dim, Number> & /*normal*/,
275 const Lambda & /*get_dirichlet_data*/) const
276 {
277 return U;
278 }
279
281
285
307 const InitialPrecomputedVector & /*piv*/,
308 const unsigned int /*i*/,
309 const state_type & /*U_i*/) const
310 {
311 return flux_contribution_type{};
312 }
313
316 const InitialPrecomputedVector & /*piv*/,
317 const unsigned int * /*js*/,
318 const state_type & /*U_j*/) const
319 {
320 return flux_contribution_type{};
321 }
322
329 const flux_contribution_type & /*flux_j*/,
330 const dealii::Tensor<1, dim, Number> & /*c_ij*/) const
331 {
332 return state_type{};
333 }
334
338 static constexpr bool have_high_order_flux = false;
339
343 const dealii::Tensor<1, dim, Number> &) const = delete;
344
346
350
352 static constexpr bool have_source_terms = false;
353
355 const unsigned int /*i*/,
356 const state_type & /*U_i*/,
357 const ScalarNumber /*tau*/) const = delete;
358
360 const unsigned int * /*js*/,
361 const state_type & /*U_j*/,
362 const ScalarNumber /*tau*/) const = delete;
363
365
369
380 template <typename ST>
381 state_type expand_state(const ST &state) const
382 {
383 return state;
384 }
385
390 state_type from_primitive_state(const state_type &primitive_state) const
391 {
392 return primitive_state;
393 }
394
400 {
401 return state;
402 }
403
409 template <typename Lambda>
411 const Lambda & /*lambda*/) const
412 {
413 return state;
414 }
415
416 }; /* HyperbolicSystemView */
417 } // namespace Skeleton
418} // 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")