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_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
CompressibleEulerMulticomponentEquations1D(; 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 e_{\text{total}} \\ \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 e_{\text{total}} +p) v_1 \\ \rho_1 v_1 \\ \rho_2 v_1 \\ \vdots \\ \rho_{n} v_1
21
\end{pmatrix}
22
23
=
24
\begin{pmatrix}
25
0 \\ 0 \\ 0 \\ 0 \\ \vdots \\ 0
26
\end{pmatrix}
27
```
28
for calorically perfect gas in one space dimension.
29
Here, ``\rho_i`` is the density of component ``i``, ``\rho=\sum_{i=1}^n\rho_i`` the sum of the individual ``\rho_i``,
30
``v_1`` the velocity, ``e_{\text{total}}`` the specific total energy, and
31
```math
32
p = (\gamma - 1) \left( \rho e - \frac{1}{2} \rho v_1^2 \right)
33
```
34
the pressure,
35
```math
36
\gamma=\frac{\sum_{i=1}^n\rho_i C_{v,i}\gamma_i}{\sum_{i=1}^n\rho_i C_{v,i}}
37
```
38
total heat capacity ratio, ``\gamma_i`` heat capacity ratio of component ``i``,
39
```math
40
C_{v,i}=\frac{R}{\gamma_i-1}
41
```
42
specific heat capacity at constant volume of component ``i``.
43
44
In case of more than one component, the specific heat ratios `gammas` and the gas constants
45
`gas_constants` should be passed as tuples, e.g., `gammas=(1.4, 1.667)`.
46
47
The remaining variables like the specific heats at constant volume `cv` or the specific heats at
48
constant pressure `cp` are then calculated considering a calorically perfect gas.
49
"""
50
struct CompressibleEulerMulticomponentEquations1D{NVARS, NCOMP, RealT <: Real} <:
51
AbstractCompressibleEulerMulticomponentEquations{1, NVARS, NCOMP}
52
gammas::SVector{NCOMP, RealT}
53
gas_constants::SVector{NCOMP, RealT}
54
cv::SVector{NCOMP, RealT}
55
cp::SVector{NCOMP, RealT}
56
57
function CompressibleEulerMulticomponentEquations1D{NVARS, NCOMP, RealT}(gammas::SVector{NCOMP,
58
RealT},
59
gas_constants::SVector{NCOMP,
60
RealT}) where {
61
NVARS,
62
NCOMP,
63
RealT <:
64
Real
65
}
66
NCOMP >= 1 ||
67
throw(DimensionMismatch("`gammas` and `gas_constants` have to be filled with at least one value"))
68
69
cv = gas_constants ./ (gammas .- 1)
70
cp = gas_constants + gas_constants ./ (gammas .- 1)
71
72
return new(gammas, gas_constants, cv, cp)
73
end
74
end
75
76
function CompressibleEulerMulticomponentEquations1D(; gammas, gas_constants)
77
_gammas = promote(gammas...)
78
_gas_constants = promote(gas_constants...)
79
RealT = promote_type(eltype(_gammas), eltype(_gas_constants),
80
typeof(gas_constants[1] / (gammas[1] - 1)))
81
__gammas = SVector(map(RealT, _gammas))
82
__gas_constants = SVector(map(RealT, _gas_constants))
83
84
NVARS = length(_gammas) + 2
85
NCOMP = length(_gammas)
86
87
return CompressibleEulerMulticomponentEquations1D{NVARS, NCOMP, RealT}(__gammas,
88
__gas_constants)
89
end
90
91
@inline function Base.real(::CompressibleEulerMulticomponentEquations1D{NVARS, NCOMP,
92
RealT}) where {
93
NVARS,
94
NCOMP,
95
RealT
96
}
97
return RealT
98
end
99
100
function varnames(::typeof(cons2cons),
101
equations::CompressibleEulerMulticomponentEquations1D)
102
cons = ("rho_v1", "rho_e_total")
103
rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))
104
return (cons..., rhos...)
105
end
106
107
function varnames(::typeof(cons2prim),
108
equations::CompressibleEulerMulticomponentEquations1D)
109
prim = ("v1", "p")
110
rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))
111
return (prim..., rhos...)
112
end
113
114
# Set initial conditions at physical location `x` for time `t`
115
116
"""
117
initial_condition_convergence_test(x, t, equations::CompressibleEulerMulticomponentEquations1D)
118
119
A smooth initial condition used for convergence tests in combination with
120
[`source_terms_convergence_test`](@ref)
121
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
122
"""
123
function initial_condition_convergence_test(x, t,
124
equations::CompressibleEulerMulticomponentEquations1D)
125
RealT = eltype(x)
126
c = 2
127
A = convert(RealT, 0.1)
128
L = 2
129
f = 1.0f0 / L
130
omega = 2 * convert(RealT, pi) * f
131
ini = c + A * sin(omega * (x[1] - t))
132
133
v1 = 1
134
135
rho = ini
136
137
# Here we compute an arbitrary number of different rhos. (one rho is double the next rho while the sum of all rhos is 1)
138
prim_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) * (1 - 2) *
139
rho / (1 -
140
2^ncomponents(equations))
141
for i in eachcomponent(equations))
142
143
prim1 = rho * v1
144
prim2 = rho^2
145
146
prim_other = SVector(prim1, prim2)
147
148
return vcat(prim_other, prim_rho)
149
end
150
151
"""
152
source_terms_convergence_test(u, x, t, equations::CompressibleEulerMulticomponentEquations1D)
153
154
Source terms used for convergence tests in combination with
155
[`initial_condition_convergence_test`](@ref)
156
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
157
158
References for the method of manufactured solutions (MMS):
159
- Kambiz Salari and Patrick Knupp (2000)
160
Code Verification by the Method of Manufactured Solutions
161
[DOI: 10.2172/759450](https://doi.org/10.2172/759450)
162
- Patrick J. Roache (2002)
163
Code Verification by the Method of Manufactured Solutions
164
[DOI: 10.1115/1.1436090](https://doi.org/10.1115/1.1436090)
165
"""
166
@inline function source_terms_convergence_test(u, x, t,
167
equations::CompressibleEulerMulticomponentEquations1D)
168
# Same settings as in `initial_condition`
169
RealT = eltype(u)
170
c = 2
171
A = convert(RealT, 0.1)
172
L = 2
173
f = 1.0f0 / L
174
omega = 2 * convert(RealT, pi) * f
175
176
gamma = totalgamma(u, equations)
177
178
x1, = x
179
si, co = sincos((t - x1) * omega)
180
tmp = (-((4 * si * A - 4 * c) + 1) * (gamma - 1) * co * A * omega) / 2
181
182
# Here we compute an arbitrary number of different rhos. (one rho is double the next rho while the sum of all rhos is 1
183
du_rho = SVector{ncomponents(equations), real(equations)}(0
184
for i in eachcomponent(equations))
185
186
du1 = tmp
187
du2 = tmp
188
189
du_other = SVector(du1, du2)
190
191
return vcat(du_other, du_rho)
192
end
193
194
"""
195
initial_condition_weak_blast_wave(x, t, equations::CompressibleEulerMulticomponentEquations1D)
196
197
A for multicomponent adapted weak blast wave adapted to multicomponent and taken from
198
- Sebastian Hennemann, Gregor J. Gassner (2020)
199
A provably entropy stable subcell shock capturing approach for high order split form DG
200
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
201
"""
202
function initial_condition_weak_blast_wave(x, t,
203
equations::CompressibleEulerMulticomponentEquations1D)
204
# From Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
205
RealT = eltype(x)
206
inicenter = SVector(0)
207
x_norm = x[1] - inicenter[1]
208
r = abs(x_norm)
209
cos_phi = x_norm > 0 ? 1 : -1
210
211
prim_rho = SVector{ncomponents(equations), real(equations)}(r > 0.5f0 ?
212
2^(i - 1) * (1 - 2) /
213
(RealT(1) -
214
2^ncomponents(equations)) :
215
2^(i - 1) * (1 - 2) *
216
RealT(1.1691) /
217
(1 -
218
2^ncomponents(equations))
219
for i in eachcomponent(equations))
220
221
v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos_phi
222
p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)
223
224
prim_other = SVector(v1, p)
225
226
return prim2cons(vcat(prim_other, prim_rho), equations)
227
end
228
229
# Calculate 1D flux for a single point
230
@inline function flux(u, orientation::Integer,
231
equations::CompressibleEulerMulticomponentEquations1D)
232
rho_v1, rho_e_total = u
233
234
rho = density(u, equations)
235
236
v1 = rho_v1 / rho
237
gamma = totalgamma(u, equations)
238
p = (gamma - 1) * (rho_e_total - 0.5f0 * rho * v1^2)
239
240
f_rho = densities(u, v1, equations)
241
f1 = rho_v1 * v1 + p
242
f2 = (rho_e_total + p) * v1
243
244
f_other = SVector(f1, f2)
245
246
return vcat(f_other, f_rho)
247
end
248
249
"""
250
flux_chandrashekar(u_ll, u_rr, orientation, equations::CompressibleEulerMulticomponentEquations1D)
251
252
Entropy conserving two-point flux by
253
- Ayoub Gouasmi, Karthik Duraisamy (2020)
254
"Formulation of Entropy-Stable schemes for the multicomponent compressible Euler equations"
255
[arXiv:1904.00972v3](https://arxiv.org/abs/1904.00972) [math.NA] 4 Feb 2020
256
"""
257
@inline function flux_chandrashekar(u_ll, u_rr, orientation::Integer,
258
equations::CompressibleEulerMulticomponentEquations1D)
259
# Unpack left and right state
260
@unpack gammas, gas_constants, cv = equations
261
rho_v1_ll, rho_e_total_ll = u_ll
262
rho_v1_rr, rho_e_total_rr = u_rr
263
rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 2],
264
u_rr[i + 2])
265
for i in eachcomponent(equations))
266
rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 2] +
267
u_rr[i + 2])
268
for i in eachcomponent(equations))
269
270
# Iterating over all partial densities
271
rho_ll = density(u_ll, equations)
272
rho_rr = density(u_rr, equations)
273
274
gamma_ll = totalgamma(u_ll, equations)
275
gamma_rr = totalgamma(u_rr, equations)
276
277
# extract velocities
278
v1_ll = rho_v1_ll / rho_ll
279
v1_rr = rho_v1_rr / rho_rr
280
v1_avg = 0.5f0 * (v1_ll + v1_rr)
281
v1_square = 0.5f0 * (v1_ll^2 + v1_rr^2)
282
v_sum = v1_avg
283
284
RealT = eltype(u_ll)
285
enth = zero(RealT)
286
help1_ll = zero(RealT)
287
help1_rr = zero(RealT)
288
289
for i in eachcomponent(equations)
290
enth += rhok_avg[i] * gas_constants[i]
291
help1_ll += u_ll[i + 2] * cv[i]
292
help1_rr += u_rr[i + 2] * cv[i]
293
end
294
295
T_ll = (rho_e_total_ll - 0.5f0 * rho_ll * (v1_ll^2)) / help1_ll
296
T_rr = (rho_e_total_rr - 0.5f0 * rho_rr * (v1_rr^2)) / help1_rr
297
T = 0.5f0 * (1 / T_ll + 1 / T_rr)
298
T_log = ln_mean(1 / T_ll, 1 / T_rr)
299
300
# Calculate fluxes depending on orientation
301
help1 = zero(RealT)
302
help2 = zero(RealT)
303
304
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg
305
for i in eachcomponent(equations))
306
for i in eachcomponent(equations)
307
help1 += f_rho[i] * cv[i]
308
help2 += f_rho[i]
309
end
310
f1 = (help2) * v1_avg + enth / T
311
f2 = (help1) / T_log - 0.5f0 * (v1_square) * (help2) + v1_avg * f1
312
313
f_other = SVector(f1, f2)
314
315
return vcat(f_other, f_rho)
316
end
317
318
"""
319
flux_ranocha(u_ll, u_rr, orientation_or_normal_direction,
320
equations::CompressibleEulerMulticomponentEquations1D)
321
322
Adaption of the entropy conserving and kinetic energy preserving two-point flux by
323
- Hendrik Ranocha (2018)
324
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
325
for Hyperbolic Balance Laws
326
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
327
See also
328
- Hendrik Ranocha (2020)
329
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
330
the Euler Equations Using Summation-by-Parts Operators
331
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
332
"""
333
@inline function flux_ranocha(u_ll, u_rr, orientation::Integer,
334
equations::CompressibleEulerMulticomponentEquations1D)
335
# Unpack left and right state
336
@unpack gammas, gas_constants, cv = equations
337
rho_v1_ll, rho_e_total_ll = u_ll
338
rho_v1_rr, rho_e_total_rr = u_rr
339
rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 2],
340
u_rr[i + 2])
341
for i in eachcomponent(equations))
342
rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 2] +
343
u_rr[i + 2])
344
for i in eachcomponent(equations))
345
346
# Iterating over all partial densities
347
rho_ll = density(u_ll, equations)
348
rho_rr = density(u_rr, equations)
349
350
# Calculating gamma
351
gamma = totalgamma(0.5f0 * (u_ll + u_rr), equations)
352
inv_gamma_minus_one = 1 / (gamma - 1)
353
354
# extract velocities
355
v1_ll = rho_v1_ll / rho_ll
356
v1_rr = rho_v1_rr / rho_rr
357
v1_avg = 0.5f0 * (v1_ll + v1_rr)
358
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr)
359
360
# density flux
361
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg
362
for i in eachcomponent(equations))
363
364
# helpful variables
365
RealT = eltype(u_ll)
366
f_rho_sum = zero(RealT)
367
help1_ll = zero(RealT)
368
help1_rr = zero(RealT)
369
enth_ll = zero(RealT)
370
enth_rr = zero(RealT)
371
for i in eachcomponent(equations)
372
enth_ll += u_ll[i + 2] * gas_constants[i]
373
enth_rr += u_rr[i + 2] * gas_constants[i]
374
f_rho_sum += f_rho[i]
375
help1_ll += u_ll[i + 2] * cv[i]
376
help1_rr += u_rr[i + 2] * cv[i]
377
end
378
379
# temperature and pressure
380
T_ll = (rho_e_total_ll - 0.5f0 * rho_ll * (v1_ll^2)) / help1_ll
381
T_rr = (rho_e_total_rr - 0.5f0 * rho_rr * (v1_rr^2)) / help1_rr
382
p_ll = T_ll * enth_ll
383
p_rr = T_rr * enth_rr
384
p_avg = 0.5f0 * (p_ll + p_rr)
385
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
386
387
# momentum and energy flux
388
f1 = f_rho_sum * v1_avg + p_avg
389
f2 = f_rho_sum * (velocity_square_avg + inv_rho_p_mean * inv_gamma_minus_one) +
390
0.5f0 * (p_ll * v1_rr + p_rr * v1_ll)
391
f_other = SVector(f1, f2)
392
393
return vcat(f_other, f_rho)
394
end
395
396
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation
397
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
398
equations::CompressibleEulerMulticomponentEquations1D)
399
rho_v1_ll, rho_e_total_ll = u_ll
400
rho_v1_rr, rho_e_total_rr = u_rr
401
402
# Calculate primitive variables and speed of sound
403
rho_ll = density(u_ll, equations)
404
rho_rr = density(u_rr, equations)
405
gamma_ll = totalgamma(u_ll, equations)
406
gamma_rr = totalgamma(u_rr, equations)
407
408
v_ll = rho_v1_ll / rho_ll
409
v_rr = rho_v1_rr / rho_rr
410
411
p_ll = (gamma_ll - 1) * (rho_e_total_ll - 0.5f0 * rho_ll * v_ll^2)
412
p_rr = (gamma_rr - 1) * (rho_e_total_rr - 0.5f0 * rho_rr * v_rr^2)
413
c_ll = sqrt(gamma_ll * p_ll / rho_ll)
414
c_rr = sqrt(gamma_rr * p_rr / rho_rr)
415
416
return max(abs(v_ll), abs(v_rr)) + max(c_ll, c_rr)
417
end
418
419
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
420
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
421
equations::CompressibleEulerMulticomponentEquations1D)
422
rho_v1_ll, rho_e_total_ll = u_ll
423
rho_v1_rr, rho_e_total_rr = u_rr
424
425
# Calculate primitive variables and speed of sound
426
rho_ll = density(u_ll, equations)
427
rho_rr = density(u_rr, equations)
428
gamma_ll = totalgamma(u_ll, equations)
429
gamma_rr = totalgamma(u_rr, equations)
430
431
v_ll = rho_v1_ll / rho_ll
432
v_rr = rho_v1_rr / rho_rr
433
434
p_ll = (gamma_ll - 1) * (rho_e_total_ll - 0.5f0 * rho_ll * v_ll^2)
435
p_rr = (gamma_rr - 1) * (rho_e_total_rr - 0.5f0 * rho_rr * v_rr^2)
436
c_ll = sqrt(gamma_ll * p_ll / rho_ll)
437
c_rr = sqrt(gamma_rr * p_rr / rho_rr)
438
439
return max(abs(v_ll) + c_ll, abs(v_rr) + c_rr)
440
end
441
442
@inline function max_abs_speeds(u,
443
equations::CompressibleEulerMulticomponentEquations1D)
444
rho_v1, rho_e_total = u
445
446
rho = density(u, equations)
447
v1 = rho_v1 / rho
448
449
gamma = totalgamma(u, equations)
450
p = (gamma - 1) * (rho_e_total - 0.5f0 * rho * (v1^2))
451
c = sqrt(gamma * p / rho)
452
453
return (abs(v1) + c,)
454
end
455
456
# Convert conservative variables to primitive
457
@inline function cons2prim(u, equations::CompressibleEulerMulticomponentEquations1D)
458
rho_v1, rho_e_total = u
459
460
prim_rho = SVector{ncomponents(equations), real(equations)}(u[i + 2]
461
for i in eachcomponent(equations))
462
463
rho = density(u, equations)
464
v1 = rho_v1 / rho
465
gamma = totalgamma(u, equations)
466
467
p = (gamma - 1) * (rho_e_total - 0.5f0 * rho * (v1^2))
468
prim_other = SVector(v1, p)
469
470
return vcat(prim_other, prim_rho)
471
end
472
473
# Convert primitive to conservative variables
474
@inline function prim2cons(prim, equations::CompressibleEulerMulticomponentEquations1D)
475
@unpack cv, gammas = equations
476
v1, p = prim
477
478
RealT = eltype(prim)
479
480
cons_rho = SVector{ncomponents(equations), RealT}(prim[i + 2]
481
for i in eachcomponent(equations))
482
rho = density(prim, equations)
483
gamma = totalgamma(prim, equations)
484
485
rho_v1 = rho * v1
486
487
rho_e_total = p / (gamma - 1) + 0.5f0 * (rho_v1 * v1)
488
489
cons_other = SVector{2, RealT}(rho_v1, rho_e_total)
490
491
return vcat(cons_other, cons_rho)
492
end
493
494
# Convert conservative variables to entropy
495
@inline function cons2entropy(u, equations::CompressibleEulerMulticomponentEquations1D)
496
@unpack cv, gammas, gas_constants = equations
497
498
rho_v1, rho_e_total = u
499
500
rho = density(u, equations)
501
502
RealT = eltype(u)
503
help1 = zero(RealT)
504
gas_constant = zero(RealT)
505
for i in eachcomponent(equations)
506
help1 += u[i + 2] * cv[i]
507
gas_constant += gas_constants[i] * (u[i + 2] / rho)
508
end
509
510
v1 = rho_v1 / rho
511
v_square = v1^2
512
gamma = totalgamma(u, equations)
513
514
p = (gamma - 1) * (rho_e_total - 0.5f0 * rho * v_square)
515
s = log(p) - gamma * log(rho) - log(gas_constant)
516
rho_p = rho / p
517
T = (rho_e_total - 0.5f0 * rho * v_square) / (help1)
518
519
entrop_rho = SVector{ncomponents(equations), real(equations)}((cv[i] *
520
(1 - log(T)) +
521
gas_constants[i] *
522
(1 + log(u[i + 2])) -
523
v1^2 / (2 * T))
524
for i in eachcomponent(equations))
525
526
w1 = gas_constant * v1 * rho_p
527
w2 = gas_constant * (-rho_p)
528
529
entrop_other = SVector(w1, w2)
530
531
return vcat(entrop_other, entrop_rho)
532
end
533
534
# Convert entropy variables to conservative variables
535
@inline function entropy2cons(w, equations::CompressibleEulerMulticomponentEquations1D)
536
@unpack gammas, gas_constants, cv, cp = equations
537
T = -1 / w[2]
538
v1 = w[1] * T
539
cons_rho = SVector{ncomponents(equations), real(equations)}(exp(1 /
540
gas_constants[i] *
541
(-cv[i] *
542
log(-w[2]) -
543
cp[i] + w[i + 2] -
544
0.5f0 * w[1]^2 /
545
w[2]))
546
for i in eachcomponent(equations))
547
548
RealT = eltype(w)
549
rho = zero(RealT)
550
help1 = zero(RealT)
551
help2 = zero(RealT)
552
p = zero(RealT)
553
for i in eachcomponent(equations)
554
rho += cons_rho[i]
555
help1 += cons_rho[i] * cv[i] * gammas[i]
556
help2 += cons_rho[i] * cv[i]
557
p += cons_rho[i] * gas_constants[i] * T
558
end
559
u1 = rho * v1
560
gamma = help1 / help2
561
u2 = p / (gamma - 1) + 0.5f0 * rho * v1^2
562
cons_other = SVector(u1, u2)
563
return vcat(cons_other, cons_rho)
564
end
565
566
@inline function total_entropy(u, equations::CompressibleEulerMulticomponentEquations1D)
567
@unpack cv, gammas, gas_constants = equations
568
rho_v1, rho_e_total = u
569
rho = density(u, equations)
570
T = temperature(u, equations)
571
572
RealT = eltype(u)
573
total_entropy = zero(RealT)
574
for i in eachcomponent(equations)
575
total_entropy -= u[i + 2] * (cv[i] * log(T) - gas_constants[i] * log(u[i + 2]))
576
end
577
578
return total_entropy
579
end
580
581
@inline function temperature(u, equations::CompressibleEulerMulticomponentEquations1D)
582
@unpack cv, gammas, gas_constants = equations
583
584
rho_v1, rho_e_total = u
585
586
rho = density(u, equations)
587
588
RealT = eltype(u)
589
help1 = zero(RealT)
590
591
for i in eachcomponent(equations)
592
help1 += u[i + 2] * cv[i]
593
end
594
595
v1 = rho_v1 / rho
596
v_square = v1^2
597
T = (rho_e_total - 0.5f0 * rho * v_square) / help1
598
599
return T
600
end
601
602
"""
603
totalgamma(u, equations::CompressibleEulerMulticomponentEquations1D)
604
605
Function that calculates the total gamma out of all partial gammas using the
606
partial density fractions as well as the partial specific heats at constant volume.
607
"""
608
@inline function totalgamma(u, equations::CompressibleEulerMulticomponentEquations1D)
609
@unpack cv, gammas = equations
610
611
RealT = eltype(u)
612
help1 = zero(RealT)
613
help2 = zero(RealT)
614
615
for i in eachcomponent(equations)
616
help1 += u[i + 2] * cv[i] * gammas[i]
617
help2 += u[i + 2] * cv[i]
618
end
619
620
return help1 / help2
621
end
622
623
@doc raw"""
624
pressure(u, equations::AbstractCompressibleEulerMulticomponentEquations)
625
626
Computes the pressure for an ideal equation of state with
627
isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`.
628
```math
629
\begin{aligned}
630
p &= (\gamma - 1) \left( E_\mathrm{tot} - E_\mathrm{kin} \right) \\
631
&= (\gamma - 1) \left( \rho e - \frac{1}{2}\rho \Vert v \Vert_2^2 \right)
632
\end{aligned}
633
```
634
"""
635
@inline function pressure(u, equations::CompressibleEulerMulticomponentEquations1D)
636
rho_v1, rho_e_total = u
637
638
rho = density(u, equations)
639
gamma = totalgamma(u, equations)
640
641
p = (gamma - 1) * (rho_e_total - 0.5f0 * (rho_v1^2) / rho)
642
643
return p
644
end
645
646
@inline function density_pressure(u,
647
equations::CompressibleEulerMulticomponentEquations1D)
648
rho_v1, rho_e_total = u
649
650
rho = density(u, equations)
651
gamma = totalgamma(u, equations)
652
rho_times_p = (gamma - 1) * (rho * rho_e_total - 0.5f0 * rho_v1^2)
653
654
return rho_times_p
655
end
656
657
@doc raw"""
658
density(u, equations::AbstractCompressibleEulerMulticomponentEquations)
659
660
Computes the total density ``\rho = \sum_{i=1}^n \rho_i`` from the conserved variables `u`.
661
"""
662
@inline function density(u, equations::CompressibleEulerMulticomponentEquations1D)
663
RealT = eltype(u)
664
rho = zero(RealT)
665
666
for i in eachcomponent(equations)
667
rho += u[i + 2]
668
end
669
670
return rho
671
end
672
673
@inline function densities(u, v, equations::CompressibleEulerMulticomponentEquations1D)
674
return SVector{ncomponents(equations), real(equations)}(u[i + 2] * v
675
for i in eachcomponent(equations))
676
end
677
678
@inline function velocity(u, equations::CompressibleEulerMulticomponentEquations1D)
679
rho = density(u, equations)
680
v1 = u[1] / rho
681
return v1
682
end
683
end # @muladd
684
685