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_3d.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
IdealGlmMhdEquations3D(gamma)
10
11
The ideal compressible GLM-MHD equations for an ideal gas with ratio of
12
specific heats `gamma` in three space dimensions.
13
"""
14
struct IdealGlmMhdEquations3D{RealT <: Real} <:
15
AbstractIdealGlmMhdEquations{3, 9}
16
gamma::RealT # ratio of specific heats
17
inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications
18
c_h::RealT # GLM cleaning speed
19
20
function IdealGlmMhdEquations3D(gamma, c_h)
21
γ, inv_gamma_minus_one, c_h = promote(gamma, inv(gamma - 1), c_h)
22
new{typeof(γ)}(γ, inv_gamma_minus_one, c_h)
23
end
24
end
25
26
function IdealGlmMhdEquations3D(gamma; initial_c_h = convert(typeof(gamma), NaN))
27
# Use `promote` to ensure that `gamma` and `initial_c_h` have the same type
28
IdealGlmMhdEquations3D(promote(gamma, initial_c_h)...)
29
end
30
31
# Outer constructor for `@reset` works correctly
32
function IdealGlmMhdEquations3D(gamma, inv_gamma_minus_one, c_h)
33
IdealGlmMhdEquations3D(gamma, c_h)
34
end
35
36
have_nonconservative_terms(::IdealGlmMhdEquations3D) = True()
37
function varnames(::typeof(cons2cons), ::IdealGlmMhdEquations3D)
38
("rho", "rho_v1", "rho_v2", "rho_v3", "rho_e", "B1", "B2", "B3", "psi")
39
end
40
function varnames(::typeof(cons2prim), ::IdealGlmMhdEquations3D)
41
("rho", "v1", "v2", "v3", "p", "B1", "B2", "B3", "psi")
42
end
43
function default_analysis_integrals(::IdealGlmMhdEquations3D)
44
(entropy_timederivative, Val(:l2_divb), Val(:linf_divb))
45
end
46
47
# Helper function to extract the magnetic field vector from the conservative variables
48
magnetic_field(u, equations::IdealGlmMhdEquations3D) = SVector(u[6], u[7], u[8])
49
50
# Set initial conditions at physical location `x` for time `t`
51
"""
52
initial_condition_constant(x, t, equations::IdealGlmMhdEquations3D)
53
54
A constant initial condition to test free-stream preservation.
55
"""
56
function initial_condition_constant(x, t, equations::IdealGlmMhdEquations3D)
57
RealT = eltype(x)
58
rho = 1
59
rho_v1 = convert(RealT, 0.1)
60
rho_v2 = -convert(RealT, 0.2)
61
rho_v3 = -0.5f0
62
rho_e = 50
63
B1 = 3
64
B2 = -convert(RealT, 1.2)
65
B3 = 0.5f0
66
psi = 0
67
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi)
68
end
69
70
"""
71
initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations3D)
72
73
An Alfvén wave as smooth initial condition used for convergence tests.
74
See for reference section 5.1 in
75
- Christoph Altmann (2012)
76
Explicit Discontinuous Galerkin Methods for Magnetohydrodynamics
77
[DOI: 10.18419/opus-3895](http://dx.doi.org/10.18419/opus-3895)
78
"""
79
function initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations3D)
80
# Alfvén wave in three space dimensions
81
# domain must be set to [-1, 1]^3, γ = 5/3
82
RealT = eltype(x)
83
p = 1
84
omega = 2 * convert(RealT, pi) # may be multiplied by frequency
85
# r: length-variable = length of computational domain
86
r = convert(RealT, 2)
87
# e: epsilon = 0.2
88
e = convert(RealT, 0.2)
89
nx = 1 / sqrt(r^2 + 1)
90
ny = r / sqrt(r^2 + 1)
91
sqr = 1
92
Va = omega / (ny * sqr)
93
phi_alv = omega / ny * (nx * (x[1] - 0.5f0 * r) + ny * (x[2] - 0.5f0 * r)) - Va * t
94
95
rho = 1
96
v1 = -e * ny * cos(phi_alv) / rho
97
v2 = e * nx * cos(phi_alv) / rho
98
v3 = e * sin(phi_alv) / rho
99
B1 = nx - rho * v1 * sqr
100
B2 = ny - rho * v2 * sqr
101
B3 = -rho * v3 * sqr
102
psi = 0
103
104
return prim2cons(SVector(rho, v1, v2, v3, p, B1, B2, B3, psi), equations)
105
end
106
107
"""
108
initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations3D)
109
110
A weak blast wave adapted from
111
- Sebastian Hennemann, Gregor J. Gassner (2020)
112
A provably entropy stable subcell shock capturing approach for high order split form DG
113
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
114
"""
115
function initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations3D)
116
# Adapted MHD version of the weak blast wave from Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
117
# Same discontinuity in the velocities but with magnetic fields
118
# Set up polar coordinates
119
RealT = eltype(x)
120
inicenter = (0, 0, 0)
121
x_norm = x[1] - inicenter[1]
122
y_norm = x[2] - inicenter[2]
123
z_norm = x[3] - inicenter[3]
124
r = sqrt(x_norm^2 + y_norm^2 + z_norm^2)
125
phi = atan(y_norm, x_norm)
126
theta = iszero(r) ? zero(RealT) : acos(z_norm / r)
127
128
# Calculate primitive variables
129
rho = r > 0.5f0 ? one(RealT) : convert(RealT, 1.1691)
130
v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos(phi) * sin(theta)
131
v2 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * sin(phi) * sin(theta)
132
v3 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos(theta)
133
p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)
134
135
return prim2cons(SVector(rho, v1, v2, v3, p, 1, 1, 1, 0), equations)
136
end
137
138
# Pre-defined source terms should be implemented as
139
# function source_terms_WHATEVER(u, x, t, equations::IdealGlmMhdEquations3D)
140
141
# Calculate 1D flux for a single point
142
@inline function flux(u, orientation::Integer, equations::IdealGlmMhdEquations3D)
143
rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
144
v1 = rho_v1 / rho
145
v2 = rho_v2 / rho
146
v3 = rho_v3 / rho
147
kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
148
mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)
149
p_over_gamma_minus_one = (rho_e - kin_en - mag_en - 0.5f0 * psi^2)
150
p = (equations.gamma - 1) * p_over_gamma_minus_one
151
if orientation == 1
152
f1 = rho_v1
153
f2 = rho_v1 * v1 + p + mag_en - B1^2
154
f3 = rho_v1 * v2 - B1 * B2
155
f4 = rho_v1 * v3 - B1 * B3
156
f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v1 -
157
B1 * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B1
158
f6 = equations.c_h * psi
159
f7 = v1 * B2 - v2 * B1
160
f8 = v1 * B3 - v3 * B1
161
f9 = equations.c_h * B1
162
elseif orientation == 2
163
f1 = rho_v2
164
f2 = rho_v2 * v1 - B2 * B1
165
f3 = rho_v2 * v2 + p + mag_en - B2^2
166
f4 = rho_v2 * v3 - B2 * B3
167
f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v2 -
168
B2 * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B2
169
f6 = v2 * B1 - v1 * B2
170
f7 = equations.c_h * psi
171
f8 = v2 * B3 - v3 * B2
172
f9 = equations.c_h * B2
173
else
174
f1 = rho_v3
175
f2 = rho_v3 * v1 - B3 * B1
176
f3 = rho_v3 * v2 - B3 * B2
177
f4 = rho_v3 * v3 + p + mag_en - B3^2
178
f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v3 -
179
B3 * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B3
180
f6 = v3 * B1 - v1 * B3
181
f7 = v3 * B2 - v2 * B3
182
f8 = equations.c_h * psi
183
f9 = equations.c_h * B3
184
end
185
186
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
187
end
188
189
# Calculate 1D flux for a single point in the normal direction
190
# Note, this directional vector is not normalized
191
@inline function flux(u, normal_direction::AbstractVector,
192
equations::IdealGlmMhdEquations3D)
193
rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
194
v1 = rho_v1 / rho
195
v2 = rho_v2 / rho
196
v3 = rho_v3 / rho
197
kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
198
mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)
199
p_over_gamma_minus_one = (rho_e - kin_en - mag_en - 0.5f0 * psi^2)
200
p = (equations.gamma - 1) * p_over_gamma_minus_one
201
202
v_normal = v1 * normal_direction[1] + v2 * normal_direction[2] +
203
v3 * normal_direction[3]
204
B_normal = B1 * normal_direction[1] + B2 * normal_direction[2] +
205
B3 * normal_direction[3]
206
rho_v_normal = rho * v_normal
207
208
f1 = rho_v_normal
209
f2 = rho_v_normal * v1 - B1 * B_normal + (p + mag_en) * normal_direction[1]
210
f3 = rho_v_normal * v2 - B2 * B_normal + (p + mag_en) * normal_direction[2]
211
f4 = rho_v_normal * v3 - B3 * B_normal + (p + mag_en) * normal_direction[3]
212
f5 = ((kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v_normal
213
-
214
B_normal * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B_normal)
215
f6 = (equations.c_h * psi * normal_direction[1] +
216
(v2 * B1 - v1 * B2) * normal_direction[2] +
217
(v3 * B1 - v1 * B3) * normal_direction[3])
218
f7 = ((v1 * B2 - v2 * B1) * normal_direction[1] +
219
equations.c_h * psi * normal_direction[2] +
220
(v3 * B2 - v2 * B3) * normal_direction[3])
221
f8 = ((v1 * B3 - v3 * B1) * normal_direction[1] +
222
(v2 * B3 - v3 * B2) * normal_direction[2] +
223
equations.c_h * psi * normal_direction[3])
224
f9 = equations.c_h * B_normal
225
226
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
227
end
228
229
"""
230
flux_nonconservative_powell(u_ll, u_rr, orientation::Integer,
231
equations::IdealGlmMhdEquations3D)
232
flux_nonconservative_powell(u_ll, u_rr,
233
normal_direction::AbstractVector,
234
equations::IdealGlmMhdEquations3D)
235
236
Non-symmetric two-point flux discretizing the nonconservative (source) term of
237
Powell and the Galilean nonconservative term associated with the GLM multiplier
238
of the [`IdealGlmMhdEquations3D`](@ref).
239
240
On curvilinear meshes, the implementation differs from the reference since we use the averaged
241
normal direction for the metrics dealiasing.
242
243
## References
244
- Marvin Bohm, Andrew R.Winters, Gregor J. Gassner, Dominik Derigs,
245
Florian Hindenlang, Joachim Saur
246
An entropy stable nodal discontinuous Galerkin method for the resistive MHD
247
equations. Part I: Theory and numerical verification
248
[DOI: 10.1016/j.jcp.2018.06.027](https://doi.org/10.1016/j.jcp.2018.06.027)
249
"""
250
@inline function flux_nonconservative_powell(u_ll, u_rr, orientation::Integer,
251
equations::IdealGlmMhdEquations3D)
252
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
253
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
254
255
v1_ll = rho_v1_ll / rho_ll
256
v2_ll = rho_v2_ll / rho_ll
257
v3_ll = rho_v3_ll / rho_ll
258
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
259
260
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
261
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2,3}, 0, 0, 0, v_{1,2,3})
262
if orientation == 1
263
f = SVector(0,
264
B1_ll * B1_rr,
265
B2_ll * B1_rr,
266
B3_ll * B1_rr,
267
v_dot_B_ll * B1_rr + v1_ll * psi_ll * psi_rr,
268
v1_ll * B1_rr,
269
v2_ll * B1_rr,
270
v3_ll * B1_rr,
271
v1_ll * psi_rr)
272
elseif orientation == 2
273
f = SVector(0,
274
B1_ll * B2_rr,
275
B2_ll * B2_rr,
276
B3_ll * B2_rr,
277
v_dot_B_ll * B2_rr + v2_ll * psi_ll * psi_rr,
278
v1_ll * B2_rr,
279
v2_ll * B2_rr,
280
v3_ll * B2_rr,
281
v2_ll * psi_rr)
282
else # orientation == 3
283
f = SVector(0,
284
B1_ll * B3_rr,
285
B2_ll * B3_rr,
286
B3_ll * B3_rr,
287
v_dot_B_ll * B3_rr + v3_ll * psi_ll * psi_rr,
288
v1_ll * B3_rr,
289
v2_ll * B3_rr,
290
v3_ll * B3_rr,
291
v3_ll * psi_rr)
292
end
293
294
return f
295
end
296
297
@inline function flux_nonconservative_powell(u_ll, u_rr,
298
normal_direction::AbstractVector,
299
equations::IdealGlmMhdEquations3D)
300
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
301
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
302
303
v1_ll = rho_v1_ll / rho_ll
304
v2_ll = rho_v2_ll / rho_ll
305
v3_ll = rho_v3_ll / rho_ll
306
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
307
308
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2] +
309
v3_ll * normal_direction[3]
310
B_dot_n_rr = B1_rr * normal_direction[1] +
311
B2_rr * normal_direction[2] +
312
B3_rr * normal_direction[3]
313
314
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
315
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2,3}, 0, 0, 0, v_{1,2,3})
316
f = SVector(0,
317
B1_ll * B_dot_n_rr,
318
B2_ll * B_dot_n_rr,
319
B3_ll * B_dot_n_rr,
320
v_dot_B_ll * B_dot_n_rr + v_dot_n_ll * psi_ll * psi_rr,
321
v1_ll * B_dot_n_rr,
322
v2_ll * B_dot_n_rr,
323
v3_ll * B_dot_n_rr,
324
v_dot_n_ll * psi_rr)
325
326
return f
327
end
328
329
"""
330
flux_derigs_etal(u_ll, u_rr, orientation, equations::IdealGlmMhdEquations3D)
331
332
Entropy conserving two-point flux by
333
- Derigs et al. (2018)
334
Ideal GLM-MHD: About the entropy consistent nine-wave magnetic field
335
divergence diminishing ideal magnetohydrodynamics equations
336
[DOI: 10.1016/j.jcp.2018.03.002](https://doi.org/10.1016/j.jcp.2018.03.002)
337
"""
338
function flux_derigs_etal(u_ll, u_rr, orientation::Integer,
339
equations::IdealGlmMhdEquations3D)
340
# Unpack left and right states to get velocities, pressure, and inverse temperature (called beta)
341
rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll, psi_ll = cons2prim(u_ll,
342
equations)
343
rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr, psi_rr = cons2prim(u_rr,
344
equations)
345
346
vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2
347
vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2
348
mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^2
349
mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^2
350
beta_ll = 0.5f0 * rho_ll / p_ll
351
beta_rr = 0.5f0 * rho_rr / p_rr
352
# for convenience store v⋅B
353
vel_dot_mag_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
354
vel_dot_mag_rr = v1_rr * B1_rr + v2_rr * B2_rr + v3_rr * B3_rr
355
356
# Compute the necessary mean values needed for either direction
357
rho_avg = 0.5f0 * (rho_ll + rho_rr)
358
rho_mean = ln_mean(rho_ll, rho_rr)
359
beta_mean = ln_mean(beta_ll, beta_rr)
360
beta_avg = 0.5f0 * (beta_ll + beta_rr)
361
v1_avg = 0.5f0 * (v1_ll + v1_rr)
362
v2_avg = 0.5f0 * (v2_ll + v2_rr)
363
v3_avg = 0.5f0 * (v3_ll + v3_rr)
364
p_mean = 0.5f0 * rho_avg / beta_avg
365
B1_avg = 0.5f0 * (B1_ll + B1_rr)
366
B2_avg = 0.5f0 * (B2_ll + B2_rr)
367
B3_avg = 0.5f0 * (B3_ll + B3_rr)
368
psi_avg = 0.5f0 * (psi_ll + psi_rr)
369
vel_norm_avg = 0.5f0 * (vel_norm_ll + vel_norm_rr)
370
mag_norm_avg = 0.5f0 * (mag_norm_ll + mag_norm_rr)
371
vel_dot_mag_avg = 0.5f0 * (vel_dot_mag_ll + vel_dot_mag_rr)
372
373
# Calculate fluxes depending on orientation with specific direction averages
374
if orientation == 1
375
f1 = rho_mean * v1_avg
376
f2 = f1 * v1_avg + p_mean + 0.5f0 * mag_norm_avg - B1_avg * B1_avg
377
f3 = f1 * v2_avg - B1_avg * B2_avg
378
f4 = f1 * v3_avg - B1_avg * B3_avg
379
f6 = equations.c_h * psi_avg
380
f7 = v1_avg * B2_avg - v2_avg * B1_avg
381
f8 = v1_avg * B3_avg - v3_avg * B1_avg
382
f9 = equations.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
f5 = (f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +
387
f2 * v1_avg + f3 * v2_avg +
388
f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg + f9 * psi_avg -
389
0.5f0 * v1_mag_avg +
390
B1_avg * vel_dot_mag_avg - equations.c_h * psi_B1_avg)
391
elseif orientation == 2
392
f1 = rho_mean * v2_avg
393
f2 = f1 * v1_avg - B2_avg * B1_avg
394
f3 = f1 * v2_avg + p_mean + 0.5f0 * mag_norm_avg - B2_avg * B2_avg
395
f4 = f1 * v3_avg - B2_avg * B3_avg
396
f6 = v2_avg * B1_avg - v1_avg * B2_avg
397
f7 = equations.c_h * psi_avg
398
f8 = v2_avg * B3_avg - v3_avg * B2_avg
399
f9 = equations.c_h * B2_avg
400
# total energy flux is complicated and involves the previous eight components
401
psi_B2_avg = 0.5f0 * (B2_ll * psi_ll + B2_rr * psi_rr)
402
v2_mag_avg = 0.5f0 * (v2_ll * mag_norm_ll + v2_rr * mag_norm_rr)
403
f5 = (f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +
404
f2 * v1_avg + f3 * v2_avg +
405
f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg + f9 * psi_avg -
406
0.5f0 * v2_mag_avg +
407
B2_avg * vel_dot_mag_avg - equations.c_h * psi_B2_avg)
408
else
409
f1 = rho_mean * v3_avg
410
f2 = f1 * v1_avg - B3_avg * B1_avg
411
f3 = f1 * v2_avg - B3_avg * B2_avg
412
f4 = f1 * v3_avg + p_mean + 0.5f0 * mag_norm_avg - B3_avg * B3_avg
413
f6 = v3_avg * B1_avg - v1_avg * B3_avg
414
f7 = v3_avg * B2_avg - v2_avg * B3_avg
415
f8 = equations.c_h * psi_avg
416
f9 = equations.c_h * B3_avg
417
# total energy flux is complicated and involves the previous eight components
418
psi_B3_avg = 0.5f0 * (B3_ll * psi_ll + B3_rr * psi_rr)
419
v3_mag_avg = 0.5f0 * (v3_ll * mag_norm_ll + v3_rr * mag_norm_rr)
420
f5 = (f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +
421
f2 * v1_avg + f3 * v2_avg +
422
f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg + f9 * psi_avg -
423
0.5f0 * v3_mag_avg +
424
B3_avg * vel_dot_mag_avg - equations.c_h * psi_B3_avg)
425
end
426
427
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
428
end
429
430
"""
431
flux_hindenlang_gassner(u_ll, u_rr, orientation_or_normal_direction,
432
equations::IdealGlmMhdEquations3D)
433
434
Entropy conserving and kinetic energy preserving two-point flux of
435
Hindenlang and Gassner (2019), extending [`flux_ranocha`](@ref) to the MHD equations.
436
437
## References
438
- Florian Hindenlang, Gregor Gassner (2019)
439
A new entropy conservative two-point flux for ideal MHD equations derived from
440
first principles.
441
Presented at HONOM 2019: European workshop on high order numerical methods
442
for evolutionary PDEs, theory and applications
443
- Hendrik Ranocha (2018)
444
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
445
for Hyperbolic Balance Laws
446
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
447
- Hendrik Ranocha (2020)
448
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
449
the Euler Equations Using Summation-by-Parts Operators
450
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
451
"""
452
@inline function flux_hindenlang_gassner(u_ll, u_rr, orientation::Integer,
453
equations::IdealGlmMhdEquations3D)
454
# Unpack left and right states
455
rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll, psi_ll = cons2prim(u_ll,
456
equations)
457
rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr, psi_rr = cons2prim(u_rr,
458
equations)
459
460
# Compute the necessary mean values needed for either direction
461
rho_mean = ln_mean(rho_ll, rho_rr)
462
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
463
# in exact arithmetic since
464
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
465
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
466
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
467
v1_avg = 0.5f0 * (v1_ll + v1_rr)
468
v2_avg = 0.5f0 * (v2_ll + v2_rr)
469
v3_avg = 0.5f0 * (v3_ll + v3_rr)
470
p_avg = 0.5f0 * (p_ll + p_rr)
471
psi_avg = 0.5f0 * (psi_ll + psi_rr)
472
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
473
magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)
474
475
# Calculate fluxes depending on orientation with specific direction averages
476
if orientation == 1
477
f1 = rho_mean * v1_avg
478
f2 = f1 * v1_avg + p_avg + magnetic_square_avg -
479
0.5f0 * (B1_ll * B1_rr + B1_rr * B1_ll)
480
f3 = f1 * v2_avg - 0.5f0 * (B1_ll * B2_rr + B1_rr * B2_ll)
481
f4 = f1 * v3_avg - 0.5f0 * (B1_ll * B3_rr + B1_rr * B3_ll)
482
#f5 below
483
f6 = equations.c_h * psi_avg
484
f7 = 0.5f0 * (v1_ll * B2_ll - v2_ll * B1_ll + v1_rr * B2_rr - v2_rr * B1_rr)
485
f8 = 0.5f0 * (v1_ll * B3_ll - v3_ll * B1_ll + v1_rr * B3_rr - v3_rr * B1_rr)
486
f9 = equations.c_h * 0.5f0 * (B1_ll + B1_rr)
487
# total energy flux is complicated and involves the previous components
488
f5 = (f1 *
489
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
490
+
491
0.5f0 * (+p_ll * v1_rr + p_rr * v1_ll
492
+ (v1_ll * B2_ll * B2_rr + v1_rr * B2_rr * B2_ll)
493
+ (v1_ll * B3_ll * B3_rr + v1_rr * B3_rr * B3_ll)
494
-
495
(v2_ll * B1_ll * B2_rr + v2_rr * B1_rr * B2_ll)
496
-
497
(v3_ll * B1_ll * B3_rr + v3_rr * B1_rr * B3_ll)
498
+
499
equations.c_h * (B1_ll * psi_rr + B1_rr * psi_ll)))
500
elseif orientation == 2
501
f1 = rho_mean * v2_avg
502
f2 = f1 * v1_avg - 0.5f0 * (B2_ll * B1_rr + B2_rr * B1_ll)
503
f3 = f1 * v2_avg + p_avg + magnetic_square_avg -
504
0.5f0 * (B2_ll * B2_rr + B2_rr * B2_ll)
505
f4 = f1 * v3_avg - 0.5f0 * (B2_ll * B3_rr + B2_rr * B3_ll)
506
#f5 below
507
f6 = 0.5f0 * (v2_ll * B1_ll - v1_ll * B2_ll + v2_rr * B1_rr - v1_rr * B2_rr)
508
f7 = equations.c_h * psi_avg
509
f8 = 0.5f0 * (v2_ll * B3_ll - v3_ll * B2_ll + v2_rr * B3_rr - v3_rr * B2_rr)
510
f9 = equations.c_h * 0.5f0 * (B2_ll + B2_rr)
511
# total energy flux is complicated and involves the previous components
512
f5 = (f1 *
513
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
514
+
515
0.5f0 * (+p_ll * v2_rr + p_rr * v2_ll
516
+ (v2_ll * B1_ll * B1_rr + v2_rr * B1_rr * B1_ll)
517
+ (v2_ll * B3_ll * B3_rr + v2_rr * B3_rr * B3_ll)
518
-
519
(v1_ll * B2_ll * B1_rr + v1_rr * B2_rr * B1_ll)
520
-
521
(v3_ll * B2_ll * B3_rr + v3_rr * B2_rr * B3_ll)
522
+
523
equations.c_h * (B2_ll * psi_rr + B2_rr * psi_ll)))
524
else # orientation == 3
525
f1 = rho_mean * v3_avg
526
f2 = f1 * v1_avg - 0.5f0 * (B3_ll * B1_rr + B3_rr * B1_ll)
527
f3 = f1 * v2_avg - 0.5f0 * (B3_ll * B2_rr + B3_rr * B2_ll)
528
f4 = f1 * v3_avg + p_avg + magnetic_square_avg -
529
0.5f0 * (B3_ll * B3_rr + B3_rr * B3_ll)
530
#f5 below
531
f6 = 0.5f0 * (v3_ll * B1_ll - v1_ll * B3_ll + v3_rr * B1_rr - v1_rr * B3_rr)
532
f7 = 0.5f0 * (v3_ll * B2_ll - v2_ll * B3_ll + v3_rr * B2_rr - v2_rr * B3_rr)
533
f8 = equations.c_h * psi_avg
534
f9 = equations.c_h * 0.5f0 * (B3_ll + B3_rr)
535
# total energy flux is complicated and involves the previous components
536
f5 = (f1 *
537
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
538
+
539
0.5f0 * (+p_ll * v3_rr + p_rr * v3_ll
540
+ (v3_ll * B1_ll * B1_rr + v3_rr * B1_rr * B1_ll)
541
+ (v3_ll * B2_ll * B2_rr + v3_rr * B2_rr * B2_ll)
542
-
543
(v1_ll * B3_ll * B1_rr + v1_rr * B3_rr * B1_ll)
544
-
545
(v2_ll * B3_ll * B2_rr + v2_rr * B3_rr * B2_ll)
546
+
547
equations.c_h * (B3_ll * psi_rr + B3_rr * psi_ll)))
548
end
549
550
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
551
end
552
553
@inline function flux_hindenlang_gassner(u_ll, u_rr, normal_direction::AbstractVector,
554
equations::IdealGlmMhdEquations3D)
555
# Unpack left and right states
556
rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll, psi_ll = cons2prim(u_ll,
557
equations)
558
rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr, psi_rr = cons2prim(u_rr,
559
equations)
560
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2] +
561
v3_ll * normal_direction[3]
562
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2] +
563
v3_rr * normal_direction[3]
564
B_dot_n_ll = B1_ll * normal_direction[1] + B2_ll * normal_direction[2] +
565
B3_ll * normal_direction[3]
566
B_dot_n_rr = B1_rr * normal_direction[1] + B2_rr * normal_direction[2] +
567
B3_rr * normal_direction[3]
568
569
# Compute the necessary mean values needed for either direction
570
rho_mean = ln_mean(rho_ll, rho_rr)
571
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
572
# in exact arithmetic since
573
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
574
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
575
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
576
v1_avg = 0.5f0 * (v1_ll + v1_rr)
577
v2_avg = 0.5f0 * (v2_ll + v2_rr)
578
v3_avg = 0.5f0 * (v3_ll + v3_rr)
579
p_avg = 0.5f0 * (p_ll + p_rr)
580
psi_avg = 0.5f0 * (psi_ll + psi_rr)
581
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
582
magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)
583
584
# Calculate fluxes depending on normal_direction
585
f1 = rho_mean * 0.5f0 * (v_dot_n_ll + v_dot_n_rr)
586
f2 = (f1 * v1_avg + (p_avg + magnetic_square_avg) * normal_direction[1]
587
-
588
0.5f0 * (B_dot_n_ll * B1_rr + B_dot_n_rr * B1_ll))
589
f3 = (f1 * v2_avg + (p_avg + magnetic_square_avg) * normal_direction[2]
590
-
591
0.5f0 * (B_dot_n_ll * B2_rr + B_dot_n_rr * B2_ll))
592
f4 = (f1 * v3_avg + (p_avg + magnetic_square_avg) * normal_direction[3]
593
-
594
0.5f0 * (B_dot_n_ll * B3_rr + B_dot_n_rr * B3_ll))
595
#f5 below
596
f6 = (equations.c_h * psi_avg * normal_direction[1]
597
+
598
0.5f0 * (v_dot_n_ll * B1_ll - v1_ll * B_dot_n_ll +
599
v_dot_n_rr * B1_rr - v1_rr * B_dot_n_rr))
600
f7 = (equations.c_h * psi_avg * normal_direction[2]
601
+
602
0.5f0 * (v_dot_n_ll * B2_ll - v2_ll * B_dot_n_ll +
603
v_dot_n_rr * B2_rr - v2_rr * B_dot_n_rr))
604
f8 = (equations.c_h * psi_avg * normal_direction[3]
605
+
606
0.5f0 * (v_dot_n_ll * B3_ll - v3_ll * B_dot_n_ll +
607
v_dot_n_rr * B3_rr - v3_rr * B_dot_n_rr))
608
f9 = equations.c_h * 0.5f0 * (B_dot_n_ll + B_dot_n_rr)
609
# total energy flux is complicated and involves the previous components
610
f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
611
+
612
0.5f0 * (+p_ll * v_dot_n_rr + p_rr * v_dot_n_ll
613
+ (v_dot_n_ll * B1_ll * B1_rr + v_dot_n_rr * B1_rr * B1_ll)
614
+ (v_dot_n_ll * B2_ll * B2_rr + v_dot_n_rr * B2_rr * B2_ll)
615
+ (v_dot_n_ll * B3_ll * B3_rr + v_dot_n_rr * B3_rr * B3_ll)
616
-
617
(v1_ll * B_dot_n_ll * B1_rr + v1_rr * B_dot_n_rr * B1_ll)
618
-
619
(v2_ll * B_dot_n_ll * B2_rr + v2_rr * B_dot_n_rr * B2_ll)
620
-
621
(v3_ll * B_dot_n_ll * B3_rr + v3_rr * B_dot_n_rr * B3_ll)
622
+
623
equations.c_h * (B_dot_n_ll * psi_rr + B_dot_n_rr * psi_ll)))
624
625
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
626
end
627
628
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation
629
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
630
equations::IdealGlmMhdEquations3D)
631
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, _ = u_ll
632
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, _ = u_rr
633
634
# Calculate the left/right velocities and fast magnetoacoustic wave speeds
635
if orientation == 1
636
v_ll = rho_v1_ll / rho_ll
637
v_rr = rho_v1_rr / rho_rr
638
elseif orientation == 2
639
v_ll = rho_v2_ll / rho_ll
640
v_rr = rho_v2_rr / rho_rr
641
else # orientation == 3
642
v_ll = rho_v3_ll / rho_ll
643
v_rr = rho_v3_rr / rho_rr
644
end
645
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
646
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
647
648
return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)
649
end
650
651
@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
652
equations::IdealGlmMhdEquations3D)
653
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, _ = u_ll
654
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, _ = u_rr
655
656
# Calculate normal velocities and fast magnetoacoustic wave speeds
657
# left
658
v1_ll = rho_v1_ll / rho_ll
659
v2_ll = rho_v2_ll / rho_ll
660
v3_ll = rho_v3_ll / rho_ll
661
v_ll = (v1_ll * normal_direction[1]
662
+ v2_ll * normal_direction[2]
663
+ v3_ll * normal_direction[3])
664
cf_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
665
# right
666
v1_rr = rho_v1_rr / rho_rr
667
v2_rr = rho_v2_rr / rho_rr
668
v3_rr = rho_v3_rr / rho_rr
669
v_rr = (v1_rr * normal_direction[1]
670
+ v2_rr * normal_direction[2]
671
+ v3_rr * normal_direction[3])
672
cf_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
673
674
# wave speeds already scaled by norm(normal_direction) in [`calc_fast_wavespeed`](@ref)
675
return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)
676
end
677
678
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
679
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
680
equations::IdealGlmMhdEquations3D)
681
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, _ = u_ll
682
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, _ = u_rr
683
684
# Calculate the left/right velocities and fast magnetoacoustic wave speeds
685
if orientation == 1
686
v_ll = rho_v1_ll / rho_ll
687
v_rr = rho_v1_rr / rho_rr
688
elseif orientation == 2
689
v_ll = rho_v2_ll / rho_ll
690
v_rr = rho_v2_rr / rho_rr
691
else # orientation == 3
692
v_ll = rho_v3_ll / rho_ll
693
v_rr = rho_v3_rr / rho_rr
694
end
695
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
696
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
697
698
return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)
699
end
700
701
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
702
@inline function max_abs_speed(u_ll, u_rr, normal_direction::AbstractVector,
703
equations::IdealGlmMhdEquations3D)
704
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, _ = u_ll
705
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, _ = u_rr
706
707
# Calculate normal velocities and fast magnetoacoustic wave speeds
708
# left
709
v1_ll = rho_v1_ll / rho_ll
710
v2_ll = rho_v2_ll / rho_ll
711
v3_ll = rho_v3_ll / rho_ll
712
v_ll = (v1_ll * normal_direction[1]
713
+ v2_ll * normal_direction[2]
714
+ v3_ll * normal_direction[3])
715
cf_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
716
# right
717
v1_rr = rho_v1_rr / rho_rr
718
v2_rr = rho_v2_rr / rho_rr
719
v3_rr = rho_v3_rr / rho_rr
720
v_rr = (v1_rr * normal_direction[1]
721
+ v2_rr * normal_direction[2]
722
+ v3_rr * normal_direction[3])
723
cf_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
724
725
# wave speeds already scaled by norm(normal_direction) in [`calc_fast_wavespeed`](@ref)
726
return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)
727
end
728
729
# Calculate estimate for minimum and maximum wave speeds for HLL-type fluxes
730
@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,
731
equations::IdealGlmMhdEquations3D)
732
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, _ = u_ll
733
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, _ = u_rr
734
735
# Calculate primitive variables and speed of sound
736
v1_ll = rho_v1_ll / rho_ll
737
v2_ll = rho_v2_ll / rho_ll
738
v3_ll = rho_v3_ll / rho_ll
739
740
v1_rr = rho_v1_rr / rho_rr
741
v2_rr = rho_v2_rr / rho_rr
742
v3_rr = rho_v3_rr / rho_rr
743
744
# Approximate the left-most and right-most eigenvalues in the Riemann fan
745
if orientation == 1 # x-direction
746
λ_min = v1_ll - calc_fast_wavespeed(u_ll, orientation, equations)
747
λ_max = v1_rr + calc_fast_wavespeed(u_rr, orientation, equations)
748
elseif orientation == 2 # y-direction
749
λ_min = v2_ll - calc_fast_wavespeed(u_ll, orientation, equations)
750
λ_max = v2_rr + calc_fast_wavespeed(u_rr, orientation, equations)
751
else # z-direction
752
λ_min = v3_ll - calc_fast_wavespeed(u_ll, orientation, equations)
753
λ_max = v3_rr + calc_fast_wavespeed(u_rr, orientation, equations)
754
end
755
756
return λ_min, λ_max
757
end
758
759
@inline function min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
760
equations::IdealGlmMhdEquations3D)
761
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, _ = u_ll
762
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, _ = u_rr
763
764
# Calculate primitive velocity variables
765
v1_ll = rho_v1_ll / rho_ll
766
v2_ll = rho_v2_ll / rho_ll
767
v3_ll = rho_v3_ll / rho_ll
768
769
v1_rr = rho_v1_rr / rho_rr
770
v2_rr = rho_v2_rr / rho_rr
771
v3_rr = rho_v3_rr / rho_rr
772
773
v_normal_ll = (v1_ll * normal_direction[1] +
774
v2_ll * normal_direction[2] +
775
v3_ll * normal_direction[3])
776
v_normal_rr = (v1_rr * normal_direction[1] +
777
v2_rr * normal_direction[2] +
778
v3_rr * normal_direction[3])
779
780
# Estimate the min/max eigenvalues in the normal direction
781
λ_min = v_normal_ll - calc_fast_wavespeed(u_ll, normal_direction, equations)
782
λ_max = v_normal_rr + calc_fast_wavespeed(u_rr, normal_direction, equations)
783
784
return λ_min, λ_max
785
end
786
787
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
788
@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,
789
equations::IdealGlmMhdEquations3D)
790
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, _ = u_ll
791
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, _ = u_rr
792
793
# Calculate primitive variables and speed of sound
794
v1_ll = rho_v1_ll / rho_ll
795
v2_ll = rho_v2_ll / rho_ll
796
v3_ll = rho_v3_ll / rho_ll
797
798
v1_rr = rho_v1_rr / rho_rr
799
v2_rr = rho_v2_rr / rho_rr
800
v3_rr = rho_v3_rr / rho_rr
801
802
# Approximate the left-most and right-most eigenvalues in the Riemann fan
803
if orientation == 1 # x-direction
804
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
805
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
806
807
λ_min = min(v1_ll - c_f_ll, v1_rr - c_f_rr)
808
λ_max = max(v1_ll + c_f_ll, v1_rr + c_f_rr)
809
elseif orientation == 2 # y-direction
810
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
811
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
812
813
λ_min = min(v2_ll - c_f_ll, v2_rr - c_f_rr)
814
λ_max = max(v2_ll + c_f_ll, v2_rr + c_f_rr)
815
else # z-direction
816
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
817
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
818
819
λ_min = min(v3_ll - c_f_ll, v3_rr - c_f_rr)
820
λ_max = max(v3_ll + c_f_ll, v3_rr + c_f_rr)
821
end
822
823
return λ_min, λ_max
824
end
825
826
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
827
@inline function min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector,
828
equations::IdealGlmMhdEquations3D)
829
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, _ = u_ll
830
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, _ = u_rr
831
832
# Calculate primitive velocity variables
833
v1_ll = rho_v1_ll / rho_ll
834
v2_ll = rho_v2_ll / rho_ll
835
v3_ll = rho_v3_ll / rho_ll
836
837
v1_rr = rho_v1_rr / rho_rr
838
v2_rr = rho_v2_rr / rho_rr
839
v3_rr = rho_v3_rr / rho_rr
840
841
v_normal_ll = (v1_ll * normal_direction[1] +
842
v2_ll * normal_direction[2] +
843
v3_ll * normal_direction[3])
844
v_normal_rr = (v1_rr * normal_direction[1] +
845
v2_rr * normal_direction[2] +
846
v3_rr * normal_direction[3])
847
848
c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
849
c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
850
851
# Estimate the min/max eigenvalues in the normal direction
852
λ_min = min(v_normal_ll - c_f_ll, v_normal_rr - c_f_rr)
853
λ_max = max(v_normal_ll + c_f_ll, v_normal_rr + c_f_rr)
854
855
return λ_min, λ_max
856
end
857
858
"""
859
min_max_speed_einfeldt(u_ll, u_rr, orientation_or_normal_direction, equations::IdealGlmMhdEquations3D)
860
861
Calculate minimum and maximum wave speeds for HLL-type fluxes as in
862
- Li (2005)
863
An HLLC Riemann solver for magneto-hydrodynamics
864
[DOI: 10.1016/j.jcp.2004.08.020](https://doi.org/10.1016/j.jcp.2004.08.020)
865
"""
866
@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,
867
equations::IdealGlmMhdEquations3D)
868
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, _ = u_ll
869
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, _ = u_rr
870
871
# Calculate primitive variables and speed of sound
872
v1_ll = rho_v1_ll / rho_ll
873
v2_ll = rho_v2_ll / rho_ll
874
v3_ll = rho_v3_ll / rho_ll
875
876
v1_rr = rho_v1_rr / rho_rr
877
v2_rr = rho_v2_rr / rho_rr
878
v3_rr = rho_v3_rr / rho_rr
879
880
# Approximate the left-most and right-most eigenvalues in the Riemann fan
881
if orientation == 1 # x-direction
882
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
883
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
884
vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)
885
λ_min = min(v1_ll - c_f_ll, vel_roe - c_f_roe)
886
λ_max = max(v1_rr + c_f_rr, vel_roe + c_f_roe)
887
elseif orientation == 2 # y-direction
888
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
889
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
890
vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)
891
λ_min = min(v2_ll - c_f_ll, vel_roe - c_f_roe)
892
λ_max = max(v2_rr + c_f_rr, vel_roe + c_f_roe)
893
else # z-direction
894
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
895
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
896
vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)
897
λ_min = min(v3_ll - c_f_ll, vel_roe - c_f_roe)
898
λ_max = max(v3_rr + c_f_rr, vel_roe + c_f_roe)
899
end
900
901
return λ_min, λ_max
902
end
903
904
@inline function min_max_speed_einfeldt(u_ll, u_rr, normal_direction::AbstractVector,
905
equations::IdealGlmMhdEquations3D)
906
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, _ = u_ll
907
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, _ = u_rr
908
909
# Calculate primitive velocity variables
910
v1_ll = rho_v1_ll / rho_ll
911
v2_ll = rho_v2_ll / rho_ll
912
v3_ll = rho_v3_ll / rho_ll
913
914
v1_rr = rho_v1_rr / rho_rr
915
v2_rr = rho_v2_rr / rho_rr
916
v3_rr = rho_v3_rr / rho_rr
917
918
v_normal_ll = (v1_ll * normal_direction[1] +
919
v2_ll * normal_direction[2] +
920
v3_ll * normal_direction[3])
921
v_normal_rr = (v1_rr * normal_direction[1] +
922
v2_rr * normal_direction[2] +
923
v3_rr * normal_direction[3])
924
925
c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
926
c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
927
v_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, normal_direction, equations)
928
929
# Estimate the min/max eigenvalues in the normal direction
930
λ_min = min(v_normal_ll - c_f_ll, v_roe - c_f_roe)
931
λ_max = max(v_normal_rr + c_f_rr, v_roe + c_f_roe)
932
933
return λ_min, λ_max
934
end
935
936
# Rotate normal vector to x-axis; normal, tangent1 and tangent2 need to be orthonormal
937
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the directions
938
# has been normalized prior to this rotation of the state vector
939
# Note, for ideal GLM-MHD only the velocities and magnetic field variables rotate
940
@inline function rotate_to_x(u, normal_vector, tangent1, tangent2,
941
equations::IdealGlmMhdEquations3D)
942
# Multiply with [ 1 0 0 0 0 0 0 0 0;
943
# 0 ― normal_vector ― 0 0 0 0 0;
944
# 0 ― tangent1 ― 0 0 0 0 0;
945
# 0 ― tangent2 ― 0 0 0 0 0;
946
# 0 0 0 0 1 0 0 0 0;
947
# 0 0 0 0 0 ― normal_vector ― 0;
948
# 0 0 0 0 0 ― tangent1 ― 0;
949
# 0 0 0 0 0 ― tangent2 ― 0;
950
# 0 0 0 0 0 0 0 0 1 ]
951
return SVector(u[1],
952
normal_vector[1] * u[2] + normal_vector[2] * u[3] +
953
normal_vector[3] * u[4],
954
tangent1[1] * u[2] + tangent1[2] * u[3] + tangent1[3] * u[4],
955
tangent2[1] * u[2] + tangent2[2] * u[3] + tangent2[3] * u[4],
956
u[5],
957
normal_vector[1] * u[6] + normal_vector[2] * u[7] +
958
normal_vector[3] * u[8],
959
tangent1[1] * u[6] + tangent1[2] * u[7] + tangent1[3] * u[8],
960
tangent2[1] * u[6] + tangent2[2] * u[7] + tangent2[3] * u[8],
961
u[9])
962
end
963
964
# Rotate x-axis to normal vector; normal, tangent1 and tangent2 need to be orthonormal
965
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the directions
966
# has been normalized prior to this back-rotation of the state vector
967
# Note, for ideal GLM-MHD only the velocities and magnetic field variables back-rotate
968
@inline function rotate_from_x(u, normal_vector, tangent1, tangent2,
969
equations::IdealGlmMhdEquations3D)
970
# Multiply with [ 1 0 0 0 0 0 0 0 0;
971
# 0 | | | 0 0 0 0 0;
972
# 0 normal_vector tangent1 tangent2 0 0 0 0 0;
973
# 0 | | | 0 0 0 0 0;
974
# 0 0 0 0 1 0 0 0 0;
975
# 0 0 0 0 0 | | | 0;
976
# 0 0 0 0 0 normal_vector tangent1 tangent2 0;
977
# 0 0 0 0 0 | | | 0;
978
# 0 0 0 0 0 0 0 0 1 ]
979
return SVector(u[1],
980
normal_vector[1] * u[2] + tangent1[1] * u[3] + tangent2[1] * u[4],
981
normal_vector[2] * u[2] + tangent1[2] * u[3] + tangent2[2] * u[4],
982
normal_vector[3] * u[2] + tangent1[3] * u[3] + tangent2[3] * u[4],
983
u[5],
984
normal_vector[1] * u[6] + tangent1[1] * u[7] + tangent2[1] * u[8],
985
normal_vector[2] * u[6] + tangent1[2] * u[7] + tangent2[2] * u[8],
986
normal_vector[3] * u[6] + tangent1[3] * u[7] + tangent2[3] * u[8],
987
u[9])
988
end
989
990
@inline function max_abs_speeds(u, equations::IdealGlmMhdEquations3D)
991
rho, rho_v1, rho_v2, rho_v3, _ = u
992
v1 = rho_v1 / rho
993
v2 = rho_v2 / rho
994
v3 = rho_v3 / rho
995
cf_x_direction = calc_fast_wavespeed(u, 1, equations)
996
cf_y_direction = calc_fast_wavespeed(u, 2, equations)
997
cf_z_direction = calc_fast_wavespeed(u, 3, equations)
998
999
return abs(v1) + cf_x_direction, abs(v2) + cf_y_direction, abs(v3) + cf_z_direction
1000
end
1001
1002
# Convert conservative variables to primitive
1003
@inline function cons2prim(u, equations::IdealGlmMhdEquations3D)
1004
rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
1005
1006
v1 = rho_v1 / rho
1007
v2 = rho_v2 / rho
1008
v3 = rho_v3 / rho
1009
p = (equations.gamma - 1) * (rho_e -
1010
0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3
1011
+ B1 * B1 + B2 * B2 + B3 * B3
1012
+ psi * psi))
1013
1014
return SVector(rho, v1, v2, v3, p, B1, B2, B3, psi)
1015
end
1016
1017
# Convert conservative variables to entropy
1018
@inline function cons2entropy(u, equations::IdealGlmMhdEquations3D)
1019
rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
1020
v1 = rho_v1 / rho
1021
v2 = rho_v2 / rho
1022
v3 = rho_v3 / rho
1023
v_square = v1^2 + v2^2 + v3^2
1024
p = (equations.gamma - 1) *
1025
(rho_e - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2) - 0.5f0 * psi^2)
1026
s = log(p) - equations.gamma * log(rho)
1027
rho_p = rho / p
1028
1029
w1 = (equations.gamma - s) * equations.inv_gamma_minus_one -
1030
0.5f0 * rho_p * v_square
1031
w2 = rho_p * v1
1032
w3 = rho_p * v2
1033
w4 = rho_p * v3
1034
w5 = -rho_p
1035
w6 = rho_p * B1
1036
w7 = rho_p * B2
1037
w8 = rho_p * B3
1038
w9 = rho_p * psi
1039
1040
return SVector(w1, w2, w3, w4, w5, w6, w7, w8, w9)
1041
end
1042
1043
# Convert primitive to conservative variables
1044
@inline function prim2cons(prim, equations::IdealGlmMhdEquations3D)
1045
rho, v1, v2, v3, p, B1, B2, B3, psi = prim
1046
1047
rho_v1 = rho * v1
1048
rho_v2 = rho * v2
1049
rho_v3 = rho * v3
1050
rho_e = p * equations.inv_gamma_minus_one +
1051
0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3) +
1052
0.5f0 * (B1^2 + B2^2 + B3^2) + 0.5f0 * psi^2
1053
1054
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi)
1055
end
1056
1057
@inline function density(u, equations::IdealGlmMhdEquations3D)
1058
return u[1]
1059
end
1060
1061
@inline function velocity(u, equations::IdealGlmMhdEquations3D)
1062
rho = u[1]
1063
v1 = u[2] / rho
1064
v2 = u[3] / rho
1065
v3 = u[4] / rho
1066
return SVector(v1, v2, v3)
1067
end
1068
1069
@inline function velocity(u, orientation::Int, equations::IdealGlmMhdEquations3D)
1070
rho = u[1]
1071
v = u[orientation + 1] / rho
1072
return v
1073
end
1074
1075
@inline function pressure(u, equations::IdealGlmMhdEquations3D)
1076
rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
1077
p = (equations.gamma - 1) * (rho_e - 0.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2) / rho
1078
-
1079
0.5f0 * (B1^2 + B2^2 + B3^2)
1080
-
1081
0.5f0 * psi^2)
1082
return p
1083
end
1084
1085
@inline function density_pressure(u, equations::IdealGlmMhdEquations3D)
1086
rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
1087
p = (equations.gamma - 1) * (rho_e - 0.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2) / rho
1088
-
1089
0.5f0 * (B1^2 + B2^2 + B3^2)
1090
-
1091
0.5f0 * psi^2)
1092
return rho * p
1093
end
1094
1095
# Compute the fastest wave speed for ideal MHD equations: c_f, the fast magnetoacoustic eigenvalue
1096
@inline function calc_fast_wavespeed(cons, orientation::Integer,
1097
equations::IdealGlmMhdEquations3D)
1098
rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = cons
1099
v1 = rho_v1 / rho
1100
v2 = rho_v2 / rho
1101
v3 = rho_v3 / rho
1102
kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
1103
mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)
1104
p = (equations.gamma - 1) * (rho_e - kin_en - mag_en - 0.5f0 * psi^2)
1105
a_square = equations.gamma * p / rho
1106
sqrt_rho = sqrt(rho)
1107
b1 = B1 / sqrt_rho
1108
b2 = B2 / sqrt_rho
1109
b3 = B3 / sqrt_rho
1110
b_square = b1 * b1 + b2 * b2 + b3 * b3
1111
if orientation == 1 # x-direction
1112
c_f = sqrt(0.5f0 * (a_square + b_square) +
1113
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b1^2))
1114
elseif orientation == 2 # y-direction
1115
c_f = sqrt(0.5f0 * (a_square + b_square) +
1116
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b2^2))
1117
else # z-direction
1118
c_f = sqrt(0.5f0 * (a_square + b_square) +
1119
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b3^2))
1120
end
1121
return c_f
1122
end
1123
1124
@inline function calc_fast_wavespeed(cons, normal_direction::AbstractVector,
1125
equations::IdealGlmMhdEquations3D)
1126
rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = cons
1127
v1 = rho_v1 / rho
1128
v2 = rho_v2 / rho
1129
v3 = rho_v3 / rho
1130
kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
1131
mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)
1132
p = (equations.gamma - 1) * (rho_e - kin_en - mag_en - 0.5f0 * psi^2)
1133
a_square = equations.gamma * p / rho
1134
sqrt_rho = sqrt(rho)
1135
b1 = B1 / sqrt_rho
1136
b2 = B2 / sqrt_rho
1137
b3 = B3 / sqrt_rho
1138
b_square = b1 * b1 + b2 * b2 + b3 * b3
1139
norm_squared = (normal_direction[1] * normal_direction[1] +
1140
normal_direction[2] * normal_direction[2] +
1141
normal_direction[3] * normal_direction[3])
1142
b_dot_n_squared = (b1 * normal_direction[1] +
1143
b2 * normal_direction[2] +
1144
b3 * normal_direction[3])^2 / norm_squared
1145
1146
c_f = sqrt((0.5f0 * (a_square + b_square) +
1147
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b_dot_n_squared)) *
1148
norm_squared)
1149
return c_f
1150
end
1151
1152
"""
1153
calc_fast_wavespeed_roe(u_ll, u_rr, orientation_or_normal_direction, equations::IdealGlmMhdEquations3D)
1154
1155
Compute the fast magnetoacoustic wave speed using Roe averages as given by
1156
- Cargo and Gallice (1997)
1157
Roe Matrices for Ideal MHD and Systematic Construction
1158
of Roe Matrices for Systems of Conservation Laws
1159
[DOI: 10.1006/jcph.1997.5773](https://doi.org/10.1006/jcph.1997.5773)
1160
"""
1161
@inline function calc_fast_wavespeed_roe(u_ll, u_rr, orientation::Integer,
1162
equations::IdealGlmMhdEquations3D)
1163
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
1164
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
1165
1166
# Calculate primitive variables
1167
v1_ll = rho_v1_ll / rho_ll
1168
v2_ll = rho_v2_ll / rho_ll
1169
v3_ll = rho_v3_ll / rho_ll
1170
kin_en_ll = 0.5f0 * (rho_v1_ll * v1_ll + rho_v2_ll * v2_ll + rho_v3_ll * v3_ll)
1171
mag_norm_ll = B1_ll * B1_ll + B2_ll * B2_ll + B3_ll * B3_ll
1172
p_ll = (equations.gamma - 1) *
1173
(rho_e_ll - kin_en_ll - 0.5f0 * mag_norm_ll - 0.5f0 * psi_ll^2)
1174
1175
v1_rr = rho_v1_rr / rho_rr
1176
v2_rr = rho_v2_rr / rho_rr
1177
v3_rr = rho_v3_rr / rho_rr
1178
kin_en_rr = 0.5f0 * (rho_v1_rr * v1_rr + rho_v2_rr * v2_rr + rho_v3_rr * v3_rr)
1179
mag_norm_rr = B1_rr * B1_rr + B2_rr * B2_rr + B3_rr * B3_rr
1180
p_rr = (equations.gamma - 1) *
1181
(rho_e_rr - kin_en_rr - 0.5f0 * mag_norm_rr - 0.5f0 * psi_rr^2)
1182
1183
# compute total pressure which is thermal + magnetic pressures
1184
p_total_ll = p_ll + 0.5f0 * mag_norm_ll
1185
p_total_rr = p_rr + 0.5f0 * mag_norm_rr
1186
1187
# compute the Roe density averages
1188
sqrt_rho_ll = sqrt(rho_ll)
1189
sqrt_rho_rr = sqrt(rho_rr)
1190
inv_sqrt_rho_add = 1 / (sqrt_rho_ll + sqrt_rho_rr)
1191
inv_sqrt_rho_prod = 1 / (sqrt_rho_ll * sqrt_rho_rr)
1192
rho_ll_roe = sqrt_rho_ll * inv_sqrt_rho_add
1193
rho_rr_roe = sqrt_rho_rr * inv_sqrt_rho_add
1194
# Roe averages
1195
# velocities and magnetic fields
1196
v1_roe = v1_ll * rho_ll_roe + v1_rr * rho_rr_roe
1197
v2_roe = v2_ll * rho_ll_roe + v2_rr * rho_rr_roe
1198
v3_roe = v3_ll * rho_ll_roe + v3_rr * rho_rr_roe
1199
B1_roe = B1_ll * rho_ll_roe + B1_rr * rho_rr_roe
1200
B2_roe = B2_ll * rho_ll_roe + B2_rr * rho_rr_roe
1201
B3_roe = B3_ll * rho_ll_roe + B3_rr * rho_rr_roe
1202
# enthalpy
1203
H_ll = (rho_e_ll + p_total_ll) / rho_ll
1204
H_rr = (rho_e_rr + p_total_rr) / rho_rr
1205
H_roe = H_ll * rho_ll_roe + H_rr * rho_rr_roe
1206
# temporary variable see equation (4.12) in Cargo and Gallice
1207
X = 0.5f0 * ((B1_ll - B1_rr)^2 + (B2_ll - B2_rr)^2 + (B3_ll - B3_rr)^2) *
1208
inv_sqrt_rho_add^2
1209
# averaged components needed to compute c_f, the fast magnetoacoustic wave speed
1210
b_square_roe = (B1_roe^2 + B2_roe^2 + B3_roe^2) * inv_sqrt_rho_prod # scaled magnectic sum
1211
a_square_roe = ((2 - equations.gamma) * X +
1212
(equations.gamma - 1) *
1213
(H_roe - 0.5f0 * (v1_roe^2 + v2_roe^2 + v3_roe^2) -
1214
b_square_roe)) # acoustic speed
1215
# finally compute the average wave speed and set the output velocity (depends on orientation)
1216
if orientation == 1 # x-direction
1217
c_a_roe = B1_roe^2 * inv_sqrt_rho_prod # (squared) Alfvén wave speed
1218
a_star_roe = sqrt((a_square_roe + b_square_roe)^2 -
1219
4 * a_square_roe * c_a_roe)
1220
c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe))
1221
vel_out_roe = v1_roe
1222
elseif orientation == 2 # y-direction
1223
c_a_roe = B2_roe^2 * inv_sqrt_rho_prod # (squared) Alfvén wave speed
1224
a_star_roe = sqrt((a_square_roe + b_square_roe)^2 -
1225
4 * a_square_roe * c_a_roe)
1226
c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe))
1227
vel_out_roe = v2_roe
1228
else # z-direction
1229
c_a_roe = B3_roe^2 * inv_sqrt_rho_prod # (squared) Alfvén wave speed
1230
a_star_roe = sqrt((a_square_roe + b_square_roe)^2 -
1231
4 * a_square_roe * c_a_roe)
1232
c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe))
1233
vel_out_roe = v3_roe
1234
end
1235
1236
return vel_out_roe, c_f_roe
1237
end
1238
1239
@inline function calc_fast_wavespeed_roe(u_ll, u_rr, normal_direction::AbstractVector,
1240
equations::IdealGlmMhdEquations3D)
1241
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
1242
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
1243
1244
# Calculate primitive variables
1245
v1_ll = rho_v1_ll / rho_ll
1246
v2_ll = rho_v2_ll / rho_ll
1247
v3_ll = rho_v3_ll / rho_ll
1248
kin_en_ll = 0.5f0 * (rho_v1_ll * v1_ll + rho_v2_ll * v2_ll + rho_v3_ll * v3_ll)
1249
mag_norm_ll = B1_ll * B1_ll + B2_ll * B2_ll + B3_ll * B3_ll
1250
p_ll = (equations.gamma - 1) *
1251
(rho_e_ll - kin_en_ll - 0.5f0 * mag_norm_ll - 0.5f0 * psi_ll^2)
1252
1253
v1_rr = rho_v1_rr / rho_rr
1254
v2_rr = rho_v2_rr / rho_rr
1255
v3_rr = rho_v3_rr / rho_rr
1256
kin_en_rr = 0.5f0 * (rho_v1_rr * v1_rr + rho_v2_rr * v2_rr + rho_v3_rr * v3_rr)
1257
mag_norm_rr = B1_rr * B1_rr + B2_rr * B2_rr + B3_rr * B3_rr
1258
p_rr = (equations.gamma - 1) *
1259
(rho_e_rr - kin_en_rr - 0.5f0 * mag_norm_rr - 0.5f0 * psi_rr^2)
1260
1261
# compute total pressure which is thermal + magnetic pressures
1262
p_total_ll = p_ll + 0.5f0 * mag_norm_ll
1263
p_total_rr = p_rr + 0.5f0 * mag_norm_rr
1264
1265
# compute the Roe density averages
1266
sqrt_rho_ll = sqrt(rho_ll)
1267
sqrt_rho_rr = sqrt(rho_rr)
1268
inv_sqrt_rho_add = 1 / (sqrt_rho_ll + sqrt_rho_rr)
1269
inv_sqrt_rho_prod = 1 / (sqrt_rho_ll * sqrt_rho_rr)
1270
rho_ll_roe = sqrt_rho_ll * inv_sqrt_rho_add
1271
rho_rr_roe = sqrt_rho_rr * inv_sqrt_rho_add
1272
# Roe averages
1273
# velocities and magnetic fields
1274
v1_roe = v1_ll * rho_ll_roe + v1_rr * rho_rr_roe
1275
v2_roe = v2_ll * rho_ll_roe + v2_rr * rho_rr_roe
1276
v3_roe = v3_ll * rho_ll_roe + v3_rr * rho_rr_roe
1277
B1_roe = B1_ll * rho_ll_roe + B1_rr * rho_rr_roe
1278
B2_roe = B2_ll * rho_ll_roe + B2_rr * rho_rr_roe
1279
B3_roe = B3_ll * rho_ll_roe + B3_rr * rho_rr_roe
1280
# enthalpy
1281
H_ll = (rho_e_ll + p_total_ll) / rho_ll
1282
H_rr = (rho_e_rr + p_total_rr) / rho_rr
1283
H_roe = H_ll * rho_ll_roe + H_rr * rho_rr_roe
1284
# temporary variable see equation (4.12) in Cargo and Gallice
1285
X = 0.5f0 * ((B1_ll - B1_rr)^2 + (B2_ll - B2_rr)^2 + (B3_ll - B3_rr)^2) *
1286
inv_sqrt_rho_add^2
1287
# averaged components needed to compute c_f, the fast magnetoacoustic wave speed
1288
b_square_roe = (B1_roe^2 + B2_roe^2 + B3_roe^2) * inv_sqrt_rho_prod # scaled magnectic sum
1289
a_square_roe = ((2 - equations.gamma) * X +
1290
(equations.gamma - 1) *
1291
(H_roe - 0.5f0 * (v1_roe^2 + v2_roe^2 + v3_roe^2) -
1292
b_square_roe)) # acoustic speed
1293
1294
# finally compute the average wave speed and set the output velocity (depends on orientation)
1295
norm_squared = (normal_direction[1] * normal_direction[1] +
1296
normal_direction[2] * normal_direction[2] +
1297
normal_direction[3] * normal_direction[3])
1298
B_roe_dot_n_squared = (B1_roe * normal_direction[1] +
1299
B2_roe * normal_direction[2] +
1300
B3_roe * normal_direction[3])^2 / norm_squared
1301
1302
c_a_roe = B_roe_dot_n_squared * inv_sqrt_rho_prod # (squared) Alfvén wave speed
1303
a_star_roe = sqrt((a_square_roe + b_square_roe)^2 - 4 * a_square_roe * c_a_roe)
1304
c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe) * norm_squared)
1305
vel_out_roe = (v1_roe * normal_direction[1] +
1306
v2_roe * normal_direction[2] +
1307
v3_roe * normal_direction[3])
1308
1309
return vel_out_roe, c_f_roe
1310
end
1311
1312
# Calculate thermodynamic entropy for a conservative state `cons`
1313
@inline function entropy_thermodynamic(cons, equations::IdealGlmMhdEquations3D)
1314
# Pressure
1315
p = (equations.gamma - 1) *
1316
(cons[5] - 0.5f0 * (cons[2]^2 + cons[3]^2 + cons[4]^2) / cons[1]
1317
-
1318
0.5f0 * (cons[6]^2 + cons[7]^2 + cons[8]^2)
1319
-
1320
0.5f0 * cons[9]^2)
1321
1322
# Thermodynamic entropy
1323
s = log(p) - equations.gamma * log(cons[1])
1324
1325
return s
1326
end
1327
1328
# Calculate mathematical entropy for a conservative state `cons`
1329
@inline function entropy_math(cons, equations::IdealGlmMhdEquations3D)
1330
S = -entropy_thermodynamic(cons, equations) * cons[1] *
1331
equations.inv_gamma_minus_one
1332
1333
return S
1334
end
1335
1336
# Default entropy is the mathematical entropy
1337
@inline entropy(cons, equations::IdealGlmMhdEquations3D) = entropy_math(cons, equations)
1338
1339
# Calculate total energy for a conservative state `cons`
1340
@inline energy_total(cons, ::IdealGlmMhdEquations3D) = cons[5]
1341
1342
# Calculate kinetic energy for a conservative state `cons`
1343
@inline function energy_kinetic(cons, equations::IdealGlmMhdEquations3D)
1344
return 0.5f0 * (cons[2]^2 + cons[3]^2 + cons[4]^2) / cons[1]
1345
end
1346
1347
# Calculate the magnetic energy for a conservative state `cons'.
1348
# OBS! For non-dinmensional form of the ideal MHD magnetic pressure ≡ magnetic energy
1349
@inline function energy_magnetic(cons, ::IdealGlmMhdEquations3D)
1350
return 0.5f0 * (cons[6]^2 + cons[7]^2 + cons[8]^2)
1351
end
1352
1353
# Calculate internal energy for a conservative state `cons`
1354
@inline function energy_internal(cons, equations::IdealGlmMhdEquations3D)
1355
return (energy_total(cons, equations)
1356
-
1357
energy_kinetic(cons, equations)
1358
-
1359
energy_magnetic(cons, equations)
1360
-
1361
cons[9]^2 / 2)
1362
end
1363
1364
# Calculate the cross helicity (\vec{v}⋅\vec{B}) for a conservative state `cons'
1365
@inline function cross_helicity(cons, ::IdealGlmMhdEquations3D)
1366
return (cons[2] * cons[6] + cons[3] * cons[7] + cons[4] * cons[8]) / cons[1]
1367
end
1368
end # @muladd
1369
1370