Path: blob/main/src/equations/ideal_glm_mhd_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"""8IdealGlmMhdEquations2D(gamma)910The ideal compressible GLM-MHD equations11```math12\frac{\partial}{\partial t}13\begin{pmatrix}14\rho \\ \rho \mathbf{v} \\ \rho e_{\text{total}} \\ \mathbf{B} \\ \psi15\end{pmatrix}16+17\nabla \cdot18\begin{pmatrix}19\rho \mathbf{v} \\ \rho (\mathbf{v} \otimes \mathbf{v}) + (p + \frac{1}{2} \Vert \mathbf{B} \Vert_2 ^2) \underline{I} - \mathbf{B} \otimes \mathbf{B} \\20\mathbf{v} (\frac{1}{2} \rho \Vert \mathbf{v} \Vert_2 ^2 + \frac{\gamma p}{\gamma - 1} + \Vert \mathbf{B} \Vert_2 ^2) - \mathbf{B} (\mathbf{v} \cdot \mathbf{B}) + c_h \psi \mathbf{B} \\21\mathbf{v} \otimes \mathbf{B} - \mathbf{B} \otimes \mathbf{v} + c_h \psi \underline{I} \\22c_h \mathbf{B}23\end{pmatrix}24+25(\nabla \cdot \mathbf{B})26\begin{pmatrix}270 \\ \mathbf{B} \\ \mathbf{v} \cdot \mathbf{B} \\ \mathbf{v} \\ 028\end{pmatrix}29+30(\nabla \psi) \cdot31\begin{pmatrix}32\mathbf{0} \\ 0 \\ \mathbf{v} \cdot \psi \\ 0 \\ \mathbf{v}33\end{pmatrix}34=35\begin{pmatrix}360 \\ \mathbf{0} \\ 0 \\ \mathbf{0} \\ 037\end{pmatrix}38```39for an ideal gas in two space dimensions.40Here, ``\mathbf{v}`` is the velocity, ``\mathbf{B}`` the magnetic field, ``c_h`` the hyperbolic divergence cleaning speed,41``\psi`` the generalized Lagrangian Multiplier (GLM),42``e_{\text{total}}`` the specific total energy, and43```math44p = (\gamma - 1) \left( \rho e_{\text{total}} - \frac{1}{2} \rho \Vert \mathbf{v} \Vert_2 ^2 - \frac{1}{2} \Vert \mathbf{B} \Vert_2 ^2 - \frac{1}{2} \psi^2 \right)45```46the pressure, ``\gamma`` the total heat capacity ratio and ``\underline{I}`` the ``3\times 3`` identity matrix.47"""48struct IdealGlmMhdEquations2D{RealT <: Real} <:49AbstractIdealGlmMhdEquations{2, 9}50gamma::RealT # ratio of specific heats51inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications52c_h::RealT # GLM cleaning speed5354function IdealGlmMhdEquations2D(gamma, c_h)55γ, inv_gamma_minus_one, c_h = promote(gamma, inv(gamma - 1), c_h)56return new{typeof(γ)}(γ, inv_gamma_minus_one, c_h)57end58end5960function IdealGlmMhdEquations2D(gamma; initial_c_h = convert(typeof(gamma), NaN))61# Use `promote` to ensure that `gamma` and `initial_c_h` have the same type62return IdealGlmMhdEquations2D(promote(gamma, initial_c_h)...)63end6465# Outer constructor for `@reset` works correctly66function IdealGlmMhdEquations2D(gamma, inv_gamma_minus_one, c_h)67return IdealGlmMhdEquations2D(gamma, c_h)68end6970have_nonconservative_terms(::IdealGlmMhdEquations2D) = True()7172function varnames(::typeof(cons2cons), ::IdealGlmMhdEquations2D)73return ("rho", "rho_v1", "rho_v2", "rho_v3", "rho_e_total", "B1", "B2", "B3", "psi")74end75function varnames(::typeof(cons2prim), ::IdealGlmMhdEquations2D)76return ("rho", "v1", "v2", "v3", "p", "B1", "B2", "B3", "psi")77end78function default_analysis_integrals(::IdealGlmMhdEquations2D)79return (entropy_timederivative, Val(:l2_divb), Val(:linf_divb))80end8182# Helper function to extract the magnetic field vector from the conservative variables83magnetic_field(u, equations::IdealGlmMhdEquations2D) = SVector(u[6], u[7], u[8])8485# Set initial conditions at physical location `x` for time `t`86"""87initial_condition_constant(x, t, equations::IdealGlmMhdEquations2D)8889A constant initial condition to test free-stream preservation.90"""91function initial_condition_constant(x, t, equations::IdealGlmMhdEquations2D)92RealT = eltype(x)93rho = 194rho_v1 = convert(RealT, 0.1)95rho_v2 = -convert(RealT, 0.2)96rho_v3 = -0.5f097rho_e_total = 5098B1 = 399B2 = -convert(RealT, 1.2)100B3 = 0.5f0101psi = 0102return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi)103end104105"""106initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations2D)107108An Alfvén wave as smooth initial condition used for convergence tests.109See for reference section 4.2 in110- Dominik Derigs, Andrew R. Winters, Gregor J. Gassner, and Stefanie Walch (2016)111A novel high-order, entropy stable, 3D AMR MHD solver with guaranteed positive pressure112[DOI: 10.1016/j.jcp.2016.04.048](https://doi.org/10.1016/j.jcp.2016.04.048)113"""114function initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations2D)115# smooth Alfvén wave test from Derigs et al. (2016)116# domain must be set to [0, 1/cos(α)] x [0, 1/sin(α)], γ = 5/3117RealT = eltype(x)118alpha = 0.25f0 * convert(RealT, pi)119x_perp = x[1] * cos(alpha) + x[2] * sin(alpha)120B_perp = convert(RealT, 0.1) * sinpi(2 * (x_perp + t))121rho = 1122v1 = -B_perp * sin(alpha)123v2 = B_perp * cos(alpha)124v3 = convert(RealT, 0.1) * cospi(2 * (x_perp + t))125p = convert(RealT, 0.1)126B1 = cos(alpha) + v1127B2 = sin(alpha) + v2128B3 = v3129psi = 0130return prim2cons(SVector(rho, v1, v2, v3, p, B1, B2, B3, psi), equations)131end132133"""134initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations2D)135136A weak blast wave adapted from137- Sebastian Hennemann, Gregor J. Gassner (2020)138A provably entropy stable subcell shock capturing approach for high order split form DG139[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)140"""141function initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations2D)142# Adapted MHD version of the weak blast wave from Hennemann & Gassner JCP paper 2020 (Sec. 6.3)143# Same discontinuity in the velocities but with magnetic fields144# Set up polar coordinates145RealT = eltype(x)146inicenter = (0, 0)147x_norm = x[1] - inicenter[1]148y_norm = x[2] - inicenter[2]149r = sqrt(x_norm^2 + y_norm^2)150phi = atan(y_norm, x_norm)151152# Calculate primitive variables153rho = r > 0.5f0 ? one(RealT) : convert(RealT, 1.1691)154v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos(phi)155v2 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * sin(phi)156p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)157158return prim2cons(SVector(rho, v1, v2, 0, p, 1, 1, 1, 0), equations)159end160161# Pre-defined source terms should be implemented as162# function source_terms_WHATEVER(u, x, t, equations::IdealGlmMhdEquations2D)163164# Calculate 1D flux for a single point165@inline function flux(u, orientation::Integer, equations::IdealGlmMhdEquations2D)166rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u167v1 = rho_v1 / rho168v2 = rho_v2 / rho169v3 = rho_v3 / rho170kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)171mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)172p_over_gamma_minus_one = (rho_e_total - kin_en - mag_en - 0.5f0 * psi^2)173p = (equations.gamma - 1) * p_over_gamma_minus_one174if orientation == 1175f1 = rho_v1176f2 = rho_v1 * v1 + p + mag_en - B1^2177f3 = rho_v1 * v2 - B1 * B2178f4 = rho_v1 * v3 - B1 * B3179f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v1 -180B1 * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B1181f6 = equations.c_h * psi182f7 = v1 * B2 - v2 * B1183f8 = v1 * B3 - v3 * B1184f9 = equations.c_h * B1185else #if orientation == 2186f1 = rho_v2187f2 = rho_v2 * v1 - B2 * B1188f3 = rho_v2 * v2 + p + mag_en - B2^2189f4 = rho_v2 * v3 - B2 * B3190f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v2 -191B2 * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B2192f6 = v2 * B1 - v1 * B2193f7 = equations.c_h * psi194f8 = v2 * B3 - v3 * B2195f9 = equations.c_h * B2196end197198return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)199end200201# Calculate 1D flux for a single point in the normal direction202# Note, this directional vector is not normalized203@inline function flux(u, normal_direction::AbstractVector,204equations::IdealGlmMhdEquations2D)205rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u206v1 = rho_v1 / rho207v2 = rho_v2 / rho208v3 = rho_v3 / rho209kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)210mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)211p_over_gamma_minus_one = (rho_e_total - kin_en - mag_en - 0.5f0 * psi^2)212p = (equations.gamma - 1) * p_over_gamma_minus_one213214v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]215B_normal = B1 * normal_direction[1] + B2 * normal_direction[2]216rho_v_normal = rho * v_normal217218f1 = rho_v_normal219f2 = rho_v_normal * v1 - B1 * B_normal + (p + mag_en) * normal_direction[1]220f3 = rho_v_normal * v2 - B2 * B_normal + (p + mag_en) * normal_direction[2]221f4 = rho_v_normal * v3 - B3 * B_normal222f5 = ((kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v_normal223-224B_normal * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B_normal)225f6 = equations.c_h * psi * normal_direction[1] +226(v2 * B1 - v1 * B2) * normal_direction[2]227f7 = equations.c_h * psi * normal_direction[2] +228(v1 * B2 - v2 * B1) * normal_direction[1]229f8 = v_normal * B3 - v3 * B_normal230f9 = equations.c_h * B_normal231232return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)233end234235"""236flux_nonconservative_powell(u_ll, u_rr, orientation::Integer,237equations::IdealGlmMhdEquations2D)238flux_nonconservative_powell(u_ll, u_rr,239normal_direction::AbstractVector,240equations::IdealGlmMhdEquations2D)241242Non-symmetric two-point flux discretizing the nonconservative (source) term of243Powell and the Galilean nonconservative term associated with the GLM multiplier244of the [`IdealGlmMhdEquations2D`](@ref).245246On curvilinear meshes, the implementation differs from the reference since we use the averaged247normal direction for the metrics dealiasing.248249## References250- Marvin Bohm, Andrew R.Winters, Gregor J. Gassner, Dominik Derigs,251Florian Hindenlang, Joachim Saur252An entropy stable nodal discontinuous Galerkin method for the resistive MHD253equations. Part I: Theory and numerical verification254[DOI: 10.1016/j.jcp.2018.06.027](https://doi.org/10.1016/j.jcp.2018.06.027)255"""256@inline function flux_nonconservative_powell(u_ll, u_rr, orientation::Integer,257equations::IdealGlmMhdEquations2D)258rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll259rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr260261v1_ll = rho_v1_ll / rho_ll262v2_ll = rho_v2_ll / rho_ll263v3_ll = rho_v3_ll / rho_ll264v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll265266# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)267# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})268if orientation == 1269f = SVector(0,270B1_ll * B1_rr,271B2_ll * B1_rr,272B3_ll * B1_rr,273v_dot_B_ll * B1_rr + v1_ll * psi_ll * psi_rr,274v1_ll * B1_rr,275v2_ll * B1_rr,276v3_ll * B1_rr,277v1_ll * psi_rr)278else # orientation == 2279f = SVector(0,280B1_ll * B2_rr,281B2_ll * B2_rr,282B3_ll * B2_rr,283v_dot_B_ll * B2_rr + v2_ll * psi_ll * psi_rr,284v1_ll * B2_rr,285v2_ll * B2_rr,286v3_ll * B2_rr,287v2_ll * psi_rr)288end289290return f291end292293@inline function flux_nonconservative_powell(u_ll, u_rr,294normal_direction::AbstractVector,295equations::IdealGlmMhdEquations2D)296rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll297rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr298299v1_ll = rho_v1_ll / rho_ll300v2_ll = rho_v2_ll / rho_ll301v3_ll = rho_v3_ll / rho_ll302v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll303304v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]305B_dot_n_rr = B1_rr * normal_direction[1] +306B2_rr * normal_direction[2]307308# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)309# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})310f = SVector(0,311B1_ll * B_dot_n_rr,312B2_ll * B_dot_n_rr,313B3_ll * B_dot_n_rr,314v_dot_B_ll * B_dot_n_rr + v_dot_n_ll * psi_ll * psi_rr,315v1_ll * B_dot_n_rr,316v2_ll * B_dot_n_rr,317v3_ll * B_dot_n_rr,318v_dot_n_ll * psi_rr)319320return f321end322323# For `VolumeIntegralSubcellLimiting` the nonconservative flux is created as a callable struct to324# enable dispatch on the type of the nonconservative term (symmetric / jump).325struct FluxNonConservativePowellLocalSymmetric <:326FluxNonConservative{NonConservativeSymmetric()}327end328329n_nonconservative_terms(::FluxNonConservativePowellLocalSymmetric) = 2330331const flux_nonconservative_powell_local_symmetric = FluxNonConservativePowellLocalSymmetric()332333"""334flux_nonconservative_powell_local_symmetric(u_ll, u_rr,335orientation::Integer,336equations::IdealGlmMhdEquations2D)337flux_nonconservative_powell_local_symmetric(u_ll, u_rr,338normal_direction::AbstractVector,339equations::IdealGlmMhdEquations2D)340341Non-symmetric two-point flux discretizing the nonconservative (source) term of342Powell and the Galilean nonconservative term associated with the GLM multiplier343of the [`IdealGlmMhdEquations2D`](@ref).344345This implementation uses a non-conservative term that can be written as the product346of local and symmetric parts. It is equivalent to the non-conservative flux of Bohm347et al. [`flux_nonconservative_powell`](@ref) for conforming meshes but it yields different348results on non-conforming meshes(!). On curvilinear meshes this formulation applies the349local normal direction compared to the averaged one used in [`flux_nonconservative_powell`](@ref).350351The two other flux functions with the same name return either the local352or symmetric portion of the non-conservative flux based on the type of the353nonconservative_type argument, employing multiple dispatch. They are used to354compute the subcell fluxes in dg_2d_subcell_limiters.jl.355356## References357- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts358Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.359"""360@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll, u_rr,361orientation::Integer,362equations::IdealGlmMhdEquations2D)363rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll364rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr365366v1_ll = rho_v1_ll / rho_ll367v2_ll = rho_v2_ll / rho_ll368v3_ll = rho_v3_ll / rho_ll369v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll370371# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)372# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})373psi_avg = (psi_ll + psi_rr) #* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code374if orientation == 1375B1_avg = (B1_ll + B1_rr) #* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code376f = SVector(0,377B1_ll * B1_avg,378B2_ll * B1_avg,379B3_ll * B1_avg,380v_dot_B_ll * B1_avg + v1_ll * psi_ll * psi_avg,381v1_ll * B1_avg,382v2_ll * B1_avg,383v3_ll * B1_avg,384v1_ll * psi_avg)385else # orientation == 2386B2_avg = (B2_ll + B2_rr) #* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code387f = SVector(0,388B1_ll * B2_avg,389B2_ll * B2_avg,390B3_ll * B2_avg,391v_dot_B_ll * B2_avg + v2_ll * psi_ll * psi_avg,392v1_ll * B2_avg,393v2_ll * B2_avg,394v3_ll * B2_avg,395v2_ll * psi_avg)396end397398return f399end400401@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll, u_rr,402normal_direction::AbstractVector,403equations::IdealGlmMhdEquations2D)404rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll405rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr406407v1_ll = rho_v1_ll / rho_ll408v2_ll = rho_v2_ll / rho_ll409v3_ll = rho_v3_ll / rho_ll410v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll411412# The factor 0.5 of the averages can be omitted since it is already applied when this413# function is called.414psi_avg = (psi_ll + psi_rr)415B1_avg = (B1_ll + B1_rr)416B2_avg = (B2_ll + B2_rr)417418B_dot_n_avg = B1_avg * normal_direction[1] + B2_avg * normal_direction[2]419v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]420421# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)422# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})423f = SVector(0,424B1_ll * B_dot_n_avg,425B2_ll * B_dot_n_avg,426B3_ll * B_dot_n_avg,427v_dot_B_ll * B_dot_n_avg + v_dot_n_ll * psi_ll * psi_avg,428v1_ll * B_dot_n_avg,429v2_ll * B_dot_n_avg,430v3_ll * B_dot_n_avg,431v_dot_n_ll * psi_avg)432433return f434end435436"""437flux_nonconservative_powell_local_symmetric(u_ll, orientation::Integer,438equations::IdealGlmMhdEquations2D,439nonconservative_type::NonConservativeLocal,440nonconservative_term::Integer)441flux_nonconservative_powell_local_symmetric(u_ll, normal_direction_ll::AbstractVector,442equations::IdealGlmMhdEquations2D,443nonconservative_type::NonConservativeLocal,444nonconservative_term::Integer)445446Local part of the Powell and GLM non-conservative terms. Needed for the calculation of447the non-conservative staggered "fluxes" for subcell limiting. See, e.g.,448- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts449Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.450This function is used to compute the subcell fluxes in dg_2d_subcell_limiters.jl.451"""452@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll,453orientation::Integer,454equations::IdealGlmMhdEquations2D,455nonconservative_type::NonConservativeLocal,456nonconservative_term::Integer)457rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll458459if nonconservative_term == 1460# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)461v1_ll = rho_v1_ll / rho_ll462v2_ll = rho_v2_ll / rho_ll463v3_ll = rho_v3_ll / rho_ll464v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll465f = SVector(0,466B1_ll,467B2_ll,468B3_ll,469v_dot_B_ll,470v1_ll,471v2_ll,472v3_ll,4730)474else #nonconservative_term ==2475# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})476if orientation == 1477v1_ll = rho_v1_ll / rho_ll478f = SVector(0,4790,4800,4810,482v1_ll * psi_ll,4830,4840,4850,486v1_ll)487else #orientation == 2488v2_ll = rho_v2_ll / rho_ll489f = SVector(0,4900,4910,4920,493v2_ll * psi_ll,4940,4950,4960,497v2_ll)498end499end500return f501end502503@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll,504normal_direction_ll::AbstractVector,505equations::IdealGlmMhdEquations2D,506nonconservative_type::NonConservativeLocal,507nonconservative_term::Integer)508rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll509510if nonconservative_term == 1511# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)512v1_ll = rho_v1_ll / rho_ll513v2_ll = rho_v2_ll / rho_ll514v3_ll = rho_v3_ll / rho_ll515v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll516517f = SVector(0,518B1_ll,519B2_ll,520B3_ll,521v_dot_B_ll,522v1_ll,523v2_ll,524v3_ll,5250)526else # nonconservative_term == 2527# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})528v1_ll = rho_v1_ll / rho_ll529v2_ll = rho_v2_ll / rho_ll530v_dot_n_ll = v1_ll * normal_direction_ll[1] + v2_ll * normal_direction_ll[2]531532f = SVector(0,5330,5340,5350,536v_dot_n_ll * psi_ll,5370,5380,5390,540v_dot_n_ll)541end542return f543end544545"""546flux_nonconservative_powell_local_symmetric(u_ll, orientation::Integer,547equations::IdealGlmMhdEquations2D,548nonconservative_type::NonConservativeSymmetric,549nonconservative_term::Integer)550flux_nonconservative_powell_local_symmetric(u_ll, normal_direction_avg::AbstractVector,551equations::IdealGlmMhdEquations2D,552nonconservative_type::NonConservativeSymmetric,553nonconservative_term::Integer)554555Symmetric part of the Powell and GLM non-conservative terms. Needed for the calculation of556the non-conservative staggered "fluxes" for subcell limiting. See, e.g.,557- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts558Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.559This function is used to compute the subcell fluxes in dg_2d_subcell_limiters.jl.560"""561@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll, u_rr,562orientation::Integer,563equations::IdealGlmMhdEquations2D,564nonconservative_type::NonConservativeSymmetric,565nonconservative_term::Integer)566rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll567rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr568569if nonconservative_term == 1570# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)571if orientation == 1572B1_avg = (B1_ll + B1_rr)#* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code573f = SVector(0,574B1_avg,575B1_avg,576B1_avg,577B1_avg,578B1_avg,579B1_avg,580B1_avg,5810)582else # orientation == 2583B2_avg = (B2_ll + B2_rr)#* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code584f = SVector(0,585B2_avg,586B2_avg,587B2_avg,588B2_avg,589B2_avg,590B2_avg,591B2_avg,5920)593end594else #nonconservative_term == 2595# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})596psi_avg = (psi_ll + psi_rr)#* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code597f = SVector(0,5980,5990,6000,601psi_avg,6020,6030,6040,605psi_avg)606end607608return f609end610611@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll, u_rr,612normal_direction_avg::AbstractVector,613equations::IdealGlmMhdEquations2D,614nonconservative_type::NonConservativeSymmetric,615nonconservative_term::Integer)616rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll617rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr618619if nonconservative_term == 1620# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)621# The factor 0.5 of the average can be omitted since it is already applied when this622# function is called.623B_dot_n_avg = (B1_ll + B1_rr) * normal_direction_avg[1] +624(B2_ll + B2_rr) * normal_direction_avg[2]625f = SVector(0,626B_dot_n_avg,627B_dot_n_avg,628B_dot_n_avg,629B_dot_n_avg,630B_dot_n_avg,631B_dot_n_avg,632B_dot_n_avg,6330)634else # nonconservative_term == 2635# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})636# The factor 0.5 of the average can be omitted since it is already applied when this637# function is called.638psi_avg = (psi_ll + psi_rr)639f = SVector(0,6400,6410,6420,643psi_avg,6440,6450,6460,647psi_avg)648end649650return f651end652653# For `VolumeIntegralSubcellLimiting` the nonconservative flux is created as a callable struct to654# enable dispatch on the type of the nonconservative term (symmetric / jump).655struct FluxNonConservativePowellLocalJump <:656FluxNonConservative{NonConservativeJump()}657end658659n_nonconservative_terms(::FluxNonConservativePowellLocalJump) = 2660661const flux_nonconservative_powell_local_jump = FluxNonConservativePowellLocalJump()662663"""664flux_nonconservative_powell_local_jump(u_ll, u_rr,665orientation::Integer,666equations::IdealGlmMhdEquations2D)667flux_nonconservative_powell_local_jump(u_ll, u_rr,668normal_direction::AbstractVector,669equations::IdealGlmMhdEquations2D)670671Non-symmetric two-point flux discretizing the nonconservative (source) term of672Powell and the Galilean nonconservative term associated with the GLM multiplier673of the [`IdealGlmMhdEquations2D`](@ref).674675This implementation uses a non-conservative term that can be written as the product676of local and jump parts.677678The two other flux functions with the same name return either the local679or jump portion of the non-conservative flux based on the type of the680nonconservative_type argument, employing multiple dispatch. They are used to681compute the subcell fluxes in dg_2d_subcell_limiters.jl.682683## References684- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts685Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.686"""687@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll, u_rr,688orientation::Integer,689equations::IdealGlmMhdEquations2D)690rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll691rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr692693v1_ll = rho_v1_ll / rho_ll694v2_ll = rho_v2_ll / rho_ll695v3_ll = rho_v3_ll / rho_ll696v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll697698# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)699# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})700psi_jump = psi_rr - psi_ll701if orientation == 1702B1_jump = B1_rr - B1_ll # The flux is already multiplied by 0.5 wherever it is used in the code703f = SVector(0,704B1_ll * B1_jump,705B2_ll * B1_jump,706B3_ll * B1_jump,707v_dot_B_ll * B1_jump + v1_ll * psi_ll * psi_jump,708v1_ll * B1_jump,709v2_ll * B1_jump,710v3_ll * B1_jump,711v1_ll * psi_jump)712else # orientation == 2713B2_jump = B2_rr - B2_ll # The flux is already multiplied by 0.5 wherever it is used in the code714f = SVector(0,715B1_ll * B2_jump,716B2_ll * B2_jump,717B3_ll * B2_jump,718v_dot_B_ll * B2_jump + v2_ll * psi_ll * psi_jump,719v1_ll * B2_jump,720v2_ll * B2_jump,721v3_ll * B2_jump,722v2_ll * psi_jump)723end724725return f726end727728@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll, u_rr,729normal_direction::AbstractVector,730equations::IdealGlmMhdEquations2D)731rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll732rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr733734v1_ll = rho_v1_ll / rho_ll735v2_ll = rho_v2_ll / rho_ll736v3_ll = rho_v3_ll / rho_ll737v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll738739psi_jump = psi_rr - psi_ll740B1_jump = B1_rr - B1_ll741B2_jump = B2_rr - B2_ll742743B_dot_n_jump = B1_jump * normal_direction[1] + B2_jump * normal_direction[2]744v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]745746# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)747# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})748f = SVector(0,749B1_ll * B_dot_n_jump,750B2_ll * B_dot_n_jump,751B3_ll * B_dot_n_jump,752v_dot_B_ll * B_dot_n_jump + v_dot_n_ll * psi_ll * psi_jump,753v1_ll * B_dot_n_jump,754v2_ll * B_dot_n_jump,755v3_ll * B_dot_n_jump,756v_dot_n_ll * psi_jump)757758return f759end760761"""762flux_nonconservative_powell_local_jump(u_ll, orientation::Integer,763equations::IdealGlmMhdEquations2D,764nonconservative_type::NonConservativeLocal,765nonconservative_term::Integer)766flux_nonconservative_powell_local_jump(u_ll, normal_direction_ll::AbstractVector,767equations::IdealGlmMhdEquations2D,768nonconservative_type::NonConservativeLocal,769nonconservative_term::Integer)770771Local part of the Powell and GLM non-conservative terms. Needed for the calculation of772the non-conservative staggered "fluxes" for subcell limiting.773This function is used to compute the subcell fluxes in dg_2d_subcell_limiters.jl.774775## References776- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts777Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.778"""779@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll,780orientation::Integer,781equations::IdealGlmMhdEquations2D,782nonconservative_type::NonConservativeLocal,783nonconservative_term::Integer)784rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll785if nonconservative_term == 1786# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)787v1_ll = rho_v1_ll / rho_ll788v2_ll = rho_v2_ll / rho_ll789v3_ll = rho_v3_ll / rho_ll790v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll791f = SVector(0,792B1_ll,793B2_ll,794B3_ll,795v_dot_B_ll,796v1_ll,797v2_ll,798v3_ll,7990)800else #nonconservative_term ==2801# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})802if orientation == 1803v1_ll = rho_v1_ll / rho_ll804f = SVector(0,8050,8060,8070,808v1_ll * psi_ll,8090,8100,8110,812v1_ll)813else #orientation == 2814v2_ll = rho_v2_ll / rho_ll815f = SVector(0,8160,8170,8180,819v2_ll * psi_ll,8200,8210,8220,823v2_ll)824end825end826return f827end828829@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll,830normal_direction_ll::AbstractVector,831equations::IdealGlmMhdEquations2D,832nonconservative_type::NonConservativeLocal,833nonconservative_term::Integer)834rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll835836if nonconservative_term == 1837# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)838v1_ll = rho_v1_ll / rho_ll839v2_ll = rho_v2_ll / rho_ll840v3_ll = rho_v3_ll / rho_ll841v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll842843f = SVector(0,844B1_ll,845B2_ll,846B3_ll,847v_dot_B_ll,848v1_ll,849v2_ll,850v3_ll,8510)852else # nonconservative_term == 2853# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})854v1_ll = rho_v1_ll / rho_ll855v2_ll = rho_v2_ll / rho_ll856v_dot_n_ll = v1_ll * normal_direction_ll[1] + v2_ll * normal_direction_ll[2]857858f = SVector(0,8590,8600,8610,862v_dot_n_ll * psi_ll,8630,8640,8650,866v_dot_n_ll)867end868return f869end870871"""872flux_nonconservative_powell_local_jump(u_ll, orientation::Integer,873equations::IdealGlmMhdEquations2D,874nonconservative_type::NonConservativeJump,875nonconservative_term::Integer)876flux_nonconservative_powell_local_jump(u_ll, normal_direction_avg::AbstractVector,877equations::IdealGlmMhdEquations2D,878nonconservative_type::NonConservativeJump,879nonconservative_term::Integer)880881Jump part of the Powell and GLM non-conservative terms. Needed for the calculation of882the non-conservative staggered "fluxes" for subcell limiting.883This function is used to compute the subcell fluxes in dg_2d_subcell_limiters.jl.884885## References886- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts887Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.888"""889@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll, u_rr,890orientation::Integer,891equations::IdealGlmMhdEquations2D,892nonconservative_type::NonConservativeJump,893nonconservative_term::Integer)894rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll895rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr896897if nonconservative_term == 1898# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)899if orientation == 1900B1_jump = B1_rr - B1_ll # The flux is already multiplied by 0.5 wherever it is used in the code901f = SVector(0,902B1_jump,903B1_jump,904B1_jump,905B1_jump,906B1_jump,907B1_jump,908B1_jump,9090)910else # orientation == 2911B2_jump = B2_rr - B2_ll # The flux is already multiplied by 0.5 wherever it is used in the code912f = SVector(0,913B2_jump,914B2_jump,915B2_jump,916B2_jump,917B2_jump,918B2_jump,919B2_jump,9200)921end922else #nonconservative_term == 2923# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})924psi_jump = psi_rr - psi_ll # The flux is already multiplied by 0.5 wherever it is used in the code925f = SVector(0,9260,9270,9280,929psi_jump,9300,9310,9320,933psi_jump)934end935936return f937end938939@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll, u_rr,940normal_direction_avg::AbstractVector,941equations::IdealGlmMhdEquations2D,942nonconservative_type::NonConservativeJump,943nonconservative_term::Integer)944rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll945rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr946947if nonconservative_term == 1948# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)949B1_jump = B1_rr - B1_ll950B2_jump = B2_rr - B2_ll951B_dot_n_jump = B1_jump * normal_direction_avg[1] +952B2_jump * normal_direction_avg[2]953f = SVector(0,954B_dot_n_jump,955B_dot_n_jump,956B_dot_n_jump,957B_dot_n_jump,958B_dot_n_jump,959B_dot_n_jump,960B_dot_n_jump,9610)962else # nonconservative_term == 2963# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})964psi_jump = (psi_rr - psi_ll)965f = SVector(0,9660,9670,9680,969psi_jump,9700,9710,9720,973psi_jump)974end975976return f977end978979"""980flux_derigs_etal(u_ll, u_rr, orientation, equations::IdealGlmMhdEquations2D)981982Entropy conserving two-point flux by983- Derigs et al. (2018)984Ideal GLM-MHD: About the entropy consistent nine-wave magnetic field985divergence diminishing ideal magnetohydrodynamics equations986[DOI: 10.1016/j.jcp.2018.03.002](https://doi.org/10.1016/j.jcp.2018.03.002)987"""988function flux_derigs_etal(u_ll, u_rr, orientation::Integer,989equations::IdealGlmMhdEquations2D)990# Unpack left and right states to get velocities, pressure, and inverse temperature (called beta)991rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll992rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr993994v1_ll = rho_v1_ll / rho_ll995v2_ll = rho_v2_ll / rho_ll996v3_ll = rho_v3_ll / rho_ll997v1_rr = rho_v1_rr / rho_rr998v2_rr = rho_v2_rr / rho_rr999v3_rr = rho_v3_rr / rho_rr1000vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^21001vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^21002mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^21003mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^21004p_ll = (equations.gamma - 1) *1005(rho_e_total_ll - 0.5f0 * rho_ll * vel_norm_ll - 0.5f0 * mag_norm_ll -10060.5f0 * psi_ll^2)1007p_rr = (equations.gamma - 1) *1008(rho_e_total_rr - 0.5f0 * rho_rr * vel_norm_rr - 0.5f0 * mag_norm_rr -10090.5f0 * psi_rr^2)1010beta_ll = 0.5f0 * rho_ll / p_ll1011beta_rr = 0.5f0 * rho_rr / p_rr1012# for convenience store v⋅B1013vel_dot_mag_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll1014vel_dot_mag_rr = v1_rr * B1_rr + v2_rr * B2_rr + v3_rr * B3_rr10151016# Compute the necessary mean values needed for either direction1017rho_avg = 0.5f0 * (rho_ll + rho_rr)1018rho_mean = ln_mean(rho_ll, rho_rr)1019beta_mean = ln_mean(beta_ll, beta_rr)1020beta_avg = 0.5f0 * (beta_ll + beta_rr)1021v1_avg = 0.5f0 * (v1_ll + v1_rr)1022v2_avg = 0.5f0 * (v2_ll + v2_rr)1023v3_avg = 0.5f0 * (v3_ll + v3_rr)1024p_mean = 0.5f0 * rho_avg / beta_avg1025B1_avg = 0.5f0 * (B1_ll + B1_rr)1026B2_avg = 0.5f0 * (B2_ll + B2_rr)1027B3_avg = 0.5f0 * (B3_ll + B3_rr)1028psi_avg = 0.5f0 * (psi_ll + psi_rr)1029vel_norm_avg = 0.5f0 * (vel_norm_ll + vel_norm_rr)1030mag_norm_avg = 0.5f0 * (mag_norm_ll + mag_norm_rr)1031vel_dot_mag_avg = 0.5f0 * (vel_dot_mag_ll + vel_dot_mag_rr)10321033# Calculate fluxes depending on orientation with specific direction averages1034if orientation == 11035f1 = rho_mean * v1_avg1036f2 = f1 * v1_avg + p_mean + 0.5f0 * mag_norm_avg - B1_avg * B1_avg1037f3 = f1 * v2_avg - B1_avg * B2_avg1038f4 = f1 * v3_avg - B1_avg * B3_avg1039f6 = equations.c_h * psi_avg1040f7 = v1_avg * B2_avg - v2_avg * B1_avg1041f8 = v1_avg * B3_avg - v3_avg * B1_avg1042f9 = equations.c_h * B1_avg1043# total energy flux is complicated and involves the previous eight components1044psi_B1_avg = 0.5f0 * (B1_ll * psi_ll + B1_rr * psi_rr)1045v1_mag_avg = 0.5f0 * (v1_ll * mag_norm_ll + v1_rr * mag_norm_rr)1046f5 = (f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +1047f2 * v1_avg + f3 * v2_avg +1048f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg + f9 * psi_avg -10490.5f0 * v1_mag_avg +1050B1_avg * vel_dot_mag_avg - equations.c_h * psi_B1_avg)1051else1052f1 = rho_mean * v2_avg1053f2 = f1 * v1_avg - B1_avg * B2_avg1054f3 = f1 * v2_avg + p_mean + 0.5f0 * mag_norm_avg - B2_avg * B2_avg1055f4 = f1 * v3_avg - B2_avg * B3_avg1056f6 = v2_avg * B1_avg - v1_avg * B2_avg1057f7 = equations.c_h * psi_avg1058f8 = v2_avg * B3_avg - v3_avg * B2_avg1059f9 = equations.c_h * B2_avg1060# total energy flux is complicated and involves the previous eight components1061psi_B2_avg = 0.5f0 * (B2_ll * psi_ll + B2_rr * psi_rr)1062v2_mag_avg = 0.5f0 * (v2_ll * mag_norm_ll + v2_rr * mag_norm_rr)1063f5 = (f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +1064f2 * v1_avg + f3 * v2_avg +1065f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg + f9 * psi_avg -10660.5f0 * v2_mag_avg +1067B2_avg * vel_dot_mag_avg - equations.c_h * psi_B2_avg)1068end10691070return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)1071end10721073"""1074flux_hindenlang_gassner(u_ll, u_rr, orientation_or_normal_direction,1075equations::IdealGlmMhdEquations2D)10761077Entropy conserving and kinetic energy preserving two-point flux of1078Hindenlang and Gassner (2019), extending [`flux_ranocha`](@ref) to the MHD equations.10791080## References1081- Florian Hindenlang, Gregor Gassner (2019)1082A new entropy conservative two-point flux for ideal MHD equations derived from1083first principles.1084Presented at HONOM 2019: European workshop on high order numerical methods1085for evolutionary PDEs, theory and applications1086- Hendrik Ranocha (2018)1087Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods1088for Hyperbolic Balance Laws1089[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)1090- Hendrik Ranocha (2020)1091Entropy Conserving and Kinetic Energy Preserving Numerical Methods for1092the Euler Equations Using Summation-by-Parts Operators1093[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)1094"""1095@inline function flux_hindenlang_gassner(u_ll, u_rr, orientation::Integer,1096equations::IdealGlmMhdEquations2D)1097# Unpack left and right states1098rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll, psi_ll = cons2prim(u_ll,1099equations)1100rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr, psi_rr = cons2prim(u_rr,1101equations)11021103# Compute the necessary mean values needed for either direction1104rho_mean = ln_mean(rho_ll, rho_rr)1105# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`1106# in exact arithmetic since1107# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)1108# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)1109inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)1110v1_avg = 0.5f0 * (v1_ll + v1_rr)1111v2_avg = 0.5f0 * (v2_ll + v2_rr)1112v3_avg = 0.5f0 * (v3_ll + v3_rr)1113p_avg = 0.5f0 * (p_ll + p_rr)1114psi_avg = 0.5f0 * (psi_ll + psi_rr)1115velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)1116magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)11171118# Calculate fluxes depending on orientation with specific direction averages1119if orientation == 11120f1 = rho_mean * v1_avg1121f2 = f1 * v1_avg + p_avg + magnetic_square_avg -11220.5f0 * (B1_ll * B1_rr + B1_rr * B1_ll)1123f3 = f1 * v2_avg - 0.5f0 * (B1_ll * B2_rr + B1_rr * B2_ll)1124f4 = f1 * v3_avg - 0.5f0 * (B1_ll * B3_rr + B1_rr * B3_ll)1125#f5 below1126f6 = equations.c_h * psi_avg1127f7 = 0.5f0 * (v1_ll * B2_ll - v2_ll * B1_ll + v1_rr * B2_rr - v2_rr * B1_rr)1128f8 = 0.5f0 * (v1_ll * B3_ll - v3_ll * B1_ll + v1_rr * B3_rr - v3_rr * B1_rr)1129f9 = equations.c_h * 0.5f0 * (B1_ll + B1_rr)1130# total energy flux is complicated and involves the previous components1131f5 = (f1 *1132(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)1133+11340.5f0 * (+p_ll * v1_rr + p_rr * v1_ll1135+ (v1_ll * B2_ll * B2_rr + v1_rr * B2_rr * B2_ll)1136+ (v1_ll * B3_ll * B3_rr + v1_rr * B3_rr * B3_ll)1137-1138(v2_ll * B1_ll * B2_rr + v2_rr * B1_rr * B2_ll)1139-1140(v3_ll * B1_ll * B3_rr + v3_rr * B1_rr * B3_ll)1141+1142equations.c_h * (B1_ll * psi_rr + B1_rr * psi_ll)))1143else # orientation == 21144f1 = rho_mean * v2_avg1145f2 = f1 * v1_avg - 0.5f0 * (B2_ll * B1_rr + B2_rr * B1_ll)1146f3 = f1 * v2_avg + p_avg + magnetic_square_avg -11470.5f0 * (B2_ll * B2_rr + B2_rr * B2_ll)1148f4 = f1 * v3_avg - 0.5f0 * (B2_ll * B3_rr + B2_rr * B3_ll)1149#f5 below1150f6 = 0.5f0 * (v2_ll * B1_ll - v1_ll * B2_ll + v2_rr * B1_rr - v1_rr * B2_rr)1151f7 = equations.c_h * psi_avg1152f8 = 0.5f0 * (v2_ll * B3_ll - v3_ll * B2_ll + v2_rr * B3_rr - v3_rr * B2_rr)1153f9 = equations.c_h * 0.5f0 * (B2_ll + B2_rr)1154# total energy flux is complicated and involves the previous components1155f5 = (f1 *1156(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)1157+11580.5f0 * (+p_ll * v2_rr + p_rr * v2_ll1159+ (v2_ll * B1_ll * B1_rr + v2_rr * B1_rr * B1_ll)1160+ (v2_ll * B3_ll * B3_rr + v2_rr * B3_rr * B3_ll)1161-1162(v1_ll * B2_ll * B1_rr + v1_rr * B2_rr * B1_ll)1163-1164(v3_ll * B2_ll * B3_rr + v3_rr * B2_rr * B3_ll)1165+1166equations.c_h * (B2_ll * psi_rr + B2_rr * psi_ll)))1167end11681169return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)1170end11711172@inline function flux_hindenlang_gassner(u_ll, u_rr, normal_direction::AbstractVector,1173equations::IdealGlmMhdEquations2D)1174# Unpack left and right states1175rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll, psi_ll = cons2prim(u_ll,1176equations)1177rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr, psi_rr = cons2prim(u_rr,1178equations)1179v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]1180v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]1181B_dot_n_ll = B1_ll * normal_direction[1] + B2_ll * normal_direction[2]1182B_dot_n_rr = B1_rr * normal_direction[1] + B2_rr * normal_direction[2]11831184# Compute the necessary mean values needed for either direction1185rho_mean = ln_mean(rho_ll, rho_rr)1186# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`1187# in exact arithmetic since1188# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)1189# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)1190inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)1191v1_avg = 0.5f0 * (v1_ll + v1_rr)1192v2_avg = 0.5f0 * (v2_ll + v2_rr)1193v3_avg = 0.5f0 * (v3_ll + v3_rr)1194p_avg = 0.5f0 * (p_ll + p_rr)1195psi_avg = 0.5f0 * (psi_ll + psi_rr)1196velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)1197magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)11981199# Calculate fluxes depending on normal_direction1200f1 = rho_mean * 0.5f0 * (v_dot_n_ll + v_dot_n_rr)1201f2 = (f1 * v1_avg + (p_avg + magnetic_square_avg) * normal_direction[1]1202-12030.5f0 * (B_dot_n_ll * B1_rr + B_dot_n_rr * B1_ll))1204f3 = (f1 * v2_avg + (p_avg + magnetic_square_avg) * normal_direction[2]1205-12060.5f0 * (B_dot_n_ll * B2_rr + B_dot_n_rr * B2_ll))1207f4 = (f1 * v3_avg1208-12090.5f0 * (B_dot_n_ll * B3_rr + B_dot_n_rr * B3_ll))1210#f5 below1211f6 = (equations.c_h * psi_avg * normal_direction[1]1212+12130.5f0 * (v_dot_n_ll * B1_ll - v1_ll * B_dot_n_ll +1214v_dot_n_rr * B1_rr - v1_rr * B_dot_n_rr))1215f7 = (equations.c_h * psi_avg * normal_direction[2]1216+12170.5f0 * (v_dot_n_ll * B2_ll - v2_ll * B_dot_n_ll +1218v_dot_n_rr * B2_rr - v2_rr * B_dot_n_rr))1219f8 = +0.5f0 * (v_dot_n_ll * B3_ll - v3_ll * B_dot_n_ll +1220v_dot_n_rr * B3_rr - v3_rr * B_dot_n_rr)1221f9 = equations.c_h * 0.5f0 * (B_dot_n_ll + B_dot_n_rr)1222# total energy flux is complicated and involves the previous components1223f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)1224+12250.5f0 * (+p_ll * v_dot_n_rr + p_rr * v_dot_n_ll1226+ (v_dot_n_ll * B1_ll * B1_rr + v_dot_n_rr * B1_rr * B1_ll)1227+ (v_dot_n_ll * B2_ll * B2_rr + v_dot_n_rr * B2_rr * B2_ll)1228+ (v_dot_n_ll * B3_ll * B3_rr + v_dot_n_rr * B3_rr * B3_ll)1229-1230(v1_ll * B_dot_n_ll * B1_rr + v1_rr * B_dot_n_rr * B1_ll)1231-1232(v2_ll * B_dot_n_ll * B2_rr + v2_rr * B_dot_n_rr * B2_ll)1233-1234(v3_ll * B_dot_n_ll * B3_rr + v3_rr * B_dot_n_rr * B3_ll)1235+1236equations.c_h * (B_dot_n_ll * psi_rr + B_dot_n_rr * psi_ll)))12371238return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)1239end12401241# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation1242@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,1243equations::IdealGlmMhdEquations2D)1244rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll1245rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr12461247# Calculate the left/right velocities and fast magnetoacoustic wave speeds1248if orientation == 11249v_ll = rho_v1_ll / rho_ll1250v_rr = rho_v1_rr / rho_rr1251else # orientation == 21252v_ll = rho_v2_ll / rho_ll1253v_rr = rho_v2_rr / rho_rr1254end1255cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)1256cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)12571258return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)1259end12601261@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,1262equations::IdealGlmMhdEquations2D)1263# return max(v_mag_ll, v_mag_rr) + max(cf_ll, cf_rr)1264rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll1265rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr12661267# Calculate normal velocities and fast magnetoacoustic wave speeds1268# left1269v1_ll = rho_v1_ll / rho_ll1270v2_ll = rho_v2_ll / rho_ll1271v_ll = (v1_ll * normal_direction[1]1272+1273v2_ll * normal_direction[2])1274cf_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)1275# right1276v1_rr = rho_v1_rr / rho_rr1277v2_rr = rho_v2_rr / rho_rr1278v_rr = (v1_rr * normal_direction[1]1279+1280v2_rr * normal_direction[2])1281cf_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)12821283# wave speeds already scaled by norm(normal_direction) in [`calc_fast_wavespeed`](@ref)1284return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)1285end12861287# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`1288@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,1289equations::IdealGlmMhdEquations2D)1290rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll1291rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr12921293# Calculate the left/right velocities and fast magnetoacoustic wave speeds1294if orientation == 11295v_ll = rho_v1_ll / rho_ll1296v_rr = rho_v1_rr / rho_rr1297else # orientation == 21298v_ll = rho_v2_ll / rho_ll1299v_rr = rho_v2_rr / rho_rr1300end1301cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)1302cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)13031304return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)1305end13061307# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`1308@inline function max_abs_speed(u_ll, u_rr, normal_direction::AbstractVector,1309equations::IdealGlmMhdEquations2D)1310# return max(v_mag_ll, v_mag_rr) + max(cf_ll, cf_rr)1311rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll1312rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr13131314# Calculate normal velocities and fast magnetoacoustic wave speeds1315# left1316v1_ll = rho_v1_ll / rho_ll1317v2_ll = rho_v2_ll / rho_ll1318v_ll = (v1_ll * normal_direction[1]1319+1320v2_ll * normal_direction[2])1321cf_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)1322# right1323v1_rr = rho_v1_rr / rho_rr1324v2_rr = rho_v2_rr / rho_rr1325v_rr = (v1_rr * normal_direction[1]1326+1327v2_rr * normal_direction[2])1328cf_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)13291330# wave speeds already scaled by norm(normal_direction) in [`calc_fast_wavespeed`](@ref)1331return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)1332end13331334# Calculate estimate for minimum and maximum wave speeds for HLL-type fluxes1335@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,1336equations::IdealGlmMhdEquations2D)1337rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll1338rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr13391340# Calculate primitive velocity variables1341v1_ll = rho_v1_ll / rho_ll1342v2_ll = rho_v2_ll / rho_ll13431344v1_rr = rho_v1_rr / rho_rr1345v2_rr = rho_v2_rr / rho_rr13461347# Approximate the left-most and right-most eigenvalues in the Riemann fan1348if orientation == 1 # x-direction1349λ_min = v1_ll - calc_fast_wavespeed(u_ll, orientation, equations)1350λ_max = v1_rr + calc_fast_wavespeed(u_rr, orientation, equations)1351else # y-direction1352λ_min = v2_ll - calc_fast_wavespeed(u_ll, orientation, equations)1353λ_max = v2_rr + calc_fast_wavespeed(u_rr, orientation, equations)1354end13551356return λ_min, λ_max1357end13581359@inline function min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,1360equations::IdealGlmMhdEquations2D)1361rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll1362rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr13631364# Calculate primitive velocity variables1365v1_ll = rho_v1_ll / rho_ll1366v2_ll = rho_v2_ll / rho_ll13671368v1_rr = rho_v1_rr / rho_rr1369v2_rr = rho_v2_rr / rho_rr13701371v_normal_ll = (v1_ll * normal_direction[1] + v2_ll * normal_direction[2])1372v_normal_rr = (v1_rr * normal_direction[1] + v2_rr * normal_direction[2])13731374c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)1375c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)13761377# Estimate the min/max eigenvalues in the normal direction1378λ_min = min(v_normal_ll - c_f_ll, v_normal_rr - c_f_rr)1379λ_max = max(v_normal_rr + c_f_rr, v_normal_rr + c_f_rr)13801381return λ_min, λ_max1382end13831384# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes1385@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,1386equations::IdealGlmMhdEquations2D)1387rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll1388rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr13891390# Calculate primitive velocity variables1391v1_ll = rho_v1_ll / rho_ll1392v2_ll = rho_v2_ll / rho_ll13931394v1_rr = rho_v1_rr / rho_rr1395v2_rr = rho_v2_rr / rho_rr13961397# Approximate the left-most and right-most eigenvalues in the Riemann fan1398if orientation == 1 # x-direction1399c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)1400c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)14011402λ_min = min(v1_ll - c_f_ll, v1_rr - c_f_rr)1403λ_max = max(v1_ll + c_f_ll, v1_rr + c_f_rr)1404else # y-direction1405c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)1406c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)14071408λ_min = min(v2_ll - c_f_ll, v2_rr - c_f_rr)1409λ_max = max(v2_ll + c_f_ll, v1_rr + c_f_rr)1410end14111412return λ_min, λ_max1413end14141415# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes1416@inline function min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector,1417equations::IdealGlmMhdEquations2D)1418rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll1419rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr14201421# Calculate primitive velocity variables1422v1_ll = rho_v1_ll / rho_ll1423v2_ll = rho_v2_ll / rho_ll14241425v1_rr = rho_v1_rr / rho_rr1426v2_rr = rho_v2_rr / rho_rr14271428v_normal_ll = (v1_ll * normal_direction[1] + v2_ll * normal_direction[2])1429v_normal_rr = (v1_rr * normal_direction[1] + v2_rr * normal_direction[2])14301431c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)1432c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)14331434# Estimate the min/max eigenvalues in the normal direction1435λ_min = min(v_normal_ll - c_f_ll, v_normal_rr - c_f_rr)1436λ_max = max(v_normal_ll + c_f_ll, v_normal_rr + c_f_rr)14371438return λ_min, λ_max1439end14401441"""1442min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer, equations::IdealGlmMhdEquations2D)14431444Calculate minimum and maximum wave speeds for HLL-type fluxes as in1445- Li (2005)1446An HLLC Riemann solver for magneto-hydrodynamics1447[DOI: 10.1016/j.jcp.2004.08.020](https://doi.org/10.1016/j.jcp.2004.08.020).1448"""1449@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,1450equations::IdealGlmMhdEquations2D)1451rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll1452rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr14531454# Calculate primitive velocity variables1455v1_ll = rho_v1_ll / rho_ll1456v2_ll = rho_v2_ll / rho_ll14571458v1_rr = rho_v1_rr / rho_rr1459v2_rr = rho_v2_rr / rho_rr14601461# Approximate the left-most and right-most eigenvalues in the Riemann fan1462if orientation == 1 # x-direction1463c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)1464c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)1465vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)1466λ_min = min(v1_ll - c_f_ll, vel_roe - c_f_roe)1467λ_max = max(v1_rr + c_f_rr, vel_roe + c_f_roe)1468else # y-direction1469c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)1470c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)1471vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)1472λ_min = min(v2_ll - c_f_ll, vel_roe - c_f_roe)1473λ_max = max(v2_rr + c_f_rr, vel_roe + c_f_roe)1474end14751476return λ_min, λ_max1477end14781479@inline function min_max_speed_einfeldt(u_ll, u_rr, normal_direction::AbstractVector,1480equations::IdealGlmMhdEquations2D)1481rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll1482rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr14831484# Calculate primitive velocity variables1485v1_ll = rho_v1_ll / rho_ll1486v2_ll = rho_v2_ll / rho_ll14871488v1_rr = rho_v1_rr / rho_rr1489v2_rr = rho_v2_rr / rho_rr14901491v_normal_ll = (v1_ll * normal_direction[1] + v2_ll * normal_direction[2])1492v_normal_rr = (v1_rr * normal_direction[1] + v2_rr * normal_direction[2])14931494c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)1495c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)1496v_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, normal_direction, equations)14971498# Estimate the min/max eigenvalues in the normal direction1499λ_min = min(v_normal_ll - c_f_ll, v_roe - c_f_roe)1500λ_max = max(v_normal_rr + c_f_rr, v_roe + c_f_roe)15011502return λ_min, λ_max1503end15041505# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction1506# has been normalized prior to this rotation of the state vector1507@inline function rotate_to_x(u, normal_vector, equations::IdealGlmMhdEquations2D)1508# cos and sin of the angle between the x-axis and the normalized normal_vector are1509# the normalized vector's x and y coordinates respectively (see unit circle).1510c = normal_vector[1]1511s = normal_vector[2]15121513# Apply the 2D rotation matrix with normal and tangent directions of the form1514# [ 1 0 0 0 0 0 0 0 0;1515# 0 n_1 n_2 0 0 0 0 0 0;1516# 0 t_1 t_2 0 0 0 0 0 0;1517# 0 0 0 1 0 0 0 0 0;1518# 0 0 0 0 1 0 0 0 0;1519# 0 0 0 0 0 n_1 n_2 0 0;1520# 0 0 0 0 0 t_1 t_2 0 0;1521# 0 0 0 0 0 0 0 1 0;1522# 0 0 0 0 0 0 0 0 1 ]1523# where t_1 = -n_2 and t_2 = n_1.1524# Note for IdealGlmMhdEquations2D only the velocities and magnetic field variables rotate15251526return SVector(u[1],1527c * u[2] + s * u[3],1528-s * u[2] + c * u[3],1529u[4],1530u[5],1531c * u[6] + s * u[7],1532-s * u[6] + c * u[7],1533u[8],1534u[9])1535end15361537# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction1538# has been normalized prior to this back-rotation of the state vector1539@inline function rotate_from_x(u, normal_vector, equations::IdealGlmMhdEquations2D)1540# cos and sin of the angle between the x-axis and the normalized normal_vector are1541# the normalized vector's x and y coordinates respectively (see unit circle).1542c = normal_vector[1]1543s = normal_vector[2]15441545# Apply the 2D back-rotation matrix with normal and tangent directions of the form1546# [ 1 0 0 0 0 0 0 0 0;1547# 0 n_1 t_1 0 0 0 0 0 0;1548# 0 n_2 t_2 0 0 0 0 0 0;1549# 0 0 0 1 0 0 0 0 0;1550# 0 0 0 0 1 0 0 0 0;1551# 0 0 0 0 0 n_1 t_1 0 0;1552# 0 0 0 0 0 n_2 t_2 0 0;1553# 0 0 0 0 0 0 0 1 0;1554# 0 0 0 0 0 0 0 0 1 ]1555# where t_1 = -n_2 and t_2 = n_1.1556# Note for IdealGlmMhdEquations2D the velocities and magnetic field variables back-rotate15571558return SVector(u[1],1559c * u[2] - s * u[3],1560s * u[2] + c * u[3],1561u[4],1562u[5],1563c * u[6] - s * u[7],1564s * u[6] + c * u[7],1565u[8],1566u[9])1567end15681569@inline function max_abs_speeds(u, equations::IdealGlmMhdEquations2D)1570rho, rho_v1, rho_v2, rho_v3, _ = u1571v1 = rho_v1 / rho1572v2 = rho_v2 / rho1573v3 = rho_v3 / rho1574cf_x_direction = calc_fast_wavespeed(u, 1, equations)1575cf_y_direction = calc_fast_wavespeed(u, 2, equations)15761577return abs(v1) + cf_x_direction, abs(v2) + cf_y_direction1578end15791580# Convert conservative variables to primitive1581@inline function cons2prim(u, equations::IdealGlmMhdEquations2D)1582rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u15831584v1 = rho_v1 / rho1585v2 = rho_v2 / rho1586v3 = rho_v3 / rho1587p = (equations.gamma - 1) * (rho_e_total -15880.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v31589+ B1 * B1 + B2 * B2 + B3 * B31590+ psi * psi))15911592return SVector(rho, v1, v2, v3, p, B1, B2, B3, psi)1593end15941595# Convert conservative variables to entropy variables1596@inline function cons2entropy(u, equations::IdealGlmMhdEquations2D)1597rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u15981599v1 = rho_v1 / rho1600v2 = rho_v2 / rho1601v3 = rho_v3 / rho1602v_square = v1^2 + v2^2 + v3^21603p = (equations.gamma - 1) *1604(rho_e_total - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2) -16050.5f0 * psi^2)1606s = log(p) - equations.gamma * log(rho)1607rho_p = rho / p16081609w1 = (equations.gamma - s) * equations.inv_gamma_minus_one -16100.5f0 * rho_p * v_square1611w2 = rho_p * v11612w3 = rho_p * v21613w4 = rho_p * v31614w5 = -rho_p1615w6 = rho_p * B11616w7 = rho_p * B21617w8 = rho_p * B31618w9 = rho_p * psi16191620return SVector(w1, w2, w3, w4, w5, w6, w7, w8, w9)1621end16221623# Convert entropy variables to conservative variables1624@inline function entropy2cons(w, equations::IdealGlmMhdEquations2D)1625w1, w2, w3, w4, w5, w6, w7, w8, w9 = w16261627v1 = -w2 / w51628v2 = -w3 / w51629v3 = -w4 / w516301631B1 = -w6 / w51632B2 = -w7 / w51633B3 = -w8 / w51634psi = -w9 / w516351636# This imitates what is done for compressible Euler 3D `entropy2cons`: we convert from1637# the entropy variables for `-rho * s / (gamma - 1)` to the entropy variables for the entropy1638# `-rho * s` used by Hughes, Franca, Mallet (1986).1639@unpack gamma = equations1640V1, V2, V3, V4, V5 = SVector(w1, w2, w3, w4, w5) * (gamma - 1)1641s = gamma - V1 + (V2^2 + V3^2 + V4^2) / (2 * V5)1642rho_iota = ((gamma - 1) / (-V5)^gamma)^(equations.inv_gamma_minus_one) *1643exp(-s * equations.inv_gamma_minus_one)1644rho = -rho_iota * V51645p = -rho / w516461647return prim2cons(SVector(rho, v1, v2, v3, p, B1, B2, B3, psi), equations)1648end16491650# Convert primitive to conservative variables1651@inline function prim2cons(prim, equations::IdealGlmMhdEquations2D)1652rho, v1, v2, v3, p, B1, B2, B3, psi = prim16531654rho_v1 = rho * v11655rho_v2 = rho * v21656rho_v3 = rho * v31657rho_e_total = p * equations.inv_gamma_minus_one +16580.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3) +16590.5f0 * (B1^2 + B2^2 + B3^2) + 0.5f0 * psi^216601661return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi)1662end16631664@inline function density(u, equations::IdealGlmMhdEquations2D)1665rho = u[1]1666return rho1667end16681669@inline function velocity(u, equations::IdealGlmMhdEquations2D)1670rho = u[1]1671v1 = u[2] / rho1672v2 = u[3] / rho1673v3 = u[4] / rho1674return SVector(v1, v2, v3)1675end16761677@inline function velocity(u, orientation::Int, equations::IdealGlmMhdEquations2D)1678rho = u[1]1679v = u[orientation + 1] / rho1680return v1681end16821683@doc raw"""1684pressure(u, equations::IdealGlmMhdEquations2D)16851686Computes the pressure for an ideal equation of state with1687isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`.1688```math1689\begin{aligned}1690p &= (\gamma - 1) \left( \rho e_{\text{total}} - \rho e_{\text{kinetic}} - \rho e_{\text{magnetic}} - \frac{1}{2}\psi^2 \right) \\1691&= (\gamma - 1) \left( \rho e - \frac{1}{2}1692\left[\rho \Vert v \Vert_2^2 + \Vert B \Vert_2^2 + \psi^2 \right] \right)1693\end{aligned}1694```1695"""1696@inline function pressure(u, equations::IdealGlmMhdEquations2D)1697rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u1698p = (equations.gamma - 1) * (rho_e_total -16990.5f0 *1700((rho_v1^2 + rho_v2^2 + rho_v3^2) / rho +1701B1^2 + B2^2 + B3^2 + psi^2))1702return p1703end17041705# Transformation from conservative variables u to d(p)/d(u)1706@inline function gradient_conservative(::typeof(pressure),1707u, equations::IdealGlmMhdEquations2D)1708rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u17091710v1 = rho_v1 / rho1711v2 = rho_v2 / rho1712v3 = rho_v3 / rho1713v_square = v1^2 + v2^2 + v3^217141715return (equations.gamma - 1) *1716SVector(0.5f0 * v_square, -v1, -v2, -v3, 1, -B1, -B2, -B3, -psi)1717end17181719@inline function density_pressure(u, equations::IdealGlmMhdEquations2D)1720rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u1721rho_times_p = (equations.gamma - 1) * (rho * rho_e_total -17220.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2 +1723rho * (B1^2 + B2^2 + B3^2 + psi^2)))1724return rho_times_p1725end17261727# Compute the fastest wave speed for ideal MHD equations: c_f, the fast magnetoacoustic eigenvalue1728@inline function calc_fast_wavespeed(cons, orientation::Integer,1729equations::IdealGlmMhdEquations2D)1730rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = cons1731v1 = rho_v1 / rho1732v2 = rho_v2 / rho1733v3 = rho_v3 / rho1734kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)1735mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)1736p = (equations.gamma - 1) * (rho_e_total - kin_en - mag_en - 0.5f0 * psi^2)1737a_square = equations.gamma * p / rho1738sqrt_rho = sqrt(rho)1739b1 = B1 / sqrt_rho1740b2 = B2 / sqrt_rho1741b3 = B3 / sqrt_rho1742b_square = b1 * b1 + b2 * b2 + b3 * b31743if orientation == 1 # x-direction1744c_f = sqrt(0.5f0 * (a_square + b_square) +17450.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b1^2))1746else1747c_f = sqrt(0.5f0 * (a_square + b_square) +17480.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b2^2))1749end1750return c_f1751end17521753@inline function calc_fast_wavespeed(cons, normal_direction::AbstractVector,1754equations::IdealGlmMhdEquations2D)1755rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = cons1756v1 = rho_v1 / rho1757v2 = rho_v2 / rho1758v3 = rho_v3 / rho1759kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)1760mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)1761p = (equations.gamma - 1) * (rho_e_total - kin_en - mag_en - 0.5f0 * psi^2)1762a_square = equations.gamma * p / rho1763sqrt_rho = sqrt(rho)1764b1 = B1 / sqrt_rho1765b2 = B2 / sqrt_rho1766b3 = B3 / sqrt_rho1767b_square = b1 * b1 + b2 * b2 + b3 * b31768norm_squared = (normal_direction[1] * normal_direction[1] +1769normal_direction[2] * normal_direction[2])1770b_dot_n_squared = (b1 * normal_direction[1] +1771b2 * normal_direction[2])^2 / norm_squared17721773c_f = sqrt((0.5f0 * (a_square + b_square) +17740.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b_dot_n_squared)) *1775norm_squared)1776return c_f1777end17781779"""1780calc_fast_wavespeed_roe(u_ll, u_rr, orientation_or_normal_direction, equations::IdealGlmMhdEquations2D)17811782Compute the fast magnetoacoustic wave speed using Roe averages1783as given by1784- Cargo and Gallice (1997)1785Roe Matrices for Ideal MHD and Systematic Construction1786of Roe Matrices for Systems of Conservation Laws1787[DOI: 10.1006/jcph.1997.5773](https://doi.org/10.1006/jcph.1997.5773)1788"""1789@inline function calc_fast_wavespeed_roe(u_ll, u_rr, orientation::Integer,1790equations::IdealGlmMhdEquations2D)1791rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll1792rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr17931794# Calculate primitive variables1795v1_ll = rho_v1_ll / rho_ll1796v2_ll = rho_v2_ll / rho_ll1797v3_ll = rho_v3_ll / rho_ll1798kin_en_ll = 0.5f0 * (rho_v1_ll * v1_ll + rho_v2_ll * v2_ll + rho_v3_ll * v3_ll)1799mag_norm_ll = B1_ll * B1_ll + B2_ll * B2_ll + B3_ll * B3_ll1800p_ll = (equations.gamma - 1) *1801(rho_e_total_ll - kin_en_ll - 0.5f0 * mag_norm_ll - 0.5f0 * psi_ll^2)18021803v1_rr = rho_v1_rr / rho_rr1804v2_rr = rho_v2_rr / rho_rr1805v3_rr = rho_v3_rr / rho_rr1806kin_en_rr = 0.5f0 * (rho_v1_rr * v1_rr + rho_v2_rr * v2_rr + rho_v3_rr * v3_rr)1807mag_norm_rr = B1_rr * B1_rr + B2_rr * B2_rr + B3_rr * B3_rr1808p_rr = (equations.gamma - 1) *1809(rho_e_total_rr - kin_en_rr - 0.5f0 * mag_norm_rr - 0.5f0 * psi_rr^2)18101811# compute total pressure which is thermal + magnetic pressures1812p_total_ll = p_ll + 0.5f0 * mag_norm_ll1813p_total_rr = p_rr + 0.5f0 * mag_norm_rr18141815# compute the Roe density averages1816sqrt_rho_ll = sqrt(rho_ll)1817sqrt_rho_rr = sqrt(rho_rr)1818inv_sqrt_rho_add = 1 / (sqrt_rho_ll + sqrt_rho_rr)1819inv_sqrt_rho_prod = 1 / (sqrt_rho_ll * sqrt_rho_rr)1820rho_ll_roe = sqrt_rho_ll * inv_sqrt_rho_add1821rho_rr_roe = sqrt_rho_rr * inv_sqrt_rho_add1822# Roe averages1823# velocities and magnetic fields1824v1_roe = v1_ll * rho_ll_roe + v1_rr * rho_rr_roe1825v2_roe = v2_ll * rho_ll_roe + v2_rr * rho_rr_roe1826v3_roe = v3_ll * rho_ll_roe + v3_rr * rho_rr_roe1827B1_roe = B1_ll * rho_ll_roe + B1_rr * rho_rr_roe1828B2_roe = B2_ll * rho_ll_roe + B2_rr * rho_rr_roe1829B3_roe = B3_ll * rho_ll_roe + B3_rr * rho_rr_roe1830# enthalpy1831H_ll = (rho_e_total_ll + p_total_ll) / rho_ll1832H_rr = (rho_e_total_rr + p_total_rr) / rho_rr1833H_roe = H_ll * rho_ll_roe + H_rr * rho_rr_roe1834# temporary variable see equation (4.12) in Cargo and Gallice1835X = 0.5f0 * ((B1_ll - B1_rr)^2 + (B2_ll - B2_rr)^2 + (B3_ll - B3_rr)^2) *1836inv_sqrt_rho_add^21837# averaged components needed to compute c_f, the fast magnetoacoustic wave speed1838b_square_roe = (B1_roe^2 + B2_roe^2 + B3_roe^2) * inv_sqrt_rho_prod # scaled magnectic sum1839a_square_roe = ((2 - equations.gamma) * X +1840(equations.gamma - 1) *1841(H_roe - 0.5f0 * (v1_roe^2 + v2_roe^2 + v3_roe^2) -1842b_square_roe)) # acoustic speed1843# finally compute the average wave speed and set the output velocity (depends on orientation)1844if orientation == 1 # x-direction1845c_a_roe = B1_roe^2 * inv_sqrt_rho_prod # (squared) Alfvén wave speed1846a_star_roe = sqrt((a_square_roe + b_square_roe)^2 -18474 * a_square_roe * c_a_roe)1848c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe))1849vel_out_roe = v1_roe1850else # y-direction1851c_a_roe = B2_roe^2 * inv_sqrt_rho_prod # (squared) Alfvén wave speed1852a_star_roe = sqrt((a_square_roe + b_square_roe)^2 -18534 * a_square_roe * c_a_roe)1854c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe))1855vel_out_roe = v2_roe1856end18571858return vel_out_roe, c_f_roe1859end18601861@inline function calc_fast_wavespeed_roe(u_ll, u_rr, normal_direction::AbstractVector,1862equations::IdealGlmMhdEquations2D)1863rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll1864rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr18651866# Calculate primitive variables1867v1_ll = rho_v1_ll / rho_ll1868v2_ll = rho_v2_ll / rho_ll1869v3_ll = rho_v3_ll / rho_ll1870kin_en_ll = 0.5f0 * (rho_v1_ll * v1_ll + rho_v2_ll * v2_ll + rho_v3_ll * v3_ll)1871mag_norm_ll = B1_ll * B1_ll + B2_ll * B2_ll + B3_ll * B3_ll1872p_ll = (equations.gamma - 1) *1873(rho_e_total_ll - kin_en_ll - 0.5f0 * mag_norm_ll - 0.5f0 * psi_ll^2)18741875v1_rr = rho_v1_rr / rho_rr1876v2_rr = rho_v2_rr / rho_rr1877v3_rr = rho_v3_rr / rho_rr1878kin_en_rr = 0.5f0 * (rho_v1_rr * v1_rr + rho_v2_rr * v2_rr + rho_v3_rr * v3_rr)1879mag_norm_rr = B1_rr * B1_rr + B2_rr * B2_rr + B3_rr * B3_rr1880p_rr = (equations.gamma - 1) *1881(rho_e_total_rr - kin_en_rr - 0.5f0 * mag_norm_rr - 0.5f0 * psi_rr^2)18821883# compute total pressure which is thermal + magnetic pressures1884p_total_ll = p_ll + 0.5f0 * mag_norm_ll1885p_total_rr = p_rr + 0.5f0 * mag_norm_rr18861887# compute the Roe density averages1888sqrt_rho_ll = sqrt(rho_ll)1889sqrt_rho_rr = sqrt(rho_rr)1890inv_sqrt_rho_add = 1 / (sqrt_rho_ll + sqrt_rho_rr)1891inv_sqrt_rho_prod = 1 / (sqrt_rho_ll * sqrt_rho_rr)1892rho_ll_roe = sqrt_rho_ll * inv_sqrt_rho_add1893rho_rr_roe = sqrt_rho_rr * inv_sqrt_rho_add1894# Roe averages1895# velocities and magnetic fields1896v1_roe = v1_ll * rho_ll_roe + v1_rr * rho_rr_roe1897v2_roe = v2_ll * rho_ll_roe + v2_rr * rho_rr_roe1898v3_roe = v3_ll * rho_ll_roe + v3_rr * rho_rr_roe1899B1_roe = B1_ll * rho_ll_roe + B1_rr * rho_rr_roe1900B2_roe = B2_ll * rho_ll_roe + B2_rr * rho_rr_roe1901B3_roe = B3_ll * rho_ll_roe + B3_rr * rho_rr_roe1902# enthalpy1903H_ll = (rho_e_total_ll + p_total_ll) / rho_ll1904H_rr = (rho_e_total_rr + p_total_rr) / rho_rr1905H_roe = H_ll * rho_ll_roe + H_rr * rho_rr_roe1906# temporary variable see equation (4.12) in Cargo and Gallice1907X = 0.5f0 * ((B1_ll - B1_rr)^2 + (B2_ll - B2_rr)^2 + (B3_ll - B3_rr)^2) *1908inv_sqrt_rho_add^21909# averaged components needed to compute c_f, the fast magnetoacoustic wave speed1910b_square_roe = (B1_roe^2 + B2_roe^2 + B3_roe^2) * inv_sqrt_rho_prod # scaled magnectic sum1911a_square_roe = ((2 - equations.gamma) * X +1912(equations.gamma - 1) *1913(H_roe - 0.5f0 * (v1_roe^2 + v2_roe^2 + v3_roe^2) -1914b_square_roe)) # acoustic speed19151916# finally compute the average wave speed and set the output velocity (depends on orientation)1917norm_squared = (normal_direction[1] * normal_direction[1] +1918normal_direction[2] * normal_direction[2])1919B_roe_dot_n_squared = (B1_roe * normal_direction[1] +1920B2_roe * normal_direction[2])^2 / norm_squared19211922c_a_roe = B_roe_dot_n_squared * inv_sqrt_rho_prod # (squared) Alfvén wave speed1923a_star_roe = sqrt((a_square_roe + b_square_roe)^2 - 4 * a_square_roe * c_a_roe)1924c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe) * norm_squared)1925vel_out_roe = (v1_roe * normal_direction[1] +1926v2_roe * normal_direction[2])19271928return vel_out_roe, c_f_roe1929end19301931# Calculate thermodynamic entropy for a conservative state `cons`1932@inline function entropy_thermodynamic(cons, equations::IdealGlmMhdEquations2D)1933# Pressure1934p = (equations.gamma - 1) *1935(cons[5] - 0.5f0 * (cons[2]^2 + cons[3]^2 + cons[4]^2) / cons[1]1936-19370.5f0 * (cons[6]^2 + cons[7]^2 + cons[8]^2)1938-19390.5f0 * cons[9]^2)19401941# Thermodynamic entropy1942s = log(p) - equations.gamma * log(cons[1])19431944return s1945end19461947# Calculate mathematical entropy for a conservative state `cons`1948@inline function entropy_math(cons, equations::IdealGlmMhdEquations2D)1949S = -entropy_thermodynamic(cons, equations) * cons[1] *1950equations.inv_gamma_minus_one19511952return S1953end19541955# Default entropy is the mathematical entropy1956@inline entropy(cons, equations::IdealGlmMhdEquations2D) = entropy_math(cons, equations)19571958# Calculate total energy for a conservative state `cons`1959@inline energy_total(cons, ::IdealGlmMhdEquations2D) = cons[5]19601961# Calculate kinetic energy for a conservative state `cons`1962@inline function energy_kinetic(cons, equations::IdealGlmMhdEquations2D)1963return 0.5f0 * (cons[2]^2 + cons[3]^2 + cons[4]^2) / cons[1]1964end19651966# Calculate the magnetic energy for a conservative state `cons`.1967# OBS! For non-dinmensional form of the ideal MHD magnetic pressure ≡ magnetic energy1968@inline function energy_magnetic(cons, ::IdealGlmMhdEquations2D)1969return 0.5f0 * (cons[6]^2 + cons[7]^2 + cons[8]^2)1970end19711972# Calculate internal energy for a conservative state `cons`1973@inline function energy_internal(cons, equations::IdealGlmMhdEquations2D)1974return (energy_total(cons, equations)1975-1976energy_kinetic(cons, equations)1977-1978energy_magnetic(cons, equations)1979-1980cons[9]^2 / 2)1981end19821983# State validation for Newton-bisection method of subcell IDP limiting1984@inline function Base.isvalid(u, equations::IdealGlmMhdEquations2D)1985if u[1] <= 0 || pressure(u, equations) <= 01986return false1987end1988return true1989end19901991# Calculate the cross helicity (\vec{v}⋅\vec{B}) for a conservative state `cons`1992@inline function cross_helicity(cons, ::IdealGlmMhdEquations2D)1993return (cons[2] * cons[6] + cons[3] * cons[7] + cons[4] * cons[8]) / cons[1]1994end1995end # @muladd199619971998