Rusanov.hpp
1 // Distributed under the MIT License.
2 // See LICENSE.txt for details.
3 
4 #pragma once
5 
6 #include <memory>
7 #include <optional>
8 
11 #include "Evolution/Systems/RadiationTransport/M1Grey/BoundaryCorrections/BoundaryCorrection.hpp"
12 #include "Evolution/Systems/RadiationTransport/M1Grey/Tags.hpp"
13 #include "NumericalAlgorithms/DiscontinuousGalerkin/Formulation.hpp"
14 #include "Options/Options.hpp"
16 #include "PointwiseFunctions/Hydro/EquationsOfState/EquationOfState.hpp"
17 #include "PointwiseFunctions/Hydro/Tags.hpp"
18 #include "Utilities/Gsl.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 
33 
34 namespace Rusanov_detail {
35 void dg_package_data_impl(
36  gsl::not_null<Scalar<DataVector>*> packaged_tilde_e,
37  gsl::not_null<tnsr::i<DataVector, 3, Frame::Inertial>*> packaged_tilde_s,
38  gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_tilde_e,
39  gsl::not_null<tnsr::i<DataVector, 3, Frame::Inertial>*>
40  packaged_normal_dot_flux_tilde_s,
41  const Scalar<DataVector>& tilde_e,
42  const tnsr::i<DataVector, 3, Frame::Inertial>& tilde_s,
43  const tnsr::I<DataVector, 3, Frame::Inertial>& flux_tilde_e,
44  const tnsr::Ij<DataVector, 3, Frame::Inertial>& flux_tilde_s,
45  const tnsr::i<DataVector, 3, Frame::Inertial>& normal_covector,
46  const tnsr::I<DataVector, 3, Frame::Inertial>& normal_vector,
47  const std::optional<tnsr::I<DataVector, 3, Frame::Inertial>>& mesh_velocity,
48  const std::optional<Scalar<DataVector>>& normal_dot_mesh_velocity) noexcept;
49 
50 void dg_boundary_terms_impl(
51  gsl::not_null<Scalar<DataVector>*> boundary_correction_tilde_e,
52  gsl::not_null<tnsr::i<DataVector, 3, Frame::Inertial>*>
53  boundary_correction_tilde_s,
54  const Scalar<DataVector>& tilde_e_int,
55  const tnsr::i<DataVector, 3, Frame::Inertial>& tilde_s_int,
56  const Scalar<DataVector>& normal_dot_flux_tilde_e_int,
57  const tnsr::i<DataVector, 3, Frame::Inertial>& normal_dot_flux_tilde_s_int,
58  const Scalar<DataVector>& tilde_e_ext,
59  const tnsr::i<DataVector, 3, Frame::Inertial>& tilde_s_ext,
60  const Scalar<DataVector>& normal_dot_flux_tilde_e_ext,
61  const tnsr::i<DataVector, 3, Frame::Inertial>& normal_dot_flux_tilde_s_ext,
62  dg::Formulation dg_formulation) noexcept;
63 } // namespace Rusanov_detail
64 
65 /*!
66  * \brief A Rusanov/local Lax-Friedrichs Riemann solver
67  *
68  * Let \f$U\f$ be the state vector of evolved variables, \f$F^i\f$ the
69  * corresponding fluxes, and \f$n_i\f$ be the outward directed unit normal to
70  * the interface. Denoting \f$F := n_i F^i\f$, the %Rusanov boundary correction
71  * is
72  *
73  * \f{align*}
74  * G_\text{Rusanov} = \frac{F_\text{int} - F_\text{ext}}{2} -
75  * \frac{\text{max}\left(\{|\lambda_\text{int}|\},
76  * \{|\lambda_\text{ext}|\}\right)}{2} \left(U_\text{ext} - U_\text{int}\right),
77  * \f}
78  *
79  * where "int" and "ext" stand for interior and exterior, and
80  * \f$\{|\lambda|\}\f$ is the set of characteristic/signal speeds. The minus
81  * sign in front of the \f$F_{\text{ext}}\f$ is necessary because the outward
82  * directed normal of the neighboring element has the opposite sign, i.e.
83  * \f$n_i^{\text{ext}}=-n_i^{\text{int}}\f$.
84  *
85  * For radiation, the max characteristic/signal speed is 1 (the speed of light).
86  * Note that the characteristic speeds of this system are *not* yet fully
87  * implemented, see the warning in the documentation for the characteristics.
88  *
89  * \note In the strong form the `dg_boundary_terms` function returns
90  * \f$G - F_\text{int}\f$
91  */
92 template <typename... NeutrinoSpecies>
93 class Rusanov final : public BoundaryCorrection<NeutrinoSpecies...> {
94  public:
95  using options = tmpl::list<>;
96  static constexpr Options::String help = {
97  "Computes the Rusanov or local Lax-Friedrichs boundary correction term "
98  "for the M1Grey radiation transport system."};
99 
100  Rusanov() = default;
101  Rusanov(const Rusanov&) = default;
102  Rusanov& operator=(const Rusanov&) = default;
103  Rusanov(Rusanov&&) = default;
104  Rusanov& operator=(Rusanov&&) = default;
105  ~Rusanov() override = default;
106 
107  /// \cond
108  explicit Rusanov(CkMigrateMessage* msg) noexcept
109  : BoundaryCorrection<NeutrinoSpecies...>(msg) {}
110  using PUP::able::register_constructor;
112  /// \endcond
113  void pup(PUP::er& p) override { // NOLINT
115  }
116 
117  std::unique_ptr<BoundaryCorrection<NeutrinoSpecies...>> get_clone()
118  const noexcept override {
119  return std::make_unique<Rusanov>(*this);
120  }
121 
122  using dg_package_field_tags = tmpl::list<
127  using dg_package_data_temporary_tags = tmpl::list<>;
128  using dg_package_data_primitive_tags = tmpl::list<>;
129  using dg_package_data_volume_tags = tmpl::list<>;
130 
131  double dg_package_data(
132  const gsl::not_null<typename Tags::TildeE<
133  Frame::Inertial, NeutrinoSpecies>::type*>... packaged_tilde_e,
134  const gsl::not_null<typename Tags::TildeS<
135  Frame::Inertial, NeutrinoSpecies>::type*>... packaged_tilde_s,
136  const gsl::not_null<typename Tags::TildeE<
138  NeutrinoSpecies>::type*>... packaged_normal_dot_flux_tilde_e,
139  const gsl::not_null<typename Tags::TildeS<
141  NeutrinoSpecies>::type*>... packaged_normal_dot_flux_tilde_s,
142 
143  const typename Tags::TildeE<Frame::Inertial,
144  NeutrinoSpecies>::type&... tilde_e,
145  const typename Tags::TildeS<Frame::Inertial,
146  NeutrinoSpecies>::type&... tilde_s,
147  const typename ::Tags::Flux<
149  Frame::Inertial>::type&... flux_tilde_e,
150  const typename ::Tags::Flux<
152  Frame::Inertial>::type&... flux_tilde_s,
153 
154  const tnsr::i<DataVector, 3, Frame::Inertial>& normal_covector,
155  const tnsr::I<DataVector, 3, Frame::Inertial>& normal_vector,
156  const std::optional<tnsr::I<DataVector, 3, Frame::Inertial>>&
157  mesh_velocity,
158  const std::optional<Scalar<DataVector>>& normal_dot_mesh_velocity)
159  const noexcept {
160  EXPAND_PACK_LEFT_TO_RIGHT(Rusanov_detail::dg_package_data_impl(
161  packaged_tilde_e, packaged_tilde_s, packaged_normal_dot_flux_tilde_e,
162  packaged_normal_dot_flux_tilde_s, tilde_e, tilde_s, flux_tilde_e,
163  flux_tilde_s, normal_covector, normal_vector, mesh_velocity,
164  normal_dot_mesh_velocity));
165  // max speed = speed of light
166  return 1.0;
167  }
168 
169  void dg_boundary_terms(
170  const gsl::not_null<typename Tags::TildeE<
172  NeutrinoSpecies>::type*>... boundary_correction_tilde_e,
173  const gsl::not_null<typename Tags::TildeS<
175  NeutrinoSpecies>::type*>... boundary_correction_tilde_s,
176  const typename Tags::TildeE<Frame::Inertial,
177  NeutrinoSpecies>::type&... tilde_e_int,
178  const typename Tags::TildeS<Frame::Inertial,
179  NeutrinoSpecies>::type&... tilde_s_int,
181  type&... normal_dot_flux_tilde_e_int,
183  type&... normal_dot_flux_tilde_s_int,
184  const typename Tags::TildeE<Frame::Inertial,
185  NeutrinoSpecies>::type&... tilde_e_ext,
186  const typename Tags::TildeS<Frame::Inertial,
187  NeutrinoSpecies>::type&... tilde_s_ext,
189  type&... normal_dot_flux_tilde_e_ext,
191  type&... normal_dot_flux_tilde_s_ext,
192  dg::Formulation dg_formulation) const noexcept {
193  EXPAND_PACK_LEFT_TO_RIGHT(Rusanov_detail::dg_boundary_terms_impl(
194  boundary_correction_tilde_e, boundary_correction_tilde_s, tilde_e_int,
195  tilde_s_int, normal_dot_flux_tilde_e_int, normal_dot_flux_tilde_s_int,
196  tilde_e_ext, tilde_s_ext, normal_dot_flux_tilde_e_ext,
197  normal_dot_flux_tilde_s_ext, dg_formulation));
198  }
199 };
200 
201 /// \cond
202 template <typename... NeutrinoSpecies>
203 // NOLINTNEXTLINE
204 PUP::able::PUP_ID Rusanov<NeutrinoSpecies...>::my_PUP_ID = 0;
205 /// \endcond
206 
207 } // namespace RadiationTransport::M1Grey::BoundaryCorrections
CharmPupable.hpp
EXPAND_PACK_LEFT_TO_RIGHT
#define EXPAND_PACK_LEFT_TO_RIGHT(...)
Expand a parameter pack evaluating the terms from left to right.
Definition: TMPL.hpp:565
Frame::Inertial
Definition: IndexType.hpp:44
Options.hpp
Tags::NormalDotFlux
Prefix indicating a boundary unit normal vector dotted into the flux.
Definition: Prefixes.hpp:96
RadiationTransport::M1Grey::BoundaryCorrections::Rusanov
A Rusanov/local Lax-Friedrichs Riemann solver.
Definition: Rusanov.hpp:93
RadiationTransport::M1Grey::Tags::TildeE
The densitized energy density of neutrinos of a given species .
Definition: Tags.hpp:31
RadiationTransport::M1Grey::Tags::TildeS
The densitized momentum density of neutrinos of a given species .
Definition: Tags.hpp:41
WRAPPED_PUPable_decl_template
#define WRAPPED_PUPable_decl_template(className)
Mark derived classes as serializable.
Definition: CharmPupable.hpp:22
RadiationTransport::M1Grey::BoundaryCorrections
Boundary corrections/numerical fluxes.
Definition: BoundaryCorrection.hpp:13
RadiationTransport::M1Grey::BoundaryCorrections::BoundaryCorrection
The base class used to create boundary corrections from input files and store them in the global cach...
Definition: BoundaryCorrection.hpp:24
DataVector
Stores a collection of function values.
Definition: DataVector.hpp:46
memory
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
optional
dg::Formulation
Formulation
The DG formulation to use.
Definition: Formulation.hpp:34
std::unique_ptr
Prefixes.hpp
gsl
Implementations from the Guideline Support Library.
Definition: ReadSpecThirdOrderPiecewisePolynomial.hpp:11
TMPL.hpp
gsl::not_null
Require a pointer to not be a nullptr
Definition: ReadSpecThirdOrderPiecewisePolynomial.hpp:13