SpECTRE Documentation Coverage Report
Current view: top level - Evolution/Systems/ForceFree/BoundaryCorrections - Rusanov.hpp Hit Total Coverage
Commit: 3c072f0ce967e2e56649d3fa12aa2a0e4fe2a42e Lines: 1 23 4.3 %
Date: 2024-04-23 20:50:18
Legend: Lines: hit not hit

          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             : #include <pup.h>
       9             : 
      10             : #include "DataStructures/DataBox/Prefixes.hpp"
      11             : #include "DataStructures/Tensor/TypeAliases.hpp"
      12             : #include "Evolution/Systems/ForceFree/BoundaryCorrections/BoundaryCorrection.hpp"
      13             : #include "Evolution/Systems/ForceFree/Tags.hpp"
      14             : #include "NumericalAlgorithms/DiscontinuousGalerkin/Formulation.hpp"
      15             : #include "Options/String.hpp"
      16             : #include "PointwiseFunctions/GeneralRelativity/TagsDeclarations.hpp"
      17             : #include "Utilities/Serialization/CharmPupable.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             : 
      31             : namespace ForceFree::BoundaryCorrections {
      32             : 
      33             : /*!
      34             :  * \brief A Rusanov/local Lax-Friedrichs Riemann solver
      35             :  *
      36             :  * Let \f$U\f$ be the evolved variables, \f$F^i\f$ the corresponding fluxes, and
      37             :  * \f$n_i\f$ be the outward directed unit normal to the interface. Denoting \f$F
      38             :  * := n_i F^i\f$, the %Rusanov boundary correction 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 characteristic/signal speed. The minus sign in
      48             :  * front of the \f$F_{\text{ext}}\f$ is necessary because the outward directed
      49             :  * normal of the neighboring element has the opposite sign, i.e.
      50             :  * \f$n_i^{\text{ext}}=-n_i^{\text{int}}\f$.
      51             :  *
      52             :  * For the GRFFE system the largest characteristic speeds \f$\lambda\f$ of our
      53             :  * interest are given as
      54             :  *
      55             :  * \f{align*}{
      56             :  *   \lambda_{\pm} = -\beta^i n_i \pm \alpha.
      57             :  * \f}
      58             :  *
      59             :  * which correspond to fast mode waves.
      60             :  *
      61             :  * \note In the strong form the `dg_boundary_terms` function returns
      62             :  * \f$G - F_\text{int}\f$
      63             :  *
      64             :  */
      65           1 : class Rusanov final : public 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 GRFFE 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* /*unused*/);
      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::TildeE, Tags::TildeB, Tags::TildePsi, Tags::TildePhi,
      95             :                  Tags::TildeQ, ::Tags::NormalDotFlux<Tags::TildeE>,
      96             :                  ::Tags::NormalDotFlux<Tags::TildeB>,
      97             :                  ::Tags::NormalDotFlux<Tags::TildePsi>,
      98             :                  ::Tags::NormalDotFlux<Tags::TildePhi>,
      99             :                  ::Tags::NormalDotFlux<Tags::TildeQ>, AbsCharSpeed>;
     100           0 :   using dg_package_data_temporary_tags =
     101             :       tmpl::list<gr::Tags::Lapse<DataVector>, gr::Tags::Shift<DataVector, 3>>;
     102           0 :   using dg_package_data_primitive_tags = tmpl::list<>;
     103           0 :   using dg_package_data_volume_tags = tmpl::list<>;
     104           0 :   using dg_boundary_terms_volume_tags = tmpl::list<>;
     105             : 
     106           0 :   static double dg_package_data(
     107             :       gsl::not_null<tnsr::I<DataVector, 3, Frame::Inertial>*> packaged_tilde_e,
     108             :       gsl::not_null<tnsr::I<DataVector, 3, Frame::Inertial>*> packaged_tilde_b,
     109             :       gsl::not_null<Scalar<DataVector>*> packaged_tilde_psi,
     110             :       gsl::not_null<Scalar<DataVector>*> packaged_tilde_phi,
     111             :       gsl::not_null<Scalar<DataVector>*> packaged_tilde_q,
     112             :       gsl::not_null<tnsr::I<DataVector, 3, Frame::Inertial>*>
     113             :           packaged_normal_dot_flux_tilde_e,
     114             :       gsl::not_null<tnsr::I<DataVector, 3, Frame::Inertial>*>
     115             :           packaged_normal_dot_flux_tilde_b,
     116             :       gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_tilde_psi,
     117             :       gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_tilde_phi,
     118             :       gsl::not_null<Scalar<DataVector>*> packaged_normal_dot_flux_tilde_q,
     119             :       gsl::not_null<Scalar<DataVector>*> packaged_abs_char_speed,
     120             : 
     121             :       const tnsr::I<DataVector, 3, Frame::Inertial>& tilde_e,
     122             :       const tnsr::I<DataVector, 3, Frame::Inertial>& tilde_b,
     123             :       const Scalar<DataVector>& tilde_psi, const Scalar<DataVector>& tilde_phi,
     124             :       const Scalar<DataVector>& tilde_q,
     125             : 
     126             :       const tnsr::IJ<DataVector, 3, Frame::Inertial>& flux_tilde_e,
     127             :       const tnsr::IJ<DataVector, 3, Frame::Inertial>& flux_tilde_b,
     128             :       const tnsr::I<DataVector, 3, Frame::Inertial>& flux_tilde_psi,
     129             :       const tnsr::I<DataVector, 3, Frame::Inertial>& flux_tilde_phi,
     130             :       const tnsr::I<DataVector, 3, Frame::Inertial>& flux_tilde_q,
     131             : 
     132             :       const Scalar<DataVector>& lapse,
     133             :       const tnsr::I<DataVector, 3, Frame::Inertial>& shift,
     134             : 
     135             :       const tnsr::i<DataVector, 3, Frame::Inertial>& normal_covector,
     136             :       const tnsr::I<DataVector, 3, Frame::Inertial>& normal_vector,
     137             :       const std::optional<tnsr::I<DataVector, 3, Frame::Inertial>>&
     138             :       /*mesh_velocity*/,
     139             :       const std::optional<Scalar<DataVector>>& normal_dot_mesh_velocity);
     140             : 
     141           0 :   static void dg_boundary_terms(
     142             :       gsl::not_null<tnsr::I<DataVector, 3, Frame::Inertial>*>
     143             :           boundary_correction_tilde_e,
     144             :       gsl::not_null<tnsr::I<DataVector, 3, Frame::Inertial>*>
     145             :           boundary_correction_tilde_b,
     146             :       gsl::not_null<Scalar<DataVector>*> boundary_correction_tilde_psi,
     147             :       gsl::not_null<Scalar<DataVector>*> boundary_correction_tilde_phi,
     148             :       gsl::not_null<Scalar<DataVector>*> boundary_correction_tilde_q,
     149             : 
     150             :       const tnsr::I<DataVector, 3, Frame::Inertial>& tilde_e_int,
     151             :       const tnsr::I<DataVector, 3, Frame::Inertial>& tilde_b_int,
     152             :       const Scalar<DataVector>& tilde_psi_int,
     153             :       const Scalar<DataVector>& tilde_phi_int,
     154             :       const Scalar<DataVector>& tilde_q_int,
     155             :       const tnsr::I<DataVector, 3, Frame::Inertial>&
     156             :           normal_dot_flux_tilde_e_int,
     157             :       const tnsr::I<DataVector, 3, Frame::Inertial>&
     158             :           normal_dot_flux_tilde_b_int,
     159             :       const Scalar<DataVector>& normal_dot_flux_tilde_psi_int,
     160             :       const Scalar<DataVector>& normal_dot_flux_tilde_phi_int,
     161             :       const Scalar<DataVector>& normal_dot_flux_tilde_q_int,
     162             :       const Scalar<DataVector>& abs_char_speed_int,
     163             : 
     164             :       const tnsr::I<DataVector, 3, Frame::Inertial>& tilde_e_ext,
     165             :       const tnsr::I<DataVector, 3, Frame::Inertial>& tilde_b_ext,
     166             :       const Scalar<DataVector>& tilde_psi_ext,
     167             :       const Scalar<DataVector>& tilde_phi_ext,
     168             :       const Scalar<DataVector>& tilde_q_ext,
     169             :       const tnsr::I<DataVector, 3, Frame::Inertial>&
     170             :           normal_dot_flux_tilde_e_ext,
     171             :       const tnsr::I<DataVector, 3, Frame::Inertial>&
     172             :           normal_dot_flux_tilde_b_ext,
     173             :       const Scalar<DataVector>& normal_dot_flux_tilde_psi_ext,
     174             :       const Scalar<DataVector>& normal_dot_flux_tilde_phi_ext,
     175             :       const Scalar<DataVector>& normal_dot_flux_tilde_q_ext,
     176             :       const Scalar<DataVector>& abs_char_speed_ext,
     177             :       dg::Formulation dg_formulation);
     178             : };
     179             : 
     180           0 : bool operator==(const Rusanov& lhs, const Rusanov& rhs);
     181           0 : bool operator!=(const Rusanov& lhs, const Rusanov& rhs);
     182             : 
     183             : }  // namespace ForceFree::BoundaryCorrections

Generated by: LCOV version 1.14