SpECTRE
v2023.01.13
|
Compute the local Lax-Friedrichs numerical flux. More...
#include <LocalLaxFriedrichs.hpp>
Classes | |
struct | MaxAbsCharSpeed |
The maximum characteristic speed modulus on one side of the interface. More... | |
Public Types | |
using | variables_tags = typename System::variables_tag::tags_list |
using | package_field_tags = tmpl::push_back< tmpl::append< db::wrap_tags_in<::Tags::NormalDotFlux, variables_tags >, variables_tags >, MaxAbsCharSpeed > |
using | package_extra_tags = tmpl::list<> |
using | argument_tags = tmpl::push_back< tmpl::append< db::wrap_tags_in<::Tags::NormalDotFlux, variables_tags >, variables_tags >, char_speeds_tag > |
using | options = tmpl::list<> |
Public Member Functions | |
void | pup (PUP::er &) |
template<class... Args> | |
void | package_data (const Args &... args) const |
template<class... Args> | |
void | operator() (const Args &... args) const |
Static Public Attributes | |
static constexpr Options::String | help |
Compute the local Lax-Friedrichs numerical flux.
Let \(U\) be the state vector of the system and \(F^i\) the corresponding volume fluxes. Let \(n_i\) be the unit normal to the interface. Denoting \(F := n_i F^i\), the local Lax-Friedrichs flux is
\begin{align*} G_\text{LLF} = \frac{F_\text{int} + F_\text{ext}}{2} + \frac{\alpha}{2}\left(U_\text{int} - U_\text{ext}\right), \end{align*}
where "int" and "ext" stand for interior and exterior, respectively, and
\begin{align*} \alpha = \text{max}\left(\{|\lambda_\text{int}|\}, \{|\lambda_\text{ext}|\}\right), \end{align*}
where \(\{|\lambda|\}\) is the set of all moduli of the characteristic speeds along a given normal.
|
staticconstexpr |