Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/equations/compressible_euler_3d.jl
2055 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
CompressibleEulerEquations3D(gamma)
10
11
The compressible Euler equations
12
```math
13
\frac{\partial}{\partial t}
14
\begin{pmatrix}
15
\rho \\ \rho v_1 \\ \rho v_2 \\ \rho v_3 \\ \rho e
16
\end{pmatrix}
17
+
18
\frac{\partial}{\partial x}
19
\begin{pmatrix}
20
\rho v_1 \\ \rho v_1^2 + p \\ \rho v_1 v_2 \\ \rho v_1 v_3 \\ ( \rho e +p) v_1
21
\end{pmatrix}
22
+
23
\frac{\partial}{\partial y}
24
\begin{pmatrix}
25
\rho v_2 \\ \rho v_1 v_2 \\ \rho v_2^2 + p \\ \rho v_1 v_3 \\ ( \rho e +p) v_2
26
\end{pmatrix}
27
+
28
\frac{\partial}{\partial z}
29
\begin{pmatrix}
30
\rho v_3 \\ \rho v_1 v_3 \\ \rho v_2 v_3 \\ \rho v_3^2 + p \\ ( \rho e +p) v_3
31
\end{pmatrix}
32
=
33
\begin{pmatrix}
34
0 \\ 0 \\ 0 \\ 0 \\ 0
35
\end{pmatrix}
36
```
37
for an ideal gas with ratio of specific heats `gamma`
38
in three space dimensions.
39
Here, ``\rho`` is the density, ``v_1``, ``v_2``, ``v_3`` the velocities, ``e`` the specific total energy **rather than** specific internal energy, and
40
```math
41
p = (\gamma - 1) \left( \rho e - \frac{1}{2} \rho (v_1^2+v_2^2+v_3^2) \right)
42
```
43
the pressure.
44
"""
45
struct CompressibleEulerEquations3D{RealT <: Real} <:
46
AbstractCompressibleEulerEquations{3, 5}
47
gamma::RealT # ratio of specific heats
48
inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications
49
50
function CompressibleEulerEquations3D(gamma)
51
γ, inv_gamma_minus_one = promote(gamma, inv(gamma - 1))
52
new{typeof(γ)}(γ, inv_gamma_minus_one)
53
end
54
end
55
56
function varnames(::typeof(cons2cons), ::CompressibleEulerEquations3D)
57
("rho", "rho_v1", "rho_v2", "rho_v3", "rho_e")
58
end
59
function varnames(::typeof(cons2prim), ::CompressibleEulerEquations3D)
60
("rho", "v1", "v2", "v3", "p")
61
end
62
63
# Set initial conditions at physical location `x` for time `t`
64
"""
65
initial_condition_constant(x, t, equations::CompressibleEulerEquations3D)
66
67
A constant initial condition to test free-stream preservation.
68
"""
69
function initial_condition_constant(x, t, equations::CompressibleEulerEquations3D)
70
RealT = eltype(x)
71
rho = 1
72
rho_v1 = convert(RealT, 0.1)
73
rho_v2 = convert(RealT, -0.2)
74
rho_v3 = convert(RealT, 0.7)
75
rho_e = 10
76
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e)
77
end
78
79
"""
80
initial_condition_convergence_test(x, t, equations::CompressibleEulerEquations3D)
81
82
A smooth initial condition used for convergence tests in combination with
83
[`source_terms_convergence_test`](@ref).
84
"""
85
function initial_condition_convergence_test(x, t,
86
equations::CompressibleEulerEquations3D)
87
RealT = eltype(x)
88
c = 2
89
A = convert(RealT, 0.1)
90
L = 2
91
f = 1.0f0 / L
92
ω = 2 * convert(RealT, pi) * f
93
ini = c + A * sin(ω * (x[1] + x[2] + x[3] - t))
94
95
rho = ini
96
rho_v1 = ini
97
rho_v2 = ini
98
rho_v3 = ini
99
rho_e = ini^2
100
101
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e)
102
end
103
104
"""
105
source_terms_convergence_test(u, x, t, equations::CompressibleEulerEquations3D)
106
107
Source terms used for convergence tests in combination with
108
[`initial_condition_convergence_test`](@ref).
109
"""
110
@inline function source_terms_convergence_test(u, x, t,
111
equations::CompressibleEulerEquations3D)
112
# Same settings as in `initial_condition`
113
RealT = eltype(u)
114
c = 2
115
A = convert(RealT, 0.1)
116
L = 2
117
f = 1.0f0 / L
118
ω = 2 * convert(RealT, pi) * f
119
γ = equations.gamma
120
121
x1, x2, x3 = x
122
si, co = sincos(ω * (x1 + x2 + x3 - t))
123
rho = c + A * si
124
rho_x = ω * A * co
125
# Note that d/dt rho = -d/dx rho = -d/dy rho = - d/dz rho.
126
127
tmp = (2 * rho - 1.5f0) * (γ - 1)
128
129
du1 = 2 * rho_x
130
du2 = rho_x * (2 + tmp)
131
du3 = du2
132
du4 = du2
133
du5 = rho_x * (4 * rho + 3 * tmp)
134
135
return SVector(du1, du2, du3, du4, du5)
136
end
137
138
"""
139
initial_condition_weak_blast_wave(x, t, equations::CompressibleEulerEquations3D)
140
141
A weak blast wave taken from
142
- Sebastian Hennemann, Gregor J. Gassner (2020)
143
A provably entropy stable subcell shock capturing approach for high order split form DG
144
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
145
"""
146
function initial_condition_weak_blast_wave(x, t,
147
equations::CompressibleEulerEquations3D)
148
# From Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
149
# Set up spherical coordinates
150
RealT = eltype(x)
151
inicenter = (0, 0, 0)
152
x_norm = x[1] - inicenter[1]
153
y_norm = x[2] - inicenter[2]
154
z_norm = x[3] - inicenter[3]
155
r = sqrt(x_norm^2 + y_norm^2 + z_norm^2)
156
phi = atan(y_norm, x_norm)
157
theta = iszero(r) ? zero(RealT) : acos(z_norm / r)
158
159
# Calculate primitive variables
160
rho = r > 0.5f0 ? one(RealT) : convert(RealT, 1.1691)
161
v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos(phi) * sin(theta)
162
v2 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * sin(phi) * sin(theta)
163
v3 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos(theta)
164
p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)
165
166
return prim2cons(SVector(rho, v1, v2, v3, p), equations)
167
end
168
169
"""
170
initial_condition_eoc_test_coupled_euler_gravity(x, t, equations::CompressibleEulerEquations3D)
171
172
Setup used for convergence tests of the Euler equations with self-gravity used in
173
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
174
A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics
175
[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)
176
in combination with [`source_terms_eoc_test_coupled_euler_gravity`](@ref)
177
or [`source_terms_eoc_test_euler`](@ref).
178
"""
179
function initial_condition_eoc_test_coupled_euler_gravity(x, t,
180
equations::CompressibleEulerEquations3D)
181
# OBS! this assumes that γ = 2 other manufactured source terms are incorrect
182
if equations.gamma != 2
183
error("adiabatic constant must be 2 for the coupling convergence test")
184
end
185
RealT = eltype(x)
186
c = 2
187
A = convert(RealT, 0.1)
188
ini = c + A * sinpi(x[1] + x[2] + x[3] - t)
189
G = 1 # gravitational constant
190
191
rho = ini
192
v1 = 1
193
v2 = 1
194
v3 = 1
195
p = ini^2 * G * 2 / (3 * convert(RealT, pi)) # "3" is the number of spatial dimensions
196
197
return prim2cons(SVector(rho, v1, v2, v3, p), equations)
198
end
199
200
"""
201
source_terms_eoc_test_coupled_euler_gravity(u, x, t, equations::CompressibleEulerEquations3D)
202
203
Setup used for convergence tests of the Euler equations with self-gravity used in
204
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
205
A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics
206
[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)
207
in combination with [`initial_condition_eoc_test_coupled_euler_gravity`](@ref).
208
"""
209
@inline function source_terms_eoc_test_coupled_euler_gravity(u, x, t,
210
equations::CompressibleEulerEquations3D)
211
# Same settings as in `initial_condition_eoc_test_coupled_euler_gravity`
212
RealT = eltype(u)
213
c = 2
214
A = convert(RealT, 0.1)
215
G = 1 # gravitational constant, must match coupling solver
216
C_grav = -4 * G / (3 * convert(RealT, pi)) # "3" is the number of spatial dimensions # 2D: -2.0*G/pi
217
218
x1, x2, x3 = x
219
si, co = sincospi(x1 + x2 + x3 - t)
220
rhox = A * convert(RealT, pi) * co
221
rho = c + A * si
222
223
# In "2 * rhox", the "2" is "number of spatial dimensions minus one"
224
du1 = 2 * rhox
225
du2 = 2 * rhox
226
du3 = 2 * rhox
227
du4 = 2 * rhox
228
du5 = 2 * rhox * (1.5f0 - C_grav * rho) # "3" in "3/2" is the number of spatial dimensions
229
230
return SVector(du1, du2, du3, du4, du5)
231
end
232
233
"""
234
source_terms_eoc_test_euler(u, x, t, equations::CompressibleEulerEquations3D)
235
236
Setup used for convergence tests of the Euler equations with self-gravity used in
237
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
238
A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics
239
[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)
240
in combination with [`initial_condition_eoc_test_coupled_euler_gravity`](@ref).
241
242
!!! note
243
This method is to be used for testing pure Euler simulations with analytic self-gravity.
244
If you intend to do coupled Euler-gravity simulations, you need to use
245
[`source_terms_eoc_test_coupled_euler_gravity`](@ref) instead.
246
"""
247
function source_terms_eoc_test_euler(u, x, t, equations::CompressibleEulerEquations3D)
248
# Same settings as in `initial_condition_eoc_test_coupled_euler_gravity`
249
RealT = eltype(u)
250
c = 2
251
A = convert(RealT, 0.1)
252
G = 1
253
C_grav = -4 * G / (3 * convert(RealT, pi)) # "3" is the number of spatial dimensions
254
255
x1, x2, x3 = x
256
si, co = sincospi(x1 + x2 + x3 - t)
257
rhox = A * convert(RealT, pi) * co
258
rho = c + A * si
259
260
du1 = rhox * 2
261
du2 = rhox * (2 - C_grav * rho)
262
du3 = rhox * (2 - C_grav * rho)
263
du4 = rhox * (2 - C_grav * rho)
264
du5 = rhox * (3 - 5 * C_grav * rho)
265
266
return SVector(du1, du2, du3, du4, du5)
267
end
268
269
"""
270
boundary_condition_slip_wall(u_inner, normal_direction, x, t, surface_flux_function,
271
equations::CompressibleEulerEquations3D)
272
273
Determine the boundary numerical surface flux for a slip wall condition.
274
Imposes a zero normal velocity at the wall.
275
Density is taken from the internal solution state and pressure is computed as an
276
exact solution of a 1D Riemann problem. Further details about this boundary state
277
are available in the paper:
278
- J. J. W. van der Vegt and H. van der Ven (2002)
279
Slip flow boundary conditions in discontinuous Galerkin discretizations of
280
the Euler equations of gas dynamics
281
[PDF](https://reports.nlr.nl/bitstream/handle/10921/692/TP-2002-300.pdf?sequence=1)
282
283
Details about the 1D pressure Riemann solution can be found in Section 6.3.3 of the book
284
- Eleuterio F. Toro (2009)
285
Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction
286
3rd edition
287
[DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)
288
289
Should be used together with [`P4estMesh`](@ref) or [`T8codeMesh`](@ref).
290
"""
291
@inline function boundary_condition_slip_wall(u_inner, normal_direction::AbstractVector,
292
x, t,
293
surface_flux_function,
294
equations::CompressibleEulerEquations3D)
295
norm_ = norm(normal_direction)
296
# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later
297
normal = normal_direction / norm_
298
299
# Some vector that can't be identical to normal_vector (unless normal_vector == 0)
300
tangent1 = SVector(normal_direction[2], normal_direction[3], -normal_direction[1])
301
# Orthogonal projection
302
tangent1 -= dot(normal, tangent1) * normal
303
tangent1 = normalize(tangent1)
304
305
# Third orthogonal vector
306
tangent2 = normalize(cross(normal_direction, tangent1))
307
308
# rotate the internal solution state
309
u_local = rotate_to_x(u_inner, normal, tangent1, tangent2, equations)
310
311
# compute the primitive variables
312
rho_local, v_normal, v_tangent1, v_tangent2, p_local = cons2prim(u_local, equations)
313
314
# Get the solution of the pressure Riemann problem
315
# See Section 6.3.3 of
316
# Eleuterio F. Toro (2009)
317
# Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction
318
# [DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)
319
if v_normal <= 0
320
sound_speed = sqrt(equations.gamma * p_local / rho_local) # local sound speed
321
p_star = p_local *
322
(1 + 0.5f0 * (equations.gamma - 1) * v_normal / sound_speed)^(2 *
323
equations.gamma *
324
equations.inv_gamma_minus_one)
325
else # v_normal > 0
326
A = 2 / ((equations.gamma + 1) * rho_local)
327
B = p_local * (equations.gamma - 1) / (equations.gamma + 1)
328
p_star = p_local +
329
0.5f0 * v_normal / A *
330
(v_normal + sqrt(v_normal^2 + 4 * A * (p_local + B)))
331
end
332
333
# For the slip wall we directly set the flux as the normal velocity is zero
334
return SVector(0,
335
p_star * normal[1],
336
p_star * normal[2],
337
p_star * normal[3],
338
0) * norm_
339
end
340
341
"""
342
boundary_condition_slip_wall(u_inner, orientation, direction, x, t,
343
surface_flux_function, equations::CompressibleEulerEquations3D)
344
345
Should be used together with [`TreeMesh`](@ref).
346
"""
347
@inline function boundary_condition_slip_wall(u_inner, orientation,
348
direction, x, t,
349
surface_flux_function,
350
equations::CompressibleEulerEquations3D)
351
# get the appropriate normal vector from the orientation
352
RealT = eltype(u_inner)
353
if orientation == 1
354
normal_direction = SVector(one(RealT), zero(RealT), zero(RealT))
355
elseif orientation == 2
356
normal_direction = SVector(zero(RealT), one(RealT), zero(RealT))
357
else # orientation == 3
358
normal_direction = SVector(zero(RealT), zero(RealT), one(RealT))
359
end
360
361
# compute and return the flux using `boundary_condition_slip_wall` routine below
362
return boundary_condition_slip_wall(u_inner, normal_direction, direction,
363
x, t, surface_flux_function, equations)
364
end
365
366
"""
367
boundary_condition_slip_wall(u_inner, normal_direction, direction, x, t,
368
surface_flux_function, equations::CompressibleEulerEquations3D)
369
370
Should be used together with [`StructuredMesh`](@ref).
371
"""
372
@inline function boundary_condition_slip_wall(u_inner, normal_direction::AbstractVector,
373
direction, x, t,
374
surface_flux_function,
375
equations::CompressibleEulerEquations3D)
376
# flip sign of normal to make it outward pointing, then flip the sign of the normal flux back
377
# to be inward pointing on the -x, -y, and -z sides due to the orientation convention used by StructuredMesh
378
if isodd(direction)
379
boundary_flux = -boundary_condition_slip_wall(u_inner, -normal_direction,
380
x, t, surface_flux_function,
381
equations)
382
else
383
boundary_flux = boundary_condition_slip_wall(u_inner, normal_direction,
384
x, t, surface_flux_function,
385
equations)
386
end
387
388
return boundary_flux
389
end
390
391
# Calculate 1D flux for a single point
392
@inline function flux(u, orientation::Integer, equations::CompressibleEulerEquations3D)
393
rho, rho_v1, rho_v2, rho_v3, rho_e = u
394
v1 = rho_v1 / rho
395
v2 = rho_v2 / rho
396
v3 = rho_v3 / rho
397
p = (equations.gamma - 1) *
398
(rho_e - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3))
399
if orientation == 1
400
f1 = rho_v1
401
f2 = rho_v1 * v1 + p
402
f3 = rho_v1 * v2
403
f4 = rho_v1 * v3
404
f5 = (rho_e + p) * v1
405
elseif orientation == 2
406
f1 = rho_v2
407
f2 = rho_v2 * v1
408
f3 = rho_v2 * v2 + p
409
f4 = rho_v2 * v3
410
f5 = (rho_e + p) * v2
411
else
412
f1 = rho_v3
413
f2 = rho_v3 * v1
414
f3 = rho_v3 * v2
415
f4 = rho_v3 * v3 + p
416
f5 = (rho_e + p) * v3
417
end
418
return SVector(f1, f2, f3, f4, f5)
419
end
420
421
@inline function flux(u, normal_direction::AbstractVector,
422
equations::CompressibleEulerEquations3D)
423
rho_e = last(u)
424
rho, v1, v2, v3, p = cons2prim(u, equations)
425
426
v_normal = v1 * normal_direction[1] + v2 * normal_direction[2] +
427
v3 * normal_direction[3]
428
rho_v_normal = rho * v_normal
429
f1 = rho_v_normal
430
f2 = rho_v_normal * v1 + p * normal_direction[1]
431
f3 = rho_v_normal * v2 + p * normal_direction[2]
432
f4 = rho_v_normal * v3 + p * normal_direction[3]
433
f5 = (rho_e + p) * v_normal
434
return SVector(f1, f2, f3, f4, f5)
435
end
436
437
"""
438
flux_shima_etal(u_ll, u_rr, orientation_or_normal_direction,
439
equations::CompressibleEulerEquations3D)
440
441
This flux is is a modification of the original kinetic energy preserving two-point flux by
442
- Yuichi Kuya, Kosuke Totani and Soshi Kawai (2018)
443
Kinetic energy and entropy preserving schemes for compressible flows
444
by split convective forms
445
[DOI: 10.1016/j.jcp.2018.08.058](https://doi.org/10.1016/j.jcp.2018.08.058)
446
447
The modification is in the energy flux to guarantee pressure equilibrium and was developed by
448
- Nao Shima, Yuichi Kuya, Yoshiharu Tamaki, Soshi Kawai (JCP 2020)
449
Preventing spurious pressure oscillations in split convective form discretizations for
450
compressible flows
451
[DOI: 10.1016/j.jcp.2020.110060](https://doi.org/10.1016/j.jcp.2020.110060)
452
"""
453
@inline function flux_shima_etal(u_ll, u_rr, orientation::Integer,
454
equations::CompressibleEulerEquations3D)
455
# Unpack left and right state
456
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
457
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
458
459
# Average each factor of products in flux
460
rho_avg = 0.5f0 * (rho_ll + rho_rr)
461
v1_avg = 0.5f0 * (v1_ll + v1_rr)
462
v2_avg = 0.5f0 * (v2_ll + v2_rr)
463
v3_avg = 0.5f0 * (v3_ll + v3_rr)
464
p_avg = 0.5f0 * (p_ll + p_rr)
465
kin_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
466
467
# Calculate fluxes depending on orientation
468
if orientation == 1
469
pv1_avg = 0.5f0 * (p_ll * v1_rr + p_rr * v1_ll)
470
f1 = rho_avg * v1_avg
471
f2 = f1 * v1_avg + p_avg
472
f3 = f1 * v2_avg
473
f4 = f1 * v3_avg
474
f5 = p_avg * v1_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv1_avg
475
elseif orientation == 2
476
pv2_avg = 0.5f0 * (p_ll * v2_rr + p_rr * v2_ll)
477
f1 = rho_avg * v2_avg
478
f2 = f1 * v1_avg
479
f3 = f1 * v2_avg + p_avg
480
f4 = f1 * v3_avg
481
f5 = p_avg * v2_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv2_avg
482
else
483
pv3_avg = 0.5f0 * (p_ll * v3_rr + p_rr * v3_ll)
484
f1 = rho_avg * v3_avg
485
f2 = f1 * v1_avg
486
f3 = f1 * v2_avg
487
f4 = f1 * v3_avg + p_avg
488
f5 = p_avg * v3_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv3_avg
489
end
490
491
return SVector(f1, f2, f3, f4, f5)
492
end
493
494
@inline function flux_shima_etal(u_ll, u_rr, normal_direction::AbstractVector,
495
equations::CompressibleEulerEquations3D)
496
# Unpack left and right state
497
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
498
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
499
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2] +
500
v3_ll * normal_direction[3]
501
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2] +
502
v3_rr * normal_direction[3]
503
504
# Average each factor of products in flux
505
rho_avg = 0.5f0 * (rho_ll + rho_rr)
506
v1_avg = 0.5f0 * (v1_ll + v1_rr)
507
v2_avg = 0.5f0 * (v2_ll + v2_rr)
508
v3_avg = 0.5f0 * (v3_ll + v3_rr)
509
v_dot_n_avg = 0.5f0 * (v_dot_n_ll + v_dot_n_rr)
510
p_avg = 0.5f0 * (p_ll + p_rr)
511
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
512
513
# Calculate fluxes depending on normal_direction
514
f1 = rho_avg * v_dot_n_avg
515
f2 = f1 * v1_avg + p_avg * normal_direction[1]
516
f3 = f1 * v2_avg + p_avg * normal_direction[2]
517
f4 = f1 * v3_avg + p_avg * normal_direction[3]
518
f5 = (f1 * velocity_square_avg +
519
p_avg * v_dot_n_avg * equations.inv_gamma_minus_one
520
+ 0.5f0 * (p_ll * v_dot_n_rr + p_rr * v_dot_n_ll))
521
522
return SVector(f1, f2, f3, f4, f5)
523
end
524
525
"""
526
flux_kennedy_gruber(u_ll, u_rr, orientation_or_normal_direction,
527
equations::CompressibleEulerEquations3D)
528
529
Kinetic energy preserving two-point flux by
530
- Kennedy and Gruber (2008)
531
Reduced aliasing formulations of the convective terms within the
532
Navier-Stokes equations for a compressible fluid
533
[DOI: 10.1016/j.jcp.2007.09.020](https://doi.org/10.1016/j.jcp.2007.09.020)
534
"""
535
@inline function flux_kennedy_gruber(u_ll, u_rr, orientation::Integer,
536
equations::CompressibleEulerEquations3D)
537
# Unpack left and right state
538
rho_e_ll = last(u_ll)
539
rho_e_rr = last(u_rr)
540
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
541
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
542
543
# Average each factor of products in flux
544
rho_avg = 0.5f0 * (rho_ll + rho_rr)
545
v1_avg = 0.5f0 * (v1_ll + v1_rr)
546
v2_avg = 0.5f0 * (v2_ll + v2_rr)
547
v3_avg = 0.5f0 * (v3_ll + v3_rr)
548
p_avg = 0.5f0 * (p_ll + p_rr)
549
e_avg = 0.5f0 * (rho_e_ll / rho_ll + rho_e_rr / rho_rr)
550
551
# Calculate fluxes depending on orientation
552
if orientation == 1
553
f1 = rho_avg * v1_avg
554
f2 = f1 * v1_avg + p_avg
555
f3 = f1 * v2_avg
556
f4 = f1 * v3_avg
557
f5 = (rho_avg * e_avg + p_avg) * v1_avg
558
elseif orientation == 2
559
f1 = rho_avg * v2_avg
560
f2 = f1 * v1_avg
561
f3 = f1 * v2_avg + p_avg
562
f4 = f1 * v3_avg
563
f5 = (rho_avg * e_avg + p_avg) * v2_avg
564
else
565
f1 = rho_avg * v3_avg
566
f2 = f1 * v1_avg
567
f3 = f1 * v2_avg
568
f4 = f1 * v3_avg + p_avg
569
f5 = (rho_avg * e_avg + p_avg) * v3_avg
570
end
571
572
return SVector(f1, f2, f3, f4, f5)
573
end
574
575
@inline function flux_kennedy_gruber(u_ll, u_rr, normal_direction::AbstractVector,
576
equations::CompressibleEulerEquations3D)
577
# Unpack left and right state
578
rho_e_ll = last(u_ll)
579
rho_e_rr = last(u_rr)
580
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
581
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
582
583
# Average each factor of products in flux
584
rho_avg = 0.5f0 * (rho_ll + rho_rr)
585
v1_avg = 0.5f0 * (v1_ll + v1_rr)
586
v2_avg = 0.5f0 * (v2_ll + v2_rr)
587
v3_avg = 0.5f0 * (v3_ll + v3_rr)
588
p_avg = 0.5f0 * (p_ll + p_rr)
589
e_avg = 0.5f0 * (rho_e_ll / rho_ll + rho_e_rr / rho_rr)
590
591
v_dot_n_avg = v1_avg * normal_direction[1] + v2_avg * normal_direction[2] +
592
v3_avg * normal_direction[3]
593
p_avg = 0.5f0 * (p_ll + p_rr)
594
e_avg = 0.5f0 * (rho_e_ll / rho_ll + rho_e_rr / rho_rr)
595
596
# Calculate fluxes depending on normal_direction
597
f1 = rho_avg * v_dot_n_avg
598
f2 = f1 * v1_avg + p_avg * normal_direction[1]
599
f3 = f1 * v2_avg + p_avg * normal_direction[2]
600
f4 = f1 * v3_avg + p_avg * normal_direction[3]
601
f5 = f1 * e_avg + p_avg * v_dot_n_avg
602
603
return SVector(f1, f2, f3, f4, f5)
604
end
605
606
"""
607
flux_chandrashekar(u_ll, u_rr, orientation_or_normal_direction, equations::CompressibleEulerEquations3D)
608
609
Entropy conserving two-point flux by
610
- Chandrashekar (2013)
611
Kinetic Energy Preserving and Entropy Stable Finite Volume Schemes
612
for Compressible Euler and Navier-Stokes Equations
613
[DOI: 10.4208/cicp.170712.010313a](https://doi.org/10.4208/cicp.170712.010313a)
614
"""
615
@inline function flux_chandrashekar(u_ll, u_rr, orientation::Integer,
616
equations::CompressibleEulerEquations3D)
617
# Unpack left and right state
618
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
619
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
620
621
beta_ll = 0.5f0 * rho_ll / p_ll
622
beta_rr = 0.5f0 * rho_rr / p_rr
623
specific_kin_ll = 0.5f0 * (v1_ll^2 + v2_ll^2 + v3_ll^2)
624
specific_kin_rr = 0.5f0 * (v1_rr^2 + v2_rr^2 + v3_rr^2)
625
626
# Compute the necessary mean values
627
rho_avg = 0.5f0 * (rho_ll + rho_rr)
628
rho_mean = ln_mean(rho_ll, rho_rr)
629
beta_mean = ln_mean(beta_ll, beta_rr)
630
beta_avg = 0.5f0 * (beta_ll + beta_rr)
631
v1_avg = 0.5f0 * (v1_ll + v1_rr)
632
v2_avg = 0.5f0 * (v2_ll + v2_rr)
633
v3_avg = 0.5f0 * (v3_ll + v3_rr)
634
p_mean = 0.5f0 * rho_avg / beta_avg
635
velocity_square_avg = specific_kin_ll + specific_kin_rr
636
637
# Calculate fluxes depending on orientation
638
if orientation == 1
639
f1 = rho_mean * v1_avg
640
f2 = f1 * v1_avg + p_mean
641
f3 = f1 * v2_avg
642
f4 = f1 * v3_avg
643
f5 = f1 * 0.5f0 *
644
(1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +
645
f2 * v1_avg + f3 * v2_avg + f4 * v3_avg
646
elseif orientation == 2
647
f1 = rho_mean * v2_avg
648
f2 = f1 * v1_avg
649
f3 = f1 * v2_avg + p_mean
650
f4 = f1 * v3_avg
651
f5 = f1 * 0.5f0 *
652
(1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +
653
f2 * v1_avg + f3 * v2_avg + f4 * v3_avg
654
else
655
f1 = rho_mean * v3_avg
656
f2 = f1 * v1_avg
657
f3 = f1 * v2_avg
658
f4 = f1 * v3_avg + p_mean
659
f5 = f1 * 0.5f0 *
660
(1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +
661
f2 * v1_avg + f3 * v2_avg + f4 * v3_avg
662
end
663
664
return SVector(f1, f2, f3, f4, f5)
665
end
666
667
@inline function flux_chandrashekar(u_ll, u_rr, normal_direction::AbstractVector,
668
equations::CompressibleEulerEquations3D)
669
# Unpack left and right state
670
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
671
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
672
673
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2] +
674
v3_ll * normal_direction[3]
675
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2] +
676
v3_rr * normal_direction[3]
677
678
beta_ll = 0.5f0 * rho_ll / p_ll
679
beta_rr = 0.5f0 * rho_rr / p_rr
680
specific_kin_ll = 0.5f0 * (v1_ll^2 + v2_ll^2 + v3_ll^2)
681
specific_kin_rr = 0.5f0 * (v1_rr^2 + v2_rr^2 + v3_rr^2)
682
683
# Compute the necessary mean values
684
rho_avg = 0.5f0 * (rho_ll + rho_rr)
685
rho_mean = ln_mean(rho_ll, rho_rr)
686
beta_mean = ln_mean(beta_ll, beta_rr)
687
beta_avg = 0.5f0 * (beta_ll + beta_rr)
688
v1_avg = 0.5f0 * (v1_ll + v1_rr)
689
v2_avg = 0.5f0 * (v2_ll + v2_rr)
690
v3_avg = 0.5f0 * (v3_ll + v3_rr)
691
p_mean = 0.5f0 * rho_avg / beta_avg
692
velocity_square_avg = specific_kin_ll + specific_kin_rr
693
694
# Multiply with average of normal velocities
695
f1 = rho_mean * 0.5f0 * (v_dot_n_ll + v_dot_n_rr)
696
f2 = f1 * v1_avg + p_mean * normal_direction[1]
697
f3 = f1 * v2_avg + p_mean * normal_direction[2]
698
f4 = f1 * v3_avg + p_mean * normal_direction[3]
699
f5 = f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +
700
f2 * v1_avg + f3 * v2_avg + f4 * v3_avg
701
702
return SVector(f1, f2, f3, f4, f5)
703
end
704
705
"""
706
flux_ranocha(u_ll, u_rr, orientation_or_normal_direction,
707
equations::CompressibleEulerEquations3D)
708
709
Entropy conserving and kinetic energy preserving two-point flux by
710
- Hendrik Ranocha (2018)
711
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
712
for Hyperbolic Balance Laws
713
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
714
See also
715
- Hendrik Ranocha (2020)
716
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
717
the Euler Equations Using Summation-by-Parts Operators
718
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
719
"""
720
@inline function flux_ranocha(u_ll, u_rr, orientation::Integer,
721
equations::CompressibleEulerEquations3D)
722
# Unpack left and right state
723
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
724
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
725
726
# Compute the necessary mean values
727
rho_mean = ln_mean(rho_ll, rho_rr)
728
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
729
# in exact arithmetic since
730
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
731
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
732
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
733
v1_avg = 0.5f0 * (v1_ll + v1_rr)
734
v2_avg = 0.5f0 * (v2_ll + v2_rr)
735
v3_avg = 0.5f0 * (v3_ll + v3_rr)
736
p_avg = 0.5f0 * (p_ll + p_rr)
737
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
738
739
# Calculate fluxes depending on orientation
740
if orientation == 1
741
f1 = rho_mean * v1_avg
742
f2 = f1 * v1_avg + p_avg
743
f3 = f1 * v2_avg
744
f4 = f1 * v3_avg
745
f5 = f1 *
746
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +
747
0.5f0 * (p_ll * v1_rr + p_rr * v1_ll)
748
elseif orientation == 2
749
f1 = rho_mean * v2_avg
750
f2 = f1 * v1_avg
751
f3 = f1 * v2_avg + p_avg
752
f4 = f1 * v3_avg
753
f5 = f1 *
754
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +
755
0.5f0 * (p_ll * v2_rr + p_rr * v2_ll)
756
else # orientation == 3
757
f1 = rho_mean * v3_avg
758
f2 = f1 * v1_avg
759
f3 = f1 * v2_avg
760
f4 = f1 * v3_avg + p_avg
761
f5 = f1 *
762
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +
763
0.5f0 * (p_ll * v3_rr + p_rr * v3_ll)
764
end
765
766
return SVector(f1, f2, f3, f4, f5)
767
end
768
769
@inline function flux_ranocha(u_ll, u_rr, normal_direction::AbstractVector,
770
equations::CompressibleEulerEquations3D)
771
# Unpack left and right state
772
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
773
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
774
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2] +
775
v3_ll * normal_direction[3]
776
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2] +
777
v3_rr * normal_direction[3]
778
779
# Compute the necessary mean values
780
rho_mean = ln_mean(rho_ll, rho_rr)
781
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
782
# in exact arithmetic since
783
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
784
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
785
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
786
v1_avg = 0.5f0 * (v1_ll + v1_rr)
787
v2_avg = 0.5f0 * (v2_ll + v2_rr)
788
v3_avg = 0.5f0 * (v3_ll + v3_rr)
789
p_avg = 0.5f0 * (p_ll + p_rr)
790
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
791
792
# Calculate fluxes depending on normal_direction
793
f1 = rho_mean * 0.5f0 * (v_dot_n_ll + v_dot_n_rr)
794
f2 = f1 * v1_avg + p_avg * normal_direction[1]
795
f3 = f1 * v2_avg + p_avg * normal_direction[2]
796
f4 = f1 * v3_avg + p_avg * normal_direction[3]
797
f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
798
+
799
0.5f0 * (p_ll * v_dot_n_rr + p_rr * v_dot_n_ll))
800
801
return SVector(f1, f2, f3, f4, f5)
802
end
803
804
"""
805
splitting_steger_warming(u, orientation::Integer,
806
equations::CompressibleEulerEquations3D)
807
splitting_steger_warming(u, which::Union{Val{:minus}, Val{:plus}}
808
orientation::Integer,
809
equations::CompressibleEulerEquations3D)
810
811
Splitting of the compressible Euler flux of Steger and Warming.
812
813
Returns a tuple of the fluxes "minus" (associated with waves going into the
814
negative axis direction) and "plus" (associated with waves going into the
815
positive axis direction). If only one of the fluxes is required, use the
816
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
817
818
!!! warning "Experimental implementation (upwind SBP)"
819
This is an experimental feature and may change in future releases.
820
821
## References
822
823
- Joseph L. Steger and R. F. Warming (1979)
824
Flux Vector Splitting of the Inviscid Gasdynamic Equations
825
With Application to Finite Difference Methods
826
[NASA Technical Memorandum](https://ntrs.nasa.gov/api/citations/19790020779/downloads/19790020779.pdf)
827
"""
828
@inline function splitting_steger_warming(u, orientation::Integer,
829
equations::CompressibleEulerEquations3D)
830
fm = splitting_steger_warming(u, Val{:minus}(), orientation, equations)
831
fp = splitting_steger_warming(u, Val{:plus}(), orientation, equations)
832
return fm, fp
833
end
834
835
@inline function splitting_steger_warming(u, ::Val{:plus}, orientation::Integer,
836
equations::CompressibleEulerEquations3D)
837
rho, rho_v1, rho_v2, rho_v3, rho_e = u
838
v1 = rho_v1 / rho
839
v2 = rho_v2 / rho
840
v3 = rho_v3 / rho
841
p = (equations.gamma - 1) *
842
(rho_e - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3))
843
a = sqrt(equations.gamma * p / rho)
844
845
if orientation == 1
846
lambda1 = v1
847
lambda2 = v1 + a
848
lambda3 = v1 - a
849
850
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
851
lambda2_p = positive_part(lambda2)
852
lambda3_p = positive_part(lambda3)
853
854
alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p
855
856
rho_2gamma = 0.5f0 * rho / equations.gamma
857
f1p = rho_2gamma * alpha_p
858
f2p = rho_2gamma * (alpha_p * v1 + a * (lambda2_p - lambda3_p))
859
f3p = rho_2gamma * alpha_p * v2
860
f4p = rho_2gamma * alpha_p * v3
861
f5p = rho_2gamma *
862
(alpha_p * 0.5f0 * (v1^2 + v2^2 + v3^2) +
863
a * v1 *
864
(lambda2_p - lambda3_p)
865
+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)
866
elseif orientation == 2
867
lambda1 = v2
868
lambda2 = v2 + a
869
lambda3 = v2 - a
870
871
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
872
lambda2_p = positive_part(lambda2)
873
lambda3_p = positive_part(lambda3)
874
875
alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p
876
877
rho_2gamma = 0.5f0 * rho / equations.gamma
878
f1p = rho_2gamma * alpha_p
879
f2p = rho_2gamma * alpha_p * v1
880
f3p = rho_2gamma * (alpha_p * v2 + a * (lambda2_p - lambda3_p))
881
f4p = rho_2gamma * alpha_p * v3
882
f5p = rho_2gamma *
883
(alpha_p * 0.5f0 * (v1^2 + v2^2 + v3^2) +
884
a * v2 *
885
(lambda2_p - lambda3_p)
886
+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)
887
else # orientation == 3
888
lambda1 = v3
889
lambda2 = v3 + a
890
lambda3 = v3 - a
891
892
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
893
lambda2_p = positive_part(lambda2)
894
lambda3_p = positive_part(lambda3)
895
896
alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p
897
898
rho_2gamma = 0.5f0 * rho / equations.gamma
899
f1p = rho_2gamma * alpha_p
900
f2p = rho_2gamma * alpha_p * v1
901
f3p = rho_2gamma * alpha_p * v2
902
f4p = rho_2gamma * (alpha_p * v3 + a * (lambda2_p - lambda3_p))
903
f5p = rho_2gamma *
904
(alpha_p * 0.5f0 * (v1^2 + v2^2 + v3^2) +
905
a * v3 *
906
(lambda2_p - lambda3_p)
907
+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)
908
end
909
return SVector(f1p, f2p, f3p, f4p, f5p)
910
end
911
912
@inline function splitting_steger_warming(u, ::Val{:minus}, orientation::Integer,
913
equations::CompressibleEulerEquations3D)
914
rho, rho_v1, rho_v2, rho_v3, rho_e = u
915
v1 = rho_v1 / rho
916
v2 = rho_v2 / rho
917
v3 = rho_v3 / rho
918
p = (equations.gamma - 1) *
919
(rho_e - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3))
920
a = sqrt(equations.gamma * p / rho)
921
922
if orientation == 1
923
lambda1 = v1
924
lambda2 = v1 + a
925
lambda3 = v1 - a
926
927
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
928
lambda2_m = negative_part(lambda2)
929
lambda3_m = negative_part(lambda3)
930
931
alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m
932
933
rho_2gamma = 0.5f0 * rho / equations.gamma
934
f1m = rho_2gamma * alpha_m
935
f2m = rho_2gamma * (alpha_m * v1 + a * (lambda2_m - lambda3_m))
936
f3m = rho_2gamma * alpha_m * v2
937
f4m = rho_2gamma * alpha_m * v3
938
f5m = rho_2gamma *
939
(alpha_m * 0.5f0 * (v1^2 + v2^2 + v3^2) +
940
a * v1 *
941
(lambda2_m - lambda3_m)
942
+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)
943
elseif orientation == 2
944
lambda1 = v2
945
lambda2 = v2 + a
946
lambda3 = v2 - a
947
948
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
949
lambda2_m = negative_part(lambda2)
950
lambda3_m = negative_part(lambda3)
951
952
alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m
953
954
rho_2gamma = 0.5f0 * rho / equations.gamma
955
f1m = rho_2gamma * alpha_m
956
f2m = rho_2gamma * alpha_m * v1
957
f3m = rho_2gamma * (alpha_m * v2 + a * (lambda2_m - lambda3_m))
958
f4m = rho_2gamma * alpha_m * v3
959
f5m = rho_2gamma *
960
(alpha_m * 0.5f0 * (v1^2 + v2^2 + v3^2) +
961
a * v2 *
962
(lambda2_m - lambda3_m)
963
+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)
964
else # orientation == 3
965
lambda1 = v3
966
lambda2 = v3 + a
967
lambda3 = v3 - a
968
969
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
970
lambda2_m = negative_part(lambda2)
971
lambda3_m = negative_part(lambda3)
972
973
alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m
974
975
rho_2gamma = 0.5f0 * rho / equations.gamma
976
f1m = rho_2gamma * alpha_m
977
f2m = rho_2gamma * alpha_m * v1
978
f3m = rho_2gamma * alpha_m * v2
979
f4m = rho_2gamma * (alpha_m * v3 + a * (lambda2_m - lambda3_m))
980
f5m = rho_2gamma *
981
(alpha_m * 0.5f0 * (v1^2 + v2^2 + v3^2) +
982
a * v3 *
983
(lambda2_m - lambda3_m)
984
+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)
985
end
986
return SVector(f1m, f2m, f3m, f4m, f5m)
987
end
988
989
"""
990
FluxLMARS(c)(u_ll, u_rr, orientation_or_normal_direction,
991
equations::CompressibleEulerEquations3D)
992
993
Low Mach number approximate Riemann solver (LMARS) for atmospheric flows using
994
an estimate `c` of the speed of sound.
995
996
References:
997
- Xi Chen et al. (2013)
998
A Control-Volume Model of the Compressible Euler Equations with a Vertical
999
Lagrangian Coordinate
1000
[DOI: 10.1175/MWR-D-12-00129.1](https://doi.org/10.1175/mwr-d-12-00129.1)
1001
"""
1002
@inline function (flux_lmars::FluxLMARS)(u_ll, u_rr, orientation::Integer,
1003
equations::CompressibleEulerEquations3D)
1004
c = flux_lmars.speed_of_sound
1005
1006
# Unpack left and right state
1007
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1008
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1009
1010
if orientation == 1
1011
v_ll = v1_ll
1012
v_rr = v1_rr
1013
elseif orientation == 2
1014
v_ll = v2_ll
1015
v_rr = v2_rr
1016
else # orientation == 3
1017
v_ll = v3_ll
1018
v_rr = v3_rr
1019
end
1020
1021
rho = 0.5f0 * (rho_ll + rho_rr)
1022
p = 0.5f0 * (p_ll + p_rr) - 0.5f0 * c * rho * (v_rr - v_ll)
1023
v = 0.5f0 * (v_ll + v_rr) - 1 / (2 * c * rho) * (p_rr - p_ll)
1024
1025
# We treat the energy term analogous to the potential temperature term in the paper by
1026
# Chen et al., i.e. we use p_ll and p_rr, and not p
1027
if v >= 0
1028
f1, f2, f3, f4, f5 = v * u_ll
1029
f5 = f5 + p_ll * v
1030
else
1031
f1, f2, f3, f4, f5 = v * u_rr
1032
f5 = f5 + p_rr * v
1033
end
1034
1035
if orientation == 1
1036
f2 += p
1037
elseif orientation == 2
1038
f3 += p
1039
else # orientation == 3
1040
f4 += p
1041
end
1042
1043
return SVector(f1, f2, f3, f4, f5)
1044
end
1045
1046
@inline function (flux_lmars::FluxLMARS)(u_ll, u_rr, normal_direction::AbstractVector,
1047
equations::CompressibleEulerEquations3D)
1048
c = flux_lmars.speed_of_sound
1049
1050
# Unpack left and right state
1051
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1052
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1053
1054
v_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2] +
1055
v3_ll * normal_direction[3]
1056
v_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2] +
1057
v3_rr * normal_direction[3]
1058
1059
# Note that this is the same as computing v_ll and v_rr with a normalized normal vector
1060
# and then multiplying v by `norm_` again, but this version is slightly faster.
1061
norm_ = norm(normal_direction)
1062
1063
rho = 0.5f0 * (rho_ll + rho_rr)
1064
p = 0.5f0 * (p_ll + p_rr) - 0.5f0 * c * rho * (v_rr - v_ll) / norm_
1065
v = 0.5f0 * (v_ll + v_rr) - 1 / (2 * c * rho) * (p_rr - p_ll) * norm_
1066
1067
# We treat the energy term analogous to the potential temperature term in the paper by
1068
# Chen et al., i.e. we use p_ll and p_rr, and not p
1069
if v >= 0
1070
f1, f2, f3, f4, f5 = v * u_ll
1071
f5 = f5 + p_ll * v
1072
else
1073
f1, f2, f3, f4, f5 = v * u_rr
1074
f5 = f5 + p_rr * v
1075
end
1076
1077
return SVector(f1,
1078
f2 + p * normal_direction[1],
1079
f3 + p * normal_direction[2],
1080
f4 + p * normal_direction[3],
1081
f5)
1082
end
1083
1084
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation as the
1085
# maximum velocity magnitude plus the maximum speed of sound
1086
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
1087
equations::CompressibleEulerEquations3D)
1088
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1089
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1090
1091
# Get the velocity value in the appropriate direction
1092
if orientation == 1
1093
v_ll = v1_ll
1094
v_rr = v1_rr
1095
elseif orientation == 2
1096
v_ll = v2_ll
1097
v_rr = v2_rr
1098
else # orientation == 3
1099
v_ll = v3_ll
1100
v_rr = v3_rr
1101
end
1102
# Calculate sound speeds
1103
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1104
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1105
1106
return max(abs(v_ll), abs(v_rr)) + max(c_ll, c_rr)
1107
end
1108
1109
@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
1110
equations::CompressibleEulerEquations3D)
1111
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1112
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1113
1114
# Calculate normal velocities and sound speed
1115
# left
1116
v_ll = (v1_ll * normal_direction[1]
1117
+ v2_ll * normal_direction[2]
1118
+ v3_ll * normal_direction[3])
1119
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1120
# right
1121
v_rr = (v1_rr * normal_direction[1]
1122
+ v2_rr * normal_direction[2]
1123
+ v3_rr * normal_direction[3])
1124
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1125
1126
return max(abs(v_ll), abs(v_rr)) + max(c_ll, c_rr) * norm(normal_direction)
1127
end
1128
1129
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
1130
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
1131
equations::CompressibleEulerEquations3D)
1132
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1133
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1134
1135
# Get the velocity value in the appropriate direction
1136
if orientation == 1
1137
v_ll = v1_ll
1138
v_rr = v1_rr
1139
elseif orientation == 2
1140
v_ll = v2_ll
1141
v_rr = v2_rr
1142
else # orientation == 3
1143
v_ll = v3_ll
1144
v_rr = v3_rr
1145
end
1146
# Calculate sound speeds
1147
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1148
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1149
1150
return max(abs(v_ll) + c_ll, abs(v_rr) + c_rr)
1151
end
1152
1153
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
1154
@inline function max_abs_speed(u_ll, u_rr, normal_direction::AbstractVector,
1155
equations::CompressibleEulerEquations3D)
1156
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1157
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1158
1159
# Calculate normal velocities and sound speeds
1160
# left
1161
v_ll = (v1_ll * normal_direction[1]
1162
+ v2_ll * normal_direction[2]
1163
+ v3_ll * normal_direction[3])
1164
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1165
# right
1166
v_rr = (v1_rr * normal_direction[1]
1167
+ v2_rr * normal_direction[2]
1168
+ v3_rr * normal_direction[3])
1169
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1170
1171
norm_ = norm(normal_direction)
1172
return max(abs(v_ll) + c_ll * norm_, abs(v_rr) + c_rr * norm_)
1173
end
1174
1175
# Calculate estimates for minimum and maximum wave speeds for HLL-type fluxes
1176
@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,
1177
equations::CompressibleEulerEquations3D)
1178
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1179
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1180
1181
if orientation == 1 # x-direction
1182
λ_min = v1_ll - sqrt(equations.gamma * p_ll / rho_ll)
1183
λ_max = v1_rr + sqrt(equations.gamma * p_rr / rho_rr)
1184
elseif orientation == 2 # y-direction
1185
λ_min = v2_ll - sqrt(equations.gamma * p_ll / rho_ll)
1186
λ_max = v2_rr + sqrt(equations.gamma * p_rr / rho_rr)
1187
else # z-direction
1188
λ_min = v3_ll - sqrt(equations.gamma * p_ll / rho_ll)
1189
λ_max = v3_rr + sqrt(equations.gamma * p_rr / rho_rr)
1190
end
1191
1192
return λ_min, λ_max
1193
end
1194
1195
@inline function min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
1196
equations::CompressibleEulerEquations3D)
1197
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1198
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1199
1200
v_normal_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2] +
1201
v3_ll * normal_direction[3]
1202
v_normal_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2] +
1203
v3_rr * normal_direction[3]
1204
1205
norm_ = norm(normal_direction)
1206
# The v_normals are already scaled by the norm
1207
λ_min = v_normal_ll - sqrt(equations.gamma * p_ll / rho_ll) * norm_
1208
λ_max = v_normal_rr + sqrt(equations.gamma * p_rr / rho_rr) * norm_
1209
1210
return λ_min, λ_max
1211
end
1212
1213
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
1214
@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,
1215
equations::CompressibleEulerEquations3D)
1216
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1217
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1218
1219
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1220
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1221
1222
if orientation == 1 # x-direction
1223
λ_min = min(v1_ll - c_ll, v1_rr - c_rr)
1224
λ_max = max(v1_ll + c_ll, v1_rr + c_rr)
1225
elseif orientation == 2 # y-direction
1226
λ_min = min(v2_ll - c_ll, v2_rr - c_rr)
1227
λ_max = max(v2_ll + c_ll, v2_rr + c_rr)
1228
else # z-direction
1229
λ_min = min(v3_ll - c_ll, v3_rr - c_rr)
1230
λ_max = max(v3_ll + c_ll, v3_rr + c_rr)
1231
end
1232
1233
return λ_min, λ_max
1234
end
1235
1236
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
1237
@inline function min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector,
1238
equations::CompressibleEulerEquations3D)
1239
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1240
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1241
1242
norm_ = norm(normal_direction)
1243
1244
c_ll = sqrt(equations.gamma * p_ll / rho_ll) * norm_
1245
c_rr = sqrt(equations.gamma * p_rr / rho_rr) * norm_
1246
1247
v_normal_ll = v1_ll * normal_direction[1] +
1248
v2_ll * normal_direction[2] +
1249
v3_ll * normal_direction[3]
1250
v_normal_rr = v1_rr * normal_direction[1] +
1251
v2_rr * normal_direction[2] +
1252
v3_rr * normal_direction[3]
1253
1254
# The v_normals are already scaled by the norm
1255
λ_min = min(v_normal_ll - c_ll, v_normal_rr - c_rr)
1256
λ_max = max(v_normal_ll + c_ll, v_normal_rr + c_rr)
1257
1258
return λ_min, λ_max
1259
end
1260
1261
# Rotate normal vector to x-axis; normal, tangent1 and tangent2 need to be orthonormal
1262
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the directions
1263
# has been normalized prior to this rotation of the state vector
1264
@inline function rotate_to_x(u, normal_vector, tangent1, tangent2,
1265
equations::CompressibleEulerEquations3D)
1266
# Multiply with [ 1 0 0 0 0;
1267
# 0 ― normal_vector ― 0;
1268
# 0 ― tangent1 ― 0;
1269
# 0 ― tangent2 ― 0;
1270
# 0 0 0 0 1 ]
1271
return SVector(u[1],
1272
normal_vector[1] * u[2] + normal_vector[2] * u[3] +
1273
normal_vector[3] * u[4],
1274
tangent1[1] * u[2] + tangent1[2] * u[3] + tangent1[3] * u[4],
1275
tangent2[1] * u[2] + tangent2[2] * u[3] + tangent2[3] * u[4],
1276
u[5])
1277
end
1278
1279
@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,
1280
normal_direction::AbstractVector,
1281
equations::CompressibleEulerEquations3D)
1282
(; gamma) = equations
1283
1284
# Step 1:
1285
# Rotate solution into the appropriate direction
1286
1287
norm_ = norm(normal_direction)
1288
# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later
1289
normal_vector = normal_direction / norm_
1290
1291
# Some vector that can't be identical to normal_vector (unless normal_vector == 0)
1292
tangent1 = SVector(normal_direction[2], normal_direction[3], -normal_direction[1])
1293
# Orthogonal projection
1294
tangent1 -= dot(normal_vector, tangent1) * normal_vector
1295
tangent1 = normalize(tangent1)
1296
1297
# Third orthogonal vector
1298
tangent2 = normalize(cross(normal_direction, tangent1))
1299
1300
u_ll_rotated = rotate_to_x(u_ll, normal_vector, tangent1, tangent2, equations)
1301
u_rr_rotated = rotate_to_x(u_rr, normal_vector, tangent1, tangent2, equations)
1302
1303
# Step 2:
1304
# Compute the averages using the rotated variables
1305
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll_rotated, equations)
1306
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr_rotated, equations)
1307
1308
b_ll = rho_ll / (2 * p_ll)
1309
b_rr = rho_rr / (2 * p_rr)
1310
1311
rho_log = ln_mean(rho_ll, rho_rr)
1312
b_log = ln_mean(b_ll, b_rr)
1313
v1_avg = 0.5f0 * (v1_ll + v1_rr)
1314
v2_avg = 0.5f0 * (v2_ll + v2_rr)
1315
v3_avg = 0.5f0 * (v3_ll + v3_rr)
1316
p_avg = 0.5f0 * (rho_ll + rho_rr) / (b_ll + b_rr)
1317
v_squared_bar = v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr
1318
h_bar = gamma / (2 * b_log * (gamma - 1)) + 0.5f0 * v_squared_bar
1319
c_bar = sqrt(gamma * p_avg / rho_log)
1320
1321
# Step 3:
1322
# Build the dissipation term as given in Appendix A of the paper
1323
# - A. R. Winters, D. Derigs, G. Gassner, S. Walch, A uniquely defined entropy stable matrix dissipation operator
1324
# for high Mach number ideal MHD and compressible Euler simulations (2017). Journal of Computational Physics.
1325
# [DOI: 10.1016/j.jcp.2016.12.006](https://doi.org/10.1016/j.jcp.2016.12.006).
1326
1327
# Get entropy variables jump in the rotated variables
1328
w_jump = cons2entropy(u_rr_rotated, equations) -
1329
cons2entropy(u_ll_rotated, equations)
1330
1331
# Entries of the diagonal scaling matrix where D = ABS(\Lambda)T
1332
lambda_1 = abs(v1_avg - c_bar) * rho_log / (2 * gamma)
1333
lambda_2 = abs(v1_avg) * rho_log * (gamma - 1) / gamma
1334
lambda_3 = abs(v1_avg) * p_avg # scaled repeated eigenvalue in the tangential direction
1335
lambda_5 = abs(v1_avg + c_bar) * rho_log / (2 * gamma)
1336
D = SVector(lambda_1, lambda_2, lambda_3, lambda_3, lambda_5)
1337
1338
# Entries of the right eigenvector matrix (others have already been precomputed)
1339
r21 = v1_avg - c_bar
1340
r25 = v1_avg + c_bar
1341
r51 = h_bar - v1_avg * c_bar
1342
r52 = 0.5f0 * v_squared_bar
1343
r55 = h_bar + v1_avg * c_bar
1344
1345
# Build R and transpose of R matrices
1346
R = @SMatrix [[1;; 1;; 0;; 0;; 1];
1347
[r21;; v1_avg;; 0;; 0;; r25];
1348
[v2_avg;; v2_avg;; 1;; 0;; v2_avg];
1349
[v3_avg;; v3_avg;; 0;; 1;; v3_avg];
1350
[r51;; r52;; v2_avg;; v3_avg;; r55]]
1351
1352
RT = @SMatrix [[1;; r21;; v2_avg;; v3_avg;; r51];
1353
[1;; v1_avg;; v2_avg;; v3_avg;; r52];
1354
[0;; 0;; 1;; 0;; v2_avg];
1355
[0;; 0;; 0;; 1;; v3_avg];
1356
[1;; r25;; v2_avg;; v3_avg;; r55]]
1357
1358
# Compute the dissipation term R * D * R^T * [[w]] from right-to-left
1359
1360
# First comes R^T * [[w]]
1361
diss = RT * w_jump
1362
# Next multiply with the eigenvalues and Barth scaling
1363
diss = D .* diss
1364
# Finally apply the remaining eigenvector matrix
1365
diss = R * diss
1366
1367
# Step 4:
1368
# Do not forget to backrotate and then return with proper normalization scaling
1369
return -0.5f0 * rotate_from_x(diss, normal_vector, tangent1, tangent2, equations) *
1370
norm_
1371
end
1372
1373
# Rotate x-axis to normal vector; normal, tangent1 and tangent2 need to be orthonormal
1374
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the directions
1375
# has been normalized prior to this back-rotation of the state vector
1376
@inline function rotate_from_x(u, normal_vector, tangent1, tangent2,
1377
equations::CompressibleEulerEquations3D)
1378
# Multiply with [ 1 0 0 0 0;
1379
# 0 | | | 0;
1380
# 0 normal_vector tangent1 tangent2 0;
1381
# 0 | | | 0;
1382
# 0 0 0 0 1 ]
1383
return SVector(u[1],
1384
normal_vector[1] * u[2] + tangent1[1] * u[3] + tangent2[1] * u[4],
1385
normal_vector[2] * u[2] + tangent1[2] * u[3] + tangent2[2] * u[4],
1386
normal_vector[3] * u[2] + tangent1[3] * u[3] + tangent2[3] * u[4],
1387
u[5])
1388
end
1389
1390
"""
1391
flux_hllc(u_ll, u_rr, orientation_or_normal_direction, equations::CompressibleEulerEquations3D)
1392
1393
Computes the HLLC flux (HLL with Contact) for compressible Euler equations developed by E.F. Toro
1394
[Lecture slides](http://www.prague-sum.com/download/2012/Toro_2-HLLC-RiemannSolver.pdf)
1395
Signal speeds: [DOI: 10.1137/S1064827593260140](https://doi.org/10.1137/S1064827593260140)
1396
"""
1397
function flux_hllc(u_ll, u_rr, orientation::Integer,
1398
equations::CompressibleEulerEquations3D)
1399
# Calculate primitive variables and speed of sound
1400
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll = u_ll
1401
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr = u_rr
1402
1403
v1_ll = rho_v1_ll / rho_ll
1404
v2_ll = rho_v2_ll / rho_ll
1405
v3_ll = rho_v3_ll / rho_ll
1406
e_ll = rho_e_ll / rho_ll
1407
p_ll = (equations.gamma - 1) *
1408
(rho_e_ll - 0.5f0 * rho_ll * (v1_ll^2 + v2_ll^2 + v3_ll^2))
1409
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1410
1411
v1_rr = rho_v1_rr / rho_rr
1412
v2_rr = rho_v2_rr / rho_rr
1413
v3_rr = rho_v3_rr / rho_rr
1414
e_rr = rho_e_rr / rho_rr
1415
p_rr = (equations.gamma - 1) *
1416
(rho_e_rr - 0.5f0 * rho_rr * (v1_rr^2 + v2_rr^2 + v3_rr^2))
1417
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1418
1419
# Obtain left and right fluxes
1420
f_ll = flux(u_ll, orientation, equations)
1421
f_rr = flux(u_rr, orientation, equations)
1422
1423
# Compute Roe averages
1424
sqrt_rho_ll = sqrt(rho_ll)
1425
sqrt_rho_rr = sqrt(rho_rr)
1426
sum_sqrt_rho = sqrt_rho_ll + sqrt_rho_rr
1427
if orientation == 1 # x-direction
1428
vel_L = v1_ll
1429
vel_R = v1_rr
1430
elseif orientation == 2 # y-direction
1431
vel_L = v2_ll
1432
vel_R = v2_rr
1433
else # z-direction
1434
vel_L = v3_ll
1435
vel_R = v3_rr
1436
end
1437
vel_roe = (sqrt_rho_ll * vel_L + sqrt_rho_rr * vel_R) / sum_sqrt_rho
1438
v1_roe = sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr
1439
v2_roe = sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr
1440
v3_roe = sqrt_rho_ll * v3_ll + sqrt_rho_rr * v3_rr
1441
vel_roe_mag = (v1_roe^2 + v2_roe^2 + v3_roe^2) / sum_sqrt_rho^2
1442
H_ll = (rho_e_ll + p_ll) / rho_ll
1443
H_rr = (rho_e_rr + p_rr) / rho_rr
1444
H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) / sum_sqrt_rho
1445
c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * vel_roe_mag))
1446
Ssl = min(vel_L - c_ll, vel_roe - c_roe)
1447
Ssr = max(vel_R + c_rr, vel_roe + c_roe)
1448
sMu_L = Ssl - vel_L
1449
sMu_R = Ssr - vel_R
1450
1451
if Ssl >= 0
1452
f1 = f_ll[1]
1453
f2 = f_ll[2]
1454
f3 = f_ll[3]
1455
f4 = f_ll[4]
1456
f5 = f_ll[5]
1457
elseif Ssr <= 0
1458
f1 = f_rr[1]
1459
f2 = f_rr[2]
1460
f3 = f_rr[3]
1461
f4 = f_rr[4]
1462
f5 = f_rr[5]
1463
else
1464
SStar = (p_rr - p_ll + rho_ll * vel_L * sMu_L - rho_rr * vel_R * sMu_R) /
1465
(rho_ll * sMu_L - rho_rr * sMu_R)
1466
if Ssl <= 0 <= SStar
1467
densStar = rho_ll * sMu_L / (Ssl - SStar)
1468
enerStar = e_ll + (SStar - vel_L) * (SStar + p_ll / (rho_ll * sMu_L))
1469
UStar1 = densStar
1470
UStar5 = densStar * enerStar
1471
if orientation == 1 # x-direction
1472
UStar2 = densStar * SStar
1473
UStar3 = densStar * v2_ll
1474
UStar4 = densStar * v3_ll
1475
elseif orientation == 2 # y-direction
1476
UStar2 = densStar * v1_ll
1477
UStar3 = densStar * SStar
1478
UStar4 = densStar * v3_ll
1479
else # z-direction
1480
UStar2 = densStar * v1_ll
1481
UStar3 = densStar * v2_ll
1482
UStar4 = densStar * SStar
1483
end
1484
f1 = f_ll[1] + Ssl * (UStar1 - rho_ll)
1485
f2 = f_ll[2] + Ssl * (UStar2 - rho_v1_ll)
1486
f3 = f_ll[3] + Ssl * (UStar3 - rho_v2_ll)
1487
f4 = f_ll[4] + Ssl * (UStar4 - rho_v3_ll)
1488
f5 = f_ll[5] + Ssl * (UStar5 - rho_e_ll)
1489
else
1490
densStar = rho_rr * sMu_R / (Ssr - SStar)
1491
enerStar = e_rr + (SStar - vel_R) * (SStar + p_rr / (rho_rr * sMu_R))
1492
UStar1 = densStar
1493
UStar5 = densStar * enerStar
1494
if orientation == 1 # x-direction
1495
UStar2 = densStar * SStar
1496
UStar3 = densStar * v2_rr
1497
UStar4 = densStar * v3_rr
1498
elseif orientation == 2 # y-direction
1499
UStar2 = densStar * v1_rr
1500
UStar3 = densStar * SStar
1501
UStar4 = densStar * v3_rr
1502
else # z-direction
1503
UStar2 = densStar * v1_rr
1504
UStar3 = densStar * v2_rr
1505
UStar4 = densStar * SStar
1506
end
1507
f1 = f_rr[1] + Ssr * (UStar1 - rho_rr)
1508
f2 = f_rr[2] + Ssr * (UStar2 - rho_v1_rr)
1509
f3 = f_rr[3] + Ssr * (UStar3 - rho_v2_rr)
1510
f4 = f_rr[4] + Ssr * (UStar4 - rho_v3_rr)
1511
f5 = f_rr[5] + Ssr * (UStar5 - rho_e_rr)
1512
end
1513
end
1514
return SVector(f1, f2, f3, f4, f5)
1515
end
1516
1517
function flux_hllc(u_ll, u_rr, normal_direction::AbstractVector,
1518
equations::CompressibleEulerEquations3D)
1519
# Calculate primitive variables and speed of sound
1520
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1521
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1522
1523
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2] +
1524
v3_ll * normal_direction[3]
1525
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2] +
1526
v3_rr * normal_direction[3]
1527
1528
norm_ = norm(normal_direction)
1529
norm_sq = norm_ * norm_
1530
inv_norm_sq = inv(norm_sq)
1531
1532
c_ll = sqrt(equations.gamma * p_ll / rho_ll) * norm_
1533
c_rr = sqrt(equations.gamma * p_rr / rho_rr) * norm_
1534
1535
# Obtain left and right fluxes
1536
f_ll = flux(u_ll, normal_direction, equations)
1537
f_rr = flux(u_rr, normal_direction, equations)
1538
1539
# Compute Roe averages
1540
sqrt_rho_ll = sqrt(rho_ll)
1541
sqrt_rho_rr = sqrt(rho_rr)
1542
sum_sqrt_rho = sqrt_rho_ll + sqrt_rho_rr
1543
1544
v1_roe = (sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr) / sum_sqrt_rho
1545
v2_roe = (sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr) / sum_sqrt_rho
1546
v3_roe = (sqrt_rho_ll * v3_ll + sqrt_rho_rr * v3_rr) / sum_sqrt_rho
1547
vel_roe = v1_roe * normal_direction[1] + v2_roe * normal_direction[2] +
1548
v3_roe * normal_direction[3]
1549
vel_roe_mag = v1_roe^2 + v2_roe^2 + v3_roe^2
1550
1551
e_ll = u_ll[5] / rho_ll
1552
e_rr = u_rr[5] / rho_rr
1553
1554
H_ll = (u_ll[5] + p_ll) / rho_ll
1555
H_rr = (u_rr[5] + p_rr) / rho_rr
1556
1557
H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) / sum_sqrt_rho
1558
c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * vel_roe_mag)) * norm_
1559
1560
Ssl = min(v_dot_n_ll - c_ll, vel_roe - c_roe)
1561
Ssr = max(v_dot_n_rr + c_rr, vel_roe + c_roe)
1562
sMu_L = Ssl - v_dot_n_ll
1563
sMu_R = Ssr - v_dot_n_rr
1564
1565
if Ssl >= 0
1566
f1 = f_ll[1]
1567
f2 = f_ll[2]
1568
f3 = f_ll[3]
1569
f4 = f_ll[4]
1570
f5 = f_ll[5]
1571
elseif Ssr <= 0
1572
f1 = f_rr[1]
1573
f2 = f_rr[2]
1574
f3 = f_rr[3]
1575
f4 = f_rr[4]
1576
f5 = f_rr[5]
1577
else
1578
SStar = (rho_ll * v_dot_n_ll * sMu_L - rho_rr * v_dot_n_rr * sMu_R +
1579
(p_rr - p_ll) * norm_sq) / (rho_ll * sMu_L - rho_rr * sMu_R)
1580
if Ssl <= 0 <= SStar
1581
densStar = rho_ll * sMu_L / (Ssl - SStar)
1582
enerStar = e_ll +
1583
(SStar - v_dot_n_ll) *
1584
(SStar * inv_norm_sq + p_ll / (rho_ll * sMu_L))
1585
UStar1 = densStar
1586
UStar2 = densStar *
1587
(v1_ll + (SStar - v_dot_n_ll) * normal_direction[1] * inv_norm_sq)
1588
UStar3 = densStar *
1589
(v2_ll + (SStar - v_dot_n_ll) * normal_direction[2] * inv_norm_sq)
1590
UStar4 = densStar *
1591
(v3_ll + (SStar - v_dot_n_ll) * normal_direction[3] * inv_norm_sq)
1592
UStar5 = densStar * enerStar
1593
f1 = f_ll[1] + Ssl * (UStar1 - u_ll[1])
1594
f2 = f_ll[2] + Ssl * (UStar2 - u_ll[2])
1595
f3 = f_ll[3] + Ssl * (UStar3 - u_ll[3])
1596
f4 = f_ll[4] + Ssl * (UStar4 - u_ll[4])
1597
f5 = f_ll[5] + Ssl * (UStar5 - u_ll[5])
1598
else
1599
densStar = rho_rr * sMu_R / (Ssr - SStar)
1600
enerStar = e_rr +
1601
(SStar - v_dot_n_rr) *
1602
(SStar * inv_norm_sq + p_rr / (rho_rr * sMu_R))
1603
UStar1 = densStar
1604
UStar2 = densStar *
1605
(v1_rr + (SStar - v_dot_n_rr) * normal_direction[1] * inv_norm_sq)
1606
UStar3 = densStar *
1607
(v2_rr + (SStar - v_dot_n_rr) * normal_direction[2] * inv_norm_sq)
1608
UStar4 = densStar *
1609
(v3_rr + (SStar - v_dot_n_rr) * normal_direction[3] * inv_norm_sq)
1610
UStar5 = densStar * enerStar
1611
f1 = f_rr[1] + Ssr * (UStar1 - u_rr[1])
1612
f2 = f_rr[2] + Ssr * (UStar2 - u_rr[2])
1613
f3 = f_rr[3] + Ssr * (UStar3 - u_rr[3])
1614
f4 = f_rr[4] + Ssr * (UStar4 - u_rr[4])
1615
f5 = f_rr[5] + Ssr * (UStar5 - u_rr[5])
1616
end
1617
end
1618
return SVector(f1, f2, f3, f4, f5)
1619
end
1620
1621
"""
1622
min_max_speed_einfeldt(u_ll, u_rr, orientation, equations::CompressibleEulerEquations3D)
1623
1624
Computes the HLLE (Harten-Lax-van Leer-Einfeldt) flux for the compressible Euler equations.
1625
Special estimates of the signal velocites and linearization of the Riemann problem developed
1626
by Einfeldt to ensure that the internal energy and density remain positive during the computation
1627
of the numerical flux.
1628
1629
- Bernd Einfeldt (1988)
1630
On Godunov-type methods for gas dynamics.
1631
[DOI: 10.1137/0725021](https://doi.org/10.1137/0725021)
1632
- Bernd Einfeldt, Claus-Dieter Munz, Philip L. Roe and Björn Sjögreen (1991)
1633
On Godunov-type methods near low densities.
1634
[DOI: 10.1016/0021-9991(91)90211-3](https://doi.org/10.1016/0021-9991(91)90211-3)
1635
"""
1636
@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,
1637
equations::CompressibleEulerEquations3D)
1638
# Calculate primitive variables, enthalpy and speed of sound
1639
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1640
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1641
1642
# `u_ll[5]` is total energy `rho_e_ll` on the left
1643
H_ll = (u_ll[5] + p_ll) / rho_ll
1644
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1645
1646
# `u_rr[5]` is total energy `rho_e_rr` on the right
1647
H_rr = (u_rr[5] + p_rr) / rho_rr
1648
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1649
1650
# Compute Roe averages
1651
sqrt_rho_ll = sqrt(rho_ll)
1652
sqrt_rho_rr = sqrt(rho_rr)
1653
inv_sum_sqrt_rho = inv(sqrt_rho_ll + sqrt_rho_rr)
1654
1655
v1_roe = (sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr) * inv_sum_sqrt_rho
1656
v2_roe = (sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr) * inv_sum_sqrt_rho
1657
v3_roe = (sqrt_rho_ll * v3_ll + sqrt_rho_rr * v3_rr) * inv_sum_sqrt_rho
1658
v_roe_mag = v1_roe^2 + v2_roe^2 + v3_roe^2
1659
1660
H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) * inv_sum_sqrt_rho
1661
c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * v_roe_mag))
1662
1663
# Compute convenience constant for positivity preservation, see
1664
# https://doi.org/10.1016/0021-9991(91)90211-3
1665
beta = sqrt(0.5f0 * (equations.gamma - 1) / equations.gamma)
1666
1667
# Estimate the edges of the Riemann fan (with positivity conservation)
1668
if orientation == 1 # x-direction
1669
SsL = min(v1_roe - c_roe, v1_ll - beta * c_ll, 0)
1670
SsR = max(v1_roe + c_roe, v1_rr + beta * c_rr, 0)
1671
elseif orientation == 2 # y-direction
1672
SsL = min(v2_roe - c_roe, v2_ll - beta * c_ll, 0)
1673
SsR = max(v2_roe + c_roe, v2_rr + beta * c_rr, 0)
1674
else # z-direction
1675
SsL = min(v3_roe - c_roe, v3_ll - beta * c_ll, 0)
1676
SsR = max(v3_roe + c_roe, v3_rr + beta * c_rr, 0)
1677
end
1678
1679
return SsL, SsR
1680
end
1681
1682
"""
1683
min_max_speed_einfeldt(u_ll, u_rr, normal_direction, equations::CompressibleEulerEquations3D)
1684
1685
Computes the HLLE (Harten-Lax-van Leer-Einfeldt) flux for the compressible Euler equations.
1686
Special estimates of the signal velocites and linearization of the Riemann problem developed
1687
by Einfeldt to ensure that the internal energy and density remain positive during the computation
1688
of the numerical flux.
1689
1690
- Bernd Einfeldt (1988)
1691
On Godunov-type methods for gas dynamics.
1692
[DOI: 10.1137/0725021](https://doi.org/10.1137/0725021)
1693
- Bernd Einfeldt, Claus-Dieter Munz, Philip L. Roe and Björn Sjögreen (1991)
1694
On Godunov-type methods near low densities.
1695
[DOI: 10.1016/0021-9991(91)90211-3](https://doi.org/10.1016/0021-9991(91)90211-3)
1696
"""
1697
@inline function min_max_speed_einfeldt(u_ll, u_rr, normal_direction::AbstractVector,
1698
equations::CompressibleEulerEquations3D)
1699
# Calculate primitive variables, enthalpy and speed of sound
1700
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = cons2prim(u_ll, equations)
1701
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = cons2prim(u_rr, equations)
1702
1703
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2] +
1704
v3_ll * normal_direction[3]
1705
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2] +
1706
v3_rr * normal_direction[3]
1707
1708
norm_ = norm(normal_direction)
1709
1710
# `u_ll[5]` is total energy `rho_e_ll` on the left
1711
H_ll = (u_ll[5] + p_ll) / rho_ll
1712
c_ll = sqrt(equations.gamma * p_ll / rho_ll) * norm_
1713
1714
# `u_rr[5]` is total energy `rho_e_rr` on the right
1715
H_rr = (u_rr[5] + p_rr) / rho_rr
1716
c_rr = sqrt(equations.gamma * p_rr / rho_rr) * norm_
1717
1718
# Compute Roe averages
1719
sqrt_rho_ll = sqrt(rho_ll)
1720
sqrt_rho_rr = sqrt(rho_rr)
1721
inv_sum_sqrt_rho = inv(sqrt_rho_ll + sqrt_rho_rr)
1722
1723
v1_roe = (sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr) * inv_sum_sqrt_rho
1724
v2_roe = (sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr) * inv_sum_sqrt_rho
1725
v3_roe = (sqrt_rho_ll * v3_ll + sqrt_rho_rr * v3_rr) * inv_sum_sqrt_rho
1726
v_roe = v1_roe * normal_direction[1] + v2_roe * normal_direction[2] +
1727
v3_roe * normal_direction[3]
1728
v_roe_mag = v1_roe^2 + v2_roe^2 + v3_roe^2
1729
1730
H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) * inv_sum_sqrt_rho
1731
c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * v_roe_mag)) * norm_
1732
1733
# Compute convenience constant for positivity preservation, see
1734
# https://doi.org/10.1016/0021-9991(91)90211-3
1735
beta = sqrt(0.5f0 * (equations.gamma - 1) / equations.gamma)
1736
1737
# Estimate the edges of the Riemann fan (with positivity conservation)
1738
SsL = min(v_roe - c_roe, v_dot_n_ll - beta * c_ll, 0)
1739
SsR = max(v_roe + c_roe, v_dot_n_rr + beta * c_rr, 0)
1740
1741
return SsL, SsR
1742
end
1743
1744
@inline function max_abs_speeds(u, equations::CompressibleEulerEquations3D)
1745
rho, v1, v2, v3, p = cons2prim(u, equations)
1746
c = sqrt(equations.gamma * p / rho)
1747
1748
return abs(v1) + c, abs(v2) + c, abs(v3) + c
1749
end
1750
1751
# Convert conservative variables to primitive
1752
@inline function cons2prim(u, equations::CompressibleEulerEquations3D)
1753
rho, rho_v1, rho_v2, rho_v3, rho_e = u
1754
1755
v1 = rho_v1 / rho
1756
v2 = rho_v2 / rho
1757
v3 = rho_v3 / rho
1758
p = (equations.gamma - 1) *
1759
(rho_e - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3))
1760
1761
return SVector(rho, v1, v2, v3, p)
1762
end
1763
1764
# Convert conservative variables to entropy
1765
@inline function cons2entropy(u, equations::CompressibleEulerEquations3D)
1766
rho, rho_v1, rho_v2, rho_v3, rho_e = u
1767
1768
v1 = rho_v1 / rho
1769
v2 = rho_v2 / rho
1770
v3 = rho_v3 / rho
1771
v_square = v1^2 + v2^2 + v3^2
1772
p = (equations.gamma - 1) * (rho_e - 0.5f0 * rho * v_square)
1773
s = log(p) - equations.gamma * log(rho)
1774
rho_p = rho / p
1775
1776
w1 = (equations.gamma - s) * equations.inv_gamma_minus_one -
1777
0.5f0 * rho_p * v_square
1778
w2 = rho_p * v1
1779
w3 = rho_p * v2
1780
w4 = rho_p * v3
1781
w5 = -rho_p
1782
1783
return SVector(w1, w2, w3, w4, w5)
1784
end
1785
1786
@inline function entropy2cons(w, equations::CompressibleEulerEquations3D)
1787
# See Hughes, Franca, Mallet (1986) A new finite element formulation for CFD
1788
# [DOI: 10.1016/0045-7825(86)90127-1](https://doi.org/10.1016/0045-7825(86)90127-1)
1789
@unpack gamma = equations
1790
1791
# convert to entropy `-rho * s` used by Hughes, France, Mallet (1986)
1792
# instead of `-rho * s / (gamma - 1)`
1793
V1, V2, V3, V4, V5 = w .* (gamma - 1)
1794
1795
# s = specific entropy, eq. (53)
1796
V_square = V2^2 + V3^2 + V4^2
1797
s = gamma - V1 + V_square / (2 * V5)
1798
1799
# eq. (52)
1800
rho_iota = ((gamma - 1) / (-V5)^gamma)^(equations.inv_gamma_minus_one) *
1801
exp(-s * equations.inv_gamma_minus_one)
1802
1803
# eq. (51)
1804
rho = -rho_iota * V5
1805
rho_v1 = rho_iota * V2
1806
rho_v2 = rho_iota * V3
1807
rho_v3 = rho_iota * V4
1808
rho_e = rho_iota * (1 - V_square / (2 * V5))
1809
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e)
1810
end
1811
1812
# Convert primitive to conservative variables
1813
@inline function prim2cons(prim, equations::CompressibleEulerEquations3D)
1814
rho, v1, v2, v3, p = prim
1815
rho_v1 = rho * v1
1816
rho_v2 = rho * v2
1817
rho_v3 = rho * v3
1818
rho_e = p * equations.inv_gamma_minus_one +
1819
0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
1820
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e)
1821
end
1822
1823
@inline function density(u, equations::CompressibleEulerEquations3D)
1824
rho = u[1]
1825
return rho
1826
end
1827
1828
@inline function velocity(u, equations::CompressibleEulerEquations3D)
1829
rho = u[1]
1830
v1 = u[2] / rho
1831
v2 = u[3] / rho
1832
v3 = u[4] / rho
1833
return SVector(v1, v2, v3)
1834
end
1835
1836
@inline function velocity(u, orientation::Int, equations::CompressibleEulerEquations3D)
1837
rho = u[1]
1838
v = u[orientation + 1] / rho
1839
return v
1840
end
1841
1842
@inline function pressure(u, equations::CompressibleEulerEquations3D)
1843
rho, rho_v1, rho_v2, rho_v3, rho_e = u
1844
p = (equations.gamma - 1) * (rho_e - 0.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2) / rho)
1845
return p
1846
end
1847
1848
@inline function density_pressure(u, equations::CompressibleEulerEquations3D)
1849
rho, rho_v1, rho_v2, rho_v3, rho_e = u
1850
rho_times_p = (equations.gamma - 1) *
1851
(rho * rho_e - 0.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2))
1852
return rho_times_p
1853
end
1854
1855
# Calculate thermodynamic entropy for a conservative state `u`
1856
@inline function entropy_thermodynamic(u, equations::CompressibleEulerEquations3D)
1857
rho, _ = u
1858
p = pressure(u, equations)
1859
1860
# Thermodynamic entropy
1861
s = log(p) - equations.gamma * log(rho)
1862
1863
return s
1864
end
1865
1866
# Calculate mathematical entropy for a conservative state `cons`
1867
@inline function entropy_math(cons, equations::CompressibleEulerEquations3D)
1868
S = -entropy_thermodynamic(cons, equations) * cons[1] *
1869
equations.inv_gamma_minus_one
1870
# Mathematical entropy
1871
1872
return S
1873
end
1874
1875
# Default entropy is the mathematical entropy
1876
@inline function entropy(cons, equations::CompressibleEulerEquations3D)
1877
entropy_math(cons, equations)
1878
end
1879
1880
# Calculate total energy for a conservative state `cons`
1881
@inline energy_total(cons, ::CompressibleEulerEquations3D) = cons[5]
1882
1883
# Calculate kinetic energy for a conservative state `cons`
1884
@inline function energy_kinetic(u, equations::CompressibleEulerEquations3D)
1885
rho, rho_v1, rho_v2, rho_v3, _ = u
1886
return 0.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2) / rho
1887
end
1888
1889
# Calculate internal energy for a conservative state `cons`
1890
@inline function energy_internal(cons, equations::CompressibleEulerEquations3D)
1891
return energy_total(cons, equations) - energy_kinetic(cons, equations)
1892
end
1893
end # @muladd
1894
1895