SpECTRE
v2024.05.11
|
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::TildeE< Frame::Inertial, NeutrinoSpecies >..., Tags::TildeS< Frame::Inertial, NeutrinoSpecies >..., ::Tags::NormalDotFlux< Tags::TildeE< Frame::Inertial, NeutrinoSpecies > >..., ::Tags::NormalDotFlux< Tags::TildeS< Frame::Inertial, NeutrinoSpecies > >... > |
using | dg_package_data_temporary_tags = tmpl::list<> |
using | dg_package_data_primitive_tags = tmpl::list<> |
using | dg_package_data_volume_tags = tmpl::list<> |
using | dg_boundary_terms_volume_tags = tmpl::list<> |
Public Types inherited from RadiationTransport::M1Grey::BoundaryCorrections::BoundaryCorrection< tmpl::list< NeutrinoSpecies... > > | |
using | creatable_classes = tmpl::list< Rusanov< tmpl::list< NeutrinoSpecies... > > > |
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< tmpl::list< NeutrinoSpecies... > > > | get_clone () const override |
double | dg_package_data (const gsl::not_null< typename Tags::TildeE< Frame::Inertial, NeutrinoSpecies >::type * >... packaged_tilde_e, const gsl::not_null< typename Tags::TildeS< Frame::Inertial, NeutrinoSpecies >::type * >... packaged_tilde_s, const gsl::not_null< typename Tags::TildeE< Frame::Inertial, NeutrinoSpecies >::type * >... packaged_normal_dot_flux_tilde_e, const gsl::not_null< typename Tags::TildeS< Frame::Inertial, NeutrinoSpecies >::type * >... packaged_normal_dot_flux_tilde_s, const typename Tags::TildeE< Frame::Inertial, NeutrinoSpecies >::type &... tilde_e, const typename Tags::TildeS< Frame::Inertial, NeutrinoSpecies >::type &... tilde_s, const typename ::Tags::Flux< Tags::TildeE< Frame::Inertial, NeutrinoSpecies >, tmpl::size_t< 3 >, Frame::Inertial >::type &... flux_tilde_e, const typename ::Tags::Flux< Tags::TildeS< Frame::Inertial, NeutrinoSpecies >, tmpl::size_t< 3 >, Frame::Inertial >::type &... flux_tilde_s, 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 > > &mesh_velocity, const std::optional< Scalar< DataVector > > &normal_dot_mesh_velocity) const |
void | dg_boundary_terms (const gsl::not_null< typename Tags::TildeE< Frame::Inertial, NeutrinoSpecies >::type * >... boundary_correction_tilde_e, const gsl::not_null< typename Tags::TildeS< Frame::Inertial, NeutrinoSpecies >::type * >... boundary_correction_tilde_s, const typename Tags::TildeE< Frame::Inertial, NeutrinoSpecies >::type &... tilde_e_int, const typename Tags::TildeS< Frame::Inertial, NeutrinoSpecies >::type &... tilde_s_int, const typename Tags::TildeE< Frame::Inertial, NeutrinoSpecies >::type &... normal_dot_flux_tilde_e_int, const typename Tags::TildeS< Frame::Inertial, NeutrinoSpecies >::type &... normal_dot_flux_tilde_s_int, const typename Tags::TildeE< Frame::Inertial, NeutrinoSpecies >::type &... tilde_e_ext, const typename Tags::TildeS< Frame::Inertial, NeutrinoSpecies >::type &... tilde_s_ext, const typename Tags::TildeE< Frame::Inertial, NeutrinoSpecies >::type &... normal_dot_flux_tilde_e_ext, const typename Tags::TildeS< Frame::Inertial, NeutrinoSpecies >::type &... normal_dot_flux_tilde_s_ext, dg::Formulation dg_formulation) const |
Public Member Functions inherited from RadiationTransport::M1Grey::BoundaryCorrections::BoundaryCorrection< tmpl::list< NeutrinoSpecies... > > | |
BoundaryCorrection (const BoundaryCorrection &)=default | |
BoundaryCorrection (BoundaryCorrection &&)=default | |
BoundaryCorrection & | operator= (const BoundaryCorrection &)=default |
BoundaryCorrection & | operator= (BoundaryCorrection &&)=default |
virtual std::unique_ptr< BoundaryCorrection< tmpl::list< NeutrinoSpecies... > > > | get_clone () const=0 |
Static Public Attributes | |
static constexpr Options::String | help |
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}}\).
For radiation, the max characteristic/signal speed is 1 (the speed of light). Note that the characteristic speeds of this system are not yet fully implemented, see the warning in the documentation for the characteristics.
dg_boundary_terms
function returns \(G - F_\text{int}\)
|
inlineoverridevirtual |
|
staticconstexpr |