Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/equations/compressible_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
CompressibleEulerEquations2D(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 e_{\text{total}}
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 e_{\text{total}} +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 e_{\text{total}} +p) v_2
26
\end{pmatrix}
27
=
28
\begin{pmatrix}
29
0 \\ 0 \\ 0 \\ 0
30
\end{pmatrix}
31
```
32
for an ideal gas with ratio of specific heats `gamma`
33
in two space dimensions.
34
Here, ``\rho`` is the density, ``v_1``, ``v_2`` the velocities, ``e_{\text{total}}`` the specific total energy, and
35
```math
36
p = (\gamma - 1) \left( \rho e_{\text{total}} - \frac{1}{2} \rho (v_1^2+v_2^2) \right)
37
```
38
the pressure.
39
"""
40
struct CompressibleEulerEquations2D{RealT <: Real} <:
41
AbstractCompressibleEulerEquations{2, 4}
42
gamma::RealT # ratio of specific heats
43
inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications
44
45
function CompressibleEulerEquations2D(gamma)
46
γ, inv_gamma_minus_one = promote(gamma, inv(gamma - 1))
47
return new{typeof(γ)}(γ, inv_gamma_minus_one)
48
end
49
end
50
51
function varnames(::typeof(cons2cons), ::CompressibleEulerEquations2D)
52
return ("rho", "rho_v1", "rho_v2", "rho_e_total")
53
end
54
varnames(::typeof(cons2prim), ::CompressibleEulerEquations2D) = ("rho", "v1", "v2", "p")
55
56
# Set initial conditions at physical location `x` for time `t`
57
"""
58
initial_condition_constant(x, t, equations::CompressibleEulerEquations2D)
59
60
A constant initial condition to test free-stream preservation.
61
"""
62
function initial_condition_constant(x, t, equations::CompressibleEulerEquations2D)
63
RealT = eltype(x)
64
rho = 1
65
rho_v1 = convert(RealT, 0.1)
66
rho_v2 = convert(RealT, -0.2)
67
rho_e_total = 10
68
return SVector(rho, rho_v1, rho_v2, rho_e_total)
69
end
70
71
"""
72
initial_condition_convergence_test(x, t, equations::CompressibleEulerEquations2D)
73
74
A smooth initial condition used for convergence tests in combination with
75
[`source_terms_convergence_test`](@ref)
76
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
77
"""
78
function initial_condition_convergence_test(x, t,
79
equations::CompressibleEulerEquations2D)
80
RealT = eltype(x)
81
c = 2
82
A = convert(RealT, 0.1)
83
L = 2
84
f = 1.0f0 / L
85
ω = 2 * convert(RealT, pi) * f
86
ini = c + A * sin(ω * (x[1] + x[2] - t))
87
88
rho = ini
89
rho_v1 = ini
90
rho_v2 = ini
91
rho_e_total = ini^2
92
93
return SVector(rho, rho_v1, rho_v2, rho_e_total)
94
end
95
96
"""
97
source_terms_convergence_test(u, x, t, equations::CompressibleEulerEquations2D)
98
99
Source terms used for convergence tests in combination with
100
[`initial_condition_convergence_test`](@ref)
101
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
102
103
References for the method of manufactured solutions (MMS):
104
- Kambiz Salari and Patrick Knupp (2000)
105
Code Verification by the Method of Manufactured Solutions
106
[DOI: 10.2172/759450](https://doi.org/10.2172/759450)
107
- Patrick J. Roache (2002)
108
Code Verification by the Method of Manufactured Solutions
109
[DOI: 10.1115/1.1436090](https://doi.org/10.1115/1.1436090)
110
"""
111
@inline function source_terms_convergence_test(u, x, t,
112
equations::CompressibleEulerEquations2D)
113
# Same settings as in `initial_condition`
114
RealT = eltype(u)
115
c = 2
116
A = convert(RealT, 0.1)
117
L = 2
118
f = 1.0f0 / L
119
ω = 2 * convert(RealT, pi) * f
120
γ = equations.gamma
121
122
x1, x2 = x
123
si, co = sincos(ω * (x1 + x2 - t))
124
rho = c + A * si
125
rho_x = ω * A * co
126
# Note that d/dt rho = -d/dx rho = -d/dy rho.
127
128
tmp = (2 * rho - 1) * (γ - 1)
129
130
du1 = rho_x
131
du2 = rho_x * (1 + tmp)
132
du3 = du2
133
du4 = 2 * rho_x * (rho + tmp)
134
135
return SVector(du1, du2, du3, du4)
136
end
137
138
"""
139
initial_condition_density_wave(x, t, equations::CompressibleEulerEquations2D)
140
141
A sine wave in the density with constant velocity and pressure; reduces the
142
compressible Euler equations to the linear advection equations.
143
This setup is the test case for stability of EC fluxes from paper
144
- Gregor J. Gassner, Magnus Svärd, Florian J. Hindenlang (2020)
145
Stability issues of entropy-stable and/or split-form high-order schemes
146
[arXiv: 2007.09026](https://arxiv.org/abs/2007.09026)
147
with the following parameters
148
- domain [-1, 1]
149
- mesh = 4x4
150
- polydeg = 5
151
"""
152
function initial_condition_density_wave(x, t, equations::CompressibleEulerEquations2D)
153
RealT = eltype(x)
154
v1 = convert(RealT, 0.1)
155
v2 = convert(RealT, 0.2)
156
rho = 1 + convert(RealT, 0.98) * sinpi(2 * (x[1] + x[2] - t * (v1 + v2)))
157
rho_v1 = rho * v1
158
rho_v2 = rho * v2
159
p = 20
160
rho_e_total = p / (equations.gamma - 1) + 0.5f0 * rho * (v1^2 + v2^2)
161
return SVector(rho, rho_v1, rho_v2, rho_e_total)
162
end
163
164
"""
165
initial_condition_weak_blast_wave(x, t, equations::CompressibleEulerEquations2D)
166
167
A weak blast wave taken from
168
- Sebastian Hennemann, Gregor J. Gassner (2020)
169
A provably entropy stable subcell shock capturing approach for high order split form DG
170
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
171
"""
172
function initial_condition_weak_blast_wave(x, t,
173
equations::CompressibleEulerEquations2D)
174
# From Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
175
# Set up polar coordinates
176
inicenter = SVector(0, 0)
177
x_norm = x[1] - inicenter[1]
178
y_norm = x[2] - inicenter[2]
179
r = sqrt(x_norm^2 + y_norm^2)
180
phi = atan(y_norm, x_norm)
181
sin_phi, cos_phi = sincos(phi)
182
183
# Calculate primitive variables
184
RealT = eltype(x)
185
rho = r > 0.5f0 ? one(RealT) : convert(RealT, 1.1691)
186
v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos_phi
187
v2 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * sin_phi
188
p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)
189
190
return prim2cons(SVector(rho, v1, v2, p), equations)
191
end
192
193
"""
194
initial_condition_eoc_test_coupled_euler_gravity(x, t, equations::CompressibleEulerEquations2D)
195
196
Setup used for convergence tests of the Euler equations with self-gravity used in
197
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
198
A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics
199
[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)
200
in combination with [`source_terms_eoc_test_coupled_euler_gravity`](@ref)
201
or [`source_terms_eoc_test_euler`](@ref).
202
"""
203
function initial_condition_eoc_test_coupled_euler_gravity(x, t,
204
equations::CompressibleEulerEquations2D)
205
# OBS! this assumes that γ = 2 other manufactured source terms are incorrect
206
if equations.gamma != 2
207
error("adiabatic constant must be 2 for the coupling convergence test")
208
end
209
RealT = eltype(x)
210
c = 2
211
A = convert(RealT, 0.1)
212
ini = c + A * sinpi(x[1] + x[2] - t)
213
G = 1 # gravitational constant
214
215
rho = ini
216
v1 = 1
217
v2 = 1
218
p = ini^2 * G / convert(RealT, pi) # * 2 / ndims, but ndims==2 here
219
220
return prim2cons(SVector(rho, v1, v2, p), equations)
221
end
222
223
"""
224
source_terms_eoc_test_coupled_euler_gravity(u, x, t, equations::CompressibleEulerEquations2D)
225
226
Setup used for convergence tests of the Euler equations with self-gravity used in
227
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
228
A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics
229
[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)
230
in combination with [`initial_condition_eoc_test_coupled_euler_gravity`](@ref).
231
"""
232
@inline function source_terms_eoc_test_coupled_euler_gravity(u, x, t,
233
equations::CompressibleEulerEquations2D)
234
# Same settings as in `initial_condition_eoc_test_coupled_euler_gravity`
235
RealT = eltype(u)
236
c = 2
237
A = convert(RealT, 0.1)
238
G = 1 # gravitational constant, must match coupling solver
239
C_grav = -2 * G / convert(RealT, pi) # 2 == 4 / ndims
240
241
x1, x2 = x
242
si, co = sincospi(x1 + x2 - t)
243
rhox = A * convert(RealT, pi) * co
244
rho = c + A * si
245
246
du1 = rhox
247
du2 = rhox
248
du3 = rhox
249
du4 = (1 - C_grav * rho) * rhox
250
251
return SVector(du1, du2, du3, du4)
252
end
253
254
"""
255
source_terms_eoc_test_euler(u, x, t, equations::CompressibleEulerEquations2D)
256
257
Setup used for convergence tests of the Euler equations with self-gravity used in
258
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
259
A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics
260
[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)
261
in combination with [`initial_condition_eoc_test_coupled_euler_gravity`](@ref).
262
"""
263
@inline function source_terms_eoc_test_euler(u, x, t,
264
equations::CompressibleEulerEquations2D)
265
# Same settings as in `initial_condition_eoc_test_coupled_euler_gravity`
266
RealT = eltype(u)
267
c = 2
268
A = convert(RealT, 0.1)
269
G = 1
270
C_grav = -2 * G / convert(RealT, pi) # 2 == 4 / ndims
271
272
x1, x2 = x
273
si, co = sincospi(x1 + x2 - t)
274
rhox = A * convert(RealT, pi) * co
275
rho = c + A * si
276
277
du1 = rhox
278
du2 = rhox * (1 - C_grav * rho)
279
du3 = rhox * (1 - C_grav * rho)
280
du4 = rhox * (1 - 3 * C_grav * rho)
281
282
return SVector(du1, du2, du3, du4)
283
end
284
285
"""
286
boundary_condition_slip_wall(u_inner, normal_direction, x, t, surface_flux_function,
287
equations::CompressibleEulerEquations2D)
288
289
Determine the boundary numerical surface flux for a slip wall condition.
290
Imposes a zero normal velocity at the wall.
291
Density is taken from the internal solution state and pressure is computed as an
292
exact solution of a 1D Riemann problem. Further details about this boundary state
293
are available in the paper:
294
- J. J. W. van der Vegt and H. van der Ven (2002)
295
Slip flow boundary conditions in discontinuous Galerkin discretizations of
296
the Euler equations of gas dynamics
297
[PDF](https://reports.nlr.nl/bitstream/handle/10921/692/TP-2002-300.pdf?sequence=1)
298
299
Details about the 1D pressure Riemann solution can be found in Section 6.3.3 of the book
300
- Eleuterio F. Toro (2009)
301
Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction
302
3rd edition
303
[DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)
304
305
Should be used together with [`UnstructuredMesh2D`](@ref), [`P4estMesh`](@ref), or [`T8codeMesh`](@ref).
306
"""
307
@inline function boundary_condition_slip_wall(u_inner, normal_direction::AbstractVector,
308
x, t,
309
surface_flux_function,
310
equations::CompressibleEulerEquations2D)
311
norm_ = norm(normal_direction)
312
# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later
313
normal = normal_direction / norm_
314
315
# rotate the internal solution state
316
u_local = rotate_to_x(u_inner, normal, equations)
317
318
# compute the primitive variables
319
rho_local, v_normal, v_tangent, p_local = cons2prim(u_local, equations)
320
321
# Get the solution of the pressure Riemann problem
322
# See Section 6.3.3 of
323
# Eleuterio F. Toro (2009)
324
# Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction
325
# [DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)
326
if v_normal <= 0
327
sound_speed = sqrt(equations.gamma * p_local / rho_local) # local sound speed
328
p_star = p_local *
329
(1 + 0.5f0 * (equations.gamma - 1) * v_normal / sound_speed)^(2 *
330
equations.gamma *
331
equations.inv_gamma_minus_one)
332
else # v_normal > 0
333
A = 2 / ((equations.gamma + 1) * rho_local)
334
B = p_local * (equations.gamma - 1) / (equations.gamma + 1)
335
p_star = p_local +
336
0.5f0 * v_normal / A *
337
(v_normal + sqrt(v_normal^2 + 4 * A * (p_local + B)))
338
end
339
340
# For the slip wall we directly set the flux as the normal velocity is zero
341
return SVector(0,
342
p_star * normal[1],
343
p_star * normal[2],
344
0) * norm_
345
end
346
347
"""
348
boundary_condition_slip_wall(u_inner, orientation, direction, x, t,
349
surface_flux_function, equations::CompressibleEulerEquations2D)
350
351
Should be used together with [`TreeMesh`](@ref).
352
"""
353
@inline function boundary_condition_slip_wall(u_inner, orientation,
354
direction, x, t,
355
surface_flux_function,
356
equations::CompressibleEulerEquations2D)
357
# get the appropriate normal vector from the orientation
358
RealT = eltype(u_inner)
359
if orientation == 1
360
normal_direction = SVector(one(RealT), zero(RealT))
361
else # orientation == 2
362
normal_direction = SVector(zero(RealT), one(RealT))
363
end
364
365
# compute and return the flux using `boundary_condition_slip_wall` routine below
366
return boundary_condition_slip_wall(u_inner, normal_direction, direction,
367
x, t, surface_flux_function, equations)
368
end
369
370
"""
371
boundary_condition_slip_wall(u_inner, normal_direction, direction, x, t,
372
surface_flux_function, equations::CompressibleEulerEquations2D)
373
374
Should be used together with [`StructuredMesh`](@ref).
375
"""
376
@inline function boundary_condition_slip_wall(u_inner, normal_direction::AbstractVector,
377
direction, x, t,
378
surface_flux_function,
379
equations::CompressibleEulerEquations2D)
380
# flip sign of normal to make it outward pointing, then flip the sign of the normal flux back
381
# to be inward pointing on the -x and -y sides due to the orientation convention used by StructuredMesh
382
if isodd(direction)
383
boundary_flux = -boundary_condition_slip_wall(u_inner, -normal_direction,
384
x, t, surface_flux_function,
385
equations)
386
else
387
boundary_flux = boundary_condition_slip_wall(u_inner, normal_direction,
388
x, t, surface_flux_function,
389
equations)
390
end
391
392
return boundary_flux
393
end
394
395
# Calculate 1D flux for a single point
396
@inline function flux(u, orientation::Integer, equations::CompressibleEulerEquations2D)
397
rho, rho_v1, rho_v2, rho_e_total = u
398
v1 = rho_v1 / rho
399
v2 = rho_v2 / rho
400
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
401
if orientation == 1
402
f1 = rho_v1
403
f2 = rho_v1 * v1 + p
404
f3 = rho_v1 * v2
405
f4 = (rho_e_total + p) * v1
406
else
407
f1 = rho_v2
408
f2 = rho_v2 * v1
409
f3 = rho_v2 * v2 + p
410
f4 = (rho_e_total + p) * v2
411
end
412
return SVector(f1, f2, f3, f4)
413
end
414
415
# Calculate 1D flux for a single point in the normal direction
416
# Note, this directional vector is not normalized
417
@inline function flux(u, normal_direction::AbstractVector,
418
equations::CompressibleEulerEquations2D)
419
rho_e_total = last(u)
420
rho, v1, v2, p = cons2prim(u, equations)
421
422
v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]
423
rho_v_normal = rho * v_normal
424
f1 = rho_v_normal
425
f2 = rho_v_normal * v1 + p * normal_direction[1]
426
f3 = rho_v_normal * v2 + p * normal_direction[2]
427
f4 = (rho_e_total + p) * v_normal
428
return SVector(f1, f2, f3, f4)
429
end
430
431
"""
432
flux_shima_etal(u_ll, u_rr, orientation_or_normal_direction,
433
equations::CompressibleEulerEquations2D)
434
435
This flux is is a modification of the original kinetic energy preserving two-point flux by
436
- Yuichi Kuya, Kosuke Totani and Soshi Kawai (2018)
437
Kinetic energy and entropy preserving schemes for compressible flows
438
by split convective forms
439
[DOI: 10.1016/j.jcp.2018.08.058](https://doi.org/10.1016/j.jcp.2018.08.058)
440
441
The modification is in the energy flux to guarantee pressure equilibrium and was developed by
442
- Nao Shima, Yuichi Kuya, Yoshiharu Tamaki, Soshi Kawai (JCP 2020)
443
Preventing spurious pressure oscillations in split convective form discretizations for
444
compressible flows
445
[DOI: 10.1016/j.jcp.2020.110060](https://doi.org/10.1016/j.jcp.2020.110060)
446
"""
447
@inline function flux_shima_etal(u_ll, u_rr, orientation::Integer,
448
equations::CompressibleEulerEquations2D)
449
# Unpack left and right state
450
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
451
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
452
453
# Average each factor of products in flux
454
rho_avg = 0.5f0 * (rho_ll + rho_rr)
455
v1_avg = 0.5f0 * (v1_ll + v1_rr)
456
v2_avg = 0.5f0 * (v2_ll + v2_rr)
457
p_avg = 0.5f0 * (p_ll + p_rr)
458
kin_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr)
459
460
# Calculate fluxes depending on orientation
461
if orientation == 1
462
pv1_avg = 0.5f0 * (p_ll * v1_rr + p_rr * v1_ll)
463
f1 = rho_avg * v1_avg
464
f2 = f1 * v1_avg + p_avg
465
f3 = f1 * v2_avg
466
f4 = p_avg * v1_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv1_avg
467
else
468
pv2_avg = 0.5f0 * (p_ll * v2_rr + p_rr * v2_ll)
469
f1 = rho_avg * v2_avg
470
f2 = f1 * v1_avg
471
f3 = f1 * v2_avg + p_avg
472
f4 = p_avg * v2_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv2_avg
473
end
474
475
return SVector(f1, f2, f3, f4)
476
end
477
478
@inline function flux_shima_etal(u_ll, u_rr, normal_direction::AbstractVector,
479
equations::CompressibleEulerEquations2D)
480
# Unpack left and right state
481
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
482
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
483
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
484
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
485
486
# Average each factor of products in flux
487
rho_avg = 0.5f0 * (rho_ll + rho_rr)
488
v1_avg = 0.5f0 * (v1_ll + v1_rr)
489
v2_avg = 0.5f0 * (v2_ll + v2_rr)
490
v_dot_n_avg = 0.5f0 * (v_dot_n_ll + v_dot_n_rr)
491
p_avg = 0.5f0 * (p_ll + p_rr)
492
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr)
493
494
# Calculate fluxes depending on normal_direction
495
f1 = rho_avg * v_dot_n_avg
496
f2 = f1 * v1_avg + p_avg * normal_direction[1]
497
f3 = f1 * v2_avg + p_avg * normal_direction[2]
498
f4 = (f1 * velocity_square_avg +
499
p_avg * v_dot_n_avg * equations.inv_gamma_minus_one
500
+ 0.5f0 * (p_ll * v_dot_n_rr + p_rr * v_dot_n_ll))
501
502
return SVector(f1, f2, f3, f4)
503
end
504
505
"""
506
flux_kennedy_gruber(u_ll, u_rr, orientation_or_normal_direction,
507
equations::CompressibleEulerEquations2D)
508
509
Kinetic energy preserving two-point flux by
510
- Kennedy and Gruber (2008)
511
Reduced aliasing formulations of the convective terms within the
512
Navier-Stokes equations for a compressible fluid
513
[DOI: 10.1016/j.jcp.2007.09.020](https://doi.org/10.1016/j.jcp.2007.09.020)
514
"""
515
@inline function flux_kennedy_gruber(u_ll, u_rr, orientation::Integer,
516
equations::CompressibleEulerEquations2D)
517
# Unpack left and right state
518
rho_e_total_ll = last(u_ll)
519
rho_e_total_rr = last(u_rr)
520
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
521
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
522
523
# Average each factor of products in flux
524
rho_avg = 0.5f0 * (rho_ll + rho_rr)
525
v1_avg = 0.5f0 * (v1_ll + v1_rr)
526
v2_avg = 0.5f0 * (v2_ll + v2_rr)
527
p_avg = 0.5f0 * (p_ll + p_rr)
528
e_avg = 0.5f0 * (rho_e_total_ll / rho_ll + rho_e_total_rr / rho_rr)
529
530
# Calculate fluxes depending on orientation
531
if orientation == 1
532
f1 = rho_avg * v1_avg
533
f2 = rho_avg * v1_avg * v1_avg + p_avg
534
f3 = rho_avg * v1_avg * v2_avg
535
f4 = (rho_avg * e_avg + p_avg) * v1_avg
536
else
537
f1 = rho_avg * v2_avg
538
f2 = rho_avg * v2_avg * v1_avg
539
f3 = rho_avg * v2_avg * v2_avg + p_avg
540
f4 = (rho_avg * e_avg + p_avg) * v2_avg
541
end
542
543
return SVector(f1, f2, f3, f4)
544
end
545
546
@inline function flux_kennedy_gruber(u_ll, u_rr, normal_direction::AbstractVector,
547
equations::CompressibleEulerEquations2D)
548
# Unpack left and right state
549
rho_e_total_ll = last(u_ll)
550
rho_e_total_rr = last(u_rr)
551
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
552
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
553
554
# Average each factor of products in flux
555
rho_avg = 0.5f0 * (rho_ll + rho_rr)
556
v1_avg = 0.5f0 * (v1_ll + v1_rr)
557
v2_avg = 0.5f0 * (v2_ll + v2_rr)
558
v_dot_n_avg = v1_avg * normal_direction[1] + v2_avg * normal_direction[2]
559
p_avg = 0.5f0 * (p_ll + p_rr)
560
e_avg = 0.5f0 * (rho_e_total_ll / rho_ll + rho_e_total_rr / rho_rr)
561
562
# Calculate fluxes depending on normal_direction
563
f1 = rho_avg * v_dot_n_avg
564
f2 = f1 * v1_avg + p_avg * normal_direction[1]
565
f3 = f1 * v2_avg + p_avg * normal_direction[2]
566
f4 = f1 * e_avg + p_avg * v_dot_n_avg
567
568
return SVector(f1, f2, f3, f4)
569
end
570
571
"""
572
flux_chandrashekar(u_ll, u_rr, orientation_or_normal_direction, equations::CompressibleEulerEquations2D)
573
574
Entropy conserving two-point flux by
575
- Chandrashekar (2013)
576
Kinetic Energy Preserving and Entropy Stable Finite Volume Schemes
577
for Compressible Euler and Navier-Stokes Equations
578
[DOI: 10.4208/cicp.170712.010313a](https://doi.org/10.4208/cicp.170712.010313a)
579
"""
580
@inline function flux_chandrashekar(u_ll, u_rr, orientation::Integer,
581
equations::CompressibleEulerEquations2D)
582
# Unpack left and right state
583
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
584
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
585
beta_ll = 0.5f0 * rho_ll / p_ll
586
beta_rr = 0.5f0 * rho_rr / p_rr
587
specific_kin_ll = 0.5f0 * (v1_ll^2 + v2_ll^2)
588
specific_kin_rr = 0.5f0 * (v1_rr^2 + v2_rr^2)
589
590
# Compute the necessary mean values
591
rho_avg = 0.5f0 * (rho_ll + rho_rr)
592
rho_mean = ln_mean(rho_ll, rho_rr)
593
beta_mean = ln_mean(beta_ll, beta_rr)
594
beta_avg = 0.5f0 * (beta_ll + beta_rr)
595
v1_avg = 0.5f0 * (v1_ll + v1_rr)
596
v2_avg = 0.5f0 * (v2_ll + v2_rr)
597
p_mean = 0.5f0 * rho_avg / beta_avg
598
velocity_square_avg = specific_kin_ll + specific_kin_rr
599
600
# Calculate fluxes depending on orientation
601
if orientation == 1
602
f1 = rho_mean * v1_avg
603
f2 = f1 * v1_avg + p_mean
604
f3 = f1 * v2_avg
605
f4 = f1 * 0.5f0 *
606
(1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +
607
f2 * v1_avg + f3 * v2_avg
608
else
609
f1 = rho_mean * v2_avg
610
f2 = f1 * v1_avg
611
f3 = f1 * v2_avg + p_mean
612
f4 = f1 * 0.5f0 *
613
(1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +
614
f2 * v1_avg + f3 * v2_avg
615
end
616
617
return SVector(f1, f2, f3, f4)
618
end
619
620
@inline function flux_chandrashekar(u_ll, u_rr, normal_direction::AbstractVector,
621
equations::CompressibleEulerEquations2D)
622
# Unpack left and right state
623
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
624
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
625
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
626
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
627
beta_ll = 0.5f0 * rho_ll / p_ll
628
beta_rr = 0.5f0 * rho_rr / p_rr
629
specific_kin_ll = 0.5f0 * (v1_ll^2 + v2_ll^2)
630
specific_kin_rr = 0.5f0 * (v1_rr^2 + v2_rr^2)
631
632
# Compute the necessary mean values
633
rho_avg = 0.5f0 * (rho_ll + rho_rr)
634
rho_mean = ln_mean(rho_ll, rho_rr)
635
beta_mean = ln_mean(beta_ll, beta_rr)
636
beta_avg = 0.5f0 * (beta_ll + beta_rr)
637
v1_avg = 0.5f0 * (v1_ll + v1_rr)
638
v2_avg = 0.5f0 * (v2_ll + v2_rr)
639
p_mean = 0.5f0 * rho_avg / beta_avg
640
velocity_square_avg = specific_kin_ll + specific_kin_rr
641
642
# Multiply with average of normal velocities
643
f1 = rho_mean * 0.5f0 * (v_dot_n_ll + v_dot_n_rr)
644
f2 = f1 * v1_avg + p_mean * normal_direction[1]
645
f3 = f1 * v2_avg + p_mean * normal_direction[2]
646
f4 = f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +
647
f2 * v1_avg + f3 * v2_avg
648
649
return SVector(f1, f2, f3, f4)
650
end
651
652
"""
653
flux_ranocha(u_ll, u_rr, orientation_or_normal_direction,
654
equations::CompressibleEulerEquations2D)
655
656
Entropy conserving and kinetic energy preserving two-point flux by
657
- Hendrik Ranocha (2018)
658
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
659
for Hyperbolic Balance Laws
660
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
661
See also
662
- Hendrik Ranocha (2020)
663
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
664
the Euler Equations Using Summation-by-Parts Operators
665
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
666
"""
667
@inline function flux_ranocha(u_ll, u_rr, orientation::Integer,
668
equations::CompressibleEulerEquations2D)
669
# Unpack left and right state
670
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
671
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
672
673
# Compute the necessary mean values
674
rho_mean = ln_mean(rho_ll, rho_rr)
675
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
676
# in exact arithmetic since
677
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
678
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
679
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
680
v1_avg = 0.5f0 * (v1_ll + v1_rr)
681
v2_avg = 0.5f0 * (v2_ll + v2_rr)
682
p_avg = 0.5f0 * (p_ll + p_rr)
683
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr)
684
685
# Calculate fluxes depending on orientation
686
if orientation == 1
687
f1 = rho_mean * v1_avg
688
f2 = f1 * v1_avg + p_avg
689
f3 = f1 * v2_avg
690
f4 = f1 *
691
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +
692
0.5f0 * (p_ll * v1_rr + p_rr * v1_ll)
693
else
694
f1 = rho_mean * v2_avg
695
f2 = f1 * v1_avg
696
f3 = f1 * v2_avg + p_avg
697
f4 = f1 *
698
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +
699
0.5f0 * (p_ll * v2_rr + p_rr * v2_ll)
700
end
701
702
return SVector(f1, f2, f3, f4)
703
end
704
705
@inline function flux_ranocha(u_ll, u_rr, normal_direction::AbstractVector,
706
equations::CompressibleEulerEquations2D)
707
# Unpack left and right state
708
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
709
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
710
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
711
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
712
713
# Compute the necessary mean values
714
rho_mean = ln_mean(rho_ll, rho_rr)
715
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
716
# in exact arithmetic since
717
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
718
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
719
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
720
v1_avg = 0.5f0 * (v1_ll + v1_rr)
721
v2_avg = 0.5f0 * (v2_ll + v2_rr)
722
p_avg = 0.5f0 * (p_ll + p_rr)
723
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr)
724
725
# Calculate fluxes depending on normal_direction
726
f1 = rho_mean * 0.5f0 * (v_dot_n_ll + v_dot_n_rr)
727
f2 = f1 * v1_avg + p_avg * normal_direction[1]
728
f3 = f1 * v2_avg + p_avg * normal_direction[2]
729
f4 = (f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
730
+
731
0.5f0 * (p_ll * v_dot_n_rr + p_rr * v_dot_n_ll))
732
733
return SVector(f1, f2, f3, f4)
734
end
735
736
"""
737
splitting_steger_warming(u, orientation::Integer,
738
equations::CompressibleEulerEquations2D)
739
splitting_steger_warming(u, which::Union{Val{:minus}, Val{:plus}}
740
orientation::Integer,
741
equations::CompressibleEulerEquations2D)
742
743
Splitting of the compressible Euler flux of Steger and Warming. For
744
curvilinear coordinates use the improved Steger-Warming-type splitting
745
[`splitting_drikakis_tsangaris`](@ref).
746
747
Returns a tuple of the fluxes "minus" (associated with waves going into the
748
negative axis direction) and "plus" (associated with waves going into the
749
positive axis direction). If only one of the fluxes is required, use the
750
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
751
752
!!! warning "Experimental implementation (upwind SBP)"
753
This is an experimental feature and may change in future releases.
754
755
## References
756
757
- Joseph L. Steger and R. F. Warming (1979)
758
Flux Vector Splitting of the Inviscid Gasdynamic Equations
759
With Application to Finite Difference Methods
760
[NASA Technical Memorandum](https://ntrs.nasa.gov/api/citations/19790020779/downloads/19790020779.pdf)
761
"""
762
@inline function splitting_steger_warming(u, orientation::Integer,
763
equations::CompressibleEulerEquations2D)
764
fm = splitting_steger_warming(u, Val{:minus}(), orientation, equations)
765
fp = splitting_steger_warming(u, Val{:plus}(), orientation, equations)
766
return fm, fp
767
end
768
769
@inline function splitting_steger_warming(u, ::Val{:plus}, orientation::Integer,
770
equations::CompressibleEulerEquations2D)
771
rho, rho_v1, rho_v2, rho_e_total = u
772
v1 = rho_v1 / rho
773
v2 = rho_v2 / rho
774
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
775
a = sqrt(equations.gamma * p / rho)
776
777
if orientation == 1
778
lambda1 = v1
779
lambda2 = v1 + a
780
lambda3 = v1 - a
781
782
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
783
lambda2_p = positive_part(lambda2)
784
lambda3_p = positive_part(lambda3)
785
786
alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p
787
788
rho_2gamma = 0.5f0 * rho / equations.gamma
789
f1p = rho_2gamma * alpha_p
790
f2p = rho_2gamma * (alpha_p * v1 + a * (lambda2_p - lambda3_p))
791
f3p = rho_2gamma * alpha_p * v2
792
f4p = rho_2gamma *
793
(alpha_p * 0.5f0 * (v1^2 + v2^2) + a * v1 * (lambda2_p - lambda3_p)
794
+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)
795
else # orientation == 2
796
lambda1 = v2
797
lambda2 = v2 + a
798
lambda3 = v2 - a
799
800
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
801
lambda2_p = positive_part(lambda2)
802
lambda3_p = positive_part(lambda3)
803
804
alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p
805
806
rho_2gamma = 0.5f0 * rho / equations.gamma
807
f1p = rho_2gamma * alpha_p
808
f2p = rho_2gamma * alpha_p * v1
809
f3p = rho_2gamma * (alpha_p * v2 + a * (lambda2_p - lambda3_p))
810
f4p = rho_2gamma *
811
(alpha_p * 0.5f0 * (v1^2 + v2^2) + a * v2 * (lambda2_p - lambda3_p)
812
+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)
813
end
814
return SVector(f1p, f2p, f3p, f4p)
815
end
816
817
@inline function splitting_steger_warming(u, ::Val{:minus}, orientation::Integer,
818
equations::CompressibleEulerEquations2D)
819
rho, rho_v1, rho_v2, rho_e_total = u
820
v1 = rho_v1 / rho
821
v2 = rho_v2 / rho
822
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
823
a = sqrt(equations.gamma * p / rho)
824
825
if orientation == 1
826
lambda1 = v1
827
lambda2 = v1 + a
828
lambda3 = v1 - a
829
830
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
831
lambda2_m = negative_part(lambda2)
832
lambda3_m = negative_part(lambda3)
833
834
alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m
835
836
rho_2gamma = 0.5f0 * rho / equations.gamma
837
f1m = rho_2gamma * alpha_m
838
f2m = rho_2gamma * (alpha_m * v1 + a * (lambda2_m - lambda3_m))
839
f3m = rho_2gamma * alpha_m * v2
840
f4m = rho_2gamma *
841
(alpha_m * 0.5f0 * (v1^2 + v2^2) + a * v1 * (lambda2_m - lambda3_m)
842
+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)
843
else # orientation == 2
844
lambda1 = v2
845
lambda2 = v2 + a
846
lambda3 = v2 - a
847
848
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
849
lambda2_m = negative_part(lambda2)
850
lambda3_m = negative_part(lambda3)
851
852
alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m
853
854
rho_2gamma = 0.5f0 * rho / equations.gamma
855
f1m = rho_2gamma * alpha_m
856
f2m = rho_2gamma * alpha_m * v1
857
f3m = rho_2gamma * (alpha_m * v2 + a * (lambda2_m - lambda3_m))
858
f4m = rho_2gamma *
859
(alpha_m * 0.5f0 * (v1^2 + v2^2) + a * v2 * (lambda2_m - lambda3_m)
860
+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)
861
end
862
return SVector(f1m, f2m, f3m, f4m)
863
end
864
865
"""
866
splitting_drikakis_tsangaris(u, orientation_or_normal_direction,
867
equations::CompressibleEulerEquations2D)
868
splitting_drikakis_tsangaris(u, which::Union{Val{:minus}, Val{:plus}}
869
orientation_or_normal_direction,
870
equations::CompressibleEulerEquations2D)
871
872
Improved variant of the Steger-Warming flux vector splitting
873
[`splitting_steger_warming`](@ref) for generalized coordinates.
874
This splitting also reformulates the energy
875
flux as in Hänel et al. to obtain conservation of the total temperature
876
for inviscid flows.
877
878
Returns a tuple of the fluxes "minus" (associated with waves going into the
879
negative axis direction) and "plus" (associated with waves going into the
880
positive axis direction). If only one of the fluxes is required, use the
881
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
882
883
!!! warning "Experimental implementation (upwind SBP)"
884
This is an experimental feature and may change in future releases.
885
886
## References
887
888
- D. Drikakis and S. Tsangaris (1993)
889
On the solution of the compressible Navier-Stokes equations using
890
improved flux vector splitting methods
891
[DOI: 10.1016/0307-904X(93)90054-K](https://doi.org/10.1016/0307-904X(93)90054-K)
892
- D. Hänel, R. Schwane and G. Seider (1987)
893
On the accuracy of upwind schemes for the solution of the Navier-Stokes equations
894
[DOI: 10.2514/6.1987-1105](https://doi.org/10.2514/6.1987-1105)
895
"""
896
@inline function splitting_drikakis_tsangaris(u, orientation_or_normal_direction,
897
equations::CompressibleEulerEquations2D)
898
fm = splitting_drikakis_tsangaris(u, Val{:minus}(), orientation_or_normal_direction,
899
equations)
900
fp = splitting_drikakis_tsangaris(u, Val{:plus}(), orientation_or_normal_direction,
901
equations)
902
return fm, fp
903
end
904
905
@inline function splitting_drikakis_tsangaris(u, ::Val{:plus}, orientation::Integer,
906
equations::CompressibleEulerEquations2D)
907
rho, rho_v1, rho_v2, rho_e_total = u
908
v1 = rho_v1 / rho
909
v2 = rho_v2 / rho
910
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
911
a = sqrt(equations.gamma * p / rho)
912
H = (rho_e_total + p) / rho
913
914
if orientation == 1
915
lambda1 = v1 + a
916
lambda2 = v1 - a
917
918
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
919
lambda2_p = positive_part(lambda2)
920
921
rhoa_2gamma = 0.5f0 * rho * a / equations.gamma
922
f1p = 0.5f0 * rho * (lambda1_p + lambda2_p)
923
f2p = f1p * v1 + rhoa_2gamma * (lambda1_p - lambda2_p)
924
f3p = f1p * v2
925
f4p = f1p * H
926
else # orientation == 2
927
lambda1 = v2 + a
928
lambda2 = v2 - a
929
930
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
931
lambda2_p = positive_part(lambda2)
932
933
rhoa_2gamma = 0.5f0 * rho * a / equations.gamma
934
f1p = 0.5f0 * rho * (lambda1_p + lambda2_p)
935
f2p = f1p * v1
936
f3p = f1p * v2 + rhoa_2gamma * (lambda1_p - lambda2_p)
937
f4p = f1p * H
938
end
939
return SVector(f1p, f2p, f3p, f4p)
940
end
941
942
@inline function splitting_drikakis_tsangaris(u, ::Val{:minus}, orientation::Integer,
943
equations::CompressibleEulerEquations2D)
944
rho, rho_v1, rho_v2, rho_e_total = u
945
v1 = rho_v1 / rho
946
v2 = rho_v2 / rho
947
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
948
a = sqrt(equations.gamma * p / rho)
949
H = (rho_e_total + p) / rho
950
951
if orientation == 1
952
lambda1 = v1 + a
953
lambda2 = v1 - a
954
955
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
956
lambda2_m = negative_part(lambda2)
957
958
rhoa_2gamma = 0.5f0 * rho * a / equations.gamma
959
f1m = 0.5f0 * rho * (lambda1_m + lambda2_m)
960
f2m = f1m * v1 + rhoa_2gamma * (lambda1_m - lambda2_m)
961
f3m = f1m * v2
962
f4m = f1m * H
963
else # orientation == 2
964
lambda1 = v2 + a
965
lambda2 = v2 - a
966
967
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
968
lambda2_m = negative_part(lambda2)
969
970
rhoa_2gamma = 0.5f0 * rho * a / equations.gamma
971
f1m = 0.5f0 * rho * (lambda1_m + lambda2_m)
972
f2m = f1m * v1
973
f3m = f1m * v2 + rhoa_2gamma * (lambda1_m - lambda2_m)
974
f4m = f1m * H
975
end
976
return SVector(f1m, f2m, f3m, f4m)
977
end
978
979
@inline function splitting_drikakis_tsangaris(u, ::Val{:plus},
980
normal_direction::AbstractVector,
981
equations::CompressibleEulerEquations2D)
982
rho, rho_v1, rho_v2, rho_e_total = u
983
v1 = rho_v1 / rho
984
v2 = rho_v2 / rho
985
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
986
a = sqrt(equations.gamma * p / rho)
987
H = (rho_e_total + p) / rho
988
989
v_n = normal_direction[1] * v1 + normal_direction[2] * v2
990
991
lambda1 = v_n + a
992
lambda2 = v_n - a
993
994
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
995
lambda2_p = positive_part(lambda2)
996
997
rhoa_2gamma = 0.5f0 * rho * a / equations.gamma
998
f1p = 0.5f0 * rho * (lambda1_p + lambda2_p)
999
f2p = f1p * v1 + rhoa_2gamma * normal_direction[1] * (lambda1_p - lambda2_p)
1000
f3p = f1p * v2 + rhoa_2gamma * normal_direction[2] * (lambda1_p - lambda2_p)
1001
f4p = f1p * H
1002
1003
return SVector(f1p, f2p, f3p, f4p)
1004
end
1005
1006
@inline function splitting_drikakis_tsangaris(u, ::Val{:minus},
1007
normal_direction::AbstractVector,
1008
equations::CompressibleEulerEquations2D)
1009
rho, rho_v1, rho_v2, rho_e_total = u
1010
v1 = rho_v1 / rho
1011
v2 = rho_v2 / rho
1012
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
1013
a = sqrt(equations.gamma * p / rho)
1014
H = (rho_e_total + p) / rho
1015
1016
v_n = normal_direction[1] * v1 + normal_direction[2] * v2
1017
1018
lambda1 = v_n + a
1019
lambda2 = v_n - a
1020
1021
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
1022
lambda2_m = negative_part(lambda2)
1023
1024
rhoa_2gamma = 0.5f0 * rho * a / equations.gamma
1025
f1m = 0.5f0 * rho * (lambda1_m + lambda2_m)
1026
f2m = f1m * v1 + rhoa_2gamma * normal_direction[1] * (lambda1_m - lambda2_m)
1027
f3m = f1m * v2 + rhoa_2gamma * normal_direction[2] * (lambda1_m - lambda2_m)
1028
f4m = f1m * H
1029
1030
return SVector(f1m, f2m, f3m, f4m)
1031
end
1032
1033
"""
1034
FluxLMARS(c)(u_ll, u_rr, orientation_or_normal_direction,
1035
equations::CompressibleEulerEquations2D)
1036
1037
Low Mach number approximate Riemann solver (LMARS) for atmospheric flows using
1038
an estimate `c` of the speed of sound.
1039
1040
References:
1041
- Xi Chen et al. (2013)
1042
A Control-Volume Model of the Compressible Euler Equations with a Vertical
1043
Lagrangian Coordinate
1044
[DOI: 10.1175/MWR-D-12-00129.1](https://doi.org/10.1175/mwr-d-12-00129.1)
1045
"""
1046
struct FluxLMARS{SpeedOfSound}
1047
# Estimate for the speed of sound
1048
speed_of_sound::SpeedOfSound
1049
end
1050
1051
@inline function (flux_lmars::FluxLMARS)(u_ll, u_rr, orientation::Integer,
1052
equations::CompressibleEulerEquations2D)
1053
c = flux_lmars.speed_of_sound
1054
1055
# Unpack left and right state
1056
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1057
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1058
1059
if orientation == 1
1060
v_ll = v1_ll
1061
v_rr = v1_rr
1062
else # orientation == 2
1063
v_ll = v2_ll
1064
v_rr = v2_rr
1065
end
1066
1067
rho = 0.5f0 * (rho_ll + rho_rr)
1068
p = 0.5f0 * (p_ll + p_rr) - 0.5f0 * c * rho * (v_rr - v_ll)
1069
v = 0.5f0 * (v_ll + v_rr) - 1 / (2 * c * rho) * (p_rr - p_ll)
1070
1071
# We treat the energy term analogous to the potential temperature term in the paper by
1072
# Chen et al., i.e. we use p_ll and p_rr, and not p
1073
if v >= 0
1074
f1, f2, f3, f4 = v * u_ll
1075
f4 = f4 + p_ll * v
1076
else
1077
f1, f2, f3, f4 = v * u_rr
1078
f4 = f4 + p_rr * v
1079
end
1080
1081
if orientation == 1
1082
f2 = f2 + p
1083
else # orientation == 2
1084
f3 = f3 + p
1085
end
1086
1087
return SVector(f1, f2, f3, f4)
1088
end
1089
1090
@inline function (flux_lmars::FluxLMARS)(u_ll, u_rr, normal_direction::AbstractVector,
1091
equations::CompressibleEulerEquations2D)
1092
c = flux_lmars.speed_of_sound
1093
1094
# Unpack left and right state
1095
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1096
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1097
1098
v_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
1099
v_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
1100
1101
# Note that this is the same as computing v_ll and v_rr with a normalized normal vector
1102
# and then multiplying v by `norm_` again, but this version is slightly faster.
1103
norm_ = norm(normal_direction)
1104
1105
rho = 0.5f0 * (rho_ll + rho_rr)
1106
p = 0.5f0 * (p_ll + p_rr) - 0.5f0 * c * rho * (v_rr - v_ll) / norm_
1107
v = 0.5f0 * (v_ll + v_rr) - 1 / (2 * c * rho) * (p_rr - p_ll) * norm_
1108
1109
# We treat the energy term analogous to the potential temperature term in the paper by
1110
# Chen et al., i.e. we use p_ll and p_rr, and not p
1111
if v >= 0
1112
f1, f2, f3, f4 = u_ll * v
1113
f4 = f4 + p_ll * v
1114
else
1115
f1, f2, f3, f4 = u_rr * v
1116
f4 = f4 + p_rr * v
1117
end
1118
1119
return SVector(f1,
1120
f2 + p * normal_direction[1],
1121
f3 + p * normal_direction[2],
1122
f4)
1123
end
1124
1125
"""
1126
splitting_vanleer_haenel(u, orientation_or_normal_direction,
1127
equations::CompressibleEulerEquations2D)
1128
splitting_vanleer_haenel(u, which::Union{Val{:minus}, Val{:plus}}
1129
orientation_or_normal_direction,
1130
equations::CompressibleEulerEquations2D)
1131
1132
Splitting of the compressible Euler flux from van Leer. This splitting further
1133
contains a reformulation due to Hänel et al. where the energy flux uses the
1134
enthalpy. The pressure splitting is independent from the splitting of the
1135
convective terms. As such there are many pressure splittings suggested across
1136
the literature. We implement the 'p4' variant suggested by Liou and Steffen as
1137
it proved the most robust in practice. For details on the curvilinear variant
1138
of this flux vector splitting see Anderson et al.
1139
1140
Returns a tuple of the fluxes "minus" (associated with waves going into the
1141
negative axis direction) and "plus" (associated with waves going into the
1142
positive axis direction). If only one of the fluxes is required, use the
1143
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
1144
1145
!!! warning "Experimental implementation (upwind SBP)"
1146
This is an experimental feature and may change in future releases.
1147
1148
## References
1149
1150
- Bram van Leer (1982)
1151
Flux-Vector Splitting for the Euler Equation
1152
[DOI: 10.1007/978-3-642-60543-7_5](https://doi.org/10.1007/978-3-642-60543-7_5)
1153
- D. Hänel, R. Schwane and G. Seider (1987)
1154
On the accuracy of upwind schemes for the solution of the Navier-Stokes equations
1155
[DOI: 10.2514/6.1987-1105](https://doi.org/10.2514/6.1987-1105)
1156
- Meng-Sing Liou and Chris J. Steffen, Jr. (1991)
1157
High-Order Polynomial Expansions (HOPE) for Flux-Vector Splitting
1158
[NASA Technical Memorandum](https://ntrs.nasa.gov/citations/19910016425)
1159
- W. Kyle Anderson, James L. Thomas, and Bram van Leer (1986)
1160
Comparison of Finite Volume Flux Vector Splittings for the Euler Equations
1161
[DOI: 10.2514/3.9465](https://doi.org/10.2514/3.9465)
1162
"""
1163
@inline function splitting_vanleer_haenel(u, orientation_or_normal_direction,
1164
equations::CompressibleEulerEquations2D)
1165
fm = splitting_vanleer_haenel(u, Val{:minus}(), orientation_or_normal_direction,
1166
equations)
1167
fp = splitting_vanleer_haenel(u, Val{:plus}(), orientation_or_normal_direction,
1168
equations)
1169
return fm, fp
1170
end
1171
1172
@inline function splitting_vanleer_haenel(u, ::Val{:plus}, orientation::Integer,
1173
equations::CompressibleEulerEquations2D)
1174
rho, rho_v1, rho_v2, rho_e_total = u
1175
v1 = rho_v1 / rho
1176
v2 = rho_v2 / rho
1177
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
1178
1179
a = sqrt(equations.gamma * p / rho)
1180
H = (rho_e_total + p) / rho
1181
1182
if orientation == 1
1183
M = v1 / a
1184
p_plus = 0.5f0 * (1 + equations.gamma * M) * p
1185
1186
f1p = 0.25f0 * rho * a * (M + 1)^2
1187
f2p = f1p * v1 + p_plus
1188
f3p = f1p * v2
1189
f4p = f1p * H
1190
else # orientation == 2
1191
M = v2 / a
1192
p_plus = 0.5f0 * (1 + equations.gamma * M) * p
1193
1194
f1p = 0.25f0 * rho * a * (M + 1)^2
1195
f2p = f1p * v1
1196
f3p = f1p * v2 + p_plus
1197
f4p = f1p * H
1198
end
1199
return SVector(f1p, f2p, f3p, f4p)
1200
end
1201
1202
@inline function splitting_vanleer_haenel(u, ::Val{:minus}, orientation::Integer,
1203
equations::CompressibleEulerEquations2D)
1204
rho, rho_v1, rho_v2, rho_e_total = u
1205
v1 = rho_v1 / rho
1206
v2 = rho_v2 / rho
1207
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
1208
1209
a = sqrt(equations.gamma * p / rho)
1210
H = (rho_e_total + p) / rho
1211
1212
if orientation == 1
1213
M = v1 / a
1214
p_minus = 0.5f0 * (1 - equations.gamma * M) * p
1215
1216
f1m = -0.25f0 * rho * a * (M - 1)^2
1217
f2m = f1m * v1 + p_minus
1218
f3m = f1m * v2
1219
f4m = f1m * H
1220
else # orientation == 2
1221
M = v2 / a
1222
p_minus = 0.5f0 * (1 - equations.gamma * M) * p
1223
1224
f1m = -0.25f0 * rho * a * (M - 1)^2
1225
f2m = f1m * v1
1226
f3m = f1m * v2 + p_minus
1227
f4m = f1m * H
1228
end
1229
return SVector(f1m, f2m, f3m, f4m)
1230
end
1231
1232
@inline function splitting_vanleer_haenel(u, ::Val{:plus},
1233
normal_direction::AbstractVector,
1234
equations::CompressibleEulerEquations2D)
1235
rho, rho_v1, rho_v2, rho_e_total = u
1236
v1 = rho_v1 / rho
1237
v2 = rho_v2 / rho
1238
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
1239
1240
a = sqrt(equations.gamma * p / rho)
1241
H = (rho_e_total + p) / rho
1242
1243
v_n = normal_direction[1] * v1 + normal_direction[2] * v2
1244
M = v_n / a
1245
p_plus = 0.5f0 * (1 + equations.gamma * M) * p
1246
1247
f1p = 0.25f0 * rho * a * (M + 1)^2
1248
f2p = f1p * v1 + normal_direction[1] * p_plus
1249
f3p = f1p * v2 + normal_direction[2] * p_plus
1250
f4p = f1p * H
1251
1252
return SVector(f1p, f2p, f3p, f4p)
1253
end
1254
1255
@inline function splitting_vanleer_haenel(u, ::Val{:minus},
1256
normal_direction::AbstractVector,
1257
equations::CompressibleEulerEquations2D)
1258
rho, rho_v1, rho_v2, rho_e_total = u
1259
v1 = rho_v1 / rho
1260
v2 = rho_v2 / rho
1261
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
1262
1263
a = sqrt(equations.gamma * p / rho)
1264
H = (rho_e_total + p) / rho
1265
1266
v_n = normal_direction[1] * v1 + normal_direction[2] * v2
1267
M = v_n / a
1268
p_minus = 0.5f0 * (1 - equations.gamma * M) * p
1269
1270
f1m = -0.25f0 * rho * a * (M - 1)^2
1271
f2m = f1m * v1 + normal_direction[1] * p_minus
1272
f3m = f1m * v2 + normal_direction[2] * p_minus
1273
f4m = f1m * H
1274
1275
return SVector(f1m, f2m, f3m, f4m)
1276
end
1277
1278
"""
1279
splitting_lax_friedrichs(u, orientation_or_normal_direction,
1280
equations::CompressibleEulerEquations2D)
1281
splitting_lax_friedrichs(u, which::Union{Val{:minus}, Val{:plus}}
1282
orientation_or_normal_direction,
1283
equations::CompressibleEulerEquations2D)
1284
1285
Naive local Lax-Friedrichs style flux splitting of the form `f⁺ = 0.5 (f + λ u)`
1286
and `f⁻ = 0.5 (f - λ u)` similar to a flux splitting one would apply, e.g.,
1287
to Burgers' equation.
1288
1289
Returns a tuple of the fluxes "minus" (associated with waves going into the
1290
negative axis direction) and "plus" (associated with waves going into the
1291
positive axis direction). If only one of the fluxes is required, use the
1292
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
1293
1294
!!! warning "Experimental implementation (upwind SBP)"
1295
This is an experimental feature and may change in future releases.
1296
"""
1297
@inline function splitting_lax_friedrichs(u, orientation_or_normal_direction,
1298
equations::CompressibleEulerEquations2D)
1299
fm = splitting_lax_friedrichs(u, Val{:minus}(), orientation_or_normal_direction,
1300
equations)
1301
fp = splitting_lax_friedrichs(u, Val{:plus}(), orientation_or_normal_direction,
1302
equations)
1303
return fm, fp
1304
end
1305
1306
@inline function splitting_lax_friedrichs(u, ::Val{:plus}, orientation::Integer,
1307
equations::CompressibleEulerEquations2D)
1308
rho, rho_v1, rho_v2, rho_e_total = u
1309
v1 = rho_v1 / rho
1310
v2 = rho_v2 / rho
1311
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
1312
1313
a = sqrt(equations.gamma * p / rho)
1314
H = (rho_e_total + p) / rho
1315
lambda = 0.5f0 * (sqrt(v1^2 + v2^2) + a)
1316
1317
if orientation == 1
1318
#lambda = 0.5 * (abs(v1) + a)
1319
f1p = 0.5f0 * rho * v1 + lambda * u[1]
1320
f2p = 0.5f0 * rho * v1 * v1 + 0.5f0 * p + lambda * u[2]
1321
f3p = 0.5f0 * rho * v1 * v2 + lambda * u[3]
1322
f4p = 0.5f0 * rho * v1 * H + lambda * u[4]
1323
else # orientation == 2
1324
#lambda = 0.5 * (abs(v2) + a)
1325
f1p = 0.5f0 * rho * v2 + lambda * u[1]
1326
f2p = 0.5f0 * rho * v2 * v1 + lambda * u[2]
1327
f3p = 0.5f0 * rho * v2 * v2 + 0.5f0 * p + lambda * u[3]
1328
f4p = 0.5f0 * rho * v2 * H + lambda * u[4]
1329
end
1330
return SVector(f1p, f2p, f3p, f4p)
1331
end
1332
1333
@inline function splitting_lax_friedrichs(u, ::Val{:minus}, orientation::Integer,
1334
equations::CompressibleEulerEquations2D)
1335
rho, rho_v1, rho_v2, rho_e_total = u
1336
v1 = rho_v1 / rho
1337
v2 = rho_v2 / rho
1338
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
1339
1340
a = sqrt(equations.gamma * p / rho)
1341
H = (rho_e_total + p) / rho
1342
lambda = 0.5f0 * (sqrt(v1^2 + v2^2) + a)
1343
1344
if orientation == 1
1345
#lambda = 0.5 * (abs(v1) + a)
1346
f1m = 0.5f0 * rho * v1 - lambda * u[1]
1347
f2m = 0.5f0 * rho * v1 * v1 + 0.5f0 * p - lambda * u[2]
1348
f3m = 0.5f0 * rho * v1 * v2 - lambda * u[3]
1349
f4m = 0.5f0 * rho * v1 * H - lambda * u[4]
1350
else # orientation == 2
1351
#lambda = 0.5 * (abs(v2) + a)
1352
f1m = 0.5f0 * rho * v2 - lambda * u[1]
1353
f2m = 0.5f0 * rho * v2 * v1 - lambda * u[2]
1354
f3m = 0.5f0 * rho * v2 * v2 + 0.5f0 * p - lambda * u[3]
1355
f4m = 0.5f0 * rho * v2 * H - lambda * u[4]
1356
end
1357
return SVector(f1m, f2m, f3m, f4m)
1358
end
1359
1360
@inline function splitting_lax_friedrichs(u, ::Val{:plus},
1361
normal_direction::AbstractVector,
1362
equations::CompressibleEulerEquations2D)
1363
rho_e_total = last(u)
1364
rho, v1, v2, p = cons2prim(u, equations)
1365
1366
a = sqrt(equations.gamma * p / rho)
1367
H = (rho_e_total + p) / rho
1368
lambda = 0.5f0 * (sqrt(v1^2 + v2^2) + a)
1369
1370
v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]
1371
rho_v_normal = rho * v_normal
1372
1373
f1p = 0.5f0 * rho_v_normal + lambda * u[1]
1374
f2p = 0.5f0 * rho_v_normal * v1 + 0.5f0 * p * normal_direction[1] + lambda * u[2]
1375
f3p = 0.5f0 * rho_v_normal * v2 + 0.5f0 * p * normal_direction[2] + lambda * u[3]
1376
f4p = 0.5f0 * rho_v_normal * H + lambda * u[4]
1377
1378
return SVector(f1p, f2p, f3p, f4p)
1379
end
1380
1381
@inline function splitting_lax_friedrichs(u, ::Val{:minus},
1382
normal_direction::AbstractVector,
1383
equations::CompressibleEulerEquations2D)
1384
rho_e_total = last(u)
1385
rho, v1, v2, p = cons2prim(u, equations)
1386
1387
a = sqrt(equations.gamma * p / rho)
1388
H = (rho_e_total + p) / rho
1389
lambda = 0.5f0 * (sqrt(v1^2 + v2^2) + a)
1390
1391
v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]
1392
rho_v_normal = rho * v_normal
1393
1394
f1m = 0.5f0 * rho_v_normal - lambda * u[1]
1395
f2m = 0.5f0 * rho_v_normal * v1 + 0.5f0 * p * normal_direction[1] - lambda * u[2]
1396
f3m = 0.5f0 * rho_v_normal * v2 + 0.5f0 * p * normal_direction[2] - lambda * u[3]
1397
f4m = 0.5f0 * rho_v_normal * H - lambda * u[4]
1398
1399
return SVector(f1m, f2m, f3m, f4m)
1400
end
1401
1402
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation as the
1403
# maximum velocity magnitude plus the maximum speed of sound
1404
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
1405
equations::CompressibleEulerEquations2D)
1406
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1407
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1408
1409
# Get the velocity value in the appropriate direction
1410
if orientation == 1
1411
v_ll = v1_ll
1412
v_rr = v1_rr
1413
else # orientation == 2
1414
v_ll = v2_ll
1415
v_rr = v2_rr
1416
end
1417
# Calculate sound speeds
1418
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1419
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1420
1421
return max(abs(v_ll), abs(v_rr)) + max(c_ll, c_rr)
1422
end
1423
1424
@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
1425
equations::CompressibleEulerEquations2D)
1426
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1427
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1428
1429
# Calculate normal velocities and sound speed
1430
# left
1431
v_ll = (v1_ll * normal_direction[1]
1432
+
1433
v2_ll * normal_direction[2])
1434
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1435
# right
1436
v_rr = (v1_rr * normal_direction[1]
1437
+
1438
v2_rr * normal_direction[2])
1439
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1440
1441
return max(abs(v_ll), abs(v_rr)) + max(c_ll, c_rr) * norm(normal_direction)
1442
end
1443
1444
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
1445
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
1446
equations::CompressibleEulerEquations2D)
1447
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1448
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1449
1450
# Get the velocity value in the appropriate direction
1451
if orientation == 1
1452
v_ll = v1_ll
1453
v_rr = v1_rr
1454
else # orientation == 2
1455
v_ll = v2_ll
1456
v_rr = v2_rr
1457
end
1458
# Calculate sound speeds
1459
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1460
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1461
1462
return max(abs(v_ll) + c_ll, abs(v_rr) + c_rr)
1463
end
1464
1465
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
1466
@inline function max_abs_speed(u_ll, u_rr, normal_direction::AbstractVector,
1467
equations::CompressibleEulerEquations2D)
1468
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1469
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1470
1471
# Calculate normal velocities and sound speeds
1472
# left
1473
v_ll = (v1_ll * normal_direction[1]
1474
+
1475
v2_ll * normal_direction[2])
1476
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1477
# right
1478
v_rr = (v1_rr * normal_direction[1]
1479
+
1480
v2_rr * normal_direction[2])
1481
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1482
1483
norm_ = norm(normal_direction)
1484
return max(abs(v_ll) + c_ll * norm_,
1485
abs(v_rr) + c_rr * norm_)
1486
end
1487
1488
# Calculate estimate for minimum and maximum wave speeds for HLL-type fluxes
1489
@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,
1490
equations::CompressibleEulerEquations2D)
1491
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1492
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1493
1494
if orientation == 1 # x-direction
1495
λ_min = v1_ll - sqrt(equations.gamma * p_ll / rho_ll)
1496
λ_max = v1_rr + sqrt(equations.gamma * p_rr / rho_rr)
1497
else # y-direction
1498
λ_min = v2_ll - sqrt(equations.gamma * p_ll / rho_ll)
1499
λ_max = v2_rr + sqrt(equations.gamma * p_rr / rho_rr)
1500
end
1501
1502
return λ_min, λ_max
1503
end
1504
1505
@inline function min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
1506
equations::CompressibleEulerEquations2D)
1507
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1508
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1509
1510
v_normal_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
1511
v_normal_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
1512
1513
norm_ = norm(normal_direction)
1514
# The v_normals are already scaled by the norm
1515
λ_min = v_normal_ll - sqrt(equations.gamma * p_ll / rho_ll) * norm_
1516
λ_max = v_normal_rr + sqrt(equations.gamma * p_rr / rho_rr) * norm_
1517
1518
return λ_min, λ_max
1519
end
1520
1521
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
1522
@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,
1523
equations::CompressibleEulerEquations2D)
1524
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1525
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1526
1527
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1528
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1529
1530
if orientation == 1 # x-direction
1531
λ_min = min(v1_ll - c_ll, v1_rr - c_rr)
1532
λ_max = max(v1_ll + c_ll, v1_rr + c_rr)
1533
else # y-direction
1534
λ_min = min(v2_ll - c_ll, v2_rr - c_rr)
1535
λ_max = max(v2_ll + c_ll, v2_rr + c_rr)
1536
end
1537
1538
return λ_min, λ_max
1539
end
1540
1541
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
1542
@inline function min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector,
1543
equations::CompressibleEulerEquations2D)
1544
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1545
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1546
1547
norm_ = norm(normal_direction)
1548
1549
c_ll = sqrt(equations.gamma * p_ll / rho_ll) * norm_
1550
c_rr = sqrt(equations.gamma * p_rr / rho_rr) * norm_
1551
1552
v_normal_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
1553
v_normal_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
1554
1555
# The v_normals are already scaled by the norm
1556
λ_min = min(v_normal_ll - c_ll, v_normal_rr - c_rr)
1557
λ_max = max(v_normal_ll + c_ll, v_normal_rr + c_rr)
1558
1559
return λ_min, λ_max
1560
end
1561
1562
@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,
1563
normal_direction::AbstractVector,
1564
equations::CompressibleEulerEquations2D)
1565
(; gamma) = equations
1566
1567
norm_ = norm(normal_direction)
1568
unit_normal_direction = normal_direction / norm_
1569
1570
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1571
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1572
1573
b_ll = rho_ll / (2 * p_ll)
1574
b_rr = rho_rr / (2 * p_rr)
1575
1576
rho_log = ln_mean(rho_ll, rho_rr)
1577
b_log = ln_mean(b_ll, b_rr)
1578
v1_avg = 0.5f0 * (v1_ll + v1_rr)
1579
v2_avg = 0.5f0 * (v2_ll + v2_rr)
1580
p_avg = 0.5f0 * (rho_ll + rho_rr) / (b_ll + b_rr) # 2 * b_avg = b_ll + b_rr
1581
v_squared_bar = v1_ll * v1_rr + v2_ll * v2_rr
1582
h_bar = gamma / (2 * b_log * (gamma - 1)) + 0.5f0 * v_squared_bar
1583
c_bar = sqrt(gamma * p_avg / rho_log)
1584
1585
v_avg_normal = dot(SVector(v1_avg, v2_avg), unit_normal_direction)
1586
1587
lambda_1 = abs(v_avg_normal - c_bar) * rho_log / (2 * gamma)
1588
lambda_2 = abs(v_avg_normal) * rho_log * (gamma - 1) / gamma
1589
lambda_3 = abs(v_avg_normal + c_bar) * rho_log / (2 * gamma)
1590
lambda_4 = abs(v_avg_normal) * p_avg
1591
1592
v1_minus_c = v1_avg - c_bar * unit_normal_direction[1]
1593
v2_minus_c = v2_avg - c_bar * unit_normal_direction[2]
1594
v1_plus_c = v1_avg + c_bar * unit_normal_direction[1]
1595
v2_plus_c = v2_avg + c_bar * unit_normal_direction[2]
1596
v1_tangential = v1_avg - v_avg_normal * unit_normal_direction[1]
1597
v2_tangential = v2_avg - v_avg_normal * unit_normal_direction[2]
1598
1599
entropy_vars_jump = cons2entropy(u_rr, equations) - cons2entropy(u_ll, equations)
1600
entropy_var_rho_jump, entropy_var_rho_v1_jump,
1601
entropy_var_rho_v2_jump, entropy_var_rho_e_total_jump = entropy_vars_jump
1602
1603
velocity_minus_c_dot_entropy_vars_jump = v1_minus_c * entropy_var_rho_v1_jump +
1604
v2_minus_c * entropy_var_rho_v2_jump
1605
velocity_plus_c_dot_entropy_vars_jump = v1_plus_c * entropy_var_rho_v1_jump +
1606
v2_plus_c * entropy_var_rho_v2_jump
1607
velocity_avg_dot_vjump = v1_avg * entropy_var_rho_v1_jump +
1608
v2_avg * entropy_var_rho_v2_jump
1609
w1 = lambda_1 * (entropy_var_rho_jump + velocity_minus_c_dot_entropy_vars_jump +
1610
(h_bar - c_bar * v_avg_normal) * entropy_var_rho_e_total_jump)
1611
w2 = lambda_2 * (entropy_var_rho_jump + velocity_avg_dot_vjump +
1612
v_squared_bar / 2 * entropy_var_rho_e_total_jump)
1613
w3 = lambda_3 * (entropy_var_rho_jump + velocity_plus_c_dot_entropy_vars_jump +
1614
(h_bar + c_bar * v_avg_normal) * entropy_var_rho_e_total_jump)
1615
1616
entropy_var_v_normal_jump = dot(SVector(entropy_var_rho_v1_jump,
1617
entropy_var_rho_v2_jump),
1618
unit_normal_direction)
1619
1620
dissipation_rho = w1 + w2 + w3
1621
1622
dissipation_rho_v1 = (w1 * v1_minus_c +
1623
w2 * v1_avg +
1624
w3 * v1_plus_c +
1625
lambda_4 * (entropy_var_rho_v1_jump -
1626
unit_normal_direction[1] * entropy_var_v_normal_jump +
1627
entropy_var_rho_e_total_jump * v1_tangential))
1628
1629
dissipation_rho_v2 = (w1 * v2_minus_c +
1630
w2 * v2_avg +
1631
w3 * v2_plus_c +
1632
lambda_4 * (entropy_var_rho_v2_jump -
1633
unit_normal_direction[2] * entropy_var_v_normal_jump +
1634
entropy_var_rho_e_total_jump * v2_tangential))
1635
1636
v_tangential_dot_entropy_vars_jump = v1_tangential * entropy_var_rho_v1_jump +
1637
v2_tangential * entropy_var_rho_v2_jump
1638
1639
dissipation_rhoe = (w1 * (h_bar - c_bar * v_avg_normal) +
1640
w2 * 0.5f0 * v_squared_bar +
1641
w3 * (h_bar + c_bar * v_avg_normal) +
1642
lambda_4 * (v_tangential_dot_entropy_vars_jump +
1643
entropy_var_rho_e_total_jump *
1644
(v1_avg^2 + v2_avg^2 - v_avg_normal^2)))
1645
1646
return -0.5f0 *
1647
SVector(dissipation_rho, dissipation_rho_v1, dissipation_rho_v2,
1648
dissipation_rhoe) * norm_
1649
end
1650
1651
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction
1652
# has been normalized prior to this rotation of the state vector
1653
@inline function rotate_to_x(u, normal_vector, equations::CompressibleEulerEquations2D)
1654
# cos and sin of the angle between the x-axis and the normalized normal_vector are
1655
# the normalized vector's x and y coordinates respectively (see unit circle).
1656
c = normal_vector[1]
1657
s = normal_vector[2]
1658
1659
# Apply the 2D rotation matrix with normal and tangent directions of the form
1660
# [ 1 0 0 0;
1661
# 0 n_1 n_2 0;
1662
# 0 t_1 t_2 0;
1663
# 0 0 0 1 ]
1664
# where t_1 = -n_2 and t_2 = n_1
1665
1666
return SVector(u[1],
1667
c * u[2] + s * u[3],
1668
-s * u[2] + c * u[3],
1669
u[4])
1670
end
1671
1672
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction
1673
# has been normalized prior to this back-rotation of the state vector
1674
@inline function rotate_from_x(u, normal_vector,
1675
equations::CompressibleEulerEquations2D)
1676
# cos and sin of the angle between the x-axis and the normalized normal_vector are
1677
# the normalized vector's x and y coordinates respectively (see unit circle).
1678
c = normal_vector[1]
1679
s = normal_vector[2]
1680
1681
# Apply the 2D back-rotation matrix with normal and tangent directions of the form
1682
# [ 1 0 0 0;
1683
# 0 n_1 t_1 0;
1684
# 0 n_2 t_2 0;
1685
# 0 0 0 1 ]
1686
# where t_1 = -n_2 and t_2 = n_1
1687
1688
return SVector(u[1],
1689
c * u[2] - s * u[3],
1690
s * u[2] + c * u[3],
1691
u[4])
1692
end
1693
1694
"""
1695
flux_hllc(u_ll, u_rr, orientation_or_normal_direction, equations::CompressibleEulerEquations2D)
1696
1697
Computes the HLLC flux (HLL with Contact) for compressible Euler equations developed by E.F. Toro
1698
[Lecture slides](http://www.prague-sum.com/download/2012/Toro_2-HLLC-RiemannSolver.pdf)
1699
Signal speeds: [DOI: 10.1137/S1064827593260140](https://doi.org/10.1137/S1064827593260140)
1700
"""
1701
function flux_hllc(u_ll, u_rr, orientation::Integer,
1702
equations::CompressibleEulerEquations2D)
1703
# Calculate primitive variables and speed of sound
1704
rho_ll, rho_v1_ll, rho_v2_ll, rho_e_total_ll = u_ll
1705
rho_rr, rho_v1_rr, rho_v2_rr, rho_e_total_rr = u_rr
1706
1707
v1_ll = rho_v1_ll / rho_ll
1708
v2_ll = rho_v2_ll / rho_ll
1709
e_ll = rho_e_total_ll / rho_ll
1710
p_ll = (equations.gamma - 1) *
1711
(rho_e_total_ll - 0.5f0 * rho_ll * (v1_ll^2 + v2_ll^2))
1712
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1713
1714
v1_rr = rho_v1_rr / rho_rr
1715
v2_rr = rho_v2_rr / rho_rr
1716
e_rr = rho_e_total_rr / rho_rr
1717
p_rr = (equations.gamma - 1) *
1718
(rho_e_total_rr - 0.5f0 * rho_rr * (v1_rr^2 + v2_rr^2))
1719
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1720
1721
# Obtain left and right fluxes
1722
f_ll = flux(u_ll, orientation, equations)
1723
f_rr = flux(u_rr, orientation, equations)
1724
1725
# Compute Roe averages
1726
sqrt_rho_ll = sqrt(rho_ll)
1727
sqrt_rho_rr = sqrt(rho_rr)
1728
sum_sqrt_rho = sqrt_rho_ll + sqrt_rho_rr
1729
if orientation == 1 # x-direction
1730
vel_L = v1_ll
1731
vel_R = v1_rr
1732
elseif orientation == 2 # y-direction
1733
vel_L = v2_ll
1734
vel_R = v2_rr
1735
end
1736
vel_roe = (sqrt_rho_ll * vel_L + sqrt_rho_rr * vel_R) / sum_sqrt_rho
1737
v1_roe = sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr
1738
v2_roe = sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr
1739
vel_roe_mag = (v1_roe^2 + v2_roe^2) / sum_sqrt_rho^2
1740
H_ll = (rho_e_total_ll + p_ll) / rho_ll
1741
H_rr = (rho_e_total_rr + p_rr) / rho_rr
1742
H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) / sum_sqrt_rho
1743
c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * vel_roe_mag))
1744
Ssl = min(vel_L - c_ll, vel_roe - c_roe)
1745
Ssr = max(vel_R + c_rr, vel_roe + c_roe)
1746
sMu_L = Ssl - vel_L
1747
sMu_R = Ssr - vel_R
1748
1749
if Ssl >= 0
1750
f1 = f_ll[1]
1751
f2 = f_ll[2]
1752
f3 = f_ll[3]
1753
f4 = f_ll[4]
1754
elseif Ssr <= 0
1755
f1 = f_rr[1]
1756
f2 = f_rr[2]
1757
f3 = f_rr[3]
1758
f4 = f_rr[4]
1759
else
1760
SStar = (p_rr - p_ll + rho_ll * vel_L * sMu_L - rho_rr * vel_R * sMu_R) /
1761
(rho_ll * sMu_L - rho_rr * sMu_R)
1762
if Ssl <= 0 <= SStar
1763
densStar = rho_ll * sMu_L / (Ssl - SStar)
1764
enerStar = e_ll + (SStar - vel_L) * (SStar + p_ll / (rho_ll * sMu_L))
1765
UStar1 = densStar
1766
UStar4 = densStar * enerStar
1767
if orientation == 1 # x-direction
1768
UStar2 = densStar * SStar
1769
UStar3 = densStar * v2_ll
1770
elseif orientation == 2 # y-direction
1771
UStar2 = densStar * v1_ll
1772
UStar3 = densStar * SStar
1773
end
1774
f1 = f_ll[1] + Ssl * (UStar1 - rho_ll)
1775
f2 = f_ll[2] + Ssl * (UStar2 - rho_v1_ll)
1776
f3 = f_ll[3] + Ssl * (UStar3 - rho_v2_ll)
1777
f4 = f_ll[4] + Ssl * (UStar4 - rho_e_total_ll)
1778
else
1779
densStar = rho_rr * sMu_R / (Ssr - SStar)
1780
enerStar = e_rr + (SStar - vel_R) * (SStar + p_rr / (rho_rr * sMu_R))
1781
UStar1 = densStar
1782
UStar4 = densStar * enerStar
1783
if orientation == 1 # x-direction
1784
UStar2 = densStar * SStar
1785
UStar3 = densStar * v2_rr
1786
elseif orientation == 2 # y-direction
1787
UStar2 = densStar * v1_rr
1788
UStar3 = densStar * SStar
1789
end
1790
f1 = f_rr[1] + Ssr * (UStar1 - rho_rr)
1791
f2 = f_rr[2] + Ssr * (UStar2 - rho_v1_rr)
1792
f3 = f_rr[3] + Ssr * (UStar3 - rho_v2_rr)
1793
f4 = f_rr[4] + Ssr * (UStar4 - rho_e_total_rr)
1794
end
1795
end
1796
return SVector(f1, f2, f3, f4)
1797
end
1798
1799
function flux_hllc(u_ll, u_rr, normal_direction::AbstractVector,
1800
equations::CompressibleEulerEquations2D)
1801
# Calculate primitive variables and speed of sound
1802
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1803
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1804
1805
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
1806
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
1807
1808
norm_ = norm(normal_direction)
1809
norm_sq = norm_ * norm_
1810
inv_norm_sq = inv(norm_sq)
1811
1812
c_ll = sqrt(equations.gamma * p_ll / rho_ll) * norm_
1813
c_rr = sqrt(equations.gamma * p_rr / rho_rr) * norm_
1814
1815
# Obtain left and right fluxes
1816
f_ll = flux(u_ll, normal_direction, equations)
1817
f_rr = flux(u_rr, normal_direction, equations)
1818
1819
# Compute Roe averages
1820
sqrt_rho_ll = sqrt(rho_ll)
1821
sqrt_rho_rr = sqrt(rho_rr)
1822
sum_sqrt_rho = sqrt_rho_ll + sqrt_rho_rr
1823
1824
v1_roe = (sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr) / sum_sqrt_rho
1825
v2_roe = (sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr) / sum_sqrt_rho
1826
vel_roe = v1_roe * normal_direction[1] + v2_roe * normal_direction[2]
1827
vel_roe_mag = v1_roe^2 + v2_roe^2
1828
1829
e_ll = u_ll[4] / rho_ll
1830
e_rr = u_rr[4] / rho_rr
1831
1832
H_ll = (u_ll[4] + p_ll) / rho_ll
1833
H_rr = (u_rr[4] + p_rr) / rho_rr
1834
1835
H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) / sum_sqrt_rho
1836
c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * vel_roe_mag)) * norm_
1837
1838
Ssl = min(v_dot_n_ll - c_ll, vel_roe - c_roe)
1839
Ssr = max(v_dot_n_rr + c_rr, vel_roe + c_roe)
1840
sMu_L = Ssl - v_dot_n_ll
1841
sMu_R = Ssr - v_dot_n_rr
1842
1843
if Ssl >= 0
1844
f1 = f_ll[1]
1845
f2 = f_ll[2]
1846
f3 = f_ll[3]
1847
f4 = f_ll[4]
1848
elseif Ssr <= 0
1849
f1 = f_rr[1]
1850
f2 = f_rr[2]
1851
f3 = f_rr[3]
1852
f4 = f_rr[4]
1853
else
1854
SStar = (rho_ll * v_dot_n_ll * sMu_L - rho_rr * v_dot_n_rr * sMu_R +
1855
(p_rr - p_ll) * norm_sq) / (rho_ll * sMu_L - rho_rr * sMu_R)
1856
if Ssl <= 0 <= SStar
1857
densStar = rho_ll * sMu_L / (Ssl - SStar)
1858
enerStar = e_ll +
1859
(SStar - v_dot_n_ll) *
1860
(SStar * inv_norm_sq + p_ll / (rho_ll * sMu_L))
1861
UStar1 = densStar
1862
UStar2 = densStar *
1863
(v1_ll + (SStar - v_dot_n_ll) * normal_direction[1] * inv_norm_sq)
1864
UStar3 = densStar *
1865
(v2_ll + (SStar - v_dot_n_ll) * normal_direction[2] * inv_norm_sq)
1866
UStar4 = densStar * enerStar
1867
f1 = f_ll[1] + Ssl * (UStar1 - u_ll[1])
1868
f2 = f_ll[2] + Ssl * (UStar2 - u_ll[2])
1869
f3 = f_ll[3] + Ssl * (UStar3 - u_ll[3])
1870
f4 = f_ll[4] + Ssl * (UStar4 - u_ll[4])
1871
else
1872
densStar = rho_rr * sMu_R / (Ssr - SStar)
1873
enerStar = e_rr +
1874
(SStar - v_dot_n_rr) *
1875
(SStar * inv_norm_sq + p_rr / (rho_rr * sMu_R))
1876
UStar1 = densStar
1877
UStar2 = densStar *
1878
(v1_rr + (SStar - v_dot_n_rr) * normal_direction[1] * inv_norm_sq)
1879
UStar3 = densStar *
1880
(v2_rr + (SStar - v_dot_n_rr) * normal_direction[2] * inv_norm_sq)
1881
UStar4 = densStar * enerStar
1882
f1 = f_rr[1] + Ssr * (UStar1 - u_rr[1])
1883
f2 = f_rr[2] + Ssr * (UStar2 - u_rr[2])
1884
f3 = f_rr[3] + Ssr * (UStar3 - u_rr[3])
1885
f4 = f_rr[4] + Ssr * (UStar4 - u_rr[4])
1886
end
1887
end
1888
return SVector(f1, f2, f3, f4)
1889
end
1890
1891
"""
1892
min_max_speed_einfeldt(u_ll, u_rr, orientation, equations::CompressibleEulerEquations2D)
1893
1894
Computes the HLLE (Harten-Lax-van Leer-Einfeldt) flux for the compressible Euler equations.
1895
Special estimates of the signal velocites and linearization of the Riemann problem developed
1896
by Einfeldt to ensure that the internal energy and density remain positive during the computation
1897
of the numerical flux.
1898
1899
- Bernd Einfeldt (1988)
1900
On Godunov-type methods for gas dynamics.
1901
[DOI: 10.1137/0725021](https://doi.org/10.1137/0725021)
1902
- Bernd Einfeldt, Claus-Dieter Munz, Philip L. Roe and Björn Sjögreen (1991)
1903
On Godunov-type methods near low densities.
1904
[DOI: 10.1016/0021-9991(91)90211-3](https://doi.org/10.1016/0021-9991(91)90211-3)
1905
"""
1906
@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,
1907
equations::CompressibleEulerEquations2D)
1908
# Calculate primitive variables, enthalpy and speed of sound
1909
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1910
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1911
1912
# `u_ll[4]` is total energy `rho_e_total_ll` on the left
1913
H_ll = (u_ll[4] + p_ll) / rho_ll
1914
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
1915
1916
# `u_rr[4]` is total energy `rho_e_total_rr` on the right
1917
H_rr = (u_rr[4] + p_rr) / rho_rr
1918
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
1919
1920
# Compute Roe averages
1921
sqrt_rho_ll = sqrt(rho_ll)
1922
sqrt_rho_rr = sqrt(rho_rr)
1923
inv_sum_sqrt_rho = inv(sqrt_rho_ll + sqrt_rho_rr)
1924
1925
v1_roe = (sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr) * inv_sum_sqrt_rho
1926
v2_roe = (sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr) * inv_sum_sqrt_rho
1927
v_roe_mag = v1_roe^2 + v2_roe^2
1928
1929
H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) * inv_sum_sqrt_rho
1930
c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * v_roe_mag))
1931
1932
# Compute convenience constant for positivity preservation, see
1933
# https://doi.org/10.1016/0021-9991(91)90211-3
1934
beta = sqrt(0.5f0 * (equations.gamma - 1) / equations.gamma)
1935
1936
# Estimate the edges of the Riemann fan (with positivity conservation)
1937
if orientation == 1 # x-direction
1938
SsL = min(v1_roe - c_roe, v1_ll - beta * c_ll, 0)
1939
SsR = max(v1_roe + c_roe, v1_rr + beta * c_rr, 0)
1940
elseif orientation == 2 # y-direction
1941
SsL = min(v2_roe - c_roe, v2_ll - beta * c_ll, 0)
1942
SsR = max(v2_roe + c_roe, v2_rr + beta * c_rr, 0)
1943
end
1944
1945
return SsL, SsR
1946
end
1947
1948
"""
1949
min_max_speed_einfeldt(u_ll, u_rr, normal_direction, equations::CompressibleEulerEquations2D)
1950
1951
Computes the HLLE (Harten-Lax-van Leer-Einfeldt) flux for the compressible Euler equations.
1952
Special estimates of the signal velocites and linearization of the Riemann problem developed
1953
by Einfeldt to ensure that the internal energy and density remain positive during the computation
1954
of the numerical flux.
1955
1956
- Bernd Einfeldt (1988)
1957
On Godunov-type methods for gas dynamics.
1958
[DOI: 10.1137/0725021](https://doi.org/10.1137/0725021)
1959
- Bernd Einfeldt, Claus-Dieter Munz, Philip L. Roe and Björn Sjögreen (1991)
1960
On Godunov-type methods near low densities.
1961
[DOI: 10.1016/0021-9991(91)90211-3](https://doi.org/10.1016/0021-9991(91)90211-3)
1962
"""
1963
@inline function min_max_speed_einfeldt(u_ll, u_rr, normal_direction::AbstractVector,
1964
equations::CompressibleEulerEquations2D)
1965
# Calculate primitive variables, enthalpy and speed of sound
1966
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
1967
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
1968
1969
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
1970
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
1971
1972
norm_ = norm(normal_direction)
1973
1974
# `u_ll[4]` is total energy `rho_e_total_ll` on the left
1975
H_ll = (u_ll[4] + p_ll) / rho_ll
1976
c_ll = sqrt(equations.gamma * p_ll / rho_ll) * norm_
1977
1978
# `u_rr[4]` is total energy `rho_e_total_rr` on the right
1979
H_rr = (u_rr[4] + p_rr) / rho_rr
1980
c_rr = sqrt(equations.gamma * p_rr / rho_rr) * norm_
1981
1982
# Compute Roe averages
1983
sqrt_rho_ll = sqrt(rho_ll)
1984
sqrt_rho_rr = sqrt(rho_rr)
1985
inv_sum_sqrt_rho = inv(sqrt_rho_ll + sqrt_rho_rr)
1986
1987
v1_roe = (sqrt_rho_ll * v1_ll + sqrt_rho_rr * v1_rr) * inv_sum_sqrt_rho
1988
v2_roe = (sqrt_rho_ll * v2_ll + sqrt_rho_rr * v2_rr) * inv_sum_sqrt_rho
1989
v_roe = v1_roe * normal_direction[1] + v2_roe * normal_direction[2]
1990
v_roe_mag = v1_roe^2 + v2_roe^2
1991
1992
H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) * inv_sum_sqrt_rho
1993
c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * v_roe_mag)) * norm_
1994
1995
# Compute convenience constant for positivity preservation, see
1996
# https://doi.org/10.1016/0021-9991(91)90211-3
1997
beta = sqrt(0.5f0 * (equations.gamma - 1) / equations.gamma)
1998
1999
# Estimate the edges of the Riemann fan (with positivity conservation)
2000
SsL = min(v_roe - c_roe, v_dot_n_ll - beta * c_ll, 0)
2001
SsR = max(v_roe + c_roe, v_dot_n_rr + beta * c_rr, 0)
2002
2003
return SsL, SsR
2004
end
2005
2006
@inline function max_abs_speeds(u, equations::CompressibleEulerEquations2D)
2007
rho, v1, v2, p = cons2prim(u, equations)
2008
c = sqrt(equations.gamma * p / rho)
2009
2010
return abs(v1) + c, abs(v2) + c
2011
end
2012
2013
# Convert conservative variables to primitive
2014
@inline function cons2prim(u, equations::CompressibleEulerEquations2D)
2015
rho, rho_v1, rho_v2, rho_e_total = u
2016
2017
v1 = rho_v1 / rho
2018
v2 = rho_v2 / rho
2019
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1 * v1 + rho_v2 * v2))
2020
2021
return SVector(rho, v1, v2, p)
2022
end
2023
2024
# Convert conservative variables to entropy
2025
@inline function cons2entropy(u, equations::CompressibleEulerEquations2D)
2026
rho, rho_v1, rho_v2, rho_e_total = u
2027
2028
v1 = rho_v1 / rho
2029
v2 = rho_v2 / rho
2030
v_square = v1^2 + v2^2
2031
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho * v_square)
2032
s = log(p) - equations.gamma * log(rho)
2033
rho_p = rho / p
2034
2035
w1 = (equations.gamma - s) * equations.inv_gamma_minus_one -
2036
0.5f0 * rho_p * v_square
2037
w2 = rho_p * v1
2038
w3 = rho_p * v2
2039
w4 = -rho_p
2040
2041
return SVector(w1, w2, w3, w4)
2042
end
2043
2044
# Transformation from conservative variables u to entropy vector ds_0/du,
2045
# using the modified specific entropy of Guermond et al. (2019): s_0 = p * rho^(-gamma) / (gamma-1).
2046
# Note: This is *not* the "conventional" specific entropy s = ln(p / rho^(gamma)).
2047
@inline function cons2entropy_guermond_etal(u, equations::CompressibleEulerEquations2D)
2048
rho, rho_v1, rho_v2, rho_e_total = u
2049
2050
v1 = rho_v1 / rho
2051
v2 = rho_v2 / rho
2052
v_square = v1^2 + v2^2
2053
inv_rho_gammap1 = (1 / rho)^(equations.gamma + 1)
2054
2055
# The derivative vector for the modified specific entropy of Guermond et al.
2056
w1 = inv_rho_gammap1 *
2057
(0.5f0 * rho * (equations.gamma + 1) * v_square -
2058
equations.gamma * rho_e_total)
2059
w2 = -rho_v1 * inv_rho_gammap1
2060
w3 = -rho_v2 * inv_rho_gammap1
2061
w4 = (1 / rho)^equations.gamma
2062
2063
return SVector(w1, w2, w3, w4)
2064
end
2065
2066
@inline function entropy2cons(w, equations::CompressibleEulerEquations2D)
2067
# See Hughes, Franca, Mallet (1986) A new finite element formulation for CFD
2068
# [DOI: 10.1016/0045-7825(86)90127-1](https://doi.org/10.1016/0045-7825(86)90127-1)
2069
@unpack gamma = equations
2070
2071
# convert to entropy `-rho * s` used by Hughes, France, Mallet (1986)
2072
# instead of `-rho * s / (gamma - 1)`
2073
V1, V2, V3, V5 = w .* (gamma - 1)
2074
2075
# s = specific entropy, eq. (53)
2076
s = gamma - V1 + (V2^2 + V3^2) / (2 * V5)
2077
2078
# eq. (52)
2079
rho_iota = ((gamma - 1) / (-V5)^gamma)^(equations.inv_gamma_minus_one) *
2080
exp(-s * equations.inv_gamma_minus_one)
2081
2082
# eq. (51)
2083
rho = -rho_iota * V5
2084
rho_v1 = rho_iota * V2
2085
rho_v2 = rho_iota * V3
2086
rho_e_total = rho_iota * (1 - (V2^2 + V3^2) / (2 * V5))
2087
return SVector(rho, rho_v1, rho_v2, rho_e_total)
2088
end
2089
2090
# Convert primitive to conservative variables
2091
@inline function prim2cons(prim, equations::CompressibleEulerEquations2D)
2092
rho, v1, v2, p = prim
2093
rho_v1 = rho * v1
2094
rho_v2 = rho * v2
2095
rho_e_total = p * equations.inv_gamma_minus_one +
2096
0.5f0 * (rho_v1 * v1 + rho_v2 * v2)
2097
return SVector(rho, rho_v1, rho_v2, rho_e_total)
2098
end
2099
2100
@inline function density(u, equations::CompressibleEulerEquations2D)
2101
rho = u[1]
2102
return rho
2103
end
2104
2105
@inline function velocity(u, equations::CompressibleEulerEquations2D)
2106
rho = u[1]
2107
v1 = u[2] / rho
2108
v2 = u[3] / rho
2109
return SVector(v1, v2)
2110
end
2111
2112
@inline function velocity(u, orientation::Int, equations::CompressibleEulerEquations2D)
2113
rho = u[1]
2114
v = u[orientation + 1] / rho
2115
return v
2116
end
2117
2118
@inline function pressure(u, equations::CompressibleEulerEquations2D)
2119
rho, rho_v1, rho_v2, rho_e_total = u
2120
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1^2 + rho_v2^2) / rho)
2121
return p
2122
end
2123
2124
# Transformation from conservative variables u to d(p)/d(u)
2125
@inline function gradient_conservative(::typeof(pressure),
2126
u, equations::CompressibleEulerEquations2D)
2127
rho, rho_v1, rho_v2, rho_e_total = u
2128
2129
v1 = rho_v1 / rho
2130
v2 = rho_v2 / rho
2131
v_square = v1^2 + v2^2
2132
2133
return (equations.gamma - 1) * SVector(0.5f0 * v_square, -v1, -v2, 1)
2134
end
2135
2136
@inline function density_pressure(u, equations::CompressibleEulerEquations2D)
2137
rho, rho_v1, rho_v2, rho_e_total = u
2138
rho_times_p = (equations.gamma - 1) *
2139
(rho * rho_e_total - 0.5f0 * (rho_v1^2 + rho_v2^2))
2140
return rho_times_p
2141
end
2142
2143
# Calculates the entropy flux in direction "orientation" and the entropy variables for a state cons
2144
# NOTE: This method seems to work currently (b82534e) but is never used anywhere. Thus it is
2145
# commented here until someone uses it or writes a test for it.
2146
# @inline function cons2entropyvars_and_flux(gamma::Float64, cons, orientation::Int)
2147
# entropy = MVector{4, Float64}(undef)
2148
# v = (cons[2] / cons[1] , cons[3] / cons[1])
2149
# v_square= v[1]*v[1]+v[2]*v[2]
2150
# p = (gamma - 1) * (cons[4] - 1/2 * (cons[2] * v[1] + cons[3] * v[2]))
2151
# rho_p = cons[1] / p
2152
# # thermodynamic entropy
2153
# s = log(p) - gamma*log(cons[1])
2154
# # mathematical entropy
2155
# S = - s*cons[1]/(gamma-1)
2156
# # entropy variables
2157
# entropy[1] = (gamma - s)/(gamma-1) - 0.5*rho_p*v_square
2158
# entropy[2] = rho_p*v[1]
2159
# entropy[3] = rho_p*v[2]
2160
# entropy[4] = -rho_p
2161
# # entropy flux
2162
# entropy_flux = S*v[orientation]
2163
# return entropy, entropy_flux
2164
# end
2165
2166
# Calculate thermodynamic entropy for a conservative state `cons`
2167
@inline function entropy_thermodynamic(cons, equations::CompressibleEulerEquations2D)
2168
# Pressure
2169
p = (equations.gamma - 1) * (cons[4] - 0.5f0 * (cons[2]^2 + cons[3]^2) / cons[1])
2170
2171
# Thermodynamic entropy
2172
s = log(p) - equations.gamma * log(cons[1])
2173
2174
return s
2175
end
2176
2177
# Calculate mathematical entropy for a conservative state `cons`
2178
@inline function entropy_math(cons, equations::CompressibleEulerEquations2D)
2179
# Mathematical entropy
2180
S = -entropy_thermodynamic(cons, equations) * cons[1] *
2181
equations.inv_gamma_minus_one
2182
2183
return S
2184
end
2185
2186
# Transformation from conservative variables u to d(s)/d(u)
2187
@inline function gradient_conservative(::typeof(entropy_math),
2188
u, equations::CompressibleEulerEquations2D)
2189
return cons2entropy(u, equations)
2190
end
2191
2192
@doc raw"""
2193
entropy_guermond_etal(u, equations::CompressibleEulerEquations2D)
2194
2195
Calculate the modified specific entropy of Guermond et al. (2019):
2196
```math
2197
s_0 = p * \rho^{-\gamma} / (\gamma-1).
2198
```
2199
Note: This is *not* the "conventional" specific entropy ``s = ln(p / \rho^\gamma)``.
2200
- Guermond at al. (2019)
2201
Invariant domain preserving discretization-independent schemes and convex limiting for hyperbolic systems.
2202
[DOI: 10.1016/j.cma.2018.11.036](https://doi.org/10.1016/j.cma.2018.11.036)
2203
"""
2204
@inline function entropy_guermond_etal(u, equations::CompressibleEulerEquations2D)
2205
rho, rho_v1, rho_v2, rho_e_total = u
2206
2207
# Modified specific entropy from Guermond et al. (2019)
2208
s = (rho_e_total - 0.5f0 * (rho_v1^2 + rho_v2^2) / rho) * (1 / rho)^equations.gamma
2209
2210
return s
2211
end
2212
2213
# Transformation from conservative variables u to d(s)/d(u)
2214
@inline function gradient_conservative(::typeof(entropy_guermond_etal),
2215
u, equations::CompressibleEulerEquations2D)
2216
return cons2entropy_guermond_etal(u, equations)
2217
end
2218
2219
# Default entropy is the mathematical entropy
2220
@inline function entropy(cons, equations::CompressibleEulerEquations2D)
2221
return entropy_math(cons, equations)
2222
end
2223
2224
# Calculate total energy for a conservative state `cons`
2225
@inline energy_total(cons, ::CompressibleEulerEquations2D) = cons[4]
2226
2227
# Calculate kinetic energy for a conservative state `cons`
2228
@inline function energy_kinetic(u, equations::CompressibleEulerEquations2D)
2229
rho, rho_v1, rho_v2, rho_e_total = u
2230
return (rho_v1^2 + rho_v2^2) / (2 * rho)
2231
end
2232
2233
# Calculate internal energy for a conservative state `cons`
2234
@inline function energy_internal(cons, equations::CompressibleEulerEquations2D)
2235
return energy_total(cons, equations) - energy_kinetic(cons, equations)
2236
end
2237
2238
@inline function entropy_potential(u, orientation::Int,
2239
equations::CompressibleEulerEquations2D)
2240
if orientation == 1
2241
return u[2]
2242
else # if orientation == 2
2243
return u[3]
2244
end
2245
end
2246
# Version for non-Cartesian meshes, i.e., everything but `TreeMesh`es.
2247
@inline function entropy_potential(u, normal_direction::AbstractVector,
2248
equations::CompressibleEulerEquations2D)
2249
return u[2] * normal_direction[1] +
2250
u[3] * normal_direction[2]
2251
end
2252
2253
# State validation for Newton-bisection method of subcell IDP limiting
2254
@inline function Base.isvalid(u, equations::CompressibleEulerEquations2D)
2255
if u[1] <= 0 || pressure(u, equations) <= 0
2256
return false
2257
end
2258
return true
2259
end
2260
end # @muladd
2261
2262