Path: blob/main/src/equations/compressible_euler_1d.jl
5586 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@doc raw"""8CompressibleEulerEquations1D(gamma)910The compressible Euler equations11```math12\frac{\partial}{\partial t}13\begin{pmatrix}14\rho \\ \rho v_1 \\ \rho e_{\text{total}}15\end{pmatrix}16+17\frac{\partial}{\partial x}18\begin{pmatrix}19\rho v_1 \\ \rho v_1^2 + p \\ (\rho e_{\text{total}} +p) v_120\end{pmatrix}21=22\begin{pmatrix}230 \\ 0 \\ 024\end{pmatrix}25```26for an ideal gas with ratio of specific heats `gamma` in one space dimension.27Here, ``\rho`` is the density, ``v_1`` the velocity, ``e_{\text{total}}`` the specific total energy, and28```math29p = (\gamma - 1) \left( \rho e_{\text{total}} - \frac{1}{2} \rho v_1^2 \right)30```31the pressure.32"""33struct CompressibleEulerEquations1D{RealT <: Real} <:34AbstractCompressibleEulerEquations{1, 3}35gamma::RealT # ratio of specific heats36inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications3738function CompressibleEulerEquations1D(gamma)39γ, inv_gamma_minus_one = promote(gamma, inv(gamma - 1))40return new{typeof(γ)}(γ, inv_gamma_minus_one)41end42end4344function varnames(::typeof(cons2cons), ::CompressibleEulerEquations1D)45return ("rho", "rho_v1", "rho_e_total")46end47varnames(::typeof(cons2prim), ::CompressibleEulerEquations1D) = ("rho", "v1", "p")4849"""50initial_condition_constant(x, t, equations::CompressibleEulerEquations1D)5152A constant initial condition to test free-stream preservation.53"""54function initial_condition_constant(x, t, equations::CompressibleEulerEquations1D)55RealT = eltype(x)56rho = 157rho_v1 = convert(RealT, 0.1)58rho_e_total = 1059return SVector(rho, rho_v1, rho_e_total)60end6162"""63initial_condition_convergence_test(x, t, equations::CompressibleEulerEquations1D)6465A smooth initial condition used for convergence tests in combination with66[`source_terms_convergence_test`](@ref)67(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).68"""69function initial_condition_convergence_test(x, t,70equations::CompressibleEulerEquations1D)71RealT = eltype(x)72c = 273A = convert(RealT, 0.1)74L = 275f = 1.0f0 / L76ω = 2 * convert(RealT, pi) * f77ini = c + A * sin(ω * (x[1] - t))7879rho = ini80rho_v1 = ini81rho_e_total = ini^28283return SVector(rho, rho_v1, rho_e_total)84end8586"""87source_terms_convergence_test(u, x, t, equations::CompressibleEulerEquations1D)8889Source terms used for convergence tests in combination with90[`initial_condition_convergence_test`](@ref)91(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).9293References for the method of manufactured solutions (MMS):94- Kambiz Salari and Patrick Knupp (2000)95Code Verification by the Method of Manufactured Solutions96[DOI: 10.2172/759450](https://doi.org/10.2172/759450)97- Patrick J. Roache (2002)98Code Verification by the Method of Manufactured Solutions99[DOI: 10.1115/1.1436090](https://doi.org/10.1115/1.1436090)100"""101@inline function source_terms_convergence_test(u, x, t,102equations::CompressibleEulerEquations1D)103# Same settings as in `initial_condition`104RealT = eltype(u)105c = 2106A = convert(RealT, 0.1)107L = 2108f = 1.0f0 / L109ω = 2 * convert(RealT, pi) * f110γ = equations.gamma111112x1, = x113114si, co = sincos(ω * (x1 - t))115rho = c + A * si116rho_x = ω * A * co117118# Note that d/dt rho = -d/dx rho.119# This yields du2 = du3 = d/dx p (derivative of pressure).120# Other terms vanish because of v = 1.121du1 = 0122du2 = rho_x * (2 * rho - 0.5f0) * (γ - 1)123du3 = du2124125return SVector(du1, du2, du3)126end127128"""129initial_condition_density_wave(x, t, equations::CompressibleEulerEquations1D)130131A sine wave in the density with constant velocity and pressure; reduces the132compressible Euler equations to the linear advection equations.133This setup is the test case for stability of EC fluxes from paper134- Gregor J. Gassner, Magnus Svärd, Florian J. Hindenlang (2020)135Stability issues of entropy-stable and/or split-form high-order schemes136[arXiv: 2007.09026](https://arxiv.org/abs/2007.09026)137with the following parameters138- domain [-1, 1]139- mesh = 4x4140- polydeg = 5141"""142function initial_condition_density_wave(x, t, equations::CompressibleEulerEquations1D)143RealT = eltype(x)144v1 = convert(RealT, 0.1)145rho = 1 + convert(RealT, 0.98) * sinpi(2 * (x[1] - t * v1))146rho_v1 = rho * v1147p = 20148rho_e_total = p / (equations.gamma - 1) + 0.5f0 * rho * v1^2149return SVector(rho, rho_v1, rho_e_total)150end151152"""153initial_condition_weak_blast_wave(x, t, equations::CompressibleEulerEquations1D)154155A weak blast wave taken from156- Sebastian Hennemann, Gregor J. Gassner (2020)157A provably entropy stable subcell shock capturing approach for high order split form DG158[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)159"""160function initial_condition_weak_blast_wave(x, t,161equations::CompressibleEulerEquations1D)162# From Hennemann & Gassner JCP paper 2020 (Sec. 6.3)163# Set up polar coordinates164RealT = eltype(x)165inicenter = SVector(0)166x_norm = x[1] - inicenter[1]167r = abs(x_norm)168# The following code is equivalent to169# phi = atan(0.0, x_norm)170# cos_phi = cos(phi)171# in 1D but faster172cos_phi = x_norm > 0 ? 1 : -1173174# Calculate primitive variables175rho = r > 0.5f0 ? one(RealT) : convert(RealT, 1.1691)176v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos_phi177p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)178179return prim2cons(SVector(rho, v1, p), equations)180end181182"""183initial_condition_eoc_test_coupled_euler_gravity(x, t, equations::CompressibleEulerEquations1D)184185One dimensional variant of the setup used for convergence tests of the Euler equations186with self-gravity from187- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)188A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics189[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)190!!! note191There is no additional source term necessary for the manufactured solution in one192spatial dimension. Thus, [`source_terms_eoc_test_coupled_euler_gravity`](@ref) is not193present there.194"""195function initial_condition_eoc_test_coupled_euler_gravity(x, t,196equations::CompressibleEulerEquations1D)197# OBS! this assumes that γ = 2 other manufactured source terms are incorrect198if equations.gamma != 2199error("adiabatic constant must be 2 for the coupling convergence test")200end201RealT = eltype(x)202c = 2203A = convert(RealT, 0.1)204ini = c + A * sinpi(x[1] - t)205G = 1 # gravitational constant206207rho = ini208v1 = 1209p = 2 * ini^2 * G / convert(RealT, pi) # * 2 / ndims, but ndims==1 here210211return prim2cons(SVector(rho, v1, p), equations)212end213214"""215boundary_condition_slip_wall(u_inner, orientation, direction, x, t,216surface_flux_function, equations::CompressibleEulerEquations1D)217Determine the boundary numerical surface flux for a slip wall condition.218Imposes a zero normal velocity at the wall.219Density is taken from the internal solution state and pressure is computed as an220exact solution of a 1D Riemann problem. Further details about this boundary state221are available in the paper:222- J. J. W. van der Vegt and H. van der Ven (2002)223Slip flow boundary conditions in discontinuous Galerkin discretizations of224the Euler equations of gas dynamics225[PDF](https://reports.nlr.nl/bitstream/handle/10921/692/TP-2002-300.pdf?sequence=1)226227Should be used together with [`TreeMesh`](@ref).228"""229@inline function boundary_condition_slip_wall(u_inner, orientation,230direction, x, t,231surface_flux_function,232equations::CompressibleEulerEquations1D)233# compute the primitive variables234rho_local, v_normal, p_local = cons2prim(u_inner, equations)235236if isodd(direction) # flip sign of normal to make it outward pointing237v_normal *= -1238end239240# Get the solution of the pressure Riemann problem241# See Section 6.3.3 of242# Eleuterio F. Toro (2009)243# Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction244# [DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)245if v_normal <= 0246sound_speed = sqrt(equations.gamma * p_local / rho_local) # local sound speed247p_star = p_local *248(1 + 0.5f0 * (equations.gamma - 1) * v_normal / sound_speed)^(2 *249equations.gamma *250equations.inv_gamma_minus_one)251else # v_normal > 0252A = 2 / ((equations.gamma + 1) * rho_local)253B = p_local * (equations.gamma - 1) / (equations.gamma + 1)254p_star = p_local +2550.5f0 * v_normal / A *256(v_normal + sqrt(v_normal^2 + 4 * A * (p_local + B)))257end258259# For the slip wall we directly set the flux as the normal velocity is zero260return SVector(0, p_star, 0)261end262263# Calculate 1D flux for a single point264@inline function flux(u, orientation::Integer, equations::CompressibleEulerEquations1D)265rho, rho_v1, rho_e_total = u266v1 = rho_v1 / rho267p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)268# Ignore orientation since it is always "1" in 1D269f1 = rho_v1270f2 = rho_v1 * v1 + p271f3 = (rho_e_total + p) * v1272return SVector(f1, f2, f3)273end274275"""276flux_shima_etal(u_ll, u_rr, orientation, equations::CompressibleEulerEquations1D)277278This flux is is a modification of the original kinetic energy preserving two-point flux by279- Yuichi Kuya, Kosuke Totani and Soshi Kawai (2018)280Kinetic energy and entropy preserving schemes for compressible flows281by split convective forms282[DOI: 10.1016/j.jcp.2018.08.058](https://doi.org/10.1016/j.jcp.2018.08.058)283284The modification is in the energy flux to guarantee pressure equilibrium and was developed by285- Nao Shima, Yuichi Kuya, Yoshiharu Tamaki, Soshi Kawai (JCP 2020)286Preventing spurious pressure oscillations in split convective form discretizations for287compressible flows288[DOI: 10.1016/j.jcp.2020.110060](https://doi.org/10.1016/j.jcp.2020.110060)289"""290@inline function flux_shima_etal(u_ll, u_rr, orientation::Integer,291equations::CompressibleEulerEquations1D)292# Unpack left and right state293rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)294rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)295296# Average each factor of products in flux297rho_avg = 0.5f0 * (rho_ll + rho_rr)298v1_avg = 0.5f0 * (v1_ll + v1_rr)299p_avg = 0.5f0 * (p_ll + p_rr)300kin_avg = 0.5f0 * (v1_ll * v1_rr)301302# Calculate fluxes303# Ignore orientation since it is always "1" in 1D304pv1_avg = 0.5f0 * (p_ll * v1_rr + p_rr * v1_ll)305f1 = rho_avg * v1_avg306f2 = f1 * v1_avg + p_avg307f3 = p_avg * v1_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv1_avg308309return SVector(f1, f2, f3)310end311312"""313flux_kennedy_gruber(u_ll, u_rr, orientation, equations::CompressibleEulerEquations1D)314315Kinetic energy preserving two-point flux by316- Kennedy and Gruber (2008)317Reduced aliasing formulations of the convective terms within the318Navier-Stokes equations for a compressible fluid319[DOI: 10.1016/j.jcp.2007.09.020](https://doi.org/10.1016/j.jcp.2007.09.020)320"""321@inline function flux_kennedy_gruber(u_ll, u_rr, orientation::Integer,322equations::CompressibleEulerEquations1D)323# Unpack left and right state324rho_e_total_ll = last(u_ll)325rho_e_total_rr = last(u_rr)326rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)327rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)328329# Average each factor of products in flux330rho_avg = 0.5f0 * (rho_ll + rho_rr)331v1_avg = 0.5f0 * (v1_ll + v1_rr)332p_avg = 0.5f0 * (p_ll + p_rr)333e_avg = 0.5f0 * (rho_e_total_ll / rho_ll + rho_e_total_rr / rho_rr)334335# Ignore orientation since it is always "1" in 1D336f1 = rho_avg * v1_avg337f2 = rho_avg * v1_avg * v1_avg + p_avg338f3 = (rho_avg * e_avg + p_avg) * v1_avg339340return SVector(f1, f2, f3)341end342343"""344flux_chandrashekar(u_ll, u_rr, orientation, equations::CompressibleEulerEquations1D)345346Entropy conserving two-point flux by347- Chandrashekar (2013)348Kinetic Energy Preserving and Entropy Stable Finite Volume Schemes349for Compressible Euler and Navier-Stokes Equations350[DOI: 10.4208/cicp.170712.010313a](https://doi.org/10.4208/cicp.170712.010313a)351"""352@inline function flux_chandrashekar(u_ll, u_rr, orientation::Integer,353equations::CompressibleEulerEquations1D)354# Unpack left and right state355rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)356rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)357beta_ll = 0.5f0 * rho_ll / p_ll358beta_rr = 0.5f0 * rho_rr / p_rr359specific_kin_ll = 0.5f0 * (v1_ll^2)360specific_kin_rr = 0.5f0 * (v1_rr^2)361362# Compute the necessary mean values363rho_avg = 0.5f0 * (rho_ll + rho_rr)364rho_mean = ln_mean(rho_ll, rho_rr)365beta_mean = ln_mean(beta_ll, beta_rr)366beta_avg = 0.5f0 * (beta_ll + beta_rr)367v1_avg = 0.5f0 * (v1_ll + v1_rr)368p_mean = 0.5f0 * rho_avg / beta_avg369velocity_square_avg = specific_kin_ll + specific_kin_rr370371# Calculate fluxes372# Ignore orientation since it is always "1" in 1D373f1 = rho_mean * v1_avg374f2 = f1 * v1_avg + p_mean375f3 = f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +376f2 * v1_avg377378return SVector(f1, f2, f3)379end380381"""382flux_ranocha(u_ll, u_rr, orientation_or_normal_direction, equations::CompressibleEulerEquations1D)383384Entropy conserving and kinetic energy preserving two-point flux by385- Hendrik Ranocha (2018)386Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods387for Hyperbolic Balance Laws388[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)389See also390- Hendrik Ranocha (2020)391Entropy Conserving and Kinetic Energy Preserving Numerical Methods for392the Euler Equations Using Summation-by-Parts Operators393[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)394"""395@inline function flux_ranocha(u_ll, u_rr, orientation::Integer,396equations::CompressibleEulerEquations1D)397# Unpack left and right state398rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)399rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)400401# Compute the necessary mean values402rho_mean = ln_mean(rho_ll, rho_rr)403# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`404# in exact arithmetic since405# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)406# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)407inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)408v1_avg = 0.5f0 * (v1_ll + v1_rr)409p_avg = 0.5f0 * (p_ll + p_rr)410velocity_square_avg = 0.5f0 * (v1_ll * v1_rr)411412# Calculate fluxes413# Ignore orientation since it is always "1" in 1D414f1 = rho_mean * v1_avg415f2 = f1 * v1_avg + p_avg416f3 = f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +4170.5f0 * (p_ll * v1_rr + p_rr * v1_ll)418419return SVector(f1, f2, f3)420end421422# While `normal_direction` isn't strictly necessary in 1D, certain solvers assume that423# the normal component is incorporated into the numerical flux.424#425# See `flux(u, normal_direction::AbstractVector, equations::AbstractEquations{1})` for a426# similar implementation.427@inline function flux_ranocha(u_ll, u_rr, normal_direction::AbstractVector,428equations::CompressibleEulerEquations1D)429return normal_direction[1] * flux_ranocha(u_ll, u_rr, 1, equations)430end431432"""433splitting_steger_warming(u, orientation::Integer,434equations::CompressibleEulerEquations1D)435splitting_steger_warming(u, which::Union{Val{:minus}, Val{:plus}}436orientation::Integer,437equations::CompressibleEulerEquations1D)438439Splitting of the compressible Euler flux of Steger and Warming.440441Returns a tuple of the fluxes "minus" (associated with waves going into the442negative axis direction) and "plus" (associated with waves going into the443positive axis direction). If only one of the fluxes is required, use the444function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.445446!!! warning "Experimental implementation (upwind SBP)"447This is an experimental feature and may change in future releases.448449## References450451- Joseph L. Steger and R. F. Warming (1979)452Flux Vector Splitting of the Inviscid Gasdynamic Equations453With Application to Finite Difference Methods454[NASA Technical Memorandum](https://ntrs.nasa.gov/api/citations/19790020779/downloads/19790020779.pdf)455"""456@inline function splitting_steger_warming(u, orientation::Integer,457equations::CompressibleEulerEquations1D)458fm = splitting_steger_warming(u, Val{:minus}(), orientation, equations)459fp = splitting_steger_warming(u, Val{:plus}(), orientation, equations)460return fm, fp461end462463@inline function splitting_steger_warming(u, ::Val{:plus}, orientation::Integer,464equations::CompressibleEulerEquations1D)465rho, rho_v1, rho_e_total = u466v1 = rho_v1 / rho467p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)468a = sqrt(equations.gamma * p / rho)469470lambda1 = v1471lambda2 = v1 + a472lambda3 = v1 - a473474lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)475lambda2_p = positive_part(lambda2)476lambda3_p = positive_part(lambda3)477478alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p479480rho_2gamma = 0.5f0 * rho / equations.gamma481f1p = rho_2gamma * alpha_p482f2p = rho_2gamma * (alpha_p * v1 + a * (lambda2_p - lambda3_p))483f3p = rho_2gamma * (alpha_p * 0.5f0 * v1^2 + a * v1 * (lambda2_p - lambda3_p)484+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)485486return SVector(f1p, f2p, f3p)487end488489@inline function splitting_steger_warming(u, ::Val{:minus}, orientation::Integer,490equations::CompressibleEulerEquations1D)491rho, rho_v1, rho_e_total = u492v1 = rho_v1 / rho493p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)494a = sqrt(equations.gamma * p / rho)495496lambda1 = v1497lambda2 = v1 + a498lambda3 = v1 - a499500lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)501lambda2_m = negative_part(lambda2)502lambda3_m = negative_part(lambda3)503504alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m505506rho_2gamma = 0.5f0 * rho / equations.gamma507f1m = rho_2gamma * alpha_m508f2m = rho_2gamma * (alpha_m * v1 + a * (lambda2_m - lambda3_m))509f3m = rho_2gamma * (alpha_m * 0.5f0 * v1^2 + a * v1 * (lambda2_m - lambda3_m)510+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)511512return SVector(f1m, f2m, f3m)513end514515"""516splitting_vanleer_haenel(u, orientation::Integer,517equations::CompressibleEulerEquations1D)518splitting_vanleer_haenel(u, which::Union{Val{:minus}, Val{:plus}}519orientation::Integer,520equations::CompressibleEulerEquations1D)521522Splitting of the compressible Euler flux from van Leer. This splitting further523contains a reformulation due to Hänel et al. where the energy flux uses the524enthalpy. The pressure splitting is independent from the splitting of the525convective terms. As such there are many pressure splittings suggested across526the literature. We implement the 'p4' variant suggested by Liou and Steffen as527it proved the most robust in practice.528529Returns a tuple of the fluxes "minus" (associated with waves going into the530negative axis direction) and "plus" (associated with waves going into the531positive axis direction). If only one of the fluxes is required, use the532function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.533534!!! warning "Experimental implementation (upwind SBP)"535This is an experimental feature and may change in future releases.536537## References538539- Bram van Leer (1982)540Flux-Vector Splitting for the Euler Equation541[DOI: 10.1007/978-3-642-60543-7_5](https://doi.org/10.1007/978-3-642-60543-7_5)542- D. Hänel, R. Schwane and G. Seider (1987)543On the accuracy of upwind schemes for the solution of the Navier-Stokes equations544[DOI: 10.2514/6.1987-1105](https://doi.org/10.2514/6.1987-1105)545- Meng-Sing Liou and Chris J. Steffen, Jr. (1991)546High-Order Polynomial Expansions (HOPE) for Flux-Vector Splitting547[NASA Technical Memorandum](https://ntrs.nasa.gov/citations/19910016425)548"""549@inline function splitting_vanleer_haenel(u, orientation::Integer,550equations::CompressibleEulerEquations1D)551fm = splitting_vanleer_haenel(u, Val{:minus}(), orientation, equations)552fp = splitting_vanleer_haenel(u, Val{:plus}(), orientation, equations)553return fm, fp554end555556@inline function splitting_vanleer_haenel(u, ::Val{:plus}, orientation::Integer,557equations::CompressibleEulerEquations1D)558rho, rho_v1, rho_e_total = u559v1 = rho_v1 / rho560p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)561562# sound speed and enthalpy563a = sqrt(equations.gamma * p / rho)564H = (rho_e_total + p) / rho565566# signed Mach number567M = v1 / a568569p_plus = 0.5f0 * (1 + equations.gamma * M) * p570571f1p = 0.25f0 * rho * a * (M + 1)^2572f2p = f1p * v1 + p_plus573f3p = f1p * H574575return SVector(f1p, f2p, f3p)576end577578@inline function splitting_vanleer_haenel(u, ::Val{:minus}, orientation::Integer,579equations::CompressibleEulerEquations1D)580rho, rho_v1, rho_e_total = u581v1 = rho_v1 / rho582p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)583584# sound speed and enthalpy585a = sqrt(equations.gamma * p / rho)586H = (rho_e_total + p) / rho587588# signed Mach number589M = v1 / a590591p_minus = 0.5f0 * (1 - equations.gamma * M) * p592593f1m = -0.25f0 * rho * a * (M - 1)^2594f2m = f1m * v1 + p_minus595f3m = f1m * H596597return SVector(f1m, f2m, f3m)598end599600# TODO: FD601# This splitting is interesting because it can handle the "el diablo" wave602# for long time runs. Computing the eigenvalues of the operator we see603# J = jacobian_ad_forward(semi);604# lamb = eigvals(J);605# maximum(real.(lamb))606# 2.1411031631522748e-6607# So the instability of this splitting is very weak. However, the 2D variant608# of this splitting on "el diablo" still crashes early. Can we learn anything609# from the design of this splitting?610"""611splitting_coirier_vanleer(u, orientation::Integer,612equations::CompressibleEulerEquations1D)613splitting_coirier_vanleer(u, which::Union{Val{:minus}, Val{:plus}}614orientation::Integer,615equations::CompressibleEulerEquations1D)616617Splitting of the compressible Euler flux from Coirier and van Leer.618The splitting has correction terms in the pressure splitting as well as619the mass and energy flux components. The motivation for these corrections620are to handle flows at the low Mach number limit.621622Returns a tuple of the fluxes "minus" (associated with waves going into the623negative axis direction) and "plus" (associated with waves going into the624positive axis direction). If only one of the fluxes is required, use the625function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.626627!!! warning "Experimental implementation (upwind SBP)"628This is an experimental feature and may change in future releases.629630## References631632- William Coirier and Bram van Leer (1991)633Numerical flux formulas for the Euler and Navier-Stokes equations.634II - Progress in flux-vector splitting635[DOI: 10.2514/6.1991-1566](https://doi.org/10.2514/6.1991-1566)636"""637@inline function splitting_coirier_vanleer(u, orientation::Integer,638equations::CompressibleEulerEquations1D)639fm = splitting_coirier_vanleer(u, Val{:minus}(), orientation, equations)640fp = splitting_coirier_vanleer(u, Val{:plus}(), orientation, equations)641return fm, fp642end643644@inline function splitting_coirier_vanleer(u, ::Val{:plus}, orientation::Integer,645equations::CompressibleEulerEquations1D)646rho, rho_v1, rho_e_total = u647v1 = rho_v1 / rho648p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)649650# sound speed and enthalpy651a = sqrt(equations.gamma * p / rho)652H = (rho_e_total + p) / rho653654# signed Mach number655M = v1 / a656657P = 2658mu = 1659nu = 0.75f0660omega = 2 # adjusted from suggested value of 1.5661662p_plus = 0.25f0 * ((M + 1)^2 * (2 - M) - nu * M * (M^2 - 1)^P) * p663664f1p = 0.25f0 * rho * a * ((M + 1)^2 - mu * (M^2 - 1)^P)665f2p = f1p * v1 + p_plus666f3p = f1p * H - omega * rho * a^3 * M^2 * (M^2 - 1)^2667668return SVector(f1p, f2p, f3p)669end670671@inline function splitting_coirier_vanleer(u, ::Val{:minus}, orientation::Integer,672equations::CompressibleEulerEquations1D)673rho, rho_v1, rho_e_total = u674v1 = rho_v1 / rho675p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)676677# sound speed and enthalpy678a = sqrt(equations.gamma * p / rho)679H = (rho_e_total + p) / rho680681# signed Mach number682M = v1 / a683684P = 2685mu = 1686nu = 0.75f0687omega = 2 # adjusted from suggested value of 1.5688689p_minus = 0.25f0 * ((M - 1)^2 * (2 + M) + nu * M * (M^2 - 1)^P) * p690691f1m = -0.25f0 * rho * a * ((M - 1)^2 - mu * (M^2 - 1)^P)692f2m = f1m * v1 + p_minus693f3m = f1m * H + omega * rho * a^3 * M^2 * (M^2 - 1)^2694695return SVector(f1m, f2m, f3m)696end697698# Calculate estimates for maximum wave speed for local Lax-Friedrichs-type dissipation as the699# maximum velocity magnitude plus the maximum speed of sound700@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,701equations::CompressibleEulerEquations1D)702rho_ll, rho_v1_ll, rho_e_total_ll = u_ll703rho_rr, rho_v1_rr, rho_e_total_rr = u_rr704705# Calculate primitive variables and speed of sound706v1_ll = rho_v1_ll / rho_ll707v_mag_ll = abs(v1_ll)708p_ll = (equations.gamma - 1) * (rho_e_total_ll - 0.5f0 * rho_ll * v_mag_ll^2)709c_ll = sqrt(equations.gamma * p_ll / rho_ll)710v1_rr = rho_v1_rr / rho_rr711v_mag_rr = abs(v1_rr)712p_rr = (equations.gamma - 1) * (rho_e_total_rr - 0.5f0 * rho_rr * v_mag_rr^2)713c_rr = sqrt(equations.gamma * p_rr / rho_rr)714715return max(v_mag_ll, v_mag_rr) + max(c_ll, c_rr)716end717718@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,719normal_direction::AbstractVector,720equations::CompressibleEulerEquations1D)721gamma = equations.gamma722723norm_ = norm(normal_direction)724unit_normal_direction = normal_direction / norm_725726rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)727rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)728729b_ll = rho_ll / (2 * p_ll)730b_rr = rho_rr / (2 * p_rr)731732rho_log = ln_mean(rho_ll, rho_rr)733b_log = ln_mean(b_ll, b_rr)734v1_avg = 0.5f0 * (v1_ll + v1_rr)735p_avg = 0.5f0 * (rho_ll + rho_rr) / (b_ll + b_rr) # 2 * b_avg = b_ll + b_rr736v1_squared_bar = v1_ll * v1_rr737h_bar = gamma / (2 * b_log * (gamma - 1)) + 0.5f0 * v1_squared_bar738c_bar = sqrt(gamma * p_avg / rho_log)739740v1_minus_c = v1_avg - c_bar * unit_normal_direction[1]741v1_plus_c = v1_avg + c_bar * unit_normal_direction[1]742v_avg_normal = v1_avg * unit_normal_direction[1]743744entropy_vars_jump = cons2entropy(u_rr, equations) - cons2entropy(u_ll, equations)745746lambda_1 = abs(v_avg_normal - c_bar) * rho_log / (2 * gamma)747lambda_2 = abs(v_avg_normal) * rho_log * (gamma - 1) / gamma748lambda_3 = abs(v_avg_normal + c_bar) * rho_log / (2 * gamma)749lambda_4 = abs(v_avg_normal) * p_avg750751entropy_var_rho_jump, entropy_var_rho_v1_jump, entropy_var_rho_e_total_jump = entropy_vars_jump752753w1 = lambda_1 * (entropy_var_rho_jump + v1_minus_c * entropy_var_rho_v1_jump +754(h_bar - c_bar * v_avg_normal) * entropy_var_rho_e_total_jump)755w2 = lambda_2 * (entropy_var_rho_jump + v1_avg * entropy_var_rho_v1_jump +756v1_squared_bar / 2 * entropy_var_rho_e_total_jump)757w3 = lambda_3 * (entropy_var_rho_jump + v1_plus_c * entropy_var_rho_v1_jump +758(h_bar + c_bar * v_avg_normal) * entropy_var_rho_e_total_jump)759760dissipation_rho = w1 + w2 + w3761762dissipation_rho_v1 = (w1 * v1_minus_c +763w2 * v1_avg +764w3 * v1_plus_c)765766dissipation_rhoe = (w1 * (h_bar - c_bar * v_avg_normal) +767w2 * 0.5f0 * v1_squared_bar +768w3 * (h_bar + c_bar * v_avg_normal) +769lambda_4 *770(entropy_var_rho_e_total_jump *771(v1_avg * v1_avg - v_avg_normal^2)))772773return -0.5f0 * SVector(dissipation_rho, dissipation_rho_v1, dissipation_rhoe) *774norm_775end776777# Calculate estimates for minimum and maximum wave speeds for HLL-type fluxes778@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,779equations::CompressibleEulerEquations1D)780rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)781rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)782783λ_min = v1_ll - sqrt(equations.gamma * p_ll / rho_ll)784λ_max = v1_rr + sqrt(equations.gamma * p_rr / rho_rr)785786return λ_min, λ_max787end788789# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`790@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,791equations::CompressibleEulerEquations1D)792rho_ll, rho_v1_ll, rho_e_total_ll = u_ll793rho_rr, rho_v1_rr, rho_e_total_rr = u_rr794795# Calculate primitive variables and speed of sound796v1_ll = rho_v1_ll / rho_ll797v_mag_ll = abs(v1_ll)798p_ll = (equations.gamma - 1) * (rho_e_total_ll - 0.5f0 * rho_ll * v_mag_ll^2)799c_ll = sqrt(equations.gamma * p_ll / rho_ll)800v1_rr = rho_v1_rr / rho_rr801v_mag_rr = abs(v1_rr)802p_rr = (equations.gamma - 1) * (rho_e_total_rr - 0.5f0 * rho_rr * v_mag_rr^2)803c_rr = sqrt(equations.gamma * p_rr / rho_rr)804805return max(v_mag_ll + c_ll, v_mag_rr + c_rr)806end807808# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes809@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,810equations::CompressibleEulerEquations1D)811rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)812rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)813814c_ll = sqrt(equations.gamma * p_ll / rho_ll)815c_rr = sqrt(equations.gamma * p_rr / rho_rr)816817λ_min = min(v1_ll - c_ll, v1_rr - c_rr)818λ_max = max(v1_ll + c_ll, v1_rr + c_rr)819820return λ_min, λ_max821end822823"""824flux_hllc(u_ll, u_rr, orientation, equations::CompressibleEulerEquations1D)825826Computes the HLLC flux (HLL with Contact) for compressible Euler equations developed by E.F. Toro827[Lecture slides](http://www.prague-sum.com/download/2012/Toro_2-HLLC-RiemannSolver.pdf)828Signal speeds: [DOI: 10.1137/S1064827593260140](https://doi.org/10.1137/S1064827593260140)829"""830function flux_hllc(u_ll, u_rr, orientation::Integer,831equations::CompressibleEulerEquations1D)832# Calculate primitive variables and speed of sound833rho_ll, rho_v1_ll, rho_e_total_ll = u_ll834rho_rr, rho_v1_rr, rho_e_total_rr = u_rr835836v1_ll = rho_v1_ll / rho_ll837e_ll = rho_e_total_ll / rho_ll838p_ll = (equations.gamma - 1) * (rho_e_total_ll - 0.5f0 * rho_ll * v1_ll^2)839c_ll = sqrt(equations.gamma * p_ll / rho_ll)840841v1_rr = rho_v1_rr / rho_rr842e_rr = rho_e_total_rr / rho_rr843p_rr = (equations.gamma - 1) * (rho_e_total_rr - 0.5f0 * rho_rr * v1_rr^2)844c_rr = sqrt(equations.gamma * p_rr / rho_rr)845846# Obtain left and right fluxes847f_ll = flux(u_ll, orientation, equations)848f_rr = flux(u_rr, orientation, equations)849850# Compute Roe averages851sqrt_rho_ll = sqrt(rho_ll)852sqrt_rho_rr = sqrt(rho_rr)853sum_sqrt_rho = sqrt_rho_ll + sqrt_rho_rr854vel_L = v1_ll855vel_R = v1_rr856vel_roe = (sqrt_rho_ll * vel_L + sqrt_rho_rr * vel_R) / sum_sqrt_rho857ekin_roe = 0.5f0 * vel_roe^2858H_ll = (rho_e_total_ll + p_ll) / rho_ll859H_rr = (rho_e_total_rr + p_rr) / rho_rr860H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) / sum_sqrt_rho861c_roe = sqrt((equations.gamma - 1) * (H_roe - ekin_roe))862863Ssl = min(vel_L - c_ll, vel_roe - c_roe)864Ssr = max(vel_R + c_rr, vel_roe + c_roe)865sMu_L = Ssl - vel_L866sMu_R = Ssr - vel_R867if Ssl >= 0868f1 = f_ll[1]869f2 = f_ll[2]870f3 = f_ll[3]871elseif Ssr <= 0872f1 = f_rr[1]873f2 = f_rr[2]874f3 = f_rr[3]875else876SStar = (p_rr - p_ll + rho_ll * vel_L * sMu_L - rho_rr * vel_R * sMu_R) /877(rho_ll * sMu_L - rho_rr * sMu_R)878if Ssl <= 0 <= SStar879densStar = rho_ll * sMu_L / (Ssl - SStar)880enerStar = e_ll + (SStar - vel_L) * (SStar + p_ll / (rho_ll * sMu_L))881UStar1 = densStar882UStar2 = densStar * SStar883UStar3 = densStar * enerStar884885f1 = f_ll[1] + Ssl * (UStar1 - rho_ll)886f2 = f_ll[2] + Ssl * (UStar2 - rho_v1_ll)887f3 = f_ll[3] + Ssl * (UStar3 - rho_e_total_ll)888else889densStar = rho_rr * sMu_R / (Ssr - SStar)890enerStar = e_rr + (SStar - vel_R) * (SStar + p_rr / (rho_rr * sMu_R))891UStar1 = densStar892UStar2 = densStar * SStar893UStar3 = densStar * enerStar894895#end896f1 = f_rr[1] + Ssr * (UStar1 - rho_rr)897f2 = f_rr[2] + Ssr * (UStar2 - rho_v1_rr)898f3 = f_rr[3] + Ssr * (UStar3 - rho_e_total_rr)899end900end901return SVector(f1, f2, f3)902end903904# While `normal_direction` isn't strictly necessary in 1D, certain solvers assume that905# the normal component is incorporated into the numerical flux.906#907# The HLLC flux along a 1D "normal" can be evaluated by scaling the velocity/momentum by908# the normal for the 1D HLLC flux, then scaling the resulting momentum flux again.909# Moreover, the 2D HLLC flux reduces to this if the normal vector is [n, 0].910function flux_hllc(u_ll, u_rr, normal_direction::AbstractVector,911equations::CompressibleEulerEquations1D)912norm_ = abs(normal_direction[1])913normal_direction_unit = normal_direction[1] * inv(norm_)914915# scale the momentum by the normal direction916f = flux_hllc(SVector(u_ll[1], normal_direction_unit * u_ll[2], u_ll[3]),917SVector(u_rr[1], normal_direction_unit * u_rr[2], u_rr[3]), 1,918equations)919920# rescale the momentum flux by the normal direction and normalize921return SVector(f[1], normal_direction_unit * f[2], f[3]) * norm_922end923924"""925min_max_speed_einfeldt(u_ll, u_rr, orientation, equations::CompressibleEulerEquations1D)926927Computes the HLLE (Harten-Lax-van Leer-Einfeldt) flux for the compressible Euler equations.928Special estimates of the signal velocites and linearization of the Riemann problem developed929by Einfeldt to ensure that the internal energy and density remain positive during the computation930of the numerical flux.931932Original publication:933- Bernd Einfeldt (1988)934On Godunov-type methods for gas dynamics.935[DOI: 10.1137/0725021](https://doi.org/10.1137/0725021)936937Compactly summarized:938- Siddhartha Mishra, Ulrik Skre Fjordholm and Rémi Abgrall939Numerical methods for conservation laws and related equations.940[Link](https://metaphor.ethz.ch/x/2019/hs/401-4671-00L/literature/mishra_hyperbolic_pdes.pdf)941"""942@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,943equations::CompressibleEulerEquations1D)944# Calculate primitive variables, enthalpy and speed of sound945rho_ll, v_ll, p_ll = cons2prim(u_ll, equations)946rho_rr, v_rr, p_rr = cons2prim(u_rr, equations)947948# `u_ll[3]` is total energy `rho_e_total_ll` on the left949H_ll = (u_ll[3] + p_ll) / rho_ll950c_ll = sqrt(equations.gamma * p_ll / rho_ll)951952# `u_rr[3]` is total energy `rho_e_total_rr` on the right953H_rr = (u_rr[3] + p_rr) / rho_rr954c_rr = sqrt(equations.gamma * p_rr / rho_rr)955956# Compute Roe averages957sqrt_rho_ll = sqrt(rho_ll)958sqrt_rho_rr = sqrt(rho_rr)959inv_sum_sqrt_rho = inv(sqrt_rho_ll + sqrt_rho_rr)960961v_roe = (sqrt_rho_ll * v_ll + sqrt_rho_rr * v_rr) * inv_sum_sqrt_rho962v_roe_mag = v_roe^2963964H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) * inv_sum_sqrt_rho965c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * v_roe_mag))966967# Compute convenience constant for positivity preservation, see968# https://doi.org/10.1016/0021-9991(91)90211-3969beta = sqrt(0.5f0 * (equations.gamma - 1) / equations.gamma)970971# Estimate the edges of the Riemann fan (with positivity conservation)972SsL = min(v_roe - c_roe, v_ll - beta * c_ll, 0)973SsR = max(v_roe + c_roe, v_rr + beta * c_rr, 0)974975return SsL, SsR976end977978@inline function max_abs_speeds(u, equations::CompressibleEulerEquations1D)979rho, rho_v1, rho_e_total = u980v1 = rho_v1 / rho981p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho * v1^2)982c = sqrt(equations.gamma * p / rho)983984return (abs(v1) + c,)985end986987# Convert conservative variables to primitive988@inline function cons2prim(u, equations::CompressibleEulerEquations1D)989rho, rho_v1, rho_e_total = u990991v1 = rho_v1 / rho992p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)993994return SVector(rho, v1, p)995end996997# Convert conservative variables to entropy998@inline function cons2entropy(u, equations::CompressibleEulerEquations1D)999rho, rho_v1, rho_e_total = u10001001v1 = rho_v1 / rho1002v_square = v1^21003p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho * v_square)1004s = log(p) - equations.gamma * log(rho)1005rho_p = rho / p10061007w1 = (equations.gamma - s) * equations.inv_gamma_minus_one -10080.5f0 * rho_p * v_square1009w2 = rho_p * v11010w3 = -rho_p10111012return SVector(w1, w2, w3)1013end10141015@inline function entropy2cons(w, equations::CompressibleEulerEquations1D)1016# See Hughes, Franca, Mallet (1986) A new finite element formulation for CFD1017# [DOI: 10.1016/0045-7825(86)90127-1](https://doi.org/10.1016/0045-7825(86)90127-1)1018@unpack gamma = equations10191020# convert to entropy `-rho * s` used by Hughes, France, Mallet (1986)1021# instead of `-rho * s / (gamma - 1)`1022V1, V2, V5 = w .* (gamma - 1)10231024# specific entropy, eq. (53)1025s = gamma - V1 + 0.5f0 * (V2^2) / V510261027# eq. (52)1028energy_internal = ((gamma - 1) / (-V5)^gamma)^(equations.inv_gamma_minus_one) *1029exp(-s * equations.inv_gamma_minus_one)10301031# eq. (51)1032rho = -V5 * energy_internal1033rho_v1 = V2 * energy_internal1034rho_e_total = (1 - 0.5f0 * (V2^2) / V5) * energy_internal1035return SVector(rho, rho_v1, rho_e_total)1036end10371038# Convert primitive to conservative variables1039@inline function prim2cons(prim, equations::CompressibleEulerEquations1D)1040rho, v1, p = prim1041rho_v1 = rho * v11042rho_e_total = p * equations.inv_gamma_minus_one + 0.5f0 * (rho_v1 * v1)1043return SVector(rho, rho_v1, rho_e_total)1044end10451046@inline function density(u, equations::CompressibleEulerEquations1D)1047rho = u[1]1048return rho1049end10501051@inline function velocity(u, orientation_or_normal,1052equations::CompressibleEulerEquations1D)1053return velocity(u, equations)1054end10551056@inline function velocity(u, equations::CompressibleEulerEquations1D)1057rho = u[1]1058v1 = u[2] / rho1059return v11060end10611062@doc raw"""1063pressure(u, equations::AbstractCompressibleEulerEquations)10641065Computes the pressure for an ideal equation of state with1066isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`.1067```math1068\begin{aligned}1069p &= (\gamma - 1) \left( E_\mathrm{tot} - E_\mathrm{kin} \right) \\1070&= (\gamma - 1) \left( \rho e - \frac{1}{2}\rho \Vert v \Vert_2^2 \right)1071\end{aligned}1072```1073"""1074@inline function pressure(u, equations::CompressibleEulerEquations1D)1075rho, rho_v1, rho_e_total = u1076p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1^2) / rho)1077return p1078end10791080@inline function density_pressure(u, equations::CompressibleEulerEquations1D)1081rho, rho_v1, rho_e_total = u1082rho_times_p = (equations.gamma - 1) * (rho * rho_e_total - 0.5f0 * (rho_v1^2))1083return rho_times_p1084end10851086@doc raw"""1087entropy_thermodynamic(cons, equations::AbstractCompressibleEulerEquations)10881089Calculate thermodynamic entropy for a conservative state `cons` as10901091```math1092s = \log(p) - \gamma \log(\rho)1093```1094"""1095@inline function entropy_thermodynamic(cons, equations::CompressibleEulerEquations1D)1096# Pressure1097p = (equations.gamma - 1) * (cons[3] - 0.5f0 * (cons[2]^2) / cons[1])10981099# Thermodynamic entropy1100s = log(p) - equations.gamma * log(cons[1])11011102return s1103end11041105@doc raw"""1106entropy_math(cons, equations::AbstractCompressibleEulerEquations)11071108Calculate mathematical entropy for a conservative state `cons` as1109```math1110S = -\frac{\rho s}{\gamma - 1}1111```1112where `s` is the thermodynamic entropy calculated by1113[`entropy_thermodynamic(cons, equations::AbstractCompressibleEulerEquations)`](@ref).1114"""1115@inline function entropy_math(cons, equations::CompressibleEulerEquations1D)1116# Mathematical entropy1117S = -entropy_thermodynamic(cons, equations) * cons[1] *1118equations.inv_gamma_minus_one11191120return S1121end11221123"""1124entropy(cons, equations::AbstractCompressibleEulerEquations)11251126Default entropy is the mathematical entropy1127[`entropy_math(cons, equations::AbstractCompressibleEulerEquations)`](@ref).1128"""1129@inline function entropy(cons, equations::CompressibleEulerEquations1D)1130return entropy_math(cons, equations)1131end11321133# Calculate total energy for a conservative state `cons`1134@inline energy_total(cons, ::CompressibleEulerEquations1D) = cons[3]11351136# Calculate kinetic energy for a conservative state `cons`1137@inline function energy_kinetic(cons, equations::CompressibleEulerEquations1D)1138return 0.5f0 * (cons[2]^2) / cons[1]1139end11401141"""1142energy_internal(cons, equations::AbstractCompressibleEulerEquations)11431144Calculate internal energy ``e_{\text{internal}}`` for a conservative state `cons` as the difference1145of total energy and kinetic energy.1146"""1147@inline function energy_internal(cons, equations::CompressibleEulerEquations1D)1148return energy_total(cons, equations) - energy_kinetic(cons, equations)1149end11501151@doc raw"""1152entropy_potential(u, orientation_or_normal_direction,1153equations::AbstractCompressibleEulerEquations)11541155Calculate the entropy potential, which for the compressible Euler equations is simply1156the momentum for the choice of mathematical [`entropy`](@ref) ``S(u) = -\frac{\rho s}{\gamma - 1}``1157with thermodynamic entropy ``s = \ln(p) - \gamma \ln(\rho)``.11581159## References1160- Eitan Tadmor (2003)1161Entropy stability theory for difference approximations of nonlinear conservation laws and related time-dependent problems1162[DOI: 10.1017/S0962492902000156](https://doi.org/10.1017/S0962492902000156)1163"""1164@inline function entropy_potential(u, orientation::Int,1165equations::CompressibleEulerEquations1D)1166return u[2]1167end1168end # @muladd116911701171