Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/equations/ideal_glm_mhd_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
IdealGlmMhdMulticomponentEquations1D(; gammas, gas_constants)
10
11
The ideal compressible multicomponent GLM-MHD equations
12
```math
13
\frac{\partial}{\partial t}
14
\begin{pmatrix}
15
\rho v_1 \\ \rho v_2 \\ \rho v_3 \\ \rho e_{\text{total}} \\ B_1 \\ B_2 \\ B_3 \\ \rho_1 \\ \vdots \\ \rho_{n}
16
\end{pmatrix}
17
+
18
\frac{\partial}{\partial x}
19
\begin{pmatrix}
20
\rho v_1^2 + p + \frac{1}{2} \Vert \mathbf{B} \Vert_2 ^2 - B_1^2 \\ \rho v_1 v_2 - B_1 B_2 \\ \rho v_1 v_3 - B_1 B_3
21
\\ (\frac{1}{2} \rho \Vert \mathbf{v} \Vert_2 ^2 + \frac{\gamma p}{\gamma - 1} + \Vert \mathbf{B} \Vert_2 ^2) v_1 - B_1 (\mathbf{v} \cdot \mathbf{B})
22
\\ 0 \\ v_1 B_2 - v_2 B_1 \\ v_1 B_3 - v_3 B_1 \\ \rho_1 v_1 \\ \vdots \\ \rho_{n} v_1
23
\end{pmatrix}
24
=
25
\begin{pmatrix}
26
0 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0 \\ \vdots \\ 0
27
\end{pmatrix}
28
```
29
for calorically perfect gases in one space dimension.
30
Here, ``\rho_i`` is the density of component ``i``, ``\rho=\sum_{i=1}^n\rho_i`` the sum of the individual ``\rho_i``,
31
``\mathbf{v}`` the velocity, ``\mathbf{B}`` the magnetic field, ``e_{\text{total}}`` the specific total energy, and
32
```math
33
p = (\gamma - 1) \left( \rho e_{\text{total}} - \frac{1}{2} \rho \Vert \mathbf{v} \Vert_2 ^2 - \frac{1}{2} \Vert \mathbf{B} \Vert_2 ^2 \right)
34
```
35
the pressure,
36
```math
37
\gamma=\frac{\sum_{i=1}^n\rho_i C_{v,i}\gamma_i}{\sum_{i=1}^n\rho_i C_{v,i}}
38
```
39
total heat capacity ratio, ``\gamma_i`` heat capacity ratio of component ``i``,
40
```math
41
C_{v,i}=\frac{R_i}{\gamma_i-1}
42
```
43
specific heat capacity at constant volume of component ``i``.
44
45
In case of more than one component, the specific heat ratios `gammas` and the gas constants
46
`gas_constants` should be passed as tuples, e.g., `gammas = (1.4, 1.667)`.
47
48
The remaining variables like the specific heats at constant volume `cv` or the specific heats at
49
constant pressure `cp` are then calculated considering a calorically perfect gas.
50
"""
51
struct IdealGlmMhdMulticomponentEquations1D{NVARS, NCOMP, RealT <: Real} <:
52
AbstractIdealGlmMhdMulticomponentEquations{1, NVARS, NCOMP}
53
gammas::SVector{NCOMP, RealT}
54
gas_constants::SVector{NCOMP, RealT}
55
cv::SVector{NCOMP, RealT}
56
cp::SVector{NCOMP, RealT}
57
58
function IdealGlmMhdMulticomponentEquations1D{NVARS, NCOMP, RealT}(gammas::SVector{NCOMP,
59
RealT},
60
gas_constants::SVector{NCOMP,
61
RealT}) where {
62
NVARS,
63
NCOMP,
64
RealT <:
65
Real
66
}
67
NCOMP >= 1 ||
68
throw(DimensionMismatch("`gammas` and `gas_constants` have to be filled with at least one value"))
69
70
cv = gas_constants ./ (gammas .- 1)
71
cp = gas_constants + gas_constants ./ (gammas .- 1)
72
73
return new(gammas, gas_constants, cv, cp)
74
end
75
end
76
77
function IdealGlmMhdMulticomponentEquations1D(; gammas, gas_constants)
78
_gammas = promote(gammas...)
79
_gas_constants = promote(gas_constants...)
80
RealT = promote_type(eltype(_gammas), eltype(_gas_constants))
81
82
__gammas = SVector(map(RealT, _gammas))
83
__gas_constants = SVector(map(RealT, _gas_constants))
84
85
NVARS = length(_gammas) + 7
86
NCOMP = length(_gammas)
87
88
return IdealGlmMhdMulticomponentEquations1D{NVARS, NCOMP, RealT}(__gammas,
89
__gas_constants)
90
end
91
92
# Outer constructor for `@reset` works correctly
93
function IdealGlmMhdMulticomponentEquations1D(gammas, gas_constants, cv, cp, c_h)
94
return IdealGlmMhdMulticomponentEquations1D(gammas = gammas,
95
gas_constants = gas_constants)
96
end
97
98
@inline function Base.real(::IdealGlmMhdMulticomponentEquations1D{NVARS, NCOMP, RealT}) where {
99
NVARS,
100
NCOMP,
101
RealT
102
}
103
return RealT
104
end
105
106
have_nonconservative_terms(::IdealGlmMhdMulticomponentEquations1D) = False()
107
108
function varnames(::typeof(cons2cons), equations::IdealGlmMhdMulticomponentEquations1D)
109
cons = ("rho_v1", "rho_v2", "rho_v3", "rho_e_total", "B1", "B2", "B3")
110
rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))
111
return (cons..., rhos...)
112
end
113
114
function varnames(::typeof(cons2prim), equations::IdealGlmMhdMulticomponentEquations1D)
115
prim = ("v1", "v2", "v3", "p", "B1", "B2", "B3")
116
rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))
117
return (prim..., rhos...)
118
end
119
120
"""
121
initial_condition_convergence_test(x, t, equations::IdealGlmMhdMulticomponentEquations1D)
122
123
An Alfvén wave as smooth initial condition used for convergence tests.
124
"""
125
function initial_condition_convergence_test(x, t,
126
equations::IdealGlmMhdMulticomponentEquations1D)
127
# smooth Alfvén wave test from Derigs et al. FLASH (2016)
128
# domain must be set to [0, 1], γ = 5/3
129
130
RealT = eltype(x)
131
rho = one(RealT)
132
prim_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) * (1 - 2) *
133
rho / (1 -
134
2^ncomponents(equations))
135
for i in eachcomponent(equations))
136
v1 = 0
137
si, co = sincospi(2 * x[1])
138
v2 = convert(RealT, 0.1) * si
139
v3 = convert(RealT, 0.1) * co
140
p = convert(RealT, 0.1)
141
B1 = 1
142
B2 = v2
143
B3 = v3
144
prim_other = SVector(v1, v2, v3, p, B1, B2, B3)
145
146
return prim2cons(vcat(prim_other, prim_rho), equations)
147
end
148
149
"""
150
initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdMulticomponentEquations1D)
151
152
A weak blast wave adapted from
153
- Sebastian Hennemann, Gregor J. Gassner (2020)
154
A provably entropy stable subcell shock capturing approach for high order split form DG
155
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
156
"""
157
function initial_condition_weak_blast_wave(x, t,
158
equations::IdealGlmMhdMulticomponentEquations1D)
159
# Adapted MHD version of the weak blast wave from Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
160
# Same discontinuity in the velocities but with magnetic fields
161
# Set up polar coordinates
162
RealT = eltype(x)
163
inicenter = (0)
164
x_norm = x[1] - inicenter[1]
165
r = sqrt(x_norm^2)
166
phi = atan(x_norm)
167
168
# Calculate primitive variables
169
# We initialize each species with a fraction of the total density `rho`, such
170
# that the sum of the densities is `rho := density(prim, equations)`. The density of
171
# a species is double the density of the next species.
172
if r > 0.5f0
173
rho = one(RealT)
174
prim_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) *
175
(1 - 2) * rho /
176
(1 -
177
2^ncomponents(equations))
178
for i in eachcomponent(equations))
179
else
180
rho = convert(RealT, 1.1691)
181
prim_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) *
182
(1 - 2) * rho /
183
(1 -
184
2^ncomponents(equations))
185
for i in eachcomponent(equations))
186
end
187
v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos(phi)
188
p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)
189
190
prim_other = SVector(v1, 0, 0, p, 1, 1, 1)
191
192
return prim2cons(vcat(prim_other, prim_rho), equations)
193
end
194
195
# Calculate 1D flux for a single point
196
@inline function flux(u, orientation::Integer,
197
equations::IdealGlmMhdMulticomponentEquations1D)
198
rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u
199
200
rho = density(u, equations)
201
202
v1 = rho_v1 / rho
203
v2 = rho_v2 / rho
204
v3 = rho_v3 / rho
205
kin_en = 0.5f0 * rho * (v1^2 + v2^2 + v3^2)
206
mag_en = 0.5f0 * (B1^2 + B2^2 + B3^2)
207
gamma = totalgamma(u, equations)
208
p = (gamma - 1) * (rho_e_total - kin_en - mag_en)
209
210
f_rho = densities(u, v1, equations)
211
f1 = rho_v1 * v1 + p + mag_en - B1^2
212
f2 = rho_v1 * v2 - B1 * B2
213
f3 = rho_v1 * v3 - B1 * B3
214
f4 = (kin_en + gamma * p / (gamma - 1) + 2 * mag_en) * v1 -
215
B1 * (v1 * B1 + v2 * B2 + v3 * B3)
216
f5 = 0
217
f6 = v1 * B2 - v2 * B1
218
f7 = v1 * B3 - v3 * B1
219
220
f_other = SVector(f1, f2, f3, f4, f5, f6, f7)
221
222
return vcat(f_other, f_rho)
223
end
224
225
"""
226
flux_derigs_etal(u_ll, u_rr, orientation, equations::IdealGlmMhdEquations1D)
227
228
Entropy conserving two-point flux adapted by
229
- Derigs et al. (2018)
230
Ideal GLM-MHD: About the entropy consistent nine-wave magnetic field
231
divergence diminishing ideal magnetohydrodynamics equations for multicomponent
232
[DOI: 10.1016/j.jcp.2018.03.002](https://doi.org/10.1016/j.jcp.2018.03.002)
233
"""
234
function flux_derigs_etal(u_ll, u_rr, orientation::Integer,
235
equations::IdealGlmMhdMulticomponentEquations1D)
236
# Unpack left and right states to get velocities, pressure, and inverse temperature (called beta)
237
rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll = u_ll
238
rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr = u_rr
239
@unpack gammas, gas_constants, cv = equations
240
241
rho_ll = density(u_ll, equations)
242
rho_rr = density(u_rr, equations)
243
244
gamma_ll = totalgamma(u_ll, equations)
245
gamma_rr = totalgamma(u_rr, equations)
246
247
rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 7],
248
u_rr[i + 7])
249
for i in eachcomponent(equations))
250
rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 7] +
251
u_rr[i + 7])
252
for i in eachcomponent(equations))
253
254
v1_ll = rho_v1_ll / rho_ll
255
v2_ll = rho_v2_ll / rho_ll
256
v3_ll = rho_v3_ll / rho_ll
257
v1_rr = rho_v1_rr / rho_rr
258
v2_rr = rho_v2_rr / rho_rr
259
v3_rr = rho_v3_rr / rho_rr
260
vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2
261
vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2
262
mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^2
263
mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^2
264
# for convenience store v⋅B
265
vel_dot_mag_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
266
vel_dot_mag_rr = v1_rr * B1_rr + v2_rr * B2_rr + v3_rr * B3_rr
267
268
# Compute the necessary mean values needed for either direction
269
v1_avg = 0.5f0 * (v1_ll + v1_rr)
270
v2_avg = 0.5f0 * (v2_ll + v2_rr)
271
v3_avg = 0.5f0 * (v3_ll + v3_rr)
272
v_sum = v1_avg + v2_avg + v3_avg
273
B1_avg = 0.5f0 * (B1_ll + B1_rr)
274
B2_avg = 0.5f0 * (B2_ll + B2_rr)
275
B3_avg = 0.5f0 * (B3_ll + B3_rr)
276
vel_norm_avg = 0.5f0 * (vel_norm_ll + vel_norm_rr)
277
mag_norm_avg = 0.5f0 * (mag_norm_ll + mag_norm_rr)
278
vel_dot_mag_avg = 0.5f0 * (vel_dot_mag_ll + vel_dot_mag_rr)
279
280
RealT = eltype(u_ll)
281
enth = zero(RealT)
282
help1_ll = zero(RealT)
283
help1_rr = zero(RealT)
284
285
for i in eachcomponent(equations)
286
enth += rhok_avg[i] * gas_constants[i]
287
help1_ll += u_ll[i + 7] * cv[i]
288
help1_rr += u_rr[i + 7] * cv[i]
289
end
290
291
T_ll = (rho_e_total_ll - 0.5f0 * rho_ll * (vel_norm_ll) - 0.5f0 * mag_norm_ll) /
292
help1_ll
293
T_rr = (rho_e_total_rr - 0.5f0 * rho_rr * (vel_norm_rr) - 0.5f0 * mag_norm_rr) /
294
help1_rr
295
T = 0.5f0 * (1 / T_ll + 1 / T_rr)
296
T_log = ln_mean(1 / T_ll, 1 / T_rr)
297
298
# Calculate fluxes depending on orientation with specific direction averages
299
help1 = zero(RealT)
300
help2 = zero(RealT)
301
302
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg
303
for i in eachcomponent(equations))
304
for i in eachcomponent(equations)
305
help1 += f_rho[i] * cv[i]
306
help2 += f_rho[i]
307
end
308
f1 = help2 * v1_avg + enth / T + 0.5f0 * mag_norm_avg - B1_avg * B1_avg
309
f2 = help2 * v2_avg - B1_avg * B2_avg
310
f3 = help2 * v3_avg - B1_avg * B3_avg
311
f5 = 0
312
f6 = v1_avg * B2_avg - v2_avg * B1_avg
313
f7 = v1_avg * B3_avg - v3_avg * B1_avg
314
315
# total energy flux is complicated and involves the previous eight components
316
v1_mag_avg = 0.5f0 * (v1_ll * mag_norm_ll + v1_rr * mag_norm_rr)
317
318
f4 = (help1 / T_log) - 0.5f0 * (vel_norm_avg) * (help2) + f1 * v1_avg +
319
f2 * v2_avg +
320
f3 * v3_avg +
321
f5 * B1_avg + f6 * B2_avg + f7 * B3_avg - 0.5f0 * v1_mag_avg +
322
B1_avg * vel_dot_mag_avg
323
324
f_other = SVector(f1, f2, f3, f4, f5, f6, f7)
325
326
return vcat(f_other, f_rho)
327
end
328
329
"""
330
flux_hindenlang_gassner(u_ll, u_rr, orientation_or_normal_direction,
331
equations::IdealGlmMhdMulticomponentEquations1D)
332
333
Adaption of the entropy conserving and kinetic energy preserving two-point flux of
334
Hindenlang (2019), extending [`flux_ranocha`](@ref) to the MHD equations.
335
## References
336
- Florian Hindenlang, Gregor Gassner (2019)
337
A new entropy conservative two-point flux for ideal MHD equations derived from
338
first principles.
339
Presented at HONOM 2019: European workshop on high order numerical methods
340
for evolutionary PDEs, theory and applications
341
- Hendrik Ranocha (2018)
342
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
343
for Hyperbolic Balance Laws
344
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
345
- Hendrik Ranocha (2020)
346
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
347
the Euler Equations Using Summation-by-Parts Operators
348
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
349
"""
350
@inline function flux_hindenlang_gassner(u_ll, u_rr, orientation::Integer,
351
equations::IdealGlmMhdMulticomponentEquations1D)
352
# Unpack left and right states
353
v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll = cons2prim(u_ll, equations)
354
v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr = cons2prim(u_rr, equations)
355
356
rho_ll = density(u_ll, equations)
357
rho_rr = density(u_rr, equations)
358
359
# Compute the necessary mean values needed for either direction
360
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
361
# in exact arithmetic since
362
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
363
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
364
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
365
v1_avg = 0.5f0 * (v1_ll + v1_rr)
366
v2_avg = 0.5f0 * (v2_ll + v2_rr)
367
v3_avg = 0.5f0 * (v3_ll + v3_rr)
368
p_avg = 0.5f0 * (p_ll + p_rr)
369
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
370
magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)
371
372
inv_gamma_minus_one = 1 / (totalgamma(0.5f0 * (u_ll + u_rr), equations) - 1)
373
374
rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 7],
375
u_rr[i + 7])
376
for i in eachcomponent(equations))
377
rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 7] +
378
u_rr[i + 7])
379
for i in eachcomponent(equations))
380
381
RealT = eltype(u_ll)
382
f1 = zero(RealT)
383
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg
384
for i in eachcomponent(equations))
385
for i in eachcomponent(equations)
386
f1 += f_rho[i]
387
end
388
389
# Calculate fluxes depending on orientation with specific direction averages
390
f2 = f1 * v1_avg + p_avg + magnetic_square_avg -
391
0.5f0 * (B1_ll * B1_rr + B1_rr * B1_ll)
392
f3 = f1 * v2_avg - 0.5f0 * (B1_ll * B2_rr + B1_rr * B2_ll)
393
f4 = f1 * v3_avg - 0.5f0 * (B1_ll * B3_rr + B1_rr * B3_ll)
394
# f5 below
395
f6 = 0
396
f7 = 0.5f0 * (v1_ll * B2_ll - v2_ll * B1_ll + v1_rr * B2_rr - v2_rr * B1_rr)
397
f8 = 0.5f0 * (v1_ll * B3_ll - v3_ll * B1_ll + v1_rr * B3_rr - v3_rr * B1_rr)
398
# total energy flux is complicated and involves the previous components
399
f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * inv_gamma_minus_one)
400
+
401
0.5f0 * (+p_ll * v1_rr + p_rr * v1_ll
402
+ (v1_ll * B2_ll * B2_rr + v1_rr * B2_rr * B2_ll)
403
+ (v1_ll * B3_ll * B3_rr + v1_rr * B3_rr * B3_ll)
404
-
405
(v2_ll * B1_ll * B2_rr + v2_rr * B1_rr * B2_ll)
406
-
407
(v3_ll * B1_ll * B3_rr + v3_rr * B1_rr * B3_ll)))
408
409
f_other = SVector(f2, f3, f4, f5, f6, f7, f8)
410
411
return vcat(f_other, f_rho)
412
end
413
414
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation
415
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
416
equations::IdealGlmMhdMulticomponentEquations1D)
417
rho_v1_ll, _ = u_ll
418
rho_v1_rr, _ = u_rr
419
420
rho_ll = density(u_ll, equations)
421
rho_rr = density(u_rr, equations)
422
423
# Calculate velocities (ignore orientation since it is always "1" in 1D)
424
# and fast magnetoacoustic wave speeds
425
# left
426
v_ll = rho_v1_ll / rho_ll
427
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
428
# right
429
v_rr = rho_v1_rr / rho_rr
430
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
431
432
return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)
433
end
434
435
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
436
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
437
equations::IdealGlmMhdMulticomponentEquations1D)
438
rho_v1_ll, _ = u_ll
439
rho_v1_rr, _ = u_rr
440
441
rho_ll = density(u_ll, equations)
442
rho_rr = density(u_rr, equations)
443
444
# Calculate velocities (ignore orientation since it is always "1" in 1D)
445
# and fast magnetoacoustic wave speeds
446
# left
447
v_ll = rho_v1_ll / rho_ll
448
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
449
# right
450
v_rr = rho_v1_rr / rho_rr
451
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
452
453
return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)
454
end
455
456
@inline function max_abs_speeds(u, equations::IdealGlmMhdMulticomponentEquations1D)
457
rho_v1, _ = u
458
459
rho = density(u, equations)
460
461
v1 = rho_v1 / rho
462
463
cf_x_direction = calc_fast_wavespeed(u, 1, equations)
464
465
return (abs(v1) + cf_x_direction,)
466
end
467
468
# Convert conservative variables to primitive
469
function cons2prim(u, equations::IdealGlmMhdMulticomponentEquations1D)
470
rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u
471
472
prim_rho = SVector{ncomponents(equations), real(equations)}(u[i + 7]
473
for i in eachcomponent(equations))
474
rho = density(u, equations)
475
476
v1 = rho_v1 / rho
477
v2 = rho_v2 / rho
478
v3 = rho_v3 / rho
479
480
gamma = totalgamma(u, equations)
481
482
p = (gamma - 1) *
483
(rho_e_total - 0.5f0 * rho * (v1^2 + v2^2 + v3^2) -
484
0.5f0 * (B1^2 + B2^2 + B3^2))
485
prim_other = SVector(v1, v2, v3, p, B1, B2, B3)
486
487
return vcat(prim_other, prim_rho)
488
end
489
490
# Convert conservative variables to entropy
491
@inline function cons2entropy(u, equations::IdealGlmMhdMulticomponentEquations1D)
492
rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u
493
@unpack cv, gammas, gas_constants = equations
494
495
rho = density(u, equations)
496
497
v1 = rho_v1 / rho
498
v2 = rho_v2 / rho
499
v3 = rho_v3 / rho
500
v_square = v1^2 + v2^2 + v3^2
501
gamma = totalgamma(u, equations)
502
p = (gamma - 1) *
503
(rho_e_total - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2))
504
s = log(p) - gamma * log(rho)
505
rho_p = rho / p
506
507
# Multicomponent stuff
508
RealT = eltype(u)
509
help1 = zero(RealT)
510
511
for i in eachcomponent(equations)
512
help1 += u[i + 7] * cv[i]
513
end
514
515
T = (rho_e_total - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2)) / (help1)
516
517
entrop_rho = SVector{ncomponents(equations), real(equations)}(-1 *
518
(cv[i] * log(T) -
519
gas_constants[i] *
520
log(u[i + 7])) +
521
gas_constants[i] +
522
cv[i] -
523
(v_square / (2 * T))
524
for i in eachcomponent(equations))
525
526
w1 = v1 / T
527
w2 = v2 / T
528
w3 = v3 / T
529
w4 = -1 / T
530
w5 = B1 / T
531
w6 = B2 / T
532
w7 = B3 / T
533
534
entrop_other = SVector(w1, w2, w3, w4, w5, w6, w7)
535
536
return vcat(entrop_other, entrop_rho)
537
end
538
539
# Convert primitive to conservative variables
540
@inline function prim2cons(prim, equations::IdealGlmMhdMulticomponentEquations1D)
541
v1, v2, v3, p, B1, B2, B3 = prim
542
543
cons_rho = SVector{ncomponents(equations), real(equations)}(prim[i + 7]
544
for i in eachcomponent(equations))
545
rho = density(prim, equations)
546
547
rho_v1 = rho * v1
548
rho_v2 = rho * v2
549
rho_v3 = rho * v3
550
551
gamma = totalgamma(prim, equations)
552
rho_e_total = p / (gamma - 1) + 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3) +
553
0.5f0 * (B1^2 + B2^2 + B3^2)
554
555
cons_other = SVector(rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3)
556
557
return vcat(cons_other, cons_rho)
558
end
559
560
@doc raw"""
561
pressure(u, equations::IdealGlmMhdMulticomponentEquations1D)
562
563
Computes the pressure for an ideal equation of state with
564
isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`.
565
```math
566
\begin{aligned}
567
p &= (\gamma - 1) \left( \rho e_{\text{total}} - \rho e_{\text{kinetic}} - \rho e_{\text{magnetic}} \right) \\
568
&= (\gamma - 1) \left( \rho e_{\text{total}} - \frac{1}{2}
569
\left[\rho \Vert v \Vert_2^2 + \Vert B \Vert_2^2 \right] \right)
570
\end{aligned}
571
```
572
"""
573
@inline function pressure(u, equations::IdealGlmMhdMulticomponentEquations1D)
574
rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u
575
rho = density(u, equations)
576
gamma = totalgamma(u, equations)
577
p = (gamma - 1) * (rho_e_total -
578
0.5f0 *
579
((rho_v1^2 + rho_v2^2 + rho_v3^2) / rho +
580
B1^2 + B2^2 + B3^2))
581
return p
582
end
583
584
@inline function density_pressure(u, equations::IdealGlmMhdMulticomponentEquations1D)
585
rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u
586
rho = density(u, equations)
587
gamma = totalgamma(u, equations)
588
rho_times_p = (gamma - 1) * (rho * rho_e_total -
589
0.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2 +
590
rho * (B1^2 + B2^2 + B3^2)))
591
return rho_times_p
592
end
593
594
# Compute the fastest wave speed for ideal MHD equations: c_f, the fast magnetoacoustic eigenvalue
595
@inline function calc_fast_wavespeed(cons, direction,
596
equations::IdealGlmMhdMulticomponentEquations1D)
597
rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = cons
598
rho = density(cons, equations)
599
v1 = rho_v1 / rho
600
v2 = rho_v2 / rho
601
v3 = rho_v3 / rho
602
v_mag = sqrt(v1^2 + v2^2 + v3^2)
603
gamma = totalgamma(cons, equations)
604
p = (gamma - 1) *
605
(rho_e_total - 0.5f0 * rho * v_mag^2 - 0.5f0 * (B1^2 + B2^2 + B3^2))
606
a_square = gamma * p / rho
607
sqrt_rho = sqrt(rho)
608
b1 = B1 / sqrt_rho
609
b2 = B2 / sqrt_rho
610
b3 = B3 / sqrt_rho
611
b_square = b1^2 + b2^2 + b3^2
612
613
c_f = sqrt(0.5f0 * (a_square + b_square) +
614
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b1^2))
615
616
return c_f
617
end
618
619
@doc raw"""
620
density(u, equations::AbstractIdealGlmMhdMulticomponentEquations)
621
622
Computes the total density ``\rho = \sum_{i=1}^n \rho_i`` from the conserved variables `u`.
623
"""
624
@inline function density(u, equations::IdealGlmMhdMulticomponentEquations1D)
625
RealT = eltype(u)
626
rho = zero(RealT)
627
628
for i in eachcomponent(equations)
629
rho += u[i + 7]
630
end
631
632
return rho
633
end
634
635
@inline function totalgamma(u, equations::IdealGlmMhdMulticomponentEquations1D)
636
@unpack cv, gammas = equations
637
638
RealT = eltype(u)
639
help1 = zero(RealT)
640
help2 = zero(RealT)
641
642
for i in eachcomponent(equations)
643
help1 += u[i + 7] * cv[i] * gammas[i]
644
help2 += u[i + 7] * cv[i]
645
end
646
647
return help1 / help2
648
end
649
650
@inline function densities(u, v, equations::IdealGlmMhdMulticomponentEquations1D)
651
return SVector{ncomponents(equations), real(equations)}(u[i + 7] * v
652
for i in eachcomponent(equations))
653
end
654
end # @muladd
655
656