Path: blob/main/src/equations/linearized_euler_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"""8LinearizedEulerEquations2D(v_mean_global, c_mean_global, rho_mean_global)910Linearized Euler equations in two space dimensions. The equations are given by11```math12\partial_t13\begin{pmatrix}14\rho' \\ v_1' \\ v_2' \\ p'15\end{pmatrix}16+17\partial_x18\begin{pmatrix}19\bar{\rho} v_1' + \bar{v_1} \rho ' \\ \bar{v_1} v_1' + \frac{p'}{\bar{\rho}} \\ \bar{v_1} v_2' \\ \bar{v_1} p' + c^2 \bar{\rho} v_1'20\end{pmatrix}21+22\partial_y23\begin{pmatrix}24\bar{\rho} v_2' + \bar{v_2} \rho ' \\ \bar{v_2} v_1' \\ \bar{v_2} v_2' + \frac{p'}{\bar{\rho}} \\ \bar{v_2} p' + c^2 \bar{\rho} v_2'25\end{pmatrix}26=27\begin{pmatrix}280 \\ 0 \\ 0 \\ 029\end{pmatrix}30```31The bar ``\bar{(\cdot)}`` indicates uniform mean flow variables and ``c`` is the speed of sound.32The unknowns are the acoustic velocities ``v' = (v_1', v_2')``, the pressure ``p'`` and the density ``\rho'``.33"""34struct LinearizedEulerEquations2D{RealT <: Real} <:35AbstractLinearizedEulerEquations{2, 4}36v_mean_global::SVector{2, RealT}37c_mean_global::RealT38rho_mean_global::RealT39end4041function LinearizedEulerEquations2D(v_mean_global::NTuple{2, <:Real},42c_mean_global::Real, rho_mean_global::Real)43if rho_mean_global < 044throw(ArgumentError("rho_mean_global must be non-negative"))45elseif c_mean_global < 046throw(ArgumentError("c_mean_global must be non-negative"))47end4849return LinearizedEulerEquations2D(SVector(v_mean_global), c_mean_global,50rho_mean_global)51end5253function LinearizedEulerEquations2D(; v_mean_global::NTuple{2, <:Real},54c_mean_global::Real, rho_mean_global::Real)55return LinearizedEulerEquations2D(v_mean_global, c_mean_global,56rho_mean_global)57end5859function varnames(::typeof(cons2cons), ::LinearizedEulerEquations2D)60return ("rho_prime", "v1_prime", "v2_prime", "p_prime")61end62function varnames(::typeof(cons2prim), ::LinearizedEulerEquations2D)63return ("rho_prime", "v1_prime", "v2_prime", "p_prime")64end6566"""67initial_condition_convergence_test(x, t, equations::LinearizedEulerEquations2D)6869A smooth initial condition used for convergence tests.70"""71function initial_condition_convergence_test(x, t, equations::LinearizedEulerEquations2D)72rho_prime = -cospi(2 * t) * (sinpi(2 * x[1]) + sinpi(2 * x[2]))73v1_prime = sinpi(2 * t) * cospi(2 * x[1])74v2_prime = sinpi(2 * t) * cospi(2 * x[2])75p_prime = rho_prime7677return SVector(rho_prime, v1_prime, v2_prime, p_prime)78end7980"""81boundary_condition_wall(u_inner, orientation, direction, x, t, surface_flux_function,82equations::LinearizedEulerEquations2D)8384Boundary conditions for a solid wall.85"""86function boundary_condition_wall(u_inner, orientation, direction, x, t,87surface_flux_function,88equations::LinearizedEulerEquations2D)89# Boundary state is equal to the inner state except for the velocity. For boundaries90# in the -x/+x direction, we multiply the velocity in the x direction by -1.91# Similarly, for boundaries in the -y/+y direction, we multiply the velocity in the92# y direction by -193if direction in (1, 2) # x direction94u_boundary = SVector(u_inner[1], -u_inner[2], u_inner[3], u_inner[4])95else # y direction96u_boundary = SVector(u_inner[1], u_inner[2], -u_inner[3], u_inner[4])97end9899# Calculate boundary flux100if iseven(direction) # u_inner is "left" of boundary, u_boundary is "right" of boundary101flux = surface_flux_function(u_inner, u_boundary, orientation, equations)102else # u_boundary is "left" of boundary, u_inner is "right" of boundary103flux = surface_flux_function(u_boundary, u_inner, orientation, equations)104end105106return flux107end108109# Calculate 1D flux for a single point110@inline function flux(u, orientation::Integer, equations::LinearizedEulerEquations2D)111@unpack v_mean_global, c_mean_global, rho_mean_global = equations112rho_prime, v1_prime, v2_prime, p_prime = u113if orientation == 1114f1 = v_mean_global[1] * rho_prime + rho_mean_global * v1_prime115f2 = v_mean_global[1] * v1_prime + p_prime / rho_mean_global116f3 = v_mean_global[1] * v2_prime117f4 = v_mean_global[1] * p_prime + c_mean_global^2 * rho_mean_global * v1_prime118else119f1 = v_mean_global[2] * rho_prime + rho_mean_global * v2_prime120f2 = v_mean_global[2] * v1_prime121f3 = v_mean_global[2] * v2_prime + p_prime / rho_mean_global122f4 = v_mean_global[2] * p_prime + c_mean_global^2 * rho_mean_global * v2_prime123end124125return SVector(f1, f2, f3, f4)126end127128# Calculate 1D flux for a single point129@inline function flux(u, normal_direction::AbstractVector,130equations::LinearizedEulerEquations2D)131@unpack v_mean_global, c_mean_global, rho_mean_global = equations132rho_prime, v1_prime, v2_prime, p_prime = u133134v_mean_normal = v_mean_global[1] * normal_direction[1] +135v_mean_global[2] * normal_direction[2]136v_prime_normal = v1_prime * normal_direction[1] + v2_prime * normal_direction[2]137138f1 = v_mean_normal * rho_prime + rho_mean_global * v_prime_normal139f2 = v_mean_normal * v1_prime + normal_direction[1] * p_prime / rho_mean_global140f3 = v_mean_normal * v2_prime + normal_direction[2] * p_prime / rho_mean_global141f4 = v_mean_normal * p_prime + c_mean_global^2 * rho_mean_global * v_prime_normal142143return SVector(f1, f2, f3, f4)144end145146"""147have_constant_speed(::LinearizedEulerEquations2D)148149Indicates whether the characteristic speeds are constant, i.e., independent of the solution.150Queried in the timestep computation [`StepsizeCallback`](@ref) and [`linear_structure`](@ref).151152# Returns153- `True()`154"""155@inline have_constant_speed(::LinearizedEulerEquations2D) = True()156157@inline function max_abs_speeds(equations::LinearizedEulerEquations2D)158@unpack v_mean_global, c_mean_global = equations159return abs(v_mean_global[1]) + c_mean_global, abs(v_mean_global[2]) + c_mean_global160end161162@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,163equations::LinearizedEulerEquations2D)164@unpack v_mean_global, c_mean_global = equations165if orientation == 1166return abs(v_mean_global[1]) + c_mean_global167else # orientation == 2168return abs(v_mean_global[2]) + c_mean_global169end170end171172@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,173equations::LinearizedEulerEquations2D)174@unpack v_mean_global, c_mean_global = equations175v_mean_normal = normal_direction[1] * v_mean_global[1] +176normal_direction[2] * v_mean_global[2]177return abs(v_mean_normal) + c_mean_global * norm(normal_direction)178end179180@doc raw"""181flux_godunov(u_ll, u_rr, orientation_or_normal_direction,182equations::LinearizedEulerEquations2D)183184An upwind flux for the linearized Euler equations based on diagonalization of the physical185flux matrix. Given the physical flux ``Au``, ``A=T \Lambda T^{-1}`` with186``\Lambda`` being a diagonal matrix that holds the eigenvalues of ``A``, decompose187``\Lambda = \Lambda^+ + \Lambda^-`` where ``\Lambda^+`` and ``\Lambda^-`` are diagonal188matrices holding the positive and negative eigenvalues of ``A``, respectively. Then for189left and right states ``u_L, u_R``, the numerical flux calculated by this function is given190by ``A^+ u_L + A^- u_R`` where ``A^{\pm} = T \Lambda^{\pm} T^{-1}``.191192The diagonalization of the flux matrix can be found in193- R. F. Warming, Richard M. Beam and B. J. Hyett (1975)194Diagonalization and simultaneous symmetrization of the gas-dynamic matrices195[DOI: 10.1090/S0025-5718-1975-0388967-5](https://doi.org/10.1090/S0025-5718-1975-0388967-5)196"""197@inline function flux_godunov(u_ll, u_rr, orientation::Integer,198equations::LinearizedEulerEquations2D)199@unpack v_mean_global, rho_mean_global, c_mean_global = equations200v1_mean = v_mean_global[1]201v2_mean = v_mean_global[2]202203rho_prime_ll, v1_prime_ll, v2_prime_ll, p_prime_ll = u_ll204rho_prime_rr, v1_prime_rr, v2_prime_rr, p_prime_rr = u_rr205206if orientation == 1207# Eigenvalues of the flux matrix208lambda1 = v1_mean209lambda2 = v1_mean - c_mean_global210lambda3 = v1_mean + c_mean_global211212lambda1_p = positive_part(lambda1)213lambda2_p = positive_part(lambda2)214lambda3_p = positive_part(lambda3)215lambda2p3_half_p = 0.5f0 * (lambda2_p + lambda3_p)216lambda3m2_half_p = 0.5f0 * (lambda3_p - lambda2_p)217218lambda1_m = negative_part(lambda1)219lambda2_m = negative_part(lambda2)220lambda3_m = negative_part(lambda3)221lambda2p3_half_m = 0.5f0 * (lambda2_m + lambda3_m)222lambda3m2_half_m = 0.5f0 * (lambda3_m - lambda2_m)223224f1p = (lambda1_p * rho_prime_ll +225lambda3m2_half_p / c_mean_global * rho_mean_global * v1_prime_ll +226(lambda2p3_half_p - lambda1_p) / c_mean_global^2 * p_prime_ll)227f2p = (lambda2p3_half_p * v1_prime_ll +228lambda3m2_half_p / c_mean_global * p_prime_ll / rho_mean_global)229f3p = lambda1_p * v2_prime_ll230f4p = (lambda3m2_half_p * c_mean_global * rho_mean_global * v1_prime_ll +231lambda2p3_half_p * p_prime_ll)232233f1m = (lambda1_m * rho_prime_rr +234lambda3m2_half_m / c_mean_global * rho_mean_global * v1_prime_rr +235(lambda2p3_half_m - lambda1_m) / c_mean_global^2 * p_prime_rr)236f2m = (lambda2p3_half_m * v1_prime_rr +237lambda3m2_half_m / c_mean_global * p_prime_rr / rho_mean_global)238f3m = lambda1_m * v2_prime_rr239f4m = (lambda3m2_half_m * c_mean_global * rho_mean_global * v1_prime_rr +240lambda2p3_half_m * p_prime_rr)241242f1 = f1p + f1m243f2 = f2p + f2m244f3 = f3p + f3m245f4 = f4p + f4m246else # orientation == 2247# Eigenvalues of the flux matrix248lambda1 = v2_mean249lambda2 = v2_mean - c_mean_global250lambda3 = v2_mean + c_mean_global251252lambda1_p = positive_part(lambda1)253lambda2_p = positive_part(lambda2)254lambda3_p = positive_part(lambda3)255lambda2p3_half_p = 0.5f0 * (lambda2_p + lambda3_p)256lambda3m2_half_p = 0.5f0 * (lambda3_p - lambda2_p)257258lambda1_m = negative_part(lambda1)259lambda2_m = negative_part(lambda2)260lambda3_m = negative_part(lambda3)261lambda2p3_half_m = 0.5f0 * (lambda2_m + lambda3_m)262lambda3m2_half_m = 0.5f0 * (lambda3_m - lambda2_m)263264f1p = (lambda1_p * rho_prime_ll +265lambda3m2_half_p / c_mean_global * rho_mean_global * v2_prime_ll +266(lambda2p3_half_p - lambda1_p) / c_mean_global^2 * p_prime_ll)267f2p = lambda1_p * v1_prime_ll268f3p = (lambda2p3_half_p * v2_prime_ll +269lambda3m2_half_p / c_mean_global * p_prime_ll / rho_mean_global)270f4p = (lambda3m2_half_p * c_mean_global * rho_mean_global * v2_prime_ll +271lambda2p3_half_p * p_prime_ll)272273f1m = (lambda1_m * rho_prime_rr +274lambda3m2_half_m / c_mean_global * rho_mean_global * v2_prime_rr +275(lambda2p3_half_m - lambda1_m) / c_mean_global^2 * p_prime_rr)276f2m = lambda1_m * v1_prime_rr277f3m = (lambda2p3_half_m * v2_prime_rr +278lambda3m2_half_m / c_mean_global * p_prime_rr / rho_mean_global)279f4m = (lambda3m2_half_m * c_mean_global * rho_mean_global * v2_prime_rr +280lambda2p3_half_m * p_prime_rr)281282f1 = f1p + f1m283f2 = f2p + f2m284f3 = f3p + f3m285f4 = f4p + f4m286end287288return SVector(f1, f2, f3, f4)289end290291@inline function flux_godunov(u_ll, u_rr, normal_direction::AbstractVector,292equations::LinearizedEulerEquations2D)293@unpack v_mean_global, rho_mean_global, c_mean_global = equations294rho_prime_ll, v1_prime_ll, v2_prime_ll, p_prime_ll = u_ll295rho_prime_rr, v1_prime_rr, v2_prime_rr, p_prime_rr = u_rr296297# Do not use `normalize` since we use `norm_` later to scale the eigenvalues298norm_ = norm(normal_direction)299normal_vector = normal_direction / norm_300301# Use normalized vector here, scaling is applied via eigenvalues of the flux matrix302v_mean_normal = v_mean_global[1] * normal_vector[1] +303v_mean_global[2] * normal_vector[2]304v_prime_normal_ll = v1_prime_ll * normal_vector[1] + v2_prime_ll * normal_vector[2]305v_prime_normal_rr = v1_prime_rr * normal_vector[1] + v2_prime_rr * normal_vector[2]306307# Eigenvalues of the flux matrix308lambda1 = v_mean_normal * norm_309lambda2 = (v_mean_normal - c_mean_global) * norm_310lambda3 = (v_mean_normal + c_mean_global) * norm_311312lambda1_p = positive_part(lambda1)313lambda2_p = positive_part(lambda2)314lambda3_p = positive_part(lambda3)315lambda2p3_half_p = 0.5f0 * (lambda2_p + lambda3_p)316lambda3m2_half_p = 0.5f0 * (lambda3_p - lambda2_p)317318lambda1_m = negative_part(lambda1)319lambda2_m = negative_part(lambda2)320lambda3_m = negative_part(lambda3)321lambda2p3_half_m = 0.5f0 * (lambda2_m + lambda3_m)322lambda3m2_half_m = 0.5f0 * (lambda3_m - lambda2_m)323324f1p = (lambda1_p * rho_prime_ll +325lambda3m2_half_p / c_mean_global * rho_mean_global * v_prime_normal_ll +326(lambda2p3_half_p - lambda1_p) / c_mean_global^2 * p_prime_ll)327f2p = (((lambda1_p * normal_vector[2]^2 +328lambda2p3_half_p * normal_vector[1]^2) * v1_prime_ll +329(lambda2p3_half_p - lambda1_p) * prod(normal_vector) * v2_prime_ll) +330lambda3m2_half_p / c_mean_global * normal_vector[1] * p_prime_ll /331rho_mean_global)332f3p = (((lambda1_p * normal_vector[1]^2 +333lambda2p3_half_p * normal_vector[2]^2) * v2_prime_ll +334(lambda2p3_half_p - lambda1_p) * prod(normal_vector) * v1_prime_ll) +335lambda3m2_half_p / c_mean_global * normal_vector[2] * p_prime_ll /336rho_mean_global)337f4p = (lambda3m2_half_p * c_mean_global * rho_mean_global * v_prime_normal_ll +338lambda2p3_half_p * p_prime_ll)339340f1m = (lambda1_m * rho_prime_rr +341lambda3m2_half_m / c_mean_global * rho_mean_global * v_prime_normal_rr +342(lambda2p3_half_m - lambda1_m) / c_mean_global^2 * p_prime_rr)343f2m = (((lambda1_m * normal_vector[2]^2 +344lambda2p3_half_m * normal_vector[1]^2) * v1_prime_rr +345(lambda2p3_half_m - lambda1_m) * prod(normal_vector) * v2_prime_rr) +346lambda3m2_half_m / c_mean_global * normal_vector[1] * p_prime_rr /347rho_mean_global)348f3m = (((lambda1_m * normal_vector[1]^2 +349lambda2p3_half_m * normal_vector[2]^2) * v2_prime_rr +350(lambda2p3_half_m - lambda1_m) * prod(normal_vector) * v1_prime_rr) +351lambda3m2_half_m / c_mean_global * normal_vector[2] * p_prime_rr /352rho_mean_global)353f4m = (lambda3m2_half_m * c_mean_global * rho_mean_global * v_prime_normal_rr +354lambda2p3_half_m * p_prime_rr)355356f1 = f1p + f1m357f2 = f2p + f2m358f3 = f3p + f3m359f4 = f4p + f4m360361return SVector(f1, f2, f3, f4)362end363364# Calculate estimate for minimum and maximum wave speeds for HLL-type fluxes365@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,366equations::LinearizedEulerEquations2D)367return min_max_speed_davis(u_ll, u_rr, orientation, equations)368end369370@inline function min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,371equations::LinearizedEulerEquations2D)372return min_max_speed_davis(u_ll, u_rr, normal_direction, equations)373end374375# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes376@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,377equations::LinearizedEulerEquations2D)378@unpack v_mean_global, c_mean_global = equations379380λ_min = v_mean_global[orientation] - c_mean_global381λ_max = v_mean_global[orientation] + c_mean_global382383return λ_min, λ_max384end385386@inline function min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector,387equations::LinearizedEulerEquations2D)388@unpack v_mean_global, c_mean_global = equations389390norm_ = norm(normal_direction)391392v_normal = v_mean_global[1] * normal_direction[1] +393v_mean_global[2] * normal_direction[2]394395# The v_normals are already scaled by the norm396λ_min = v_normal - c_mean_global * norm_397λ_max = v_normal + c_mean_global * norm_398399return λ_min, λ_max400end401402# Convert conservative variables to primitive403@inline cons2prim(u, equations::LinearizedEulerEquations2D) = u404@inline cons2entropy(u, ::LinearizedEulerEquations2D) = u405end # muladd406407408