SpECTRE  v2025.03.17
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
grmhd::ValenciaDivClean::BoundaryCorrections::Rusanov Class Referencefinal

A Rusanov/local Lax-Friedrichs Riemann solver. More...

#include <Rusanov.hpp>

Public Types

using options = implementation defined
 
using dg_package_field_tags = implementation defined
 
using dg_package_data_temporary_tags = implementation defined
 
using dg_package_data_primitive_tags = implementation defined
 
using dg_package_data_volume_tags = implementation defined
 
using dg_boundary_terms_volume_tags = implementation defined
 
- Public Types inherited from grmhd::ValenciaDivClean::BoundaryCorrections::BoundaryCorrection
using creatable_classes = implementation defined
 

Public Member Functions

 Rusanov (const Rusanov &)=default
 
Rusanovoperator= (const Rusanov &)=default
 
 Rusanov (Rusanov &&)=default
 
Rusanovoperator= (Rusanov &&)=default
 
void pup (PUP::er &p) override
 
std::unique_ptr< BoundaryCorrectionget_clone () const override
 
- Public Member Functions inherited from grmhd::ValenciaDivClean::BoundaryCorrections::BoundaryCorrection
 BoundaryCorrection (const BoundaryCorrection &)=default
 
BoundaryCorrectionoperator= (const BoundaryCorrection &)=default
 
 BoundaryCorrection (BoundaryCorrection &&)=default
 
BoundaryCorrectionoperator= (BoundaryCorrection &&)=default
 
virtual std::unique_ptr< BoundaryCorrectionget_clone () const =0
 

Static Public Member Functions

static double dg_package_data (gsl::not_null< Scalar< DataVector > * > packaged_tilde_d, gsl::not_null< Scalar< DataVector > * > packaged_tilde_ye, gsl::not_null< Scalar< DataVector > * > packaged_tilde_tau, gsl::not_null< tnsr::i< DataVector, 3, Frame::Inertial > * > packaged_tilde_s, gsl::not_null< tnsr::I< DataVector, 3, Frame::Inertial > * > packaged_tilde_b, gsl::not_null< Scalar< DataVector > * > packaged_tilde_phi, gsl::not_null< Scalar< DataVector > * > packaged_normal_dot_flux_tilde_d, gsl::not_null< Scalar< DataVector > * > packaged_normal_dot_flux_tilde_ye, gsl::not_null< Scalar< DataVector > * > packaged_normal_dot_flux_tilde_tau, gsl::not_null< tnsr::i< DataVector, 3, Frame::Inertial > * > packaged_normal_dot_flux_tilde_s, gsl::not_null< tnsr::I< DataVector, 3, Frame::Inertial > * > packaged_normal_dot_flux_tilde_b, gsl::not_null< Scalar< DataVector > * > packaged_normal_dot_flux_tilde_phi, gsl::not_null< Scalar< DataVector > * > packaged_abs_char_speed, const Scalar< DataVector > &tilde_d, const Scalar< DataVector > &tilde_ye, const Scalar< DataVector > &tilde_tau, const tnsr::i< DataVector, 3, Frame::Inertial > &tilde_s, const tnsr::I< DataVector, 3, Frame::Inertial > &tilde_b, const Scalar< DataVector > &tilde_phi, const tnsr::I< DataVector, 3, Frame::Inertial > &flux_tilde_d, const tnsr::I< DataVector, 3, Frame::Inertial > &flux_tilde_ye, const tnsr::I< DataVector, 3, Frame::Inertial > &flux_tilde_tau, const tnsr::Ij< DataVector, 3, Frame::Inertial > &flux_tilde_s, const tnsr::IJ< DataVector, 3, Frame::Inertial > &flux_tilde_b, const tnsr::I< DataVector, 3, Frame::Inertial > &flux_tilde_phi, const Scalar< DataVector > &lapse, const tnsr::I< DataVector, 3, Frame::Inertial > &shift, const tnsr::i< DataVector, 3, Frame::Inertial > &, const Scalar< DataVector > &, const Scalar< DataVector > &, const Scalar< DataVector > &, const tnsr::I< DataVector, 3, Frame::Inertial > &, const tnsr::i< DataVector, 3, Frame::Inertial > &normal_covector, const tnsr::I< DataVector, 3, Frame::Inertial > &normal_vector, const std::optional< tnsr::I< DataVector, 3, Frame::Inertial > > &, const std::optional< Scalar< DataVector > > &normal_dot_mesh_velocity, const EquationsOfState::EquationOfState< true, 3 > &)
 
