Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/equations/linearized_euler_2d.jl
5586 views
1
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).
2
# Since these FMAs can increase the performance of many numerical algorithms,
3
# we need to opt-in explicitly.
4
# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.
5
@muladd begin
6
#! format: noindent
7
8
@doc raw"""
9
LinearizedEulerEquations2D(v_mean_global, c_mean_global, rho_mean_global)
10
11
Linearized Euler equations in two space dimensions. The equations are given by
12
```math
13
\partial_t
14
\begin{pmatrix}
15
\rho' \\ v_1' \\ v_2' \\ p'
16
\end{pmatrix}
17
+
18
\partial_x
19
\begin{pmatrix}
20
\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'
21
\end{pmatrix}
22
+
23
\partial_y
24
\begin{pmatrix}
25
\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'
26
\end{pmatrix}
27
=
28
\begin{pmatrix}
29
0 \\ 0 \\ 0 \\ 0
30
\end{pmatrix}
31
```
32
The bar ``\bar{(\cdot)}`` indicates uniform mean flow variables and ``c`` is the speed of sound.
33
The unknowns are the acoustic velocities ``v' = (v_1', v_2')``, the pressure ``p'`` and the density ``\rho'``.
34
"""
35
struct LinearizedEulerEquations2D{RealT <: Real} <:
36
AbstractLinearizedEulerEquations{2, 4}
37
v_mean_global::SVector{2, RealT}
38
c_mean_global::RealT
39
rho_mean_global::RealT
40
end
41
42
function LinearizedEulerEquations2D(v_mean_global::NTuple{2, <:Real},
43
c_mean_global::Real, rho_mean_global::Real)
44
if rho_mean_global < 0
45
throw(ArgumentError("rho_mean_global must be non-negative"))
46
elseif c_mean_global < 0
47
throw(ArgumentError("c_mean_global must be non-negative"))
48
end
49
50
return LinearizedEulerEquations2D(SVector(v_mean_global), c_mean_global,
51
rho_mean_global)
52
end
53
54
function LinearizedEulerEquations2D(; v_mean_global::NTuple{2, <:Real},
55
c_mean_global::Real, rho_mean_global::Real)
56
return LinearizedEulerEquations2D(v_mean_global, c_mean_global,
57
rho_mean_global)
58
end
59
60
function varnames(::typeof(cons2cons), ::LinearizedEulerEquations2D)
61
return ("rho_prime", "v1_prime", "v2_prime", "p_prime")
62
end
63
function varnames(::typeof(cons2prim), ::LinearizedEulerEquations2D)
64
return ("rho_prime", "v1_prime", "v2_prime", "p_prime")
65
end
66
67
"""
68
initial_condition_convergence_test(x, t, equations::LinearizedEulerEquations2D)
69
70
A smooth initial condition used for convergence tests.
71
"""
72
function initial_condition_convergence_test(x, t, equations::LinearizedEulerEquations2D)
73
rho_prime = -cospi(2 * t) * (sinpi(2 * x[1]) + sinpi(2 * x[2]))
74
v1_prime = sinpi(2 * t) * cospi(2 * x[1])
75
v2_prime = sinpi(2 * t) * cospi(2 * x[2])
76
p_prime = rho_prime
77
78
return SVector(rho_prime, v1_prime, v2_prime, p_prime)
79
end
80
81
"""
82
boundary_condition_wall(u_inner, orientation, direction, x, t, surface_flux_function,
83
equations::LinearizedEulerEquations2D)
84
85
Boundary conditions for a solid wall.
86
"""
87
function boundary_condition_wall(u_inner, orientation, direction, x, t,
88
surface_flux_function,
89
equations::LinearizedEulerEquations2D)
90
# Boundary state is equal to the inner state except for the velocity. For boundaries
91
# in the -x/+x direction, we multiply the velocity in the x direction by -1.
92
# Similarly, for boundaries in the -y/+y direction, we multiply the velocity in the
93
# y direction by -1
94
if direction in (1, 2) # x direction
95
u_boundary = SVector(u_inner[1], -u_inner[2], u_inner[3], u_inner[4])
96
else # y direction
97
u_boundary = SVector(u_inner[1], u_inner[2], -u_inner[3], u_inner[4])
98
end
99
100
# Calculate boundary flux
101
if iseven(direction) # u_inner is "left" of boundary, u_boundary is "right" of boundary
102
flux = surface_flux_function(u_inner, u_boundary, orientation, equations)
103
else # u_boundary is "left" of boundary, u_inner is "right" of boundary
104
flux = surface_flux_function(u_boundary, u_inner, orientation, equations)
105
end
106
107
return flux
108
end
109
110
# Calculate 1D flux for a single point
111
@inline function flux(u, orientation::Integer, equations::LinearizedEulerEquations2D)
112
@unpack v_mean_global, c_mean_global, rho_mean_global = equations
113
rho_prime, v1_prime, v2_prime, p_prime = u
114
if orientation == 1
115
f1 = v_mean_global[1] * rho_prime + rho_mean_global * v1_prime
116
f2 = v_mean_global[1] * v1_prime + p_prime / rho_mean_global
117
f3 = v_mean_global[1] * v2_prime
118
f4 = v_mean_global[1] * p_prime + c_mean_global^2 * rho_mean_global * v1_prime
119
else
120
f1 = v_mean_global[2] * rho_prime + rho_mean_global * v2_prime
121
f2 = v_mean_global[2] * v1_prime
122
f3 = v_mean_global[2] * v2_prime + p_prime / rho_mean_global
123
f4 = v_mean_global[2] * p_prime + c_mean_global^2 * rho_mean_global * v2_prime
124
end
125
126
return SVector(f1, f2, f3, f4)
127
end
128
129
# Calculate 1D flux for a single point
130
@inline function flux(u, normal_direction::AbstractVector,
131
equations::LinearizedEulerEquations2D)
132
@unpack v_mean_global, c_mean_global, rho_mean_global = equations
133
rho_prime, v1_prime, v2_prime, p_prime = u
134
135
v_mean_normal = v_mean_global[1] * normal_direction[1] +
136
v_mean_global[2] * normal_direction[2]
137
v_prime_normal = v1_prime * normal_direction[1] + v2_prime * normal_direction[2]
138
139
f1 = v_mean_normal * rho_prime + rho_mean_global * v_prime_normal
140
f2 = v_mean_normal * v1_prime + normal_direction[1] * p_prime / rho_mean_global
141
f3 = v_mean_normal * v2_prime + normal_direction[2] * p_prime / rho_mean_global
142
f4 = v_mean_normal * p_prime + c_mean_global^2 * rho_mean_global * v_prime_normal
143
144
return SVector(f1, f2, f3, f4)
145
end
146
147
"""
148
have_constant_speed(::LinearizedEulerEquations2D)
149
150
Indicates whether the characteristic speeds are constant, i.e., independent of the solution.
151
Queried in the timestep computation [`StepsizeCallback`](@ref) and [`linear_structure`](@ref).
152
153
# Returns
154
- `True()`
155
"""
156
@inline have_constant_speed(::LinearizedEulerEquations2D) = True()
157
158
@inline function max_abs_speeds(equations::LinearizedEulerEquations2D)
159
@unpack v_mean_global, c_mean_global = equations
160
return abs(v_mean_global[1]) + c_mean_global, abs(v_mean_global[2]) + c_mean_global
161
end
162
163
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
164
equations::LinearizedEulerEquations2D)
165
@unpack v_mean_global, c_mean_global = equations
166
if orientation == 1
167
return abs(v_mean_global[1]) + c_mean_global
168
else # orientation == 2
169
return abs(v_mean_global[2]) + c_mean_global
170
end
171
end
172
173
@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
174
equations::LinearizedEulerEquations2D)
175
@unpack v_mean_global, c_mean_global = equations
176
v_mean_normal = normal_direction[1] * v_mean_global[1] +
177
normal_direction[2] * v_mean_global[2]
178
return abs(v_mean_normal) + c_mean_global * norm(normal_direction)
179
end
180
181
@doc raw"""
182
flux_godunov(u_ll, u_rr, orientation_or_normal_direction,
183
equations::LinearizedEulerEquations2D)
184
185
An upwind flux for the linearized Euler equations based on diagonalization of the physical
186
flux matrix. Given the physical flux ``Au``, ``A=T \Lambda T^{-1}`` with
187
``\Lambda`` being a diagonal matrix that holds the eigenvalues of ``A``, decompose
188
``\Lambda = \Lambda^+ + \Lambda^-`` where ``\Lambda^+`` and ``\Lambda^-`` are diagonal
189
matrices holding the positive and negative eigenvalues of ``A``, respectively. Then for
190
left and right states ``u_L, u_R``, the numerical flux calculated by this function is given
191
by ``A^+ u_L + A^- u_R`` where ``A^{\pm} = T \Lambda^{\pm} T^{-1}``.
192
193
The diagonalization of the flux matrix can be found in
194
- R. F. Warming, Richard M. Beam and B. J. Hyett (1975)
195
Diagonalization and simultaneous symmetrization of the gas-dynamic matrices
196
[DOI: 10.1090/S0025-5718-1975-0388967-5](https://doi.org/10.1090/S0025-5718-1975-0388967-5)
197
"""
198
@inline function flux_godunov(u_ll, u_rr, orientation::Integer,
199
equations::LinearizedEulerEquations2D)
200
@unpack v_mean_global, rho_mean_global, c_mean_global = equations
201
v1_mean = v_mean_global[1]
202
v2_mean = v_mean_global[2]
203
204
rho_prime_ll, v1_prime_ll, v2_prime_ll, p_prime_ll = u_ll
205
rho_prime_rr, v1_prime_rr, v2_prime_rr, p_prime_rr = u_rr
206
207
if orientation == 1
208
# Eigenvalues of the flux matrix
209
lambda1 = v1_mean
210
lambda2 = v1_mean - c_mean_global
211
lambda3 = v1_mean + c_mean_global
212
213
lambda1_p = positive_part(lambda1)
214
lambda2_p = positive_part(lambda2)
215
lambda3_p = positive_part(lambda3)
216
lambda2p3_half_p = 0.5f0 * (lambda2_p + lambda3_p)
217
lambda3m2_half_p = 0.5f0 * (lambda3_p - lambda2_p)
218
219
lambda1_m = negative_part(lambda1)
220
lambda2_m = negative_part(lambda2)
221
lambda3_m = negative_part(lambda3)
222
lambda2p3_half_m = 0.5f0 * (lambda2_m + lambda3_m)
223
lambda3m2_half_m = 0.5f0 * (lambda3_m - lambda2_m)
224
225
f1p = (lambda1_p * rho_prime_ll +
226
lambda3m2_half_p / c_mean_global * rho_mean_global * v1_prime_ll +
227
(lambda2p3_half_p - lambda1_p) / c_mean_global^2 * p_prime_ll)
228
f2p = (lambda2p3_half_p * v1_prime_ll +
229
lambda3m2_half_p / c_mean_global * p_prime_ll / rho_mean_global)
230
f3p = lambda1_p * v2_prime_ll
231
f4p = (lambda3m2_half_p * c_mean_global * rho_mean_global * v1_prime_ll +
232
lambda2p3_half_p * p_prime_ll)
233
234
f1m = (lambda1_m * rho_prime_rr +
235
lambda3m2_half_m / c_mean_global * rho_mean_global * v1_prime_rr +
236
(lambda2p3_half_m - lambda1_m) / c_mean_global^2 * p_prime_rr)
237
f2m = (lambda2p3_half_m * v1_prime_rr +
238
lambda3m2_half_m / c_mean_global * p_prime_rr / rho_mean_global)
239
f3m = lambda1_m * v2_prime_rr
240
f4m = (lambda3m2_half_m * c_mean_global * rho_mean_global * v1_prime_rr +
241
lambda2p3_half_m * p_prime_rr)
242
243
f1 = f1p + f1m
244
f2 = f2p + f2m
245
f3 = f3p + f3m
246
f4 = f4p + f4m
247
else # orientation == 2
248
# Eigenvalues of the flux matrix
249
lambda1 = v2_mean
250
lambda2 = v2_mean - c_mean_global
251
lambda3 = v2_mean + c_mean_global
252
253
lambda1_p = positive_part(lambda1)
254
lambda2_p = positive_part(lambda2)
255
lambda3_p = positive_part(lambda3)
256
lambda2p3_half_p = 0.5f0 * (lambda2_p + lambda3_p)
257
lambda3m2_half_p = 0.5f0 * (lambda3_p - lambda2_p)
258
259
lambda1_m = negative_part(lambda1)
260
lambda2_m = negative_part(lambda2)
261
lambda3_m = negative_part(lambda3)
262
lambda2p3_half_m = 0.5f0 * (lambda2_m + lambda3_m)
263
lambda3m2_half_m = 0.5f0 * (lambda3_m - lambda2_m)
264
265
f1p = (lambda1_p * rho_prime_ll +
266
lambda3m2_half_p / c_mean_global * rho_mean_global * v2_prime_ll +
267
(lambda2p3_half_p - lambda1_p) / c_mean_global^2 * p_prime_ll)
268
f2p = lambda1_p * v1_prime_ll
269
f3p = (lambda2p3_half_p * v2_prime_ll +
270
lambda3m2_half_p / c_mean_global * p_prime_ll / rho_mean_global)
271
f4p = (lambda3m2_half_p * c_mean_global * rho_mean_global * v2_prime_ll +
272
lambda2p3_half_p * p_prime_ll)
273
274
f1m = (lambda1_m * rho_prime_rr +
275
lambda3m2_half_m / c_mean_global * rho_mean_global * v2_prime_rr +
276
(lambda2p3_half_m - lambda1_m) / c_mean_global^2 * p_prime_rr)
277
f2m = lambda1_m * v1_prime_rr
278
f3m = (lambda2p3_half_m * v2_prime_rr +
279
lambda3m2_half_m / c_mean_global * p_prime_rr / rho_mean_global)
280
f4m = (lambda3m2_half_m * c_mean_global * rho_mean_global * v2_prime_rr +
281
lambda2p3_half_m * p_prime_rr)
282
283
f1 = f1p + f1m
284
f2 = f2p + f2m
285
f3 = f3p + f3m
286
f4 = f4p + f4m
287
end
288
289
return SVector(f1, f2, f3, f4)
290
end
291
292
@inline function flux_godunov(u_ll, u_rr, normal_direction::AbstractVector,
293
equations::LinearizedEulerEquations2D)
294
@unpack v_mean_global, rho_mean_global, c_mean_global = equations
295
rho_prime_ll, v1_prime_ll, v2_prime_ll, p_prime_ll = u_ll
296
rho_prime_rr, v1_prime_rr, v2_prime_rr, p_prime_rr = u_rr
297
298
# Do not use `normalize` since we use `norm_` later to scale the eigenvalues
299
norm_ = norm(normal_direction)
300
normal_vector = normal_direction / norm_
301
302
# Use normalized vector here, scaling is applied via eigenvalues of the flux matrix
303
v_mean_normal = v_mean_global[1] * normal_vector[1] +
304
v_mean_global[2] * normal_vector[2]
305
v_prime_normal_ll = v1_prime_ll * normal_vector[1] + v2_prime_ll * normal_vector[2]
306
v_prime_normal_rr = v1_prime_rr * normal_vector[1] + v2_prime_rr * normal_vector[2]
307
308
# Eigenvalues of the flux matrix
309
lambda1 = v_mean_normal * norm_
310
lambda2 = (v_mean_normal - c_mean_global) * norm_
311
lambda3 = (v_mean_normal + c_mean_global) * norm_
312
313
lambda1_p = positive_part(lambda1)
314
lambda2_p = positive_part(lambda2)
315
lambda3_p = positive_part(lambda3)
316
lambda2p3_half_p = 0.5f0 * (lambda2_p + lambda3_p)
317
lambda3m2_half_p = 0.5f0 * (lambda3_p - lambda2_p)
318
319
lambda1_m = negative_part(lambda1)
320
lambda2_m = negative_part(lambda2)
321
lambda3_m = negative_part(lambda3)
322
lambda2p3_half_m = 0.5f0 * (lambda2_m + lambda3_m)
323
lambda3m2_half_m = 0.5f0 * (lambda3_m - lambda2_m)
324
325
f1p = (lambda1_p * rho_prime_ll +
326
lambda3m2_half_p / c_mean_global * rho_mean_global * v_prime_normal_ll +
327
(lambda2p3_half_p - lambda1_p) / c_mean_global^2 * p_prime_ll)
328
f2p = (((lambda1_p * normal_vector[2]^2 +
329
lambda2p3_half_p * normal_vector[1]^2) * v1_prime_ll +
330
(lambda2p3_half_p - lambda1_p) * prod(normal_vector) * v2_prime_ll) +
331
lambda3m2_half_p / c_mean_global * normal_vector[1] * p_prime_ll /
332
rho_mean_global)
333
f3p = (((lambda1_p * normal_vector[1]^2 +
334
lambda2p3_half_p * normal_vector[2]^2) * v2_prime_ll +
335
(lambda2p3_half_p - lambda1_p) * prod(normal_vector) * v1_prime_ll) +
336
lambda3m2_half_p / c_mean_global * normal_vector[2] * p_prime_ll /
337
rho_mean_global)
338
f4p = (lambda3m2_half_p * c_mean_global * rho_mean_global * v_prime_normal_ll +
339
lambda2p3_half_p * p_prime_ll)
340
341
f1m = (lambda1_m * rho_prime_rr +
342
lambda3m2_half_m / c_mean_global * rho_mean_global * v_prime_normal_rr +
343
(lambda2p3_half_m - lambda1_m) / c_mean_global^2 * p_prime_rr)
344
f2m = (((lambda1_m * normal_vector[2]^2 +
345
lambda2p3_half_m * normal_vector[1]^2) * v1_prime_rr +
346
(lambda2p3_half_m - lambda1_m) * prod(normal_vector) * v2_prime_rr) +
347
lambda3m2_half_m / c_mean_global * normal_vector[1] * p_prime_rr /
348
rho_mean_global)
349
f3m = (((lambda1_m * normal_vector[1]^2 +
350
lambda2p3_half_m * normal_vector[2]^2) * v2_prime_rr +
351
(lambda2p3_half_m - lambda1_m) * prod(normal_vector) * v1_prime_rr) +
352
lambda3m2_half_m / c_mean_global * normal_vector[2] * p_prime_rr /
353
rho_mean_global)
354
f4m = (lambda3m2_half_m * c_mean_global * rho_mean_global * v_prime_normal_rr +
355
lambda2p3_half_m * p_prime_rr)
356
357
f1 = f1p + f1m
358
f2 = f2p + f2m
359
f3 = f3p + f3m
360
f4 = f4p + f4m
361
362
return SVector(f1, f2, f3, f4)
363
end
364
365
# Calculate estimate for minimum and maximum wave speeds for HLL-type fluxes
366
@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,
367
equations::LinearizedEulerEquations2D)
368
return min_max_speed_davis(u_ll, u_rr, orientation, equations)
369
end
370
371
@inline function min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
372
equations::LinearizedEulerEquations2D)
373
return min_max_speed_davis(u_ll, u_rr, normal_direction, equations)
374
end
375
376
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
377
@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,
378
equations::LinearizedEulerEquations2D)
379
@unpack v_mean_global, c_mean_global = equations
380
381
λ_min = v_mean_global[orientation] - c_mean_global
382
λ_max = v_mean_global[orientation] + c_mean_global
383
384
return λ_min, λ_max
385
end
386
387
@inline function min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector,
388
equations::LinearizedEulerEquations2D)
389
@unpack v_mean_global, c_mean_global = equations
390
391
norm_ = norm(normal_direction)
392
393
v_normal = v_mean_global[1] * normal_direction[1] +
394
v_mean_global[2] * normal_direction[2]
395
396
# The v_normals are already scaled by the norm
397
λ_min = v_normal - c_mean_global * norm_
398
λ_max = v_normal + c_mean_global * norm_
399
400
return λ_min, λ_max
401
end
402
403
# Convert conservative variables to primitive
404
@inline cons2prim(u, equations::LinearizedEulerEquations2D) = u
405
@inline cons2entropy(u, ::LinearizedEulerEquations2D) = u
406
end # muladd
407
408