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/GrMhd/ValenciaDivClean/BoundaryCorrections/BoundaryCorrection.hpp"
12 #include "Evolution/Systems/GrMhd/ValenciaDivClean/Tags.hpp"
13 #include "NumericalAlgorithms/DiscontinuousGalerkin/Formulation.hpp"
14 #include "Options/Options.hpp"
16 #include "PointwiseFunctions/GeneralRelativity/Tags.hpp"
17 #include "Utilities/Gsl.hpp"
18 #include "Utilities/TMPL.hpp"
19 
20 /// \cond
21 class DataVector;
22 namespace gsl {
23 template <typename T>
24 class not_null;
25 } // namespace gsl
26 namespace PUP {
27 class er;
28 } // namespace PUP
29 /// \endcond
30 
32 /*!
33  * \brief A Rusanov/local Lax-Friedrichs Riemann solver
34  *
35  * Let \f$U\f$ be the state vector of evolved variables, \f$F^i\f$ the
36  * corresponding fluxes, and \f$n_i\f$ be the outward directed unit normal to
37  * the interface. Denoting \f$F := n_i F^i\f$, the %Rusanov boundary correction
38  * is
39  *
40  * \f{align*}
41  * G_\text{Rusanov} = \frac{F_\text{int} - F_\text{ext}}{2} -
42  * \frac{\text{max}\left(\{|\lambda_\text{int}|\},
43  * \{|\lambda_\text{ext}|\}\right)}{2} \left(U_\text{ext} - U_\text{int}\right),
44  * \f}
45  *
46  * where "int" and "ext" stand for interior and exterior, and
47  * \f$\{|\lambda|\}\f$ is the set of characteristic/signal speeds. The minus
48  * sign in front of the \f$F_{\text{ext}}\f$ is necessary because the outward
49  * directed normal of the neighboring element has the opposite sign, i.e.
50  * \f$n_i^{\text{ext}}=-n_i^{\text{int}}\f$. The characteristic/signal speeds
51  * are given in the documentation for
52  * grmhd::ValenciaDivClean::characteristic_speeds(). Since the fluid is
53  * travelling slower than the speed of light, the speeds we are interested in
54  * are
55  *
56  * \f{align*}{
57  * \lambda_{\pm}&=\pm\alpha-\beta^i n_i,
58  * \f}
59  *
60  * which correspond to the divergence cleaning field.
61  *
62  * \note
63  * - In the strong form the `dg_boundary_terms` function returns
64  * \f$G - F_\text{int}\f$
65  * - It may be possible to use the slower speeds for the magnetic field and
66  * fluid part of the system in order to make the flux less dissipative for
67  * those variables.
68  */
69 class Rusanov final : public BoundaryCorrection {
70  private:
71  struct AbsCharSpeed : db::SimpleTag {
72  using type = Scalar<DataVector>;
73  };
74 
75  public:
76  using options = tmpl::list<>;
77  static constexpr Options::String help = {
78  "Computes the Rusanov or local Lax-Friedrichs boundary correction term "
79  "for the GRMHD system."};
80 
81  Rusanov() = default;
82  Rusanov(const Rusanov&) = default;
83  Rusanov& operator=(const Rusanov&) = default;
84  Rusanov(Rusanov&&) = default;
85  Rusanov& operator=(Rusanov&&) = default;
86  ~Rusanov() override = default;
87 
88  /// \cond
89  explicit Rusanov(CkMigrateMessage* /*unused*/) noexcept;
90  using PUP::able::register_constructor;
92  /// \endcond
93  void pup(PUP::er& p) override; // NOLINT
94 
95  std::unique_ptr<BoundaryCorrection> get_clone() const noexcept override;
96 
97  using dg_package_field_tags =
98  tmpl::list<Tags::TildeD, Tags::TildeTau, Tags::TildeS<Frame::Inertial>,
105  using dg_package_data_temporary_tags =
106  tmpl::list<gr::Tags::Lapse<DataVector>,
108  using dg_package_data_primitive_tags = tmpl::list<>;
109  using dg_package_data_volume_tags = tmpl::list<>;
110 
111  static double dg_package_data(
112  gsl::not_null<Scalar<DataVector>*> packaged_tilde_d,
113  gsl::not_null<Scalar<DataVector>*> packaged_tilde_tau,
114  gsl::not_null<tnsr::i<DataVector, 3, Frame::Inertial>*> packaged_tilde_s,
115  gsl::not_null<tnsr::I<DataVector, 3, Frame::Inertial>*> packaged_tilde_b,
116  gsl::not_null<Scalar<DataVector>*> packaged_tilde_phi,
117  gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_tilde_d,
118  gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_tilde_tau,
119  gsl::not_null<tnsr::i<DataVector, 3, Frame::Inertial>*>
120  packaged_normal_dot_flux_tilde_s,
121  gsl::not_null<tnsr::I<DataVector, 3, Frame::Inertial>*>
122  packaged_normal_dot_flux_tilde_b,
123  gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_tilde_phi,
124  gsl::not_null<Scalar<DataVector>*> packaged_abs_char_speed,
125 
126  const Scalar<DataVector>& tilde_d, const Scalar<DataVector>& tilde_tau,
127  const tnsr::i<DataVector, 3, Frame::Inertial>& tilde_s,
128  const tnsr::I<DataVector, 3, Frame::Inertial>& tilde_b,
129  const Scalar<DataVector>& tilde_phi,
130 
131  const tnsr::I<DataVector, 3, Frame::Inertial>& flux_tilde_d,
132  const tnsr::I<DataVector, 3, Frame::Inertial>& flux_tilde_tau,
133  const tnsr::Ij<DataVector, 3, Frame::Inertial>& flux_tilde_s,
134  const tnsr::IJ<DataVector, 3, Frame::Inertial>& flux_tilde_b,
135  const tnsr::I<DataVector, 3, Frame::Inertial>& flux_tilde_phi,
136 
137  const Scalar<DataVector>& lapse,
138  const tnsr::I<DataVector, 3, Frame::Inertial>& shift,
139 
140  const tnsr::i<DataVector, 3, Frame::Inertial>& normal_covector,
141  const tnsr::I<DataVector, 3, Frame::Inertial>& normal_vector,
142  const std::optional<tnsr::I<DataVector, 3, Frame::Inertial>>&
143  /*mesh_velocity*/,
145  normal_dot_mesh_velocity) noexcept;
146 
147  static void dg_boundary_terms(
148  gsl::not_null<Scalar<DataVector>*> boundary_correction_tilde_d,
149  gsl::not_null<Scalar<DataVector>*> boundary_correction_tilde_tau,
150  gsl::not_null<tnsr::i<DataVector, 3, Frame::Inertial>*>
151  boundary_correction_tilde_s,
152  gsl::not_null<tnsr::I<DataVector, 3, Frame::Inertial>*>
153  boundary_correction_tilde_b,
154  gsl::not_null<Scalar<DataVector>*> boundary_correction_tilde_phi,
155  const Scalar<DataVector>& tilde_d_int,
156  const Scalar<DataVector>& tilde_tau_int,
157  const tnsr::i<DataVector, 3, Frame::Inertial>& tilde_s_int,
158  const tnsr::I<DataVector, 3, Frame::Inertial>& tilde_b_int,
159  const Scalar<DataVector>& tilde_phi_int,
160  const Scalar<DataVector>& normal_dot_flux_tilde_d_int,
161  const Scalar<DataVector>& normal_dot_flux_tilde_tau_int,
162  const tnsr::i<DataVector, 3, Frame::Inertial>&
163  normal_dot_flux_tilde_s_int,
164  const tnsr::I<DataVector, 3, Frame::Inertial>&
165  normal_dot_flux_tilde_b_int,
166  const Scalar<DataVector>& normal_dot_flux_tilde_phi_int,
167  const Scalar<DataVector>& abs_char_speed_int,
168  const Scalar<DataVector>& tilde_d_ext,
169  const Scalar<DataVector>& tilde_tau_ext,
170  const tnsr::i<DataVector, 3, Frame::Inertial>& tilde_s_ext,
171  const tnsr::I<DataVector, 3, Frame::Inertial>& tilde_b_ext,
172  const Scalar<DataVector>& tilde_phi_ext,
173  const Scalar<DataVector>& normal_dot_flux_tilde_d_ext,
174  const Scalar<DataVector>& normal_dot_flux_tilde_tau_ext,
175  const tnsr::i<DataVector, 3, Frame::Inertial>&
176  normal_dot_flux_tilde_s_ext,
177  const tnsr::I<DataVector, 3, Frame::Inertial>&
178  normal_dot_flux_tilde_b_ext,
179  const Scalar<DataVector>& normal_dot_flux_tilde_phi_ext,
180  const Scalar<DataVector>& abs_char_speed_ext,
181  dg::Formulation dg_formulation) noexcept;
182 };
183 } // namespace grmhd::ValenciaDivClean::BoundaryCorrections
CharmPupable.hpp
Options.hpp
db::SimpleTag
Mark a struct as a simple tag by inheriting from this.
Definition: Tag.hpp:36
Tags::NormalDotFlux
Prefix indicating a boundary unit normal vector dotted into the flux.
Definition: Prefixes.hpp:96
gr::lapse
Scalar< DataType > lapse(const tnsr::I< DataType, SpatialDim, Frame > &shift, const tnsr::aa< DataType, SpatialDim, Frame > &spacetime_metric) noexcept
Compute lapse from shift and spacetime metric.
WRAPPED_PUPable_decl_template
#define WRAPPED_PUPable_decl_template(className)
Mark derived classes as serializable.
Definition: CharmPupable.hpp:22
DataVector
Stores a collection of function values.
Definition: DataVector.hpp:46
memory
gr::Tags::Shift
Definition: Tags.hpp:48
grmhd::ValenciaDivClean::Tags::TildeB
The densitized magnetic field .
Definition: Tags.hpp:49
gr::shift
tnsr::I< DataType, SpatialDim, Frame > shift(const tnsr::aa< DataType, SpatialDim, Frame > &spacetime_metric, const tnsr::II< DataType, SpatialDim, Frame > &inverse_spatial_metric) noexcept
Compute shift from spacetime metric and inverse spatial metric.
grmhd::ValenciaDivClean::BoundaryCorrections
Boundary corrections/numerical fluxes.
Definition: BoundaryCorrection.hpp:13
Scalar
Tensor< T, Symmetry<>, index_list<> > Scalar
Definition: TypeAliases.hpp:21
grmhd::ValenciaDivClean::BoundaryCorrections::Rusanov
A Rusanov/local Lax-Friedrichs Riemann solver.
Definition: Rusanov.hpp:69
Gsl.hpp
TypeAliases.hpp
Options::String
const char *const String
The string used in option structs.
Definition: Options.hpp:32
optional
grmhd::ValenciaDivClean::Tags::TildePhi
The densitized divergence-cleaning field .
Definition: Tags.hpp:55
dg::Formulation
Formulation
The DG formulation to use.
Definition: Formulation.hpp:34
std::unique_ptr
Prefixes.hpp
grmhd::ValenciaDivClean::BoundaryCorrections::BoundaryCorrection
The base class used to make boundary corrections factory createable so they can be specified in the i...
Definition: BoundaryCorrection.hpp:22
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