10#include <compile_time_options.h>
14#include <deal.II/base/vectorization.h>
19 namespace ShallowWater
26 template <
int dim,
typename Number =
double>
34 HyperbolicSystem::problem_dimension<dim>;
45 HyperbolicSystem::n_precomputed_values<dim>;
88 : hyperbolic_system(hyperbolic_system)
89 , precomputed_values(precomputed_values)
103 void add(
const unsigned int *js,
105 const dealii::Tensor<1, dim, Number> &c_ij);
110 Number
alpha(
const Number h_i)
const;
129 Number pressure_i = 0.;
141 template <
int dim,
typename Number>
142 DEAL_II_ALWAYS_INLINE
inline void
145 h_i = hyperbolic_system.water_depth(U_i);
148 const auto &[mathematical_entropy] =
149 precomputed_values.template get_tensor<Number, precomputed_type>(i);
151 eta_i = mathematical_entropy;
153 d_eta_i = hyperbolic_system.mathematical_entropy_derivative(U_i);
154 f_i = hyperbolic_system.f(U_i);
155 pressure_i = hyperbolic_system.pressure(U_i);
162 template <
int dim,
typename Number>
163 DEAL_II_ALWAYS_INLINE
inline void
166 const dealii::Tensor<1, dim, Number> &c_ij)
170 const auto &[eta_j] =
171 precomputed_values.template get_tensor<Number, precomputed_type>(js);
173 const auto velocity_j = hyperbolic_system.momentum(U_j) *
174 hyperbolic_system.inverse_water_depth_sharp(U_j);
175 const auto f_j = hyperbolic_system.f(U_j);
176 const auto pressure_j = hyperbolic_system.pressure(U_j);
178 left += (eta_j + pressure_j) * (velocity_j * c_ij);
180 for (
unsigned int k = 0; k < problem_dimension; ++k)
181 right[k] += (f_j[k] - f_i[k]) * c_ij;
185 template <
int dim,
typename Number>
186 DEAL_II_ALWAYS_INLINE
inline Number
192 for (
unsigned int k = 0; k < problem_dimension; ++k) {
193 my_sum += d_eta_i[k] * right[k];
196 Number numerator = std::abs(left - my_sum);
197 Number denominator = std::abs(left) + std::abs(my_sum);
199 constexpr ScalarNumber eps = std::numeric_limits<ScalarNumber>::epsilon();
200 const auto quotient = std::abs(numerator + eps) /
201 (denominator + hd_i * std::abs(eta_i) + eps);
206 return std::min(Number(1.), evc_alpha_0_ * quotient);
dealii::Tensor< 1, problem_dimension< dim >, dealii::Tensor< 1, dim, Number > > flux_type
std::array< Number, n_precomputed_values< dim > > precomputed_type
dealii::Tensor< 1, problem_dimension< dim >, Number > state_type
static constexpr unsigned int problem_dimension
static constexpr unsigned int n_precomputed_values
Number alpha(const Number h_i) const
HyperbolicSystem::flux_type< dim, Number > flux_type
void reset(const unsigned int i, const state_type &U_i)
HyperbolicSystem::precomputed_type< dim, Number > precomputed_type
Indicator(const HyperbolicSystem &hyperbolic_system, const MultiComponentVector< ScalarNumber, n_precomputed_values > &precomputed_values)
HyperbolicSystem::state_type< dim, Number > state_type
typename get_value_type< Number >::type ScalarNumber
void add(const unsigned int *js, const state_type &U_j, const dealii::Tensor< 1, dim, Number > &c_ij)