Path: blob/main/src/equations/compressible_euler_2d.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"""8CompressibleEulerEquations2D(gamma)910The compressible Euler equations11```math12\frac{\partial}{\partial t}13\begin{pmatrix}14\rho \\ \rho v_1 \\ \rho v_2 \\ \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 v_1 v_2 \\ (\rho e_{\text{total}} +p) v_120\end{pmatrix}21+22\frac{\partial}{\partial y}23\begin{pmatrix}24\rho v_2 \\ \rho v_1 v_2 \\ \rho v_2^2 + p \\ (\rho e_{\text{total}} +p) v_225\end{pmatrix}26=27\begin{pmatrix}280 \\ 0 \\ 0 \\ 029\end{pmatrix}30```31for an ideal gas with ratio of specific heats `gamma`32in two space dimensions.33Here, ``\rho`` is the density, ``v_1``, ``v_2`` the velocities, ``e_{\text{total}}`` the specific total energy, and34```math35p = (\gamma - 1) \left( \rho e_{\text{total}} - \frac{1}{2} \rho (v_1^2+v_2^2) \right)36```37the pressure.38"""39struct CompressibleEulerEquations2D{RealT <: Real} <:40AbstractCompressibleEulerEquations{2, 4}41gamma::RealT # ratio of specific heats42inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications4344function CompressibleEulerEquations2D(gamma)45γ, inv_gamma_minus_one = promote(gamma, inv(gamma - 1))46return new{typeof(γ)}(γ, inv_gamma_minus_one)47end48end4950function varnames(::typeof(cons2cons), ::CompressibleEulerEquations2D)51return ("rho", "rho_v1", "rho_v2", "rho_e_total")52end53varnames(::typeof(cons2prim), ::CompressibleEulerEquations2D) = ("rho", "v1", "v2", "p")5455# Set initial conditions at physical location `x` for time `t`56"""57initial_condition_constant(x, t, equations::CompressibleEulerEquations2D)5859A constant initial condition to test free-stream preservation.60"""61function initial_condition_constant(x, t, equations::CompressibleEulerEquations2D)62RealT = eltype(x)63rho = 164rho_v1 = convert(RealT, 0.1)65rho_v2 = convert(RealT, -0.2)66rho_e_total = 1067return SVector(rho, rho_v1, rho_v2, rho_e_total)68end6970"""71initial_condition_convergence_test(x, t, equations::CompressibleEulerEquations2D)7273A smooth initial condition used for convergence tests in combination with74[`source_terms_convergence_test`](@ref)75(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).76"""77function initial_condition_convergence_test(x, t,78equations::CompressibleEulerEquations2D)79RealT = eltype(x)80c = 281A = convert(RealT, 0.1)82L = 283f = 1.0f0 / L84ω = 2 * convert(RealT, pi) * f85ini = c + A * sin(ω * (x[1] + x[2] - t))8687rho = ini88rho_v1 = ini89rho_v2 = ini90rho_e_total = ini^29192return SVector(rho, rho_v1, rho_v2, rho_e_total)93end9495"""96source_terms_convergence_test(u, x, t, equations::CompressibleEulerEquations2D)9798Source terms used for convergence tests in combination with99[`initial_condition_convergence_test`](@ref)100(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).101102References for the method of manufactured solutions (MMS):103- Kambiz Salari and Patrick Knupp (2000)104Code Verification by the Method of Manufactured Solutions105[DOI: 10.2172/759450](https://doi.org/10.2172/759450)106- Patrick J. Roache (2002)107Code Verification by the Method of Manufactured Solutions108[DOI: 10.1115/1.1436090](https://doi.org/10.1115/1.1436090)109"""110@inline function source_terms_convergence_test(u, x, t,111equations::CompressibleEulerEquations2D)112# Same settings as in `initial_condition`113RealT = eltype(u)114c = 2115A = convert(RealT, 0.1)116L = 2117f = 1.0f0 / L118ω = 2 * convert(RealT, pi) * f119γ = equations.gamma120121x1, x2 = x122si, co = sincos(ω * (x1 + x2 - t))123rho = c + A * si124rho_x = ω * A * co125# Note that d/dt rho = -d/dx rho = -d/dy rho.126127tmp = (2 * rho - 1) * (γ - 1)128129du1 = rho_x130du2 = rho_x * (1 + tmp)131du3 = du2132du4 = 2 * rho_x * (rho + tmp)133134return SVector(du1, du2, du3, du4)135end136137"""138initial_condition_density_wave(x, t, equations::CompressibleEulerEquations2D)139140A sine wave in the density with constant velocity and pressure; reduces the141compressible Euler equations to the linear advection equations.142This setup is the test case for stability of EC fluxes from paper143- Gregor J. Gassner, Magnus Svärd, Florian J. Hindenlang (2020)144Stability issues of entropy-stable and/or split-form high-order schemes145[arXiv: 2007.09026](https://arxiv.org/abs/2007.09026)146with the following parameters147- domain [-1, 1]148- mesh = 4x4149- polydeg = 5150"""151function initial_condition_density_wave(x, t, equations::CompressibleEulerEquations2D)152RealT = eltype(x)153v1 = convert(RealT, 0.1)154v2 = convert(RealT, 0.2)155rho = 1 + convert(RealT, 0.98) * sinpi(2 * (x[1] + x[2] - t * (v1 + v2)))156rho_v1 = rho * v1157rho_v2 = rho * v2158p = 20159rho_e_total = p / (equations.gamma - 1) + 0.5f0 * rho * (v1^2 + v2^2)160return SVector(rho, rho_v1, rho_v2, rho_e_total)161end162163"""164initial_condition_weak_blast_wave(x, t, equations::CompressibleEulerEquations2D)165166A weak blast wave taken from167- Sebastian Hennemann, Gregor J. Gassner (2020)168A provably entropy stable subcell shock capturing approach for high order split form DG169[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)170"""171function initial_condition_weak_blast_wave(x, t,172equations::CompressibleEulerEquations2D)173# From Hennemann & Gassner JCP paper 2020 (Sec. 6.3)174# Set up polar coordinates175inicenter = SVector(0, 0)176x_norm = x[1] - inicenter[1]177y_norm = x[2] - inicenter[2]178r = sqrt(x_norm^2 + y_norm^2)179phi = atan(y_norm, x_norm)180sin_phi, cos_phi = sincos(phi)181182# Calculate primitive variables183RealT = eltype(x)184rho = r > 0.5f0 ? one(RealT) : convert(RealT, 1.1691)185v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos_phi186v2 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * sin_phi187p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)188189return prim2cons(SVector(rho, v1, v2, p), equations)190end191192"""193initial_condition_eoc_test_coupled_euler_gravity(x, t, equations::CompressibleEulerEquations2D)194195Setup used for convergence tests of the Euler equations with self-gravity used in196- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)197A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics198[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)199in combination with [`source_terms_eoc_test_coupled_euler_gravity`](@ref)200or [`source_terms_eoc_test_euler`](@ref).201"""202function initial_condition_eoc_test_coupled_euler_gravity(x, t,203equations::CompressibleEulerEquations2D)204# OBS! this assumes that γ = 2 other manufactured source terms are incorrect205if equations.gamma != 2206error("adiabatic constant must be 2 for the coupling convergence test")207end208RealT = eltype(x)209c = 2210A = convert(RealT, 0.1)211ini = c + A * sinpi(x[1] + x[2] - t)212G = 1 # gravitational constant213214rho = ini215v1 = 1216v2 = 1217p = ini^2 * G / convert(RealT, pi) # * 2 / ndims, but ndims==2 here218219return prim2cons(SVector(rho, v1, v2, p), equations)220end221222"""223source_terms_eoc_test_coupled_euler_gravity(u, x, t, equations::CompressibleEulerEquations2D)224225Setup used for convergence tests of the Euler equations with self-gravity used in226- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)227A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics228[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)229in combination with [`initial_condition_eoc_test_coupled_euler_gravity`](@ref).230"""231@inline function source_terms_eoc_test_coupled_euler_gravity(u, x, t,232equations::CompressibleEulerEquations2D)233# Same settings as in `initial_condition_eoc_test_coupled_euler_gravity`234RealT = eltype(u)235c = 2236A = convert(RealT, 0.1)237G = 1 # gravitational constant, must match coupling solver238C_grav = -2 * G / convert(RealT, pi) # 2 == 4 / ndims239240x1, x2 = x241si, co = sincospi(x1 + x2 - t)242rhox = A * convert(RealT, pi) * co243rho = c + A * si244245du1 = rhox246du2 = rhox247du3 = rhox248du4 = (1 - C_grav * rho) * rhox249250return SVector(du1, du2, du3, du4)251end252253"""254source_terms_eoc_test_euler(u, x, t, equations::CompressibleEulerEquations2D)255256Setup used for convergence tests of the Euler equations with self-gravity used in257- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)258A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics259[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)260in combination with [`initial_condition_eoc_test_coupled_euler_gravity`](@ref).261"""262@inline function source_terms_eoc_test_euler(u, x, t,263equations::CompressibleEulerEquations2D)264# Same settings as in `initial_condition_eoc_test_coupled_euler_gravity`265RealT = eltype(u)266c = 2267A = convert(RealT, 0.1)268G = 1269C_grav = -2 * G / convert(RealT, pi) # 2 == 4 / ndims270271x1, x2 = x272si, co = sincospi(x1 + x2 - t)273rhox = A * convert(RealT, pi) * co274rho = c + A * si275276du1 = rhox277du2 = rhox * (1 - C_grav * rho)278du3 = rhox * (1 - C_grav * rho)279du4 = rhox * (1 - 3 * C_grav * rho)280281return SVector(du1, du2, du3, du4)282end283284"""285boundary_condition_slip_wall(u_inner, normal_direction, x, t, surface_flux_function,286equations::CompressibleEulerEquations2D)287288Determine the boundary numerical surface flux for a slip wall condition.289Imposes a zero normal velocity at the wall.290Density is taken from the internal solution state and pressure is computed as an291exact solution of a 1D Riemann problem. Further details about this boundary state292are available in the paper:293- J. J. W. van der Vegt and H. van der Ven (2002)294Slip flow boundary conditions in discontinuous Galerkin discretizations of295the Euler equations of gas dynamics296[PDF](https://reports.nlr.nl/bitstream/handle/10921/692/TP-2002-300.pdf?sequence=1)297298Details about the 1D pressure Riemann solution can be found in Section 6.3.3 of the book299- Eleuterio F. Toro (2009)300Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction3013rd edition302[DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)303304Should be used together with [`UnstructuredMesh2D`](@ref), [`P4estMesh`](@ref), or [`T8codeMesh`](@ref).305"""306@inline function boundary_condition_slip_wall(u_inner, normal_direction::AbstractVector,307x, t,308surface_flux_function,309equations::CompressibleEulerEquations2D)310norm_ = norm(normal_direction)311# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later312normal = normal_direction / norm_313314# rotate the internal solution state315u_local = rotate_to_x(u_inner, normal, equations)316317# compute the primitive variables318rho_local, v_normal, v_tangent, p_local = cons2prim(u_local, equations)319320# Get the solution of the pressure Riemann problem321# See Section 6.3.3 of322# Eleuterio F. Toro (2009)323# Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction324# [DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)325if v_normal <= 0326sound_speed = sqrt(equations.gamma * p_local / rho_local) # local sound speed327p_star = p_local *328(1 + 0.5f0 * (equations.gamma - 1) * v_normal / sound_speed)^(2 *329equations.gamma *330equations.inv_gamma_minus_one)331else # v_normal > 0332A = 2 / ((equations.gamma + 1) * rho_local)333B = p_local * (equations.gamma - 1) / (equations.gamma + 1)334p_star = p_local +3350.5f0 * v_normal / A *336(v_normal + sqrt(v_normal^2 + 4 * A * (p_local + B)))337end338339# For the slip wall we directly set the flux as the normal velocity is zero340return SVector(0,341p_star * normal[1],342p_star * normal[2],3430) * norm_344end345346"""347boundary_condition_slip_wall(u_inner, orientation, direction, x, t,348surface_flux_function, equations::CompressibleEulerEquations2D)349350Should be used together with [`TreeMesh`](@ref).351"""352@inline function boundary_condition_slip_wall(u_inner, orientation,353direction, x, t,354surface_flux_function,355equations::CompressibleEulerEquations2D)356# get the appropriate normal vector from the orientation357RealT = eltype(u_inner)358if orientation == 1359normal_direction = SVector(one(RealT), zero(RealT))360else # orientation == 2361normal_direction = SVector(zero(RealT), one(RealT))362end363364# compute and return the flux using `boundary_condition_slip_wall` routine below365return boundary_condition_slip_wall(u_inner, normal_direction, direction,366x, t, surface_flux_function, equations)367end368369"""370boundary_condition_slip_wall(u_inner, normal_direction, direction, x, t,371surface_flux_function, equations::CompressibleEulerEquations2D)372373Should be used together with [`StructuredMesh`](@ref).374"""375@inline function boundary_condition_slip_wall(u_inner, normal_direction::AbstractVector,376direction, x, t,377surface_flux_function,378equations::CompressibleEulerEquations2D)379# flip sign of normal to make it outward pointing, then flip the sign of the normal flux back380# to be inward pointing on the -x and -y sides due to the orientation convention used by StructuredMesh381if isodd(direction)382boundary_flux = -boundary_condition_slip_wall(u_inner, -normal_direction,383x, t, surface_flux_function,384equations)385else386boundary_flux = boundary_condition_slip_wall(u_inner, normal_direction,387x, t, surface_flux_function,388equations)389end390391return boundary_flux392end393394# Calculate 1D flux for a single point395@inline function flux(u, orientation::Integer, equations::CompressibleEulerEquations2D)396rho, rho_v1, rho_v2, rho_e_total = u397v1 = rho_v1 / rho398v2 = rho_v2 / rho399p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))400if orientation == 1401f1 = rho_v1402f2 = rho_v1 * v1 + p403f3 = rho_v1 * v2404f4 = (rho_e_total + p) * v1405else406f1 = rho_v2407f2 = rho_v2 * v1408f3 = rho_v2 * v2 + p409f4 = (rho_e_total + p) * v2410end411return SVector(f1, f2, f3, f4)412end413414# Calculate 1D flux for a single point in the normal direction415# Note, this directional vector is not normalized416@inline function flux(u, normal_direction::AbstractVector,417equations::CompressibleEulerEquations2D)418rho_e_total = last(u)419rho, v1, v2, p = cons2prim(u, equations)420421v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]422rho_v_normal = rho * v_normal423f1 = rho_v_normal424f2 = rho_v_normal * v1 + p * normal_direction[1]425f3 = rho_v_normal * v2 + p * normal_direction[2]426f4 = (rho_e_total + p) * v_normal427return SVector(f1, f2, f3, f4)428end429430"""431flux_shima_etal(u_ll, u_rr, orientation_or_normal_direction,432equations::CompressibleEulerEquations2D)433434This flux is is a modification of the original kinetic energy preserving two-point flux by435- Yuichi Kuya, Kosuke Totani and Soshi Kawai (2018)436Kinetic energy and entropy preserving schemes for compressible flows437by split convective forms438[DOI: 10.1016/j.jcp.2018.08.058](https://doi.org/10.1016/j.jcp.2018.08.058)439440The modification is in the energy flux to guarantee pressure equilibrium and was developed by441- Nao Shima, Yuichi Kuya, Yoshiharu Tamaki, Soshi Kawai (JCP 2020)442Preventing spurious pressure oscillations in split convective form discretizations for443compressible flows444[DOI: 10.1016/j.jcp.2020.110060](https://doi.org/10.1016/j.jcp.2020.110060)445"""446@inline function flux_shima_etal(u_ll, u_rr, orientation::Integer,447equations::CompressibleEulerEquations2D)448# Unpack left and right state449rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)450rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)451452# Average each factor of products in flux453rho_avg = 0.5f0 * (rho_ll + rho_rr)454v1_avg = 0.5f0 * (v1_ll + v1_rr)455v2_avg = 0.5f0 * (v2_ll + v2_rr)456p_avg = 0.5f0 * (p_ll + p_rr)457kin_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr)458459# Calculate fluxes depending on orientation460if orientation == 1461pv1_avg = 0.5f0 * (p_ll * v1_rr + p_rr * v1_ll)462f1 = rho_avg * v1_avg463f2 = f1 * v1_avg + p_avg464f3 = f1 * v2_avg465f4 = p_avg * v1_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv1_avg466else467pv2_avg = 0.5f0 * (p_ll * v2_rr + p_rr * v2_ll)468f1 = rho_avg * v2_avg469f2 = f1 * v1_avg470f3 = f1 * v2_avg + p_avg471f4 = p_avg * v2_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv2_avg472end473474return SVector(f1, f2, f3, f4)475end476477@inline function flux_shima_etal(u_ll, u_rr, normal_direction::AbstractVector,478equations::CompressibleEulerEquations2D)479# Unpack left and right state480rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)481rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)482v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]483v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]484485# Average each factor of products in flux486rho_avg = 0.5f0 * (rho_ll + rho_rr)487v1_avg = 0.5f0 * (v1_ll + v1_rr)488v2_avg = 0.5f0 * (v2_ll + v2_rr)489v_dot_n_avg = 0.5f0 * (v_dot_n_ll + v_dot_n_rr)490p_avg = 0.5f0 * (p_ll + p_rr)491velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr)492493# Calculate fluxes depending on normal_direction494f1 = rho_avg * v_dot_n_avg495f2 = f1 * v1_avg + p_avg * normal_direction[1]496f3 = f1 * v2_avg + p_avg * normal_direction[2]497f4 = (f1 * velocity_square_avg +498p_avg * v_dot_n_avg * equations.inv_gamma_minus_one499+ 0.5f0 * (p_ll * v_dot_n_rr + p_rr * v_dot_n_ll))500501return SVector(f1, f2, f3, f4)502end503504"""505flux_kennedy_gruber(u_ll, u_rr, orientation_or_normal_direction,506equations::CompressibleEulerEquations2D)507508Kinetic energy preserving two-point flux by509- Kennedy and Gruber (2008)510Reduced aliasing formulations of the convective terms within the511Navier-Stokes equations for a compressible fluid512[DOI: 10.1016/j.jcp.2007.09.020](https://doi.org/10.1016/j.jcp.2007.09.020)513"""514@inline function flux_kennedy_gruber(u_ll, u_rr, orientation::Integer,515equations::CompressibleEulerEquations2D)516# Unpack left and right state517rho_e_total_ll = last(u_ll)518rho_e_total_rr = last(u_rr)519rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)520rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)521522# Average each factor of products in flux523rho_avg = 0.5f0 * (rho_ll + rho_rr)524v1_avg = 0.5f0 * (v1_ll + v1_rr)525v2_avg = 0.5f0 * (v2_ll + v2_rr)526p_avg = 0.5f0 * (p_ll + p_rr)527e_avg = 0.5f0 * (rho_e_total_ll / rho_ll + rho_e_total_rr / rho_rr)528529# Calculate fluxes depending on orientation530if orientation == 1531f1 = rho_avg * v1_avg532f2 = rho_avg * v1_avg * v1_avg + p_avg533f3 = rho_avg * v1_avg * v2_avg534f4 = (rho_avg * e_avg + p_avg) * v1_avg535else536f1 = rho_avg * v2_avg537f2 = rho_avg * v2_avg * v1_avg538f3 = rho_avg * v2_avg * v2_avg + p_avg539f4 = (rho_avg * e_avg + p_avg) * v2_avg540end541542return SVector(f1, f2, f3, f4)543end544545@inline function flux_kennedy_gruber(u_ll, u_rr, normal_direction::AbstractVector,546equations::CompressibleEulerEquations2D)547# Unpack left and right state548rho_e_total_ll = last(u_ll)549rho_e_total_rr = last(u_rr)550rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)551rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)552553# Average each factor of products in flux554rho_avg = 0.5f0 * (rho_ll + rho_rr)555v1_avg = 0.5f0 * (v1_ll + v1_rr)556v2_avg = 0.5f0 * (v2_ll + v2_rr)557v_dot_n_avg = v1_avg * normal_direction[1] + v2_avg * normal_direction[2]558p_avg = 0.5f0 * (p_ll + p_rr)559e_avg = 0.5f0 * (rho_e_total_ll / rho_ll + rho_e_total_rr / rho_rr)560561# Calculate fluxes depending on normal_direction562f1 = rho_avg * v_dot_n_avg563f2 = f1 * v1_avg + p_avg * normal_direction[1]564f3 = f1 * v2_avg + p_avg * normal_direction[2]565f4 = f1 * e_avg + p_avg * v_dot_n_avg566567return SVector(f1, f2, f3, f4)568end569570"""571flux_chandrashekar(u_ll, u_rr, orientation_or_normal_direction, equations::CompressibleEulerEquations2D)572573Entropy conserving two-point flux by574- Chandrashekar (2013)575Kinetic Energy Preserving and Entropy Stable Finite Volume Schemes576for Compressible Euler and Navier-Stokes Equations577[DOI: 10.4208/cicp.170712.010313a](https://doi.org/10.4208/cicp.170712.010313a)578"""579@inline function flux_chandrashekar(u_ll, u_rr, orientation::Integer,580equations::CompressibleEulerEquations2D)581# Unpack left and right state582rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)583rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)584beta_ll = 0.5f0 * rho_ll / p_ll585beta_rr = 0.5f0 * rho_rr / p_rr586specific_kin_ll = 0.5f0 * (v1_ll^2 + v2_ll^2)587specific_kin_rr = 0.5f0 * (v1_rr^2 + v2_rr^2)588589# Compute the necessary mean values590rho_avg = 0.5f0 * (rho_ll + rho_rr)591rho_mean = ln_mean(rho_ll, rho_rr)592beta_mean = ln_mean(beta_ll, beta_rr)593beta_avg = 0.5f0 * (beta_ll + beta_rr)594v1_avg = 0.5f0 * (v1_ll + v1_rr)595v2_avg = 0.5f0 * (v2_ll + v2_rr)596p_mean = 0.5f0 * rho_avg / beta_avg597velocity_square_avg = specific_kin_ll + specific_kin_rr598599# Calculate fluxes depending on orientation600if orientation == 1601f1 = rho_mean * v1_avg602f2 = f1 * v1_avg + p_mean603f3 = f1 * v2_avg604f4 = f1 * 0.5f0 *605(1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +606f2 * v1_avg + f3 * v2_avg607else608f1 = rho_mean * v2_avg609f2 = f1 * v1_avg610f3 = f1 * v2_avg + p_mean611f4 = f1 * 0.5f0 *612(1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +613f2 * v1_avg + f3 * v2_avg614end615616return SVector(f1, f2, f3, f4)617end618619@inline function flux_chandrashekar(u_ll, u_rr, normal_direction::AbstractVector,620equations::CompressibleEulerEquations2D)621# Unpack left and right state622rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)623rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)624v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]625v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]626beta_ll = 0.5f0 * rho_ll / p_ll627beta_rr = 0.5f0 * rho_rr / p_rr628specific_kin_ll = 0.5f0 * (v1_ll^2 + v2_ll^2)629specific_kin_rr = 0.5f0 * (v1_rr^2 + v2_rr^2)630631# Compute the necessary mean values632rho_avg = 0.5f0 * (rho_ll + rho_rr)633rho_mean = ln_mean(rho_ll, rho_rr)634beta_mean = ln_mean(beta_ll, beta_rr)635beta_avg = 0.5f0 * (beta_ll + beta_rr)636v1_avg = 0.5f0 * (v1_ll + v1_rr)637v2_avg = 0.5f0 * (v2_ll + v2_rr)638p_mean = 0.5f0 * rho_avg / beta_avg639velocity_square_avg = specific_kin_ll + specific_kin_rr640641# Multiply with average of normal velocities642f1 = rho_mean * 0.5f0 * (v_dot_n_ll + v_dot_n_rr)643f2 = f1 * v1_avg + p_mean * normal_direction[1]644f3 = f1 * v2_avg + p_mean * normal_direction[2]645f4 = f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +646f2 * v1_avg + f3 * v2_avg647648return SVector(f1, f2, f3, f4)649end650651"""652flux_ranocha(u_ll, u_rr, orientation_or_normal_direction,653equations::CompressibleEulerEquations2D)654655Entropy conserving and kinetic energy preserving two-point flux by656- Hendrik Ranocha (2018)657Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods658for Hyperbolic Balance Laws659[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)660See also661- Hendrik Ranocha (2020)662Entropy Conserving and Kinetic Energy Preserving Numerical Methods for663the Euler Equations Using Summation-by-Parts Operators664[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)665"""666@inline function flux_ranocha(u_ll, u_rr, orientation::Integer,667equations::CompressibleEulerEquations2D)668# Unpack left and right state669rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)670rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)671672# Compute the necessary mean values673rho_mean = ln_mean(rho_ll, rho_rr)674# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`675# in exact arithmetic since676# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)677# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)678inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)679v1_avg = 0.5f0 * (v1_ll + v1_rr)680v2_avg = 0.5f0 * (v2_ll + v2_rr)681p_avg = 0.5f0 * (p_ll + p_rr)682velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr)683684# Calculate fluxes depending on orientation685if orientation == 1686f1 = rho_mean * v1_avg687f2 = f1 * v1_avg + p_avg688f3 = f1 * v2_avg689f4 = f1 *690(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +6910.5f0 * (p_ll * v1_rr + p_rr * v1_ll)692else693f1 = rho_mean * v2_avg694f2 = f1 * v1_avg695f3 = f1 * v2_avg + p_avg696f4 = f1 *697(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +6980.5f0 * (p_ll * v2_rr + p_rr * v2_ll)699end700701return SVector(f1, f2, f3, f4)702end703704@inline function flux_ranocha(u_ll, u_rr, normal_direction::AbstractVector,705equations::CompressibleEulerEquations2D)706# Unpack left and right state707rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)708rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)709v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]710v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]711712# Compute the necessary mean values713rho_mean = ln_mean(rho_ll, rho_rr)714# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`715# in exact arithmetic since716# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)717# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)718inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)719v1_avg = 0.5f0 * (v1_ll + v1_rr)720v2_avg = 0.5f0 * (v2_ll + v2_rr)721p_avg = 0.5f0 * (p_ll + p_rr)722velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr)723724# Calculate fluxes depending on normal_direction725f1 = rho_mean * 0.5f0 * (v_dot_n_ll + v_dot_n_rr)726f2 = f1 * v1_avg + p_avg * normal_direction[1]727f3 = f1 * v2_avg + p_avg * normal_direction[2]728f4 = (f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)729+7300.5f0 * (p_ll * v_dot_n_rr + p_rr * v_dot_n_ll))731732return SVector(f1, f2, f3, f4)733end734735"""736splitting_steger_warming(u, orientation::Integer,737equations::CompressibleEulerEquations2D)738splitting_steger_warming(u, which::Union{Val{:minus}, Val{:plus}}739orientation::Integer,740equations::CompressibleEulerEquations2D)741742Splitting of the compressible Euler flux of Steger and Warming. For743curvilinear coordinates use the improved Steger-Warming-type splitting744[`splitting_drikakis_tsangaris`](@ref).745746Returns a tuple of the fluxes "minus" (associated with waves going into the747negative axis direction) and "plus" (associated with waves going into the748positive axis direction). If only one of the fluxes is required, use the749function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.750751!!! warning "Experimental implementation (upwind SBP)"752This is an experimental feature and may change in future releases.753754## References755756- Joseph L. Steger and R. F. Warming (1979)757Flux Vector Splitting of the Inviscid Gasdynamic Equations758With Application to Finite Difference Methods759[NASA Technical Memorandum](https://ntrs.nasa.gov/api/citations/19790020779/downloads/19790020779.pdf)760"""761@inline function splitting_steger_warming(u, orientation::Integer,762equations::CompressibleEulerEquations2D)763fm = splitting_steger_warming(u, Val{:minus}(), orientation, equations)764fp = splitting_steger_warming(u, Val{:plus}(), orientation, equations)765return fm, fp766end767768@inline function splitting_steger_warming(u, ::Val{:plus}, orientation::Integer,769equations::CompressibleEulerEquations2D)770rho, rho_v1, rho_v2, rho_e_total = u771v1 = rho_v1 / rho772v2 = rho_v2 / rho773p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))774a = sqrt(equations.gamma * p / rho)775776if orientation == 1777lambda1 = v1778lambda2 = v1 + a779lambda3 = v1 - a780781lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)782lambda2_p = positive_part(lambda2)783lambda3_p = positive_part(lambda3)784785alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p786787rho_2gamma = 0.5f0 * rho / equations.gamma788f1p = rho_2gamma * alpha_p789f2p = rho_2gamma * (alpha_p * v1 + a * (lambda2_p - lambda3_p))790f3p = rho_2gamma * alpha_p * v2791f4p = rho_2gamma *792(alpha_p * 0.5f0 * (v1^2 + v2^2) + a * v1 * (lambda2_p - lambda3_p)793+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)794else # orientation == 2795lambda1 = v2796lambda2 = v2 + a797lambda3 = v2 - a798799lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)800lambda2_p = positive_part(lambda2)801lambda3_p = positive_part(lambda3)802803alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p804805rho_2gamma = 0.5f0 * rho / equations.gamma806f1p = rho_2gamma * alpha_p807f2p = rho_2gamma * alpha_p * v1808f3p = rho_2gamma * (alpha_p * v2 + a * (lambda2_p - lambda3_p))809f4p = rho_2gamma *810(alpha_p * 0.5f0 * (v1^2 + v2^2) + a * v2 * (lambda2_p - lambda3_p)811+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)812end813return SVector(f1p, f2p, f3p, f4p)814end815816@inline function splitting_steger_warming(u, ::Val{:minus}, orientation::Integer,817equations::CompressibleEulerEquations2D)818rho, rho_v1, rho_v2, rho_e_total = u819v1 = rho_v1 / rho820v2 = rho_v2 / rho821p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))822a = sqrt(equations.gamma * p / rho)823824if orientation == 1825lambda1 = v1826lambda2 = v1 + a827lambda3 = v1 - a828829lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)830lambda2_m = negative_part(lambda2)831lambda3_m = negative_part(lambda3)832833alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m834835rho_2gamma = 0.5f0 * rho / equations.gamma836f1m = rho_2gamma * alpha_m837f2m = rho_2gamma * (alpha_m * v1 + a * (lambda2_m - lambda3_m))838f3m = rho_2gamma * alpha_m * v2839f4m = rho_2gamma *840(alpha_m * 0.5f0 * (v1^2 + v2^2) + a * v1 * (lambda2_m - lambda3_m)841+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)842else # orientation == 2843lambda1 = v2844lambda2 = v2 + a845lambda3 = v2 - a846847lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)848lambda2_m = negative_part(lambda2)849lambda3_m = negative_part(lambda3)850851alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m852853rho_2gamma = 0.5f0 * rho / equations.gamma854f1m = rho_2gamma * alpha_m855f2m = rho_2gamma * alpha_m * v1856f3m = rho_2gamma * (alpha_m * v2 + a * (lambda2_m - lambda3_m))857f4m = rho_2gamma *858(alpha_m * 0.5f0 * (v1^2 + v2^2) + a * v2 * (lambda2_m - lambda3_m)859+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)860end861return SVector(f1m, f2m, f3m, f4m)862end863864"""865splitting_drikakis_tsangaris(u, orientation_or_normal_direction,866equations::CompressibleEulerEquations2D)867splitting_drikakis_tsangaris(u, which::Union{Val{:minus}, Val{:plus}}868orientation_or_normal_direction,869equations::CompressibleEulerEquations2D)870871Improved variant of the Steger-Warming flux vector splitting872[`splitting_steger_warming`](@ref) for generalized coordinates.873This splitting also reformulates the energy874flux as in Hänel et al. to obtain conservation of the total temperature875for inviscid flows.876877Returns a tuple of the fluxes "minus" (associated with waves going into the878negative axis direction) and "plus" (associated with waves going into the879positive axis direction). If only one of the fluxes is required, use the880function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.881882!!! warning "Experimental implementation (upwind SBP)"883This is an experimental feature and may change in future releases.884885## References886887- D. Drikakis and S. Tsangaris (1993)888On the solution of the compressible Navier-Stokes equations using889improved flux vector splitting methods890[DOI: 10.1016/0307-904X(93)90054-K](https://doi.org/10.1016/0307-904X(93)90054-K)891- D. Hänel, R. Schwane and G. Seider (1987)892On the accuracy of upwind schemes for the solution of the Navier-Stokes equations893[DOI: 10.2514/6.1987-1105](https://doi.org/10.2514/6.1987-1105)894"""895@inline function splitting_drikakis_tsangaris(u, orientation_or_normal_direction,896equations::CompressibleEulerEquations2D)897fm = splitting_drikakis_tsangaris(u, Val{:minus}(), orientation_or_normal_direction,898equations)899fp = splitting_drikakis_tsangaris(u, Val{:plus}(), orientation_or_normal_direction,900equations)901return fm, fp902end903904@inline function splitting_drikakis_tsangaris(u, ::Val{:plus}, orientation::Integer,905equations::CompressibleEulerEquations2D)906rho, rho_v1, rho_v2, rho_e_total = u907v1 = rho_v1 / rho908v2 = rho_v2 / rho909p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))910a = sqrt(equations.gamma * p / rho)911H = (rho_e_total + p) / rho912913if orientation == 1914lambda1 = v1 + a915lambda2 = v1 - a916917lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)918lambda2_p = positive_part(lambda2)919920rhoa_2gamma = 0.5f0 * rho * a / equations.gamma921f1p = 0.5f0 * rho * (lambda1_p + lambda2_p)922f2p = f1p * v1 + rhoa_2gamma * (lambda1_p - lambda2_p)923f3p = f1p * v2924f4p = f1p * H925else # orientation == 2926lambda1 = v2 + a927lambda2 = v2 - a928929lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)930lambda2_p = positive_part(lambda2)931932rhoa_2gamma = 0.5f0 * rho * a / equations.gamma933f1p = 0.5f0 * rho * (lambda1_p + lambda2_p)934f2p = f1p * v1935f3p = f1p * v2 + rhoa_2gamma * (lambda1_p - lambda2_p)936f4p = f1p * H937end938return SVector(f1p, f2p, f3p, f4p)939end940941@inline function splitting_drikakis_tsangaris(u, ::Val{:minus}, orientation::Integer,942equations::CompressibleEulerEquations2D)943rho, rho_v1, rho_v2, rho_e_total = u944v1 = rho_v1 / rho945v2 = rho_v2 / rho946p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))947a = sqrt(equations.gamma * p / rho)948H = (rho_e_total + p) / rho949950if orientation == 1951lambda1 = v1 + a952lambda2 = v1 - a953954lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)955lambda2_m = negative_part(lambda2)956957rhoa_2gamma = 0.5f0 * rho * a / equations.gamma958f1m = 0.5f0 * rho * (lambda1_m + lambda2_m)959f2m = f1m * v1 + rhoa_2gamma * (lambda1_m - lambda2_m)960f3m = f1m * v2961f4m = f1m * H962else # orientation == 2963lambda1 = v2 + a964lambda2 = v2 - a965966lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)967lambda2_m = negative_part(lambda2)968969rhoa_2gamma = 0.5f0 * rho * a / equations.gamma970f1m = 0.5f0 * rho * (lambda1_m + lambda2_m)971f2m = f1m * v1972f3m = f1m * v2 + rhoa_2gamma * (lambda1_m - lambda2_m)973f4m = f1m * H974end975return SVector(f1m, f2m, f3m, f4m)976end977978@inline function splitting_drikakis_tsangaris(u, ::Val{:plus},979normal_direction::AbstractVector,980equations::CompressibleEulerEquations2D)981rho, rho_v1, rho_v2, rho_e_total = u982v1 = rho_v1 / rho983v2 = rho_v2 / rho984p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))985a = sqrt(equations.gamma * p / rho)986H = (rho_e_total + p) / rho987988v_n = normal_direction[1] * v1 + normal_direction[2] * v2989990lambda1 = v_n + a991lambda2 = v_n - a992993lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)994lambda2_p = positive_part(lambda2)995996rhoa_2gamma = 0.5f0 * rho * a / equations.gamma997f1p = 0.5f0 * rho * (lambda1_p + lambda2_p)998f2p = f1p * v1 + rhoa_2gamma * normal_direction[1] * (lambda1_p - lambda2_p)999f3p = f1p * v2 + rhoa_2gamma * normal_direction[2] * (lambda1_p - lambda2_p)1000f4p = f1p * H10011002return SVector(f1p, f2p, f3p, f4p)1003end10041005@inline function splitting_drikakis_tsangaris(u, ::Val{:minus},1006normal_direction::AbstractVector,1007equations::CompressibleEulerEquations2D)1008rho, rho_v1, rho_v2, rho_e_total = u1009v1 = rho_v1 / rho1010v2 = rho_v2 / rho1011p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))1012a = sqrt(equations.gamma * p / rho)1013H = (rho_e_total + p) / rho10141015v_n = normal_direction[1] * v1 + normal_direction[2] * v210161017lambda1 = v_n + a1018lambda2 = v_n - a10191020lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)1021lambda2_m = negative_part(lambda2)10221023rhoa_2gamma = 0.5f0 * rho * a / equations.gamma1024f1m = 0.5f0 * rho * (lambda1_m + lambda2_m)1025f2m = f1m * v1 + rhoa_2gamma * normal_direction[1] * (lambda1_m - lambda2_m)1026f3m = f1m * v2 + rhoa_2gamma * normal_direction[2] * (lambda1_m - lambda2_m)1027f4m = f1m * H10281029return SVector(f1m, f2m, f3m, f4m)1030end10311032"""1033FluxLMARS(c)(u_ll, u_rr, orientation_or_normal_direction,1034equations::CompressibleEulerEquations2D)10351036Low Mach number approximate Riemann solver (LMARS) for atmospheric flows using1037an estimate `c` of the speed of sound.10381039References:1040- Xi Chen et al. (2013)1041A Control-Volume Model of the Compressible Euler Equations with a Vertical1042Lagrangian Coordinate1043[DOI: 10.1175/MWR-D-12-00129.1](https://doi.org/10.1175/mwr-d-12-00129.1)1044"""1045struct FluxLMARS{SpeedOfSound}1046# Estimate for the speed of sound1047speed_of_sound::SpeedOfSound1048end10491050@inline function (flux_lmars::FluxLMARS)(u_ll, u_rr, orientation::Integer,1051equations::CompressibleEulerEquations2D)1052c = flux_lmars.speed_of_sound10531054# Unpack left and right state1055rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1056rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)10571058if orientation == 11059v_ll = v1_ll1060v_rr = v1_rr1061else # orientation == 21062v_ll = v2_ll1063v_rr = v2_rr1064end10651066rho = 0.5f0 * (rho_ll + rho_rr)1067p = 0.5f0 * (p_ll + p_rr) - 0.5f0 * c * rho * (v_rr - v_ll)1068v = 0.5f0 * (v_ll + v_rr) - 1 / (2 * c * rho) * (p_rr - p_ll)10691070# We treat the energy term analogous to the potential temperature term in the paper by1071# Chen et al., i.e. we use p_ll and p_rr, and not p1072if v >= 01073f1, f2, f3, f4 = v * u_ll1074f4 = f4 + p_ll * v1075else1076f1, f2, f3, f4 = v * u_rr1077f4 = f4 + p_rr * v1078end10791080if orientation == 11081f2 = f2 + p1082else # orientation == 21083f3 = f3 + p1084end10851086return SVector(f1, f2, f3, f4)1087end10881089@inline function (flux_lmars::FluxLMARS)(u_ll, u_rr, normal_direction::AbstractVector,1090equations::CompressibleEulerEquations2D)1091c = flux_lmars.speed_of_sound10921093# Unpack left and right state1094rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1095rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)10961097v_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]1098v_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]10991100# Note that this is the same as computing v_ll and v_rr with a normalized normal vector1101# and then multiplying v by `norm_` again, but this version is slightly faster.1102norm_ = norm(normal_direction)11031104rho = 0.5f0 * (rho_ll + rho_rr)1105p = 0.5f0 * (p_ll + p_rr) - 0.5f0 * c * rho * (v_rr - v_ll) / norm_1106v = 0.5f0 * (v_ll + v_rr) - 1 / (2 * c * rho) * (p_rr - p_ll) * norm_11071108# We treat the energy term analogous to the potential temperature term in the paper by1109# Chen et al., i.e. we use p_ll and p_rr, and not p1110if v >= 01111f1, f2, f3, f4 = u_ll * v1112f4 = f4 + p_ll * v1113else1114f1, f2, f3, f4 = u_rr * v1115f4 = f4 + p_rr * v1116end11171118return SVector(f1,1119f2 + p * normal_direction[1],1120f3 + p * normal_direction[2],1121f4)1122end11231124"""1125splitting_vanleer_haenel(u, orientation_or_normal_direction,1126equations::CompressibleEulerEquations2D)1127splitting_vanleer_haenel(u, which::Union{Val{:minus}, Val{:plus}}1128orientation_or_normal_direction,1129equations::CompressibleEulerEquations2D)11301131Splitting of the compressible Euler flux from van Leer. This splitting further1132contains a reformulation due to Hänel et al. where the energy flux uses the1133enthalpy. The pressure splitting is independent from the splitting of the1134convective terms. As such there are many pressure splittings suggested across1135the literature. We implement the 'p4' variant suggested by Liou and Steffen as1136it proved the most robust in practice. For details on the curvilinear variant1137of this flux vector splitting see Anderson et al.11381139Returns a tuple of the fluxes "minus" (associated with waves going into the1140negative axis direction) and "plus" (associated with waves going into the1141positive axis direction). If only one of the fluxes is required, use the1142function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.11431144!!! warning "Experimental implementation (upwind SBP)"1145This is an experimental feature and may change in future releases.11461147## References11481149- Bram van Leer (1982)1150Flux-Vector Splitting for the Euler Equation1151[DOI: 10.1007/978-3-642-60543-7_5](https://doi.org/10.1007/978-3-642-60543-7_5)1152- D. Hänel, R. Schwane and G. Seider (1987)1153On the accuracy of upwind schemes for the solution of the Navier-Stokes equations1154[DOI: 10.2514/6.1987-1105](https://doi.org/10.2514/6.1987-1105)1155- Meng-Sing Liou and Chris J. Steffen, Jr. (1991)1156High-Order Polynomial Expansions (HOPE) for Flux-Vector Splitting1157[NASA Technical Memorandum](https://ntrs.nasa.gov/citations/19910016425)1158- W. Kyle Anderson, James L. Thomas, and Bram van Leer (1986)1159Comparison of Finite Volume Flux Vector Splittings for the Euler Equations1160[DOI: 10.2514/3.9465](https://doi.org/10.2514/3.9465)1161"""1162@inline function splitting_vanleer_haenel(u, orientation_or_normal_direction,1163equations::CompressibleEulerEquations2D)1164fm = splitting_vanleer_haenel(u, Val{:minus}(), orientation_or_normal_direction,1165equations)1166fp = splitting_vanleer_haenel(u, Val{:plus}(), orientation_or_normal_direction,1167equations)1168return fm, fp1169end11701171@inline function splitting_vanleer_haenel(u, ::Val{:plus}, orientation::Integer,1172equations::CompressibleEulerEquations2D)1173rho, rho_v1, rho_v2, rho_e_total = u1174v1 = rho_v1 / rho1175v2 = rho_v2 / rho1176p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))11771178a = sqrt(equations.gamma * p / rho)1179H = (rho_e_total + p) / rho11801181if orientation == 11182M = v1 / a1183p_plus = 0.5f0 * (1 + equations.gamma * M) * p11841185f1p = 0.25f0 * rho * a * (M + 1)^21186f2p = f1p * v1 + p_plus1187f3p = f1p * v21188f4p = f1p * H1189else # orientation == 21190M = v2 / a1191p_plus = 0.5f0 * (1 + equations.gamma * M) * p11921193f1p = 0.25f0 * rho * a * (M + 1)^21194f2p = f1p * v11195f3p = f1p * v2 + p_plus1196f4p = f1p * H1197end1198return SVector(f1p, f2p, f3p, f4p)1199end12001201@inline function splitting_vanleer_haenel(u, ::Val{:minus}, orientation::Integer,1202equations::CompressibleEulerEquations2D)1203rho, rho_v1, rho_v2, rho_e_total = u1204v1 = rho_v1 / rho1205v2 = rho_v2 / rho1206p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))12071208a = sqrt(equations.gamma * p / rho)1209H = (rho_e_total + p) / rho12101211if orientation == 11212M = v1 / a1213p_minus = 0.5f0 * (1 - equations.gamma * M) * p12141215f1m = -0.25f0 * rho * a * (M - 1)^21216f2m = f1m * v1 + p_minus1217f3m = f1m * v21218f4m = f1m * H1219else # orientation == 21220M = v2 / a1221p_minus = 0.5f0 * (1 - equations.gamma * M) * p12221223f1m = -0.25f0 * rho * a * (M - 1)^21224f2m = f1m * v11225f3m = f1m * v2 + p_minus1226f4m = f1m * H1227end1228return SVector(f1m, f2m, f3m, f4m)1229end12301231@inline function splitting_vanleer_haenel(u, ::Val{:plus},1232normal_direction::AbstractVector,1233equations::CompressibleEulerEquations2D)1234rho, rho_v1, rho_v2, rho_e_total = u1235v1 = rho_v1 / rho1236v2 = rho_v2 / rho1237p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))12381239a = sqrt(equations.gamma * p / rho)1240H = (rho_e_total + p) / rho12411242v_n = normal_direction[1] * v1 + normal_direction[2] * v21243M = v_n / a1244p_plus = 0.5f0 * (1 + equations.gamma * M) * p12451246f1p = 0.25f0 * rho * a * (M + 1)^21247f2p = f1p * v1 + normal_direction[1] * p_plus1248f3p = f1p * v2 + normal_direction[2] * p_plus1249f4p = f1p * H12501251return SVector(f1p, f2p, f3p, f4p)1252end12531254@inline function splitting_vanleer_haenel(u, ::Val{:minus},1255normal_direction::AbstractVector,1256equations::CompressibleEulerEquations2D)1257rho, rho_v1, rho_v2, rho_e_total = u1258v1 = rho_v1 / rho1259v2 = rho_v2 / rho1260p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))12611262a = sqrt(equations.gamma * p / rho)1263H = (rho_e_total + p) / rho12641265v_n = normal_direction[1] * v1 + normal_direction[2] * v21266M = v_n / a1267p_minus = 0.5f0 * (1 - equations.gamma * M) * p12681269f1m = -0.25f0 * rho * a * (M - 1)^21270f2m = f1m * v1 + normal_direction[1] * p_minus1271f3m = f1m * v2 + normal_direction[2] * p_minus1272f4m = f1m * H12731274return SVector(f1m, f2m, f3m, f4m)1275end12761277"""1278splitting_lax_friedrichs(u, orientation_or_normal_direction,1279equations::CompressibleEulerEquations2D)1280splitting_lax_friedrichs(u, which::Union{Val{:minus}, Val{:plus}}1281orientation_or_normal_direction,1282equations::CompressibleEulerEquations2D)12831284Naive local Lax-Friedrichs style flux splitting of the form `f⁺ = 0.5 (f + λ u)`1285and `f⁻ = 0.5 (f - λ u)` similar to a flux splitting one would apply, e.g.,1286to Burgers' equation.12871288Returns a tuple of the fluxes "minus" (associated with waves going into the1289negative axis direction) and "plus" (associated with waves going into the1290positive axis direction). If only one of the fluxes is required, use the1291function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.12921293!!! warning "Experimental implementation (upwind SBP)"1294This is an experimental feature and may change in future releases.1295"""1296@inline function splitting_lax_friedrichs(u, orientation_or_normal_direction,1297equations::CompressibleEulerEquations2D)1298fm = splitting_lax_friedrichs(u, Val{:minus}(), orientation_or_normal_direction,1299equations)1300fp = splitting_lax_friedrichs(u, Val{:plus}(), orientation_or_normal_direction,1301equations)1302return fm, fp1303end13041305@inline function splitting_lax_friedrichs(u, ::Val{:plus}, orientation::Integer,1306equations::CompressibleEulerEquations2D)1307rho, rho_v1, rho_v2, rho_e_total = u1308v1 = rho_v1 / rho1309v2 = rho_v2 / rho1310p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))13111312a = sqrt(equations.gamma * p / rho)1313H = (rho_e_total + p) / rho1314lambda = 0.5f0 * (sqrt(v1^2 + v2^2) + a)13151316if orientation == 11317#lambda = 0.5 * (abs(v1) + a)1318f1p = 0.5f0 * rho * v1 + lambda * u[1]1319f2p = 0.5f0 * rho * v1 * v1 + 0.5f0 * p + lambda * u[2]1320f3p = 0.5f0 * rho * v1 * v2 + lambda * u[3]1321f4p = 0.5f0 * rho * v1 * H + lambda * u[4]1322else # orientation == 21323#lambda = 0.5 * (abs(v2) + a)1324f1p = 0.5f0 * rho * v2 + lambda * u[1]1325f2p = 0.5f0 * rho * v2 * v1 + lambda * u[2]1326f3p = 0.5f0 * rho * v2 * v2 + 0.5f0 * p + lambda * u[3]1327f4p = 0.5f0 * rho * v2 * H + lambda * u[4]1328end1329return SVector(f1p, f2p, f3p, f4p)1330end13311332@inline function splitting_lax_friedrichs(u, ::Val{:minus}, orientation::Integer,1333equations::CompressibleEulerEquations2D)1334rho, rho_v1, rho_v2, rho_e_total = u1335v1 = rho_v1 / rho1336v2 = rho_v2 / rho1337p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))13381339a = sqrt(equations.gamma * p / rho)1340H = (rho_e_total + p) / rho1341lambda = 0.5f0 * (sqrt(v1^2 + v2^2) + a)13421343if orientation == 11344#lambda = 0.5 * (abs(v1) + a)1345f1m = 0.5f0 * rho * v1 - lambda * u[1]1346f2m = 0.5f0 * rho * v1 * v1 + 0.5f0 * p - lambda * u[2]1347f3m = 0.5f0 * rho * v1 * v2 - lambda * u[3]1348f4m = 0.5f0 * rho * v1 * H - lambda * u[4]1349else # orientation == 21350#lambda = 0.5 * (abs(v2) + a)1351f1m = 0.5f0 * rho * v2 - lambda * u[1]1352f2m = 0.5f0 * rho * v2 * v1 - lambda * u[2]1353f3m = 0.5f0 * rho * v2 * v2 + 0.5f0 * p - lambda * u[3]1354f4m = 0.5f0 * rho * v2 * H - lambda * u[4]1355end1356return SVector(f1m, f2m, f3m, f4m)1357end13581359@inline function splitting_lax_friedrichs(u, ::Val{:plus},1360normal_direction::AbstractVector,1361equations::CompressibleEulerEquations2D)1362rho_e_total = last(u)1363rho, v1, v2, p = cons2prim(u, equations)13641365a = sqrt(equations.gamma * p / rho)1366H = (rho_e_total + p) / rho1367lambda = 0.5f0 * (sqrt(v1^2 + v2^2) + a)13681369v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]1370rho_v_normal = rho * v_normal13711372f1p = 0.5f0 * rho_v_normal + lambda * u[1]1373f2p = 0.5f0 * rho_v_normal * v1 + 0.5f0 * p * normal_direction[1] + lambda * u[2]1374f3p = 0.5f0 * rho_v_normal * v2 + 0.5f0 * p * normal_direction[2] + lambda * u[3]1375f4p = 0.5f0 * rho_v_normal * H + lambda * u[4]13761377return SVector(f1p, f2p, f3p, f4p)1378end13791380@inline function splitting_lax_friedrichs(u, ::Val{:minus},1381normal_direction::AbstractVector,1382equations::CompressibleEulerEquations2D)1383rho_e_total = last(u)1384rho, v1, v2, p = cons2prim(u, equations)13851386a = sqrt(equations.gamma * p / rho)1387H = (rho_e_total + p) / rho1388lambda = 0.5f0 * (sqrt(v1^2 + v2^2) + a)13891390v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]1391rho_v_normal = rho * v_normal13921393f1m = 0.5f0 * rho_v_normal - lambda * u[1]1394f2m = 0.5f0 * rho_v_normal * v1 + 0.5f0 * p * normal_direction[1] - lambda * u[2]1395f3m = 0.5f0 * rho_v_normal * v2 + 0.5f0 * p * normal_direction[2] - lambda * u[3]1396f4m = 0.5f0 * rho_v_normal * H - lambda * u[4]13971398return SVector(f1m, f2m, f3m, f4m)1399end14001401# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation as the1402# maximum velocity magnitude plus the maximum speed of sound1403@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,1404equations::CompressibleEulerEquations2D)1405rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1406rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)14071408# Get the velocity value in the appropriate direction1409if orientation == 11410v_ll = v1_ll1411v_rr = v1_rr1412else # orientation == 21413v_ll = v2_ll1414v_rr = v2_rr1415end1416# Calculate sound speeds1417c_ll = sqrt(equations.gamma * p_ll / rho_ll)1418c_rr = sqrt(equations.gamma * p_rr / rho_rr)14191420return max(abs(v_ll), abs(v_rr)) + max(c_ll, c_rr)1421end14221423@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,1424equations::CompressibleEulerEquations2D)1425rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1426rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)14271428# Calculate normal velocities and sound speed1429# left1430v_ll = (v1_ll * normal_direction[1]1431+1432v2_ll * normal_direction[2])1433c_ll = sqrt(equations.gamma * p_ll / rho_ll)1434# right1435v_rr = (v1_rr * normal_direction[1]1436+1437v2_rr * normal_direction[2])1438c_rr = sqrt(equations.gamma * p_rr / rho_rr)14391440return max(abs(v_ll), abs(v_rr)) + max(c_ll, c_rr) * norm(normal_direction)1441end14421443# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`1444@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,1445equations::CompressibleEulerEquations2D)1446rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1447rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)14481449# Get the velocity value in the appropriate direction1450if orientation == 11451v_ll = v1_ll1452v_rr = v1_rr1453else # orientation == 21454v_ll = v2_ll1455v_rr = v2_rr1456end1457# Calculate sound speeds1458c_ll = sqrt(equations.gamma * p_ll / rho_ll)1459c_rr = sqrt(equations.gamma * p_rr / rho_rr)14601461return max(abs(v_ll) + c_ll, abs(v_rr) + c_rr)1462end14631464# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`1465@inline function max_abs_speed(u_ll, u_rr, normal_direction::AbstractVector,1466equations::CompressibleEulerEquations2D)1467rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1468rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)14691470# Calculate normal velocities and sound speeds1471# left1472v_ll = (v1_ll * normal_direction[1]1473+1474v2_ll * normal_direction[2])1475c_ll = sqrt(equations.gamma * p_ll / rho_ll)1476# right1477v_rr = (v1_rr * normal_direction[1]1478+1479v2_rr * normal_direction[2])1480c_rr = sqrt(equations.gamma * p_rr / rho_rr)14811482norm_ = norm(normal_direction)1483return max(abs(v_ll) + c_ll * norm_,1484abs(v_rr) + c_rr * norm_)1485end14861487# Calculate estimate for minimum and maximum wave speeds for HLL-type fluxes1488@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,1489equations::CompressibleEulerEquations2D)1490rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1491rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)14921493if orientation == 1 # x-direction1494λ_min = v1_ll - sqrt(equations.gamma * p_ll / rho_ll)1495λ_max = v1_rr + sqrt(equations.gamma * p_rr / rho_rr)1496else # y-direction1497λ_min = v2_ll - sqrt(equations.gamma * p_ll / rho_ll)1498λ_max = v2_rr + sqrt(equations.gamma * p_rr / rho_rr)1499end15001501return λ_min, λ_max1502end15031504@inline function min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,1505equations::CompressibleEulerEquations2D)1506rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1507rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)15081509v_normal_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]1510v_normal_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]15111512norm_ = norm(normal_direction)1513# The v_normals are already scaled by the norm1514λ_min = v_normal_ll - sqrt(equations.gamma * p_ll / rho_ll) * norm_1515λ_max = v_normal_rr + sqrt(equations.gamma * p_rr / rho_rr) * norm_15161517return λ_min, λ_max1518end15191520# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes1521@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,1522equations::CompressibleEulerEquations2D)1523rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1524rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)15251526c_ll = sqrt(equations.gamma * p_ll / rho_ll)1527c_rr = sqrt(equations.gamma * p_rr / rho_rr)15281529if orientation == 1 # x-direction1530λ_min = min(v1_ll - c_ll, v1_rr - c_rr)1531λ_max = max(v1_ll + c_ll, v1_rr + c_rr)1532else # y-direction1533λ_min = min(v2_ll - c_ll, v2_rr - c_rr)1534λ_max = max(v2_ll + c_ll, v2_rr + c_rr)1535end15361537return λ_min, λ_max1538end15391540# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes1541@inline function min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector,1542equations::CompressibleEulerEquations2D)1543rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1544rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)15451546norm_ = norm(normal_direction)15471548c_ll = sqrt(equations.gamma * p_ll / rho_ll) * norm_1549c_rr = sqrt(equations.gamma * p_rr / rho_rr) * norm_15501551v_normal_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]1552v_normal_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]15531554# The v_normals are already scaled by the norm1555λ_min = min(v_normal_ll - c_ll, v_normal_rr - c_rr)1556λ_max = max(v_normal_ll + c_ll, v_normal_rr + c_rr)15571558return λ_min, λ_max1559end15601561@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,1562normal_direction::AbstractVector,1563equations::CompressibleEulerEquations2D)1564(; gamma) = equations15651566norm_ = norm(normal_direction)1567unit_normal_direction = normal_direction / norm_15681569rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1570rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)15711572b_ll = rho_ll / (2 * p_ll)1573b_rr = rho_rr / (2 * p_rr)15741575rho_log = ln_mean(rho_ll, rho_rr)1576b_log = ln_mean(b_ll, b_rr)1577v1_avg = 0.5f0 * (v1_ll + v1_rr)1578v2_avg = 0.5f0 * (v2_ll + v2_rr)1579p_avg = 0.5f0 * (rho_ll + rho_rr) / (b_ll + b_rr) # 2 * b_avg = b_ll + b_rr1580v_squared_bar = v1_ll * v1_rr + v2_ll * v2_rr1581h_bar = gamma / (2 * b_log * (gamma - 1)) + 0.5f0 * v_squared_bar1582c_bar = sqrt(gamma * p_avg / rho_log)15831584v_avg_normal = dot(SVector(v1_avg, v2_avg), unit_normal_direction)15851586lambda_1 = abs(v_avg_normal - c_bar) * rho_log / (2 * gamma)1587lambda_2 = abs(v_avg_normal) * rho_log * (gamma - 1) / gamma1588lambda_3 = abs(v_avg_normal + c_bar) * rho_log / (2 * gamma)1589lambda_4 = abs(v_avg_normal) * p_avg15901591v1_minus_c = v1_avg - c_bar * unit_normal_direction[1]1592v2_minus_c = v2_avg - c_bar * unit_normal_direction[2]1593v1_plus_c = v1_avg + c_bar * unit_normal_direction[1]1594v2_plus_c = v2_avg + c_bar * unit_normal_direction[2]1595v1_tangential = v1_avg - v_avg_normal * unit_normal_direction[1]1596v2_tangential = v2_avg - v_avg_normal * unit_normal_direction[2]15971598entropy_vars_jump = cons2entropy(u_rr, equations) - cons2entropy(u_ll, equations)1599entropy_var_rho_jump, entropy_var_rho_v1_jump,1600entropy_var_rho_v2_jump, entropy_var_rho_e_total_jump = entropy_vars_jump16011602velocity_minus_c_dot_entropy_vars_jump = v1_minus_c * entropy_var_rho_v1_jump +1603v2_minus_c * entropy_var_rho_v2_jump1604velocity_plus_c_dot_entropy_vars_jump = v1_plus_c * entropy_var_rho_v1_jump +1605v2_plus_c * entropy_var_rho_v2_jump1606velocity_avg_dot_vjump = v1_avg * entropy_var_rho_v1_jump +1607v2_avg * entropy_var_rho_v2_jump1608w1 = lambda_1 * (entropy_var_rho_jump + velocity_minus_c_dot_entropy_vars_jump +1609(h_bar - c_bar * v_avg_normal) * entropy_var_rho_e_total_jump)1610w2 = lambda_2 * (entropy_var_rho_jump + velocity_avg_dot_vjump +1611v_squared_bar / 2 * entropy_var_rho_e_total_jump)1612w3 = lambda_3 * (entropy_var_rho_jump + velocity_plus_c_dot_entropy_vars_jump +1613(h_bar + c_bar * v_avg_normal) * entropy_var_rho_e_total_jump)16141615entropy_var_v_normal_jump = dot(SVector(entropy_var_rho_v1_jump,1616entropy_var_rho_v2_jump),1617unit_normal_direction)16181619dissipation_rho = w1 + w2 + w316201621dissipation_rho_v1 = (w1 * v1_minus_c +1622w2 * v1_avg +1623w3 * v1_plus_c +1624lambda_4 * (entropy_var_rho_v1_jump -1625unit_normal_direction[1] * entropy_var_v_normal_jump +1626entropy_var_rho_e_total_jump * v1_tangential))16271628dissipation_rho_v2 = (w1 * v2_minus_c +1629w2 * v2_avg +1630w3 * v2_plus_c +1631lambda_4 * (entropy_var_rho_v2_jump -1632unit_normal_direction[2] * entropy_var_v_normal_jump +1633entropy_var_rho_e_total_jump * v2_tangential))16341635v_tangential_dot_entropy_vars_jump = v1_tangential * entropy_var_rho_v1_jump +1636v2_tangential * entropy_var_rho_v2_jump16371638dissipation_rhoe = (w1 * (h_bar - c_bar * v_avg_normal) +1639w2 * 0.5f0 * v_squared_bar +1640w3 * (h_bar + c_bar * v_avg_normal) +1641lambda_4 * (v_tangential_dot_entropy_vars_jump +1642entropy_var_rho_e_total_jump *1643(v1_avg^2 + v2_avg^2 - v_avg_normal^2)))16441645return -0.5f0 *1646SVector(dissipation_rho, dissipation_rho_v1, dissipation_rho_v2,1647dissipation_rhoe) * norm_1648end16491650# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction1651# has been normalized prior to this rotation of the state vector1652@inline function rotate_to_x(u, normal_vector, equations::CompressibleEulerEquations2D)1653# cos and sin of the angle between the x-axis and the normalized normal_vector are1654# the normalized vector's x and y coordinates respectively (see unit circle).1655c = normal_vector[1]1656s = normal_vector[2]16571658# Apply the 2D rotation matrix with normal and tangent directions of the form1659# [ 1 0 0 0;1660# 0 n_1 n_2 0;1661# 0 t_1 t_2 0;1662# 0 0 0 1 ]1663# where t_1 = -n_2 and t_2 = n_116641665return SVector(u[1],1666c * u[2] + s * u[3],1667-s * u[2] + c * u[3],1668u[4])1669end16701671# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction1672# has been normalized prior to this back-rotation of the state vector1673@inline function rotate_from_x(u, normal_vector,1674equations::CompressibleEulerEquations2D)1675# cos and sin of the angle between the x-axis and the normalized normal_vector are1676# the normalized vector's x and y coordinates respectively (see unit circle).1677c = normal_vector[1]1678s = normal_vector[2]16791680# Apply the 2D back-rotation matrix with normal and tangent directions of the form1681# [ 1 0 0 0;1682# 0 n_1 t_1 0;1683# 0 n_2 t_2 0;1684# 0 0 0 1 ]1685# where t_1 = -n_2 and t_2 = n_116861687return SVector(u[1],1688c * u[2] - s * u[3],1689s * u[2] + c * u[3],1690u[4])1691end16921693"""1694flux_hllc(u_ll, u_rr, orientation_or_normal_direction, equations::CompressibleEulerEquations2D)16951696Computes the HLLC flux (HLL with Contact) for compressible Euler equations developed by E.F. Toro1697[Lecture slides](http://www.prague-sum.com/download/2012/Toro_2-HLLC-RiemannSolver.pdf)1698Signal speeds: [DOI: 10.1137/S1064827593260140](https://doi.org/10.1137/S1064827593260140)1699"""1700function flux_hllc(u_ll, u_rr, orientation::Integer,1701equations::CompressibleEulerEquations2D)1702# Calculate primitive variables and speed of sound1703rho_ll, rho_v1_ll, rho_v2_ll, rho_e_total_ll = u_ll1704rho_rr, rho_v1_rr, rho_v2_rr, rho_e_total_rr = u_rr17051706v1_ll = rho_v1_ll / rho_ll1707v2_ll = rho_v2_ll / rho_ll1708e_ll = rho_e_total_ll / rho_ll1709p_ll = (equations.gamma - 1) *1710(rho_e_total_ll - 0.5f0 * rho_ll * (v1_ll^2 + v2_ll^2))1711c_ll = sqrt(equations.gamma * p_ll / rho_ll)17121713v1_rr = rho_v1_rr / rho_rr1714v2_rr = rho_v2_rr / rho_rr1715e_rr = rho_e_total_rr / rho_rr1716p_rr = (equations.gamma - 1) *1717(rho_e_total_rr - 0.5f0 * rho_rr * (v1_rr^2 + v2_rr^2))1718c_rr = sqrt(equations.gamma * p_rr / rho_rr)17191720# Obtain left and right fluxes1721f_ll = flux(u_ll, orientation, equations)1722f_rr = flux(u_rr, orientation, equations)17231724# Compute Roe averages1725sqrt_rho_ll = sqrt(rho_ll)1726sqrt_rho_rr = sqrt(rho_rr)1727sum_sqrt_rho = sqrt_rho_ll + sqrt_rho_rr1728if orientation == 1 # x-direction1729vel_L = v1_ll1730vel_R = v1_rr1731elseif orientation == 2 # y-direction1732vel_L = v2_ll1733vel_R = v2_rr1734end1735vel_roe = (sqrt_rho_ll * vel_L + sqrt_rho_rr * vel_R) / sum_sqrt_rho1736v1_roe = sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr1737v2_roe = sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr1738vel_roe_mag = (v1_roe^2 + v2_roe^2) / sum_sqrt_rho^21739H_ll = (rho_e_total_ll + p_ll) / rho_ll1740H_rr = (rho_e_total_rr + p_rr) / rho_rr1741H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) / sum_sqrt_rho1742c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * vel_roe_mag))1743Ssl = min(vel_L - c_ll, vel_roe - c_roe)1744Ssr = max(vel_R + c_rr, vel_roe + c_roe)1745sMu_L = Ssl - vel_L1746sMu_R = Ssr - vel_R17471748if Ssl >= 01749f1 = f_ll[1]1750f2 = f_ll[2]1751f3 = f_ll[3]1752f4 = f_ll[4]1753elseif Ssr <= 01754f1 = f_rr[1]1755f2 = f_rr[2]1756f3 = f_rr[3]1757f4 = f_rr[4]1758else1759SStar = (p_rr - p_ll + rho_ll * vel_L * sMu_L - rho_rr * vel_R * sMu_R) /1760(rho_ll * sMu_L - rho_rr * sMu_R)1761if Ssl <= 0 <= SStar1762densStar = rho_ll * sMu_L / (Ssl - SStar)1763enerStar = e_ll + (SStar - vel_L) * (SStar + p_ll / (rho_ll * sMu_L))1764UStar1 = densStar1765UStar4 = densStar * enerStar1766if orientation == 1 # x-direction1767UStar2 = densStar * SStar1768UStar3 = densStar * v2_ll1769elseif orientation == 2 # y-direction1770UStar2 = densStar * v1_ll1771UStar3 = densStar * SStar1772end1773f1 = f_ll[1] + Ssl * (UStar1 - rho_ll)1774f2 = f_ll[2] + Ssl * (UStar2 - rho_v1_ll)1775f3 = f_ll[3] + Ssl * (UStar3 - rho_v2_ll)1776f4 = f_ll[4] + Ssl * (UStar4 - rho_e_total_ll)1777else1778densStar = rho_rr * sMu_R / (Ssr - SStar)1779enerStar = e_rr + (SStar - vel_R) * (SStar + p_rr / (rho_rr * sMu_R))1780UStar1 = densStar1781UStar4 = densStar * enerStar1782if orientation == 1 # x-direction1783UStar2 = densStar * SStar1784UStar3 = densStar * v2_rr1785elseif orientation == 2 # y-direction1786UStar2 = densStar * v1_rr1787UStar3 = densStar * SStar1788end1789f1 = f_rr[1] + Ssr * (UStar1 - rho_rr)1790f2 = f_rr[2] + Ssr * (UStar2 - rho_v1_rr)1791f3 = f_rr[3] + Ssr * (UStar3 - rho_v2_rr)1792f4 = f_rr[4] + Ssr * (UStar4 - rho_e_total_rr)1793end1794end1795return SVector(f1, f2, f3, f4)1796end17971798function flux_hllc(u_ll, u_rr, normal_direction::AbstractVector,1799equations::CompressibleEulerEquations2D)1800# Calculate primitive variables and speed of sound1801rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1802rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)18031804v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]1805v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]18061807norm_ = norm(normal_direction)1808norm_sq = norm_ * norm_1809inv_norm_sq = inv(norm_sq)18101811c_ll = sqrt(equations.gamma * p_ll / rho_ll) * norm_1812c_rr = sqrt(equations.gamma * p_rr / rho_rr) * norm_18131814# Obtain left and right fluxes1815f_ll = flux(u_ll, normal_direction, equations)1816f_rr = flux(u_rr, normal_direction, equations)18171818# Compute Roe averages1819sqrt_rho_ll = sqrt(rho_ll)1820sqrt_rho_rr = sqrt(rho_rr)1821sum_sqrt_rho = sqrt_rho_ll + sqrt_rho_rr18221823v1_roe = (sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr) / sum_sqrt_rho1824v2_roe = (sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr) / sum_sqrt_rho1825vel_roe = v1_roe * normal_direction[1] + v2_roe * normal_direction[2]1826vel_roe_mag = v1_roe^2 + v2_roe^218271828e_ll = u_ll[4] / rho_ll1829e_rr = u_rr[4] / rho_rr18301831H_ll = (u_ll[4] + p_ll) / rho_ll1832H_rr = (u_rr[4] + p_rr) / rho_rr18331834H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) / sum_sqrt_rho1835c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * vel_roe_mag)) * norm_18361837Ssl = min(v_dot_n_ll - c_ll, vel_roe - c_roe)1838Ssr = max(v_dot_n_rr + c_rr, vel_roe + c_roe)1839sMu_L = Ssl - v_dot_n_ll1840sMu_R = Ssr - v_dot_n_rr18411842if Ssl >= 01843f1 = f_ll[1]1844f2 = f_ll[2]1845f3 = f_ll[3]1846f4 = f_ll[4]1847elseif Ssr <= 01848f1 = f_rr[1]1849f2 = f_rr[2]1850f3 = f_rr[3]1851f4 = f_rr[4]1852else1853SStar = (rho_ll * v_dot_n_ll * sMu_L - rho_rr * v_dot_n_rr * sMu_R +1854(p_rr - p_ll) * norm_sq) / (rho_ll * sMu_L - rho_rr * sMu_R)1855if Ssl <= 0 <= SStar1856densStar = rho_ll * sMu_L / (Ssl - SStar)1857enerStar = e_ll +1858(SStar - v_dot_n_ll) *1859(SStar * inv_norm_sq + p_ll / (rho_ll * sMu_L))1860UStar1 = densStar1861UStar2 = densStar *1862(v1_ll + (SStar - v_dot_n_ll) * normal_direction[1] * inv_norm_sq)1863UStar3 = densStar *1864(v2_ll + (SStar - v_dot_n_ll) * normal_direction[2] * inv_norm_sq)1865UStar4 = densStar * enerStar1866f1 = f_ll[1] + Ssl * (UStar1 - u_ll[1])1867f2 = f_ll[2] + Ssl * (UStar2 - u_ll[2])1868f3 = f_ll[3] + Ssl * (UStar3 - u_ll[3])1869f4 = f_ll[4] + Ssl * (UStar4 - u_ll[4])1870else1871densStar = rho_rr * sMu_R / (Ssr - SStar)1872enerStar = e_rr +1873(SStar - v_dot_n_rr) *1874(SStar * inv_norm_sq + p_rr / (rho_rr * sMu_R))1875UStar1 = densStar1876UStar2 = densStar *1877(v1_rr + (SStar - v_dot_n_rr) * normal_direction[1] * inv_norm_sq)1878UStar3 = densStar *1879(v2_rr + (SStar - v_dot_n_rr) * normal_direction[2] * inv_norm_sq)1880UStar4 = densStar * enerStar1881f1 = f_rr[1] + Ssr * (UStar1 - u_rr[1])1882f2 = f_rr[2] + Ssr * (UStar2 - u_rr[2])1883f3 = f_rr[3] + Ssr * (UStar3 - u_rr[3])1884f4 = f_rr[4] + Ssr * (UStar4 - u_rr[4])1885end1886end1887return SVector(f1, f2, f3, f4)1888end18891890"""1891min_max_speed_einfeldt(u_ll, u_rr, orientation, equations::CompressibleEulerEquations2D)18921893Computes the HLLE (Harten-Lax-van Leer-Einfeldt) flux for the compressible Euler equations.1894Special estimates of the signal velocites and linearization of the Riemann problem developed1895by Einfeldt to ensure that the internal energy and density remain positive during the computation1896of the numerical flux.18971898- Bernd Einfeldt (1988)1899On Godunov-type methods for gas dynamics.1900[DOI: 10.1137/0725021](https://doi.org/10.1137/0725021)1901- Bernd Einfeldt, Claus-Dieter Munz, Philip L. Roe and Björn Sjögreen (1991)1902On Godunov-type methods near low densities.1903[DOI: 10.1016/0021-9991(91)90211-3](https://doi.org/10.1016/0021-9991(91)90211-3)1904"""1905@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,1906equations::CompressibleEulerEquations2D)1907# Calculate primitive variables, enthalpy and speed of sound1908rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1909rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)19101911# `u_ll[4]` is total energy `rho_e_total_ll` on the left1912H_ll = (u_ll[4] + p_ll) / rho_ll1913c_ll = sqrt(equations.gamma * p_ll / rho_ll)19141915# `u_rr[4]` is total energy `rho_e_total_rr` on the right1916H_rr = (u_rr[4] + p_rr) / rho_rr1917c_rr = sqrt(equations.gamma * p_rr / rho_rr)19181919# Compute Roe averages1920sqrt_rho_ll = sqrt(rho_ll)1921sqrt_rho_rr = sqrt(rho_rr)1922inv_sum_sqrt_rho = inv(sqrt_rho_ll + sqrt_rho_rr)19231924v1_roe = (sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr) * inv_sum_sqrt_rho1925v2_roe = (sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr) * inv_sum_sqrt_rho1926v_roe_mag = v1_roe^2 + v2_roe^219271928H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) * inv_sum_sqrt_rho1929c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * v_roe_mag))19301931# Compute convenience constant for positivity preservation, see1932# https://doi.org/10.1016/0021-9991(91)90211-31933beta = sqrt(0.5f0 * (equations.gamma - 1) / equations.gamma)19341935# Estimate the edges of the Riemann fan (with positivity conservation)1936if orientation == 1 # x-direction1937SsL = min(v1_roe - c_roe, v1_ll - beta * c_ll, 0)1938SsR = max(v1_roe + c_roe, v1_rr + beta * c_rr, 0)1939elseif orientation == 2 # y-direction1940SsL = min(v2_roe - c_roe, v2_ll - beta * c_ll, 0)1941SsR = max(v2_roe + c_roe, v2_rr + beta * c_rr, 0)1942end19431944return SsL, SsR1945end19461947"""1948min_max_speed_einfeldt(u_ll, u_rr, normal_direction, equations::CompressibleEulerEquations2D)19491950Computes the HLLE (Harten-Lax-van Leer-Einfeldt) flux for the compressible Euler equations.1951Special estimates of the signal velocites and linearization of the Riemann problem developed1952by Einfeldt to ensure that the internal energy and density remain positive during the computation1953of the numerical flux.19541955- Bernd Einfeldt (1988)1956On Godunov-type methods for gas dynamics.1957[DOI: 10.1137/0725021](https://doi.org/10.1137/0725021)1958- Bernd Einfeldt, Claus-Dieter Munz, Philip L. Roe and Björn Sjögreen (1991)1959On Godunov-type methods near low densities.1960[DOI: 10.1016/0021-9991(91)90211-3](https://doi.org/10.1016/0021-9991(91)90211-3)1961"""1962@inline function min_max_speed_einfeldt(u_ll, u_rr, normal_direction::AbstractVector,1963equations::CompressibleEulerEquations2D)1964# Calculate primitive variables, enthalpy and speed of sound1965rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)1966rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)19671968v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]1969v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]19701971norm_ = norm(normal_direction)19721973# `u_ll[4]` is total energy `rho_e_total_ll` on the left1974H_ll = (u_ll[4] + p_ll) / rho_ll1975c_ll = sqrt(equations.gamma * p_ll / rho_ll) * norm_19761977# `u_rr[4]` is total energy `rho_e_total_rr` on the right1978H_rr = (u_rr[4] + p_rr) / rho_rr1979c_rr = sqrt(equations.gamma * p_rr / rho_rr) * norm_19801981# Compute Roe averages1982sqrt_rho_ll = sqrt(rho_ll)1983sqrt_rho_rr = sqrt(rho_rr)1984inv_sum_sqrt_rho = inv(sqrt_rho_ll + sqrt_rho_rr)19851986v1_roe = (sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr) * inv_sum_sqrt_rho1987v2_roe = (sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr) * inv_sum_sqrt_rho1988v_roe = v1_roe * normal_direction[1] + v2_roe * normal_direction[2]1989v_roe_mag = v1_roe^2 + v2_roe^219901991H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) * inv_sum_sqrt_rho1992c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * v_roe_mag)) * norm_19931994# Compute convenience constant for positivity preservation, see1995# https://doi.org/10.1016/0021-9991(91)90211-31996beta = sqrt(0.5f0 * (equations.gamma - 1) / equations.gamma)19971998# Estimate the edges of the Riemann fan (with positivity conservation)1999SsL = min(v_roe - c_roe, v_dot_n_ll - beta * c_ll, 0)2000SsR = max(v_roe + c_roe, v_dot_n_rr + beta * c_rr, 0)20012002return SsL, SsR2003end20042005@inline function max_abs_speeds(u, equations::CompressibleEulerEquations2D)2006rho, v1, v2, p = cons2prim(u, equations)2007c = sqrt(equations.gamma * p / rho)20082009return abs(v1) + c, abs(v2) + c2010end20112012# Convert conservative variables to primitive2013@inline function cons2prim(u, equations::CompressibleEulerEquations2D)2014rho, rho_v1, rho_v2, rho_e_total = u20152016v1 = rho_v1 / rho2017v2 = rho_v2 / rho2018p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))20192020return SVector(rho, v1, v2, p)2021end20222023# Convert conservative variables to entropy2024@inline function cons2entropy(u, equations::CompressibleEulerEquations2D)2025rho, rho_v1, rho_v2, rho_e_total = u20262027v1 = rho_v1 / rho2028v2 = rho_v2 / rho2029v_square = v1^2 + v2^22030p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho * v_square)2031s = log(p) - equations.gamma * log(rho)2032rho_p = rho / p20332034w1 = (equations.gamma - s) * equations.inv_gamma_minus_one -20350.5f0 * rho_p * v_square2036w2 = rho_p * v12037w3 = rho_p * v22038w4 = -rho_p20392040return SVector(w1, w2, w3, w4)2041end20422043# Transformation from conservative variables u to entropy vector ds_0/du,2044# using the modified specific entropy of Guermond et al. (2019): s_0 = p * rho^(-gamma) / (gamma-1).2045# Note: This is *not* the "conventional" specific entropy s = ln(p / rho^(gamma)).2046@inline function cons2entropy_guermond_etal(u, equations::CompressibleEulerEquations2D)2047rho, rho_v1, rho_v2, rho_e_total = u20482049v1 = rho_v1 / rho2050v2 = rho_v2 / rho2051v_square = v1^2 + v2^22052inv_rho_gammap1 = (1 / rho)^(equations.gamma + 1)20532054# The derivative vector for the modified specific entropy of Guermond et al.2055w1 = inv_rho_gammap1 *2056(0.5f0 * rho * (equations.gamma + 1) * v_square -2057equations.gamma * rho_e_total)2058w2 = -rho_v1 * inv_rho_gammap12059w3 = -rho_v2 * inv_rho_gammap12060w4 = (1 / rho)^equations.gamma20612062return SVector(w1, w2, w3, w4)2063end20642065@inline function entropy2cons(w, equations::CompressibleEulerEquations2D)2066# See Hughes, Franca, Mallet (1986) A new finite element formulation for CFD2067# [DOI: 10.1016/0045-7825(86)90127-1](https://doi.org/10.1016/0045-7825(86)90127-1)2068@unpack gamma = equations20692070# convert to entropy `-rho * s` used by Hughes, France, Mallet (1986)2071# instead of `-rho * s / (gamma - 1)`2072V1, V2, V3, V5 = w .* (gamma - 1)20732074# s = specific entropy, eq. (53)2075s = gamma - V1 + (V2^2 + V3^2) / (2 * V5)20762077# eq. (52)2078rho_iota = ((gamma - 1) / (-V5)^gamma)^(equations.inv_gamma_minus_one) *2079exp(-s * equations.inv_gamma_minus_one)20802081# eq. (51)2082rho = -rho_iota * V52083rho_v1 = rho_iota * V22084rho_v2 = rho_iota * V32085rho_e_total = rho_iota * (1 - (V2^2 + V3^2) / (2 * V5))2086return SVector(rho, rho_v1, rho_v2, rho_e_total)2087end20882089# Convert primitive to conservative variables2090@inline function prim2cons(prim, equations::CompressibleEulerEquations2D)2091rho, v1, v2, p = prim2092rho_v1 = rho * v12093rho_v2 = rho * v22094rho_e_total = p * equations.inv_gamma_minus_one +20950.5f0 * (rho_v1 * v1 + rho_v2 * v2)2096return SVector(rho, rho_v1, rho_v2, rho_e_total)2097end20982099@inline function density(u, equations::CompressibleEulerEquations2D)2100rho = u[1]2101return rho2102end21032104@inline function velocity(u, equations::CompressibleEulerEquations2D)2105rho = u[1]2106v1 = u[2] / rho2107v2 = u[3] / rho2108return SVector(v1, v2)2109end21102111@inline function velocity(u, orientation::Int, equations::CompressibleEulerEquations2D)2112rho = u[1]2113v = u[orientation + 1] / rho2114return v2115end21162117@inline function pressure(u, equations::CompressibleEulerEquations2D)2118rho, rho_v1, rho_v2, rho_e_total = u2119p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1^2 + rho_v2^2) / rho)2120return p2121end21222123# Transformation from conservative variables u to d(p)/d(u)2124@inline function gradient_conservative(::typeof(pressure),2125u, equations::CompressibleEulerEquations2D)2126rho, rho_v1, rho_v2, rho_e_total = u21272128v1 = rho_v1 / rho2129v2 = rho_v2 / rho2130v_square = v1^2 + v2^221312132return (equations.gamma - 1) * SVector(0.5f0 * v_square, -v1, -v2, 1)2133end21342135@inline function density_pressure(u, equations::CompressibleEulerEquations2D)2136rho, rho_v1, rho_v2, rho_e_total = u2137rho_times_p = (equations.gamma - 1) *2138(rho * rho_e_total - 0.5f0 * (rho_v1^2 + rho_v2^2))2139return rho_times_p2140end21412142# Calculates the entropy flux in direction "orientation" and the entropy variables for a state cons2143# NOTE: This method seems to work currently (b82534e) but is never used anywhere. Thus it is2144# commented here until someone uses it or writes a test for it.2145# @inline function cons2entropyvars_and_flux(gamma::Float64, cons, orientation::Int)2146# entropy = MVector{4, Float64}(undef)2147# v = (cons[2] / cons[1] , cons[3] / cons[1])2148# v_square= v[1]*v[1]+v[2]*v[2]2149# p = (gamma - 1) * (cons[4] - 1/2 * (cons[2] * v[1] + cons[3] * v[2]))2150# rho_p = cons[1] / p2151# # thermodynamic entropy2152# s = log(p) - gamma*log(cons[1])2153# # mathematical entropy2154# S = - s*cons[1]/(gamma-1)2155# # entropy variables2156# entropy[1] = (gamma - s)/(gamma-1) - 0.5*rho_p*v_square2157# entropy[2] = rho_p*v[1]2158# entropy[3] = rho_p*v[2]2159# entropy[4] = -rho_p2160# # entropy flux2161# entropy_flux = S*v[orientation]2162# return entropy, entropy_flux2163# end21642165# Calculate thermodynamic entropy for a conservative state `cons`2166@inline function entropy_thermodynamic(cons, equations::CompressibleEulerEquations2D)2167# Pressure2168p = (equations.gamma - 1) * (cons[4] - 0.5f0 * (cons[2]^2 + cons[3]^2) / cons[1])21692170# Thermodynamic entropy2171s = log(p) - equations.gamma * log(cons[1])21722173return s2174end21752176# Calculate mathematical entropy for a conservative state `cons`2177@inline function entropy_math(cons, equations::CompressibleEulerEquations2D)2178# Mathematical entropy2179S = -entropy_thermodynamic(cons, equations) * cons[1] *2180equations.inv_gamma_minus_one21812182return S2183end21842185# Transformation from conservative variables u to d(s)/d(u)2186@inline function gradient_conservative(::typeof(entropy_math),2187u, equations::CompressibleEulerEquations2D)2188return cons2entropy(u, equations)2189end21902191@doc raw"""2192entropy_guermond_etal(u, equations::CompressibleEulerEquations2D)21932194Calculate the modified specific entropy of Guermond et al. (2019):2195```math2196s_0 = p * \rho^{-\gamma} / (\gamma-1).2197```2198Note: This is *not* the "conventional" specific entropy ``s = ln(p / \rho^\gamma)``.2199- Guermond at al. (2019)2200Invariant domain preserving discretization-independent schemes and convex limiting for hyperbolic systems.2201[DOI: 10.1016/j.cma.2018.11.036](https://doi.org/10.1016/j.cma.2018.11.036)2202"""2203@inline function entropy_guermond_etal(u, equations::CompressibleEulerEquations2D)2204rho, rho_v1, rho_v2, rho_e_total = u22052206# Modified specific entropy from Guermond et al. (2019)2207s = (rho_e_total - 0.5f0 * (rho_v1^2 + rho_v2^2) / rho) * (1 / rho)^equations.gamma22082209return s2210end22112212# Transformation from conservative variables u to d(s)/d(u)2213@inline function gradient_conservative(::typeof(entropy_guermond_etal),2214u, equations::CompressibleEulerEquations2D)2215return cons2entropy_guermond_etal(u, equations)2216end22172218# Default entropy is the mathematical entropy2219@inline function entropy(cons, equations::CompressibleEulerEquations2D)2220return entropy_math(cons, equations)2221end22222223# Calculate total energy for a conservative state `cons`2224@inline energy_total(cons, ::CompressibleEulerEquations2D) = cons[4]22252226# Calculate kinetic energy for a conservative state `cons`2227@inline function energy_kinetic(u, equations::CompressibleEulerEquations2D)2228rho, rho_v1, rho_v2, rho_e_total = u2229return (rho_v1^2 + rho_v2^2) / (2 * rho)2230end22312232# Calculate internal energy for a conservative state `cons`2233@inline function energy_internal(cons, equations::CompressibleEulerEquations2D)2234return energy_total(cons, equations) - energy_kinetic(cons, equations)2235end22362237@inline function entropy_potential(u, orientation::Int,2238equations::CompressibleEulerEquations2D)2239if orientation == 12240return u[2]2241else # if orientation == 22242return u[3]2243end2244end2245# Version for non-Cartesian meshes, i.e., everything but `TreeMesh`es.2246@inline function entropy_potential(u, normal_direction::AbstractVector,2247equations::CompressibleEulerEquations2D)2248return u[2] * normal_direction[1] +2249u[3] * normal_direction[2]2250end22512252# State validation for Newton-bisection method of subcell IDP limiting2253@inline function Base.isvalid(u, equations::CompressibleEulerEquations2D)2254if u[1] <= 0 || pressure(u, equations) <= 02255return false2256end2257return true2258end2259end # @muladd226022612262