Line data Source code
1 0 : // Distributed under the MIT License.
2 : // See LICENSE.txt for details.
3 :
4 : #pragma once
5 :
6 : #include <memory>
7 : #include <optional>
8 :
9 : #include "DataStructures/DataBox/Prefixes.hpp"
10 : #include "DataStructures/Tensor/TypeAliases.hpp"
11 : #include "Evolution/BoundaryCorrection.hpp"
12 : #include "Evolution/Systems/NewtonianEuler/Tags.hpp"
13 : #include "NumericalAlgorithms/DiscontinuousGalerkin/Formulation.hpp"
14 : #include "Options/String.hpp"
15 : #include "PointwiseFunctions/Hydro/EquationsOfState/EquationOfState.hpp"
16 : #include "PointwiseFunctions/Hydro/Tags.hpp"
17 : #include "Utilities/Gsl.hpp"
18 : #include "Utilities/Serialization/CharmPupable.hpp"
19 : #include "Utilities/TMPL.hpp"
20 :
21 : /// \cond
22 : class DataVector;
23 : namespace gsl {
24 : template <typename T>
25 : class not_null;
26 : } // namespace gsl
27 : namespace PUP {
28 : class er;
29 : } // namespace PUP
30 : /// \endcond
31 :
32 : namespace NewtonianEuler::BoundaryCorrections {
33 : /*!
34 : * \brief A Rusanov/local Lax-Friedrichs Riemann solver
35 : *
36 : * Let \f$U\f$ be the state vector of evolved variables, \f$F^i\f$ the
37 : * corresponding fluxes, and \f$n_i\f$ be the outward directed unit normal to
38 : * the interface. Denoting \f$F := n_i F^i\f$, the %Rusanov boundary correction
39 : * is
40 : *
41 : * \f{align*}
42 : * G_\text{Rusanov} = \frac{F_\text{int} - F_\text{ext}}{2} -
43 : * \frac{\text{max}\left(\{|\lambda_\text{int}|\},
44 : * \{|\lambda_\text{ext}|\}\right)}{2} \left(U_\text{ext} - U_\text{int}\right),
45 : * \f}
46 : *
47 : * where "int" and "ext" stand for interior and exterior, and
48 : * \f$\{|\lambda|\}\f$ is the set of characteristic/signal speeds. The minus
49 : * sign in front of the \f$F_{\text{ext}}\f$ is necessary because the outward
50 : * directed normal of the neighboring element has the opposite sign, i.e.
51 : * \f$n_i^{\text{ext}}=-n_i^{\text{int}}\f$. The characteristic/signal speeds
52 : * are given by:
53 : *
54 : * \f{align*}
55 : * \lambda_{\pm}&=v^i n_i \pm c_s, \\
56 : * \lambda_v&=v^i n_i
57 : * \f}
58 : *
59 : * where \f$v^i\f$ is the spatial velocity and \f$c_s\f$ the sound speed.
60 : *
61 : * \note In the strong form the `dg_boundary_terms` function returns
62 : * \f$G - F_\text{int}\f$
63 : */
64 : template <size_t Dim>
65 1 : class Rusanov final : public evolution::BoundaryCorrection {
66 : private:
67 0 : struct AbsCharSpeed : db::SimpleTag {
68 0 : using type = Scalar<DataVector>;
69 : };
70 :
71 : public:
72 0 : using options = tmpl::list<>;
73 0 : static constexpr Options::String help = {
74 : "Computes the Rusanov or local Lax-Friedrichs boundary correction term "
75 : "for the Newtonian Euler/hydrodynamics system."};
76 :
77 0 : Rusanov() = default;
78 0 : Rusanov(const Rusanov&) = default;
79 0 : Rusanov& operator=(const Rusanov&) = default;
80 0 : Rusanov(Rusanov&&) = default;
81 0 : Rusanov& operator=(Rusanov&&) = default;
82 0 : ~Rusanov() override = default;
83 :
84 : /// \cond
85 : explicit Rusanov(CkMigrateMessage* msg);
86 : using PUP::able::register_constructor;
87 : WRAPPED_PUPable_decl_template(Rusanov); // NOLINT
88 : /// \endcond
89 0 : void pup(PUP::er& p) override; // NOLINT
90 :
91 0 : std::unique_ptr<BoundaryCorrection> get_clone() const override;
92 :
93 0 : using dg_package_field_tags =
94 : tmpl::list<Tags::MassDensityCons, Tags::MomentumDensity<Dim>,
95 : Tags::EnergyDensity,
96 : ::Tags::NormalDotFlux<Tags::MassDensityCons>,
97 : ::Tags::NormalDotFlux<Tags::MomentumDensity<Dim>>,
98 : ::Tags::NormalDotFlux<Tags::EnergyDensity>, AbsCharSpeed>;
99 0 : using dg_package_data_temporary_tags = tmpl::list<>;
100 0 : using dg_package_data_primitive_tags =
101 : tmpl::list<hydro::Tags::SpatialVelocity<DataVector, Dim>,
102 : hydro::Tags::SpecificInternalEnergy<DataVector>>;
103 0 : using dg_package_data_volume_tags =
104 : tmpl::list<hydro::Tags::EquationOfState<false, 2>>;
105 0 : using dg_boundary_terms_volume_tags = tmpl::list<>;
106 :
107 0 : double dg_package_data(
108 : gsl::not_null<Scalar<DataVector>*> packaged_mass_density,
109 : gsl::not_null<tnsr::I<DataVector, Dim, Frame::Inertial>*>
110 : packaged_momentum_density,
111 : gsl::not_null<Scalar<DataVector>*> packaged_energy_density,
112 : gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_mass_density,
113 : gsl::not_null<tnsr::I<DataVector, Dim, Frame::Inertial>*>
114 : packaged_normal_dot_flux_momentum_density,
115 : gsl::not_null<Scalar<DataVector>*>
116 : packaged_normal_dot_flux_energy_density,
117 : gsl::not_null<Scalar<DataVector>*> packaged_abs_char_speed,
118 :
119 : const Scalar<DataVector>& mass_density,
120 : const tnsr::I<DataVector, Dim, Frame::Inertial>& momentum_density,
121 : const Scalar<DataVector>& energy_density,
122 :
123 : const tnsr::I<DataVector, Dim, Frame::Inertial>& flux_mass_density,
124 : const tnsr::IJ<DataVector, Dim, Frame::Inertial>& flux_momentum_density,
125 : const tnsr::I<DataVector, Dim, Frame::Inertial>& flux_energy_density,
126 :
127 : const tnsr::I<DataVector, Dim, Frame::Inertial>& velocity,
128 : const Scalar<DataVector>& specific_internal_energy,
129 :
130 : const tnsr::i<DataVector, Dim, Frame::Inertial>& normal_covector,
131 : const std::optional<tnsr::I<DataVector, Dim, Frame::Inertial>>&
132 : /*mesh_velocity*/,
133 : const std::optional<Scalar<DataVector>>& normal_dot_mesh_velocity,
134 : const EquationsOfState::EquationOfState<false, 2>& equation_of_state)
135 : const;
136 :
137 0 : void dg_boundary_terms(
138 : gsl::not_null<Scalar<DataVector>*> boundary_correction_mass_density,
139 : gsl::not_null<tnsr::I<DataVector, Dim, Frame::Inertial>*>
140 : boundary_correction_momentum_density,
141 : gsl::not_null<Scalar<DataVector>*> boundary_correction_energy_density,
142 : const Scalar<DataVector>& mass_density_int,
143 : const tnsr::I<DataVector, Dim, Frame::Inertial>& momentum_density_int,
144 : const Scalar<DataVector>& energy_density_int,
145 : const Scalar<DataVector>& normal_dot_flux_mass_density_int,
146 : const tnsr::I<DataVector, Dim, Frame::Inertial>&
147 : normal_dot_flux_momentum_density_int,
148 : const Scalar<DataVector>& normal_dot_flux_energy_density_int,
149 : const Scalar<DataVector>& abs_char_speed_int,
150 : const Scalar<DataVector>& mass_density_ext,
151 : const tnsr::I<DataVector, Dim, Frame::Inertial>& momentum_density_ext,
152 : const Scalar<DataVector>& energy_density_ext,
153 : const Scalar<DataVector>& normal_dot_flux_mass_density_ext,
154 : const tnsr::I<DataVector, Dim, Frame::Inertial>&
155 : normal_dot_flux_momentum_density_ext,
156 : const Scalar<DataVector>& normal_dot_flux_energy_density_ext,
157 : const Scalar<DataVector>& abs_char_speed_ext,
158 : dg::Formulation dg_formulation) const;
159 : };
160 : } // namespace NewtonianEuler::BoundaryCorrections
|