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_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
IdealGlmMhdMulticomponentEquations2D
10
11
The ideal compressible multicomponent GLM-MHD equations in two space dimensions.
12
"""
13
struct IdealGlmMhdMulticomponentEquations2D{NVARS, NCOMP, RealT <: Real} <:
14
AbstractIdealGlmMhdMulticomponentEquations{2, NVARS, NCOMP}
15
gammas::SVector{NCOMP, RealT}
16
gas_constants::SVector{NCOMP, RealT}
17
cv::SVector{NCOMP, RealT}
18
cp::SVector{NCOMP, RealT}
19
c_h::RealT # GLM cleaning speed
20
21
function IdealGlmMhdMulticomponentEquations2D{NVARS, NCOMP, RealT}(gammas::SVector{NCOMP,
22
RealT},
23
gas_constants::SVector{NCOMP,
24
RealT},
25
c_h::RealT) where {
26
NVARS,
27
NCOMP,
28
RealT <:
29
Real
30
}
31
NCOMP >= 1 ||
32
throw(DimensionMismatch("`gammas` and `gas_constants` have to be filled with at least one value"))
33
34
cv = gas_constants ./ (gammas .- 1)
35
cp = gas_constants + gas_constants ./ (gammas .- 1)
36
37
new(gammas, gas_constants, cv, cp, c_h)
38
end
39
end
40
41
function IdealGlmMhdMulticomponentEquations2D(; gammas, gas_constants)
42
_gammas = promote(gammas...)
43
_gas_constants = promote(gas_constants...)
44
RealT = promote_type(eltype(_gammas), eltype(_gas_constants))
45
__gammas = SVector(map(RealT, _gammas))
46
__gas_constants = SVector(map(RealT, _gas_constants))
47
48
c_h = convert(RealT, NaN)
49
50
NVARS = length(_gammas) + 8
51
NCOMP = length(_gammas)
52
53
return IdealGlmMhdMulticomponentEquations2D{NVARS, NCOMP, RealT}(__gammas,
54
__gas_constants,
55
c_h)
56
end
57
58
# Outer constructor for `@reset` works correctly
59
function IdealGlmMhdMulticomponentEquations2D(gammas, gas_constants, cv, cp, c_h)
60
_gammas = promote(gammas...)
61
_gas_constants = promote(gas_constants...)
62
RealT = promote_type(eltype(_gammas), eltype(_gas_constants))
63
__gammas = SVector(map(RealT, _gammas))
64
__gas_constants = SVector(map(RealT, _gas_constants))
65
66
c_h = convert(RealT, c_h)
67
68
NVARS = length(_gammas) + 8
69
NCOMP = length(_gammas)
70
71
return IdealGlmMhdMulticomponentEquations2D{NVARS, NCOMP, RealT}(__gammas,
72
__gas_constants,
73
c_h)
74
end
75
76
@inline function Base.real(::IdealGlmMhdMulticomponentEquations2D{NVARS, NCOMP, RealT}) where {
77
NVARS,
78
NCOMP,
79
RealT
80
}
81
RealT
82
end
83
84
have_nonconservative_terms(::IdealGlmMhdMulticomponentEquations2D) = True()
85
86
function varnames(::typeof(cons2cons), equations::IdealGlmMhdMulticomponentEquations2D)
87
cons = ("rho_v1", "rho_v2", "rho_v3", "rho_e", "B1", "B2", "B3", "psi")
88
rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))
89
return (cons..., rhos...)
90
end
91
92
function varnames(::typeof(cons2prim), equations::IdealGlmMhdMulticomponentEquations2D)
93
prim = ("v1", "v2", "v3", "p", "B1", "B2", "B3", "psi")
94
rhos = ntuple(n -> "rho" * string(n), Val(ncomponents(equations)))
95
return (prim..., rhos...)
96
end
97
98
function default_analysis_integrals(::IdealGlmMhdMulticomponentEquations2D)
99
(entropy_timederivative, Val(:l2_divb), Val(:linf_divb))
100
end
101
102
# Helper function to extract the magnetic field vector from the conservative variables
103
magnetic_field(u, equations::IdealGlmMhdMulticomponentEquations2D) = SVector(u[5], u[6],
104
u[7])
105
106
"""
107
initial_condition_convergence_test(x, t, equations::IdealGlmMhdMulticomponentEquations2D)
108
109
An Alfvén wave as smooth initial condition used for convergence tests.
110
"""
111
function initial_condition_convergence_test(x, t,
112
equations::IdealGlmMhdMulticomponentEquations2D)
113
# smooth Alfvén wave test from Derigs et al. FLASH (2016)
114
# domain must be set to [0, 1/cos(α)] x [0, 1/sin(α)], γ = 5/3
115
RealT = eltype(x)
116
alpha = 0.25f0 * convert(RealT, pi)
117
x_perp = x[1] * cos(alpha) + x[2] * sin(alpha)
118
B_perp = convert(RealT, 0.1) * sinpi(2 * x_perp)
119
rho = one(RealT)
120
prim_rho = SVector{ncomponents(equations), real(equations)}(2^(i - 1) * (1 - 2) *
121
rho / (1 -
122
2^ncomponents(equations))
123
for i in eachcomponent(equations))
124
125
v1 = -B_perp * sin(alpha)
126
v2 = B_perp * cos(alpha)
127
v3 = convert(RealT, 0.1) * cospi(2 * x_perp)
128
p = convert(RealT, 0.1)
129
B1 = cos(alpha) + v1
130
B2 = sin(alpha) + v2
131
B3 = v3
132
psi = 0
133
prim_other = SVector(v1, v2, v3, p, B1, B2, B3, psi)
134
135
return prim2cons(vcat(prim_other, prim_rho), equations)
136
end
137
138
"""
139
initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdMulticomponentEquations2D)
140
141
A weak blast wave adapted from
142
- Sebastian Hennemann, Gregor J. Gassner (2020)
143
A provably entropy stable subcell shock capturing approach for high order split form DG
144
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
145
"""
146
function initial_condition_weak_blast_wave(x, t,
147
equations::IdealGlmMhdMulticomponentEquations2D)
148
# Adapted MHD version of the weak blast wave from Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
149
# Same discontinuity in the velocities but with magnetic fields
150
# Set up polar coordinates
151
RealT = eltype(x)
152
inicenter = SVector(0, 0)
153
x_norm = x[1] - inicenter[1]
154
y_norm = x[2] - inicenter[2]
155
r = sqrt(x_norm^2 + y_norm^2)
156
phi = atan(y_norm, x_norm)
157
sin_phi, cos_phi = sincos(phi)
158
159
# We initialize each species with a fraction of the total density `rho`, such
160
# that the sum of the densities is `rho := density(prim, equations)`. The density of
161
# a species is double the density of the next species.
162
prim_rho = SVector{ncomponents(equations), real(equations)}(r > 0.5f0 ?
163
2^(i - 1) * (1 - 2) /
164
(RealT(1) -
165
2^ncomponents(equations)) :
166
2^(i - 1) * (1 - 2) *
167
RealT(1.1691) /
168
(1 -
169
2^ncomponents(equations))
170
for i in eachcomponent(equations))
171
172
v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos_phi
173
v2 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * sin_phi
174
p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)
175
176
prim_other = SVector(v1, v2, 0, p, 1, 1, 1, 0)
177
178
return prim2cons(vcat(prim_other, prim_rho), equations)
179
end
180
181
# Calculate 1D flux for a single point
182
@inline function flux(u, orientation::Integer,
183
equations::IdealGlmMhdMulticomponentEquations2D)
184
rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
185
@unpack c_h = equations
186
187
rho = density(u, equations)
188
189
v1 = rho_v1 / rho
190
v2 = rho_v2 / rho
191
v3 = rho_v3 / rho
192
kin_en = 0.5f0 * rho * (v1^2 + v2^2 + v3^2)
193
mag_en = 0.5f0 * (B1^2 + B2^2 + B3^2)
194
gamma = totalgamma(u, equations)
195
p = (gamma - 1) * (rho_e - kin_en - mag_en - 0.5f0 * psi^2)
196
197
if orientation == 1
198
f_rho = densities(u, v1, equations)
199
f1 = rho_v1 * v1 + p + mag_en - B1^2
200
f2 = rho_v1 * v2 - B1 * B2
201
f3 = rho_v1 * v3 - B1 * B3
202
f4 = (kin_en + gamma * p / (gamma - 1) + 2 * mag_en) * v1 -
203
B1 * (v1 * B1 + v2 * B2 + v3 * B3) + c_h * psi * B1
204
f5 = c_h * psi
205
f6 = v1 * B2 - v2 * B1
206
f7 = v1 * B3 - v3 * B1
207
f8 = c_h * B1
208
else # orientation == 2
209
f_rho = densities(u, v2, equations)
210
f1 = rho_v2 * v1 - B1 * B2
211
f2 = rho_v2 * v2 + p + mag_en - B2^2
212
f3 = rho_v2 * v3 - B2 * B3
213
f4 = (kin_en + gamma * p / (gamma - 1) + 2 * mag_en) * v2 -
214
B2 * (v1 * B1 + v2 * B2 + v3 * B3) + c_h * psi * B2
215
f5 = v2 * B1 - v1 * B2
216
f6 = c_h * psi
217
f7 = v2 * B3 - v3 * B2
218
f8 = c_h * B2
219
end
220
221
f_other = SVector(f1, f2, f3, f4, f5, f6, f7, f8)
222
223
return vcat(f_other, f_rho)
224
end
225
226
"""
227
flux_nonconservative_powell(u_ll, u_rr, orientation::Integer,
228
equations::IdealGlmMhdMulticomponentEquations2D)
229
230
Non-symmetric two-point flux discretizing the nonconservative (source) term of
231
Powell and the Galilean nonconservative term associated with the GLM multiplier
232
of the [`IdealGlmMhdMulticomponentEquations2D`](@ref).
233
234
## References
235
- Marvin Bohm, Andrew R.Winters, Gregor J. Gassner, Dominik Derigs,
236
Florian Hindenlang, Joachim Saur
237
An entropy stable nodal discontinuous Galerkin method for the resistive MHD
238
equations. Part I: Theory and numerical verification
239
[DOI: 10.1016/j.jcp.2018.06.027](https://doi.org/10.1016/j.jcp.2018.06.027)
240
"""
241
@inline function flux_nonconservative_powell(u_ll, u_rr, orientation::Integer,
242
equations::IdealGlmMhdMulticomponentEquations2D)
243
rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
244
rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
245
246
rho_ll = density(u_ll, equations)
247
248
v1_ll = rho_v1_ll / rho_ll
249
v2_ll = rho_v2_ll / rho_ll
250
v3_ll = rho_v3_ll / rho_ll
251
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
252
253
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
254
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
255
# Note that the order of conserved variables is changed compared to the
256
# standard GLM MHD equations, i.e., the densities are moved to the end
257
# Here, we compute the non-density components at first and append zero density
258
# components afterwards
259
zero_densities = SVector{ncomponents(equations), real(equations)}(ntuple(_ -> zero(real(equations)),
260
Val(ncomponents(equations))))
261
if orientation == 1
262
f = SVector(B1_ll * B1_rr,
263
B2_ll * B1_rr,
264
B3_ll * B1_rr,
265
v_dot_B_ll * B1_rr + v1_ll * psi_ll * psi_rr,
266
v1_ll * B1_rr,
267
v2_ll * B1_rr,
268
v3_ll * B1_rr,
269
v1_ll * psi_rr)
270
else # orientation == 2
271
f = SVector(B1_ll * B2_rr,
272
B2_ll * B2_rr,
273
B3_ll * B2_rr,
274
v_dot_B_ll * B2_rr + v2_ll * psi_ll * psi_rr,
275
v1_ll * B2_rr,
276
v2_ll * B2_rr,
277
v3_ll * B2_rr,
278
v2_ll * psi_rr)
279
end
280
281
return vcat(f, zero_densities)
282
end
283
284
"""
285
flux_derigs_etal(u_ll, u_rr, orientation, equations::IdealGlmMhdMulticomponentEquations2D)
286
287
Entropy conserving two-point flux adapted by
288
- Derigs et al. (2018)
289
Ideal GLM-MHD: About the entropy consistent nine-wave magnetic field
290
divergence diminishing ideal magnetohydrodynamics equations for multicomponent
291
[DOI: 10.1016/j.jcp.2018.03.002](https://doi.org/10.1016/j.jcp.2018.03.002)
292
"""
293
function flux_derigs_etal(u_ll, u_rr, orientation::Integer,
294
equations::IdealGlmMhdMulticomponentEquations2D)
295
# Unpack left and right states to get velocities, pressure, and inverse temperature (called beta)
296
rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
297
rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
298
@unpack gammas, gas_constants, cv, c_h = equations
299
300
rho_ll = density(u_ll, equations)
301
rho_rr = density(u_rr, equations)
302
303
gamma_ll = totalgamma(u_ll, equations)
304
gamma_rr = totalgamma(u_rr, equations)
305
306
rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 8],
307
u_rr[i + 8])
308
for i in eachcomponent(equations))
309
rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 8] +
310
u_rr[i + 8])
311
for i in eachcomponent(equations))
312
313
v1_ll = rho_v1_ll / rho_ll
314
v2_ll = rho_v2_ll / rho_ll
315
v3_ll = rho_v3_ll / rho_ll
316
v1_rr = rho_v1_rr / rho_rr
317
v2_rr = rho_v2_rr / rho_rr
318
v3_rr = rho_v3_rr / rho_rr
319
v1_sq = 0.5f0 * (v1_ll^2 + v1_rr^2)
320
v2_sq = 0.5f0 * (v2_ll^2 + v2_rr^2)
321
v3_sq = 0.5f0 * (v3_ll^2 + v3_rr^2)
322
v_sq = v1_sq + v2_sq + v3_sq
323
B1_sq = 0.5f0 * (B1_ll^2 + B1_rr^2)
324
B2_sq = 0.5f0 * (B2_ll^2 + B2_rr^2)
325
B3_sq = 0.5f0 * (B3_ll^2 + B3_rr^2)
326
B_sq = B1_sq + B2_sq + B3_sq
327
vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2
328
vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2
329
mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^2
330
mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^2
331
# for convenience store v⋅B
332
vel_dot_mag_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
333
vel_dot_mag_rr = v1_rr * B1_rr + v2_rr * B2_rr + v3_rr * B3_rr
334
335
# Compute the necessary mean values needed for either direction
336
v1_avg = 0.5f0 * (v1_ll + v1_rr)
337
v2_avg = 0.5f0 * (v2_ll + v2_rr)
338
v3_avg = 0.5f0 * (v3_ll + v3_rr)
339
v_sum = v1_avg + v2_avg + v3_avg
340
B1_avg = 0.5f0 * (B1_ll + B1_rr)
341
B2_avg = 0.5f0 * (B2_ll + B2_rr)
342
B3_avg = 0.5f0 * (B3_ll + B3_rr)
343
psi_avg = 0.5f0 * (psi_ll + psi_rr)
344
vel_norm_avg = 0.5f0 * (vel_norm_ll + vel_norm_rr)
345
mag_norm_avg = 0.5f0 * (mag_norm_ll + mag_norm_rr)
346
vel_dot_mag_avg = 0.5f0 * (vel_dot_mag_ll + vel_dot_mag_rr)
347
348
RealT = eltype(u_ll)
349
enth = zero(RealT)
350
help1_ll = zero(RealT)
351
help1_rr = zero(RealT)
352
353
for i in eachcomponent(equations)
354
enth += rhok_avg[i] * gas_constants[i]
355
help1_ll += u_ll[i + 8] * cv[i]
356
help1_rr += u_rr[i + 8] * cv[i]
357
end
358
359
T_ll = (rho_e_ll - 0.5f0 * rho_ll * (vel_norm_ll) - 0.5f0 * mag_norm_ll -
360
0.5f0 * psi_ll^2) / help1_ll
361
T_rr = (rho_e_rr - 0.5f0 * rho_rr * (vel_norm_rr) - 0.5f0 * mag_norm_rr -
362
0.5f0 * psi_rr^2) / help1_rr
363
T = 0.5f0 * (1 / T_ll + 1 / T_rr)
364
T_log = ln_mean(1 / T_ll, 1 / T_rr)
365
366
# Calculate fluxes depending on orientation with specific direction averages
367
help1 = zero(RealT)
368
help2 = zero(RealT)
369
if orientation == 1
370
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg
371
for i in eachcomponent(equations))
372
for i in eachcomponent(equations)
373
help1 += f_rho[i] * cv[i]
374
help2 += f_rho[i]
375
end
376
f1 = help2 * v1_avg + enth / T + 0.5f0 * mag_norm_avg - B1_avg * B1_avg
377
f2 = help2 * v2_avg - B1_avg * B2_avg
378
f3 = help2 * v3_avg - B1_avg * B3_avg
379
f5 = c_h * psi_avg
380
f6 = v1_avg * B2_avg - v2_avg * B1_avg
381
f7 = v1_avg * B3_avg - v3_avg * B1_avg
382
f8 = c_h * B1_avg
383
# total energy flux is complicated and involves the previous eight components
384
psi_B1_avg = 0.5f0 * (B1_ll * psi_ll + B1_rr * psi_rr)
385
v1_mag_avg = 0.5f0 * (v1_ll * mag_norm_ll + v1_rr * mag_norm_rr)
386
387
f4 = (help1 / T_log) - 0.5f0 * (vel_norm_avg) * (help2) + f1 * v1_avg +
388
f2 * v2_avg + f3 * v3_avg +
389
f5 * B1_avg + f6 * B2_avg + f7 * B3_avg + f8 * psi_avg -
390
0.5f0 * v1_mag_avg +
391
B1_avg * vel_dot_mag_avg - c_h * psi_B1_avg
392
393
else
394
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v2_avg
395
for i in eachcomponent(equations))
396
for i in eachcomponent(equations)
397
help1 += f_rho[i] * cv[i]
398
help2 += f_rho[i]
399
end
400
f1 = help2 * v1_avg - B1_avg * B2_avg
401
f2 = help2 * v2_avg + enth / T + 0.5f0 * mag_norm_avg - B2_avg * B2_avg
402
f3 = help2 * v3_avg - B2_avg * B3_avg
403
f5 = v2_avg * B1_avg - v1_avg * B2_avg
404
f6 = c_h * psi_avg
405
f7 = v2_avg * B3_avg - v3_avg * B2_avg
406
f8 = c_h * B2_avg
407
408
# total energy flux is complicated and involves the previous eight components
409
psi_B2_avg = 0.5f0 * (B2_ll * psi_ll + B2_rr * psi_rr)
410
v2_mag_avg = 0.5f0 * (v2_ll * mag_norm_ll + v2_rr * mag_norm_rr)
411
412
f4 = (help1 / T_log) - 0.5f0 * (vel_norm_avg) * (help2) + f1 * v1_avg +
413
f2 * v2_avg + f3 * v3_avg +
414
f5 * B1_avg + f6 * B2_avg + f7 * B3_avg + f8 * psi_avg -
415
0.5f0 * v2_mag_avg +
416
B2_avg * vel_dot_mag_avg - c_h * psi_B2_avg
417
end
418
419
f_other = SVector(f1, f2, f3, f4, f5, f6, f7, f8)
420
421
return vcat(f_other, f_rho)
422
end
423
424
"""
425
flux_hindenlang_gassner(u_ll, u_rr, orientation_or_normal_direction,
426
equations::IdealGlmMhdMulticomponentEquations2D)
427
428
Adaption of the entropy conserving and kinetic energy preserving two-point flux of
429
Hindenlang (2019), extending [`flux_ranocha`](@ref) to the MHD equations.
430
## References
431
- Florian Hindenlang, Gregor Gassner (2019)
432
A new entropy conservative two-point flux for ideal MHD equations derived from
433
first principles.
434
Presented at HONOM 2019: European workshop on high order numerical methods
435
for evolutionary PDEs, theory and applications
436
- Hendrik Ranocha (2018)
437
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
438
for Hyperbolic Balance Laws
439
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
440
- Hendrik Ranocha (2020)
441
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
442
the Euler Equations Using Summation-by-Parts Operators
443
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
444
"""
445
@inline function flux_hindenlang_gassner(u_ll, u_rr, orientation::Integer,
446
equations::IdealGlmMhdMulticomponentEquations2D)
447
# Unpack left and right states
448
v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll, psi_ll = cons2prim(u_ll, equations)
449
v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr, psi_rr = cons2prim(u_rr, equations)
450
451
rho_ll = density(u_ll, equations)
452
rho_rr = density(u_rr, equations)
453
454
# Compute the necessary mean values needed for either direction
455
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
456
# in exact arithmetic since
457
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
458
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
459
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
460
v1_avg = 0.5f0 * (v1_ll + v1_rr)
461
v2_avg = 0.5f0 * (v2_ll + v2_rr)
462
v3_avg = 0.5f0 * (v3_ll + v3_rr)
463
p_avg = 0.5f0 * (p_ll + p_rr)
464
psi_avg = 0.5f0 * (psi_ll + psi_rr)
465
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
466
magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)
467
468
inv_gamma_minus_one = 1 / (totalgamma(0.5f0 * (u_ll + u_rr), equations) - 1)
469
470
rhok_mean = SVector{ncomponents(equations), real(equations)}(ln_mean(u_ll[i + 8],
471
u_rr[i + 8])
472
for i in eachcomponent(equations))
473
rhok_avg = SVector{ncomponents(equations), real(equations)}(0.5f0 * (u_ll[i + 8] +
474
u_rr[i + 8])
475
for i in eachcomponent(equations))
476
477
RealT = eltype(u_ll)
478
if orientation == 1
479
f1 = zero(RealT)
480
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v1_avg
481
for i in eachcomponent(equations))
482
for i in eachcomponent(equations)
483
f1 += f_rho[i]
484
end
485
486
# Calculate fluxes depending on orientation with specific direction averages
487
f2 = f1 * v1_avg + p_avg + magnetic_square_avg -
488
0.5f0 * (B1_ll * B1_rr + B1_rr * B1_ll)
489
f3 = f1 * v2_avg - 0.5f0 * (B1_ll * B2_rr + B1_rr * B2_ll)
490
f4 = f1 * v3_avg - 0.5f0 * (B1_ll * B3_rr + B1_rr * B3_ll)
491
# f5 below
492
f6 = f6 = equations.c_h * psi_avg
493
f7 = 0.5f0 * (v1_ll * B2_ll - v2_ll * B1_ll + v1_rr * B2_rr - v2_rr * B1_rr)
494
f8 = 0.5f0 * (v1_ll * B3_ll - v3_ll * B1_ll + v1_rr * B3_rr - v3_rr * B1_rr)
495
f9 = equations.c_h * 0.5f0 * (B1_ll + B1_rr)
496
# total energy flux is complicated and involves the previous components
497
f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * inv_gamma_minus_one)
498
+
499
0.5f0 * (+p_ll * v1_rr + p_rr * v1_ll
500
+ (v1_ll * B2_ll * B2_rr + v1_rr * B2_rr * B2_ll)
501
+ (v1_ll * B3_ll * B3_rr + v1_rr * B3_rr * B3_ll)
502
-
503
(v2_ll * B1_ll * B2_rr + v2_rr * B1_rr * B2_ll)
504
-
505
(v3_ll * B1_ll * B3_rr + v3_rr * B1_rr * B3_ll)
506
+
507
equations.c_h * (B1_ll * psi_rr + B1_rr * psi_ll)))
508
else
509
f1 = zero(RealT)
510
f_rho = SVector{ncomponents(equations), real(equations)}(rhok_mean[i] * v2_avg
511
for i in eachcomponent(equations))
512
for i in eachcomponent(equations)
513
f1 += f_rho[i]
514
end
515
516
# Calculate fluxes depending on orientation with specific direction averages
517
f2 = f1 * v1_avg - 0.5f0 * (B2_ll * B1_rr + B2_rr * B1_ll)
518
f3 = f1 * v2_avg + p_avg + magnetic_square_avg -
519
0.5f0 * (B2_ll * B2_rr + B2_rr * B2_ll)
520
f4 = f1 * v3_avg - 0.5f0 * (B2_ll * B3_rr + B2_rr * B3_ll)
521
#f5 below
522
f6 = 0.5f0 * (v2_ll * B1_ll - v1_ll * B2_ll + v2_rr * B1_rr - v1_rr * B2_rr)
523
f7 = equations.c_h * psi_avg
524
f8 = 0.5f0 * (v2_ll * B3_ll - v3_ll * B2_ll + v2_rr * B3_rr - v3_rr * B2_rr)
525
f9 = equations.c_h * 0.5f0 * (B2_ll + B2_rr)
526
# total energy flux is complicated and involves the previous components
527
f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * inv_gamma_minus_one)
528
+
529
0.5f0 * (+p_ll * v2_rr + p_rr * v2_ll
530
+ (v2_ll * B1_ll * B1_rr + v2_rr * B1_rr * B1_ll)
531
+ (v2_ll * B3_ll * B3_rr + v2_rr * B3_rr * B3_ll)
532
-
533
(v1_ll * B2_ll * B1_rr + v1_rr * B2_rr * B1_ll)
534
-
535
(v3_ll * B2_ll * B3_rr + v3_rr * B2_rr * B3_ll)
536
+
537
equations.c_h * (B2_ll * psi_rr + B2_rr * psi_ll)))
538
end
539
540
f_other = SVector(f2, f3, f4, f5, f6, f7, f8, f9)
541
542
return vcat(f_other, f_rho)
543
end
544
545
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation
546
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
547
equations::IdealGlmMhdMulticomponentEquations2D)
548
rho_v1_ll, rho_v2_ll, _ = u_ll
549
rho_v1_rr, rho_v2_rr, _ = u_rr
550
551
rho_ll = density(u_ll, equations)
552
rho_rr = density(u_rr, equations)
553
554
# Calculate velocities and fast magnetoacoustic wave speeds
555
if orientation == 1
556
v_ll = rho_v1_ll / rho_ll
557
v_rr = rho_v1_rr / rho_rr
558
else # orientation == 2
559
v_ll = rho_v2_ll / rho_ll
560
v_rr = rho_v2_rr / rho_rr
561
end
562
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
563
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
564
565
return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)
566
end
567
568
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
569
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
570
equations::IdealGlmMhdMulticomponentEquations2D)
571
rho_v1_ll, rho_v2_ll, _ = u_ll
572
rho_v1_rr, rho_v2_rr, _ = u_rr
573
574
rho_ll = density(u_ll, equations)
575
rho_rr = density(u_rr, equations)
576
577
# Calculate velocities and fast magnetoacoustic wave speeds
578
if orientation == 1
579
v_ll = rho_v1_ll / rho_ll
580
v_rr = rho_v1_rr / rho_rr
581
else # orientation == 2
582
v_ll = rho_v2_ll / rho_ll
583
v_rr = rho_v2_rr / rho_rr
584
end
585
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
586
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
587
588
return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)
589
end
590
591
@inline function max_abs_speeds(u, equations::IdealGlmMhdMulticomponentEquations2D)
592
rho_v1, rho_v2, _ = u
593
594
rho = density(u, equations)
595
596
v1 = rho_v1 / rho
597
v2 = rho_v2 / rho
598
599
cf_x_direction = calc_fast_wavespeed(u, 1, equations)
600
cf_y_direction = calc_fast_wavespeed(u, 2, equations)
601
602
return (abs(v1) + cf_x_direction, abs(v2) + cf_y_direction)
603
end
604
605
@inline function density_pressure(u, equations::IdealGlmMhdMulticomponentEquations2D)
606
rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
607
rho = density(u, equations)
608
gamma = totalgamma(u, equations)
609
p = (gamma - 1) * (rho_e - 0.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2) / rho
610
-
611
0.5f0 * (B1^2 + B2^2 + B3^2)
612
-
613
0.5f0 * psi^2)
614
return rho * p
615
end
616
617
# Convert conservative variables to primitive
618
function cons2prim(u, equations::IdealGlmMhdMulticomponentEquations2D)
619
rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
620
621
prim_rho = SVector{ncomponents(equations), real(equations)}(u[i + 8]
622
for i in eachcomponent(equations))
623
rho = density(u, equations)
624
625
v1 = rho_v1 / rho
626
v2 = rho_v2 / rho
627
v3 = rho_v3 / rho
628
629
gamma = totalgamma(u, equations)
630
631
p = (gamma - 1) *
632
(rho_e - 0.5f0 * rho * (v1^2 + v2^2 + v3^2) - 0.5f0 * (B1^2 + B2^2 + B3^2) -
633
0.5f0 * psi^2)
634
prim_other = SVector(v1, v2, v3, p, B1, B2, B3, psi)
635
636
return vcat(prim_other, prim_rho)
637
end
638
639
# Convert conservative variables to entropy
640
@inline function cons2entropy(u, equations::IdealGlmMhdMulticomponentEquations2D)
641
rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
642
@unpack cv, gammas, gas_constants = equations
643
644
rho = density(u, equations)
645
646
v1 = rho_v1 / rho
647
v2 = rho_v2 / rho
648
v3 = rho_v3 / rho
649
v_square = v1^2 + v2^2 + v3^2
650
gamma = totalgamma(u, equations)
651
p = (gamma - 1) *
652
(rho_e - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2) - 0.5f0 * psi^2)
653
s = log(p) - gamma * log(rho)
654
rho_p = rho / p
655
656
# Multicomponent stuff
657
help1 = zero(v1)
658
659
for i in eachcomponent(equations)
660
help1 += u[i + 8] * cv[i]
661
end
662
663
T = (rho_e - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2) - 0.5f0 * psi^2) /
664
(help1)
665
666
entrop_rho = SVector{ncomponents(equations), real(equations)}(-1 *
667
(cv[i] * log(T) -
668
gas_constants[i] *
669
log(u[i + 8])) +
670
gas_constants[i] +
671
cv[i] -
672
(v_square / (2 * T))
673
for i in eachcomponent(equations))
674
675
w1 = v1 / T
676
w2 = v2 / T
677
w3 = v3 / T
678
w4 = -1 / T
679
w5 = B1 / T
680
w6 = B2 / T
681
w7 = B3 / T
682
w8 = psi / T
683
684
entrop_other = SVector(w1, w2, w3, w4, w5, w6, w7, w8)
685
686
return vcat(entrop_other, entrop_rho)
687
end
688
689
# Convert primitive to conservative variables
690
@inline function prim2cons(prim, equations::IdealGlmMhdMulticomponentEquations2D)
691
v1, v2, v3, p, B1, B2, B3, psi = prim
692
693
cons_rho = SVector{ncomponents(equations), real(equations)}(prim[i + 8]
694
for i in eachcomponent(equations))
695
rho = density(prim, equations)
696
697
rho_v1 = rho * v1
698
rho_v2 = rho * v2
699
rho_v3 = rho * v3
700
701
gamma = totalgamma(prim, equations)
702
rho_e = p / (gamma - 1) + 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3) +
703
0.5f0 * (B1^2 + B2^2 + B3^2) + 0.5f0 * psi^2
704
705
cons_other = SVector(rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3,
706
psi)
707
708
return vcat(cons_other, cons_rho)
709
end
710
711
# Compute the fastest wave speed for ideal MHD equations: c_f, the fast magnetoacoustic eigenvalue
712
@inline function calc_fast_wavespeed(cons, direction,
713
equations::IdealGlmMhdMulticomponentEquations2D)
714
rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = cons
715
rho = density(cons, equations)
716
v1 = rho_v1 / rho
717
v2 = rho_v2 / rho
718
v3 = rho_v3 / rho
719
v_mag = sqrt(v1^2 + v2^2 + v3^2)
720
gamma = totalgamma(cons, equations)
721
p = (gamma - 1) *
722
(rho_e - 0.5f0 * rho * v_mag^2 - 0.5f0 * (B1^2 + B2^2 + B3^2) - 0.5f0 * psi^2)
723
a_square = gamma * p / rho
724
sqrt_rho = sqrt(rho)
725
b1 = B1 / sqrt_rho
726
b2 = B2 / sqrt_rho
727
b3 = B3 / sqrt_rho
728
b_square = b1^2 + b2^2 + b3^2
729
if direction == 1 # x-direction
730
c_f = sqrt(0.5f0 * (a_square + b_square) +
731
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b1^2))
732
else
733
c_f = sqrt(0.5f0 * (a_square + b_square) +
734
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b2^2))
735
end
736
return c_f
737
end
738
739
@inline function density(u, equations::IdealGlmMhdMulticomponentEquations2D)
740
RealT = eltype(u)
741
rho = zero(RealT)
742
743
for i in eachcomponent(equations)
744
rho += u[i + 8]
745
end
746
747
return rho
748
end
749
750
@inline function totalgamma(u, equations::IdealGlmMhdMulticomponentEquations2D)
751
@unpack cv, gammas = equations
752
753
RealT = eltype(u)
754
help1 = zero(RealT)
755
help2 = zero(RealT)
756
757
for i in eachcomponent(equations)
758
help1 += u[i + 8] * cv[i] * gammas[i]
759
help2 += u[i + 8] * cv[i]
760
end
761
762
return help1 / help2
763
end
764
765
@inline function densities(u, v, equations::IdealGlmMhdMulticomponentEquations2D)
766
return SVector{ncomponents(equations), real(equations)}(u[i + 8] * v
767
for i in eachcomponent(equations))
768
end
769
end # @muladd
770
771