|
SpECTRE
v2025.08.19
|
A Rusanov/local Lax-Friedrichs Riemann solver. More...
#include <Rusanov.hpp>
Public Types | |
| using | options = tmpl::list<> |
| using | dg_package_field_tags = tmpl::list< Tags::U, ::Tags::NormalDotFlux< Tags::U >, AbsCharSpeed > |
| using | dg_package_data_temporary_tags = tmpl::list<> |
| using | dg_package_data_volume_tags = tmpl::list<> |
| using | dg_boundary_terms_volume_tags = tmpl::list<> |
Public Member Functions | |
| Rusanov (const Rusanov &)=default | |
| Rusanov & | operator= (const Rusanov &)=default |
| Rusanov (Rusanov &&)=default | |
| Rusanov & | operator= (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 | |
| BoundaryCorrection & | operator= (const BoundaryCorrection &)=default |
| BoundaryCorrection (BoundaryCorrection &&)=default | |
| BoundaryCorrection & | operator= (BoundaryCorrection &&)=default |
| BoundaryCorrection (CkMigrateMessage *msg) | |
| virtual std::unique_ptr< BoundaryCorrection > | get_clone () const =0 |
Static Public Member Functions | |
| static double | dg_package_data (gsl::not_null< Scalar< DataVector > * > packaged_u, gsl::not_null< Scalar< DataVector > * > packaged_normal_dot_flux, gsl::not_null< Scalar< DataVector > * > packaged_abs_char_speed, const Scalar< DataVector > &u, const tnsr::I< DataVector, 1, Frame::Inertial > &flux_u, const tnsr::i< DataVector, 1, Frame::Inertial > &normal_covector, const std::optional< tnsr::I< DataVector, 1, Frame::Inertial > > &mesh_velocity, const std::optional< Scalar< DataVector > > &normal_dot_mesh_velocity) |
| static void | dg_boundary_terms (gsl::not_null< Scalar< DataVector > * > boundary_correction_u, const Scalar< DataVector > &u_int, const Scalar< DataVector > &normal_dot_flux_u_int, const Scalar< DataVector > &abs_char_speed_int, const Scalar< DataVector > &u_ext, const Scalar< DataVector > &normal_dot_flux_u_ext, const Scalar< DataVector > &abs_char_speed_ext, dg::Formulation dg_formulation) |
Static Public Attributes | |
| static constexpr Options::String | help |
A Rusanov/local Lax-Friedrichs Riemann solver.
Let \(U\) be the evolved variable, \(F^i\) the flux, 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(|U_\text{int}|, |U_\text{ext}|\right)}{2} \left(U_\text{ext} - U_\text{int}\right), \end{align*}
where "int" and "ext" stand for interior and exterior. 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}}\).
dg_boundary_terms function returns \(G - F_\text{int}\)
|
overridevirtual |
Implements evolution::BoundaryCorrection.
|
staticconstexpr |