Hllc.hpp
1 // Distributed under the MIT License.
2 // See LICENSE.txt for details.
3 
4 #pragma once
5 
6 #include <cstddef>
7 #include <string>
8 
9 #include "DataStructures/DataBox/Tag.hpp"
11 #include "Domain/FaceNormal.hpp"
12 #include "Evolution/Systems/NewtonianEuler/Characteristics.hpp"
13 #include "Evolution/Systems/NewtonianEuler/TagsDeclarations.hpp"
14 #include "NumericalAlgorithms/DiscontinuousGalerkin/Protocols.hpp"
15 #include "Options/Options.hpp"
16 #include "Utilities/Gsl.hpp"
17 #include "Utilities/ProtocolHelpers.hpp"
18 #include "Utilities/TMPL.hpp"
19 
20 /// \cond
21 namespace PUP {
22 class er;
23 } // namespace PUP
24 
25 namespace Tags {
26 template <typename>
27 struct NormalDotFlux;
28 template <typename>
29 struct Normalized;
30 } // namespace Tags
31 
32 class DataVector;
33 template <typename>
34 class Variables;
35 /// \endcond
36 
37 namespace NewtonianEuler {
38 namespace NumericalFluxes {
39 
40 /*!
41  * \brief Compute the HLLC numerical flux
42  *
43  * Class implementing the HLLC flux for the Newtonian Euler equations,
44  * originally introduced by E. F. Toro, M. Spruce and W. Speares \cite Toro1994.
45  * Let \f$F^k = [F^k(\rho), F^k(S^i), F^k(e)]^T\f$ be the fluxes of
46  * the conservative quantities \f$U = [\rho, S^i, e]^T\f$
47  * (see NewtonianEuler::ComputeFluxes). Denoting \f$v_n = n_k v^k\f$ and
48  * \f$F = n_kF^k\f$, where \f$v^k\f$ is the velocity and
49  * \f$n_k\f$ is the interface unit normal, the HLLC flux is
50  *
51  * \f{align*}
52  * G_\text{HLLC} =
53  * \begin{cases}
54  * F_\text{int}, & S_\text{min} > 0, \\
55  * G_{*\text{int}}, & S_\text{min} \leq 0\quad\text{and}\quad S_* > 0, \\
56  * G_{*\text{ext}}, & S_* \leq 0\quad\text{and}\quad S_\text{max} > 0, \\
57  * F_\text{ext}, & S_\text{max} \leq 0,
58  * \end{cases}
59  * \f}
60  *
61  * where
62  *
63  * \f{align*}
64  * G_{*\text{int}} &= \frac{S_*\left(F_\text{int}
65  * - S_\text{min}U_\text{int}\right)
66  * - S_\text{min} p_{*\text{int}}D_*}{S_* - S_\text{min}},\\
67  * G_{*\text{ext}} &= \frac{S_*\left(F_\text{ext}
68  * - S_\text{max}U_\text{ext}\right)
69  * - S_\text{max} p_{*\text{ext}}D_*}{S_* - S_\text{max}},
70  * \f}
71  *
72  * with
73  *
74  * \f{align*}
75  * p_{*\text{int}} &\equiv p_\text{int} + \rho_\text{int}
76  * \left[(v_n)_\text{int} - S_\text{min}\right]
77  * \left[(v_n)_\text{int} - S_*\right], \\
78  * p_{*\text{ext}} &\equiv p_\text{ext} + \rho_\text{ext}
79  * \left[(v_n)_\text{ext} - S_\text{max}\right]
80  * \left[(v_n)_\text{ext} - S_*\right], \\
81  * S_* &\equiv \frac{p_\text{int} + F(\rho)_\text{int}\left[(v_n)_\text{int}
82  * - S_\text{min}\right] - p_\text{ext}
83  * - F(\rho)_\text{ext}\left[(v_n)_\text{ext} -
84  * S_\text{max}\right]}{F(\rho)_\text{int} - \rho_\text{int}S_\text{min}
85  * - F(\rho)_\text{ext} + \rho_\text{ext}S_\text{max}},\\
86  * D_* &\equiv \left[\begin{array}{c}
87  * 0\\ n^i\\ S_*
88  * \end{array}\right],
89  * \f}
90  *
91  * and \f$S_\text{min}\f$ and \f$S_\text{max}\f$ are estimates of the minimum
92  * and maximum signal speeds bounding the ingoing and outgoing
93  * wavespeeds that arise when solving the Riemann problem. One requires
94  * \f$S_\text{min} \leq S_* \leq S_\text{max}\f$. As estimates, we use
95  *
96  * \f{align*}
97  * S_\text{min} &=
98  * \text{min}\left(\{\lambda_\text{int}\},\{\lambda_\text{ext}\}, 0\right)\\
99  * S_\text{max} &=
100  * \text{max}\left(\{\lambda_\text{int}\},\{\lambda_\text{ext}\}, 0\right),
101  * \f}
102  *
103  * where \f$\{\lambda\}\f$ is the set of all the characteristic speeds along a
104  * given normal. This way, the definition of \f$G_\text{HLLC}\f$ simplifies to
105  *
106  * \f{align*}
107  * G_\text{HLLC} =
108  * \begin{cases}
109  * \dfrac{S_*\left(F_\text{int} - S_\text{min}U_\text{int}\right)
110  * - S_\text{min} p_{*\text{int}}D_*}{S_* - S_\text{min}}, & S_* > 0, \\
111  * \dfrac{S_*\left(F_\text{ext} - S_\text{max}U_\text{ext}\right)
112  * - S_\text{max} p_{*\text{ext}}D_*}{S_* - S_\text{max}}, & S_* \leq 0. \\
113  * \end{cases}
114  * \f}
115  *
116  * \warning The HLLC flux implemented here does not incorporate any cure for
117  * the Carbuncle phenomenon or other shock instabilities reported in the
118  * literature. Prefer using another numerical flux in more than 1-d.
119  */
120 template <size_t Dim, typename Frame>
121 struct Hllc : tt::ConformsTo<dg::protocols::NumericalFlux> {
123 
124  /// Estimate for one of the signal speeds
126  using type = Scalar<DataVector>;
127  };
128 
129  /// Estimate for the other signal speed
131  using type = Scalar<DataVector>;
132  };
133 
134  /// The normal component of the velocity
136  using type = Scalar<DataVector>;
137  };
138 
139  using options = tmpl::list<>;
140  static constexpr Options::String help = {
141  "Compute the HLLC flux for the Newtonian Euler system."};
142 
143  // clang-tidy: non-const reference
144  void pup(PUP::er& /*p*/) noexcept {} // NOLINT
145 
146  using variables_tags =
147  tmpl::list<Tags::MassDensityCons, Tags::MomentumDensity<Dim>,
149 
150  using package_field_tags = tmpl::list<
157  NormalVelocity, LargestIngoingSpeed, LargestOutgoingSpeed>;
158  using package_extra_tags = tmpl::list<>;
159 
160  using argument_tags = tmpl::list<
166  char_speeds_tag,
168 
169  void package_data(
170  gsl::not_null<Scalar<DataVector>*> packaged_n_dot_f_mass_density,
171  gsl::not_null<tnsr::I<DataVector, Dim, Frame>*>
172  packaged_n_dot_f_momentum_density,
173  gsl::not_null<Scalar<DataVector>*> packaged_n_dot_f_energy_density,
174  gsl::not_null<Scalar<DataVector>*> packaged_mass_density,
175  gsl::not_null<tnsr::I<DataVector, Dim, Frame>*> packaged_momentum_density,
176  gsl::not_null<Scalar<DataVector>*> packaged_energy_density,
177  gsl::not_null<Scalar<DataVector>*> packaged_pressure,
178  gsl::not_null<tnsr::i<DataVector, Dim, Frame>*> packaged_face_normal,
179  gsl::not_null<Scalar<DataVector>*> packaged_normal_velocity,
180  gsl::not_null<Scalar<DataVector>*> packaged_largest_ingoing_speed,
181  gsl::not_null<Scalar<DataVector>*> packaged_largest_outgoing_speed,
182  const Scalar<DataVector>& normal_dot_flux_mass_density,
183  const tnsr::I<DataVector, Dim, Frame>& normal_dot_flux_momentum_density,
184  const Scalar<DataVector>& normal_dot_flux_energy_density,
185  const Scalar<DataVector>& mass_density,
186  const tnsr::I<DataVector, Dim, Frame>& momentum_density,
187  const Scalar<DataVector>& energy_density,
188  const tnsr::I<DataVector, Dim, Frame>& velocity,
189  const Scalar<DataVector>& pressure,
190  const typename char_speeds_tag::type& characteristic_speeds,
191  const tnsr::i<DataVector, Dim, Frame>& interface_unit_normal) const
192  noexcept;
193 
194  void operator()(
195  gsl::not_null<Scalar<DataVector>*> normal_dot_numerical_flux_mass_density,
196  gsl::not_null<tnsr::I<DataVector, Dim, Frame>*>
197  normal_dot_numerical_flux_momentum_density,
199  normal_dot_numerical_flux_energy_density,
200  const Scalar<DataVector>& normal_dot_flux_mass_density_int,
201  const tnsr::I<DataVector, Dim, Frame>&
202  normal_dot_flux_momentum_density_int,
203  const Scalar<DataVector>& normal_dot_flux_energy_density_int,
204  const Scalar<DataVector>& mass_density_int,
205  const tnsr::I<DataVector, Dim, Frame>& momentum_density_int,
206  const Scalar<DataVector>& energy_density_int,
207  const Scalar<DataVector>& pressure_int,
208  const tnsr::i<DataVector, Dim, Frame>& interface_unit_normal,
209  const Scalar<DataVector>& normal_velocity_int,
210  const Scalar<DataVector>& largest_ingoing_speed_int,
211  const Scalar<DataVector>& largest_outgoing_speed_int,
212  const Scalar<DataVector>& minus_normal_dot_flux_mass_density_ext,
213  const tnsr::I<DataVector, Dim, Frame>&
214  minus_normal_dot_flux_momentum_density_ext,
215  const Scalar<DataVector>& minus_normal_dot_flux_energy_density_ext,
216  const Scalar<DataVector>& mass_density_ext,
217  const tnsr::I<DataVector, Dim, Frame>& momentum_density_ext,
218  const Scalar<DataVector>& energy_density_ext,
219  const Scalar<DataVector>& pressure_ext,
220  const tnsr::i<DataVector, Dim, Frame>& minus_interface_unit_normal,
221  const Scalar<DataVector>& minus_normal_velocity_ext,
222  // names are inverted w.r.t interior data. See package_data()
223  const Scalar<DataVector>& minus_largest_outgoing_speed_ext,
224  const Scalar<DataVector>& minus_largest_ingoing_speed_ext) const noexcept;
225 };
226 
227 } // namespace NumericalFluxes
228 } // namespace NewtonianEuler
NewtonianEuler::Tags::CharacteristicSpeeds
The characteristic speeds.
Definition: Tags.hpp:84
NewtonianEuler::NumericalFluxes::Hllc::NormalVelocity
The normal component of the velocity.
Definition: Hllc.hpp:135
Options.hpp
FaceNormal.hpp
db::SimpleTag
Tags for the DataBox inherit from this type.
Definition: Tag.hpp:23
Tags::NormalDotFlux
Prefix indicating a boundary unit normal vector dotted into the flux.
Definition: Prefixes.hpp:96
NewtonianEuler
Items related to evolving the Newtonian Euler system.
Definition: EvolveNewtonianEulerFwd.hpp:8
NewtonianEuler::Tags::MassDensityCons
The mass density of the fluid (as a conservative variable).
Definition: Tags.hpp:31
NewtonianEuler::NumericalFluxes::Hllc::LargestIngoingSpeed
Estimate for one of the signal speeds.
Definition: Hllc.hpp:125
Tags::Normalized
Definition: Magnitude.hpp:138
cstddef
NewtonianEuler::NumericalFluxes::Hllc::LargestOutgoingSpeed
Estimate for the other signal speed.
Definition: Hllc.hpp:130
NewtonianEuler::characteristic_speeds
void characteristic_speeds(gsl::not_null< std::array< DataVector, Dim+2 > * > char_speeds, const tnsr::I< DataVector, Dim > &velocity, const Scalar< DataVector > &sound_speed, const tnsr::i< DataVector, Dim > &normal) noexcept
Compute the characteristic speeds of NewtonianEuler system.
DataVector
Stores a collection of function values.
Definition: DataVector.hpp:42
NewtonianEuler::Tags::Pressure
The fluid pressure.
Definition: Tags.hpp:66
NewtonianEuler::Tags::EnergyDensity
The energy density of the fluid.
Definition: Tags.hpp:45
NewtonianEuler::Tags::MomentumDensity
The momentum density of the fluid.
Definition: Tags.hpp:37
Scalar
Tensor< T, Symmetry<>, index_list<> > Scalar
Definition: TypeAliases.hpp:21
Gsl.hpp
TypeAliases.hpp
Options::String
const char *const String
The string used in option structs.
Definition: Options.hpp:32
NewtonianEuler::Tags::Velocity
The macroscopic or flow velocity of the fluid.
Definition: Tags.hpp:51
NewtonianEuler::NumericalFluxes::Hllc
Compute the HLLC numerical flux.
Definition: Hllc.hpp:121
tt::ConformsTo
Indicate a class conforms to the Protocol.
Definition: ProtocolHelpers.hpp:22
TMPL.hpp
gsl::not_null
Require a pointer to not be a nullptr
Definition: Gsl.hpp:183
string