11#include <deal.II/base/config.h>
19#include <deal.II/base/parameter_acceptor.h>
20#include <deal.II/base/tensor.h>
26 namespace ScalarConservation
28 template <
int dim,
typename Number>
29 class HyperbolicSystemView;
43 static inline std::string
problem_name =
"Scalar conservation equation";
56 template <
int dim,
typename Number>
71 std::shared_ptr<Flux> selected_flux_;
73 template <
int dim,
typename Number>
84 template <
int dim,
typename Number>
93 : hyperbolic_system_(hyperbolic_system)
100 template <
int dim2,
typename Number2>
116 DEAL_II_ALWAYS_INLINE
inline const std::string &
flux()
const
118 return hyperbolic_system_.flux_;
124 const auto &
flux = hyperbolic_system_.selected_flux_;
137 DEAL_II_ALWAYS_INLINE
inline dealii::Tensor<1, dim, Number>
143 DEAL_II_ALWAYS_INLINE
inline dealii::Tensor<1, dim, Number>
170 using state_type = dealii::Tensor<1, problem_dimension, Number>;
176 dealii::Tensor<1, problem_dimension, dealii::Tensor<1, dim, Number>>;
188 std::array<std::string, problem_dimension>{
"u"};
195 std::array<std::string, problem_dimension>{
"u"};
211 []() -> std::array<std::string, n_precomputed_values> {
212 if constexpr (dim == 1)
214 else if constexpr (dim == 2)
215 return {
"f_1",
"f_2",
"df_1",
"df_2"};
216 else if constexpr (dim == 3)
217 return {
"f_1",
"f_2",
"f_3",
"df_1",
"df_2",
"df_3"};
230 std::array<Number, n_initial_precomputed_values>;
236 std::array<std::string, n_initial_precomputed_values>{};
242 StateVector<ScalarNumber, problem_dimension, n_precomputed_values>;
280 template <
typename DISPATCH,
typename SPARSITY>
282 const DISPATCH &dispatch_check,
283 const SPARSITY &sparsity_simd,
286 unsigned int right)
const;
355 template <
typename Lambda>
359 const dealii::Tensor<1, dim, Number> &normal,
360 const Lambda &get_dirichlet_data)
const;
372 dealii::Tensor<1, dim, Number>
405 const unsigned int i,
411 const unsigned int *js,
421 const dealii::Tensor<1, dim, Number> &c_ij)
const;
429 const dealii::Tensor<1, dim, Number> &c_ij)
const =
delete;
441 const unsigned int i,
446 const unsigned int *js,
466 template <
typename ST>
478 return primitive_state;
495 template <
typename Lambda>
497 const Lambda & )
const
513 : ParameterAcceptor(subsection)
516 add_parameter(
"flux",
518 "The scalar flux. Valid names are given by any of the "
519 "subsections defined below");
527 const auto populate_functions = [
this]() {
528 bool initialized =
false;
529 for (
auto &it : flux_list_)
532 if (it->name() == flux_) {
534 it->parse_parameters_call_back();
535 problem_name =
"Scalar conservation equation (" + it->name() +
536 ": " + it->flux_formula() +
")";
541 AssertThrow(initialized,
543 "Could not find a flux description with name \"" +
547 ParameterAcceptor::parse_parameters_call_back.connect(populate_functions);
548 populate_functions();
552 template <
int dim,
typename Number>
553 DEAL_II_ALWAYS_INLINE
inline dealii::Tensor<1, dim, Number>
556 const auto &flux = hyperbolic_system_.selected_flux_;
557 dealii::Tensor<1, dim, Number> result;
560 for (
unsigned int k = 0; k < dim; ++k) {
561 if constexpr (std::is_same_v<ScalarNumber, Number>) {
562 result[k] = flux->value(u, k);
564 for (
unsigned int s = 0; s < Number::size(); ++s) {
565 result[k][s] = flux->value(u[s], k);
574 template <
int dim,
typename Number>
575 DEAL_II_ALWAYS_INLINE
inline dealii::Tensor<1, dim, Number>
577 const Number &u)
const
579 const auto &flux = hyperbolic_system_.selected_flux_;
580 dealii::Tensor<1, dim, Number> result;
583 for (
unsigned int k = 0; k < dim; ++k) {
584 if constexpr (std::is_same_v<ScalarNumber, Number>) {
585 result[k] = flux->gradient(u, k);
587 for (
unsigned int s = 0; s < Number::size(); ++s) {
588 result[k][s] = flux->gradient(u[s], k);
597 template <
int dim,
typename Number>
598 template <
typename DISPATCH,
typename SPARSITY>
599 DEAL_II_ALWAYS_INLINE
inline void
601 unsigned int cycle [[maybe_unused]],
602 const DISPATCH &dispatch_check,
603 const SPARSITY &sparsity_simd,
606 unsigned int right)
const
608 Assert(cycle == 0, dealii::ExcInternalError());
610 const auto &U = std::get<0>(state_vector);
611 auto &precomputed = std::get<1>(state_vector);
615 unsigned int stride_size = get_stride_size<Number>;
618 for (
unsigned int i = left; i < right; i += stride_size) {
621 const unsigned int row_length = sparsity_simd.row_length(i);
627 const auto U_i = U.template get_tensor<Number>(i);
628 const auto u_i = state(U_i);
630 const auto f_i = flux_function(u_i);
631 const auto df_i = flux_gradient_function(u_i);
635 for (
unsigned int k = 0; k < n_precomputed_values / 2; ++k) {
637 prec_i[dim + k] = df_i[k];
640 precomputed.template write_tensor<Number>(prec_i, i);
645 template <
int dim,
typename Number>
646 DEAL_II_ALWAYS_INLINE
inline Number
653 template <
int dim,
typename Number>
654 DEAL_II_ALWAYS_INLINE
inline Number
661 template <
int dim,
typename Number>
662 DEAL_II_ALWAYS_INLINE
inline Number
664 const Number &u)
const
670 template <
int dim,
typename Number>
671 DEAL_II_ALWAYS_INLINE
inline Number
673 const Number &u)
const
675 return std::abs(k - u);
679 template <
int dim,
typename Number>
680 DEAL_II_ALWAYS_INLINE
inline Number
682 const Number &k,
const Number &u)
const
684 constexpr auto gte = dealii::SIMDComparison::greater_than_or_equal;
686 return dealii::compare_and_apply_mask<gte>(u, k, Number(1.), Number(-1.));
690 template <
int dim,
typename Number>
691 template <
typename Lambda>
692 DEAL_II_ALWAYS_INLINE
inline auto
694 dealii::types::boundary_id
id,
696 const dealii::Tensor<1, dim, Number> & ,
697 const Lambda &get_dirichlet_data)
const ->
state_type
702 result = get_dirichlet_data();
707 "Invalid boundary ID »Boundary::dirichlet_momentum«, "
708 "enforcing Dirichlet boundary conditions on a momentum "
709 "is not possible for scalar conservation equations."));
714 dealii::ExcMessage(
"Invalid boundary ID »Boundary::slip«, slip "
715 "boundary conditions are unavailable for scalar "
716 "conservation equations."));
722 dealii::ExcMessage(
"Invalid boundary ID »Boundary::no_slip«, "
723 "no-slip boundary conditions are unavailable "
724 "for scalar conservation equations."));
730 dealii::ExcMessage(
"Invalid boundary ID »Boundary::dynamic«, "
731 "dynamic boundary conditions are unavailable "
732 "for scalar conservation equations."));
736 AssertThrow(
false, dealii::ExcNotImplemented());
743 template <
int dim,
typename Number>
744 DEAL_II_ALWAYS_INLINE
inline dealii::Tensor<1, dim, Number>
748 dealii::Tensor<1, dim, Number> result;
750 if constexpr (dim == 1) {
751 const auto &[f, df] = precomputed;
754 }
else if constexpr (dim == 2) {
755 const auto &[f_1, f_2, df_1, df_2] = precomputed;
759 }
else if constexpr (dim == 3) {
760 const auto &[f_1, f_2, f_3, df_1, df_2, df_3] = precomputed;
770 template <
int dim,
typename Number>
771 DEAL_II_ALWAYS_INLINE
inline dealii::Tensor<1, dim, Number>
775 dealii::Tensor<1, dim, Number> result;
777 if constexpr (dim == 1) {
778 const auto &[f, df] = precomputed;
781 }
else if constexpr (dim == 2) {
782 const auto &[f_1, f_2, df_1, df_2] = precomputed;
786 }
else if constexpr (dim == 3) {
787 const auto &[f_1, f_2, f_3, df_1, df_2, df_3] = precomputed;
797 template <
int dim,
typename Number>
798 DEAL_II_ALWAYS_INLINE
inline auto
802 const unsigned int i,
807 result[0] = construct_flux_tensor(
808 pv.template get_tensor<Number, precomputed_type>(i));
813 template <
int dim,
typename Number>
814 DEAL_II_ALWAYS_INLINE
inline auto
818 const unsigned int *js,
823 result[0] = construct_flux_tensor(
824 pv.template get_tensor<Number, precomputed_type>(js));
829 template <
int dim,
typename Number>
830 DEAL_II_ALWAYS_INLINE
inline auto
834 const dealii::Tensor<1, dim, Number> &c_ij)
const ->
state_type
dealii::Tensor< 1, problem_dimension, Number > state_type
std::array< Number, n_precomputed_values > precomputed_type
Vectors::StateVector< ScalarNumber, problem_dimension, n_precomputed_values > StateVector
static const auto initial_precomputed_names
flux_contribution_type flux_contribution(const PrecomputedVector &pv, const InitialPrecomputedVector &, const unsigned int i, const state_type &) const
state_type apply_galilei_transform(const state_type &state, const Lambda &) const
static constexpr unsigned int n_precomputation_cycles
Number kruzkov_entropy_derivative(const Number &k, const Number &u) const
DEAL_II_ALWAYS_INLINE const std::string & flux() const
static constexpr unsigned int problem_dimension
HyperbolicSystemView(const HyperbolicSystem &hyperbolic_system)
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > flux_function(const Number &u) const
dealii::Tensor< 1, problem_dimension, dealii::Tensor< 1, dim, Number > > flux_type
std::array< Number, n_initial_precomputed_values > initial_precomputed_type
state_type from_primitive_state(const state_type &primitive_state) const
std::array< Number, n_precomputed_values > precomputed_type
static const auto component_names
Number square_entropy_derivative(const Number &u) const
dealii::Tensor< 1, dim, Number > construct_flux_gradient_tensor(const precomputed_type &precomputed_state) const
state_type to_primitive_state(const state_type &state) const
DEAL_II_ALWAYS_INLINE ScalarNumber derivative_approximation_delta() const
Number square_entropy(const Number &u) const
static const auto primitive_component_names
void precomputation_loop(unsigned int cycle, const DISPATCH &dispatch_check, const SPARSITY &sparsity_simd, StateVector &state_vector, unsigned int left, unsigned int right) const
static constexpr bool have_high_order_flux
state_type nodal_source(const PrecomputedVector &pv, const unsigned int *js, const state_type &U_j, const ScalarNumber tau) const =delete
typename get_value_type< Number >::type ScalarNumber
static Number state(const state_type &U)
Number kruzkov_entropy(const Number &k, const Number &u) const
state_type nodal_source(const PrecomputedVector &pv, const unsigned int i, const state_type &U_i, const ScalarNumber tau) const =delete
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > flux_gradient_function(const Number &u) const
static constexpr unsigned int n_initial_precomputed_values
static constexpr unsigned int n_precomputed_values
dealii::Tensor< 1, problem_dimension, Number > state_type
bool is_admissible(const state_type &) const
state_type expand_state(const ST &state) const
static const auto precomputed_names
state_type high_order_flux_divergence(const flux_contribution_type &, const flux_contribution_type &, const dealii::Tensor< 1, dim, Number > &c_ij) const =delete
flux_type flux_contribution_type
Vectors::StateVector< ScalarNumber, problem_dimension, n_precomputed_values > StateVector
dealii::Tensor< 1, dim, Number > construct_flux_tensor(const precomputed_type &precomputed_state) const
state_type apply_boundary_conditions(const dealii::types::boundary_id id, const state_type &U, const dealii::Tensor< 1, dim, Number > &normal, const Lambda &get_dirichlet_data) const
state_type flux_divergence(const flux_contribution_type &flux_i, const flux_contribution_type &flux_j, const dealii::Tensor< 1, dim, Number > &c_ij) const
static constexpr bool have_source_terms
static std::string problem_name
HyperbolicSystem(const std::string &subsection="/HyperbolicSystem")
std::set< std::shared_ptr< Flux > > flux_list_type
void populate_flux_list(flux_list_type &flux_list, const std::string &subsection)
DEAL_II_ALWAYS_INLINE FT add(const FT &flux_left_ij, const FT &flux_right_ij)
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, problem_dim, T > contract(const FT &flux_ij, const TT &c_ij)