Path: blob/main/src/equations/ideal_glm_mhd_1d.jl
5586 views
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).1# Since these FMAs can increase the performance of many numerical algorithms,2# we need to opt-in explicitly.3# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.4@muladd begin5#! format: noindent67@doc raw"""8IdealGlmMhdEquations1D(gamma)910The ideal compressible GLM-MHD equations11```math12\frac{\partial}{\partial t}13\begin{pmatrix}14\rho \\ \rho v_1 \\ \rho v_2 \\ \rho v_3 \\ \rho e_{\text{total}} \\ B_1 \\ B_2 \\ B_315\end{pmatrix}16+17\frac{\partial}{\partial x}18\begin{pmatrix}19\rho v_1 \\ \rho v_1^2 + p + \frac{1}{2} \Vert \mathbf{B} \Vert_2 ^2 - B_1^2 \\ \rho v_1 v_2 - B_1 B_2 \\ \rho v_1 v_3 - B_1 B_320\\ (\frac{1}{2} \rho \Vert \mathbf{v} \Vert_2 ^2 + \frac{\gamma p}{\gamma - 1} + \Vert \mathbf{B} \Vert_2 ^2) v_1 - B_1 (\mathbf{v} \cdot \mathbf{B})21\\ 0 \\ v_1 B_2 - v_2 B_1 \\ v_1 B_3 - v_3 B_122\end{pmatrix}23=24\begin{pmatrix}250 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0 \\ 026\end{pmatrix}27```28for an ideal gas in one space dimension.29Here, ``\mathbf{v}`` is the velocity, ``\mathbf{B}`` the magnetic field, ``e_{\text{total}}`` the specific total energy, and30```math31p = (\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 \right)32```33the pressure and ``\gamma`` the heat capacity ratio.3435!!! note36There is no divergence cleaning variable `psi` because the divergence-free constraint37is satisfied trivially in one spatial dimension.38"""39struct IdealGlmMhdEquations1D{RealT <: Real} <: AbstractIdealGlmMhdEquations{1, 8}40gamma::RealT # ratio of specific heats41inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications4243function IdealGlmMhdEquations1D(gamma)44γ, inv_gamma_minus_one = promote(gamma, inv(gamma - 1))45return new{typeof(γ)}(γ, inv_gamma_minus_one)46end47end4849have_nonconservative_terms(::IdealGlmMhdEquations1D) = False()50function varnames(::typeof(cons2cons), ::IdealGlmMhdEquations1D)51return ("rho", "rho_v1", "rho_v2", "rho_v3", "rho_e_total", "B1", "B2", "B3")52end53function varnames(::typeof(cons2prim), ::IdealGlmMhdEquations1D)54return ("rho", "v1", "v2", "v3", "p", "B1", "B2", "B3")55end56function default_analysis_integrals(::IdealGlmMhdEquations1D)57return (entropy_timederivative, Val(:l2_divb), Val(:linf_divb))58end5960"""61initial_condition_constant(x, t, equations::IdealGlmMhdEquations1D)6263A constant initial condition to test free-stream preservation.64"""65function initial_condition_constant(x, t, equations::IdealGlmMhdEquations1D)66RealT = eltype(x)67rho = 168rho_v1 = convert(RealT, 0.1)69rho_v2 = -convert(RealT, 0.2)70rho_v3 = -0.5f071rho_e_total = 5072B1 = 373B2 = -convert(RealT, 1.2)74B3 = 0.5f075return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3)76end7778"""79initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations1D)8081An Alfvén wave as smooth initial condition used for convergence tests.82See for reference section 4.2 in83- Dominik Derigs, Andrew R. Winters, Gregor J. Gassner, and Stefanie Walch (2016)84A novel high-order, entropy stable, 3D AMR MHD solver with guaranteed positive pressure85[DOI: 10.1016/j.jcp.2016.04.048](https://doi.org/10.1016/j.jcp.2016.04.048)86"""87function initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations1D)88# smooth Alfvén wave test from Derigs et al. (2016)89# domain must be set to [0, 1], γ = 5/390RealT = eltype(x)91rho = 192v1 = 093si, co = sincospi(2 * (x[1] + t)) # Adding `t` makes non-integer time valid94v2 = convert(RealT, 0.1) * si95v3 = convert(RealT, 0.1) * co96p = convert(RealT, 0.1)97B1 = 198B2 = v299B3 = v3100return prim2cons(SVector(rho, v1, v2, v3, p, B1, B2, B3), equations)101end102103"""104initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations1D)105106A weak blast wave adapted from107- Sebastian Hennemann, Gregor J. Gassner (2020)108A provably entropy stable subcell shock capturing approach for high order split form DG109[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)110"""111function initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations1D)112# Adapted MHD version of the weak blast wave from Hennemann & Gassner JCP paper 2020 (Sec. 6.3)113# Same discontinuity in the velocities but with magnetic fields114# Set up polar coordinates115RealT = eltype(x)116inicenter = (0,)117x_norm = x[1] - inicenter[1]118r = sqrt(x_norm^2)119phi = atan(x_norm)120121# Calculate primitive variables122rho = r > 0.5f0 ? one(RealT) : convert(RealT, 1.1691)123v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos(phi)124p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)125126return prim2cons(SVector(rho, v1, 0, 0, p, 1, 1, 1, 0), equations)127end128129# Calculate 1D flux for a single point130@inline function flux(u, orientation::Integer, equations::IdealGlmMhdEquations1D)131rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u132v1 = rho_v1 / rho133v2 = rho_v2 / rho134v3 = rho_v3 / rho135kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)136mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)137p_over_gamma_minus_one = (rho_e_total - kin_en - mag_en)138p = (equations.gamma - 1) * p_over_gamma_minus_one139140# Ignore orientation since it is always "1" in 1D141f1 = rho_v1142f2 = rho_v1 * v1 + p + mag_en - B1^2143f3 = rho_v1 * v2 - B1 * B2144f4 = rho_v1 * v3 - B1 * B3145f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v1 -146B1 * (v1 * B1 + v2 * B2 + v3 * B3)147f6 = 0148f7 = v1 * B2 - v2 * B1149f8 = v1 * B3 - v3 * B1150151return SVector(f1, f2, f3, f4, f5, f6, f7, f8)152end153154"""155flux_derigs_etal(u_ll, u_rr, orientation, equations::IdealGlmMhdEquations1D)156157Entropy conserving two-point flux by158- Derigs et al. (2018)159Ideal GLM-MHD: About the entropy consistent nine-wave magnetic field160divergence diminishing ideal magnetohydrodynamics equations161[DOI: 10.1016/j.jcp.2018.03.002](https://doi.org/10.1016/j.jcp.2018.03.002)162"""163function flux_derigs_etal(u_ll, u_rr, orientation::Integer,164equations::IdealGlmMhdEquations1D)165# Unpack left and right states to get velocities, pressure, and inverse temperature (called beta)166rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll = u_ll167rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr = u_rr168169v1_ll = rho_v1_ll / rho_ll170v2_ll = rho_v2_ll / rho_ll171v3_ll = rho_v3_ll / rho_ll172v1_rr = rho_v1_rr / rho_rr173v2_rr = rho_v2_rr / rho_rr174v3_rr = rho_v3_rr / rho_rr175vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2176vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2177mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^2178mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^2179p_ll = (equations.gamma - 1) *180(rho_e_total_ll - 0.5f0 * rho_ll * vel_norm_ll - 0.5f0 * mag_norm_ll)181p_rr = (equations.gamma - 1) *182(rho_e_total_rr - 0.5f0 * rho_rr * vel_norm_rr - 0.5f0 * mag_norm_rr)183beta_ll = 0.5f0 * rho_ll / p_ll184beta_rr = 0.5f0 * rho_rr / p_rr185# for convenience store v⋅B186vel_dot_mag_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll187vel_dot_mag_rr = v1_rr * B1_rr + v2_rr * B2_rr + v3_rr * B3_rr188189# Compute the necessary mean values needed for either direction190rho_avg = 0.5f0 * (rho_ll + rho_rr)191rho_mean = ln_mean(rho_ll, rho_rr)192beta_mean = ln_mean(beta_ll, beta_rr)193beta_avg = 0.5f0 * (beta_ll + beta_rr)194v1_avg = 0.5f0 * (v1_ll + v1_rr)195v2_avg = 0.5f0 * (v2_ll + v2_rr)196v3_avg = 0.5f0 * (v3_ll + v3_rr)197p_mean = 0.5f0 * rho_avg / beta_avg198B1_avg = 0.5f0 * (B1_ll + B1_rr)199B2_avg = 0.5f0 * (B2_ll + B2_rr)200B3_avg = 0.5f0 * (B3_ll + B3_rr)201vel_norm_avg = 0.5f0 * (vel_norm_ll + vel_norm_rr)202mag_norm_avg = 0.5f0 * (mag_norm_ll + mag_norm_rr)203vel_dot_mag_avg = 0.5f0 * (vel_dot_mag_ll + vel_dot_mag_rr)204205# Ignore orientation since it is always "1" in 1D206f1 = rho_mean * v1_avg207f2 = f1 * v1_avg + p_mean + 0.5f0 * mag_norm_avg - B1_avg * B1_avg208f3 = f1 * v2_avg - B1_avg * B2_avg209f4 = f1 * v3_avg - B1_avg * B3_avg210f6 = 0211f7 = v1_avg * B2_avg - v2_avg * B1_avg212f8 = v1_avg * B3_avg - v3_avg * B1_avg213# total energy flux is complicated and involves the previous eight components214v1_mag_avg = 0.5f0 * (v1_ll * mag_norm_ll + v1_rr * mag_norm_rr)215f5 = (f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +216f2 * v1_avg + f3 * v2_avg +217f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg - 0.5f0 * v1_mag_avg +218B1_avg * vel_dot_mag_avg)219220return SVector(f1, f2, f3, f4, f5, f6, f7, f8)221end222223"""224flux_hindenlang_gassner(u_ll, u_rr, orientation_or_normal_direction,225equations::IdealGlmMhdEquations1D)226227Entropy conserving and kinetic energy preserving two-point flux of228Hindenlang and Gassner (2019), extending [`flux_ranocha`](@ref) to the MHD equations.229230## References231- Florian Hindenlang, Gregor Gassner (2019)232A new entropy conservative two-point flux for ideal MHD equations derived from233first principles.234Presented at HONOM 2019: European workshop on high order numerical methods235for evolutionary PDEs, theory and applications236- Hendrik Ranocha (2018)237Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods238for Hyperbolic Balance Laws239[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)240- Hendrik Ranocha (2020)241Entropy Conserving and Kinetic Energy Preserving Numerical Methods for242the Euler Equations Using Summation-by-Parts Operators243[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)244"""245@inline function flux_hindenlang_gassner(u_ll, u_rr, orientation::Integer,246equations::IdealGlmMhdEquations1D)247# Unpack left and right states248rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll = cons2prim(u_ll, equations)249rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr = cons2prim(u_rr, equations)250251# Compute the necessary mean values needed for either direction252rho_mean = ln_mean(rho_ll, rho_rr)253# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`254# in exact arithmetic since255# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)256# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)257inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)258v1_avg = 0.5f0 * (v1_ll + v1_rr)259v2_avg = 0.5f0 * (v2_ll + v2_rr)260v3_avg = 0.5f0 * (v3_ll + v3_rr)261p_avg = 0.5f0 * (p_ll + p_rr)262velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)263magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)264265# Calculate fluxes depending on orientation with specific direction averages266f1 = rho_mean * v1_avg267f2 = f1 * v1_avg + p_avg + magnetic_square_avg -2680.5f0 * (B1_ll * B1_rr + B1_rr * B1_ll)269f3 = f1 * v2_avg - 0.5f0 * (B1_ll * B2_rr + B1_rr * B2_ll)270f4 = f1 * v3_avg - 0.5f0 * (B1_ll * B3_rr + B1_rr * B3_ll)271#f5 below272f6 = 0273f7 = 0.5f0 * (v1_ll * B2_ll - v2_ll * B1_ll + v1_rr * B2_rr - v2_rr * B1_rr)274f8 = 0.5f0 * (v1_ll * B3_ll - v3_ll * B1_ll + v1_rr * B3_rr - v3_rr * B1_rr)275# total energy flux is complicated and involves the previous components276f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)277+2780.5f0 * (+p_ll * v1_rr + p_rr * v1_ll279+ (v1_ll * B2_ll * B2_rr + v1_rr * B2_rr * B2_ll)280+ (v1_ll * B3_ll * B3_rr + v1_rr * B3_rr * B3_ll)281-282(v2_ll * B1_ll * B2_rr + v2_rr * B1_rr * B2_ll)283-284(v3_ll * B1_ll * B3_rr + v3_rr * B1_rr * B3_ll)))285286return SVector(f1, f2, f3, f4, f5, f6, f7, f8)287end288289"""290flux_hllc(u_ll, u_rr, orientation, equations::IdealGlmMhdEquations1D)291292- Li (2005)293An HLLC Riemann solver for magneto-hydrodynamics294[DOI: 10.1016/j.jcp.2004.08.020](https://doi.org/10.1016/j.jcp.2004.08.020).295"""296function flux_hllc(u_ll, u_rr, orientation::Integer,297equations::IdealGlmMhdEquations1D)298# Unpack left and right states299rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll = cons2prim(u_ll, equations)300rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr = cons2prim(u_rr, equations)301302# Total pressure, i.e., thermal + magnetic pressures (eq. (12))303p_tot_ll = p_ll + 0.5f0 * (B1_ll^2 + B2_ll^2 + B3_ll^2)304p_tot_rr = p_rr + 0.5f0 * (B1_rr^2 + B2_rr^2 + B3_rr^2)305306# Conserved variables307rho_v1_ll = u_ll[2]308rho_v2_ll = u_ll[3]309rho_v3_ll = u_ll[4]310311rho_v1_rr = u_rr[2]312rho_v2_rr = u_rr[3]313rho_v3_rr = u_rr[4]314315# Obtain left and right fluxes316f_ll = flux(u_ll, orientation, equations)317f_rr = flux(u_rr, orientation, equations)318319SsL, SsR = min_max_speed_einfeldt(u_ll, u_rr, orientation, equations)320sMu_L = SsL - v1_ll321sMu_R = SsR - v1_rr322if SsL >= 0323f1 = f_ll[1]324f2 = f_ll[2]325f3 = f_ll[3]326f4 = f_ll[4]327f5 = f_ll[5]328f6 = f_ll[6]329f7 = f_ll[7]330f8 = f_ll[8]331elseif SsR <= 0332f1 = f_rr[1]333f2 = f_rr[2]334f3 = f_rr[3]335f4 = f_rr[4]336f5 = f_rr[5]337f6 = f_rr[6]338f7 = f_rr[7]339f8 = f_rr[8]340else341# Compute the "HLLC-speed", eq. (14) from paper mentioned above342#=343SStar = (rho_rr * v1_rr * sMu_R - rho_ll * v1_ll * sMu_L + p_tot_ll - p_tot_rr - B1_ll^2 + B1_rr^2 ) /344(rho_rr * sMu_R - rho_ll * sMu_L)345=#346# Simplification for 1D: B1 is constant347SStar = (rho_rr * v1_rr * sMu_R - rho_ll * v1_ll * sMu_L + p_tot_ll - p_tot_rr) /348(rho_rr * sMu_R - rho_ll * sMu_L)349350Sdiff = SsR - SsL351352# Compute HLL values for vStar, BStar353# These correspond to eq. (28) and (30) from the referenced paper354# and the classic HLL intermediate state given by (2)355rho_HLL = (SsR * rho_rr - SsL * rho_ll - (f_rr[1] - f_ll[1])) / Sdiff356357v1Star = (SsR * rho_v1_rr - SsL * rho_v1_ll - (f_rr[2] - f_ll[2])) /358(Sdiff * rho_HLL)359v2Star = (SsR * rho_v2_rr - SsL * rho_v2_ll - (f_rr[3] - f_ll[3])) /360(Sdiff * rho_HLL)361v3Star = (SsR * rho_v3_rr - SsL * rho_v3_ll - (f_rr[4] - f_ll[4])) /362(Sdiff * rho_HLL)363364#B1Star = (SsR * B1_rr - SsL * B1_ll - (f_rr[6] - f_ll[6])) / Sdiff365# 1D B1 = constant => B1_ll = B1_rr = B1Star366B1Star = B1_ll367368B2Star = (SsR * B2_rr - SsL * B2_ll - (f_rr[7] - f_ll[7])) / Sdiff369B3Star = (SsR * B3_rr - SsL * B3_ll - (f_rr[8] - f_ll[8])) / Sdiff370if SsL <= SStar371SdiffStar = SsL - SStar372373densStar = rho_ll * sMu_L / SdiffStar # (19)374375mom_1_Star = densStar * SStar # (20)376mom_2_Star = densStar * v2_ll -377(B1Star * B2Star - B1_ll * B2_ll) / SdiffStar # (21)378mom_3_Star = densStar * v3_ll -379(B1Star * B3Star - B1_ll * B3_ll) / SdiffStar # (22)380381#p_tot_Star = rho_ll * sMu_L * (SStar - v1_ll) + p_tot_ll - B1_ll^2 + B1Star^2 # (17)382# 1D B1 = constant => B1_ll = B1_rr = B1Star383p_tot_Star = rho_ll * sMu_L * (SStar - v1_ll) + p_tot_ll # (17)384385enerStar = u_ll[5] * sMu_L / SdiffStar +386(p_tot_Star * SStar - p_tot_ll * v1_ll - (B1Star *387(B1Star * v1Star + B2Star * v2Star + B3Star * v3Star) -388B1_ll * (B1_ll * v1_ll + B2_ll * v2_ll + B3_ll * v3_ll))) /389SdiffStar # (23)390391# Classic HLLC update (32)392f1 = f_ll[1] + SsL * (densStar - u_ll[1])393f2 = f_ll[2] + SsL * (mom_1_Star - u_ll[2])394f3 = f_ll[3] + SsL * (mom_2_Star - u_ll[3])395f4 = f_ll[4] + SsL * (mom_3_Star - u_ll[4])396f5 = f_ll[5] + SsL * (enerStar - u_ll[5])397f6 = f_ll[6] + SsL * (B1Star - u_ll[6])398f7 = f_ll[7] + SsL * (B2Star - u_ll[7])399f8 = f_ll[8] + SsL * (B3Star - u_ll[8])400else # SStar <= Ssr401SdiffStar = SsR - SStar402403densStar = rho_rr * sMu_R / SdiffStar # (19)404405mom_1_Star = densStar * SStar # (20)406mom_2_Star = densStar * v2_rr -407(B1Star * B2Star - B1_rr * B2_rr) / SdiffStar # (21)408mom_3_Star = densStar * v3_rr -409(B1Star * B3Star - B1_rr * B3_rr) / SdiffStar # (22)410411#p_tot_Star = rho_rr * sMu_R * (SStar - v1_rr) + p_tot_rr - B1_rr^2 + B1Star^2 # (17)412# 1D B1 = constant => B1_ll = B1_rr = B1Star413p_tot_Star = rho_rr * sMu_R * (SStar - v1_rr) + p_tot_rr # (17)414415enerStar = u_rr[5] * sMu_R / SdiffStar +416(p_tot_Star * SStar - p_tot_rr * v1_rr - (B1Star *417(B1Star * v1Star + B2Star * v2Star + B3Star * v3Star) -418B1_rr * (B1_rr * v1_rr + B2_rr * v2_rr + B3_rr * v3_rr))) /419SdiffStar # (23)420421# Classic HLLC update (32)422f1 = f_rr[1] + SsR * (densStar - u_rr[1])423f2 = f_rr[2] + SsR * (mom_1_Star - u_rr[2])424f3 = f_rr[3] + SsR * (mom_2_Star - u_rr[3])425f4 = f_rr[4] + SsR * (mom_3_Star - u_rr[4])426f5 = f_rr[5] + SsR * (enerStar - u_rr[5])427f6 = f_rr[6] + SsR * (B1Star - u_rr[6])428f7 = f_rr[7] + SsR * (B2Star - u_rr[7])429f8 = f_rr[8] + SsR * (B3Star - u_rr[8])430end431end432return SVector(f1, f2, f3, f4, f5, f6, f7, f8)433end434435# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation436@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,437equations::IdealGlmMhdEquations1D)438rho_ll, rho_v1_ll, _ = u_ll439rho_rr, rho_v1_rr, _ = u_rr440441# Calculate velocities (ignore orientation since it is always "1" in 1D)442# and fast magnetoacoustic wave speeds443# left444v_ll = rho_v1_ll / rho_ll445cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)446# right447v_rr = rho_v1_rr / rho_rr448cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)449450return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)451end452453# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`454@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,455equations::IdealGlmMhdEquations1D)456rho_ll, rho_v1_ll, _ = u_ll457rho_rr, rho_v1_rr, _ = u_rr458459# Calculate velocities (ignore orientation since it is always "1" in 1D)460# and fast magnetoacoustic wave speeds461# left462v_ll = rho_v1_ll / rho_ll463cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)464# right465v_rr = rho_v1_rr / rho_rr466cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)467468return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)469end470471# Calculate estimates for minimum and maximum wave speeds for HLL-type fluxes472@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,473equations::IdealGlmMhdEquations1D)474rho_ll, rho_v1_ll, _ = u_ll475rho_rr, rho_v1_rr, _ = u_rr476477# Calculate primitive variables478v1_ll = rho_v1_ll / rho_ll479v1_rr = rho_v1_rr / rho_rr480481λ_min = v1_ll - calc_fast_wavespeed(u_ll, orientation, equations)482λ_max = v1_rr + calc_fast_wavespeed(u_rr, orientation, equations)483484return λ_min, λ_max485end486487# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes488@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,489equations::IdealGlmMhdEquations1D)490rho_ll, rho_v1_ll, _ = u_ll491rho_rr, rho_v1_rr, _ = u_rr492493# Calculate primitive variables494v1_ll = rho_v1_ll / rho_ll495v1_rr = rho_v1_rr / rho_rr496497# Approximate the left-most and right-most eigenvalues in the Riemann fan498c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)499c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)500501λ_min = min(v1_ll - c_f_ll, v1_rr - c_f_rr)502λ_max = max(v1_ll + c_f_ll, v1_rr + c_f_rr)503504return λ_min, λ_max505end506507"""508min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer, equations::IdealGlmMhdEquations1D)509510Calculate minimum and maximum wave speeds for HLL-type fluxes as in511- Li (2005)512An HLLC Riemann solver for magneto-hydrodynamics513[DOI: 10.1016/j.jcp.2004.08.020](https://doi.org/10.1016/j.jcp.2004.08.020).514"""515@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,516equations::IdealGlmMhdEquations1D)517rho_ll, rho_v1_ll, _ = u_ll518rho_rr, rho_v1_rr, _ = u_rr519520# Calculate primitive variables521v1_ll = rho_v1_ll / rho_ll522v1_rr = rho_v1_rr / rho_rr523524# Approximate the left-most and right-most eigenvalues in the Riemann fan525c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)526c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)527vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)528λ_min = min(v1_ll - c_f_ll, vel_roe - c_f_roe)529λ_max = max(v1_rr + c_f_rr, vel_roe + c_f_roe)530531return λ_min, λ_max532end533534@inline function max_abs_speeds(u, equations::IdealGlmMhdEquations1D)535rho, rho_v1, _ = u536v1 = rho_v1 / rho537cf_x_direction = calc_fast_wavespeed(u, 1, equations)538539return abs(v1) + cf_x_direction540end541542# Convert conservative variables to primitive543@inline function cons2prim(u, equations::IdealGlmMhdEquations1D)544rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u545546v1 = rho_v1 / rho547v2 = rho_v2 / rho548v3 = rho_v3 / rho549p = (equations.gamma - 1) * (rho_e_total -5500.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3551+ B1 * B1 + B2 * B2 + B3 * B3))552553return SVector(rho, v1, v2, v3, p, B1, B2, B3)554end555556# Convert conservative variables to entropy557@inline function cons2entropy(u, equations::IdealGlmMhdEquations1D)558rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u559560v1 = rho_v1 / rho561v2 = rho_v2 / rho562v3 = rho_v3 / rho563v_square = v1^2 + v2^2 + v3^2564p = (equations.gamma - 1) *565(rho_e_total - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2))566s = log(p) - equations.gamma * log(rho)567rho_p = rho / p568569w1 = (equations.gamma - s) / (equations.gamma - 1) - 0.5f0 * rho_p * v_square570w2 = rho_p * v1571w3 = rho_p * v2572w4 = rho_p * v3573w5 = -rho_p574w6 = rho_p * B1575w7 = rho_p * B2576w8 = rho_p * B3577578return SVector(w1, w2, w3, w4, w5, w6, w7, w8)579end580581# Convert primitive to conservative variables582@inline function prim2cons(prim, equations::IdealGlmMhdEquations1D)583rho, v1, v2, v3, p, B1, B2, B3 = prim584585rho_v1 = rho * v1586rho_v2 = rho * v2587rho_v3 = rho * v3588rho_e_total = p / (equations.gamma - 1) +5890.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3) +5900.5f0 * (B1^2 + B2^2 + B3^2)591592return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3)593end594595@inline function density(u, equations::IdealGlmMhdEquations1D)596rho = u[1]597return rho598end599600@inline function velocity(u, equations::IdealGlmMhdEquations1D)601rho = u[1]602v1 = u[2] / rho603v2 = u[3] / rho604v3 = u[4] / rho605return SVector(v1, v2, v3)606end607608@inline function velocity(u, orientation::Int, equations::IdealGlmMhdEquations1D)609rho = u[1]610v = u[orientation + 1] / rho611return v612end613614@doc raw"""615pressure(u, equations::IdealGlmMhdEquations1D)616617Computes the pressure for an ideal equation of state with618isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`.619```math620\begin{aligned}621p &= (\gamma - 1) \left( \rho e_{\text{total}} - \rho e_{\text{kinetic}} - \rho e_{\text{magnetic}} \right) \\622&= (\gamma - 1) \left( \rho e - \frac{1}{2}623\left[\rho \Vert v \Vert_2^2 + \Vert B \Vert_2^2 \right] \right)624\end{aligned}625```626"""627@inline function pressure(u, equations::IdealGlmMhdEquations1D)628rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u629p = (equations.gamma - 1) * (rho_e_total -6300.5f0 *631((rho_v1^2 + rho_v2^2 + rho_v3^2) / rho +632B1^2 + B2^2 + B3^2))633return p634end635636@inline function density_pressure(u, equations::IdealGlmMhdEquations1D)637rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u638rho_times_p = (equations.gamma - 1) * (rho * rho_e_total -6390.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2 +640rho * (B1^2 + B2^2 + B3^2)))641return rho_times_p642end643644# Compute the fastest wave speed for ideal MHD equations: c_f, the fast magnetoacoustic eigenvalue645@inline function calc_fast_wavespeed(cons, direction, equations::IdealGlmMhdEquations1D)646rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = cons647v1 = rho_v1 / rho648v2 = rho_v2 / rho649v3 = rho_v3 / rho650v_mag = sqrt(v1^2 + v2^2 + v3^2)651p = (equations.gamma - 1) *652(rho_e_total - 0.5f0 * rho * v_mag^2 - 0.5f0 * (B1^2 + B2^2 + B3^2))653a_square = equations.gamma * p / rho654sqrt_rho = sqrt(rho)655b1 = B1 / sqrt_rho656b2 = B2 / sqrt_rho657b3 = B3 / sqrt_rho658b_square = b1^2 + b2^2 + b3^2659660c_f = sqrt(0.5f0 * (a_square + b_square) +6610.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b1^2))662return c_f663end664665"""666calc_fast_wavespeed_roe(u_ll, u_rr, direction, equations::IdealGlmMhdEquations1D)667668Compute the fast magnetoacoustic wave speed using Roe averages669as given by670- Cargo and Gallice (1997)671Roe Matrices for Ideal MHD and Systematic Construction672of Roe Matrices for Systems of Conservation Laws673[DOI: 10.1006/jcph.1997.5773](https://doi.org/10.1006/jcph.1997.5773)674"""675@inline function calc_fast_wavespeed_roe(u_ll, u_rr, direction,676equations::IdealGlmMhdEquations1D)677rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll = u_ll678rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr = u_rr679680# Calculate primitive variables681v1_ll = rho_v1_ll / rho_ll682v2_ll = rho_v2_ll / rho_ll683v3_ll = rho_v3_ll / rho_ll684vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2685mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^2686p_ll = (equations.gamma - 1) *687(rho_e_total_ll - 0.5f0 * rho_ll * vel_norm_ll - 0.5f0 * mag_norm_ll)688689v1_rr = rho_v1_rr / rho_rr690v2_rr = rho_v2_rr / rho_rr691v3_rr = rho_v3_rr / rho_rr692vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2693mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^2694p_rr = (equations.gamma - 1) *695(rho_e_total_rr - 0.5f0 * rho_rr * vel_norm_rr - 0.5f0 * mag_norm_rr)696697# compute total pressure which is thermal + magnetic pressures698p_total_ll = p_ll + 0.5f0 * mag_norm_ll699p_total_rr = p_rr + 0.5f0 * mag_norm_rr700701# compute the Roe density averages702sqrt_rho_ll = sqrt(rho_ll)703sqrt_rho_rr = sqrt(rho_rr)704inv_sqrt_rho_add = 1 / (sqrt_rho_ll + sqrt_rho_rr)705inv_sqrt_rho_prod = 1 / (sqrt_rho_ll * sqrt_rho_rr)706rho_ll_roe = sqrt_rho_ll * inv_sqrt_rho_add707rho_rr_roe = sqrt_rho_rr * inv_sqrt_rho_add708# Roe averages709# velocities and magnetic fields710v1_roe = v1_ll * rho_ll_roe + v1_rr * rho_rr_roe711v2_roe = v2_ll * rho_ll_roe + v2_rr * rho_rr_roe712v3_roe = v3_ll * rho_ll_roe + v3_rr * rho_rr_roe713B1_roe = B1_ll * rho_ll_roe + B1_rr * rho_rr_roe714B2_roe = B2_ll * rho_ll_roe + B2_rr * rho_rr_roe715B3_roe = B3_ll * rho_ll_roe + B3_rr * rho_rr_roe716# enthalpy717H_ll = (rho_e_total_ll + p_total_ll) / rho_ll718H_rr = (rho_e_total_rr + p_total_rr) / rho_rr719H_roe = H_ll * rho_ll_roe + H_rr * rho_rr_roe720# temporary variable see equations (4.12) in Cargo and Gallice721X = 0.5f0 * ((B1_ll - B1_rr)^2 + (B2_ll - B2_rr)^2 + (B3_ll - B3_rr)^2) *722inv_sqrt_rho_add^2723# averaged components needed to compute c_f, the fast magnetoacoustic wave speed724b_square_roe = (B1_roe^2 + B2_roe^2 + B3_roe^2) * inv_sqrt_rho_prod # scaled magnectic sum725a_square_roe = ((2 - equations.gamma) * X +726(equations.gamma - 1) *727(H_roe - 0.5f0 * (v1_roe^2 + v2_roe^2 + v3_roe^2) -728b_square_roe)) # acoustic speed729# finally compute the average wave speed and set the output velocity730# Ignore orientation since it is always "1" in 1D731c_a_roe = B1_roe^2 * inv_sqrt_rho_prod # (squared) Alfvén wave speed732a_star_roe = sqrt((a_square_roe + b_square_roe)^2 - 4 * a_square_roe * c_a_roe)733c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe))734735return v1_roe, c_f_roe736end737738@doc raw"""739entropy_thermodynamic(cons, equations::AbstractIdealGlmMhdEquations)740741Calculate thermodynamic entropy for a conservative state `cons` as742743```math744s = \log(p) - \gamma \log(\rho)745```746"""747@inline function entropy_thermodynamic(cons, equations::IdealGlmMhdEquations1D)748# Pressure749p = (equations.gamma - 1) *750(cons[5] - 0.5f0 * (cons[2]^2 + cons[3]^2 + cons[4]^2) / cons[1]751-7520.5f0 * (cons[6]^2 + cons[7]^2 + cons[8]^2))753754# Thermodynamic entropy755s = log(p) - equations.gamma * log(cons[1])756757return s758end759760@doc raw"""761entropy_math(cons, equations::AbstractIdealGlmMhdEquations)762763Calculate mathematical entropy for a conservative state `cons` as764```math765S = -\frac{\rho s}{\gamma - 1}766```767where `s` is the thermodynamic entropy calculated by768[`entropy_thermodynamic(cons, equations::AbstractIdealGlmMhdEquations)`](@ref).769"""770@inline function entropy_math(cons, equations::IdealGlmMhdEquations1D)771S = -entropy_thermodynamic(cons, equations) * cons[1] / (equations.gamma - 1)772773return S774end775776"""777entropy(cons, equations::AbstractIdealGlmMhdEquations)778779Default entropy is the mathematical entropy780[`entropy_math(cons, equations::AbstractIdealGlmMhdEquations)`](@ref).781"""782@inline entropy(cons, equations::IdealGlmMhdEquations1D) = entropy_math(cons, equations)783784# Calculate total energy for a conservative state `cons`785@inline energy_total(cons, ::IdealGlmMhdEquations1D) = cons[5]786787# Calculate kinetic energy for a conservative state `cons`788@inline function energy_kinetic(cons, equations::IdealGlmMhdEquations1D)789return 0.5f0 * (cons[2]^2 + cons[3]^2 + cons[4]^2) / cons[1]790end791792# Calculate the magnetic energy for a conservative state `cons`.793# OBS! For non-dinmensional form of the ideal MHD magnetic pressure ≡ magnetic energy794@inline function energy_magnetic(cons, ::IdealGlmMhdEquations1D)795return 0.5f0 * (cons[6]^2 + cons[7]^2 + cons[8]^2)796end797798"""799energy_internal(cons, equations::AbstractIdealGlmMhdEquations)800801Calculate internal energy for a conservative state `cons` as the difference802between total energy and kinetic + magnetic energies.803"""804@inline function energy_internal(cons, equations::IdealGlmMhdEquations1D)805return (energy_total(cons, equations)806-807energy_kinetic(cons, equations)808-809energy_magnetic(cons, equations))810end811812# Calculate the cross helicity (\vec{v}⋅\vec{B}) for a conservative state `cons`813@inline function cross_helicity(cons, ::IdealGlmMhdEquations1D)814return (cons[2] * cons[6] + cons[3] * cons[7] + cons[4] * cons[8]) / cons[1]815end816end # @muladd817818819