static void dg_boundary_terms (gsl::not_null< Scalar< DataVector > * > boundary_correction_tilde_d, gsl::not_null< Scalar< DataVector > * > boundary_correction_tilde_ye, gsl::not_null< Scalar< DataVector > * > boundary_correction_tilde_tau, gsl::not_null< tnsr::i< DataVector, 3, Frame::Inertial > * > boundary_correction_tilde_s, gsl::not_null< tnsr::I< DataVector, 3, Frame::Inertial > * > boundary_correction_tilde_b, gsl::not_null< Scalar< DataVector > * > boundary_correction_tilde_phi, const Scalar< DataVector > &tilde_d_int, const Scalar< DataVector > &tilde_ye_int, const Scalar< DataVector > &tilde_tau_int, const tnsr::i< DataVector, 3, Frame::Inertial > &tilde_s_int, const tnsr::I< DataVector, 3, Frame::Inertial > &tilde_b_int, const Scalar< DataVector > &tilde_phi_int, const Scalar< DataVector > &normal_dot_flux_tilde_d_int, const Scalar< DataVector > &normal_dot_flux_tilde_ye_int, const Scalar< DataVector > &normal_dot_flux_tilde_tau_int, const tnsr::i< DataVector, 3, Frame::Inertial > &normal_dot_flux_tilde_s_int, const tnsr::I< DataVector, 3, Frame::Inertial > &normal_dot_flux_tilde_b_int, const Scalar< DataVector > &normal_dot_flux_tilde_phi_int, const Scalar< DataVector > &abs_char_speed_int, const Scalar< DataVector > &tilde_d_ext, const Scalar< DataVector > &tilde_ye_ext, const Scalar< DataVector > &tilde_tau_ext, const tnsr::i< DataVector, 3, Frame::Inertial > &tilde_s_ext, const tnsr::I< DataVector, 3, Frame::Inertial > &tilde_b_ext, const Scalar< DataVector > &tilde_phi_ext, const Scalar< DataVector > &normal_dot_flux_tilde_d_ext, const Scalar< DataVector > &normal_dot_flux_tilde_ye_ext, const Scalar< DataVector > &normal_dot_flux_tilde_tau_ext, const tnsr::i< DataVector, 3, Frame::Inertial > &normal_dot_flux_tilde_s_ext, const tnsr::I< DataVector, 3, Frame::Inertial > &normal_dot_flux_tilde_b_ext, const Scalar< DataVector > &normal_dot_flux_tilde_phi_ext, const Scalar< DataVector > &abs_char_speed_ext, dg::Formulation dg_formulation)
 

Static Public Attributes

static constexpr Options::String help
 

Detailed Description

A Rusanov/local Lax-Friedrichs Riemann solver.

Let U be the state vector of evolved variables, Fi the corresponding fluxes, and ni be the outward directed unit normal to the interface. Denoting F:=niFi, the Rusanov boundary correction is

GRusanov=FintFext2max({|λint|},{|λext|})2(UextUint),

where "int" and "ext" stand for interior and exterior, and {|λ|} is the set of characteristic/signal speeds. The minus sign in front of the Fext is necessary because the outward directed normal of the neighboring element has the opposite sign, i.e. niext=niint. The characteristic/signal speeds are given in the documentation for grmhd::ValenciaDivClean::characteristic_speeds(). Since the fluid is travelling slower than the speed of light, the speeds we are interested in are

λ±=±αβini,

which correspond to the divergence cleaning field.

Note
  • In the strong form the dg_boundary_terms function returns GFint
  • It may be possible to use the slower speeds for the magnetic field and fluid part of the system in order to make the flux less dissipative for those variables.

Member Function Documentation

◆ get_clone()

std::unique_ptr< BoundaryCorrection > grmhd::ValenciaDivClean::BoundaryCorrections::Rusanov::get_clone ( ) const
overridevirtual

Member Data Documentation

◆ help

constexpr Options::String grmhd::ValenciaDivClean::BoundaryCorrections::Rusanov::help
staticconstexpr
Initial value:
= {
"Computes the Rusanov or local Lax-Friedrichs boundary correction term "
"for the GRMHD system."}

The documentation for this class was generated from the following file: