Path: blob/main/src/equations/numerical_fluxes.jl
2802 views
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).1# Since these FMAs can increase the performance of many numerical algorithms,2# we need to opt-in explicitly.3# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.4@muladd begin5#! format: noindent67# This file contains general numerical fluxes that are not specific to certain equations89"""10flux_central(u_ll, u_rr, orientation_or_normal_direction, equations::AbstractEquations)1112The classical central numerical flux `f((u_ll) + f(u_rr)) / 2`. When this flux is13used as volume flux, the discretization is equivalent to the classical weak form14DG method (except floating point errors).15"""16@inline function flux_central(u_ll, u_rr, orientation_or_normal_direction,17equations::AbstractEquations)18# Calculate regular 1D fluxes19f_ll = flux(u_ll, orientation_or_normal_direction, equations)20f_rr = flux(u_rr, orientation_or_normal_direction, equations)2122# Average regular fluxes23return 0.5f0 * (f_ll + f_rr)24end2526"""27FluxPlusDissipation(numerical_flux, dissipation)2829Combine a `numerical_flux` with a `dissipation` operator to create a new numerical flux.30"""31struct FluxPlusDissipation{NumericalFlux, Dissipation}32numerical_flux::NumericalFlux33dissipation::Dissipation34end3536@inline function (numflux::FluxPlusDissipation)(u_ll, u_rr,37orientation_or_normal_direction,38equations)39@unpack numerical_flux, dissipation = numflux4041return (numerical_flux(u_ll, u_rr, orientation_or_normal_direction, equations)42+43dissipation(u_ll, u_rr, orientation_or_normal_direction, equations))44end4546function Base.show(io::IO, f::FluxPlusDissipation)47print(io, "FluxPlusDissipation(", f.numerical_flux, ", ", f.dissipation, ")")48return nothing49end5051"""52FluxRotated(numerical_flux)5354Compute a `numerical_flux` flux in direction of a normal vector by rotating the solution,55computing the numerical flux in x-direction, and rotating the calculated flux back.5657Requires a rotationally invariant equation with equation-specific functions58[`rotate_to_x`](@ref) and [`rotate_from_x`](@ref).59"""60struct FluxRotated{NumericalFlux}61numerical_flux::NumericalFlux62end6364# Rotated surface flux computation (2D version)65@inline function (flux_rotated::FluxRotated)(u,66normal_direction::AbstractVector,67equations::AbstractEquations{2})68@unpack numerical_flux = flux_rotated6970norm_ = norm(normal_direction)71# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later72normal_vector = normal_direction / norm_7374u_rotated = rotate_to_x(u, normal_vector, equations)7576f = numerical_flux(u_rotated, 1, equations)7778return rotate_from_x(f, normal_vector, equations) * norm_79end8081# Rotated surface flux computation (2D version)82@inline function (flux_rotated::FluxRotated)(u_ll, u_rr,83normal_direction::AbstractVector,84equations::AbstractEquations{2})85@unpack numerical_flux = flux_rotated8687norm_ = norm(normal_direction)88# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later89normal_vector = normal_direction / norm_9091u_ll_rotated = rotate_to_x(u_ll, normal_vector, equations)92u_rr_rotated = rotate_to_x(u_rr, normal_vector, equations)9394f = numerical_flux(u_ll_rotated, u_rr_rotated, 1, equations)9596return rotate_from_x(f, normal_vector, equations) * norm_97end9899# Rotated surface flux computation (3D version)100@inline function (flux_rotated::FluxRotated)(u_ll, u_rr,101normal_direction::AbstractVector,102equations::AbstractEquations{3})103@unpack numerical_flux = flux_rotated104105# Storing these vectors could increase the performance by 20 percent106norm_ = norm(normal_direction)107# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later108normal_vector = normal_direction / norm_109110# Some vector that can't be identical to normal_vector (unless normal_vector == 0)111tangent1 = SVector(normal_direction[2], normal_direction[3], -normal_direction[1])112# Orthogonal projection113tangent1 -= dot(normal_vector, tangent1) * normal_vector114tangent1 = normalize(tangent1)115116# Third orthogonal vector117tangent2 = normalize(cross(normal_direction, tangent1))118119u_ll_rotated = rotate_to_x(u_ll, normal_vector, tangent1, tangent2, equations)120u_rr_rotated = rotate_to_x(u_rr, normal_vector, tangent1, tangent2, equations)121122f = numerical_flux(u_ll_rotated, u_rr_rotated, 1, equations)123124return rotate_from_x(f, normal_vector, tangent1, tangent2, equations) * norm_125end126127Base.show(io::IO, f::FluxRotated) = print(io, "FluxRotated(", f.numerical_flux, ")")128129"""130DissipationGlobalLaxFriedrichs(λ)131132Create a global Lax-Friedrichs dissipation operator with dissipation coefficient `λ`.133"""134struct DissipationGlobalLaxFriedrichs{RealT}135λ::RealT136end137138@inline function (dissipation::DissipationGlobalLaxFriedrichs)(u_ll, u_rr,139orientation::Integer,140equations)141@unpack λ = dissipation142return -λ / 2 * (u_rr - u_ll)143end144145@inline function (dissipation::DissipationGlobalLaxFriedrichs)(u_ll, u_rr,146normal_direction::AbstractVector,147equations)148@unpack λ = dissipation149return -λ / 2 * norm(normal_direction) * (u_rr - u_ll)150end151152function Base.show(io::IO, d::DissipationGlobalLaxFriedrichs)153print(io, "DissipationGlobalLaxFriedrichs(", d.λ, ")")154return nothing155end156157"""158DissipationLocalLaxFriedrichs(max_abs_speed=max_abs_speed)159160Create a local Lax-Friedrichs dissipation operator where the maximum absolute wave speed161is estimated as162`max_abs_speed(u_ll, u_rr, orientation_or_normal_direction, equations)`,163defaulting to [`max_abs_speed`](@ref).164"""165struct DissipationLocalLaxFriedrichs{MaxAbsSpeed}166max_abs_speed::MaxAbsSpeed167end168169DissipationLocalLaxFriedrichs() = DissipationLocalLaxFriedrichs(max_abs_speed)170171@inline function (dissipation::DissipationLocalLaxFriedrichs)(u_ll, u_rr,172orientation_or_normal_direction,173equations)174λ = dissipation.max_abs_speed(u_ll, u_rr, orientation_or_normal_direction,175equations)176return -0.5f0 * λ * (u_rr - u_ll)177end178179function Base.show(io::IO, d::DissipationLocalLaxFriedrichs)180print(io, "DissipationLocalLaxFriedrichs(", d.max_abs_speed, ")")181return nothing182end183184"""185max_abs_speed_naive(u_ll, u_rr, orientation::Integer, equations)186max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector, equations)187188Simple and fast estimate of the maximal wave speed of the Riemann problem with left and right states189`u_ll, u_rr`, based only on the local wave speeds associated to `u_ll` and `u_rr`.190191For non-integer arguments `normal_direction` in one dimension, `max_abs_speed_naive` returns192`abs(normal_direction[1]) * max_abs_speed_naive(u_ll, u_rr, 1, equations)`.193194Slightly more diffusive/overestimating than [`max_abs_speed`](@ref).195"""196function max_abs_speed_naive end197198# for non-integer `orientation_or_normal` arguments.199@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,200equations::AbstractEquations{1})201return abs(normal_direction[1]) * max_abs_speed_naive(u_ll, u_rr, 1, equations)202end203204"""205max_abs_speed(u_ll, u_rr, orientation::Integer, equations)206max_abs_speed(u_ll, u_rr, normal_direction::AbstractVector, equations)207208Simple and fast estimate of the maximal wave speed of the Riemann problem with left and right states209`u_ll, u_rr`, based only on the local wave speeds associated to `u_ll` and `u_rr`.210Less diffusive, i.e., overestimating than [`max_abs_speed_naive`](@ref).211212In particular, `max_abs_speed(u, u, i, equations)` gives the same result as `max_abs_speeds(u, equations)[i]`,213i.e., the wave speeds used in `max_dt` which computes the maximum stable time step size through the214[`StepsizeCallback`](@ref).215216For non-integer arguments `normal_direction` in one dimension, `max_abs_speed_naive` returns217`abs(normal_direction[1]) * max_abs_speed_naive(u_ll, u_rr, 1, equations)`.218219Defaults to [`min_max_speed_naive`](@ref) if no specialized version for the 'equations` at hand is available.220"""221@inline function max_abs_speed(u_ll, u_rr,222orientation_or_normal_direction,223equations::AbstractEquations)224# Use naive version as "backup" if no specialized version for the equations at hand is available225return max_abs_speed_naive(u_ll, u_rr, orientation_or_normal_direction, equations)226end227228const FluxLaxFriedrichs{MaxAbsSpeed} = FluxPlusDissipation{typeof(flux_central),229DissipationLocalLaxFriedrichs{MaxAbsSpeed}}230"""231FluxLaxFriedrichs(max_abs_speed=max_abs_speed)232233Local Lax-Friedrichs (Rusanov) flux with maximum wave speed estimate provided by234`max_abs_speed`, cf. [`DissipationLocalLaxFriedrichs`](@ref) and235[`max_abs_speed_naive`](@ref).236"""237function FluxLaxFriedrichs(max_abs_speed = max_abs_speed)238return FluxPlusDissipation(flux_central,239DissipationLocalLaxFriedrichs(max_abs_speed))240end241242function Base.show(io::IO, f::FluxLaxFriedrichs)243print(io, "FluxLaxFriedrichs(", f.dissipation.max_abs_speed, ")")244return nothing245end246247"""248flux_lax_friedrichs249250See [`FluxLaxFriedrichs`](@ref).251"""252const flux_lax_friedrichs = FluxLaxFriedrichs()253254@doc raw"""255DissipationLaxFriedrichsEntropyVariables(max_abs_speed=max_abs_speed)256257Create a local Lax-Friedrichs-type dissipation operator that is provably entropy stable. This operator258must be used together with an entropy-conservative two-point flux function (e.g., `flux_ec`) to yield259an entropy-stable surface flux. The surface flux function can be initialized as:260```julia261flux_es = FluxPlusDissipation(flux_ec, DissipationLaxFriedrichsEntropyVariables())262```263264In particular, the numerical flux has the form265```math266f^{\mathrm{ES}} = f^{\mathrm{EC}} - \frac{1}{2} \lambda_{\mathrm{max}} H (w_r - w_l),267```268where ``f^{\mathrm{EC}}`` is the entropy-conservative two-point flux function (computed with, e.g., `flux_ec`), ``\lambda_{\mathrm{max}}``269is the maximum wave speed estimated as `max_abs_speed(u_l, u_r, orientation_or_normal_direction, equations)`,270defaulting to [`max_abs_speed`](@ref), ``H`` is a symmetric positive-definite dissipation matrix that271depends on the left and right states `u_l` and `u_r`, and ``(w_r - w_l)`` is the jump in entropy variables.272Ideally, ``H (w_r - w_l) = (u_r - u_l)``, such that the dissipation operator is consistent with the local273Lax-Friedrichs dissipation.274275The entropy-stable dissipation operator is computed with the function276`function (dissipation::DissipationLaxFriedrichsEntropyVariables)(u_l, u_r, orientation_or_normal_direction, equations)`,277which must be specialized for each equation.278279For the derivation of the dissipation matrix for the multi-ion GLM-MHD equations, see:280- A. Rueda-Ramírez, A. Sikstel, G. Gassner, An Entropy-Stable Discontinuous Galerkin Discretization281of the Ideal Multi-Ion Magnetohydrodynamics System (2024). Journal of Computational Physics.282[DOI: 10.1016/j.jcp.2024.113655](https://doi.org/10.1016/j.jcp.2024.113655).283"""284struct DissipationLaxFriedrichsEntropyVariables{MaxAbsSpeed}285max_abs_speed::MaxAbsSpeed286end287288DissipationLaxFriedrichsEntropyVariables() = DissipationLaxFriedrichsEntropyVariables(max_abs_speed)289290function Base.show(io::IO, d::DissipationLaxFriedrichsEntropyVariables)291print(io, "DissipationLaxFriedrichsEntropyVariables(", d.max_abs_speed, ")")292return nothing293end294295@doc raw"""296DissipationMatrixWintersEtal()297298Creates the Roe-like entropy-stable matrix dissipation operator from Winters et al. (2017). This operator299must be used together with an entropy-conservative two-point flux function300(e.g., [`flux_ranocha`](@ref) or [`flux_chandrashekar`](@ref)) to yield301an entropy-stable surface flux. The surface flux function can be initialized as:302```julia303flux_es = FluxPlusDissipation(flux_ec, DissipationMatrixWintersEtal())304```305The 1D and 2D implementations are adapted from the [Atum.jl library](https://github.com/mwarusz/Atum.jl/blob/c7ed44f2b7972ac726ef345da7b98b0bda60e2a3/src/balancelaws/euler.jl#L198).306The 3D implementation is adapted from the [FLUXO library](https://github.com/project-fluxo/fluxo)307308For the derivation of the matrix dissipation operator, see:309- A. R. Winters, D. Derigs, G. Gassner, S. Walch, A uniquely defined entropy stable matrix dissipation operator310for high Mach number ideal MHD and compressible Euler simulations (2017). Journal of Computational Physics.311[DOI: 10.1016/j.jcp.2016.12.006](https://doi.org/10.1016/j.jcp.2016.12.006).312"""313struct DissipationMatrixWintersEtal end314315@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,316orientation::Integer,317equations::AbstractEquations{1})318return dissipation(u_ll, u_rr, SVector(1), equations)319end320321@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,322orientation::Integer,323equations::AbstractEquations{2})324if orientation == 1325return dissipation(u_ll, u_rr, SVector(1, 0), equations)326else # orientation == 2327return dissipation(u_ll, u_rr, SVector(0, 1), equations)328end329end330331@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,332orientation::Integer,333equations::AbstractEquations{3})334if orientation == 1335return dissipation(u_ll, u_rr, SVector(1, 0, 0), equations)336elseif orientation == 2337return dissipation(u_ll, u_rr, SVector(0, 1, 0), equations)338else # orientation == 3339return dissipation(u_ll, u_rr, SVector(0, 0, 1), equations)340end341end342343"""344FluxHLL(min_max_speed=min_max_speed_davis)345346Create an HLL (Harten, Lax, van Leer) numerical flux where the minimum and maximum347wave speeds are estimated as348`λ_min, λ_max = min_max_speed(u_ll, u_rr, orientation_or_normal_direction, equations)`,349defaulting to [`min_max_speed_davis`](@ref).350Original paper:351- Amiram Harten, Peter D. Lax, Bram van Leer (1983)352On Upstream Differencing and Godunov-Type Schemes for Hyperbolic Conservation Laws353[DOI: 10.1137/1025002](https://doi.org/10.1137/1025002)354"""355struct FluxHLL{MinMaxSpeed}356min_max_speed::MinMaxSpeed357end358359FluxHLL() = FluxHLL(min_max_speed_davis)360361"""362min_max_speed_naive(u_ll, u_rr, orientation::Integer, equations)363min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector, equations)364365Simple and fast estimate(!) of the minimal and maximal wave speed of the Riemann problem with366left and right states `u_ll, u_rr`, usually based only on the local wave speeds associated to367`u_ll` and `u_rr`.368Slightly more diffusive than [`min_max_speed_davis`](@ref).369- Amiram Harten, Peter D. Lax, Bram van Leer (1983)370On Upstream Differencing and Godunov-Type Schemes for Hyperbolic Conservation Laws371[DOI: 10.1137/1025002](https://doi.org/10.1137/1025002)372373See eq. (10.37) from374- Eleuterio F. Toro (2009)375Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction376[DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)377378See also [`FluxHLL`](@ref), [`min_max_speed_davis`](@ref), [`min_max_speed_einfeldt`](@ref).379"""380function min_max_speed_naive end381382"""383min_max_speed_davis(u_ll, u_rr, orientation::Integer, equations)384min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector, equations)385386Simple and fast estimates of the minimal and maximal wave speed of the Riemann problem with387left and right states `u_ll, u_rr`, usually based only on the local wave speeds associated to388`u_ll` and `u_rr`.389390- S.F. Davis (1988)391Simplified Second-Order Godunov-Type Methods392[DOI: 10.1137/0909030](https://doi.org/10.1137/0909030)393394See eq. (10.38) from395- Eleuterio F. Toro (2009)396Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction397[DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)398See also [`FluxHLL`](@ref), [`min_max_speed_naive`](@ref), [`min_max_speed_einfeldt`](@ref).399"""400function min_max_speed_davis end401402"""403min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer, equations)404min_max_speed_einfeldt(u_ll, u_rr, normal_direction::AbstractVector, equations)405406More advanced mininmal and maximal wave speed computation based on407- Bernd Einfeldt (1988)408On Godunov-type methods for gas dynamics.409[DOI: 10.1137/0725021](https://doi.org/10.1137/0725021)410- Bernd Einfeldt, Claus-Dieter Munz, Philip L. Roe and Björn Sjögreen (1991)411On Godunov-type methods near low densities.412[DOI: 10.1016/0021-9991(91)90211-3](https://doi.org/10.1016/0021-9991(91)90211-3)413414originally developed for the compressible Euler equations.415A compact representation can be found in [this lecture notes, eq. (9.28)](https://metaphor.ethz.ch/x/2019/hs/401-4671-00L/literature/mishra_hyperbolic_pdes.pdf).416417See also [`FluxHLL`](@ref), [`min_max_speed_naive`](@ref), [`min_max_speed_davis`](@ref).418"""419function min_max_speed_einfeldt end420421@inline function (numflux::FluxHLL)(u_ll, u_rr, orientation_or_normal_direction,422equations)423λ_min, λ_max = numflux.min_max_speed(u_ll, u_rr, orientation_or_normal_direction,424equations)425426if λ_min >= 0 && λ_max >= 0427return flux(u_ll, orientation_or_normal_direction, equations)428elseif λ_max <= 0 && λ_min <= 0429return flux(u_rr, orientation_or_normal_direction, equations)430else431f_ll = flux(u_ll, orientation_or_normal_direction, equations)432f_rr = flux(u_rr, orientation_or_normal_direction, equations)433inv_λ_max_minus_λ_min = inv(λ_max - λ_min)434factor_ll = λ_max * inv_λ_max_minus_λ_min435factor_rr = λ_min * inv_λ_max_minus_λ_min436factor_diss = λ_min * λ_max * inv_λ_max_minus_λ_min437return factor_ll * f_ll - factor_rr * f_rr + factor_diss * (u_rr - u_ll)438end439end440441Base.show(io::IO, numflux::FluxHLL) = print(io, "FluxHLL(", numflux.min_max_speed, ")")442443"""444flux_hll445446See [`FluxHLL`](@ref).447"""448const flux_hll = FluxHLL()449450"""451flux_hlle452453See [`min_max_speed_einfeldt`](@ref).454This is a [`FluxHLL`](@ref)-type two-wave solver with special estimates of the wave speeds.455"""456const flux_hlle = FluxHLL(min_max_speed_einfeldt)457458"""459flux_shima_etal_turbo(u_ll, u_rr, orientation_or_normal_direction, equations)460461Equivalent to [`flux_shima_etal`](@ref) except that it may use specialized462methods, e.g., when used with [`VolumeIntegralFluxDifferencing`](@ref).463These specialized methods may enable better use of SIMD instructions to464increase runtime efficiency on modern hardware.465"""466@inline function flux_shima_etal_turbo(u_ll, u_rr, orientation_or_normal_direction,467equations)468return flux_shima_etal(u_ll, u_rr, orientation_or_normal_direction, equations)469end470471"""472flux_ranocha_turbo(u_ll, u_rr, orientation_or_normal_direction, equations)473474Equivalent to [`flux_ranocha`](@ref) except that it may use specialized475methods, e.g., when used with [`VolumeIntegralFluxDifferencing`](@ref).476These specialized methods may enable better use of SIMD instructions to477increase runtime efficiency on modern hardware.478"""479@inline function flux_ranocha_turbo(u_ll, u_rr, orientation_or_normal_direction,480equations)481return flux_ranocha(u_ll, u_rr, orientation_or_normal_direction, equations)482end483484"""485flux_wintermeyer_etal(u_ll, u_rr, orientation_or_normal_direction, equations)486487Total energy conservative (mathematical entropy for shallow water equations) split form.488When the bottom topography is nonzero this scheme will be well-balanced when used as a `volume_flux`.489For the `surface_flux` either [`flux_wintermeyer_etal`](@ref) or [`flux_fjordholm_etal`](@ref) can490be used to ensure well-balancedness and entropy conservation.491492!!! note493This function is defined in Trixi.jl to have a common interface for the494methods implemented in the subpackages [TrixiAtmo.jl](https://github.com/trixi-framework/TrixiAtmo.jl)495and [TrixiShallowWater.jl](https://github.com/trixi-framework/TrixiShallowWater.jl).496"""497function flux_wintermeyer_etal end498499"""500flux_nonconservative_wintermeyer_etal(u_ll, u_rr, orientation_or_normal_direction, equations)501502Non-symmetric two-point volume flux discretizing the nonconservative (source) term503that contains the gradient of the bottom topography for the shallow water equations.504505Gives entropy conservation and well-balancedness on both the volume and surface when combined with506[`flux_wintermeyer_etal`](@ref).507508!!! note509This function is defined in Trixi.jl to have a common interface for the510methods implemented in the subpackages [TrixiAtmo.jl](https://github.com/trixi-framework/TrixiAtmo.jl)511and [TrixiShallowWater.jl](https://github.com/trixi-framework/TrixiShallowWater.jl).512"""513function flux_nonconservative_wintermeyer_etal end514515"""516flux_fjordholm_etal(u_ll, u_rr, orientation_or_normal_direction, equations)517518Total energy conservative (mathematical entropy for shallow water equations). When the bottom topography519is nonzero this should only be used as a surface flux otherwise the scheme will not be well-balanced.520For well-balancedness in the volume flux use [`flux_wintermeyer_etal`](@ref).521522!!! note523This function is defined in Trixi.jl to have a common interface for the524methods implemented in the subpackages [TrixiAtmo.jl](https://github.com/trixi-framework/TrixiAtmo.jl)525and [TrixiShallowWater.jl](https://github.com/trixi-framework/TrixiShallowWater.jl).526"""527function flux_fjordholm_etal end528529"""530flux_nonconservative_fjordholm_etal(u_ll, u_rr, orientation_or_normal_direction, equations)531532Non-symmetric two-point surface flux discretizing the nonconservative (source) term of533that contains the gradient of the bottom topography for the shallow water equations.534535This flux can be used together with [`flux_fjordholm_etal`](@ref) at interfaces to ensure entropy536conservation and well-balancedness.537538!!! note539This function is defined in Trixi.jl to have a common interface for the540methods implemented in the subpackages [TrixiAtmo.jl](https://github.com/trixi-framework/TrixiAtmo.jl)541and [TrixiShallowWater.jl](https://github.com/trixi-framework/TrixiShallowWater.jl).542"""543function flux_nonconservative_fjordholm_etal end544545"""546FluxUpwind(splitting)547548A numerical flux `f(u_left, u_right) = f⁺(u_left) + f⁻(u_right)` based on549flux vector splitting.550551The [`SurfaceIntegralUpwind`](@ref) with a given `splitting` is equivalent to552the [`SurfaceIntegralStrongForm`](@ref) with `FluxUpwind(splitting)`553as numerical flux (up to floating point differences). Note, that554[`SurfaceIntegralUpwind`](@ref) is only available on [`TreeMesh`](@ref).555556!!! warning "Experimental implementation (upwind SBP)"557This is an experimental feature and may change in future releases.558"""559struct FluxUpwind{Splitting}560splitting::Splitting561end562563@inline function (numflux::FluxUpwind)(u_ll, u_rr, orientation::Int, equations)564@unpack splitting = numflux565fm = splitting(u_rr, Val{:minus}(), orientation, equations)566fp = splitting(u_ll, Val{:plus}(), orientation, equations)567return fm + fp568end569570@inline function (numflux::FluxUpwind)(u_ll, u_rr,571normal_direction::AbstractVector,572equations::AbstractEquations{2})573@unpack splitting = numflux574f_tilde_m = splitting(u_rr, Val{:minus}(), normal_direction, equations)575f_tilde_p = splitting(u_ll, Val{:plus}(), normal_direction, equations)576return f_tilde_m + f_tilde_p577end578579Base.show(io::IO, f::FluxUpwind) = print(io, "FluxUpwind(", f.splitting, ")")580end # @muladd581582583