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_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
IdealGlmMhdEquations1D(gamma)
10
11
The ideal compressible GLM-MHD equations
12
```math
13
\frac{\partial}{\partial t}
14
\begin{pmatrix}
15
\rho \\ \rho v_1 \\ \rho v_2 \\ \rho v_3 \\ \rho e_{\text{total}} \\ B_1 \\ B_2 \\ B_3
16
\end{pmatrix}
17
+
18
\frac{\partial}{\partial x}
19
\begin{pmatrix}
20
\rho v_1 \\ \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
23
\end{pmatrix}
24
=
25
\begin{pmatrix}
26
0 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0 \\ 0
27
\end{pmatrix}
28
```
29
for an ideal gas in one space dimension.
30
Here, ``\mathbf{v}`` is the velocity, ``\mathbf{B}`` the magnetic field, ``e_{\text{total}}`` the specific total energy, and
31
```math
32
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)
33
```
34
the pressure and ``\gamma`` the heat capacity ratio.
35
36
!!! note
37
There is no divergence cleaning variable `psi` because the divergence-free constraint
38
is satisfied trivially in one spatial dimension.
39
"""
40
struct IdealGlmMhdEquations1D{RealT <: Real} <: AbstractIdealGlmMhdEquations{1, 8}
41
gamma::RealT # ratio of specific heats
42
inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications
43
44
function IdealGlmMhdEquations1D(gamma)
45
γ, inv_gamma_minus_one = promote(gamma, inv(gamma - 1))
46
return new{typeof(γ)}(γ, inv_gamma_minus_one)
47
end
48
end
49
50
have_nonconservative_terms(::IdealGlmMhdEquations1D) = False()
51
function varnames(::typeof(cons2cons), ::IdealGlmMhdEquations1D)
52
return ("rho", "rho_v1", "rho_v2", "rho_v3", "rho_e_total", "B1", "B2", "B3")
53
end
54
function varnames(::typeof(cons2prim), ::IdealGlmMhdEquations1D)
55
return ("rho", "v1", "v2", "v3", "p", "B1", "B2", "B3")
56
end
57
function default_analysis_integrals(::IdealGlmMhdEquations1D)
58
return (entropy_timederivative, Val(:l2_divb), Val(:linf_divb))
59
end
60
61
"""
62
initial_condition_constant(x, t, equations::IdealGlmMhdEquations1D)
63
64
A constant initial condition to test free-stream preservation.
65
"""
66
function initial_condition_constant(x, t, equations::IdealGlmMhdEquations1D)
67
RealT = eltype(x)
68
rho = 1
69
rho_v1 = convert(RealT, 0.1)
70
rho_v2 = -convert(RealT, 0.2)
71
rho_v3 = -0.5f0
72
rho_e_total = 50
73
B1 = 3
74
B2 = -convert(RealT, 1.2)
75
B3 = 0.5f0
76
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3)
77
end
78
79
"""
80
initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations1D)
81
82
An Alfvén wave as smooth initial condition used for convergence tests.
83
See for reference section 4.2 in
84
- Dominik Derigs, Andrew R. Winters, Gregor J. Gassner, and Stefanie Walch (2016)
85
A novel high-order, entropy stable, 3D AMR MHD solver with guaranteed positive pressure
86
[DOI: 10.1016/j.jcp.2016.04.048](https://doi.org/10.1016/j.jcp.2016.04.048)
87
"""
88
function initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations1D)
89
# smooth Alfvén wave test from Derigs et al. (2016)
90
# domain must be set to [0, 1], γ = 5/3
91
RealT = eltype(x)
92
rho = 1
93
v1 = 0
94
si, co = sincospi(2 * (x[1] + t)) # Adding `t` makes non-integer time valid
95
v2 = convert(RealT, 0.1) * si
96
v3 = convert(RealT, 0.1) * co
97
p = convert(RealT, 0.1)
98
B1 = 1
99
B2 = v2
100
B3 = v3
101
return prim2cons(SVector(rho, v1, v2, v3, p, B1, B2, B3), equations)
102
end
103
104
"""
105
initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations1D)
106
107
A weak blast wave adapted from
108
- Sebastian Hennemann, Gregor J. Gassner (2020)
109
A provably entropy stable subcell shock capturing approach for high order split form DG
110
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
111
"""
112
function initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations1D)
113
# Adapted MHD version of the weak blast wave from Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
114
# Same discontinuity in the velocities but with magnetic fields
115
# Set up polar coordinates
116
RealT = eltype(x)
117
inicenter = (0,)
118
x_norm = x[1] - inicenter[1]
119
r = sqrt(x_norm^2)
120
phi = atan(x_norm)
121
122
# Calculate primitive variables
123
rho = r > 0.5f0 ? one(RealT) : convert(RealT, 1.1691)
124
v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos(phi)
125
p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)
126
127
return prim2cons(SVector(rho, v1, 0, 0, p, 1, 1, 1, 0), equations)
128
end
129
130
# Calculate 1D flux for a single point
131
@inline function flux(u, orientation::Integer, equations::IdealGlmMhdEquations1D)
132
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u
133
v1 = rho_v1 / rho
134
v2 = rho_v2 / rho
135
v3 = rho_v3 / rho
136
kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
137
mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)
138
p_over_gamma_minus_one = (rho_e_total - kin_en - mag_en)
139
p = (equations.gamma - 1) * p_over_gamma_minus_one
140
141
# Ignore orientation since it is always "1" in 1D
142
f1 = rho_v1
143
f2 = rho_v1 * v1 + p + mag_en - B1^2
144
f3 = rho_v1 * v2 - B1 * B2
145
f4 = rho_v1 * v3 - B1 * B3
146
f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v1 -
147
B1 * (v1 * B1 + v2 * B2 + v3 * B3)
148
f6 = 0
149
f7 = v1 * B2 - v2 * B1
150
f8 = v1 * B3 - v3 * B1
151
152
return SVector(f1, f2, f3, f4, f5, f6, f7, f8)
153
end
154
155
"""
156
flux_derigs_etal(u_ll, u_rr, orientation, equations::IdealGlmMhdEquations1D)
157
158
Entropy conserving two-point flux by
159
- Derigs et al. (2018)
160
Ideal GLM-MHD: About the entropy consistent nine-wave magnetic field
161
divergence diminishing ideal magnetohydrodynamics equations
162
[DOI: 10.1016/j.jcp.2018.03.002](https://doi.org/10.1016/j.jcp.2018.03.002)
163
"""
164
function flux_derigs_etal(u_ll, u_rr, orientation::Integer,
165
equations::IdealGlmMhdEquations1D)
166
# Unpack left and right states to get velocities, pressure, and inverse temperature (called beta)
167
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll = u_ll
168
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr = u_rr
169
170
v1_ll = rho_v1_ll / rho_ll
171
v2_ll = rho_v2_ll / rho_ll
172
v3_ll = rho_v3_ll / rho_ll
173
v1_rr = rho_v1_rr / rho_rr
174
v2_rr = rho_v2_rr / rho_rr
175
v3_rr = rho_v3_rr / rho_rr
176
vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2
177
vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2
178
mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^2
179
mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^2
180
p_ll = (equations.gamma - 1) *
181
(rho_e_total_ll - 0.5f0 * rho_ll * vel_norm_ll - 0.5f0 * mag_norm_ll)
182
p_rr = (equations.gamma - 1) *
183
(rho_e_total_rr - 0.5f0 * rho_rr * vel_norm_rr - 0.5f0 * mag_norm_rr)
184
beta_ll = 0.5f0 * rho_ll / p_ll
185
beta_rr = 0.5f0 * rho_rr / p_rr
186
# for convenience store v⋅B
187
vel_dot_mag_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
188
vel_dot_mag_rr = v1_rr * B1_rr + v2_rr * B2_rr + v3_rr * B3_rr
189
190
# Compute the necessary mean values needed for either direction
191
rho_avg = 0.5f0 * (rho_ll + rho_rr)
192
rho_mean = ln_mean(rho_ll, rho_rr)
193
beta_mean = ln_mean(beta_ll, beta_rr)
194
beta_avg = 0.5f0 * (beta_ll + beta_rr)
195
v1_avg = 0.5f0 * (v1_ll + v1_rr)
196
v2_avg = 0.5f0 * (v2_ll + v2_rr)
197
v3_avg = 0.5f0 * (v3_ll + v3_rr)
198
p_mean = 0.5f0 * rho_avg / beta_avg
199
B1_avg = 0.5f0 * (B1_ll + B1_rr)
200
B2_avg = 0.5f0 * (B2_ll + B2_rr)
201
B3_avg = 0.5f0 * (B3_ll + B3_rr)
202
vel_norm_avg = 0.5f0 * (vel_norm_ll + vel_norm_rr)
203
mag_norm_avg = 0.5f0 * (mag_norm_ll + mag_norm_rr)
204
vel_dot_mag_avg = 0.5f0 * (vel_dot_mag_ll + vel_dot_mag_rr)
205
206
# Ignore orientation since it is always "1" in 1D
207
f1 = rho_mean * v1_avg
208
f2 = f1 * v1_avg + p_mean + 0.5f0 * mag_norm_avg - B1_avg * B1_avg
209
f3 = f1 * v2_avg - B1_avg * B2_avg
210
f4 = f1 * v3_avg - B1_avg * B3_avg
211
f6 = 0
212
f7 = v1_avg * B2_avg - v2_avg * B1_avg
213
f8 = v1_avg * B3_avg - v3_avg * B1_avg
214
# total energy flux is complicated and involves the previous eight components
215
v1_mag_avg = 0.5f0 * (v1_ll * mag_norm_ll + v1_rr * mag_norm_rr)
216
f5 = (f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +
217
f2 * v1_avg + f3 * v2_avg +
218
f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg - 0.5f0 * v1_mag_avg +
219
B1_avg * vel_dot_mag_avg)
220
221
return SVector(f1, f2, f3, f4, f5, f6, f7, f8)
222
end
223
224
"""
225
flux_hindenlang_gassner(u_ll, u_rr, orientation_or_normal_direction,
226
equations::IdealGlmMhdEquations1D)
227
228
Entropy conserving and kinetic energy preserving two-point flux of
229
Hindenlang and Gassner (2019), extending [`flux_ranocha`](@ref) to the MHD equations.
230
231
## References
232
- Florian Hindenlang, Gregor Gassner (2019)
233
A new entropy conservative two-point flux for ideal MHD equations derived from
234
first principles.
235
Presented at HONOM 2019: European workshop on high order numerical methods
236
for evolutionary PDEs, theory and applications
237
- Hendrik Ranocha (2018)
238
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
239
for Hyperbolic Balance Laws
240
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
241
- Hendrik Ranocha (2020)
242
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
243
the Euler Equations Using Summation-by-Parts Operators
244
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
245
"""
246
@inline function flux_hindenlang_gassner(u_ll, u_rr, orientation::Integer,
247
equations::IdealGlmMhdEquations1D)
248
# Unpack left and right states
249
rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll = cons2prim(u_ll, equations)
250
rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr = cons2prim(u_rr, equations)
251
252
# Compute the necessary mean values needed for either direction
253
rho_mean = ln_mean(rho_ll, rho_rr)
254
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
255
# in exact arithmetic since
256
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
257
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
258
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
259
v1_avg = 0.5f0 * (v1_ll + v1_rr)
260
v2_avg = 0.5f0 * (v2_ll + v2_rr)
261
v3_avg = 0.5f0 * (v3_ll + v3_rr)
262
p_avg = 0.5f0 * (p_ll + p_rr)
263
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
264
magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)
265
266
# Calculate fluxes depending on orientation with specific direction averages
267
f1 = rho_mean * v1_avg
268
f2 = f1 * v1_avg + p_avg + magnetic_square_avg -
269
0.5f0 * (B1_ll * B1_rr + B1_rr * B1_ll)
270
f3 = f1 * v2_avg - 0.5f0 * (B1_ll * B2_rr + B1_rr * B2_ll)
271
f4 = f1 * v3_avg - 0.5f0 * (B1_ll * B3_rr + B1_rr * B3_ll)
272
#f5 below
273
f6 = 0
274
f7 = 0.5f0 * (v1_ll * B2_ll - v2_ll * B1_ll + v1_rr * B2_rr - v2_rr * B1_rr)
275
f8 = 0.5f0 * (v1_ll * B3_ll - v3_ll * B1_ll + v1_rr * B3_rr - v3_rr * B1_rr)
276
# total energy flux is complicated and involves the previous components
277
f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
278
+
279
0.5f0 * (+p_ll * v1_rr + p_rr * v1_ll
280
+ (v1_ll * B2_ll * B2_rr + v1_rr * B2_rr * B2_ll)
281
+ (v1_ll * B3_ll * B3_rr + v1_rr * B3_rr * B3_ll)
282
-
283
(v2_ll * B1_ll * B2_rr + v2_rr * B1_rr * B2_ll)
284
-
285
(v3_ll * B1_ll * B3_rr + v3_rr * B1_rr * B3_ll)))
286
287
return SVector(f1, f2, f3, f4, f5, f6, f7, f8)
288
end
289
290
"""
291
flux_hllc(u_ll, u_rr, orientation, equations::IdealGlmMhdEquations1D)
292
293
- Li (2005)
294
An HLLC Riemann solver for magneto-hydrodynamics
295
[DOI: 10.1016/j.jcp.2004.08.020](https://doi.org/10.1016/j.jcp.2004.08.020).
296
"""
297
function flux_hllc(u_ll, u_rr, orientation::Integer,
298
equations::IdealGlmMhdEquations1D)
299
# Unpack left and right states
300
rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll = cons2prim(u_ll, equations)
301
rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr = cons2prim(u_rr, equations)
302
303
# Total pressure, i.e., thermal + magnetic pressures (eq. (12))
304
p_tot_ll = p_ll + 0.5f0 * (B1_ll^2 + B2_ll^2 + B3_ll^2)
305
p_tot_rr = p_rr + 0.5f0 * (B1_rr^2 + B2_rr^2 + B3_rr^2)
306
307
# Conserved variables
308
rho_v1_ll = u_ll[2]
309
rho_v2_ll = u_ll[3]
310
rho_v3_ll = u_ll[4]
311
312
rho_v1_rr = u_rr[2]
313
rho_v2_rr = u_rr[3]
314
rho_v3_rr = u_rr[4]
315
316
# Obtain left and right fluxes
317
f_ll = flux(u_ll, orientation, equations)
318
f_rr = flux(u_rr, orientation, equations)
319
320
SsL, SsR = min_max_speed_einfeldt(u_ll, u_rr, orientation, equations)
321
sMu_L = SsL - v1_ll
322
sMu_R = SsR - v1_rr
323
if SsL >= 0
324
f1 = f_ll[1]
325
f2 = f_ll[2]
326
f3 = f_ll[3]
327
f4 = f_ll[4]
328
f5 = f_ll[5]
329
f6 = f_ll[6]
330
f7 = f_ll[7]
331
f8 = f_ll[8]
332
elseif SsR <= 0
333
f1 = f_rr[1]
334
f2 = f_rr[2]
335
f3 = f_rr[3]
336
f4 = f_rr[4]
337
f5 = f_rr[5]
338
f6 = f_rr[6]
339
f7 = f_rr[7]
340
f8 = f_rr[8]
341
else
342
# Compute the "HLLC-speed", eq. (14) from paper mentioned above
343
#=
344
SStar = (rho_rr * v1_rr * sMu_R - rho_ll * v1_ll * sMu_L + p_tot_ll - p_tot_rr - B1_ll^2 + B1_rr^2 ) /
345
(rho_rr * sMu_R - rho_ll * sMu_L)
346
=#
347
# Simplification for 1D: B1 is constant
348
SStar = (rho_rr * v1_rr * sMu_R - rho_ll * v1_ll * sMu_L + p_tot_ll - p_tot_rr) /
349
(rho_rr * sMu_R - rho_ll * sMu_L)
350
351
Sdiff = SsR - SsL
352
353
# Compute HLL values for vStar, BStar
354
# These correspond to eq. (28) and (30) from the referenced paper
355
# and the classic HLL intermediate state given by (2)
356
rho_HLL = (SsR * rho_rr - SsL * rho_ll - (f_rr[1] - f_ll[1])) / Sdiff
357
358
v1Star = (SsR * rho_v1_rr - SsL * rho_v1_ll - (f_rr[2] - f_ll[2])) /
359
(Sdiff * rho_HLL)
360
v2Star = (SsR * rho_v2_rr - SsL * rho_v2_ll - (f_rr[3] - f_ll[3])) /
361
(Sdiff * rho_HLL)
362
v3Star = (SsR * rho_v3_rr - SsL * rho_v3_ll - (f_rr[4] - f_ll[4])) /
363
(Sdiff * rho_HLL)
364
365
#B1Star = (SsR * B1_rr - SsL * B1_ll - (f_rr[6] - f_ll[6])) / Sdiff
366
# 1D B1 = constant => B1_ll = B1_rr = B1Star
367
B1Star = B1_ll
368
369
B2Star = (SsR * B2_rr - SsL * B2_ll - (f_rr[7] - f_ll[7])) / Sdiff
370
B3Star = (SsR * B3_rr - SsL * B3_ll - (f_rr[8] - f_ll[8])) / Sdiff
371
if SsL <= SStar
372
SdiffStar = SsL - SStar
373
374
densStar = rho_ll * sMu_L / SdiffStar # (19)
375
376
mom_1_Star = densStar * SStar # (20)
377
mom_2_Star = densStar * v2_ll -
378
(B1Star * B2Star - B1_ll * B2_ll) / SdiffStar # (21)
379
mom_3_Star = densStar * v3_ll -
380
(B1Star * B3Star - B1_ll * B3_ll) / SdiffStar # (22)
381
382
#p_tot_Star = rho_ll * sMu_L * (SStar - v1_ll) + p_tot_ll - B1_ll^2 + B1Star^2 # (17)
383
# 1D B1 = constant => B1_ll = B1_rr = B1Star
384
p_tot_Star = rho_ll * sMu_L * (SStar - v1_ll) + p_tot_ll # (17)
385
386
enerStar = u_ll[5] * sMu_L / SdiffStar +
387
(p_tot_Star * SStar - p_tot_ll * v1_ll - (B1Star *
388
(B1Star * v1Star + B2Star * v2Star + B3Star * v3Star) -
389
B1_ll * (B1_ll * v1_ll + B2_ll * v2_ll + B3_ll * v3_ll))) /
390
SdiffStar # (23)
391
392
# Classic HLLC update (32)
393
f1 = f_ll[1] + SsL * (densStar - u_ll[1])
394
f2 = f_ll[2] + SsL * (mom_1_Star - u_ll[2])
395
f3 = f_ll[3] + SsL * (mom_2_Star - u_ll[3])
396
f4 = f_ll[4] + SsL * (mom_3_Star - u_ll[4])
397
f5 = f_ll[5] + SsL * (enerStar - u_ll[5])
398
f6 = f_ll[6] + SsL * (B1Star - u_ll[6])
399
f7 = f_ll[7] + SsL * (B2Star - u_ll[7])
400
f8 = f_ll[8] + SsL * (B3Star - u_ll[8])
401
else # SStar <= Ssr
402
SdiffStar = SsR - SStar
403
404
densStar = rho_rr * sMu_R / SdiffStar # (19)
405
406
mom_1_Star = densStar * SStar # (20)
407
mom_2_Star = densStar * v2_rr -
408
(B1Star * B2Star - B1_rr * B2_rr) / SdiffStar # (21)
409
mom_3_Star = densStar * v3_rr -
410
(B1Star * B3Star - B1_rr * B3_rr) / SdiffStar # (22)
411
412
#p_tot_Star = rho_rr * sMu_R * (SStar - v1_rr) + p_tot_rr - B1_rr^2 + B1Star^2 # (17)
413
# 1D B1 = constant => B1_ll = B1_rr = B1Star
414
p_tot_Star = rho_rr * sMu_R * (SStar - v1_rr) + p_tot_rr # (17)
415
416
enerStar = u_rr[5] * sMu_R / SdiffStar +
417
(p_tot_Star * SStar - p_tot_rr * v1_rr - (B1Star *
418
(B1Star * v1Star + B2Star * v2Star + B3Star * v3Star) -
419
B1_rr * (B1_rr * v1_rr + B2_rr * v2_rr + B3_rr * v3_rr))) /
420
SdiffStar # (23)
421
422
# Classic HLLC update (32)
423
f1 = f_rr[1] + SsR * (densStar - u_rr[1])
424
f2 = f_rr[2] + SsR * (mom_1_Star - u_rr[2])
425
f3 = f_rr[3] + SsR * (mom_2_Star - u_rr[3])
426
f4 = f_rr[4] + SsR * (mom_3_Star - u_rr[4])
427
f5 = f_rr[5] + SsR * (enerStar - u_rr[5])
428
f6 = f_rr[6] + SsR * (B1Star - u_rr[6])
429
f7 = f_rr[7] + SsR * (B2Star - u_rr[7])
430
f8 = f_rr[8] + SsR * (B3Star - u_rr[8])
431
end
432
end
433
return SVector(f1, f2, f3, f4, f5, f6, f7, f8)
434
end
435
436
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation
437
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
438
equations::IdealGlmMhdEquations1D)
439
rho_ll, rho_v1_ll, _ = u_ll
440
rho_rr, rho_v1_rr, _ = u_rr
441
442
# Calculate velocities (ignore orientation since it is always "1" in 1D)
443
# and fast magnetoacoustic wave speeds
444
# left
445
v_ll = rho_v1_ll / rho_ll
446
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
447
# right
448
v_rr = rho_v1_rr / rho_rr
449
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
450
451
return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)
452
end
453
454
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
455
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
456
equations::IdealGlmMhdEquations1D)
457
rho_ll, rho_v1_ll, _ = u_ll
458
rho_rr, rho_v1_rr, _ = u_rr
459
460
# Calculate velocities (ignore orientation since it is always "1" in 1D)
461
# and fast magnetoacoustic wave speeds
462
# left
463
v_ll = rho_v1_ll / rho_ll
464
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
465
# right
466
v_rr = rho_v1_rr / rho_rr
467
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
468
469
return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)
470
end
471
472
# Calculate estimates for minimum and maximum wave speeds for HLL-type fluxes
473
@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,
474
equations::IdealGlmMhdEquations1D)
475
rho_ll, rho_v1_ll, _ = u_ll
476
rho_rr, rho_v1_rr, _ = u_rr
477
478
# Calculate primitive variables
479
v1_ll = rho_v1_ll / rho_ll
480
v1_rr = rho_v1_rr / rho_rr
481
482
λ_min = v1_ll - calc_fast_wavespeed(u_ll, orientation, equations)
483
λ_max = v1_rr + calc_fast_wavespeed(u_rr, orientation, equations)
484
485
return λ_min, λ_max
486
end
487
488
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
489
@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,
490
equations::IdealGlmMhdEquations1D)
491
rho_ll, rho_v1_ll, _ = u_ll
492
rho_rr, rho_v1_rr, _ = u_rr
493
494
# Calculate primitive variables
495
v1_ll = rho_v1_ll / rho_ll
496
v1_rr = rho_v1_rr / rho_rr
497
498
# Approximate the left-most and right-most eigenvalues in the Riemann fan
499
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
500
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
501
502
λ_min = min(v1_ll - c_f_ll, v1_rr - c_f_rr)
503
λ_max = max(v1_ll + c_f_ll, v1_rr + c_f_rr)
504
505
return λ_min, λ_max
506
end
507
508
"""
509
min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer, equations::IdealGlmMhdEquations1D)
510
511
Calculate minimum and maximum wave speeds for HLL-type fluxes as in
512
- Li (2005)
513
An HLLC Riemann solver for magneto-hydrodynamics
514
[DOI: 10.1016/j.jcp.2004.08.020](https://doi.org/10.1016/j.jcp.2004.08.020).
515
"""
516
@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,
517
equations::IdealGlmMhdEquations1D)
518
rho_ll, rho_v1_ll, _ = u_ll
519
rho_rr, rho_v1_rr, _ = u_rr
520
521
# Calculate primitive variables
522
v1_ll = rho_v1_ll / rho_ll
523
v1_rr = rho_v1_rr / rho_rr
524
525
# Approximate the left-most and right-most eigenvalues in the Riemann fan
526
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
527
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
528
vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)
529
λ_min = min(v1_ll - c_f_ll, vel_roe - c_f_roe)
530
λ_max = max(v1_rr + c_f_rr, vel_roe + c_f_roe)
531
532
return λ_min, λ_max
533
end
534
535
@inline function max_abs_speeds(u, equations::IdealGlmMhdEquations1D)
536
rho, rho_v1, _ = u
537
v1 = rho_v1 / rho
538
cf_x_direction = calc_fast_wavespeed(u, 1, equations)
539
540
return abs(v1) + cf_x_direction
541
end
542
543
# Convert conservative variables to primitive
544
@inline function cons2prim(u, equations::IdealGlmMhdEquations1D)
545
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u
546
547
v1 = rho_v1 / rho
548
v2 = rho_v2 / rho
549
v3 = rho_v3 / rho
550
p = (equations.gamma - 1) * (rho_e_total -
551
0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3
552
+ B1 * B1 + B2 * B2 + B3 * B3))
553
554
return SVector(rho, v1, v2, v3, p, B1, B2, B3)
555
end
556
557
# Convert conservative variables to entropy
558
@inline function cons2entropy(u, equations::IdealGlmMhdEquations1D)
559
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u
560
561
v1 = rho_v1 / rho
562
v2 = rho_v2 / rho
563
v3 = rho_v3 / rho
564
v_square = v1^2 + v2^2 + v3^2
565
p = (equations.gamma - 1) *
566
(rho_e_total - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2))
567
s = log(p) - equations.gamma * log(rho)
568
rho_p = rho / p
569
570
w1 = (equations.gamma - s) / (equations.gamma - 1) - 0.5f0 * rho_p * v_square
571
w2 = rho_p * v1
572
w3 = rho_p * v2
573
w4 = rho_p * v3
574
w5 = -rho_p
575
w6 = rho_p * B1
576
w7 = rho_p * B2
577
w8 = rho_p * B3
578
579
return SVector(w1, w2, w3, w4, w5, w6, w7, w8)
580
end
581
582
# Convert primitive to conservative variables
583
@inline function prim2cons(prim, equations::IdealGlmMhdEquations1D)
584
rho, v1, v2, v3, p, B1, B2, B3 = prim
585
586
rho_v1 = rho * v1
587
rho_v2 = rho * v2
588
rho_v3 = rho * v3
589
rho_e_total = p / (equations.gamma - 1) +
590
0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3) +
591
0.5f0 * (B1^2 + B2^2 + B3^2)
592
593
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3)
594
end
595
596
@inline function density(u, equations::IdealGlmMhdEquations1D)
597
rho = u[1]
598
return rho
599
end
600
601
@inline function velocity(u, equations::IdealGlmMhdEquations1D)
602
rho = u[1]
603
v1 = u[2] / rho
604
v2 = u[3] / rho
605
v3 = u[4] / rho
606
return SVector(v1, v2, v3)
607
end
608
609
@inline function velocity(u, orientation::Int, equations::IdealGlmMhdEquations1D)
610
rho = u[1]
611
v = u[orientation + 1] / rho
612
return v
613
end
614
615
@doc raw"""
616
pressure(u, equations::IdealGlmMhdEquations1D)
617
618
Computes the pressure for an ideal equation of state with
619
isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`.
620
```math
621
\begin{aligned}
622
p &= (\gamma - 1) \left( \rho e_{\text{total}} - \rho e_{\text{kinetic}} - \rho e_{\text{magnetic}} \right) \\
623
&= (\gamma - 1) \left( \rho e - \frac{1}{2}
624
\left[\rho \Vert v \Vert_2^2 + \Vert B \Vert_2^2 \right] \right)
625
\end{aligned}
626
```
627
"""
628
@inline function pressure(u, equations::IdealGlmMhdEquations1D)
629
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u
630
p = (equations.gamma - 1) * (rho_e_total -
631
0.5f0 *
632
((rho_v1^2 + rho_v2^2 + rho_v3^2) / rho +
633
B1^2 + B2^2 + B3^2))
634
return p
635
end
636
637
@inline function density_pressure(u, equations::IdealGlmMhdEquations1D)
638
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = u
639
rho_times_p = (equations.gamma - 1) * (rho * rho_e_total -
640
0.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2 +
641
rho * (B1^2 + B2^2 + B3^2)))
642
return rho_times_p
643
end
644
645
# Compute the fastest wave speed for ideal MHD equations: c_f, the fast magnetoacoustic eigenvalue
646
@inline function calc_fast_wavespeed(cons, direction, equations::IdealGlmMhdEquations1D)
647
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3 = cons
648
v1 = rho_v1 / rho
649
v2 = rho_v2 / rho
650
v3 = rho_v3 / rho
651
v_mag = sqrt(v1^2 + v2^2 + v3^2)
652
p = (equations.gamma - 1) *
653
(rho_e_total - 0.5f0 * rho * v_mag^2 - 0.5f0 * (B1^2 + B2^2 + B3^2))
654
a_square = equations.gamma * p / rho
655
sqrt_rho = sqrt(rho)
656
b1 = B1 / sqrt_rho
657
b2 = B2 / sqrt_rho
658
b3 = B3 / sqrt_rho
659
b_square = b1^2 + b2^2 + b3^2
660
661
c_f = sqrt(0.5f0 * (a_square + b_square) +
662
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b1^2))
663
return c_f
664
end
665
666
"""
667
calc_fast_wavespeed_roe(u_ll, u_rr, direction, equations::IdealGlmMhdEquations1D)
668
669
Compute the fast magnetoacoustic wave speed using Roe averages
670
as given by
671
- Cargo and Gallice (1997)
672
Roe Matrices for Ideal MHD and Systematic Construction
673
of Roe Matrices for Systems of Conservation Laws
674
[DOI: 10.1006/jcph.1997.5773](https://doi.org/10.1006/jcph.1997.5773)
675
"""
676
@inline function calc_fast_wavespeed_roe(u_ll, u_rr, direction,
677
equations::IdealGlmMhdEquations1D)
678
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll = u_ll
679
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr = u_rr
680
681
# Calculate primitive variables
682
v1_ll = rho_v1_ll / rho_ll
683
v2_ll = rho_v2_ll / rho_ll
684
v3_ll = rho_v3_ll / rho_ll
685
vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2
686
mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^2
687
p_ll = (equations.gamma - 1) *
688
(rho_e_total_ll - 0.5f0 * rho_ll * vel_norm_ll - 0.5f0 * mag_norm_ll)
689
690
v1_rr = rho_v1_rr / rho_rr
691
v2_rr = rho_v2_rr / rho_rr
692
v3_rr = rho_v3_rr / rho_rr
693
vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2
694
mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^2
695
p_rr = (equations.gamma - 1) *
696
(rho_e_total_rr - 0.5f0 * rho_rr * vel_norm_rr - 0.5f0 * mag_norm_rr)
697
698
# compute total pressure which is thermal + magnetic pressures
699
p_total_ll = p_ll + 0.5f0 * mag_norm_ll
700
p_total_rr = p_rr + 0.5f0 * mag_norm_rr
701
702
# compute the Roe density averages
703
sqrt_rho_ll = sqrt(rho_ll)
704
sqrt_rho_rr = sqrt(rho_rr)
705
inv_sqrt_rho_add = 1 / (sqrt_rho_ll + sqrt_rho_rr)
706
inv_sqrt_rho_prod = 1 / (sqrt_rho_ll * sqrt_rho_rr)
707
rho_ll_roe = sqrt_rho_ll * inv_sqrt_rho_add
708
rho_rr_roe = sqrt_rho_rr * inv_sqrt_rho_add
709
# Roe averages
710
# velocities and magnetic fields
711
v1_roe = v1_ll * rho_ll_roe + v1_rr * rho_rr_roe
712
v2_roe = v2_ll * rho_ll_roe + v2_rr * rho_rr_roe
713
v3_roe = v3_ll * rho_ll_roe + v3_rr * rho_rr_roe
714
B1_roe = B1_ll * rho_ll_roe + B1_rr * rho_rr_roe
715
B2_roe = B2_ll * rho_ll_roe + B2_rr * rho_rr_roe
716
B3_roe = B3_ll * rho_ll_roe + B3_rr * rho_rr_roe
717
# enthalpy
718
H_ll = (rho_e_total_ll + p_total_ll) / rho_ll
719
H_rr = (rho_e_total_rr + p_total_rr) / rho_rr
720
H_roe = H_ll * rho_ll_roe + H_rr * rho_rr_roe
721
# temporary variable see equations (4.12) in Cargo and Gallice
722
X = 0.5f0 * ((B1_ll - B1_rr)^2 + (B2_ll - B2_rr)^2 + (B3_ll - B3_rr)^2) *
723
inv_sqrt_rho_add^2
724
# averaged components needed to compute c_f, the fast magnetoacoustic wave speed
725
b_square_roe = (B1_roe^2 + B2_roe^2 + B3_roe^2) * inv_sqrt_rho_prod # scaled magnectic sum
726
a_square_roe = ((2 - equations.gamma) * X +
727
(equations.gamma - 1) *
728
(H_roe - 0.5f0 * (v1_roe^2 + v2_roe^2 + v3_roe^2) -
729
b_square_roe)) # acoustic speed
730
# finally compute the average wave speed and set the output velocity
731
# Ignore orientation since it is always "1" in 1D
732
c_a_roe = B1_roe^2 * inv_sqrt_rho_prod # (squared) Alfvén wave speed
733
a_star_roe = sqrt((a_square_roe + b_square_roe)^2 - 4 * a_square_roe * c_a_roe)
734
c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe))
735
736
return v1_roe, c_f_roe
737
end
738
739
@doc raw"""
740
entropy_thermodynamic(cons, equations::AbstractIdealGlmMhdEquations)
741
742
Calculate thermodynamic entropy for a conservative state `cons` as
743
744
```math
745
s = \log(p) - \gamma \log(\rho)
746
```
747
"""
748
@inline function entropy_thermodynamic(cons, equations::IdealGlmMhdEquations1D)
749
# Pressure
750
p = (equations.gamma - 1) *
751
(cons[5] - 0.5f0 * (cons[2]^2 + cons[3]^2 + cons[4]^2) / cons[1]
752
-
753
0.5f0 * (cons[6]^2 + cons[7]^2 + cons[8]^2))
754
755
# Thermodynamic entropy
756
s = log(p) - equations.gamma * log(cons[1])
757
758
return s
759
end
760
761
@doc raw"""
762
entropy_math(cons, equations::AbstractIdealGlmMhdEquations)
763
764
Calculate mathematical entropy for a conservative state `cons` as
765
```math
766
S = -\frac{\rho s}{\gamma - 1}
767
```
768
where `s` is the thermodynamic entropy calculated by
769
[`entropy_thermodynamic(cons, equations::AbstractIdealGlmMhdEquations)`](@ref).
770
"""
771
@inline function entropy_math(cons, equations::IdealGlmMhdEquations1D)
772
S = -entropy_thermodynamic(cons, equations) * cons[1] / (equations.gamma - 1)
773
774
return S
775
end
776
777
"""
778
entropy(cons, equations::AbstractIdealGlmMhdEquations)
779
780
Default entropy is the mathematical entropy
781
[`entropy_math(cons, equations::AbstractIdealGlmMhdEquations)`](@ref).
782
"""
783
@inline entropy(cons, equations::IdealGlmMhdEquations1D) = entropy_math(cons, equations)
784
785
# Calculate total energy for a conservative state `cons`
786
@inline energy_total(cons, ::IdealGlmMhdEquations1D) = cons[5]
787
788
# Calculate kinetic energy for a conservative state `cons`
789
@inline function energy_kinetic(cons, equations::IdealGlmMhdEquations1D)
790
return 0.5f0 * (cons[2]^2 + cons[3]^2 + cons[4]^2) / cons[1]
791
end
792
793
# Calculate the magnetic energy for a conservative state `cons`.
794
# OBS! For non-dinmensional form of the ideal MHD magnetic pressure ≡ magnetic energy
795
@inline function energy_magnetic(cons, ::IdealGlmMhdEquations1D)
796
return 0.5f0 * (cons[6]^2 + cons[7]^2 + cons[8]^2)
797
end
798
799
"""
800
energy_internal(cons, equations::AbstractIdealGlmMhdEquations)
801
802
Calculate internal energy for a conservative state `cons` as the difference
803
between total energy and kinetic + magnetic energies.
804
"""
805
@inline function energy_internal(cons, equations::IdealGlmMhdEquations1D)
806
return (energy_total(cons, equations)
807
-
808
energy_kinetic(cons, equations)
809
-
810
energy_magnetic(cons, equations))
811
end
812
813
# Calculate the cross helicity (\vec{v}⋅\vec{B}) for a conservative state `cons`
814
@inline function cross_helicity(cons, ::IdealGlmMhdEquations1D)
815
return (cons[2] * cons[6] + cons[3] * cons[7] + cons[4] * cons[8]) / cons[1]
816
end
817
end # @muladd
818
819