Path: blob/main/src/equations/compressible_euler_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"""8CompressibleEulerMulticomponentEquations1D(; gammas, gas_constants)910Multicomponent version of the compressible Euler equations11```math12\frac{\partial}{\partial t}13\begin{pmatrix}14\rho v_1 \\ \rho e_{\text{total}} \\ \rho_1 \\ \rho_2 \\ \vdots \\ \rho_{n}15\end{pmatrix}16+17\frac{\partial}{\partial x}18\begin{pmatrix}19\rho v_1^2 + p \\ (\rho e_{\text{total}} +p) v_1 \\ \rho_1 v_1 \\ \rho_2 v_1 \\ \vdots \\ \rho_{n} v_120\end{pmatrix}2122=23\begin{pmatrix}240 \\ 0 \\ 0 \\ 0 \\ \vdots \\ 025\end{pmatrix}26```27for calorically perfect gas in one space dimension.28Here, ``\rho_i`` is the density of component ``i``, ``\rho=\sum_{i=1}^n\rho_i`` the sum of the individual ``\rho_i``,29``v_1`` the velocity, ``e_{\text{total}}`` the specific total energy, and30```math31p = (\gamma - 1) \left( \rho e - \frac{1}{2} \rho v_1^2 \right)32```33the pressure,34```math35\gamma=\frac{\sum_{i=1}^n\rho_i C_{v,i}\gamma_i}{\sum_{i=1}^n\rho_i C_{v,i}}36```37total heat capacity ratio, ``\gamma_i`` heat capacity ratio of component ``i``,38```math39C_{v,i}=\frac{R}{\gamma_i-1}40```41specific heat capacity at constant volume of component ``i``.4243In case of more than one component, the specific heat ratios `gammas` and the gas constants44`gas_constants` should be passed as tuples, e.g., `gammas=(1.4, 1.667)`.4546The remaining variables like the specific heats at constant volume `cv` or the specific heats at47constant pressure `cp` are then calculated considering a calorically perfect gas.48"""49struct CompressibleEulerMulticomponentEquations1D{NVARS, NCOMP, RealT <: Real} <:50AbstractCompressibleEulerMulticomponentEquations{1, NVARS, NCOMP}51gammas::SVector{NCOMP, RealT}52gas_constants::SVector{NCOMP, RealT}53cv::SVector{NCOMP, RealT}54cp::SVector{NCOMP, RealT}5556function CompressibleEulerMulticomponentEquations1D{NVARS, NCOMP, RealT}(gammas::SVector{NCOMP,57RealT},58gas_constants::SVector{NCOMP,59RealT}) where {60NVARS,61NCOMP,62RealT <:63Real64}65NCOMP >= 1 ||66throw(DimensionMismatch("`gammas` and `gas_constants` have to be filled with at least one value"))6768cv = gas_constants ./ (gammas .- 1)69cp = gas_constants + gas_constants ./ (gammas .- 1)7071return new(gammas, gas_constants, cv, cp)72end73end7475function CompressibleEulerMulticomponentEquations1D(; gammas, gas_constants)76_gammas = promote(gammas...)77_gas_constants = promote(gas_constants...)78RealT = promote_type(eltype(_gammas), eltype(_gas_constants),79typeof(gas_constants[1] / (gammas[1] - 1)))80__gammas = SVector(map(RealT, _gammas))81__gas_constants = SVector(map(RealT, _gas_constants))8283NVARS = length(_gammas) + 284NCOMP = length(_gammas)8586return CompressibleEulerMulticomponentEquations1D{NVARS, NCOMP, RealT}(__gammas,87__gas_constants)88end8990@inline function Base.real(::CompressibleEulerMulticomponentEquations1D{NVARS, NCOMP,91RealT}) where {92NVARS,93NCOMP,94RealT95}96return RealT97end9899function varnames(::typeof(cons2cons),100equations::CompressibleEulerMulticomponentEquations1D)101cons = ("rho_v1", "rho_e_total")102rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))103return (cons..., rhos...)104end105106function varnames(::typeof(cons2prim),107equations::CompressibleEulerMulticomponentEquations1D)108prim = ("v1", "p")109rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))110return (prim..., rhos...)111end112113# Set initial conditions at physical location `x` for time `t`114115"""116initial_condition_convergence_test(x, t, equations::CompressibleEulerMulticomponentEquations1D)117118A smooth initial condition used for convergence tests in combination with119[`source_terms_convergence_test`](@ref)120(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).121"""122function initial_condition_convergence_test(x, t,123equations::CompressibleEulerMulticomponentEquations1D)124RealT = eltype(x)125c = 2126A = convert(RealT, 0.1)127L = 2128f = 1.0f0 / L129omega = 2 * convert(RealT, pi) * f130ini = c + A * sin(omega * (x[1] - t))131132v1 = 1133134rho = ini135136# Here we compute an arbitrary number of different rhos. (one rho is double the next rho while the sum of all rhos is 1)137prim_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) * (1 - 2) *138rho / (1 -1392^ncomponents(equations))140for i in eachcomponent(equations))141142prim1 = rho * v1143prim2 = rho^2144145prim_other = SVector(prim1, prim2)146147return vcat(prim_other, prim_rho)148end149150"""151source_terms_convergence_test(u, x, t, equations::CompressibleEulerMulticomponentEquations1D)152153Source terms used for convergence tests in combination with154[`initial_condition_convergence_test`](@ref)155(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).156157References for the method of manufactured solutions (MMS):158- Kambiz Salari and Patrick Knupp (2000)159Code Verification by the Method of Manufactured Solutions160[DOI: 10.2172/759450](https://doi.org/10.2172/759450)161- Patrick J. Roache (2002)162Code Verification by the Method of Manufactured Solutions163[DOI: 10.1115/1.1436090](https://doi.org/10.1115/1.1436090)164"""165@inline function source_terms_convergence_test(u, x, t,166equations::CompressibleEulerMulticomponentEquations1D)167# Same settings as in `initial_condition`168RealT = eltype(u)169c = 2170A = convert(RealT, 0.1)171L = 2172f = 1.0f0 / L173omega = 2 * convert(RealT, pi) * f174175gamma = totalgamma(u, equations)176177x1, = x178si, co = sincos((t - x1) * omega)179tmp = (-((4 * si * A - 4 * c) + 1) * (gamma - 1) * co * A * omega) / 2180181# Here we compute an arbitrary number of different rhos. (one rho is double the next rho while the sum of all rhos is 1182du_rho = SVector{ncomponents(equations), real(equations)}(0183for i in eachcomponent(equations))184185du1 = tmp186du2 = tmp187188du_other = SVector(du1, du2)189190return vcat(du_other, du_rho)191end192193"""194initial_condition_weak_blast_wave(x, t, equations::CompressibleEulerMulticomponentEquations1D)195196A for multicomponent adapted weak blast wave adapted to multicomponent and taken from197- Sebastian Hennemann, Gregor J. Gassner (2020)198A provably entropy stable subcell shock capturing approach for high order split form DG199[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)200"""201function initial_condition_weak_blast_wave(x, t,202equations::CompressibleEulerMulticomponentEquations1D)203# From Hennemann & Gassner JCP paper 2020 (Sec. 6.3)204RealT = eltype(x)205inicenter = SVector(0)206x_norm = x[1] - inicenter[1]207r = abs(x_norm)208cos_phi = x_norm > 0 ? 1 : -1209210prim_rho = SVector{ncomponents(equations), real(equations)}(r > 0.5f0 ?2112^(i - 1) * (1 - 2) /212(RealT(1) -2132^ncomponents(equations)) :2142^(i - 1) * (1 - 2) *215RealT(1.1691) /216(1 -2172^ncomponents(equations))218for i in eachcomponent(equations))219220v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos_phi221p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)222223prim_other = SVector(v1, p)224225return prim2cons(vcat(prim_other, prim_rho), equations)226end227228# Calculate 1D flux for a single point229@inline function flux(u, orientation::Integer,230equations::CompressibleEulerMulticomponentEquations1D)231rho_v1, rho_e_total = u232233rho = density(u, equations)234235v1 = rho_v1 / rho236gamma = totalgamma(u, equations)237p = (gamma - 1) * (rho_e_total - 0.5f0 * rho * v1^2)238239f_rho = densities(u, v1, equations)240f1 = rho_v1 * v1 + p241f2 = (rho_e_total + p) * v1242243f_other = SVector(f1, f2)244245return vcat(f_other, f_rho)246end247248"""249flux_chandrashekar(u_ll, u_rr, orientation, equations::CompressibleEulerMulticomponentEquations1D)250251Entropy conserving two-point flux by252- Ayoub Gouasmi, Karthik Duraisamy (2020)253"Formulation of Entropy-Stable schemes for the multicomponent compressible Euler equations"254[arXiv:1904.00972v3](https://arxiv.org/abs/1904.00972) [math.NA] 4 Feb 2020255"""256@inline function flux_chandrashekar(u_ll, u_rr, orientation::Integer,257equations::CompressibleEulerMulticomponentEquations1D)258# Unpack left and right state259@unpack gammas, gas_constants, cv = equations260rho_v1_ll, rho_e_total_ll = u_ll261rho_v1_rr, rho_e_total_rr = u_rr262rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 2],263u_rr[i + 2])264for i in eachcomponent(equations))265rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 2] +266u_rr[i + 2])267for i in eachcomponent(equations))268269# Iterating over all partial densities270rho_ll = density(u_ll, equations)271rho_rr = density(u_rr, equations)272273gamma_ll = totalgamma(u_ll, equations)274gamma_rr = totalgamma(u_rr, equations)275276# extract velocities277v1_ll = rho_v1_ll / rho_ll278v1_rr = rho_v1_rr / rho_rr279v1_avg = 0.5f0 * (v1_ll + v1_rr)280v1_square = 0.5f0 * (v1_ll^2 + v1_rr^2)281v_sum = v1_avg282283RealT = eltype(u_ll)284enth = zero(RealT)285help1_ll = zero(RealT)286help1_rr = zero(RealT)287288for i in eachcomponent(equations)289enth += rhok_avg[i] * gas_constants[i]290help1_ll += u_ll[i + 2] * cv[i]291help1_rr += u_rr[i + 2] * cv[i]292end293294T_ll = (rho_e_total_ll - 0.5f0 * rho_ll * (v1_ll^2)) / help1_ll295T_rr = (rho_e_total_rr - 0.5f0 * rho_rr * (v1_rr^2)) / help1_rr296T = 0.5f0 * (1 / T_ll + 1 / T_rr)297T_log = ln_mean(1 / T_ll, 1 / T_rr)298299# Calculate fluxes depending on orientation300help1 = zero(RealT)301help2 = zero(RealT)302303f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg304for i in eachcomponent(equations))305for i in eachcomponent(equations)306help1 += f_rho[i] * cv[i]307help2 += f_rho[i]308end309f1 = (help2) * v1_avg + enth / T310f2 = (help1) / T_log - 0.5f0 * (v1_square) * (help2) + v1_avg * f1311312f_other = SVector(f1, f2)313314return vcat(f_other, f_rho)315end316317"""318flux_ranocha(u_ll, u_rr, orientation_or_normal_direction,319equations::CompressibleEulerMulticomponentEquations1D)320321Adaption of the entropy conserving and kinetic energy preserving two-point flux by322- Hendrik Ranocha (2018)323Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods324for Hyperbolic Balance Laws325[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)326See also327- Hendrik Ranocha (2020)328Entropy Conserving and Kinetic Energy Preserving Numerical Methods for329the Euler Equations Using Summation-by-Parts Operators330[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)331"""332@inline function flux_ranocha(u_ll, u_rr, orientation::Integer,333equations::CompressibleEulerMulticomponentEquations1D)334# Unpack left and right state335@unpack gammas, gas_constants, cv = equations336rho_v1_ll, rho_e_total_ll = u_ll337rho_v1_rr, rho_e_total_rr = u_rr338rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 2],339u_rr[i + 2])340for i in eachcomponent(equations))341rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 2] +342u_rr[i + 2])343for i in eachcomponent(equations))344345# Iterating over all partial densities346rho_ll = density(u_ll, equations)347rho_rr = density(u_rr, equations)348349# Calculating gamma350gamma = totalgamma(0.5f0 * (u_ll + u_rr), equations)351inv_gamma_minus_one = 1 / (gamma - 1)352353# extract velocities354v1_ll = rho_v1_ll / rho_ll355v1_rr = rho_v1_rr / rho_rr356v1_avg = 0.5f0 * (v1_ll + v1_rr)357velocity_square_avg = 0.5f0 * (v1_ll * v1_rr)358359# density flux360f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg361for i in eachcomponent(equations))362363# helpful variables364RealT = eltype(u_ll)365f_rho_sum = zero(RealT)366help1_ll = zero(RealT)367help1_rr = zero(RealT)368enth_ll = zero(RealT)369enth_rr = zero(RealT)370for i in eachcomponent(equations)371enth_ll += u_ll[i + 2] * gas_constants[i]372enth_rr += u_rr[i + 2] * gas_constants[i]373f_rho_sum += f_rho[i]374help1_ll += u_ll[i + 2] * cv[i]375help1_rr += u_rr[i + 2] * cv[i]376end377378# temperature and pressure379T_ll = (rho_e_total_ll - 0.5f0 * rho_ll * (v1_ll^2)) / help1_ll380T_rr = (rho_e_total_rr - 0.5f0 * rho_rr * (v1_rr^2)) / help1_rr381p_ll = T_ll * enth_ll382p_rr = T_rr * enth_rr383p_avg = 0.5f0 * (p_ll + p_rr)384inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)385386# momentum and energy flux387f1 = f_rho_sum * v1_avg + p_avg388f2 = f_rho_sum * (velocity_square_avg + inv_rho_p_mean * inv_gamma_minus_one) +3890.5f0 * (p_ll * v1_rr + p_rr * v1_ll)390f_other = SVector(f1, f2)391392return vcat(f_other, f_rho)393end394395# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation396@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,397equations::CompressibleEulerMulticomponentEquations1D)398rho_v1_ll, rho_e_total_ll = u_ll399rho_v1_rr, rho_e_total_rr = u_rr400401# Calculate primitive variables and speed of sound402rho_ll = density(u_ll, equations)403rho_rr = density(u_rr, equations)404gamma_ll = totalgamma(u_ll, equations)405gamma_rr = totalgamma(u_rr, equations)406407v_ll = rho_v1_ll / rho_ll408v_rr = rho_v1_rr / rho_rr409410p_ll = (gamma_ll - 1) * (rho_e_total_ll - 0.5f0 * rho_ll * v_ll^2)411p_rr = (gamma_rr - 1) * (rho_e_total_rr - 0.5f0 * rho_rr * v_rr^2)412c_ll = sqrt(gamma_ll * p_ll / rho_ll)413c_rr = sqrt(gamma_rr * p_rr / rho_rr)414415return max(abs(v_ll), abs(v_rr)) + max(c_ll, c_rr)416end417418# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`419@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,420equations::CompressibleEulerMulticomponentEquations1D)421rho_v1_ll, rho_e_total_ll = u_ll422rho_v1_rr, rho_e_total_rr = u_rr423424# Calculate primitive variables and speed of sound425rho_ll = density(u_ll, equations)426rho_rr = density(u_rr, equations)427gamma_ll = totalgamma(u_ll, equations)428gamma_rr = totalgamma(u_rr, equations)429430v_ll = rho_v1_ll / rho_ll431v_rr = rho_v1_rr / rho_rr432433p_ll = (gamma_ll - 1) * (rho_e_total_ll - 0.5f0 * rho_ll * v_ll^2)434p_rr = (gamma_rr - 1) * (rho_e_total_rr - 0.5f0 * rho_rr * v_rr^2)435c_ll = sqrt(gamma_ll * p_ll / rho_ll)436c_rr = sqrt(gamma_rr * p_rr / rho_rr)437438return max(abs(v_ll) + c_ll, abs(v_rr) + c_rr)439end440441@inline function max_abs_speeds(u,442equations::CompressibleEulerMulticomponentEquations1D)443rho_v1, rho_e_total = u444445rho = density(u, equations)446v1 = rho_v1 / rho447448gamma = totalgamma(u, equations)449p = (gamma - 1) * (rho_e_total - 0.5f0 * rho * (v1^2))450c = sqrt(gamma * p / rho)451452return (abs(v1) + c,)453end454455# Convert conservative variables to primitive456@inline function cons2prim(u, equations::CompressibleEulerMulticomponentEquations1D)457rho_v1, rho_e_total = u458459prim_rho = SVector{ncomponents(equations), real(equations)}(u[i + 2]460for i in eachcomponent(equations))461462rho = density(u, equations)463v1 = rho_v1 / rho464gamma = totalgamma(u, equations)465466p = (gamma - 1) * (rho_e_total - 0.5f0 * rho * (v1^2))467prim_other = SVector(v1, p)468469return vcat(prim_other, prim_rho)470end471472# Convert primitive to conservative variables473@inline function prim2cons(prim, equations::CompressibleEulerMulticomponentEquations1D)474@unpack cv, gammas = equations475v1, p = prim476477RealT = eltype(prim)478479cons_rho = SVector{ncomponents(equations), RealT}(prim[i + 2]480for i in eachcomponent(equations))481rho = density(prim, equations)482gamma = totalgamma(prim, equations)483484rho_v1 = rho * v1485486rho_e_total = p / (gamma - 1) + 0.5f0 * (rho_v1 * v1)487488cons_other = SVector{2, RealT}(rho_v1, rho_e_total)489490return vcat(cons_other, cons_rho)491end492493# Convert conservative variables to entropy494@inline function cons2entropy(u, equations::CompressibleEulerMulticomponentEquations1D)495@unpack cv, gammas, gas_constants = equations496497rho_v1, rho_e_total = u498499rho = density(u, equations)500501RealT = eltype(u)502help1 = zero(RealT)503gas_constant = zero(RealT)504for i in eachcomponent(equations)505help1 += u[i + 2] * cv[i]506gas_constant += gas_constants[i] * (u[i + 2] / rho)507end508509v1 = rho_v1 / rho510v_square = v1^2511gamma = totalgamma(u, equations)512513p = (gamma - 1) * (rho_e_total - 0.5f0 * rho * v_square)514s = log(p) - gamma * log(rho) - log(gas_constant)515rho_p = rho / p516T = (rho_e_total - 0.5f0 * rho * v_square) / (help1)517518entrop_rho = SVector{ncomponents(equations), real(equations)}((cv[i] *519(1 - log(T)) +520gas_constants[i] *521(1 + log(u[i + 2])) -522v1^2 / (2 * T))523for i in eachcomponent(equations))524525w1 = gas_constant * v1 * rho_p526w2 = gas_constant * (-rho_p)527528entrop_other = SVector(w1, w2)529530return vcat(entrop_other, entrop_rho)531end532533# Convert entropy variables to conservative variables534@inline function entropy2cons(w, equations::CompressibleEulerMulticomponentEquations1D)535@unpack gammas, gas_constants, cv, cp = equations536T = -1 / w[2]537v1 = w[1] * T538cons_rho = SVector{ncomponents(equations), real(equations)}(exp(1 /539gas_constants[i] *540(-cv[i] *541log(-w[2]) -542cp[i] + w[i + 2] -5430.5f0 * w[1]^2 /544w[2]))545for i in eachcomponent(equations))546547RealT = eltype(w)548rho = zero(RealT)549help1 = zero(RealT)550help2 = zero(RealT)551p = zero(RealT)552for i in eachcomponent(equations)553rho += cons_rho[i]554help1 += cons_rho[i] * cv[i] * gammas[i]555help2 += cons_rho[i] * cv[i]556p += cons_rho[i] * gas_constants[i] * T557end558u1 = rho * v1559gamma = help1 / help2560u2 = p / (gamma - 1) + 0.5f0 * rho * v1^2561cons_other = SVector(u1, u2)562return vcat(cons_other, cons_rho)563end564565@inline function total_entropy(u, equations::CompressibleEulerMulticomponentEquations1D)566@unpack cv, gammas, gas_constants = equations567rho_v1, rho_e_total = u568rho = density(u, equations)569T = temperature(u, equations)570571RealT = eltype(u)572total_entropy = zero(RealT)573for i in eachcomponent(equations)574total_entropy -= u[i + 2] * (cv[i] * log(T) - gas_constants[i] * log(u[i + 2]))575end576577return total_entropy578end579580@inline function temperature(u, equations::CompressibleEulerMulticomponentEquations1D)581@unpack cv, gammas, gas_constants = equations582583rho_v1, rho_e_total = u584585rho = density(u, equations)586587RealT = eltype(u)588help1 = zero(RealT)589590for i in eachcomponent(equations)591help1 += u[i + 2] * cv[i]592end593594v1 = rho_v1 / rho595v_square = v1^2596T = (rho_e_total - 0.5f0 * rho * v_square) / help1597598return T599end600601"""602totalgamma(u, equations::CompressibleEulerMulticomponentEquations1D)603604Function that calculates the total gamma out of all partial gammas using the605partial density fractions as well as the partial specific heats at constant volume.606"""607@inline function totalgamma(u, equations::CompressibleEulerMulticomponentEquations1D)608@unpack cv, gammas = equations609610RealT = eltype(u)611help1 = zero(RealT)612help2 = zero(RealT)613614for i in eachcomponent(equations)615help1 += u[i + 2] * cv[i] * gammas[i]616help2 += u[i + 2] * cv[i]617end618619return help1 / help2620end621622@doc raw"""623pressure(u, equations::AbstractCompressibleEulerMulticomponentEquations)624625Computes the pressure for an ideal equation of state with626isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`.627```math628\begin{aligned}629p &= (\gamma - 1) \left( E_\mathrm{tot} - E_\mathrm{kin} \right) \\630&= (\gamma - 1) \left( \rho e - \frac{1}{2}\rho \Vert v \Vert_2^2 \right)631\end{aligned}632```633"""634@inline function pressure(u, equations::CompressibleEulerMulticomponentEquations1D)635rho_v1, rho_e_total = u636637rho = density(u, equations)638gamma = totalgamma(u, equations)639640p = (gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1^2) / rho)641642return p643end644645@inline function density_pressure(u,646equations::CompressibleEulerMulticomponentEquations1D)647rho_v1, rho_e_total = u648649rho = density(u, equations)650gamma = totalgamma(u, equations)651rho_times_p = (gamma - 1) * (rho * rho_e_total - 0.5f0 * rho_v1^2)652653return rho_times_p654end655656@doc raw"""657density(u, equations::AbstractCompressibleEulerMulticomponentEquations)658659Computes the total density ``\rho = \sum_{i=1}^n \rho_i`` from the conserved variables `u`.660"""661@inline function density(u, equations::CompressibleEulerMulticomponentEquations1D)662RealT = eltype(u)663rho = zero(RealT)664665for i in eachcomponent(equations)666rho += u[i + 2]667end668669return rho670end671672@inline function densities(u, v, equations::CompressibleEulerMulticomponentEquations1D)673return SVector{ncomponents(equations), real(equations)}(u[i + 2] * v674for i in eachcomponent(equations))675end676677@inline function velocity(u, equations::CompressibleEulerMulticomponentEquations1D)678rho = density(u, equations)679v1 = u[1] / rho680return v1681end682end # @muladd683684685