SpECTRE  v2026.04.01
Loading...
Searching...
No Matches
grmhd::ValenciaDivClean::BoundaryCorrections::Rusanov Class Referencefinal

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

#include <Rusanov.hpp>

Public Types

using options = tmpl::list<>
using dg_package_field_tags
using dg_package_data_temporary_tags
using dg_package_data_primitive_tags
using dg_package_data_volume_tags
using dg_boundary_terms_volume_tags = tmpl::list<>

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< BoundaryCorrection > get_clone () const override
Public Member Functions inherited from evolution::BoundaryCorrection
 BoundaryCorrection (const BoundaryCorrection &)=default
BoundaryCorrectionoperator= (const BoundaryCorrection &)=default
 BoundaryCorrection (BoundaryCorrection &&)=default
BoundaryCorrectionoperator= (BoundaryCorrection &&)=default
 BoundaryCorrection (CkMigrateMessage *msg)

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, \(F^i\) the corresponding fluxes, and \(n_i\) be the outward directed unit normal to the interface. Denoting \(F := n_i F^i\), the Rusanov boundary correction is

\begin{align*}G_\text{Rusanov} = \frac{F_\text{int} - F_\text{ext}}{2} - \frac{\text{max}\left(\{|\lambda_\text{int}|\}, \{|\lambda_\text{ext}|\}\right)}{2} \left(U_\text{ext} - U_\text{int}\right), \end{align*}

where "int" and "ext" stand for interior and exterior, and \(\{|\lambda|\}\) is the set of characteristic/signal speeds. The minus sign in front of the \(F_{\text{ext}}\) is necessary because the outward directed normal of the neighboring element has the opposite sign, i.e. \(n_i^{\text{ext}}=-n_i^{\text{int}}\). 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

\begin{align*} \lambda_{\pm}&=\pm\alpha-\beta^i n_i, \end{align*}

which correspond to the divergence cleaning field.

Note
  • In the strong form the dg_boundary_terms function returns \(G - F_\text{int}\)
  • 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 Typedef Documentation

◆ dg_package_data_primitive_tags

using grmhd::ValenciaDivClean::BoundaryCorrections::Rusanov::dg_package_data_primitive_tags
Initial value:
tmpl::list<hydro::Tags::RestMassDensity<DataVector>,
The electron fraction .
Definition Tags.hpp:105
The spatial velocity of the fluid, where . Here is the spatial part of the 4-velocity of the fluid,...
Definition Tags.hpp:257
The temperature of the fluid.
Definition Tags.hpp:292

◆ dg_package_data_temporary_tags

using grmhd::ValenciaDivClean::BoundaryCorrections::Rusanov::dg_package_data_temporary_tags
Initial value:
tmpl::list<
Definition Tags.hpp:65
Definition Tags.hpp:61
The spatial velocity one-form , where is raised and lowered with the spatial metric.
Definition Tags.hpp:265

◆ dg_package_data_volume_tags

using grmhd::ValenciaDivClean::BoundaryCorrections::Rusanov::dg_package_data_volume_tags
Initial value:
tmpl::list<hydro::Tags::GrmhdEquationOfState>

◆ dg_package_field_tags

using grmhd::ValenciaDivClean::BoundaryCorrections::Rusanov::dg_package_field_tags
Initial value:
The densitized divergence cleaning field associated with the magnetic field.
Definition Tags.hpp:93
Prefix indicating a boundary unit normal vector dotted into the flux.
Definition Prefixes.hpp:168
The densitized magnetic field .
Definition Tags.hpp:56
The densitized rest-mass density .
Definition Tags.hpp:32
The densitized momentum density .
Definition Tags.hpp:49
The densitized energy density .
Definition Tags.hpp:43
The densitized electron number density times the baryon mass .
Definition Tags.hpp:38

Member Function Documentation

◆ get_clone()

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

Member Data Documentation

◆ help

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:
  • src/Evolution/Systems/GrMhd/ValenciaDivClean/BoundaryCorrections/Rusanov.hpp