SpECTRE  v2026.04.01
Loading...
Searching...
No Matches
NewtonianEuler::BoundaryCorrections::Rusanov< Dim > Class Template 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 = tmpl::list<>
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
double dg_package_data (gsl::not_null< Scalar< DataVector > * > packaged_mass_density, gsl::not_null< tnsr::I< DataVector, Dim, Frame::Inertial > * > packaged_momentum_density, gsl::not_null< Scalar< DataVector > * > packaged_energy_density, gsl::not_null< Scalar< DataVector > * > packaged_normal_dot_flux_mass_density, gsl::not_null< tnsr::I< DataVector, Dim, Frame::Inertial > * > packaged_normal_dot_flux_momentum_density, gsl::not_null< Scalar< DataVector > * > packaged_normal_dot_flux_energy_density, gsl::not_null< Scalar< DataVector > * > packaged_abs_char_speed, const Scalar< DataVector > &mass_density, const tnsr::I< DataVector, Dim, Frame::Inertial > &momentum_density, const Scalar< DataVector > &energy_density, const tnsr::I< DataVector, Dim, Frame::Inertial > &flux_mass_density, const tnsr::IJ< DataVector, Dim, Frame::Inertial > &flux_momentum_density, const tnsr::I< DataVector, Dim, Frame::Inertial > &flux_energy_density, const tnsr::I< DataVector, Dim, Frame::Inertial > &velocity, const Scalar< DataVector > &specific_internal_energy, const tnsr::i< DataVector, Dim, Frame::Inertial > &normal_covector, const std::optional< tnsr::I< DataVector, Dim, Frame::Inertial > > &, const std::optional< Scalar< DataVector > > &normal_dot_mesh_velocity, const EquationsOfState::EquationOfState< false, 2 > &equation_of_state) const
void dg_boundary_terms (gsl::not_null< Scalar< DataVector > * > boundary_correction_mass_density, gsl::not_null< tnsr::I< DataVector, Dim, Frame::Inertial > * > boundary_correction_momentum_density, gsl::not_null< Scalar< DataVector > * > boundary_correction_energy_density, const Scalar< DataVector > &mass_density_int, const tnsr::I< DataVector, Dim, Frame::Inertial > &momentum_density_int, const Scalar< DataVector > &energy_density_int, const Scalar< DataVector > &normal_dot_flux_mass_density_int, const tnsr::I< DataVector, Dim, Frame::Inertial > &normal_dot_flux_momentum_density_int, const Scalar< DataVector > &normal_dot_flux_energy_density_int, const Scalar< DataVector > &abs_char_speed_int, const Scalar< DataVector > &mass_density_ext, const tnsr::I< DataVector, Dim, Frame::Inertial > &momentum_density_ext, const Scalar< DataVector > &energy_density_ext, const Scalar< DataVector > &normal_dot_flux_mass_density_ext, const tnsr::I< DataVector, Dim, Frame::Inertial > &normal_dot_flux_momentum_density_ext, const Scalar< DataVector > &normal_dot_flux_energy_density_ext, const Scalar< DataVector > &abs_char_speed_ext, dg::Formulation dg_formulation) const
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 Attributes

static constexpr Options::String help

Detailed Description

template<size_t Dim>
class NewtonianEuler::BoundaryCorrections::Rusanov< Dim >

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 by:

\begin{align*}\lambda_{\pm}&=v^i n_i \pm c_s, \\ \lambda_v&=v^i n_i \end{align*}

where \(v^i\) is the spatial velocity and \(c_s\) the sound speed.

Note
In the strong form the dg_boundary_terms function returns \(G - F_\text{int}\)

Member Typedef Documentation

◆ dg_package_data_primitive_tags

template<size_t Dim>
using NewtonianEuler::BoundaryCorrections::Rusanov< Dim >::dg_package_data_primitive_tags
Initial value:
tmpl::list<hydro::Tags::SpatialVelocity<DataVector, Dim>,
The specific internal energy .
Definition Tags.hpp:286

◆ dg_package_data_volume_tags

template<size_t Dim>
using NewtonianEuler::BoundaryCorrections::Rusanov< Dim >::dg_package_data_volume_tags
Initial value:
tmpl::list<hydro::Tags::EquationOfState<false, 2>>

◆ dg_package_field_tags

template<size_t Dim>
using NewtonianEuler::BoundaryCorrections::Rusanov< Dim >::dg_package_field_tags
Initial value:
tmpl::list<Tags::MassDensityCons, Tags::MomentumDensity<Dim>,
The energy density of the fluid.
Definition Tags.hpp:46
Prefix indicating a boundary unit normal vector dotted into the flux.
Definition Prefixes.hpp:168

Member Function Documentation

◆ get_clone()

template<size_t Dim>
std::unique_ptr< BoundaryCorrection > NewtonianEuler::BoundaryCorrections::Rusanov< Dim >::get_clone ( ) const
overridevirtual

Member Data Documentation

◆ help

template<size_t Dim>
Options::String NewtonianEuler::BoundaryCorrections::Rusanov< Dim >::help
staticconstexpr
Initial value:
= {
"Computes the Rusanov or local Lax-Friedrichs boundary correction term "
"for the Newtonian Euler/hydrodynamics system."}

The documentation for this class was generated from the following file:
  • src/Evolution/Systems/NewtonianEuler/BoundaryCorrections/Rusanov.hpp