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/Systems/RelativisticEuler/Valencia/BoundaryCorrections/BoundaryCorrection.hpp"
12 : #include "Evolution/Systems/RelativisticEuler/Valencia/Tags.hpp"
13 : #include "NumericalAlgorithms/DiscontinuousGalerkin/Formulation.hpp"
14 : #include "Options/String.hpp"
15 : #include "PointwiseFunctions/GeneralRelativity/Tags.hpp"
16 : #include "PointwiseFunctions/Hydro/EquationsOfState/EquationOfState.hpp"
17 : #include "PointwiseFunctions/Hydro/Tags.hpp"
18 : #include "Utilities/Gsl.hpp"
19 : #include "Utilities/Serialization/CharmPupable.hpp"
20 : #include "Utilities/TMPL.hpp"
21 :
22 : /// \cond
23 : class DataVector;
24 : namespace gsl {
25 : template <typename T>
26 : class not_null;
27 : } // namespace gsl
28 : namespace PUP {
29 : class er;
30 : } // namespace PUP
31 : /// \endcond
32 :
33 : namespace RelativisticEuler::Valencia::BoundaryCorrections {
34 : /*!
35 : * \brief A Rusanov/local Lax-Friedrichs Riemann solver
36 : *
37 : * Let \f$U\f$ be the state vector of evolved variables, \f$F^i\f$ the
38 : * corresponding fluxes, and \f$n_i\f$ be the outward directed unit normal to
39 : * the interface. Denoting \f$F := n_i F^i\f$, the %Rusanov boundary correction
40 : * is
41 : *
42 : * \f{align*}
43 : * G_\text{Rusanov} = \frac{F_\text{int} - F_\text{ext}}{2} -
44 : * \frac{\text{max}\left(\{|\lambda_\text{int}|\},
45 : * \{|\lambda_\text{ext}|\}\right)}{2} \left(U_\text{ext} - U_\text{int}\right),
46 : * \f}
47 : *
48 : * where "int" and "ext" stand for interior and exterior, and
49 : * \f$\{|\lambda|\}\f$ is the set of characteristic/signal speeds. The minus
50 : * sign in front of the \f$F_{\text{ext}}\f$ is necessary because the outward
51 : * directed normal of the neighboring element has the opposite sign, i.e.
52 : * \f$n_i^{\text{ext}}=-n_i^{\text{int}}\f$. The characteristic/signal speeds
53 : * expressions are listed in the documentation of `characteristic_speeds`.
54 : *
55 : * \note In the strong form the `dg_boundary_terms` function returns
56 : * \f$G - F_\text{int}\f$
57 : */
58 : template <size_t Dim>
59 1 : class Rusanov final : public BoundaryCorrection<Dim> {
60 : private:
61 0 : struct AbsCharSpeed : db::SimpleTag {
62 0 : using type = Scalar<DataVector>;
63 : };
64 :
65 : public:
66 0 : using options = tmpl::list<>;
67 0 : static constexpr Options::String help = {
68 : "Computes the Rusanov or local Lax-Friedrichs boundary correction term "
69 : "for the Valencia formulation of the relativistic Euler/hydrodynamics "
70 : "system."};
71 :
72 0 : Rusanov() = default;
73 0 : Rusanov(const Rusanov&) = default;
74 0 : Rusanov& operator=(const Rusanov&) = default;
75 0 : Rusanov(Rusanov&&) = default;
76 0 : Rusanov& operator=(Rusanov&&) = default;
77 0 : ~Rusanov() override = default;
78 :
79 : /// \cond
80 : explicit Rusanov(CkMigrateMessage* msg);
81 : using PUP::able::register_constructor;
82 : WRAPPED_PUPable_decl_template(Rusanov); // NOLINT
83 : /// \endcond
84 0 : void pup(PUP::er& p) override; // NOLINT
85 :
86 0 : std::unique_ptr<BoundaryCorrection<Dim>> get_clone() const override;
87 :
88 0 : using dg_package_field_tags =
89 : tmpl::list<Tags::TildeD, Tags::TildeTau, Tags::TildeS<Dim>,
90 : ::Tags::NormalDotFlux<Tags::TildeD>,
91 : ::Tags::NormalDotFlux<Tags::TildeTau>,
92 : ::Tags::NormalDotFlux<Tags::TildeS<Dim>>, AbsCharSpeed>;
93 0 : using dg_package_data_temporary_tags =
94 : tmpl::list<gr::Tags::Lapse<DataVector>, gr::Tags::Shift<DataVector, Dim>,
95 : gr::Tags::SpatialMetric<DataVector, Dim>>;
96 0 : using dg_package_data_primitive_tags =
97 : tmpl::list<hydro::Tags::RestMassDensity<DataVector>,
98 : hydro::Tags::SpecificInternalEnergy<DataVector>,
99 : hydro::Tags::SpecificEnthalpy<DataVector>,
100 : hydro::Tags::SpatialVelocity<DataVector, Dim>>;
101 0 : using dg_package_data_volume_tags =
102 : tmpl::list<hydro::Tags::EquationOfStateBase>;
103 0 : using dg_boundary_terms_volume_tags = tmpl::list<>;
104 :
105 : template <size_t ThermodynamicDim>
106 0 : double dg_package_data(
107 : gsl::not_null<Scalar<DataVector>*> packaged_tilde_d,
108 : gsl::not_null<Scalar<DataVector>*> packaged_tilde_tau,
109 : gsl::not_null<tnsr::i<DataVector, Dim, Frame::Inertial>*>
110 : packaged_tilde_s,
111 : gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_tilde_d,
112 : gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_tilde_tau,
113 : gsl::not_null<tnsr::i<DataVector, Dim, Frame::Inertial>*>
114 : packaged_normal_dot_flux_tilde_s,
115 : gsl::not_null<Scalar<DataVector>*> packaged_abs_char_speed,
116 :
117 : const Scalar<DataVector>& tilde_d, const Scalar<DataVector>& tilde_tau,
118 : const tnsr::i<DataVector, Dim, Frame::Inertial>& tilde_s,
119 :
120 : const tnsr::I<DataVector, Dim, Frame::Inertial>& flux_tilde_d,
121 : const tnsr::I<DataVector, Dim, Frame::Inertial>& flux_tilde_tau,
122 : const tnsr::Ij<DataVector, Dim, Frame::Inertial>& flux_tilde_s,
123 :
124 : const Scalar<DataVector>& lapse,
125 : const tnsr::I<DataVector, Dim, Frame::Inertial>& shift,
126 : const tnsr::ii<DataVector, Dim, Frame::Inertial>& spatial_metric,
127 : const Scalar<DataVector>& rest_mass_density,
128 : const Scalar<DataVector>& specific_internal_energy,
129 : const Scalar<DataVector>& specific_enthalpy,
130 : const tnsr::I<DataVector, Dim, Frame::Inertial>& spatial_velocity,
131 :
132 : const tnsr::i<DataVector, Dim, Frame::Inertial>& normal_covector,
133 : const tnsr::I<DataVector, Dim, Frame::Inertial>& normal_vector,
134 : const std::optional<tnsr::I<DataVector, Dim, Frame::Inertial>>&
135 : /*mesh_velocity*/,
136 : const std::optional<Scalar<DataVector>>& normal_dot_mesh_velocity,
137 : const EquationsOfState::EquationOfState<true, ThermodynamicDim>&
138 : equation_of_state) const;
139 :
140 0 : void dg_boundary_terms(
141 : gsl::not_null<Scalar<DataVector>*> boundary_correction_tilde_d,
142 : gsl::not_null<Scalar<DataVector>*> boundary_correction_tilde_tau,
143 : gsl::not_null<tnsr::i<DataVector, Dim, Frame::Inertial>*>
144 : boundary_correction_tilde_s,
145 : const Scalar<DataVector>& tilde_d_int,
146 : const Scalar<DataVector>& tilde_tau_int,
147 : const tnsr::i<DataVector, Dim, Frame::Inertial>& tilde_s_int,
148 : const Scalar<DataVector>& normal_dot_flux_tilde_d_int,
149 : const Scalar<DataVector>& normal_dot_flux_tilde_tau_int,
150 : const tnsr::i<DataVector, Dim, Frame::Inertial>&
151 : normal_dot_flux_tilde_s_int,
152 : const Scalar<DataVector>& abs_char_speed_int,
153 : const Scalar<DataVector>& tilde_d_ext,
154 : const Scalar<DataVector>& tilde_tau_ext,
155 : const tnsr::i<DataVector, Dim, Frame::Inertial>& tilde_s_ext,
156 : const Scalar<DataVector>& normal_dot_flux_tilde_d_ext,
157 : const Scalar<DataVector>& normal_dot_flux_tilde_tau_ext,
158 : const tnsr::i<DataVector, Dim, Frame::Inertial>&
159 : normal_dot_flux_tilde_s_ext,
160 : const Scalar<DataVector>& abs_char_speed_ext,
161 : dg::Formulation dg_formulation) const;
162 : };
163 : } // namespace RelativisticEuler::Valencia::BoundaryCorrections
|