Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/equations/compressible_euler_1d.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
CompressibleEulerEquations1D(gamma)
10
11
The compressible Euler equations
12
```math
13
\frac{\partial}{\partial t}
14
\begin{pmatrix}
15
\rho \\ \rho v_1 \\ \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 e_{\text{total}} +p) v_1
21
\end{pmatrix}
22
=
23
\begin{pmatrix}
24
0 \\ 0 \\ 0
25
\end{pmatrix}
26
```
27
for an ideal gas with ratio of specific heats `gamma` in one space dimension.
28
Here, ``\rho`` is the density, ``v_1`` the velocity, ``e_{\text{total}}`` the specific total energy, and
29
```math
30
p = (\gamma - 1) \left( \rho e_{\text{total}} - \frac{1}{2} \rho v_1^2 \right)
31
```
32
the pressure.
33
"""
34
struct CompressibleEulerEquations1D{RealT <: Real} <:
35
AbstractCompressibleEulerEquations{1, 3}
36
gamma::RealT # ratio of specific heats
37
inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications
38
39
function CompressibleEulerEquations1D(gamma)
40
γ, inv_gamma_minus_one = promote(gamma, inv(gamma - 1))
41
return new{typeof(γ)}(γ, inv_gamma_minus_one)
42
end
43
end
44
45
function varnames(::typeof(cons2cons), ::CompressibleEulerEquations1D)
46
return ("rho", "rho_v1", "rho_e_total")
47
end
48
varnames(::typeof(cons2prim), ::CompressibleEulerEquations1D) = ("rho", "v1", "p")
49
50
"""
51
initial_condition_constant(x, t, equations::CompressibleEulerEquations1D)
52
53
A constant initial condition to test free-stream preservation.
54
"""
55
function initial_condition_constant(x, t, equations::CompressibleEulerEquations1D)
56
RealT = eltype(x)
57
rho = 1
58
rho_v1 = convert(RealT, 0.1)
59
rho_e_total = 10
60
return SVector(rho, rho_v1, rho_e_total)
61
end
62
63
"""
64
initial_condition_convergence_test(x, t, equations::CompressibleEulerEquations1D)
65
66
A smooth initial condition used for convergence tests in combination with
67
[`source_terms_convergence_test`](@ref)
68
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
69
"""
70
function initial_condition_convergence_test(x, t,
71
equations::CompressibleEulerEquations1D)
72
RealT = eltype(x)
73
c = 2
74
A = convert(RealT, 0.1)
75
L = 2
76
f = 1.0f0 / L
77
ω = 2 * convert(RealT, pi) * f
78
ini = c + A * sin(ω * (x[1] - t))
79
80
rho = ini
81
rho_v1 = ini
82
rho_e_total = ini^2
83
84
return SVector(rho, rho_v1, rho_e_total)
85
end
86
87
"""
88
source_terms_convergence_test(u, x, t, equations::CompressibleEulerEquations1D)
89
90
Source terms used for convergence tests in combination with
91
[`initial_condition_convergence_test`](@ref)
92
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
93
94
References for the method of manufactured solutions (MMS):
95
- Kambiz Salari and Patrick Knupp (2000)
96
Code Verification by the Method of Manufactured Solutions
97
[DOI: 10.2172/759450](https://doi.org/10.2172/759450)
98
- Patrick J. Roache (2002)
99
Code Verification by the Method of Manufactured Solutions
100
[DOI: 10.1115/1.1436090](https://doi.org/10.1115/1.1436090)
101
"""
102
@inline function source_terms_convergence_test(u, x, t,
103
equations::CompressibleEulerEquations1D)
104
# Same settings as in `initial_condition`
105
RealT = eltype(u)
106
c = 2
107
A = convert(RealT, 0.1)
108
L = 2
109
f = 1.0f0 / L
110
ω = 2 * convert(RealT, pi) * f
111
γ = equations.gamma
112
113
x1, = x
114
115
si, co = sincos(ω * (x1 - t))
116
rho = c + A * si
117
rho_x = ω * A * co
118
119
# Note that d/dt rho = -d/dx rho.
120
# This yields du2 = du3 = d/dx p (derivative of pressure).
121
# Other terms vanish because of v = 1.
122
du1 = 0
123
du2 = rho_x * (2 * rho - 0.5f0) * (γ - 1)
124
du3 = du2
125
126
return SVector(du1, du2, du3)
127
end
128
129
"""
130
initial_condition_density_wave(x, t, equations::CompressibleEulerEquations1D)
131
132
A sine wave in the density with constant velocity and pressure; reduces the
133
compressible Euler equations to the linear advection equations.
134
This setup is the test case for stability of EC fluxes from paper
135
- Gregor J. Gassner, Magnus Svärd, Florian J. Hindenlang (2020)
136
Stability issues of entropy-stable and/or split-form high-order schemes
137
[arXiv: 2007.09026](https://arxiv.org/abs/2007.09026)
138
with the following parameters
139
- domain [-1, 1]
140
- mesh = 4x4
141
- polydeg = 5
142
"""
143
function initial_condition_density_wave(x, t, equations::CompressibleEulerEquations1D)
144
RealT = eltype(x)
145
v1 = convert(RealT, 0.1)
146
rho = 1 + convert(RealT, 0.98) * sinpi(2 * (x[1] - t * v1))
147
rho_v1 = rho * v1
148
p = 20
149
rho_e_total = p / (equations.gamma - 1) + 0.5f0 * rho * v1^2
150
return SVector(rho, rho_v1, rho_e_total)
151
end
152
153
"""
154
initial_condition_weak_blast_wave(x, t, equations::CompressibleEulerEquations1D)
155
156
A weak blast wave taken from
157
- Sebastian Hennemann, Gregor J. Gassner (2020)
158
A provably entropy stable subcell shock capturing approach for high order split form DG
159
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
160
"""
161
function initial_condition_weak_blast_wave(x, t,
162
equations::CompressibleEulerEquations1D)
163
# From Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
164
# Set up polar coordinates
165
RealT = eltype(x)
166
inicenter = SVector(0)
167
x_norm = x[1] - inicenter[1]
168
r = abs(x_norm)
169
# The following code is equivalent to
170
# phi = atan(0.0, x_norm)
171
# cos_phi = cos(phi)
172
# in 1D but faster
173
cos_phi = x_norm > 0 ? 1 : -1
174
175
# Calculate primitive variables
176
rho = r > 0.5f0 ? one(RealT) : convert(RealT, 1.1691)
177
v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos_phi
178
p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)
179
180
return prim2cons(SVector(rho, v1, p), equations)
181
end
182
183
"""
184
initial_condition_eoc_test_coupled_euler_gravity(x, t, equations::CompressibleEulerEquations1D)
185
186
One dimensional variant of the setup used for convergence tests of the Euler equations
187
with self-gravity from
188
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
189
A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics
190
[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)
191
!!! note
192
There is no additional source term necessary for the manufactured solution in one
193
spatial dimension. Thus, [`source_terms_eoc_test_coupled_euler_gravity`](@ref) is not
194
present there.
195
"""
196
function initial_condition_eoc_test_coupled_euler_gravity(x, t,
197
equations::CompressibleEulerEquations1D)
198
# OBS! this assumes that γ = 2 other manufactured source terms are incorrect
199
if equations.gamma != 2
200
error("adiabatic constant must be 2 for the coupling convergence test")
201
end
202
RealT = eltype(x)
203
c = 2
204
A = convert(RealT, 0.1)
205
ini = c + A * sinpi(x[1] - t)
206
G = 1 # gravitational constant
207
208
rho = ini
209
v1 = 1
210
p = 2 * ini^2 * G / convert(RealT, pi) # * 2 / ndims, but ndims==1 here
211
212
return prim2cons(SVector(rho, v1, p), equations)
213
end
214
215
"""
216
boundary_condition_slip_wall(u_inner, orientation, direction, x, t,
217
surface_flux_function, equations::CompressibleEulerEquations1D)
218
Determine the boundary numerical surface flux for a slip wall condition.
219
Imposes a zero normal velocity at the wall.
220
Density is taken from the internal solution state and pressure is computed as an
221
exact solution of a 1D Riemann problem. Further details about this boundary state
222
are available in the paper:
223
- J. J. W. van der Vegt and H. van der Ven (2002)
224
Slip flow boundary conditions in discontinuous Galerkin discretizations of
225
the Euler equations of gas dynamics
226
[PDF](https://reports.nlr.nl/bitstream/handle/10921/692/TP-2002-300.pdf?sequence=1)
227
228
Should be used together with [`TreeMesh`](@ref).
229
"""
230
@inline function boundary_condition_slip_wall(u_inner, orientation,
231
direction, x, t,
232
surface_flux_function,
233
equations::CompressibleEulerEquations1D)
234
# compute the primitive variables
235
rho_local, v_normal, p_local = cons2prim(u_inner, equations)
236
237
if isodd(direction) # flip sign of normal to make it outward pointing
238
v_normal *= -1
239
end
240
241
# Get the solution of the pressure Riemann problem
242
# See Section 6.3.3 of
243
# Eleuterio F. Toro (2009)
244
# Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction
245
# [DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)
246
if v_normal <= 0
247
sound_speed = sqrt(equations.gamma * p_local / rho_local) # local sound speed
248
p_star = p_local *
249
(1 + 0.5f0 * (equations.gamma - 1) * v_normal / sound_speed)^(2 *
250
equations.gamma *
251
equations.inv_gamma_minus_one)
252
else # v_normal > 0
253
A = 2 / ((equations.gamma + 1) * rho_local)
254
B = p_local * (equations.gamma - 1) / (equations.gamma + 1)
255
p_star = p_local +
256
0.5f0 * v_normal / A *
257
(v_normal + sqrt(v_normal^2 + 4 * A * (p_local + B)))
258
end
259
260
# For the slip wall we directly set the flux as the normal velocity is zero
261
return SVector(0, p_star, 0)
262
end
263
264
# Calculate 1D flux for a single point
265
@inline function flux(u, orientation::Integer, equations::CompressibleEulerEquations1D)
266
rho, rho_v1, rho_e_total = u
267
v1 = rho_v1 / rho
268
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)
269
# Ignore orientation since it is always "1" in 1D
270
f1 = rho_v1
271
f2 = rho_v1 * v1 + p
272
f3 = (rho_e_total + p) * v1
273
return SVector(f1, f2, f3)
274
end
275
276
"""
277
flux_shima_etal(u_ll, u_rr, orientation, equations::CompressibleEulerEquations1D)
278
279
This flux is is a modification of the original kinetic energy preserving two-point flux by
280
- Yuichi Kuya, Kosuke Totani and Soshi Kawai (2018)
281
Kinetic energy and entropy preserving schemes for compressible flows
282
by split convective forms
283
[DOI: 10.1016/j.jcp.2018.08.058](https://doi.org/10.1016/j.jcp.2018.08.058)
284
285
The modification is in the energy flux to guarantee pressure equilibrium and was developed by
286
- Nao Shima, Yuichi Kuya, Yoshiharu Tamaki, Soshi Kawai (JCP 2020)
287
Preventing spurious pressure oscillations in split convective form discretizations for
288
compressible flows
289
[DOI: 10.1016/j.jcp.2020.110060](https://doi.org/10.1016/j.jcp.2020.110060)
290
"""
291
@inline function flux_shima_etal(u_ll, u_rr, orientation::Integer,
292
equations::CompressibleEulerEquations1D)
293
# Unpack left and right state
294
rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)
295
rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)
296
297
# Average each factor of products in flux
298
rho_avg = 0.5f0 * (rho_ll + rho_rr)
299
v1_avg = 0.5f0 * (v1_ll + v1_rr)
300
p_avg = 0.5f0 * (p_ll + p_rr)
301
kin_avg = 0.5f0 * (v1_ll * v1_rr)
302
303
# Calculate fluxes
304
# Ignore orientation since it is always "1" in 1D
305
pv1_avg = 0.5f0 * (p_ll * v1_rr + p_rr * v1_ll)
306
f1 = rho_avg * v1_avg
307
f2 = f1 * v1_avg + p_avg
308
f3 = p_avg * v1_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv1_avg
309
310
return SVector(f1, f2, f3)
311
end
312
313
"""
314
flux_kennedy_gruber(u_ll, u_rr, orientation, equations::CompressibleEulerEquations1D)
315
316
Kinetic energy preserving two-point flux by
317
- Kennedy and Gruber (2008)
318
Reduced aliasing formulations of the convective terms within the
319
Navier-Stokes equations for a compressible fluid
320
[DOI: 10.1016/j.jcp.2007.09.020](https://doi.org/10.1016/j.jcp.2007.09.020)
321
"""
322
@inline function flux_kennedy_gruber(u_ll, u_rr, orientation::Integer,
323
equations::CompressibleEulerEquations1D)
324
# Unpack left and right state
325
rho_e_total_ll = last(u_ll)
326
rho_e_total_rr = last(u_rr)
327
rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)
328
rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)
329
330
# Average each factor of products in flux
331
rho_avg = 0.5f0 * (rho_ll + rho_rr)
332
v1_avg = 0.5f0 * (v1_ll + v1_rr)
333
p_avg = 0.5f0 * (p_ll + p_rr)
334
e_avg = 0.5f0 * (rho_e_total_ll / rho_ll + rho_e_total_rr / rho_rr)
335
336
# Ignore orientation since it is always "1" in 1D
337
f1 = rho_avg * v1_avg
338
f2 = rho_avg * v1_avg * v1_avg + p_avg
339
f3 = (rho_avg * e_avg + p_avg) * v1_avg
340
341
return SVector(f1, f2, f3)
342
end
343
344
"""
345
flux_chandrashekar(u_ll, u_rr, orientation, equations::CompressibleEulerEquations1D)
346
347
Entropy conserving two-point flux by
348
- Chandrashekar (2013)
349
Kinetic Energy Preserving and Entropy Stable Finite Volume Schemes
350
for Compressible Euler and Navier-Stokes Equations
351
[DOI: 10.4208/cicp.170712.010313a](https://doi.org/10.4208/cicp.170712.010313a)
352
"""
353
@inline function flux_chandrashekar(u_ll, u_rr, orientation::Integer,
354
equations::CompressibleEulerEquations1D)
355
# Unpack left and right state
356
rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)
357
rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)
358
beta_ll = 0.5f0 * rho_ll / p_ll
359
beta_rr = 0.5f0 * rho_rr / p_rr
360
specific_kin_ll = 0.5f0 * (v1_ll^2)
361
specific_kin_rr = 0.5f0 * (v1_rr^2)
362
363
# Compute the necessary mean values
364
rho_avg = 0.5f0 * (rho_ll + rho_rr)
365
rho_mean = ln_mean(rho_ll, rho_rr)
366
beta_mean = ln_mean(beta_ll, beta_rr)
367
beta_avg = 0.5f0 * (beta_ll + beta_rr)
368
v1_avg = 0.5f0 * (v1_ll + v1_rr)
369
p_mean = 0.5f0 * rho_avg / beta_avg
370
velocity_square_avg = specific_kin_ll + specific_kin_rr
371
372
# Calculate fluxes
373
# Ignore orientation since it is always "1" in 1D
374
f1 = rho_mean * v1_avg
375
f2 = f1 * v1_avg + p_mean
376
f3 = f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +
377
f2 * v1_avg
378
379
return SVector(f1, f2, f3)
380
end
381
382
"""
383
flux_ranocha(u_ll, u_rr, orientation_or_normal_direction, equations::CompressibleEulerEquations1D)
384
385
Entropy conserving and kinetic energy preserving two-point flux by
386
- Hendrik Ranocha (2018)
387
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
388
for Hyperbolic Balance Laws
389
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
390
See also
391
- Hendrik Ranocha (2020)
392
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
393
the Euler Equations Using Summation-by-Parts Operators
394
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
395
"""
396
@inline function flux_ranocha(u_ll, u_rr, orientation::Integer,
397
equations::CompressibleEulerEquations1D)
398
# Unpack left and right state
399
rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)
400
rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)
401
402
# Compute the necessary mean values
403
rho_mean = ln_mean(rho_ll, rho_rr)
404
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
405
# in exact arithmetic since
406
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
407
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
408
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
409
v1_avg = 0.5f0 * (v1_ll + v1_rr)
410
p_avg = 0.5f0 * (p_ll + p_rr)
411
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr)
412
413
# Calculate fluxes
414
# Ignore orientation since it is always "1" in 1D
415
f1 = rho_mean * v1_avg
416
f2 = f1 * v1_avg + p_avg
417
f3 = f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +
418
0.5f0 * (p_ll * v1_rr + p_rr * v1_ll)
419
420
return SVector(f1, f2, f3)
421
end
422
423
# While `normal_direction` isn't strictly necessary in 1D, certain solvers assume that
424
# the normal component is incorporated into the numerical flux.
425
#
426
# See `flux(u, normal_direction::AbstractVector, equations::AbstractEquations{1})` for a
427
# similar implementation.
428
@inline function flux_ranocha(u_ll, u_rr, normal_direction::AbstractVector,
429
equations::CompressibleEulerEquations1D)
430
return normal_direction[1] * flux_ranocha(u_ll, u_rr, 1, equations)
431
end
432
433
"""
434
splitting_steger_warming(u, orientation::Integer,
435
equations::CompressibleEulerEquations1D)
436
splitting_steger_warming(u, which::Union{Val{:minus}, Val{:plus}}
437
orientation::Integer,
438
equations::CompressibleEulerEquations1D)
439
440
Splitting of the compressible Euler flux of Steger and Warming.
441
442
Returns a tuple of the fluxes "minus" (associated with waves going into the
443
negative axis direction) and "plus" (associated with waves going into the
444
positive axis direction). If only one of the fluxes is required, use the
445
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
446
447
!!! warning "Experimental implementation (upwind SBP)"
448
This is an experimental feature and may change in future releases.
449
450
## References
451
452
- Joseph L. Steger and R. F. Warming (1979)
453
Flux Vector Splitting of the Inviscid Gasdynamic Equations
454
With Application to Finite Difference Methods
455
[NASA Technical Memorandum](https://ntrs.nasa.gov/api/citations/19790020779/downloads/19790020779.pdf)
456
"""
457
@inline function splitting_steger_warming(u, orientation::Integer,
458
equations::CompressibleEulerEquations1D)
459
fm = splitting_steger_warming(u, Val{:minus}(), orientation, equations)
460
fp = splitting_steger_warming(u, Val{:plus}(), orientation, equations)
461
return fm, fp
462
end
463
464
@inline function splitting_steger_warming(u, ::Val{:plus}, orientation::Integer,
465
equations::CompressibleEulerEquations1D)
466
rho, rho_v1, rho_e_total = u
467
v1 = rho_v1 / rho
468
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)
469
a = sqrt(equations.gamma * p / rho)
470
471
lambda1 = v1
472
lambda2 = v1 + a
473
lambda3 = v1 - a
474
475
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
476
lambda2_p = positive_part(lambda2)
477
lambda3_p = positive_part(lambda3)
478
479
alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p
480
481
rho_2gamma = 0.5f0 * rho / equations.gamma
482
f1p = rho_2gamma * alpha_p
483
f2p = rho_2gamma * (alpha_p * v1 + a * (lambda2_p - lambda3_p))
484
f3p = rho_2gamma * (alpha_p * 0.5f0 * v1^2 + a * v1 * (lambda2_p - lambda3_p)
485
+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)
486
487
return SVector(f1p, f2p, f3p)
488
end
489
490
@inline function splitting_steger_warming(u, ::Val{:minus}, orientation::Integer,
491
equations::CompressibleEulerEquations1D)
492
rho, rho_v1, rho_e_total = u
493
v1 = rho_v1 / rho
494
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)
495
a = sqrt(equations.gamma * p / rho)
496
497
lambda1 = v1
498
lambda2 = v1 + a
499
lambda3 = v1 - a
500
501
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
502
lambda2_m = negative_part(lambda2)
503
lambda3_m = negative_part(lambda3)
504
505
alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m
506
507
rho_2gamma = 0.5f0 * rho / equations.gamma
508
f1m = rho_2gamma * alpha_m
509
f2m = rho_2gamma * (alpha_m * v1 + a * (lambda2_m - lambda3_m))
510
f3m = rho_2gamma * (alpha_m * 0.5f0 * v1^2 + a * v1 * (lambda2_m - lambda3_m)
511
+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)
512
513
return SVector(f1m, f2m, f3m)
514
end
515
516
"""
517
splitting_vanleer_haenel(u, orientation::Integer,
518
equations::CompressibleEulerEquations1D)
519
splitting_vanleer_haenel(u, which::Union{Val{:minus}, Val{:plus}}
520
orientation::Integer,
521
equations::CompressibleEulerEquations1D)
522
523
Splitting of the compressible Euler flux from van Leer. This splitting further
524
contains a reformulation due to Hänel et al. where the energy flux uses the
525
enthalpy. The pressure splitting is independent from the splitting of the
526
convective terms. As such there are many pressure splittings suggested across
527
the literature. We implement the 'p4' variant suggested by Liou and Steffen as
528
it proved the most robust in practice.
529
530
Returns a tuple of the fluxes "minus" (associated with waves going into the
531
negative axis direction) and "plus" (associated with waves going into the
532
positive axis direction). If only one of the fluxes is required, use the
533
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
534
535
!!! warning "Experimental implementation (upwind SBP)"
536
This is an experimental feature and may change in future releases.
537
538
## References
539
540
- Bram van Leer (1982)
541
Flux-Vector Splitting for the Euler Equation
542
[DOI: 10.1007/978-3-642-60543-7_5](https://doi.org/10.1007/978-3-642-60543-7_5)
543
- D. Hänel, R. Schwane and G. Seider (1987)
544
On the accuracy of upwind schemes for the solution of the Navier-Stokes equations
545
[DOI: 10.2514/6.1987-1105](https://doi.org/10.2514/6.1987-1105)
546
- Meng-Sing Liou and Chris J. Steffen, Jr. (1991)
547
High-Order Polynomial Expansions (HOPE) for Flux-Vector Splitting
548
[NASA Technical Memorandum](https://ntrs.nasa.gov/citations/19910016425)
549
"""
550
@inline function splitting_vanleer_haenel(u, orientation::Integer,
551
equations::CompressibleEulerEquations1D)
552
fm = splitting_vanleer_haenel(u, Val{:minus}(), orientation, equations)
553
fp = splitting_vanleer_haenel(u, Val{:plus}(), orientation, equations)
554
return fm, fp
555
end
556
557
@inline function splitting_vanleer_haenel(u, ::Val{:plus}, orientation::Integer,
558
equations::CompressibleEulerEquations1D)
559
rho, rho_v1, rho_e_total = u
560
v1 = rho_v1 / rho
561
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)
562
563
# sound speed and enthalpy
564
a = sqrt(equations.gamma * p / rho)
565
H = (rho_e_total + p) / rho
566
567
# signed Mach number
568
M = v1 / a
569
570
p_plus = 0.5f0 * (1 + equations.gamma * M) * p
571
572
f1p = 0.25f0 * rho * a * (M + 1)^2
573
f2p = f1p * v1 + p_plus
574
f3p = f1p * H
575
576
return SVector(f1p, f2p, f3p)
577
end
578
579
@inline function splitting_vanleer_haenel(u, ::Val{:minus}, orientation::Integer,
580
equations::CompressibleEulerEquations1D)
581
rho, rho_v1, rho_e_total = u
582
v1 = rho_v1 / rho
583
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)
584
585
# sound speed and enthalpy
586
a = sqrt(equations.gamma * p / rho)
587
H = (rho_e_total + p) / rho
588
589
# signed Mach number
590
M = v1 / a
591
592
p_minus = 0.5f0 * (1 - equations.gamma * M) * p
593
594
f1m = -0.25f0 * rho * a * (M - 1)^2
595
f2m = f1m * v1 + p_minus
596
f3m = f1m * H
597
598
return SVector(f1m, f2m, f3m)
599
end
600
601
# TODO: FD
602
# This splitting is interesting because it can handle the "el diablo" wave
603
# for long time runs. Computing the eigenvalues of the operator we see
604
# J = jacobian_ad_forward(semi);
605
# lamb = eigvals(J);
606
# maximum(real.(lamb))
607
# 2.1411031631522748e-6
608
# So the instability of this splitting is very weak. However, the 2D variant
609
# of this splitting on "el diablo" still crashes early. Can we learn anything
610
# from the design of this splitting?
611
"""
612
splitting_coirier_vanleer(u, orientation::Integer,
613
equations::CompressibleEulerEquations1D)
614
splitting_coirier_vanleer(u, which::Union{Val{:minus}, Val{:plus}}
615
orientation::Integer,
616
equations::CompressibleEulerEquations1D)
617
618
Splitting of the compressible Euler flux from Coirier and van Leer.
619
The splitting has correction terms in the pressure splitting as well as
620
the mass and energy flux components. The motivation for these corrections
621
are to handle flows at the low Mach number limit.
622
623
Returns a tuple of the fluxes "minus" (associated with waves going into the
624
negative axis direction) and "plus" (associated with waves going into the
625
positive axis direction). If only one of the fluxes is required, use the
626
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
627
628
!!! warning "Experimental implementation (upwind SBP)"
629
This is an experimental feature and may change in future releases.
630
631
## References
632
633
- William Coirier and Bram van Leer (1991)
634
Numerical flux formulas for the Euler and Navier-Stokes equations.
635
II - Progress in flux-vector splitting
636
[DOI: 10.2514/6.1991-1566](https://doi.org/10.2514/6.1991-1566)
637
"""
638
@inline function splitting_coirier_vanleer(u, orientation::Integer,
639
equations::CompressibleEulerEquations1D)
640
fm = splitting_coirier_vanleer(u, Val{:minus}(), orientation, equations)
641
fp = splitting_coirier_vanleer(u, Val{:plus}(), orientation, equations)
642
return fm, fp
643
end
644
645
@inline function splitting_coirier_vanleer(u, ::Val{:plus}, orientation::Integer,
646
equations::CompressibleEulerEquations1D)
647
rho, rho_v1, rho_e_total = u
648
v1 = rho_v1 / rho
649
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)
650
651
# sound speed and enthalpy
652
a = sqrt(equations.gamma * p / rho)
653
H = (rho_e_total + p) / rho
654
655
# signed Mach number
656
M = v1 / a
657
658
P = 2
659
mu = 1
660
nu = 0.75f0
661
omega = 2 # adjusted from suggested value of 1.5
662
663
p_plus = 0.25f0 * ((M + 1)^2 * (2 - M) - nu * M * (M^2 - 1)^P) * p
664
665
f1p = 0.25f0 * rho * a * ((M + 1)^2 - mu * (M^2 - 1)^P)
666
f2p = f1p * v1 + p_plus
667
f3p = f1p * H - omega * rho * a^3 * M^2 * (M^2 - 1)^2
668
669
return SVector(f1p, f2p, f3p)
670
end
671
672
@inline function splitting_coirier_vanleer(u, ::Val{:minus}, orientation::Integer,
673
equations::CompressibleEulerEquations1D)
674
rho, rho_v1, rho_e_total = u
675
v1 = rho_v1 / rho
676
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)
677
678
# sound speed and enthalpy
679
a = sqrt(equations.gamma * p / rho)
680
H = (rho_e_total + p) / rho
681
682
# signed Mach number
683
M = v1 / a
684
685
P = 2
686
mu = 1
687
nu = 0.75f0
688
omega = 2 # adjusted from suggested value of 1.5
689
690
p_minus = 0.25f0 * ((M - 1)^2 * (2 + M) + nu * M * (M^2 - 1)^P) * p
691
692
f1m = -0.25f0 * rho * a * ((M - 1)^2 - mu * (M^2 - 1)^P)
693
f2m = f1m * v1 + p_minus
694
f3m = f1m * H + omega * rho * a^3 * M^2 * (M^2 - 1)^2
695
696
return SVector(f1m, f2m, f3m)
697
end
698
699
# Calculate estimates for maximum wave speed for local Lax-Friedrichs-type dissipation as the
700
# maximum velocity magnitude plus the maximum speed of sound
701
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
702
equations::CompressibleEulerEquations1D)
703
rho_ll, rho_v1_ll, rho_e_total_ll = u_ll
704
rho_rr, rho_v1_rr, rho_e_total_rr = u_rr
705
706
# Calculate primitive variables and speed of sound
707
v1_ll = rho_v1_ll / rho_ll
708
v_mag_ll = abs(v1_ll)
709
p_ll = (equations.gamma - 1) * (rho_e_total_ll - 0.5f0 * rho_ll * v_mag_ll^2)
710
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
711
v1_rr = rho_v1_rr / rho_rr
712
v_mag_rr = abs(v1_rr)
713
p_rr = (equations.gamma - 1) * (rho_e_total_rr - 0.5f0 * rho_rr * v_mag_rr^2)
714
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
715
716
return max(v_mag_ll, v_mag_rr) + max(c_ll, c_rr)
717
end
718
719
@inline function (dissipation::DissipationMatrixWintersEtal)(u_ll, u_rr,
720
normal_direction::AbstractVector,
721
equations::CompressibleEulerEquations1D)
722
gamma = equations.gamma
723
724
norm_ = norm(normal_direction)
725
unit_normal_direction = normal_direction / norm_
726
727
rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)
728
rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)
729
730
b_ll = rho_ll / (2 * p_ll)
731
b_rr = rho_rr / (2 * p_rr)
732
733
rho_log = ln_mean(rho_ll, rho_rr)
734
b_log = ln_mean(b_ll, b_rr)
735
v1_avg = 0.5f0 * (v1_ll + v1_rr)
736
p_avg = 0.5f0 * (rho_ll + rho_rr) / (b_ll + b_rr) # 2 * b_avg = b_ll + b_rr
737
v1_squared_bar = v1_ll * v1_rr
738
h_bar = gamma / (2 * b_log * (gamma - 1)) + 0.5f0 * v1_squared_bar
739
c_bar = sqrt(gamma * p_avg / rho_log)
740
741
v1_minus_c = v1_avg - c_bar * unit_normal_direction[1]
742
v1_plus_c = v1_avg + c_bar * unit_normal_direction[1]
743
v_avg_normal = v1_avg * unit_normal_direction[1]
744
745
entropy_vars_jump = cons2entropy(u_rr, equations) - cons2entropy(u_ll, equations)
746
747
lambda_1 = abs(v_avg_normal - c_bar) * rho_log / (2 * gamma)
748
lambda_2 = abs(v_avg_normal) * rho_log * (gamma - 1) / gamma
749
lambda_3 = abs(v_avg_normal + c_bar) * rho_log / (2 * gamma)
750
lambda_4 = abs(v_avg_normal) * p_avg
751
752
entropy_var_rho_jump, entropy_var_rho_v1_jump, entropy_var_rho_e_total_jump = entropy_vars_jump
753
754
w1 = lambda_1 * (entropy_var_rho_jump + v1_minus_c * entropy_var_rho_v1_jump +
755
(h_bar - c_bar * v_avg_normal) * entropy_var_rho_e_total_jump)
756
w2 = lambda_2 * (entropy_var_rho_jump + v1_avg * entropy_var_rho_v1_jump +
757
v1_squared_bar / 2 * entropy_var_rho_e_total_jump)
758
w3 = lambda_3 * (entropy_var_rho_jump + v1_plus_c * entropy_var_rho_v1_jump +
759
(h_bar + c_bar * v_avg_normal) * entropy_var_rho_e_total_jump)
760
761
dissipation_rho = w1 + w2 + w3
762
763
dissipation_rho_v1 = (w1 * v1_minus_c +
764
w2 * v1_avg +
765
w3 * v1_plus_c)
766
767
dissipation_rhoe = (w1 * (h_bar - c_bar * v_avg_normal) +
768
w2 * 0.5f0 * v1_squared_bar +
769
w3 * (h_bar + c_bar * v_avg_normal) +
770
lambda_4 *
771
(entropy_var_rho_e_total_jump *
772
(v1_avg * v1_avg - v_avg_normal^2)))
773
774
return -0.5f0 * SVector(dissipation_rho, dissipation_rho_v1, dissipation_rhoe) *
775
norm_
776
end
777
778
# Calculate estimates for minimum and maximum wave speeds for HLL-type fluxes
779
@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,
780
equations::CompressibleEulerEquations1D)
781
rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)
782
rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)
783
784
λ_min = v1_ll - sqrt(equations.gamma * p_ll / rho_ll)
785
λ_max = v1_rr + sqrt(equations.gamma * p_rr / rho_rr)
786
787
return λ_min, λ_max
788
end
789
790
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
791
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
792
equations::CompressibleEulerEquations1D)
793
rho_ll, rho_v1_ll, rho_e_total_ll = u_ll
794
rho_rr, rho_v1_rr, rho_e_total_rr = u_rr
795
796
# Calculate primitive variables and speed of sound
797
v1_ll = rho_v1_ll / rho_ll
798
v_mag_ll = abs(v1_ll)
799
p_ll = (equations.gamma - 1) * (rho_e_total_ll - 0.5f0 * rho_ll * v_mag_ll^2)
800
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
801
v1_rr = rho_v1_rr / rho_rr
802
v_mag_rr = abs(v1_rr)
803
p_rr = (equations.gamma - 1) * (rho_e_total_rr - 0.5f0 * rho_rr * v_mag_rr^2)
804
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
805
806
return max(v_mag_ll + c_ll, v_mag_rr + c_rr)
807
end
808
809
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
810
@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,
811
equations::CompressibleEulerEquations1D)
812
rho_ll, v1_ll, p_ll = cons2prim(u_ll, equations)
813
rho_rr, v1_rr, p_rr = cons2prim(u_rr, equations)
814
815
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
816
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
817
818
λ_min = min(v1_ll - c_ll, v1_rr - c_rr)
819
λ_max = max(v1_ll + c_ll, v1_rr + c_rr)
820
821
return λ_min, λ_max
822
end
823
824
"""
825
flux_hllc(u_ll, u_rr, orientation, equations::CompressibleEulerEquations1D)
826
827
Computes the HLLC flux (HLL with Contact) for compressible Euler equations developed by E.F. Toro
828
[Lecture slides](http://www.prague-sum.com/download/2012/Toro_2-HLLC-RiemannSolver.pdf)
829
Signal speeds: [DOI: 10.1137/S1064827593260140](https://doi.org/10.1137/S1064827593260140)
830
"""
831
function flux_hllc(u_ll, u_rr, orientation::Integer,
832
equations::CompressibleEulerEquations1D)
833
# Calculate primitive variables and speed of sound
834
rho_ll, rho_v1_ll, rho_e_total_ll = u_ll
835
rho_rr, rho_v1_rr, rho_e_total_rr = u_rr
836
837
v1_ll = rho_v1_ll / rho_ll
838
e_ll = rho_e_total_ll / rho_ll
839
p_ll = (equations.gamma - 1) * (rho_e_total_ll - 0.5f0 * rho_ll * v1_ll^2)
840
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
841
842
v1_rr = rho_v1_rr / rho_rr
843
e_rr = rho_e_total_rr / rho_rr
844
p_rr = (equations.gamma - 1) * (rho_e_total_rr - 0.5f0 * rho_rr * v1_rr^2)
845
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
846
847
# Obtain left and right fluxes
848
f_ll = flux(u_ll, orientation, equations)
849
f_rr = flux(u_rr, orientation, equations)
850
851
# Compute Roe averages
852
sqrt_rho_ll = sqrt(rho_ll)
853
sqrt_rho_rr = sqrt(rho_rr)
854
sum_sqrt_rho = sqrt_rho_ll + sqrt_rho_rr
855
vel_L = v1_ll
856
vel_R = v1_rr
857
vel_roe = (sqrt_rho_ll * vel_L + sqrt_rho_rr * vel_R) / sum_sqrt_rho
858
ekin_roe = 0.5f0 * vel_roe^2
859
H_ll = (rho_e_total_ll + p_ll) / rho_ll
860
H_rr = (rho_e_total_rr + p_rr) / rho_rr
861
H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) / sum_sqrt_rho
862
c_roe = sqrt((equations.gamma - 1) * (H_roe - ekin_roe))
863
864
Ssl = min(vel_L - c_ll, vel_roe - c_roe)
865
Ssr = max(vel_R + c_rr, vel_roe + c_roe)
866
sMu_L = Ssl - vel_L
867
sMu_R = Ssr - vel_R
868
if Ssl >= 0
869
f1 = f_ll[1]
870
f2 = f_ll[2]
871
f3 = f_ll[3]
872
elseif Ssr <= 0
873
f1 = f_rr[1]
874
f2 = f_rr[2]
875
f3 = f_rr[3]
876
else
877
SStar = (p_rr - p_ll + rho_ll * vel_L * sMu_L - rho_rr * vel_R * sMu_R) /
878
(rho_ll * sMu_L - rho_rr * sMu_R)
879
if Ssl <= 0 <= SStar
880
densStar = rho_ll * sMu_L / (Ssl - SStar)
881
enerStar = e_ll + (SStar - vel_L) * (SStar + p_ll / (rho_ll * sMu_L))
882
UStar1 = densStar
883
UStar2 = densStar * SStar
884
UStar3 = densStar * enerStar
885
886
f1 = f_ll[1] + Ssl * (UStar1 - rho_ll)
887
f2 = f_ll[2] + Ssl * (UStar2 - rho_v1_ll)
888
f3 = f_ll[3] + Ssl * (UStar3 - rho_e_total_ll)
889
else
890
densStar = rho_rr * sMu_R / (Ssr - SStar)
891
enerStar = e_rr + (SStar - vel_R) * (SStar + p_rr / (rho_rr * sMu_R))
892
UStar1 = densStar
893
UStar2 = densStar * SStar
894
UStar3 = densStar * enerStar
895
896
#end
897
f1 = f_rr[1] + Ssr * (UStar1 - rho_rr)
898
f2 = f_rr[2] + Ssr * (UStar2 - rho_v1_rr)
899
f3 = f_rr[3] + Ssr * (UStar3 - rho_e_total_rr)
900
end
901
end
902
return SVector(f1, f2, f3)
903
end
904
905
# While `normal_direction` isn't strictly necessary in 1D, certain solvers assume that
906
# the normal component is incorporated into the numerical flux.
907
#
908
# The HLLC flux along a 1D "normal" can be evaluated by scaling the velocity/momentum by
909
# the normal for the 1D HLLC flux, then scaling the resulting momentum flux again.
910
# Moreover, the 2D HLLC flux reduces to this if the normal vector is [n, 0].
911
function flux_hllc(u_ll, u_rr, normal_direction::AbstractVector,
912
equations::CompressibleEulerEquations1D)
913
norm_ = abs(normal_direction[1])
914
normal_direction_unit = normal_direction[1] * inv(norm_)
915
916
# scale the momentum by the normal direction
917
f = flux_hllc(SVector(u_ll[1], normal_direction_unit * u_ll[2], u_ll[3]),
918
SVector(u_rr[1], normal_direction_unit * u_rr[2], u_rr[3]), 1,
919
equations)
920
921
# rescale the momentum flux by the normal direction and normalize
922
return SVector(f[1], normal_direction_unit * f[2], f[3]) * norm_
923
end
924
925
"""
926
min_max_speed_einfeldt(u_ll, u_rr, orientation, equations::CompressibleEulerEquations1D)
927
928
Computes the HLLE (Harten-Lax-van Leer-Einfeldt) flux for the compressible Euler equations.
929
Special estimates of the signal velocites and linearization of the Riemann problem developed
930
by Einfeldt to ensure that the internal energy and density remain positive during the computation
931
of the numerical flux.
932
933
Original publication:
934
- Bernd Einfeldt (1988)
935
On Godunov-type methods for gas dynamics.
936
[DOI: 10.1137/0725021](https://doi.org/10.1137/0725021)
937
938
Compactly summarized:
939
- Siddhartha Mishra, Ulrik Skre Fjordholm and Rémi Abgrall
940
Numerical methods for conservation laws and related equations.
941
[Link](https://metaphor.ethz.ch/x/2019/hs/401-4671-00L/literature/mishra_hyperbolic_pdes.pdf)
942
"""
943
@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,
944
equations::CompressibleEulerEquations1D)
945
# Calculate primitive variables, enthalpy and speed of sound
946
rho_ll, v_ll, p_ll = cons2prim(u_ll, equations)
947
rho_rr, v_rr, p_rr = cons2prim(u_rr, equations)
948
949
# `u_ll[3]` is total energy `rho_e_total_ll` on the left
950
H_ll = (u_ll[3] + p_ll) / rho_ll
951
c_ll = sqrt(equations.gamma * p_ll / rho_ll)
952
953
# `u_rr[3]` is total energy `rho_e_total_rr` on the right
954
H_rr = (u_rr[3] + p_rr) / rho_rr
955
c_rr = sqrt(equations.gamma * p_rr / rho_rr)
956
957
# Compute Roe averages
958
sqrt_rho_ll = sqrt(rho_ll)
959
sqrt_rho_rr = sqrt(rho_rr)
960
inv_sum_sqrt_rho = inv(sqrt_rho_ll + sqrt_rho_rr)
961
962
v_roe = (sqrt_rho_ll * v_ll + sqrt_rho_rr * v_rr) * inv_sum_sqrt_rho
963
v_roe_mag = v_roe^2
964
965
H_roe = (sqrt_rho_ll * H_ll + sqrt_rho_rr * H_rr) * inv_sum_sqrt_rho
966
c_roe = sqrt((equations.gamma - 1) * (H_roe - 0.5f0 * v_roe_mag))
967
968
# Compute convenience constant for positivity preservation, see
969
# https://doi.org/10.1016/0021-9991(91)90211-3
970
beta = sqrt(0.5f0 * (equations.gamma - 1) / equations.gamma)
971
972
# Estimate the edges of the Riemann fan (with positivity conservation)
973
SsL = min(v_roe - c_roe, v_ll - beta * c_ll, 0)
974
SsR = max(v_roe + c_roe, v_rr + beta * c_rr, 0)
975
976
return SsL, SsR
977
end
978
979
@inline function max_abs_speeds(u, equations::CompressibleEulerEquations1D)
980
rho, rho_v1, rho_e_total = u
981
v1 = rho_v1 / rho
982
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho * v1^2)
983
c = sqrt(equations.gamma * p / rho)
984
985
return (abs(v1) + c,)
986
end
987
988
# Convert conservative variables to primitive
989
@inline function cons2prim(u, equations::CompressibleEulerEquations1D)
990
rho, rho_v1, rho_e_total = u
991
992
v1 = rho_v1 / rho
993
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho_v1 * v1)
994
995
return SVector(rho, v1, p)
996
end
997
998
# Convert conservative variables to entropy
999
@inline function cons2entropy(u, equations::CompressibleEulerEquations1D)
1000
rho, rho_v1, rho_e_total = u
1001
1002
v1 = rho_v1 / rho
1003
v_square = v1^2
1004
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * rho * v_square)
1005
s = log(p) - equations.gamma * log(rho)
1006
rho_p = rho / p
1007
1008
w1 = (equations.gamma - s) * equations.inv_gamma_minus_one -
1009
0.5f0 * rho_p * v_square
1010
w2 = rho_p * v1
1011
w3 = -rho_p
1012
1013
return SVector(w1, w2, w3)
1014
end
1015
1016
@inline function entropy2cons(w, equations::CompressibleEulerEquations1D)
1017
# See Hughes, Franca, Mallet (1986) A new finite element formulation for CFD
1018
# [DOI: 10.1016/0045-7825(86)90127-1](https://doi.org/10.1016/0045-7825(86)90127-1)
1019
@unpack gamma = equations
1020
1021
# convert to entropy `-rho * s` used by Hughes, France, Mallet (1986)
1022
# instead of `-rho * s / (gamma - 1)`
1023
V1, V2, V5 = w .* (gamma - 1)
1024
1025
# specific entropy, eq. (53)
1026
s = gamma - V1 + 0.5f0 * (V2^2) / V5
1027
1028
# eq. (52)
1029
energy_internal = ((gamma - 1) / (-V5)^gamma)^(equations.inv_gamma_minus_one) *
1030
exp(-s * equations.inv_gamma_minus_one)
1031
1032
# eq. (51)
1033
rho = -V5 * energy_internal
1034
rho_v1 = V2 * energy_internal
1035
rho_e_total = (1 - 0.5f0 * (V2^2) / V5) * energy_internal
1036
return SVector(rho, rho_v1, rho_e_total)
1037
end
1038
1039
# Convert primitive to conservative variables
1040
@inline function prim2cons(prim, equations::CompressibleEulerEquations1D)
1041
rho, v1, p = prim
1042
rho_v1 = rho * v1
1043
rho_e_total = p * equations.inv_gamma_minus_one + 0.5f0 * (rho_v1 * v1)
1044
return SVector(rho, rho_v1, rho_e_total)
1045
end
1046
1047
@inline function density(u, equations::CompressibleEulerEquations1D)
1048
rho = u[1]
1049
return rho
1050
end
1051
1052
@inline function velocity(u, orientation_or_normal,
1053
equations::CompressibleEulerEquations1D)
1054
return velocity(u, equations)
1055
end
1056
1057
@inline function velocity(u, equations::CompressibleEulerEquations1D)
1058
rho = u[1]
1059
v1 = u[2] / rho
1060
return v1
1061
end
1062
1063
@doc raw"""
1064
pressure(u, equations::AbstractCompressibleEulerEquations)
1065
1066
Computes the pressure for an ideal equation of state with
1067
isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`.
1068
```math
1069
\begin{aligned}
1070
p &= (\gamma - 1) \left( E_\mathrm{tot} - E_\mathrm{kin} \right) \\
1071
&= (\gamma - 1) \left( \rho e - \frac{1}{2}\rho \Vert v \Vert_2^2 \right)
1072
\end{aligned}
1073
```
1074
"""
1075
@inline function pressure(u, equations::CompressibleEulerEquations1D)
1076
rho, rho_v1, rho_e_total = u
1077
p = (equations.gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1^2) / rho)
1078
return p
1079
end
1080
1081
@inline function density_pressure(u, equations::CompressibleEulerEquations1D)
1082
rho, rho_v1, rho_e_total = u
1083
rho_times_p = (equations.gamma - 1) * (rho * rho_e_total - 0.5f0 * (rho_v1^2))
1084
return rho_times_p
1085
end
1086
1087
@doc raw"""
1088
entropy_thermodynamic(cons, equations::AbstractCompressibleEulerEquations)
1089
1090
Calculate thermodynamic entropy for a conservative state `cons` as
1091
1092
```math
1093
s = \log(p) - \gamma \log(\rho)
1094
```
1095
"""
1096
@inline function entropy_thermodynamic(cons, equations::CompressibleEulerEquations1D)
1097
# Pressure
1098
p = (equations.gamma - 1) * (cons[3] - 0.5f0 * (cons[2]^2) / cons[1])
1099
1100
# Thermodynamic entropy
1101
s = log(p) - equations.gamma * log(cons[1])
1102
1103
return s
1104
end
1105
1106
@doc raw"""
1107
entropy_math(cons, equations::AbstractCompressibleEulerEquations)
1108
1109
Calculate mathematical entropy for a conservative state `cons` as
1110
```math
1111
S = -\frac{\rho s}{\gamma - 1}
1112
```
1113
where `s` is the thermodynamic entropy calculated by
1114
[`entropy_thermodynamic(cons, equations::AbstractCompressibleEulerEquations)`](@ref).
1115
"""
1116
@inline function entropy_math(cons, equations::CompressibleEulerEquations1D)
1117
# Mathematical entropy
1118
S = -entropy_thermodynamic(cons, equations) * cons[1] *
1119
equations.inv_gamma_minus_one
1120
1121
return S
1122
end
1123
1124
"""
1125
entropy(cons, equations::AbstractCompressibleEulerEquations)
1126
1127
Default entropy is the mathematical entropy
1128
[`entropy_math(cons, equations::AbstractCompressibleEulerEquations)`](@ref).
1129
"""
1130
@inline function entropy(cons, equations::CompressibleEulerEquations1D)
1131
return entropy_math(cons, equations)
1132
end
1133
1134
# Calculate total energy for a conservative state `cons`
1135
@inline energy_total(cons, ::CompressibleEulerEquations1D) = cons[3]
1136
1137
# Calculate kinetic energy for a conservative state `cons`
1138
@inline function energy_kinetic(cons, equations::CompressibleEulerEquations1D)
1139
return 0.5f0 * (cons[2]^2) / cons[1]
1140
end
1141
1142
"""
1143
energy_internal(cons, equations::AbstractCompressibleEulerEquations)
1144
1145
Calculate internal energy ``e_{\text{internal}}`` for a conservative state `cons` as the difference
1146
of total energy and kinetic energy.
1147
"""
1148
@inline function energy_internal(cons, equations::CompressibleEulerEquations1D)
1149
return energy_total(cons, equations) - energy_kinetic(cons, equations)
1150
end
1151
1152
@doc raw"""
1153
entropy_potential(u, orientation_or_normal_direction,
1154
equations::AbstractCompressibleEulerEquations)
1155
1156
Calculate the entropy potential, which for the compressible Euler equations is simply
1157
the momentum for the choice of mathematical [`entropy`](@ref) ``S(u) = -\frac{\rho s}{\gamma - 1}``
1158
with thermodynamic entropy ``s = \ln(p) - \gamma \ln(\rho)``.
1159
1160
## References
1161
- Eitan Tadmor (2003)
1162
Entropy stability theory for difference approximations of nonlinear conservation laws and related time-dependent problems
1163
[DOI: 10.1017/S0962492902000156](https://doi.org/10.1017/S0962492902000156)
1164
"""
1165
@inline function entropy_potential(u, orientation::Int,
1166
equations::CompressibleEulerEquations1D)
1167
return u[2]
1168
end
1169
end # @muladd
1170
1171