Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/equations/compressible_euler_multicomponent_2d.jl
2055 views
1
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).
2
# Since these FMAs can increase the performance of many numerical algorithms,
3
# we need to opt-in explicitly.
4
# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.
5
@muladd begin
6
#! format: noindent
7
8
@doc raw"""
9
CompressibleEulerMulticomponentEquations2D(; gammas, gas_constants)
10
11
Multicomponent version of the compressible Euler equations
12
```math
13
\frac{\partial}{\partial t}
14
\begin{pmatrix}
15
\rho v_1 \\ \rho v_2 \\ \rho e \\ \rho_1 \\ \rho_2 \\ \vdots \\ \rho_{n}
16
\end{pmatrix}
17
+
18
\frac{\partial}{\partial x}
19
\begin{pmatrix}
20
\rho v_1^2 + p \\ \rho v_1 v_2 \\ ( \rho e +p) v_1 \\ \rho_1 v_1 \\ \rho_2 v_1 \\ \vdots \\ \rho_{n} v_1
21
\end{pmatrix}
22
+
23
\frac{\partial}{\partial y}
24
\begin{pmatrix}
25
\rho v_1 v_2 \\ \rho v_2^2 + p \\ ( \rho e +p) v_2 \\ \rho_1 v_2 \\ \rho_2 v_2 \\ \vdots \\ \rho_{n} v_2
26
\end{pmatrix}
27
=
28
\begin{pmatrix}
29
0 \\ 0 \\ 0 \\ 0 \\ 0 \\ \vdots \\ 0
30
\end{pmatrix}
31
```
32
for calorically perfect gas in two space dimensions.
33
Here, ``\rho_i`` is the density of component ``i``, ``\rho=\sum_{i=1}^n\rho_i`` the sum of the individual ``\rho_i``,
34
``v_1``, ``v_2`` the velocities, ``e`` the specific total energy **rather than** specific internal energy, and
35
```math
36
p = (\gamma - 1) \left( \rho e - \frac{1}{2} \rho (v_1^2 + v_2^2) \right)
37
```
38
the pressure,
39
```math
40
\gamma=\frac{\sum_{i=1}^n\rho_i C_{v,i}\gamma_i}{\sum_{i=1}^n\rho_i C_{v,i}}
41
```
42
total heat capacity ratio, ``\gamma_i`` heat capacity ratio of component ``i``,
43
```math
44
C_{v,i}=\frac{R}{\gamma_i-1}
45
```
46
specific heat capacity at constant volume of component ``i``.
47
48
In case of more than one component, the specific heat ratios `gammas` and the gas constants
49
`gas_constants` in [kJ/(kg*K)] should be passed as tuples, e.g., `gammas=(1.4, 1.667)`.
50
51
The remaining variables like the specific heats at constant volume `cv` or the specific heats at
52
constant pressure `cp` are then calculated considering a calorically perfect gas.
53
"""
54
struct CompressibleEulerMulticomponentEquations2D{NVARS, NCOMP, RealT <: Real} <:
55
AbstractCompressibleEulerMulticomponentEquations{2, NVARS, NCOMP}
56
gammas::SVector{NCOMP, RealT}
57
gas_constants::SVector{NCOMP, RealT}
58
cv::SVector{NCOMP, RealT}
59
cp::SVector{NCOMP, RealT}
60
61
function CompressibleEulerMulticomponentEquations2D{NVARS, NCOMP, RealT}(gammas::SVector{NCOMP,
62
RealT},
63
gas_constants::SVector{NCOMP,
64
RealT}) where {
65
NVARS,
66
NCOMP,
67
RealT <:
68
Real
69
}
70
NCOMP >= 1 ||
71
throw(DimensionMismatch("`gammas` and `gas_constants` have to be filled with at least one value"))
72
73
cv = gas_constants ./ (gammas .- 1)
74
cp = gas_constants + gas_constants ./ (gammas .- 1)
75
76
new(gammas, gas_constants, cv, cp)
77
end
78
end
79
80
function CompressibleEulerMulticomponentEquations2D(; gammas, gas_constants)
81
_gammas = promote(gammas...)
82
_gas_constants = promote(gas_constants...)
83
RealT = promote_type(eltype(_gammas), eltype(_gas_constants),
84
typeof(gas_constants[1] / (gammas[1] - 1)))
85
__gammas = SVector(map(RealT, _gammas))
86
__gas_constants = SVector(map(RealT, _gas_constants))
87
88
NVARS = length(_gammas) + 3
89
NCOMP = length(_gammas)
90
91
return CompressibleEulerMulticomponentEquations2D{NVARS, NCOMP, RealT}(__gammas,
92
__gas_constants)
93
end
94
95
@inline function Base.real(::CompressibleEulerMulticomponentEquations2D{NVARS, NCOMP,
96
RealT}) where {
97
NVARS,
98
NCOMP,
99
RealT
100
}
101
RealT
102
end
103
104
function varnames(::typeof(cons2cons),
105
equations::CompressibleEulerMulticomponentEquations2D)
106
cons = ("rho_v1", "rho_v2", "rho_e")
107
rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))
108
return (cons..., rhos...)
109
end
110
111
function varnames(::typeof(cons2prim),
112
equations::CompressibleEulerMulticomponentEquations2D)
113
prim = ("v1", "v2", "p")
114
rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))
115
return (prim..., rhos...)
116
end
117
118
# Set initial conditions at physical location `x` for time `t`
119
120
"""
121
initial_condition_convergence_test(x, t, equations::CompressibleEulerMulticomponentEquations2D)
122
123
A smooth initial condition used for convergence tests in combination with
124
[`source_terms_convergence_test`](@ref)
125
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
126
"""
127
function initial_condition_convergence_test(x, t,
128
equations::CompressibleEulerMulticomponentEquations2D)
129
RealT = eltype(x)
130
c = 2
131
A = convert(RealT, 0.1)
132
L = 2
133
f = 1.0f0 / L
134
omega = 2 * convert(RealT, pi) * f
135
ini = c + A * sin(omega * (x[1] + x[2] - t))
136
137
v1 = 1
138
v2 = 1
139
140
rho = ini
141
142
# Here we compute an arbitrary number of different rhos. (one rho is double the next rho while the sum of all rhos is 1)
143
prim_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) * (1 - 2) *
144
rho / (1 -
145
2^ncomponents(equations))
146
for i in eachcomponent(equations))
147
148
prim1 = rho * v1
149
prim2 = rho * v2
150
prim3 = rho^2
151
152
prim_other = SVector(prim1, prim2, prim3)
153
154
return vcat(prim_other, prim_rho)
155
end
156
157
"""
158
source_terms_convergence_test(u, x, t, equations::CompressibleEulerMulticomponentEquations2D)
159
160
Source terms used for convergence tests in combination with
161
[`initial_condition_convergence_test`](@ref)
162
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
163
"""
164
@inline function source_terms_convergence_test(u, x, t,
165
equations::CompressibleEulerMulticomponentEquations2D)
166
# Same settings as in `initial_condition`
167
RealT = eltype(u)
168
c = 2
169
A = convert(RealT, 0.1)
170
L = 2
171
f = 1.0f0 / L
172
omega = 2 * convert(RealT, pi) * f
173
174
gamma = totalgamma(u, equations)
175
176
x1, x2 = x
177
si, co = sincos((x1 + x2 - t) * omega)
178
tmp1 = co * A * omega
179
tmp2 = si * A
180
tmp3 = gamma - 1
181
tmp4 = (2 * c - 1) * tmp3
182
tmp5 = (2 * tmp2 * gamma - 2 * tmp2 + tmp4 + 1) * tmp1
183
tmp6 = tmp2 + c
184
185
# Here we compute an arbitrary number of different rhos. (one rho is double the next rho while the sum of all rhos is 1
186
du_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) * (1 - 2) *
187
tmp1 / (1 -
188
2^ncomponents(equations))
189
for i in eachcomponent(equations))
190
191
du1 = tmp5
192
du2 = tmp5
193
du3 = 2 * ((tmp6 - 1) * tmp3 + tmp6 * gamma) * tmp1
194
195
du_other = SVector(du1, du2, du3)
196
197
return vcat(du_other, du_rho)
198
end
199
200
"""
201
initial_condition_weak_blast_wave(x, t, equations::CompressibleEulerMulticomponentEquations2D)
202
203
A for multicomponent adapted weak blast wave taken from
204
- Sebastian Hennemann, Gregor J. Gassner (2020)
205
A provably entropy stable subcell shock capturing approach for high order split form DG
206
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
207
"""
208
function initial_condition_weak_blast_wave(x, t,
209
equations::CompressibleEulerMulticomponentEquations2D)
210
# From Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
211
# Set up polar coordinates
212
RealT = eltype(x)
213
inicenter = SVector(0, 0)
214
x_norm = x[1] - inicenter[1]
215
y_norm = x[2] - inicenter[2]
216
r = sqrt(x_norm^2 + y_norm^2)
217
phi = atan(y_norm, x_norm)
218
sin_phi, cos_phi = sincos(phi)
219
220
prim_rho = SVector{ncomponents(equations), real(equations)}(r > 0.5f0 ?
221
2^(i - 1) * (1 - 2) /
222
(RealT(1) -
223
2^ncomponents(equations)) :
224
2^(i - 1) * (1 - 2) *
225
RealT(1.1691) /
226
(1 -
227
2^ncomponents(equations))
228
for i in eachcomponent(equations))
229
230
v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos_phi
231
v2 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * sin_phi
232
p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)
233
234
prim_other = SVector(v1, v2, p)
235
236
return prim2cons(vcat(prim_other, prim_rho), equations)
237
end
238
239
# Calculate 1D flux for a single point
240
@inline function flux(u, orientation::Integer,
241
equations::CompressibleEulerMulticomponentEquations2D)
242
rho_v1, rho_v2, rho_e = u
243
244
rho = density(u, equations)
245
246
v1 = rho_v1 / rho
247
v2 = rho_v2 / rho
248
gamma = totalgamma(u, equations)
249
p = (gamma - 1) * (rho_e - 0.5f0 * rho * (v1^2 + v2^2))
250
251
if orientation == 1
252
f_rho = densities(u, v1, equations)
253
f1 = rho_v1 * v1 + p
254
f2 = rho_v2 * v1
255
f3 = (rho_e + p) * v1
256
else
257
f_rho = densities(u, v2, equations)
258
f1 = rho_v1 * v2
259
f2 = rho_v2 * v2 + p
260
f3 = (rho_e + p) * v2
261
end
262
263
f_other = SVector(f1, f2, f3)
264
265
return vcat(f_other, f_rho)
266
end
267
268
# Calculate 1D flux for a single point
269
@inline function flux(u, normal_direction::AbstractVector,
270
equations::CompressibleEulerMulticomponentEquations2D)
271
rho_v1, rho_v2, rho_e = u
272
273
rho = density(u, equations)
274
275
v1 = rho_v1 / rho
276
v2 = rho_v2 / rho
277
v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]
278
gamma = totalgamma(u, equations)
279
p = (gamma - 1) * (rho_e - 0.5f0 * rho * (v1^2 + v2^2))
280
281
f_rho = densities(u, v_normal, equations)
282
f1 = rho_v1 * v_normal + p * normal_direction[1]
283
f2 = rho_v2 * v_normal + p * normal_direction[2]
284
f3 = (rho_e + p) * v_normal
285
286
f_other = SVector(f1, f2, f3)
287
288
return vcat(f_other, f_rho)
289
end
290
291
"""
292
flux_chandrashekar(u_ll, u_rr, orientation, equations::CompressibleEulerMulticomponentEquations2D)
293
294
Adaption of the entropy conserving two-point flux by
295
- Ayoub Gouasmi, Karthik Duraisamy (2020)
296
"Formulation of Entropy-Stable schemes for the multicomponent compressible Euler equations"
297
[arXiv:1904.00972v3](https://arxiv.org/abs/1904.00972) [math.NA] 4 Feb 2020
298
"""
299
@inline function flux_chandrashekar(u_ll, u_rr, orientation::Integer,
300
equations::CompressibleEulerMulticomponentEquations2D)
301
# Unpack left and right state
302
@unpack gammas, gas_constants, cv = equations
303
rho_v1_ll, rho_v2_ll, rho_e_ll = u_ll
304
rho_v1_rr, rho_v2_rr, rho_e_rr = u_rr
305
rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 3],
306
u_rr[i + 3])
307
for i in eachcomponent(equations))
308
rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 3] +
309
u_rr[i + 3])
310
for i in eachcomponent(equations))
311
312
# Iterating over all partial densities
313
rho_ll = density(u_ll, equations)
314
rho_rr = density(u_rr, equations)
315
316
# extract velocities
317
v1_ll = rho_v1_ll / rho_ll
318
v2_ll = rho_v2_ll / rho_ll
319
v1_rr = rho_v1_rr / rho_rr
320
v2_rr = rho_v2_rr / rho_rr
321
v1_avg = 0.5f0 * (v1_ll + v1_rr)
322
v2_avg = 0.5f0 * (v2_ll + v2_rr)
323
v1_square = 0.5f0 * (v1_ll^2 + v1_rr^2)
324
v2_square = 0.5f0 * (v2_ll^2 + v2_rr^2)
325
v_sum = v1_avg + v2_avg
326
327
RealT = eltype(u_ll)
328
enth = zero(RealT)
329
help1_ll = zero(RealT)
330
help1_rr = zero(RealT)
331
332
for i in eachcomponent(equations)
333
enth += rhok_avg[i] * gas_constants[i]
334
help1_ll += u_ll[i + 3] * cv[i]
335
help1_rr += u_rr[i + 3] * cv[i]
336
end
337
338
T_ll = (rho_e_ll - 0.5f0 * rho_ll * (v1_ll^2 + v2_ll^2)) / help1_ll
339
T_rr = (rho_e_rr - 0.5f0 * rho_rr * (v1_rr^2 + v2_rr^2)) / help1_rr
340
T = 0.5f0 * (1 / T_ll + 1 / T_rr)
341
T_log = ln_mean(1 / T_ll, 1 / T_rr)
342
343
# Calculate fluxes depending on orientation
344
help1 = zero(RealT)
345
help2 = zero(RealT)
346
if orientation == 1
347
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg
348
for i in eachcomponent(equations))
349
for i in eachcomponent(equations)
350
help1 += f_rho[i] * cv[i]
351
help2 += f_rho[i]
352
end
353
f1 = (help2) * v1_avg + enth / T
354
f2 = (help2) * v2_avg
355
f3 = (help1) / T_log - 0.5f0 * (v1_square + v2_square) * (help2) + v1_avg * f1 +
356
v2_avg * f2
357
else
358
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v2_avg
359
for i in eachcomponent(equations))
360
for i in eachcomponent(equations)
361
help1 += f_rho[i] * cv[i]
362
help2 += f_rho[i]
363
end
364
f1 = (help2) * v1_avg
365
f2 = (help2) * v2_avg + enth / T
366
f3 = (help1) / T_log - 0.5f0 * (v1_square + v2_square) * (help2) + v1_avg * f1 +
367
v2_avg * f2
368
end
369
f_other = SVector(f1, f2, f3)
370
371
return vcat(f_other, f_rho)
372
end
373
374
"""
375
flux_ranocha(u_ll, u_rr, orientation_or_normal_direction,
376
equations::CompressibleEulerMulticomponentEquations2D)
377
378
Adaption of the entropy conserving and kinetic energy preserving two-point flux by
379
- Hendrik Ranocha (2018)
380
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
381
for Hyperbolic Balance Laws
382
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
383
See also
384
- Hendrik Ranocha (2020)
385
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
386
the Euler Equations Using Summation-by-Parts Operators
387
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
388
"""
389
@inline function flux_ranocha(u_ll, u_rr, orientation::Integer,
390
equations::CompressibleEulerMulticomponentEquations2D)
391
# Unpack left and right state
392
@unpack gammas, gas_constants, cv = equations
393
rho_v1_ll, rho_v2_ll, rho_e_ll = u_ll
394
rho_v1_rr, rho_v2_rr, rho_e_rr = u_rr
395
rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 3],
396
u_rr[i + 3])
397
for i in eachcomponent(equations))
398
rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 3] +
399
u_rr[i + 3])
400
for i in eachcomponent(equations))
401
402
# Iterating over all partial densities
403
rho_ll = density(u_ll, equations)
404
rho_rr = density(u_rr, equations)
405
406
# Calculating gamma
407
gamma = totalgamma(0.5f0 * (u_ll + u_rr), equations)
408
inv_gamma_minus_one = 1 / (gamma - 1)
409
410
# extract velocities
411
v1_ll = rho_v1_ll / rho_ll
412
v1_rr = rho_v1_rr / rho_rr
413
v1_avg = 0.5f0 * (v1_ll + v1_rr)
414
v2_ll = rho_v2_ll / rho_ll
415
v2_rr = rho_v2_rr / rho_rr
416
v2_avg = 0.5f0 * (v2_ll + v2_rr)
417
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr)
418
419
# helpful variables
420
RealT = eltype(u_ll)
421
help1_ll = zero(RealT)
422
help1_rr = zero(RealT)
423
enth_ll = zero(RealT)
424
enth_rr = zero(RealT)
425
for i in eachcomponent(equations)
426
enth_ll += u_ll[i + 3] * gas_constants[i]
427
enth_rr += u_rr[i + 3] * gas_constants[i]
428
help1_ll += u_ll[i + 3] * cv[i]
429
help1_rr += u_rr[i + 3] * cv[i]
430
end
431
432
# temperature and pressure
433
T_ll = (rho_e_ll - 0.5f0 * rho_ll * (v1_ll^2 + v2_ll^2)) / help1_ll
434
T_rr = (rho_e_rr - 0.5f0 * rho_rr * (v1_rr^2 + v2_rr^2)) / help1_rr
435
p_ll = T_ll * enth_ll
436
p_rr = T_rr * enth_rr
437
p_avg = 0.5f0 * (p_ll + p_rr)
438
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
439
440
f_rho_sum = zero(RealT)
441
if orientation == 1
442
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg
443
for i in eachcomponent(equations))
444
for i in eachcomponent(equations)
445
f_rho_sum += f_rho[i]
446
end
447
f1 = f_rho_sum * v1_avg + p_avg
448
f2 = f_rho_sum * v2_avg
449
f3 = f_rho_sum * (velocity_square_avg + inv_rho_p_mean * inv_gamma_minus_one) +
450
0.5f0 * (p_ll * v1_rr + p_rr * v1_ll)
451
else
452
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v2_avg
453
for i in eachcomponent(equations))
454
for i in eachcomponent(equations)
455
f_rho_sum += f_rho[i]
456
end
457
f1 = f_rho_sum * v1_avg
458
f2 = f_rho_sum * v2_avg + p_avg
459
f3 = f_rho_sum * (velocity_square_avg + inv_rho_p_mean * inv_gamma_minus_one) +
460
0.5f0 * (p_ll * v2_rr + p_rr * v2_ll)
461
end
462
463
# momentum and energy flux
464
f_other = SVector(f1, f2, f3)
465
466
return vcat(f_other, f_rho)
467
end
468
469
@inline function flux_ranocha(u_ll, u_rr, normal_direction::AbstractVector,
470
equations::CompressibleEulerMulticomponentEquations2D)
471
# Unpack left and right state
472
@unpack gammas, gas_constants, cv = equations
473
rho_v1_ll, rho_v2_ll, rho_e_ll = u_ll
474
rho_v1_rr, rho_v2_rr, rho_e_rr = u_rr
475
rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 3],
476
u_rr[i + 3])
477
for i in eachcomponent(equations))
478
rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 3] +
479
u_rr[i + 3])
480
for i in eachcomponent(equations))
481
482
# Iterating over all partial densities
483
rho_ll = density(u_ll, equations)
484
rho_rr = density(u_rr, equations)
485
486
# Calculating gamma
487
gamma = totalgamma(0.5f0 * (u_ll + u_rr), equations)
488
inv_gamma_minus_one = 1 / (gamma - 1)
489
490
# extract velocities
491
v1_ll = rho_v1_ll / rho_ll
492
v1_rr = rho_v1_rr / rho_rr
493
v1_avg = 0.5f0 * (v1_ll + v1_rr)
494
v2_ll = rho_v2_ll / rho_ll
495
v2_rr = rho_v2_rr / rho_rr
496
v2_avg = 0.5f0 * (v2_ll + v2_rr)
497
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr)
498
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
499
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
500
501
# helpful variables
502
RealT = eltype(u_ll)
503
help1_ll = zero(RealT)
504
help1_rr = zero(RealT)
505
enth_ll = zero(RealT)
506
enth_rr = zero(RealT)
507
for i in eachcomponent(equations)
508
enth_ll += u_ll[i + 3] * gas_constants[i]
509
enth_rr += u_rr[i + 3] * gas_constants[i]
510
help1_ll += u_ll[i + 3] * cv[i]
511
help1_rr += u_rr[i + 3] * cv[i]
512
end
513
514
# temperature and pressure
515
T_ll = (rho_e_ll - 0.5f0 * rho_ll * (v1_ll^2 + v2_ll^2)) / help1_ll
516
T_rr = (rho_e_rr - 0.5f0 * rho_rr * (v1_rr^2 + v2_rr^2)) / help1_rr
517
p_ll = T_ll * enth_ll
518
p_rr = T_rr * enth_rr
519
p_avg = 0.5f0 * (p_ll + p_rr)
520
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
521
522
f_rho_sum = zero(RealT)
523
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * 0.5f0 *
524
(v_dot_n_ll + v_dot_n_rr)
525
for i in eachcomponent(equations))
526
for i in eachcomponent(equations)
527
f_rho_sum += f_rho[i]
528
end
529
f1 = f_rho_sum * v1_avg + p_avg * normal_direction[1]
530
f2 = f_rho_sum * v2_avg + p_avg * normal_direction[2]
531
f3 = f_rho_sum * (velocity_square_avg + inv_rho_p_mean * inv_gamma_minus_one) +
532
0.5f0 * (p_ll * v_dot_n_rr + p_rr * v_dot_n_ll)
533
534
# momentum and energy flux
535
f_other = SVector(f1, f2, f3)
536
537
return vcat(f_other, f_rho)
538
end
539
540
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation
541
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
542
equations::CompressibleEulerMulticomponentEquations2D)
543
rho_v1_ll, rho_v2_ll, rho_e_ll = u_ll
544
rho_v1_rr, rho_v2_rr, rho_e_rr = u_rr
545
546
# Get the density and gas gamma
547
rho_ll = density(u_ll, equations)
548
rho_rr = density(u_rr, equations)
549
gamma_ll = totalgamma(u_ll, equations)
550
gamma_rr = totalgamma(u_rr, equations)
551
552
# Get the velocities based on direction
553
if orientation == 1
554
v_ll = rho_v1_ll / rho_ll
555
v_rr = rho_v1_rr / rho_rr
556
else # orientation == 2
557
v_ll = rho_v2_ll / rho_ll
558
v_rr = rho_v2_rr / rho_rr
559
end
560
561
# Compute the sound speeds on the left and right
562
p_ll = (gamma_ll - 1) * (rho_e_ll - 0.5f0 * (rho_v1_ll^2 + rho_v2_ll^2) / rho_ll)
563
c_ll = sqrt(gamma_ll * p_ll / rho_ll)
564
p_rr = (gamma_rr - 1) * (rho_e_rr - 0.5f0 * (rho_v1_rr^2 + rho_v2_rr^2) / rho_rr)
565
c_rr = sqrt(gamma_rr * p_rr / rho_rr)
566
567
return max(abs(v_ll), abs(v_rr)) + max(c_ll, c_rr)
568
end
569
570
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
571
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
572
equations::CompressibleEulerMulticomponentEquations2D)
573
rho_v1_ll, rho_v2_ll, rho_e_ll = u_ll
574
rho_v1_rr, rho_v2_rr, rho_e_rr = u_rr
575
576
# Get the density and gas gamma
577
rho_ll = density(u_ll, equations)
578
rho_rr = density(u_rr, equations)
579
gamma_ll = totalgamma(u_ll, equations)
580
gamma_rr = totalgamma(u_rr, equations)
581
582
# Get the velocities based on direction
583
if orientation == 1
584
v_ll = rho_v1_ll / rho_ll
585
v_rr = rho_v1_rr / rho_rr
586
else # orientation == 2
587
v_ll = rho_v2_ll / rho_ll
588
v_rr = rho_v2_rr / rho_rr
589
end
590
591
# Compute the sound speeds on the left and right
592
p_ll = (gamma_ll - 1) * (rho_e_ll - 0.5f0 * (rho_v1_ll^2 + rho_v2_ll^2) / rho_ll)
593
c_ll = sqrt(gamma_ll * p_ll / rho_ll)
594
p_rr = (gamma_rr - 1) * (rho_e_rr - 0.5f0 * (rho_v1_rr^2 + rho_v2_rr^2) / rho_rr)
595
c_rr = sqrt(gamma_rr * p_rr / rho_rr)
596
597
return max(abs(v_ll) + c_ll, abs(v_rr) + c_rr)
598
end
599
600
@inline function max_abs_speeds(u,
601
equations::CompressibleEulerMulticomponentEquations2D)
602
rho_v1, rho_v2, rho_e = u
603
604
rho = density(u, equations)
605
v1 = rho_v1 / rho
606
v2 = rho_v2 / rho
607
608
gamma = totalgamma(u, equations)
609
p = (gamma - 1) * (rho_e - 0.5f0 * rho * (v1^2 + v2^2))
610
c = sqrt(gamma * p / rho)
611
612
return (abs(v1) + c, abs(v2) + c)
613
end
614
615
@inline function rotate_to_x(u, normal_vector,
616
equations::CompressibleEulerMulticomponentEquations2D)
617
# cos and sin of the angle between the x-axis and the normalized normal_vector are
618
# the normalized vector's x and y coordinates respectively (see unit circle).
619
c = normal_vector[1]
620
s = normal_vector[2]
621
622
# Apply the 2D rotation matrix with normal and tangent directions of the form
623
# [ n_1 n_2 0 0;
624
# t_1 t_2 0 0;
625
# 0 0 1 0
626
# 0 0 0 1]
627
# where t_1 = -n_2 and t_2 = n_1
628
629
densities = @view u[4:end]
630
return SVector(c * u[1] + s * u[2],
631
-s * u[1] + c * u[2],
632
u[3],
633
densities...)
634
end
635
636
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction
637
# has been normalized prior to this back-rotation of the state vector
638
@inline function rotate_from_x(u, normal_vector,
639
equations::CompressibleEulerMulticomponentEquations2D)
640
# cos and sin of the angle between the x-axis and the normalized normal_vector are
641
# the normalized vector's x and y coordinates respectively (see unit circle).
642
c = normal_vector[1]
643
s = normal_vector[2]
644
645
# Apply the 2D back-rotation matrix with normal and tangent directions of the form
646
# [ n_1 t_1 0 0;
647
# n_2 t_2 0 0;
648
# 0 0 1 0;
649
# 0 0 0 1 ]
650
# where t_1 = -n_2 and t_2 = n_1
651
652
densities = @view u[4:end]
653
return SVector(c * u[1] - s * u[2],
654
s * u[1] + c * u[2],
655
u[3],
656
densities...)
657
end
658
659
# Convert conservative variables to primitive
660
@inline function cons2prim(u, equations::CompressibleEulerMulticomponentEquations2D)
661
rho_v1, rho_v2, rho_e = u
662
663
prim_rho = SVector{ncomponents(equations), real(equations)}(u[i + 3]
664
for i in eachcomponent(equations))
665
666
rho = density(u, equations)
667
v1 = rho_v1 / rho
668
v2 = rho_v2 / rho
669
gamma = totalgamma(u, equations)
670
p = (gamma - 1) * (rho_e - 0.5f0 * rho * (v1^2 + v2^2))
671
prim_other = SVector(v1, v2, p)
672
673
return vcat(prim_other, prim_rho)
674
end
675
676
# Convert conservative variables to entropy
677
@inline function cons2entropy(u, equations::CompressibleEulerMulticomponentEquations2D)
678
@unpack cv, gammas, gas_constants = equations
679
rho_v1, rho_v2, rho_e = u
680
681
rho = density(u, equations)
682
683
# Multicomponent stuff
684
RealT = eltype(u)
685
help1 = zero(RealT)
686
gas_constant = zero(RealT)
687
for i in eachcomponent(equations)
688
help1 += u[i + 3] * cv[i]
689
gas_constant += gas_constants[i] * (u[i + 3] / rho)
690
end
691
692
v1 = rho_v1 / rho
693
v2 = rho_v2 / rho
694
v_square = v1^2 + v2^2
695
gamma = totalgamma(u, equations)
696
697
p = (gamma - 1) * (rho_e - 0.5f0 * rho * v_square)
698
s = log(p) - gamma * log(rho) - log(gas_constant)
699
rho_p = rho / p
700
T = (rho_e - 0.5f0 * rho * v_square) / (help1)
701
702
entrop_rho = SVector{ncomponents(equations), real(equations)}((cv[i] *
703
(1 - log(T)) +
704
gas_constants[i] *
705
(1 + log(u[i + 3])) -
706
v_square / (2 * T))
707
for i in eachcomponent(equations))
708
709
w1 = gas_constant * v1 * rho_p
710
w2 = gas_constant * v2 * rho_p
711
w3 = gas_constant * (-rho_p)
712
713
entrop_other = SVector(w1, w2, w3)
714
715
return vcat(entrop_other, entrop_rho)
716
end
717
718
# Convert entropy variables to conservative variables
719
@inline function entropy2cons(w, equations::CompressibleEulerMulticomponentEquations2D)
720
@unpack gammas, gas_constants, cp, cv = equations
721
T = -1 / w[3]
722
v1 = w[1] * T
723
v2 = w[2] * T
724
v_squared = v1^2 + v2^2
725
cons_rho = SVector{ncomponents(equations), real(equations)}(exp((w[i + 3] -
726
cv[i] *
727
(1 - log(T)) +
728
v_squared /
729
(2 * T)) /
730
gas_constants[i] -
731
1)
732
for i in eachcomponent(equations))
733
734
RealT = eltype(w)
735
rho = zero(RealT)
736
help1 = zero(RealT)
737
help2 = zero(RealT)
738
p = zero(RealT)
739
for i in eachcomponent(equations)
740
rho += cons_rho[i]
741
help1 += cons_rho[i] * cv[i] * gammas[i]
742
help2 += cons_rho[i] * cv[i]
743
p += cons_rho[i] * gas_constants[i] * T
744
end
745
u1 = rho * v1
746
u2 = rho * v2
747
gamma = help1 / help2
748
u3 = p / (gamma - 1) + 0.5f0 * rho * v_squared
749
cons_other = SVector(u1, u2, u3)
750
return vcat(cons_other, cons_rho)
751
end
752
753
# Convert primitive to conservative variables
754
@inline function prim2cons(prim, equations::CompressibleEulerMulticomponentEquations2D)
755
@unpack cv, gammas = equations
756
v1, v2, p = prim
757
758
cons_rho = SVector{ncomponents(equations), real(equations)}(prim[i + 3]
759
for i in eachcomponent(equations))
760
rho = density(prim, equations)
761
gamma = totalgamma(prim, equations)
762
763
rho_v1 = rho * v1
764
rho_v2 = rho * v2
765
rho_e = p / (gamma - 1) + 0.5f0 * (rho_v1 * v1 + rho_v2 * v2)
766
767
cons_other = SVector(rho_v1, rho_v2, rho_e)
768
769
return vcat(cons_other, cons_rho)
770
end
771
772
@inline function total_entropy(u, equations::CompressibleEulerMulticomponentEquations2D)
773
@unpack cv, gammas, gas_constants = equations
774
rho = density(u, equations)
775
T = temperature(u, equations)
776
777
RealT = eltype(u)
778
total_entropy = zero(RealT)
779
for i in eachcomponent(equations)
780
total_entropy -= u[i + 3] * (cv[i] * log(T) - gas_constants[i] * log(u[i + 3]))
781
end
782
783
return total_entropy
784
end
785
786
@inline function temperature(u, equations::CompressibleEulerMulticomponentEquations2D)
787
@unpack cv, gammas, gas_constants = equations
788
789
rho_v1, rho_v2, rho_e = u
790
791
rho = density(u, equations)
792
RealT = eltype(u)
793
help1 = zero(RealT)
794
795
for i in eachcomponent(equations)
796
help1 += u[i + 3] * cv[i]
797
end
798
799
v1 = rho_v1 / rho
800
v2 = rho_v2 / rho
801
v_square = v1^2 + v2^2
802
T = (rho_e - 0.5f0 * rho * v_square) / help1
803
804
return T
805
end
806
807
"""
808
totalgamma(u, equations::CompressibleEulerMulticomponentEquations2D)
809
810
Function that calculates the total gamma out of all partial gammas using the
811
partial density fractions as well as the partial specific heats at constant volume.
812
"""
813
@inline function totalgamma(u, equations::CompressibleEulerMulticomponentEquations2D)
814
@unpack cv, gammas = equations
815
816
RealT = eltype(u)
817
help1 = zero(RealT)
818
help2 = zero(RealT)
819
820
for i in eachcomponent(equations)
821
help1 += u[i + 3] * cv[i] * gammas[i]
822
help2 += u[i + 3] * cv[i]
823
end
824
825
return help1 / help2
826
end
827
828
@inline function density_pressure(u,
829
equations::CompressibleEulerMulticomponentEquations2D)
830
rho_v1, rho_v2, rho_e = u
831
832
rho = density(u, equations)
833
gamma = totalgamma(u, equations)
834
rho_times_p = (gamma - 1) * (rho * rho_e - 0.5f0 * (rho_v1^2 + rho_v2^2))
835
836
return rho_times_p
837
end
838
839
@inline function density(u, equations::CompressibleEulerMulticomponentEquations2D)
840
RealT = eltype(u)
841
rho = zero(RealT)
842
843
for i in eachcomponent(equations)
844
rho += u[i + 3]
845
end
846
847
return rho
848
end
849
850
@inline function densities(u, v, equations::CompressibleEulerMulticomponentEquations2D)
851
return SVector{ncomponents(equations), real(equations)}(u[i + 3] * v
852
for i in eachcomponent(equations))
853
end
854
855
@inline function velocity(u, equations::CompressibleEulerMulticomponentEquations2D)
856
rho = density(u, equations)
857
v1 = u[1] / rho
858
v2 = u[2] / rho
859
return SVector(v1, v2)
860
end
861
862
@inline function velocity(u, orientation::Int,
863
equations::CompressibleEulerMulticomponentEquations2D)
864
rho = density(u, equations)
865
v = u[orientation] / rho
866
return v
867
end
868
end # @muladd
869
870