Path: blob/main/src/equations/ideal_glm_mhd_multicomponent_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"""8IdealGlmMhdMulticomponentEquations1D(; gammas, gas_constants)910The ideal compressible multicomponent GLM-MHD equations11```math12\frac{\partial}{\partial t}13\begin{pmatrix}14\rho v_1 \\ \rho v_2 \\ \rho v_3 \\ \rho e_{\text{total}} \\ B_1 \\ B_2 \\ B_3 \\ \rho_1 \\ \vdots \\ \rho_{n}15\end{pmatrix}16+17\frac{\partial}{\partial x}18\begin{pmatrix}19\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_1 \\ \rho_1 v_1 \\ \vdots \\ \rho_{n} v_122\end{pmatrix}23=24\begin{pmatrix}250 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0 \\ \vdots \\ 026\end{pmatrix}27```28for calorically perfect gases in one space dimension.29Here, ``\rho_i`` is the density of component ``i``, ``\rho=\sum_{i=1}^n\rho_i`` the sum of the individual ``\rho_i``,30``\mathbf{v}`` the velocity, ``\mathbf{B}`` the magnetic field, ``e_{\text{total}}`` the specific total energy, and31```math32p = (\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)33```34the pressure,35```math36\gamma=\frac{\sum_{i=1}^n\rho_i C_{v,i}\gamma_i}{\sum_{i=1}^n\rho_i C_{v,i}}37```38total heat capacity ratio, ``\gamma_i`` heat capacity ratio of component ``i``,39```math40C_{v,i}=\frac{R_i}{\gamma_i-1}41```42specific heat capacity at constant volume of component ``i``.4344In case of more than one component, the specific heat ratios `gammas` and the gas constants45`gas_constants` should be passed as tuples, e.g., `gammas = (1.4, 1.667)`.4647The remaining variables like the specific heats at constant volume `cv` or the specific heats at48constant pressure `cp` are then calculated considering a calorically perfect gas.49"""50struct IdealGlmMhdMulticomponentEquations1D{NVARS, NCOMP, RealT <: Real} <:51AbstractIdealGlmMhdMulticomponentEquations{1, NVARS, NCOMP}52gammas::SVector{NCOMP, RealT}53gas_constants::SVector{NCOMP, RealT}54cv::SVector{NCOMP, RealT}55cp::SVector{NCOMP, RealT}5657function IdealGlmMhdMulticomponentEquations1D{NVARS, NCOMP, RealT}(gammas::SVector{NCOMP,58RealT},59gas_constants::SVector{NCOMP,60RealT}) where {61NVARS,62NCOMP,63RealT <:64Real65}66NCOMP >= 1 ||67throw(DimensionMismatch("`gammas` and `gas_constants` have to be filled with at least one value"))6869cv = gas_constants ./ (gammas .- 1)70cp = gas_constants + gas_constants ./ (gammas .- 1)7172return new(gammas, gas_constants, cv, cp)73end74end7576function IdealGlmMhdMulticomponentEquations1D(; gammas, gas_constants)77_gammas = promote(gammas...)78_gas_constants = promote(gas_constants...)79RealT = promote_type(eltype(_gammas), eltype(_gas_constants))8081__gammas = SVector(map(RealT, _gammas))82__gas_constants = SVector(map(RealT, _gas_constants))8384NVARS = length(_gammas) + 785NCOMP = length(_gammas)8687return IdealGlmMhdMulticomponentEquations1D{NVARS, NCOMP, RealT}(__gammas,88__gas_constants)89end9091# Outer constructor for `@reset` works correctly92function IdealGlmMhdMulticomponentEquations1D(gammas, gas_constants, cv, cp, c_h)93return IdealGlmMhdMulticomponentEquations1D(gammas = gammas,94gas_constants = gas_constants)95end9697@inline function Base.real(::IdealGlmMhdMulticomponentEquations1D{NVARS, NCOMP, RealT}) where {98NVARS,99NCOMP,100RealT101}102return RealT103end104105have_nonconservative_terms(::IdealGlmMhdMulticomponentEquations1D) = False()106107function varnames(::typeof(cons2cons), equations::IdealGlmMhdMulticomponentEquations1D)108cons = ("rho_v1", "rho_v2", "rho_v3", "rho_e_total", "B1", "B2", "B3")109rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))110return (cons..., rhos...)111end112113function varnames(::typeof(cons2prim), equations::IdealGlmMhdMulticomponentEquations1D)114prim = ("v1", "v2", "v3", "p", "B1", "B2", "B3")115rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))116return (prim..., rhos...)117end118119"""120initial_condition_convergence_test(x, t, equations::IdealGlmMhdMulticomponentEquations1D)121122An Alfvén wave as smooth initial condition used for convergence tests.123"""124function initial_condition_convergence_test(x, t,125equations::IdealGlmMhdMulticomponentEquations1D)126# smooth Alfvén wave test from Derigs et al. FLASH (2016)127# domain must be set to [0, 1], γ = 5/3128129RealT = eltype(x)130rho = one(RealT)131prim_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) * (1 - 2) *132rho / (1 -1332^ncomponents(equations))134for i in eachcomponent(equations))135v1 = 0136si, co = sincospi(2 * x[1])137v2 = convert(RealT, 0.1) * si138v3 = convert(RealT, 0.1) * co139p = convert(RealT, 0.1)140B1 = 1141B2 = v2142B3 = v3143prim_other = SVector(v1, v2, v3, p, B1, B2, B3)144145return prim2cons(vcat(prim_other, prim_rho), equations)146end147148"""149initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdMulticomponentEquations1D)150151A weak blast wave adapted from152- Sebastian Hennemann, Gregor J. Gassner (2020)153A provably entropy stable subcell shock capturing approach for high order split form DG154[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)155"""156function initial_condition_weak_blast_wave(x, t,157equations::IdealGlmMhdMulticomponentEquations1D)158# Adapted MHD version of the weak blast wave from Hennemann & Gassner JCP paper 2020 (Sec. 6.3)159# Same discontinuity in the velocities but with magnetic fields160# Set up polar coordinates161RealT = eltype(x)162inicenter = (0)163x_norm = x[1] - inicenter[1]164r = sqrt(x_norm^2)165phi = atan(x_norm)166167# Calculate primitive variables168# We initialize each species with a fraction of the total density `rho`, such169# that the sum of the densities is `rho := density(prim, equations)`. The density of170# a species is double the density of the next species.171if r > 0.5f0172rho = one(RealT)173prim_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) *174(1 - 2) * rho /175(1 -1762^ncomponents(equations))177for i in eachcomponent(equations))178else179rho = convert(RealT, 1.1691)180prim_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) *181(1 - 2) * rho /182(1 -1832^ncomponents(equations))184for i in eachcomponent(equations))185end186v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos(phi)187p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)188189prim_other = SVector(v1, 0, 0, p, 1, 1, 1)190191return prim2cons(vcat(prim_other, prim_rho), equations)192end193194# Calculate 1D flux for a single point195@inline function flux(u, orientation::Integer,196equations::IdealGlmMhdMulticomponentEquations1D)197rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u198199rho = density(u, equations)200201v1 = rho_v1 / rho202v2 = rho_v2 / rho203v3 = rho_v3 / rho204kin_en = 0.5f0 * rho * (v1^2 + v2^2 + v3^2)205mag_en = 0.5f0 * (B1^2 + B2^2 + B3^2)206gamma = totalgamma(u, equations)207p = (gamma - 1) * (rho_e_total - kin_en - mag_en)208209f_rho = densities(u, v1, equations)210f1 = rho_v1 * v1 + p + mag_en - B1^2211f2 = rho_v1 * v2 - B1 * B2212f3 = rho_v1 * v3 - B1 * B3213f4 = (kin_en + gamma * p / (gamma - 1) + 2 * mag_en) * v1 -214B1 * (v1 * B1 + v2 * B2 + v3 * B3)215f5 = 0216f6 = v1 * B2 - v2 * B1217f7 = v1 * B3 - v3 * B1218219f_other = SVector(f1, f2, f3, f4, f5, f6, f7)220221return vcat(f_other, f_rho)222end223224"""225flux_derigs_etal(u_ll, u_rr, orientation, equations::IdealGlmMhdEquations1D)226227Entropy conserving two-point flux adapted by228- Derigs et al. (2018)229Ideal GLM-MHD: About the entropy consistent nine-wave magnetic field230divergence diminishing ideal magnetohydrodynamics equations for multicomponent231[DOI: 10.1016/j.jcp.2018.03.002](https://doi.org/10.1016/j.jcp.2018.03.002)232"""233function flux_derigs_etal(u_ll, u_rr, orientation::Integer,234equations::IdealGlmMhdMulticomponentEquations1D)235# Unpack left and right states to get velocities, pressure, and inverse temperature (called beta)236rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll = u_ll237rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr = u_rr238@unpack gammas, gas_constants, cv = equations239240rho_ll = density(u_ll, equations)241rho_rr = density(u_rr, equations)242243gamma_ll = totalgamma(u_ll, equations)244gamma_rr = totalgamma(u_rr, equations)245246rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 7],247u_rr[i + 7])248for i in eachcomponent(equations))249rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 7] +250u_rr[i + 7])251for i in eachcomponent(equations))252253v1_ll = rho_v1_ll / rho_ll254v2_ll = rho_v2_ll / rho_ll255v3_ll = rho_v3_ll / rho_ll256v1_rr = rho_v1_rr / rho_rr257v2_rr = rho_v2_rr / rho_rr258v3_rr = rho_v3_rr / rho_rr259vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2260vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2261mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^2262mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^2263# for convenience store v⋅B264vel_dot_mag_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll265vel_dot_mag_rr = v1_rr * B1_rr + v2_rr * B2_rr + v3_rr * B3_rr266267# Compute the necessary mean values needed for either direction268v1_avg = 0.5f0 * (v1_ll + v1_rr)269v2_avg = 0.5f0 * (v2_ll + v2_rr)270v3_avg = 0.5f0 * (v3_ll + v3_rr)271v_sum = v1_avg + v2_avg + v3_avg272B1_avg = 0.5f0 * (B1_ll + B1_rr)273B2_avg = 0.5f0 * (B2_ll + B2_rr)274B3_avg = 0.5f0 * (B3_ll + B3_rr)275vel_norm_avg = 0.5f0 * (vel_norm_ll + vel_norm_rr)276mag_norm_avg = 0.5f0 * (mag_norm_ll + mag_norm_rr)277vel_dot_mag_avg = 0.5f0 * (vel_dot_mag_ll + vel_dot_mag_rr)278279RealT = eltype(u_ll)280enth = zero(RealT)281help1_ll = zero(RealT)282help1_rr = zero(RealT)283284for i in eachcomponent(equations)285enth += rhok_avg[i] * gas_constants[i]286help1_ll += u_ll[i + 7] * cv[i]287help1_rr += u_rr[i + 7] * cv[i]288end289290T_ll = (rho_e_total_ll - 0.5f0 * rho_ll * (vel_norm_ll) - 0.5f0 * mag_norm_ll) /291help1_ll292T_rr = (rho_e_total_rr - 0.5f0 * rho_rr * (vel_norm_rr) - 0.5f0 * mag_norm_rr) /293help1_rr294T = 0.5f0 * (1 / T_ll + 1 / T_rr)295T_log = ln_mean(1 / T_ll, 1 / T_rr)296297# Calculate fluxes depending on orientation with specific direction averages298help1 = zero(RealT)299help2 = zero(RealT)300301f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg302for i in eachcomponent(equations))303for i in eachcomponent(equations)304help1 += f_rho[i] * cv[i]305help2 += f_rho[i]306end307f1 = help2 * v1_avg + enth / T + 0.5f0 * mag_norm_avg - B1_avg * B1_avg308f2 = help2 * v2_avg - B1_avg * B2_avg309f3 = help2 * v3_avg - B1_avg * B3_avg310f5 = 0311f6 = v1_avg * B2_avg - v2_avg * B1_avg312f7 = v1_avg * B3_avg - v3_avg * B1_avg313314# total energy flux is complicated and involves the previous eight components315v1_mag_avg = 0.5f0 * (v1_ll * mag_norm_ll + v1_rr * mag_norm_rr)316317f4 = (help1 / T_log) - 0.5f0 * (vel_norm_avg) * (help2) + f1 * v1_avg +318f2 * v2_avg +319f3 * v3_avg +320f5 * B1_avg + f6 * B2_avg + f7 * B3_avg - 0.5f0 * v1_mag_avg +321B1_avg * vel_dot_mag_avg322323f_other = SVector(f1, f2, f3, f4, f5, f6, f7)324325return vcat(f_other, f_rho)326end327328"""329flux_hindenlang_gassner(u_ll, u_rr, orientation_or_normal_direction,330equations::IdealGlmMhdMulticomponentEquations1D)331332Adaption of the entropy conserving and kinetic energy preserving two-point flux of333Hindenlang (2019), extending [`flux_ranocha`](@ref) to the MHD equations.334## References335- Florian Hindenlang, Gregor Gassner (2019)336A new entropy conservative two-point flux for ideal MHD equations derived from337first principles.338Presented at HONOM 2019: European workshop on high order numerical methods339for evolutionary PDEs, theory and applications340- Hendrik Ranocha (2018)341Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods342for Hyperbolic Balance Laws343[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)344- Hendrik Ranocha (2020)345Entropy Conserving and Kinetic Energy Preserving Numerical Methods for346the Euler Equations Using Summation-by-Parts Operators347[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)348"""349@inline function flux_hindenlang_gassner(u_ll, u_rr, orientation::Integer,350equations::IdealGlmMhdMulticomponentEquations1D)351# Unpack left and right states352v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll = cons2prim(u_ll, equations)353v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr = cons2prim(u_rr, equations)354355rho_ll = density(u_ll, equations)356rho_rr = density(u_rr, equations)357358# Compute the necessary mean values needed for either direction359# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`360# in exact arithmetic since361# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)362# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)363inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)364v1_avg = 0.5f0 * (v1_ll + v1_rr)365v2_avg = 0.5f0 * (v2_ll + v2_rr)366v3_avg = 0.5f0 * (v3_ll + v3_rr)367p_avg = 0.5f0 * (p_ll + p_rr)368velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)369magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)370371inv_gamma_minus_one = 1 / (totalgamma(0.5f0 * (u_ll + u_rr), equations) - 1)372373rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 7],374u_rr[i + 7])375for i in eachcomponent(equations))376rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 7] +377u_rr[i + 7])378for i in eachcomponent(equations))379380RealT = eltype(u_ll)381f1 = zero(RealT)382f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg383for i in eachcomponent(equations))384for i in eachcomponent(equations)385f1 += f_rho[i]386end387388# Calculate fluxes depending on orientation with specific direction averages389f2 = f1 * v1_avg + p_avg + magnetic_square_avg -3900.5f0 * (B1_ll * B1_rr + B1_rr * B1_ll)391f3 = f1 * v2_avg - 0.5f0 * (B1_ll * B2_rr + B1_rr * B2_ll)392f4 = f1 * v3_avg - 0.5f0 * (B1_ll * B3_rr + B1_rr * B3_ll)393# f5 below394f6 = 0395f7 = 0.5f0 * (v1_ll * B2_ll - v2_ll * B1_ll + v1_rr * B2_rr - v2_rr * B1_rr)396f8 = 0.5f0 * (v1_ll * B3_ll - v3_ll * B1_ll + v1_rr * B3_rr - v3_rr * B1_rr)397# total energy flux is complicated and involves the previous components398f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * inv_gamma_minus_one)399+4000.5f0 * (+p_ll * v1_rr + p_rr * v1_ll401+ (v1_ll * B2_ll * B2_rr + v1_rr * B2_rr * B2_ll)402+ (v1_ll * B3_ll * B3_rr + v1_rr * B3_rr * B3_ll)403-404(v2_ll * B1_ll * B2_rr + v2_rr * B1_rr * B2_ll)405-406(v3_ll * B1_ll * B3_rr + v3_rr * B1_rr * B3_ll)))407408f_other = SVector(f2, f3, f4, f5, f6, f7, f8)409410return vcat(f_other, f_rho)411end412413# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation414@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,415equations::IdealGlmMhdMulticomponentEquations1D)416rho_v1_ll, _ = u_ll417rho_v1_rr, _ = u_rr418419rho_ll = density(u_ll, equations)420rho_rr = density(u_rr, equations)421422# Calculate velocities (ignore orientation since it is always "1" in 1D)423# and fast magnetoacoustic wave speeds424# left425v_ll = rho_v1_ll / rho_ll426cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)427# right428v_rr = rho_v1_rr / rho_rr429cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)430431return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)432end433434# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`435@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,436equations::IdealGlmMhdMulticomponentEquations1D)437rho_v1_ll, _ = u_ll438rho_v1_rr, _ = u_rr439440rho_ll = density(u_ll, equations)441rho_rr = density(u_rr, equations)442443# Calculate velocities (ignore orientation since it is always "1" in 1D)444# and fast magnetoacoustic wave speeds445# left446v_ll = rho_v1_ll / rho_ll447cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)448# right449v_rr = rho_v1_rr / rho_rr450cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)451452return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)453end454455@inline function max_abs_speeds(u, equations::IdealGlmMhdMulticomponentEquations1D)456rho_v1, _ = u457458rho = density(u, equations)459460v1 = rho_v1 / rho461462cf_x_direction = calc_fast_wavespeed(u, 1, equations)463464return (abs(v1) + cf_x_direction,)465end466467# Convert conservative variables to primitive468function cons2prim(u, equations::IdealGlmMhdMulticomponentEquations1D)469rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u470471prim_rho = SVector{ncomponents(equations), real(equations)}(u[i + 7]472for i in eachcomponent(equations))473rho = density(u, equations)474475v1 = rho_v1 / rho476v2 = rho_v2 / rho477v3 = rho_v3 / rho478479gamma = totalgamma(u, equations)480481p = (gamma - 1) *482(rho_e_total - 0.5f0 * rho * (v1^2 + v2^2 + v3^2) -4830.5f0 * (B1^2 + B2^2 + B3^2))484prim_other = SVector(v1, v2, v3, p, B1, B2, B3)485486return vcat(prim_other, prim_rho)487end488489# Convert conservative variables to entropy490@inline function cons2entropy(u, equations::IdealGlmMhdMulticomponentEquations1D)491rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u492@unpack cv, gammas, gas_constants = equations493494rho = density(u, equations)495496v1 = rho_v1 / rho497v2 = rho_v2 / rho498v3 = rho_v3 / rho499v_square = v1^2 + v2^2 + v3^2500gamma = totalgamma(u, equations)501p = (gamma - 1) *502(rho_e_total - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2))503s = log(p) - gamma * log(rho)504rho_p = rho / p505506# Multicomponent stuff507RealT = eltype(u)508help1 = zero(RealT)509510for i in eachcomponent(equations)511help1 += u[i + 7] * cv[i]512end513514T = (rho_e_total - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2)) / (help1)515516entrop_rho = SVector{ncomponents(equations), real(equations)}(-1 *517(cv[i] * log(T) -518gas_constants[i] *519log(u[i + 7])) +520gas_constants[i] +521cv[i] -522(v_square / (2 * T))523for i in eachcomponent(equations))524525w1 = v1 / T526w2 = v2 / T527w3 = v3 / T528w4 = -1 / T529w5 = B1 / T530w6 = B2 / T531w7 = B3 / T532533entrop_other = SVector(w1, w2, w3, w4, w5, w6, w7)534535return vcat(entrop_other, entrop_rho)536end537538# Convert primitive to conservative variables539@inline function prim2cons(prim, equations::IdealGlmMhdMulticomponentEquations1D)540v1, v2, v3, p, B1, B2, B3 = prim541542cons_rho = SVector{ncomponents(equations), real(equations)}(prim[i + 7]543for i in eachcomponent(equations))544rho = density(prim, equations)545546rho_v1 = rho * v1547rho_v2 = rho * v2548rho_v3 = rho * v3549550gamma = totalgamma(prim, equations)551rho_e_total = p / (gamma - 1) + 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3) +5520.5f0 * (B1^2 + B2^2 + B3^2)553554cons_other = SVector(rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3)555556return vcat(cons_other, cons_rho)557end558559@doc raw"""560pressure(u, equations::IdealGlmMhdMulticomponentEquations1D)561562Computes the pressure for an ideal equation of state with563isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`.564```math565\begin{aligned}566p &= (\gamma - 1) \left( \rho e_{\text{total}} - \rho e_{\text{kinetic}} - \rho e_{\text{magnetic}} \right) \\567&= (\gamma - 1) \left( \rho e_{\text{total}} - \frac{1}{2}568\left[\rho \Vert v \Vert_2^2 + \Vert B \Vert_2^2 \right] \right)569\end{aligned}570```571"""572@inline function pressure(u, equations::IdealGlmMhdMulticomponentEquations1D)573rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u574rho = density(u, equations)575gamma = totalgamma(u, equations)576p = (gamma - 1) * (rho_e_total -5770.5f0 *578((rho_v1^2 + rho_v2^2 + rho_v3^2) / rho +579B1^2 + B2^2 + B3^2))580return p581end582583@inline function density_pressure(u, equations::IdealGlmMhdMulticomponentEquations1D)584rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u585rho = density(u, equations)586gamma = totalgamma(u, equations)587rho_times_p = (gamma - 1) * (rho * rho_e_total -5880.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2 +589rho * (B1^2 + B2^2 + B3^2)))590return rho_times_p591end592593# Compute the fastest wave speed for ideal MHD equations: c_f, the fast magnetoacoustic eigenvalue594@inline function calc_fast_wavespeed(cons, direction,595equations::IdealGlmMhdMulticomponentEquations1D)596rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = cons597rho = density(cons, equations)598v1 = rho_v1 / rho599v2 = rho_v2 / rho600v3 = rho_v3 / rho601v_mag = sqrt(v1^2 + v2^2 + v3^2)602gamma = totalgamma(cons, equations)603p = (gamma - 1) *604(rho_e_total - 0.5f0 * rho * v_mag^2 - 0.5f0 * (B1^2 + B2^2 + B3^2))605a_square = gamma * p / rho606sqrt_rho = sqrt(rho)607b1 = B1 / sqrt_rho608b2 = B2 / sqrt_rho609b3 = B3 / sqrt_rho610b_square = b1^2 + b2^2 + b3^2611612c_f = sqrt(0.5f0 * (a_square + b_square) +6130.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b1^2))614615return c_f616end617618@doc raw"""619density(u, equations::AbstractIdealGlmMhdMulticomponentEquations)620621Computes the total density ``\rho = \sum_{i=1}^n \rho_i`` from the conserved variables `u`.622"""623@inline function density(u, equations::IdealGlmMhdMulticomponentEquations1D)624RealT = eltype(u)625rho = zero(RealT)626627for i in eachcomponent(equations)628rho += u[i + 7]629end630631return rho632end633634@inline function totalgamma(u, equations::IdealGlmMhdMulticomponentEquations1D)635@unpack cv, gammas = equations636637RealT = eltype(u)638help1 = zero(RealT)639help2 = zero(RealT)640641for i in eachcomponent(equations)642help1 += u[i + 7] * cv[i] * gammas[i]643help2 += u[i + 7] * cv[i]644end645646return help1 / help2647end648649@inline function densities(u, v, equations::IdealGlmMhdMulticomponentEquations1D)650return SVector{ncomponents(equations), real(equations)}(u[i + 7] * v651for i in eachcomponent(equations))652end653end # @muladd654655656