Path: blob/main/src/equations/compressible_euler_quasi_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"""8CompressibleEulerEquationsQuasi1D(gamma)910The quasi-1d compressible Euler equations (see Chan et al. [DOI: 10.48550/arXiv.2307.12089](https://doi.org/10.48550/arXiv.2307.12089) for details)11```math12\frac{\partial}{\partial t}13\begin{pmatrix}14a \rho \\ a \rho v_1 \\ a e_{\text{total}}15\end{pmatrix}16+17\frac{\partial}{\partial x}18\begin{pmatrix}19a \rho v_1 \\ a \rho v_1^2 \\ a v_1 (e_{\text{total}} +p)20\end{pmatrix}21+22a \frac{\partial}{\partial x}23\begin{pmatrix}240 \\ p \\ 025\end{pmatrix}26=27\begin{pmatrix}280 \\ 0 \\ 029\end{pmatrix}30```31for an ideal gas with ratio of specific heats `gamma` in one space dimension.32Here, ``\rho`` is the density, ``v_1`` the velocity, ``e_{\text{total}}`` the specific total energy,33``a`` the (possibly) variable nozzle width, and34```math35p = (\gamma - 1) \left( e_{\text{total}} - \frac{1}{2} \rho v_1^2 \right)36```37the pressure.3839The nozzle width function ``a(x)`` is set inside the initial condition routine40for a particular problem setup. To test the conservative form of the compressible Euler equations one can set the41nozzle width variable ``a`` to one.4243In addition to the unknowns, Trixi.jl currently stores the nozzle width values at the approximation points44despite being fixed in time.45This affects the implementation and use of these equations in various ways:46* The flux values corresponding to the nozzle width must be zero.47* The nozzle width values must be included when defining initial conditions, boundary conditions or48source terms.49* [`AnalysisCallback`](@ref) analyzes this variable.50* Trixi.jl's visualization tools will visualize the nozzle width by default.51"""52struct CompressibleEulerEquationsQuasi1D{RealT <: Real} <:53AbstractCompressibleEulerEquations{1, 4}54gamma::RealT # ratio of specific heats55inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications5657function CompressibleEulerEquationsQuasi1D(gamma)58γ, inv_gamma_minus_one = promote(gamma, inv(gamma - 1))59return new{typeof(γ)}(γ, inv_gamma_minus_one)60end61end6263have_nonconservative_terms(::CompressibleEulerEquationsQuasi1D) = True()64function varnames(::typeof(cons2cons), ::CompressibleEulerEquationsQuasi1D)65return ("a_rho", "a_rho_v1", "a_e_total", "a")66end67function varnames(::typeof(cons2prim), ::CompressibleEulerEquationsQuasi1D)68return ("rho", "v1", "p", "a")69end7071"""72initial_condition_convergence_test(x, t, equations::CompressibleEulerEquationsQuasi1D)7374A smooth initial condition used for convergence tests in combination with75[`source_terms_convergence_test`](@ref)76(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).77"""78function initial_condition_convergence_test(x, t,79equations::CompressibleEulerEquationsQuasi1D)80RealT = eltype(x)81c = 282A = convert(RealT, 0.1)83L = 284f = 1.0f0 / L85ω = 2 * convert(RealT, pi) * f86ini = c + A * sin(ω * (x[1] - t))8788rho = ini89v1 = 190e_total = ini^2 / rho91p = (equations.gamma - 1) * (e_total - 0.5f0 * rho * v1^2)92a = 1.5f0 - 0.5f0 * cospi(x[1])9394return prim2cons(SVector(rho, v1, p, a), equations)95end9697"""98source_terms_convergence_test(u, x, t, equations::CompressibleEulerEquationsQuasi1D)99100Source terms used for convergence tests in combination with101[`initial_condition_convergence_test`](@ref)102(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).103104This manufactured solution source term is specifically designed for the nozzle width105```math106a(x) = 1.5 - 0.5 \\cos(x \\pi)107```108as defined in [`initial_condition_convergence_test`](@ref).109110References for the method of manufactured solutions (MMS):111- Kambiz Salari and Patrick Knupp (2000)112Code Verification by the Method of Manufactured Solutions113[DOI: 10.2172/759450](https://doi.org/10.2172/759450)114- Patrick J. Roache (2002)115Code Verification by the Method of Manufactured Solutions116[DOI: 10.1115/1.1436090](https://doi.org/10.1115/1.1436090)117"""118@inline function source_terms_convergence_test(u, x, t,119equations::CompressibleEulerEquationsQuasi1D)120# Same settings as in `initial_condition_convergence_test`.121# Derivatives calculated with ForwardDiff.jl122RealT = eltype(u)123c = 2124A = convert(RealT, 0.1)125L = 2126f = 1.0f0 / L127ω = 2 * convert(RealT, pi) * f128x1, = x129ini(x1, t) = c + A * sin(ω * (x1 - t))130131rho(x1, t) = ini(x1, t)132v1(x1, t) = 1133e_total(x1, t) = ini(x1, t)^2 / rho(x1, t)134p1(x1, t) = (equations.gamma - 1) *135(e_total(x1, t) - 0.5f0 * rho(x1, t) * v1(x1, t)^2)136a(x1, t) = 1.5f0 - 0.5f0 * cospi(x1)137138arho(x1, t) = a(x1, t) * rho(x1, t)139arhou(x1, t) = arho(x1, t) * v1(x1, t)140aE(x1, t) = a(x1, t) * e_total(x1, t)141142darho_dt(x1, t) = ForwardDiff.derivative(t -> arho(x1, t), t)143darhou_dx(x1, t) = ForwardDiff.derivative(x1 -> arhou(x1, t), x1)144145arhouu(x1, t) = arhou(x1, t) * v1(x1, t)146darhou_dt(x1, t) = ForwardDiff.derivative(t -> arhou(x1, t), t)147darhouu_dx(x1, t) = ForwardDiff.derivative(x1 -> arhouu(x1, t), x1)148dp1_dx(x1, t) = ForwardDiff.derivative(x1 -> p1(x1, t), x1)149150auEp(x1, t) = a(x1, t) * v1(x1, t) * (e_total(x1, t) + p1(x1, t))151daE_dt(x1, t) = ForwardDiff.derivative(t -> aE(x1, t), t)152dauEp_dx(x1, t) = ForwardDiff.derivative(x1 -> auEp(x1, t), x1)153154du1 = darho_dt(x1, t) + darhou_dx(x1, t)155du2 = darhou_dt(x1, t) + darhouu_dx(x1, t) + a(x1, t) * dp1_dx(x1, t)156du3 = daE_dt(x1, t) + dauEp_dx(x1, t)157158return SVector(du1, du2, du3, 0)159end160161# Calculate 1D flux for a single point162@inline function flux(u, orientation::Integer,163equations::CompressibleEulerEquationsQuasi1D)164a_rho, a_rho_v1, a_e_total, a = u165rho, v1, p, a = cons2prim(u, equations)166e = a_e_total / a167168# Ignore orientation since it is always "1" in 1D169f1 = a_rho_v1170f2 = a_rho_v1 * v1171f3 = a * v1 * (e + p)172173return SVector(f1, f2, f3, 0)174end175176"""177flux_nonconservative_chan_etal(u_ll, u_rr, orientation::Integer,178equations::CompressibleEulerEquationsQuasi1D)179flux_nonconservative_chan_etal(u_ll, u_rr, normal_direction,180equations::CompressibleEulerEquationsQuasi1D)181flux_nonconservative_chan_etal(u_ll, u_rr, normal_ll, normal_rr,182equations::CompressibleEulerEquationsQuasi1D)183184Non-symmetric two-point volume flux discretizing the nonconservative (source) term185that contains the gradient of the pressure [`CompressibleEulerEquationsQuasi1D`](@ref)186and the nozzle width.187188Further details are available in the paper:189- Jesse Chan, Khemraj Shukla, Xinhui Wu, Ruofeng Liu, Prani Nalluri (2023)190High order entropy stable schemes for the quasi-one-dimensional191shallow water and compressible Euler equations192[DOI: 10.48550/arXiv.2307.12089](https://doi.org/10.48550/arXiv.2307.12089)193"""194@inline function flux_nonconservative_chan_etal(u_ll, u_rr, orientation::Integer,195equations::CompressibleEulerEquationsQuasi1D)196#Variables197_, _, p_ll, a_ll = cons2prim(u_ll, equations)198_, _, p_rr, _ = cons2prim(u_rr, equations)199200# For flux differencing using non-conservative terms, we return the201# non-conservative flux scaled by 2. This cancels with a factor of 0.5202# in the arithmetic average of {p}.203p_avg = p_ll + p_rr204205return SVector(0, a_ll * p_avg, 0, 0)206end207208# While `normal_direction` isn't strictly necessary in 1D, certain solvers assume that209# the normal component is incorporated into the numerical flux.210#211# See `flux(u, normal_direction::AbstractVector, equations::AbstractEquations{1})` for a212# similar implementation.213@inline function flux_nonconservative_chan_etal(u_ll, u_rr,214normal_direction::AbstractVector,215equations::CompressibleEulerEquationsQuasi1D)216return normal_direction[1] *217flux_nonconservative_chan_etal(u_ll, u_rr, 1, equations)218end219220@inline function flux_nonconservative_chan_etal(u_ll, u_rr,221normal_ll::AbstractVector,222normal_rr::AbstractVector,223equations::CompressibleEulerEquationsQuasi1D)224# normal_ll should be equal to normal_rr in 1D225return flux_nonconservative_chan_etal(u_ll, u_rr, normal_ll, equations)226end227228"""229@inline function flux_chan_etal(u_ll, u_rr, orientation::Integer,230equations::CompressibleEulerEquationsQuasi1D)231232Conservative (symmetric) part of the entropy conservative flux for quasi 1D compressible Euler equations split form.233This flux is a generalization of [`flux_ranocha`](@ref) for [`CompressibleEulerEquations1D`](@ref).234Further details are available in the paper:235- Jesse Chan, Khemraj Shukla, Xinhui Wu, Ruofeng Liu, Prani Nalluri (2023)236High order entropy stable schemes for the quasi-one-dimensional237shallow water and compressible Euler equations238[DOI: 10.48550/arXiv.2307.12089](https://doi.org/10.48550/arXiv.2307.12089)239"""240@inline function flux_chan_etal(u_ll, u_rr, orientation::Integer,241equations::CompressibleEulerEquationsQuasi1D)242# Unpack left and right state243rho_ll, v1_ll, p_ll, a_ll = cons2prim(u_ll, equations)244rho_rr, v1_rr, p_rr, a_rr = cons2prim(u_rr, equations)245246# Compute the necessary mean values247rho_mean = ln_mean(rho_ll, rho_rr)248# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`249# in exact arithmetic since250# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)251# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)252inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)253v1_avg = 0.5f0 * (v1_ll + v1_rr)254a_v1_avg = 0.5f0 * (a_ll * v1_ll + a_rr * v1_rr)255p_avg = 0.5f0 * (p_ll + p_rr)256velocity_square_avg = 0.5f0 * (v1_ll * v1_rr)257258# Calculate fluxes259# Ignore orientation since it is always "1" in 1D260f1 = rho_mean * a_v1_avg261f2 = rho_mean * a_v1_avg * v1_avg262f3 = f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +2630.5f0 * (p_ll * a_rr * v1_rr + p_rr * a_ll * v1_ll)264265return SVector(f1, f2, f3, 0)266end267268# While `normal_direction` isn't strictly necessary in 1D, certain solvers assume that269# the normal component is incorporated into the numerical flux.270#271# See `flux(u, normal_direction::AbstractVector, equations::AbstractEquations{1})` for a272# similar implementation.273@inline function flux_chan_etal(u_ll, u_rr, normal_direction::AbstractVector,274equations::CompressibleEulerEquationsQuasi1D)275return normal_direction[1] * flux_chan_etal(u_ll, u_rr, 1, equations)276end277278# Calculate estimates for maximum wave speed for local Lax-Friedrichs-type dissipation as the279# maximum velocity magnitude plus the maximum speed of sound280@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,281equations::CompressibleEulerEquationsQuasi1D)282a_rho_ll, a_rho_v1_ll, a_e_total_ll, a_ll = u_ll283a_rho_rr, a_rho_v1_rr, a_e_total_rr, a_rr = u_rr284285# Calculate primitive variables and speed of sound286rho_ll = a_rho_ll / a_ll287e_ll = a_e_total_ll / a_ll288v1_ll = a_rho_v1_ll / a_rho_ll289v_mag_ll = abs(v1_ll)290p_ll = (equations.gamma - 1) * (e_ll - 0.5f0 * rho_ll * v_mag_ll^2)291c_ll = sqrt(equations.gamma * p_ll / rho_ll)292rho_rr = a_rho_rr / a_rr293e_rr = a_e_total_rr / a_rr294v1_rr = a_rho_v1_rr / a_rho_rr295v_mag_rr = abs(v1_rr)296p_rr = (equations.gamma - 1) * (e_rr - 0.5f0 * rho_rr * v_mag_rr^2)297c_rr = sqrt(equations.gamma * p_rr / rho_rr)298299return max(v_mag_ll, v_mag_rr) + max(c_ll, c_rr)300end301302# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`303@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,304equations::CompressibleEulerEquationsQuasi1D)305a_rho_ll, a_rho_v1_ll, a_e_total_ll, a_ll = u_ll306a_rho_rr, a_rho_v1_rr, a_e_total_rr, a_rr = u_rr307308# Calculate primitive variables and speed of sound309rho_ll = a_rho_ll / a_ll310e_ll = a_e_total_ll / a_ll311v1_ll = a_rho_v1_ll / a_rho_ll312v_mag_ll = abs(v1_ll)313p_ll = (equations.gamma - 1) * (e_ll - 0.5f0 * rho_ll * v_mag_ll^2)314c_ll = sqrt(equations.gamma * p_ll / rho_ll)315rho_rr = a_rho_rr / a_rr316e_rr = a_e_total_rr / a_rr317v1_rr = a_rho_v1_rr / a_rho_rr318v_mag_rr = abs(v1_rr)319p_rr = (equations.gamma - 1) * (e_rr - 0.5f0 * rho_rr * v_mag_rr^2)320c_rr = sqrt(equations.gamma * p_rr / rho_rr)321322return max(v_mag_ll + c_ll, v_mag_rr + c_rr)323end324325@inline function max_abs_speeds(u, equations::CompressibleEulerEquationsQuasi1D)326a_rho, a_rho_v1, a_e_total, a = u327rho = a_rho / a328v1 = a_rho_v1 / a_rho329e = a_e_total / a330p = (equations.gamma - 1) * (e - 0.5f0 * rho * v1^2)331c = sqrt(equations.gamma * p / rho)332333return (abs(v1) + c,)334end335336# Convert conservative variables to primitive. We use the convention that the primitive337# variables for the quasi-1D equations are `(rho, v1, p)` (i.e., the same as the primitive338# variables for `CompressibleEulerEquations1D`)339@inline function cons2prim(u, equations::CompressibleEulerEquationsQuasi1D)340a_rho, a_rho_v1, a_e_total, a = u341q = cons2prim(SVector(a_rho, a_rho_v1, a_e_total) / a,342CompressibleEulerEquations1D(equations.gamma))343344return SVector(q[1], q[2], q[3], a)345end346347"""348entropy(u, equations::CompressibleEulerEquationsQuasi1D)349350The entropy for the quasi-1D compressible Euler equations is the351[`entropy(cons, equations::AbstractCompressibleEulerEquations)`](@ref) for the352(1D) compressible Euler equations scaled by the channel width `a`.353"""354@inline function entropy(u, equations::CompressibleEulerEquationsQuasi1D)355a_rho, a_rho_v1, a_e_total, a = u356return a * entropy(SVector(a_rho, a_rho_v1, a_e_total) / a,357CompressibleEulerEquations1D(equations.gamma))358end359360# Convert conservative variables to entropy. The entropy variables for the361# quasi-1D compressible Euler equations are identical to the entropy variables362# for the standard Euler equations for an appropriate definition of `entropy`.363@inline function cons2entropy(u, equations::CompressibleEulerEquationsQuasi1D)364a_rho, a_rho_v1, a_e_total, a = u365w = cons2entropy(SVector(a_rho, a_rho_v1, a_e_total) / a,366CompressibleEulerEquations1D(equations.gamma))367368# we follow the convention for other spatially-varying equations and return the spatially369# varying coefficient `a` as the final entropy variable.370return SVector(w[1], w[2], w[3], a)371end372373@inline function entropy2cons(w, equations::CompressibleEulerEquationsQuasi1D)374w_rho, w_rho_v1, w_rho_e_total, a = w375u = entropy2cons(SVector(w_rho, w_rho_v1, w_rho_e_total),376CompressibleEulerEquations1D(equations.gamma))377return SVector(a * u[1], a * u[2], a * u[3], a)378end379380# Convert primitive to conservative variables381@inline function prim2cons(u, equations::CompressibleEulerEquationsQuasi1D)382rho, v1, p, a = u383q = prim2cons(u, CompressibleEulerEquations1D(equations.gamma))384385return SVector(a * q[1], a * q[2], a * q[3], a)386end387388@inline function density(u, equations::CompressibleEulerEquationsQuasi1D)389a_rho, _, _, a = u390rho = a_rho / a391return rho392end393394@doc raw"""395pressure(u, equations::CompressibleEulerEquationsQuasi1D)396397Computes the pressure for an ideal equation of state with398isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`,399see [`pressure(u, equations::CompressibleEulerEquations1D)`](@ref).400"""401@inline function pressure(u, equations::CompressibleEulerEquationsQuasi1D)402a_rho, a_rho_v1, a_e_total, a = u403return pressure(SVector(a_rho, a_rho_v1, a_e_total) / a,404CompressibleEulerEquations1D(equations.gamma))405end406407@inline function density_pressure(u, equations::CompressibleEulerEquationsQuasi1D)408a_rho, a_rho_v1, a_e_total, a = u409return density_pressure(SVector(a_rho, a_rho_v1, a_e_total) / a,410CompressibleEulerEquations1D(equations.gamma))411end412end # @muladd413414415