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_2d.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
IdealGlmMhdEquations2D(gamma)
10
11
The ideal compressible GLM-MHD equations
12
```math
13
\frac{\partial}{\partial t}
14
\begin{pmatrix}
15
\rho \\ \rho \mathbf{v} \\ \rho e_{\text{total}} \\ \mathbf{B} \\ \psi
16
\end{pmatrix}
17
+
18
\nabla \cdot
19
\begin{pmatrix}
20
\rho \mathbf{v} \\ \rho (\mathbf{v} \otimes \mathbf{v}) + (p + \frac{1}{2} \Vert \mathbf{B} \Vert_2 ^2) \underline{I} - \mathbf{B} \otimes \mathbf{B} \\
21
\mathbf{v} (\frac{1}{2} \rho \Vert \mathbf{v} \Vert_2 ^2 + \frac{\gamma p}{\gamma - 1} + \Vert \mathbf{B} \Vert_2 ^2) - \mathbf{B} (\mathbf{v} \cdot \mathbf{B}) + c_h \psi \mathbf{B} \\
22
\mathbf{v} \otimes \mathbf{B} - \mathbf{B} \otimes \mathbf{v} + c_h \psi \underline{I} \\
23
c_h \mathbf{B}
24
\end{pmatrix}
25
+
26
(\nabla \cdot \mathbf{B})
27
\begin{pmatrix}
28
0 \\ \mathbf{B} \\ \mathbf{v} \cdot \mathbf{B} \\ \mathbf{v} \\ 0
29
\end{pmatrix}
30
+
31
(\nabla \psi) \cdot
32
\begin{pmatrix}
33
\mathbf{0} \\ 0 \\ \mathbf{v} \cdot \psi \\ 0 \\ \mathbf{v}
34
\end{pmatrix}
35
=
36
\begin{pmatrix}
37
0 \\ \mathbf{0} \\ 0 \\ \mathbf{0} \\ 0
38
\end{pmatrix}
39
```
40
for an ideal gas in two space dimensions.
41
Here, ``\mathbf{v}`` is the velocity, ``\mathbf{B}`` the magnetic field, ``c_h`` the hyperbolic divergence cleaning speed,
42
``\psi`` the generalized Lagrangian Multiplier (GLM),
43
``e_{\text{total}}`` the specific total energy, and
44
```math
45
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 - \frac{1}{2} \psi^2 \right)
46
```
47
the pressure, ``\gamma`` the total heat capacity ratio and ``\underline{I}`` the ``3\times 3`` identity matrix.
48
"""
49
struct IdealGlmMhdEquations2D{RealT <: Real} <:
50
AbstractIdealGlmMhdEquations{2, 9}
51
gamma::RealT # ratio of specific heats
52
inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications
53
c_h::RealT # GLM cleaning speed
54
55
function IdealGlmMhdEquations2D(gamma, c_h)
56
γ, inv_gamma_minus_one, c_h = promote(gamma, inv(gamma - 1), c_h)
57
return new{typeof(γ)}(γ, inv_gamma_minus_one, c_h)
58
end
59
end
60
61
function IdealGlmMhdEquations2D(gamma; initial_c_h = convert(typeof(gamma), NaN))
62
# Use `promote` to ensure that `gamma` and `initial_c_h` have the same type
63
return IdealGlmMhdEquations2D(promote(gamma, initial_c_h)...)
64
end
65
66
# Outer constructor for `@reset` works correctly
67
function IdealGlmMhdEquations2D(gamma, inv_gamma_minus_one, c_h)
68
return IdealGlmMhdEquations2D(gamma, c_h)
69
end
70
71
have_nonconservative_terms(::IdealGlmMhdEquations2D) = True()
72
73
function varnames(::typeof(cons2cons), ::IdealGlmMhdEquations2D)
74
return ("rho", "rho_v1", "rho_v2", "rho_v3", "rho_e_total", "B1", "B2", "B3", "psi")
75
end
76
function varnames(::typeof(cons2prim), ::IdealGlmMhdEquations2D)
77
return ("rho", "v1", "v2", "v3", "p", "B1", "B2", "B3", "psi")
78
end
79
function default_analysis_integrals(::IdealGlmMhdEquations2D)
80
return (entropy_timederivative, Val(:l2_divb), Val(:linf_divb))
81
end
82
83
# Helper function to extract the magnetic field vector from the conservative variables
84
magnetic_field(u, equations::IdealGlmMhdEquations2D) = SVector(u[6], u[7], u[8])
85
86
# Set initial conditions at physical location `x` for time `t`
87
"""
88
initial_condition_constant(x, t, equations::IdealGlmMhdEquations2D)
89
90
A constant initial condition to test free-stream preservation.
91
"""
92
function initial_condition_constant(x, t, equations::IdealGlmMhdEquations2D)
93
RealT = eltype(x)
94
rho = 1
95
rho_v1 = convert(RealT, 0.1)
96
rho_v2 = -convert(RealT, 0.2)
97
rho_v3 = -0.5f0
98
rho_e_total = 50
99
B1 = 3
100
B2 = -convert(RealT, 1.2)
101
B3 = 0.5f0
102
psi = 0
103
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi)
104
end
105
106
"""
107
initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations2D)
108
109
An Alfvén wave as smooth initial condition used for convergence tests.
110
See for reference section 4.2 in
111
- Dominik Derigs, Andrew R. Winters, Gregor J. Gassner, and Stefanie Walch (2016)
112
A novel high-order, entropy stable, 3D AMR MHD solver with guaranteed positive pressure
113
[DOI: 10.1016/j.jcp.2016.04.048](https://doi.org/10.1016/j.jcp.2016.04.048)
114
"""
115
function initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations2D)
116
# smooth Alfvén wave test from Derigs et al. (2016)
117
# domain must be set to [0, 1/cos(α)] x [0, 1/sin(α)], γ = 5/3
118
RealT = eltype(x)
119
alpha = 0.25f0 * convert(RealT, pi)
120
x_perp = x[1] * cos(alpha) + x[2] * sin(alpha)
121
B_perp = convert(RealT, 0.1) * sinpi(2 * (x_perp + t))
122
rho = 1
123
v1 = -B_perp * sin(alpha)
124
v2 = B_perp * cos(alpha)
125
v3 = convert(RealT, 0.1) * cospi(2 * (x_perp + t))
126
p = convert(RealT, 0.1)
127
B1 = cos(alpha) + v1
128
B2 = sin(alpha) + v2
129
B3 = v3
130
psi = 0
131
return prim2cons(SVector(rho, v1, v2, v3, p, B1, B2, B3, psi), equations)
132
end
133
134
"""
135
initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations2D)
136
137
A weak blast wave adapted from
138
- Sebastian Hennemann, Gregor J. Gassner (2020)
139
A provably entropy stable subcell shock capturing approach for high order split form DG
140
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
141
"""
142
function initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations2D)
143
# Adapted MHD version of the weak blast wave from Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
144
# Same discontinuity in the velocities but with magnetic fields
145
# Set up polar coordinates
146
RealT = eltype(x)
147
inicenter = (0, 0)
148
x_norm = x[1] - inicenter[1]
149
y_norm = x[2] - inicenter[2]
150
r = sqrt(x_norm^2 + y_norm^2)
151
phi = atan(y_norm, x_norm)
152
153
# Calculate primitive variables
154
rho = r > 0.5f0 ? one(RealT) : convert(RealT, 1.1691)
155
v1 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * cos(phi)
156
v2 = r > 0.5f0 ? zero(RealT) : convert(RealT, 0.1882) * sin(phi)
157
p = r > 0.5f0 ? one(RealT) : convert(RealT, 1.245)
158
159
return prim2cons(SVector(rho, v1, v2, 0, p, 1, 1, 1, 0), equations)
160
end
161
162
# Pre-defined source terms should be implemented as
163
# function source_terms_WHATEVER(u, x, t, equations::IdealGlmMhdEquations2D)
164
165
# Calculate 1D flux for a single point
166
@inline function flux(u, orientation::Integer, equations::IdealGlmMhdEquations2D)
167
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u
168
v1 = rho_v1 / rho
169
v2 = rho_v2 / rho
170
v3 = rho_v3 / rho
171
kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
172
mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)
173
p_over_gamma_minus_one = (rho_e_total - kin_en - mag_en - 0.5f0 * psi^2)
174
p = (equations.gamma - 1) * p_over_gamma_minus_one
175
if orientation == 1
176
f1 = rho_v1
177
f2 = rho_v1 * v1 + p + mag_en - B1^2
178
f3 = rho_v1 * v2 - B1 * B2
179
f4 = rho_v1 * v3 - B1 * B3
180
f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v1 -
181
B1 * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B1
182
f6 = equations.c_h * psi
183
f7 = v1 * B2 - v2 * B1
184
f8 = v1 * B3 - v3 * B1
185
f9 = equations.c_h * B1
186
else #if orientation == 2
187
f1 = rho_v2
188
f2 = rho_v2 * v1 - B2 * B1
189
f3 = rho_v2 * v2 + p + mag_en - B2^2
190
f4 = rho_v2 * v3 - B2 * B3
191
f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v2 -
192
B2 * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B2
193
f6 = v2 * B1 - v1 * B2
194
f7 = equations.c_h * psi
195
f8 = v2 * B3 - v3 * B2
196
f9 = equations.c_h * B2
197
end
198
199
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
200
end
201
202
# Calculate 1D flux for a single point in the normal direction
203
# Note, this directional vector is not normalized
204
@inline function flux(u, normal_direction::AbstractVector,
205
equations::IdealGlmMhdEquations2D)
206
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u
207
v1 = rho_v1 / rho
208
v2 = rho_v2 / rho
209
v3 = rho_v3 / rho
210
kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
211
mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)
212
p_over_gamma_minus_one = (rho_e_total - kin_en - mag_en - 0.5f0 * psi^2)
213
p = (equations.gamma - 1) * p_over_gamma_minus_one
214
215
v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]
216
B_normal = B1 * normal_direction[1] + B2 * normal_direction[2]
217
rho_v_normal = rho * v_normal
218
219
f1 = rho_v_normal
220
f2 = rho_v_normal * v1 - B1 * B_normal + (p + mag_en) * normal_direction[1]
221
f3 = rho_v_normal * v2 - B2 * B_normal + (p + mag_en) * normal_direction[2]
222
f4 = rho_v_normal * v3 - B3 * B_normal
223
f5 = ((kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v_normal
224
-
225
B_normal * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B_normal)
226
f6 = equations.c_h * psi * normal_direction[1] +
227
(v2 * B1 - v1 * B2) * normal_direction[2]
228
f7 = equations.c_h * psi * normal_direction[2] +
229
(v1 * B2 - v2 * B1) * normal_direction[1]
230
f8 = v_normal * B3 - v3 * B_normal
231
f9 = equations.c_h * B_normal
232
233
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
234
end
235
236
"""
237
flux_nonconservative_powell(u_ll, u_rr, orientation::Integer,
238
equations::IdealGlmMhdEquations2D)
239
flux_nonconservative_powell(u_ll, u_rr,
240
normal_direction::AbstractVector,
241
equations::IdealGlmMhdEquations2D)
242
243
Non-symmetric two-point flux discretizing the nonconservative (source) term of
244
Powell and the Galilean nonconservative term associated with the GLM multiplier
245
of the [`IdealGlmMhdEquations2D`](@ref).
246
247
On curvilinear meshes, the implementation differs from the reference since we use the averaged
248
normal direction for the metrics dealiasing.
249
250
## References
251
- Marvin Bohm, Andrew R.Winters, Gregor J. Gassner, Dominik Derigs,
252
Florian Hindenlang, Joachim Saur
253
An entropy stable nodal discontinuous Galerkin method for the resistive MHD
254
equations. Part I: Theory and numerical verification
255
[DOI: 10.1016/j.jcp.2018.06.027](https://doi.org/10.1016/j.jcp.2018.06.027)
256
"""
257
@inline function flux_nonconservative_powell(u_ll, u_rr, orientation::Integer,
258
equations::IdealGlmMhdEquations2D)
259
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
260
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
261
262
v1_ll = rho_v1_ll / rho_ll
263
v2_ll = rho_v2_ll / rho_ll
264
v3_ll = rho_v3_ll / rho_ll
265
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
266
267
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
268
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
269
if orientation == 1
270
f = SVector(0,
271
B1_ll * B1_rr,
272
B2_ll * B1_rr,
273
B3_ll * B1_rr,
274
v_dot_B_ll * B1_rr + v1_ll * psi_ll * psi_rr,
275
v1_ll * B1_rr,
276
v2_ll * B1_rr,
277
v3_ll * B1_rr,
278
v1_ll * psi_rr)
279
else # orientation == 2
280
f = SVector(0,
281
B1_ll * B2_rr,
282
B2_ll * B2_rr,
283
B3_ll * B2_rr,
284
v_dot_B_ll * B2_rr + v2_ll * psi_ll * psi_rr,
285
v1_ll * B2_rr,
286
v2_ll * B2_rr,
287
v3_ll * B2_rr,
288
v2_ll * psi_rr)
289
end
290
291
return f
292
end
293
294
@inline function flux_nonconservative_powell(u_ll, u_rr,
295
normal_direction::AbstractVector,
296
equations::IdealGlmMhdEquations2D)
297
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
298
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
299
300
v1_ll = rho_v1_ll / rho_ll
301
v2_ll = rho_v2_ll / rho_ll
302
v3_ll = rho_v3_ll / rho_ll
303
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
304
305
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
306
B_dot_n_rr = B1_rr * normal_direction[1] +
307
B2_rr * normal_direction[2]
308
309
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
310
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
311
f = SVector(0,
312
B1_ll * B_dot_n_rr,
313
B2_ll * B_dot_n_rr,
314
B3_ll * B_dot_n_rr,
315
v_dot_B_ll * B_dot_n_rr + v_dot_n_ll * psi_ll * psi_rr,
316
v1_ll * B_dot_n_rr,
317
v2_ll * B_dot_n_rr,
318
v3_ll * B_dot_n_rr,
319
v_dot_n_ll * psi_rr)
320
321
return f
322
end
323
324
# For `VolumeIntegralSubcellLimiting` the nonconservative flux is created as a callable struct to
325
# enable dispatch on the type of the nonconservative term (symmetric / jump).
326
struct FluxNonConservativePowellLocalSymmetric <:
327
FluxNonConservative{NonConservativeSymmetric()}
328
end
329
330
n_nonconservative_terms(::FluxNonConservativePowellLocalSymmetric) = 2
331
332
const flux_nonconservative_powell_local_symmetric = FluxNonConservativePowellLocalSymmetric()
333
334
"""
335
flux_nonconservative_powell_local_symmetric(u_ll, u_rr,
336
orientation::Integer,
337
equations::IdealGlmMhdEquations2D)
338
flux_nonconservative_powell_local_symmetric(u_ll, u_rr,
339
normal_direction::AbstractVector,
340
equations::IdealGlmMhdEquations2D)
341
342
Non-symmetric two-point flux discretizing the nonconservative (source) term of
343
Powell and the Galilean nonconservative term associated with the GLM multiplier
344
of the [`IdealGlmMhdEquations2D`](@ref).
345
346
This implementation uses a non-conservative term that can be written as the product
347
of local and symmetric parts. It is equivalent to the non-conservative flux of Bohm
348
et al. [`flux_nonconservative_powell`](@ref) for conforming meshes but it yields different
349
results on non-conforming meshes(!). On curvilinear meshes this formulation applies the
350
local normal direction compared to the averaged one used in [`flux_nonconservative_powell`](@ref).
351
352
The two other flux functions with the same name return either the local
353
or symmetric portion of the non-conservative flux based on the type of the
354
nonconservative_type argument, employing multiple dispatch. They are used to
355
compute the subcell fluxes in dg_2d_subcell_limiters.jl.
356
357
## References
358
- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts
359
Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.
360
"""
361
@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll, u_rr,
362
orientation::Integer,
363
equations::IdealGlmMhdEquations2D)
364
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
365
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
366
367
v1_ll = rho_v1_ll / rho_ll
368
v2_ll = rho_v2_ll / rho_ll
369
v3_ll = rho_v3_ll / rho_ll
370
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
371
372
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
373
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
374
psi_avg = (psi_ll + psi_rr) #* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
375
if orientation == 1
376
B1_avg = (B1_ll + B1_rr) #* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
377
f = SVector(0,
378
B1_ll * B1_avg,
379
B2_ll * B1_avg,
380
B3_ll * B1_avg,
381
v_dot_B_ll * B1_avg + v1_ll * psi_ll * psi_avg,
382
v1_ll * B1_avg,
383
v2_ll * B1_avg,
384
v3_ll * B1_avg,
385
v1_ll * psi_avg)
386
else # orientation == 2
387
B2_avg = (B2_ll + B2_rr) #* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
388
f = SVector(0,
389
B1_ll * B2_avg,
390
B2_ll * B2_avg,
391
B3_ll * B2_avg,
392
v_dot_B_ll * B2_avg + v2_ll * psi_ll * psi_avg,
393
v1_ll * B2_avg,
394
v2_ll * B2_avg,
395
v3_ll * B2_avg,
396
v2_ll * psi_avg)
397
end
398
399
return f
400
end
401
402
@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll, u_rr,
403
normal_direction::AbstractVector,
404
equations::IdealGlmMhdEquations2D)
405
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
406
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
407
408
v1_ll = rho_v1_ll / rho_ll
409
v2_ll = rho_v2_ll / rho_ll
410
v3_ll = rho_v3_ll / rho_ll
411
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
412
413
# The factor 0.5 of the averages can be omitted since it is already applied when this
414
# function is called.
415
psi_avg = (psi_ll + psi_rr)
416
B1_avg = (B1_ll + B1_rr)
417
B2_avg = (B2_ll + B2_rr)
418
419
B_dot_n_avg = B1_avg * normal_direction[1] + B2_avg * normal_direction[2]
420
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
421
422
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
423
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
424
f = SVector(0,
425
B1_ll * B_dot_n_avg,
426
B2_ll * B_dot_n_avg,
427
B3_ll * B_dot_n_avg,
428
v_dot_B_ll * B_dot_n_avg + v_dot_n_ll * psi_ll * psi_avg,
429
v1_ll * B_dot_n_avg,
430
v2_ll * B_dot_n_avg,
431
v3_ll * B_dot_n_avg,
432
v_dot_n_ll * psi_avg)
433
434
return f
435
end
436
437
"""
438
flux_nonconservative_powell_local_symmetric(u_ll, orientation::Integer,
439
equations::IdealGlmMhdEquations2D,
440
nonconservative_type::NonConservativeLocal,
441
nonconservative_term::Integer)
442
flux_nonconservative_powell_local_symmetric(u_ll, normal_direction_ll::AbstractVector,
443
equations::IdealGlmMhdEquations2D,
444
nonconservative_type::NonConservativeLocal,
445
nonconservative_term::Integer)
446
447
Local part of the Powell and GLM non-conservative terms. Needed for the calculation of
448
the non-conservative staggered "fluxes" for subcell limiting. See, e.g.,
449
- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts
450
Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.
451
This function is used to compute the subcell fluxes in dg_2d_subcell_limiters.jl.
452
"""
453
@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll,
454
orientation::Integer,
455
equations::IdealGlmMhdEquations2D,
456
nonconservative_type::NonConservativeLocal,
457
nonconservative_term::Integer)
458
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
459
460
if nonconservative_term == 1
461
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
462
v1_ll = rho_v1_ll / rho_ll
463
v2_ll = rho_v2_ll / rho_ll
464
v3_ll = rho_v3_ll / rho_ll
465
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
466
f = SVector(0,
467
B1_ll,
468
B2_ll,
469
B3_ll,
470
v_dot_B_ll,
471
v1_ll,
472
v2_ll,
473
v3_ll,
474
0)
475
else #nonconservative_term ==2
476
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
477
if orientation == 1
478
v1_ll = rho_v1_ll / rho_ll
479
f = SVector(0,
480
0,
481
0,
482
0,
483
v1_ll * psi_ll,
484
0,
485
0,
486
0,
487
v1_ll)
488
else #orientation == 2
489
v2_ll = rho_v2_ll / rho_ll
490
f = SVector(0,
491
0,
492
0,
493
0,
494
v2_ll * psi_ll,
495
0,
496
0,
497
0,
498
v2_ll)
499
end
500
end
501
return f
502
end
503
504
@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll,
505
normal_direction_ll::AbstractVector,
506
equations::IdealGlmMhdEquations2D,
507
nonconservative_type::NonConservativeLocal,
508
nonconservative_term::Integer)
509
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
510
511
if nonconservative_term == 1
512
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
513
v1_ll = rho_v1_ll / rho_ll
514
v2_ll = rho_v2_ll / rho_ll
515
v3_ll = rho_v3_ll / rho_ll
516
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
517
518
f = SVector(0,
519
B1_ll,
520
B2_ll,
521
B3_ll,
522
v_dot_B_ll,
523
v1_ll,
524
v2_ll,
525
v3_ll,
526
0)
527
else # nonconservative_term == 2
528
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
529
v1_ll = rho_v1_ll / rho_ll
530
v2_ll = rho_v2_ll / rho_ll
531
v_dot_n_ll = v1_ll * normal_direction_ll[1] + v2_ll * normal_direction_ll[2]
532
533
f = SVector(0,
534
0,
535
0,
536
0,
537
v_dot_n_ll * psi_ll,
538
0,
539
0,
540
0,
541
v_dot_n_ll)
542
end
543
return f
544
end
545
546
"""
547
flux_nonconservative_powell_local_symmetric(u_ll, orientation::Integer,
548
equations::IdealGlmMhdEquations2D,
549
nonconservative_type::NonConservativeSymmetric,
550
nonconservative_term::Integer)
551
flux_nonconservative_powell_local_symmetric(u_ll, normal_direction_avg::AbstractVector,
552
equations::IdealGlmMhdEquations2D,
553
nonconservative_type::NonConservativeSymmetric,
554
nonconservative_term::Integer)
555
556
Symmetric part of the Powell and GLM non-conservative terms. Needed for the calculation of
557
the non-conservative staggered "fluxes" for subcell limiting. See, e.g.,
558
- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts
559
Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.
560
This function is used to compute the subcell fluxes in dg_2d_subcell_limiters.jl.
561
"""
562
@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll, u_rr,
563
orientation::Integer,
564
equations::IdealGlmMhdEquations2D,
565
nonconservative_type::NonConservativeSymmetric,
566
nonconservative_term::Integer)
567
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
568
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
569
570
if nonconservative_term == 1
571
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
572
if orientation == 1
573
B1_avg = (B1_ll + B1_rr)#* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
574
f = SVector(0,
575
B1_avg,
576
B1_avg,
577
B1_avg,
578
B1_avg,
579
B1_avg,
580
B1_avg,
581
B1_avg,
582
0)
583
else # orientation == 2
584
B2_avg = (B2_ll + B2_rr)#* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
585
f = SVector(0,
586
B2_avg,
587
B2_avg,
588
B2_avg,
589
B2_avg,
590
B2_avg,
591
B2_avg,
592
B2_avg,
593
0)
594
end
595
else #nonconservative_term == 2
596
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
597
psi_avg = (psi_ll + psi_rr)#* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
598
f = SVector(0,
599
0,
600
0,
601
0,
602
psi_avg,
603
0,
604
0,
605
0,
606
psi_avg)
607
end
608
609
return f
610
end
611
612
@inline function (noncons_flux::FluxNonConservativePowellLocalSymmetric)(u_ll, u_rr,
613
normal_direction_avg::AbstractVector,
614
equations::IdealGlmMhdEquations2D,
615
nonconservative_type::NonConservativeSymmetric,
616
nonconservative_term::Integer)
617
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
618
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
619
620
if nonconservative_term == 1
621
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
622
# The factor 0.5 of the average can be omitted since it is already applied when this
623
# function is called.
624
B_dot_n_avg = (B1_ll + B1_rr) * normal_direction_avg[1] +
625
(B2_ll + B2_rr) * normal_direction_avg[2]
626
f = SVector(0,
627
B_dot_n_avg,
628
B_dot_n_avg,
629
B_dot_n_avg,
630
B_dot_n_avg,
631
B_dot_n_avg,
632
B_dot_n_avg,
633
B_dot_n_avg,
634
0)
635
else # nonconservative_term == 2
636
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
637
# The factor 0.5 of the average can be omitted since it is already applied when this
638
# function is called.
639
psi_avg = (psi_ll + psi_rr)
640
f = SVector(0,
641
0,
642
0,
643
0,
644
psi_avg,
645
0,
646
0,
647
0,
648
psi_avg)
649
end
650
651
return f
652
end
653
654
# For `VolumeIntegralSubcellLimiting` the nonconservative flux is created as a callable struct to
655
# enable dispatch on the type of the nonconservative term (symmetric / jump).
656
struct FluxNonConservativePowellLocalJump <:
657
FluxNonConservative{NonConservativeJump()}
658
end
659
660
n_nonconservative_terms(::FluxNonConservativePowellLocalJump) = 2
661
662
const flux_nonconservative_powell_local_jump = FluxNonConservativePowellLocalJump()
663
664
"""
665
flux_nonconservative_powell_local_jump(u_ll, u_rr,
666
orientation::Integer,
667
equations::IdealGlmMhdEquations2D)
668
flux_nonconservative_powell_local_jump(u_ll, u_rr,
669
normal_direction::AbstractVector,
670
equations::IdealGlmMhdEquations2D)
671
672
Non-symmetric two-point flux discretizing the nonconservative (source) term of
673
Powell and the Galilean nonconservative term associated with the GLM multiplier
674
of the [`IdealGlmMhdEquations2D`](@ref).
675
676
This implementation uses a non-conservative term that can be written as the product
677
of local and jump parts.
678
679
The two other flux functions with the same name return either the local
680
or jump portion of the non-conservative flux based on the type of the
681
nonconservative_type argument, employing multiple dispatch. They are used to
682
compute the subcell fluxes in dg_2d_subcell_limiters.jl.
683
684
## References
685
- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts
686
Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.
687
"""
688
@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll, u_rr,
689
orientation::Integer,
690
equations::IdealGlmMhdEquations2D)
691
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
692
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
693
694
v1_ll = rho_v1_ll / rho_ll
695
v2_ll = rho_v2_ll / rho_ll
696
v3_ll = rho_v3_ll / rho_ll
697
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
698
699
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
700
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
701
psi_jump = psi_rr - psi_ll
702
if orientation == 1
703
B1_jump = B1_rr - B1_ll # The flux is already multiplied by 0.5 wherever it is used in the code
704
f = SVector(0,
705
B1_ll * B1_jump,
706
B2_ll * B1_jump,
707
B3_ll * B1_jump,
708
v_dot_B_ll * B1_jump + v1_ll * psi_ll * psi_jump,
709
v1_ll * B1_jump,
710
v2_ll * B1_jump,
711
v3_ll * B1_jump,
712
v1_ll * psi_jump)
713
else # orientation == 2
714
B2_jump = B2_rr - B2_ll # The flux is already multiplied by 0.5 wherever it is used in the code
715
f = SVector(0,
716
B1_ll * B2_jump,
717
B2_ll * B2_jump,
718
B3_ll * B2_jump,
719
v_dot_B_ll * B2_jump + v2_ll * psi_ll * psi_jump,
720
v1_ll * B2_jump,
721
v2_ll * B2_jump,
722
v3_ll * B2_jump,
723
v2_ll * psi_jump)
724
end
725
726
return f
727
end
728
729
@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll, u_rr,
730
normal_direction::AbstractVector,
731
equations::IdealGlmMhdEquations2D)
732
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
733
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
734
735
v1_ll = rho_v1_ll / rho_ll
736
v2_ll = rho_v2_ll / rho_ll
737
v3_ll = rho_v3_ll / rho_ll
738
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
739
740
psi_jump = psi_rr - psi_ll
741
B1_jump = B1_rr - B1_ll
742
B2_jump = B2_rr - B2_ll
743
744
B_dot_n_jump = B1_jump * normal_direction[1] + B2_jump * normal_direction[2]
745
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
746
747
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
748
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
749
f = SVector(0,
750
B1_ll * B_dot_n_jump,
751
B2_ll * B_dot_n_jump,
752
B3_ll * B_dot_n_jump,
753
v_dot_B_ll * B_dot_n_jump + v_dot_n_ll * psi_ll * psi_jump,
754
v1_ll * B_dot_n_jump,
755
v2_ll * B_dot_n_jump,
756
v3_ll * B_dot_n_jump,
757
v_dot_n_ll * psi_jump)
758
759
return f
760
end
761
762
"""
763
flux_nonconservative_powell_local_jump(u_ll, orientation::Integer,
764
equations::IdealGlmMhdEquations2D,
765
nonconservative_type::NonConservativeLocal,
766
nonconservative_term::Integer)
767
flux_nonconservative_powell_local_jump(u_ll, normal_direction_ll::AbstractVector,
768
equations::IdealGlmMhdEquations2D,
769
nonconservative_type::NonConservativeLocal,
770
nonconservative_term::Integer)
771
772
Local part of the Powell and GLM non-conservative terms. Needed for the calculation of
773
the non-conservative staggered "fluxes" for subcell limiting.
774
This function is used to compute the subcell fluxes in dg_2d_subcell_limiters.jl.
775
776
## References
777
- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts
778
Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.
779
"""
780
@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll,
781
orientation::Integer,
782
equations::IdealGlmMhdEquations2D,
783
nonconservative_type::NonConservativeLocal,
784
nonconservative_term::Integer)
785
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
786
if nonconservative_term == 1
787
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
788
v1_ll = rho_v1_ll / rho_ll
789
v2_ll = rho_v2_ll / rho_ll
790
v3_ll = rho_v3_ll / rho_ll
791
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
792
f = SVector(0,
793
B1_ll,
794
B2_ll,
795
B3_ll,
796
v_dot_B_ll,
797
v1_ll,
798
v2_ll,
799
v3_ll,
800
0)
801
else #nonconservative_term ==2
802
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
803
if orientation == 1
804
v1_ll = rho_v1_ll / rho_ll
805
f = SVector(0,
806
0,
807
0,
808
0,
809
v1_ll * psi_ll,
810
0,
811
0,
812
0,
813
v1_ll)
814
else #orientation == 2
815
v2_ll = rho_v2_ll / rho_ll
816
f = SVector(0,
817
0,
818
0,
819
0,
820
v2_ll * psi_ll,
821
0,
822
0,
823
0,
824
v2_ll)
825
end
826
end
827
return f
828
end
829
830
@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll,
831
normal_direction_ll::AbstractVector,
832
equations::IdealGlmMhdEquations2D,
833
nonconservative_type::NonConservativeLocal,
834
nonconservative_term::Integer)
835
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
836
837
if nonconservative_term == 1
838
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
839
v1_ll = rho_v1_ll / rho_ll
840
v2_ll = rho_v2_ll / rho_ll
841
v3_ll = rho_v3_ll / rho_ll
842
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
843
844
f = SVector(0,
845
B1_ll,
846
B2_ll,
847
B3_ll,
848
v_dot_B_ll,
849
v1_ll,
850
v2_ll,
851
v3_ll,
852
0)
853
else # nonconservative_term == 2
854
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
855
v1_ll = rho_v1_ll / rho_ll
856
v2_ll = rho_v2_ll / rho_ll
857
v_dot_n_ll = v1_ll * normal_direction_ll[1] + v2_ll * normal_direction_ll[2]
858
859
f = SVector(0,
860
0,
861
0,
862
0,
863
v_dot_n_ll * psi_ll,
864
0,
865
0,
866
0,
867
v_dot_n_ll)
868
end
869
return f
870
end
871
872
"""
873
flux_nonconservative_powell_local_jump(u_ll, orientation::Integer,
874
equations::IdealGlmMhdEquations2D,
875
nonconservative_type::NonConservativeJump,
876
nonconservative_term::Integer)
877
flux_nonconservative_powell_local_jump(u_ll, normal_direction_avg::AbstractVector,
878
equations::IdealGlmMhdEquations2D,
879
nonconservative_type::NonConservativeJump,
880
nonconservative_term::Integer)
881
882
Jump part of the Powell and GLM non-conservative terms. Needed for the calculation of
883
the non-conservative staggered "fluxes" for subcell limiting.
884
This function is used to compute the subcell fluxes in dg_2d_subcell_limiters.jl.
885
886
## References
887
- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts
888
Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.
889
"""
890
@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll, u_rr,
891
orientation::Integer,
892
equations::IdealGlmMhdEquations2D,
893
nonconservative_type::NonConservativeJump,
894
nonconservative_term::Integer)
895
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
896
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
897
898
if nonconservative_term == 1
899
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
900
if orientation == 1
901
B1_jump = B1_rr - B1_ll # The flux is already multiplied by 0.5 wherever it is used in the code
902
f = SVector(0,
903
B1_jump,
904
B1_jump,
905
B1_jump,
906
B1_jump,
907
B1_jump,
908
B1_jump,
909
B1_jump,
910
0)
911
else # orientation == 2
912
B2_jump = B2_rr - B2_ll # The flux is already multiplied by 0.5 wherever it is used in the code
913
f = SVector(0,
914
B2_jump,
915
B2_jump,
916
B2_jump,
917
B2_jump,
918
B2_jump,
919
B2_jump,
920
B2_jump,
921
0)
922
end
923
else #nonconservative_term == 2
924
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
925
psi_jump = psi_rr - psi_ll # The flux is already multiplied by 0.5 wherever it is used in the code
926
f = SVector(0,
927
0,
928
0,
929
0,
930
psi_jump,
931
0,
932
0,
933
0,
934
psi_jump)
935
end
936
937
return f
938
end
939
940
@inline function (noncons_flux::FluxNonConservativePowellLocalJump)(u_ll, u_rr,
941
normal_direction_avg::AbstractVector,
942
equations::IdealGlmMhdEquations2D,
943
nonconservative_type::NonConservativeJump,
944
nonconservative_term::Integer)
945
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
946
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
947
948
if nonconservative_term == 1
949
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
950
B1_jump = B1_rr - B1_ll
951
B2_jump = B2_rr - B2_ll
952
B_dot_n_jump = B1_jump * normal_direction_avg[1] +
953
B2_jump * normal_direction_avg[2]
954
f = SVector(0,
955
B_dot_n_jump,
956
B_dot_n_jump,
957
B_dot_n_jump,
958
B_dot_n_jump,
959
B_dot_n_jump,
960
B_dot_n_jump,
961
B_dot_n_jump,
962
0)
963
else # nonconservative_term == 2
964
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
965
psi_jump = (psi_rr - psi_ll)
966
f = SVector(0,
967
0,
968
0,
969
0,
970
psi_jump,
971
0,
972
0,
973
0,
974
psi_jump)
975
end
976
977
return f
978
end
979
980
"""
981
flux_derigs_etal(u_ll, u_rr, orientation, equations::IdealGlmMhdEquations2D)
982
983
Entropy conserving two-point flux by
984
- Derigs et al. (2018)
985
Ideal GLM-MHD: About the entropy consistent nine-wave magnetic field
986
divergence diminishing ideal magnetohydrodynamics equations
987
[DOI: 10.1016/j.jcp.2018.03.002](https://doi.org/10.1016/j.jcp.2018.03.002)
988
"""
989
function flux_derigs_etal(u_ll, u_rr, orientation::Integer,
990
equations::IdealGlmMhdEquations2D)
991
# Unpack left and right states to get velocities, pressure, and inverse temperature (called beta)
992
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
993
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
994
995
v1_ll = rho_v1_ll / rho_ll
996
v2_ll = rho_v2_ll / rho_ll
997
v3_ll = rho_v3_ll / rho_ll
998
v1_rr = rho_v1_rr / rho_rr
999
v2_rr = rho_v2_rr / rho_rr
1000
v3_rr = rho_v3_rr / rho_rr
1001
vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2
1002
vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2
1003
mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^2
1004
mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^2
1005
p_ll = (equations.gamma - 1) *
1006
(rho_e_total_ll - 0.5f0 * rho_ll * vel_norm_ll - 0.5f0 * mag_norm_ll -
1007
0.5f0 * psi_ll^2)
1008
p_rr = (equations.gamma - 1) *
1009
(rho_e_total_rr - 0.5f0 * rho_rr * vel_norm_rr - 0.5f0 * mag_norm_rr -
1010
0.5f0 * psi_rr^2)
1011
beta_ll = 0.5f0 * rho_ll / p_ll
1012
beta_rr = 0.5f0 * rho_rr / p_rr
1013
# for convenience store v⋅B
1014
vel_dot_mag_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
1015
vel_dot_mag_rr = v1_rr * B1_rr + v2_rr * B2_rr + v3_rr * B3_rr
1016
1017
# Compute the necessary mean values needed for either direction
1018
rho_avg = 0.5f0 * (rho_ll + rho_rr)
1019
rho_mean = ln_mean(rho_ll, rho_rr)
1020
beta_mean = ln_mean(beta_ll, beta_rr)
1021
beta_avg = 0.5f0 * (beta_ll + beta_rr)
1022
v1_avg = 0.5f0 * (v1_ll + v1_rr)
1023
v2_avg = 0.5f0 * (v2_ll + v2_rr)
1024
v3_avg = 0.5f0 * (v3_ll + v3_rr)
1025
p_mean = 0.5f0 * rho_avg / beta_avg
1026
B1_avg = 0.5f0 * (B1_ll + B1_rr)
1027
B2_avg = 0.5f0 * (B2_ll + B2_rr)
1028
B3_avg = 0.5f0 * (B3_ll + B3_rr)
1029
psi_avg = 0.5f0 * (psi_ll + psi_rr)
1030
vel_norm_avg = 0.5f0 * (vel_norm_ll + vel_norm_rr)
1031
mag_norm_avg = 0.5f0 * (mag_norm_ll + mag_norm_rr)
1032
vel_dot_mag_avg = 0.5f0 * (vel_dot_mag_ll + vel_dot_mag_rr)
1033
1034
# Calculate fluxes depending on orientation with specific direction averages
1035
if orientation == 1
1036
f1 = rho_mean * v1_avg
1037
f2 = f1 * v1_avg + p_mean + 0.5f0 * mag_norm_avg - B1_avg * B1_avg
1038
f3 = f1 * v2_avg - B1_avg * B2_avg
1039
f4 = f1 * v3_avg - B1_avg * B3_avg
1040
f6 = equations.c_h * psi_avg
1041
f7 = v1_avg * B2_avg - v2_avg * B1_avg
1042
f8 = v1_avg * B3_avg - v3_avg * B1_avg
1043
f9 = equations.c_h * B1_avg
1044
# total energy flux is complicated and involves the previous eight components
1045
psi_B1_avg = 0.5f0 * (B1_ll * psi_ll + B1_rr * psi_rr)
1046
v1_mag_avg = 0.5f0 * (v1_ll * mag_norm_ll + v1_rr * mag_norm_rr)
1047
f5 = (f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +
1048
f2 * v1_avg + f3 * v2_avg +
1049
f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg + f9 * psi_avg -
1050
0.5f0 * v1_mag_avg +
1051
B1_avg * vel_dot_mag_avg - equations.c_h * psi_B1_avg)
1052
else
1053
f1 = rho_mean * v2_avg
1054
f2 = f1 * v1_avg - B1_avg * B2_avg
1055
f3 = f1 * v2_avg + p_mean + 0.5f0 * mag_norm_avg - B2_avg * B2_avg
1056
f4 = f1 * v3_avg - B2_avg * B3_avg
1057
f6 = v2_avg * B1_avg - v1_avg * B2_avg
1058
f7 = equations.c_h * psi_avg
1059
f8 = v2_avg * B3_avg - v3_avg * B2_avg
1060
f9 = equations.c_h * B2_avg
1061
# total energy flux is complicated and involves the previous eight components
1062
psi_B2_avg = 0.5f0 * (B2_ll * psi_ll + B2_rr * psi_rr)
1063
v2_mag_avg = 0.5f0 * (v2_ll * mag_norm_ll + v2_rr * mag_norm_rr)
1064
f5 = (f1 * 0.5f0 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +
1065
f2 * v1_avg + f3 * v2_avg +
1066
f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg + f9 * psi_avg -
1067
0.5f0 * v2_mag_avg +
1068
B2_avg * vel_dot_mag_avg - equations.c_h * psi_B2_avg)
1069
end
1070
1071
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
1072
end
1073
1074
"""
1075
flux_hindenlang_gassner(u_ll, u_rr, orientation_or_normal_direction,
1076
equations::IdealGlmMhdEquations2D)
1077
1078
Entropy conserving and kinetic energy preserving two-point flux of
1079
Hindenlang and Gassner (2019), extending [`flux_ranocha`](@ref) to the MHD equations.
1080
1081
## References
1082
- Florian Hindenlang, Gregor Gassner (2019)
1083
A new entropy conservative two-point flux for ideal MHD equations derived from
1084
first principles.
1085
Presented at HONOM 2019: European workshop on high order numerical methods
1086
for evolutionary PDEs, theory and applications
1087
- Hendrik Ranocha (2018)
1088
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
1089
for Hyperbolic Balance Laws
1090
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
1091
- Hendrik Ranocha (2020)
1092
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
1093
the Euler Equations Using Summation-by-Parts Operators
1094
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
1095
"""
1096
@inline function flux_hindenlang_gassner(u_ll, u_rr, orientation::Integer,
1097
equations::IdealGlmMhdEquations2D)
1098
# Unpack left and right states
1099
rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll, psi_ll = cons2prim(u_ll,
1100
equations)
1101
rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr, psi_rr = cons2prim(u_rr,
1102
equations)
1103
1104
# Compute the necessary mean values needed for either direction
1105
rho_mean = ln_mean(rho_ll, rho_rr)
1106
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
1107
# in exact arithmetic since
1108
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
1109
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
1110
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
1111
v1_avg = 0.5f0 * (v1_ll + v1_rr)
1112
v2_avg = 0.5f0 * (v2_ll + v2_rr)
1113
v3_avg = 0.5f0 * (v3_ll + v3_rr)
1114
p_avg = 0.5f0 * (p_ll + p_rr)
1115
psi_avg = 0.5f0 * (psi_ll + psi_rr)
1116
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
1117
magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)
1118
1119
# Calculate fluxes depending on orientation with specific direction averages
1120
if orientation == 1
1121
f1 = rho_mean * v1_avg
1122
f2 = f1 * v1_avg + p_avg + magnetic_square_avg -
1123
0.5f0 * (B1_ll * B1_rr + B1_rr * B1_ll)
1124
f3 = f1 * v2_avg - 0.5f0 * (B1_ll * B2_rr + B1_rr * B2_ll)
1125
f4 = f1 * v3_avg - 0.5f0 * (B1_ll * B3_rr + B1_rr * B3_ll)
1126
#f5 below
1127
f6 = equations.c_h * psi_avg
1128
f7 = 0.5f0 * (v1_ll * B2_ll - v2_ll * B1_ll + v1_rr * B2_rr - v2_rr * B1_rr)
1129
f8 = 0.5f0 * (v1_ll * B3_ll - v3_ll * B1_ll + v1_rr * B3_rr - v3_rr * B1_rr)
1130
f9 = equations.c_h * 0.5f0 * (B1_ll + B1_rr)
1131
# total energy flux is complicated and involves the previous components
1132
f5 = (f1 *
1133
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
1134
+
1135
0.5f0 * (+p_ll * v1_rr + p_rr * v1_ll
1136
+ (v1_ll * B2_ll * B2_rr + v1_rr * B2_rr * B2_ll)
1137
+ (v1_ll * B3_ll * B3_rr + v1_rr * B3_rr * B3_ll)
1138
-
1139
(v2_ll * B1_ll * B2_rr + v2_rr * B1_rr * B2_ll)
1140
-
1141
(v3_ll * B1_ll * B3_rr + v3_rr * B1_rr * B3_ll)
1142
+
1143
equations.c_h * (B1_ll * psi_rr + B1_rr * psi_ll)))
1144
else # orientation == 2
1145
f1 = rho_mean * v2_avg
1146
f2 = f1 * v1_avg - 0.5f0 * (B2_ll * B1_rr + B2_rr * B1_ll)
1147
f3 = f1 * v2_avg + p_avg + magnetic_square_avg -
1148
0.5f0 * (B2_ll * B2_rr + B2_rr * B2_ll)
1149
f4 = f1 * v3_avg - 0.5f0 * (B2_ll * B3_rr + B2_rr * B3_ll)
1150
#f5 below
1151
f6 = 0.5f0 * (v2_ll * B1_ll - v1_ll * B2_ll + v2_rr * B1_rr - v1_rr * B2_rr)
1152
f7 = equations.c_h * psi_avg
1153
f8 = 0.5f0 * (v2_ll * B3_ll - v3_ll * B2_ll + v2_rr * B3_rr - v3_rr * B2_rr)
1154
f9 = equations.c_h * 0.5f0 * (B2_ll + B2_rr)
1155
# total energy flux is complicated and involves the previous components
1156
f5 = (f1 *
1157
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
1158
+
1159
0.5f0 * (+p_ll * v2_rr + p_rr * v2_ll
1160
+ (v2_ll * B1_ll * B1_rr + v2_rr * B1_rr * B1_ll)
1161
+ (v2_ll * B3_ll * B3_rr + v2_rr * B3_rr * B3_ll)
1162
-
1163
(v1_ll * B2_ll * B1_rr + v1_rr * B2_rr * B1_ll)
1164
-
1165
(v3_ll * B2_ll * B3_rr + v3_rr * B2_rr * B3_ll)
1166
+
1167
equations.c_h * (B2_ll * psi_rr + B2_rr * psi_ll)))
1168
end
1169
1170
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
1171
end
1172
1173
@inline function flux_hindenlang_gassner(u_ll, u_rr, normal_direction::AbstractVector,
1174
equations::IdealGlmMhdEquations2D)
1175
# Unpack left and right states
1176
rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll, psi_ll = cons2prim(u_ll,
1177
equations)
1178
rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr, psi_rr = cons2prim(u_rr,
1179
equations)
1180
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
1181
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
1182
B_dot_n_ll = B1_ll * normal_direction[1] + B2_ll * normal_direction[2]
1183
B_dot_n_rr = B1_rr * normal_direction[1] + B2_rr * normal_direction[2]
1184
1185
# Compute the necessary mean values needed for either direction
1186
rho_mean = ln_mean(rho_ll, rho_rr)
1187
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
1188
# in exact arithmetic since
1189
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
1190
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
1191
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
1192
v1_avg = 0.5f0 * (v1_ll + v1_rr)
1193
v2_avg = 0.5f0 * (v2_ll + v2_rr)
1194
v3_avg = 0.5f0 * (v3_ll + v3_rr)
1195
p_avg = 0.5f0 * (p_ll + p_rr)
1196
psi_avg = 0.5f0 * (psi_ll + psi_rr)
1197
velocity_square_avg = 0.5f0 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
1198
magnetic_square_avg = 0.5f0 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)
1199
1200
# Calculate fluxes depending on normal_direction
1201
f1 = rho_mean * 0.5f0 * (v_dot_n_ll + v_dot_n_rr)
1202
f2 = (f1 * v1_avg + (p_avg + magnetic_square_avg) * normal_direction[1]
1203
-
1204
0.5f0 * (B_dot_n_ll * B1_rr + B_dot_n_rr * B1_ll))
1205
f3 = (f1 * v2_avg + (p_avg + magnetic_square_avg) * normal_direction[2]
1206
-
1207
0.5f0 * (B_dot_n_ll * B2_rr + B_dot_n_rr * B2_ll))
1208
f4 = (f1 * v3_avg
1209
-
1210
0.5f0 * (B_dot_n_ll * B3_rr + B_dot_n_rr * B3_ll))
1211
#f5 below
1212
f6 = (equations.c_h * psi_avg * normal_direction[1]
1213
+
1214
0.5f0 * (v_dot_n_ll * B1_ll - v1_ll * B_dot_n_ll +
1215
v_dot_n_rr * B1_rr - v1_rr * B_dot_n_rr))
1216
f7 = (equations.c_h * psi_avg * normal_direction[2]
1217
+
1218
0.5f0 * (v_dot_n_ll * B2_ll - v2_ll * B_dot_n_ll +
1219
v_dot_n_rr * B2_rr - v2_rr * B_dot_n_rr))
1220
f8 = +0.5f0 * (v_dot_n_ll * B3_ll - v3_ll * B_dot_n_ll +
1221
v_dot_n_rr * B3_rr - v3_rr * B_dot_n_rr)
1222
f9 = equations.c_h * 0.5f0 * (B_dot_n_ll + B_dot_n_rr)
1223
# total energy flux is complicated and involves the previous components
1224
f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
1225
+
1226
0.5f0 * (+p_ll * v_dot_n_rr + p_rr * v_dot_n_ll
1227
+ (v_dot_n_ll * B1_ll * B1_rr + v_dot_n_rr * B1_rr * B1_ll)
1228
+ (v_dot_n_ll * B2_ll * B2_rr + v_dot_n_rr * B2_rr * B2_ll)
1229
+ (v_dot_n_ll * B3_ll * B3_rr + v_dot_n_rr * B3_rr * B3_ll)
1230
-
1231
(v1_ll * B_dot_n_ll * B1_rr + v1_rr * B_dot_n_rr * B1_ll)
1232
-
1233
(v2_ll * B_dot_n_ll * B2_rr + v2_rr * B_dot_n_rr * B2_ll)
1234
-
1235
(v3_ll * B_dot_n_ll * B3_rr + v3_rr * B_dot_n_rr * B3_ll)
1236
+
1237
equations.c_h * (B_dot_n_ll * psi_rr + B_dot_n_rr * psi_ll)))
1238
1239
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
1240
end
1241
1242
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation
1243
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
1244
equations::IdealGlmMhdEquations2D)
1245
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
1246
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
1247
1248
# Calculate the left/right velocities and fast magnetoacoustic wave speeds
1249
if orientation == 1
1250
v_ll = rho_v1_ll / rho_ll
1251
v_rr = rho_v1_rr / rho_rr
1252
else # orientation == 2
1253
v_ll = rho_v2_ll / rho_ll
1254
v_rr = rho_v2_rr / rho_rr
1255
end
1256
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
1257
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
1258
1259
return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)
1260
end
1261
1262
@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
1263
equations::IdealGlmMhdEquations2D)
1264
# return max(v_mag_ll, v_mag_rr) + max(cf_ll, cf_rr)
1265
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
1266
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
1267
1268
# Calculate normal velocities and fast magnetoacoustic wave speeds
1269
# left
1270
v1_ll = rho_v1_ll / rho_ll
1271
v2_ll = rho_v2_ll / rho_ll
1272
v_ll = (v1_ll * normal_direction[1]
1273
+
1274
v2_ll * normal_direction[2])
1275
cf_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
1276
# right
1277
v1_rr = rho_v1_rr / rho_rr
1278
v2_rr = rho_v2_rr / rho_rr
1279
v_rr = (v1_rr * normal_direction[1]
1280
+
1281
v2_rr * normal_direction[2])
1282
cf_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
1283
1284
# wave speeds already scaled by norm(normal_direction) in [`calc_fast_wavespeed`](@ref)
1285
return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)
1286
end
1287
1288
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
1289
@inline function max_abs_speed(u_ll, u_rr, orientation::Integer,
1290
equations::IdealGlmMhdEquations2D)
1291
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
1292
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
1293
1294
# Calculate the left/right velocities and fast magnetoacoustic wave speeds
1295
if orientation == 1
1296
v_ll = rho_v1_ll / rho_ll
1297
v_rr = rho_v1_rr / rho_rr
1298
else # orientation == 2
1299
v_ll = rho_v2_ll / rho_ll
1300
v_rr = rho_v2_rr / rho_rr
1301
end
1302
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
1303
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
1304
1305
return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)
1306
end
1307
1308
# Less "cautious", i.e., less overestimating `λ_max` compared to `max_abs_speed_naive`
1309
@inline function max_abs_speed(u_ll, u_rr, normal_direction::AbstractVector,
1310
equations::IdealGlmMhdEquations2D)
1311
# return max(v_mag_ll, v_mag_rr) + max(cf_ll, cf_rr)
1312
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
1313
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
1314
1315
# Calculate normal velocities and fast magnetoacoustic wave speeds
1316
# left
1317
v1_ll = rho_v1_ll / rho_ll
1318
v2_ll = rho_v2_ll / rho_ll
1319
v_ll = (v1_ll * normal_direction[1]
1320
+
1321
v2_ll * normal_direction[2])
1322
cf_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
1323
# right
1324
v1_rr = rho_v1_rr / rho_rr
1325
v2_rr = rho_v2_rr / rho_rr
1326
v_rr = (v1_rr * normal_direction[1]
1327
+
1328
v2_rr * normal_direction[2])
1329
cf_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
1330
1331
# wave speeds already scaled by norm(normal_direction) in [`calc_fast_wavespeed`](@ref)
1332
return max(abs(v_ll) + cf_ll, abs(v_rr) + cf_rr)
1333
end
1334
1335
# Calculate estimate for minimum and maximum wave speeds for HLL-type fluxes
1336
@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,
1337
equations::IdealGlmMhdEquations2D)
1338
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
1339
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
1340
1341
# Calculate primitive velocity variables
1342
v1_ll = rho_v1_ll / rho_ll
1343
v2_ll = rho_v2_ll / rho_ll
1344
1345
v1_rr = rho_v1_rr / rho_rr
1346
v2_rr = rho_v2_rr / rho_rr
1347
1348
# Approximate the left-most and right-most eigenvalues in the Riemann fan
1349
if orientation == 1 # x-direction
1350
λ_min = v1_ll - calc_fast_wavespeed(u_ll, orientation, equations)
1351
λ_max = v1_rr + calc_fast_wavespeed(u_rr, orientation, equations)
1352
else # y-direction
1353
λ_min = v2_ll - calc_fast_wavespeed(u_ll, orientation, equations)
1354
λ_max = v2_rr + calc_fast_wavespeed(u_rr, orientation, equations)
1355
end
1356
1357
return λ_min, λ_max
1358
end
1359
1360
@inline function min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
1361
equations::IdealGlmMhdEquations2D)
1362
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
1363
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
1364
1365
# Calculate primitive velocity variables
1366
v1_ll = rho_v1_ll / rho_ll
1367
v2_ll = rho_v2_ll / rho_ll
1368
1369
v1_rr = rho_v1_rr / rho_rr
1370
v2_rr = rho_v2_rr / rho_rr
1371
1372
v_normal_ll = (v1_ll * normal_direction[1] + v2_ll * normal_direction[2])
1373
v_normal_rr = (v1_rr * normal_direction[1] + v2_rr * normal_direction[2])
1374
1375
c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
1376
c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
1377
1378
# Estimate the min/max eigenvalues in the normal direction
1379
λ_min = min(v_normal_ll - c_f_ll, v_normal_rr - c_f_rr)
1380
λ_max = max(v_normal_rr + c_f_rr, v_normal_rr + c_f_rr)
1381
1382
return λ_min, λ_max
1383
end
1384
1385
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
1386
@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,
1387
equations::IdealGlmMhdEquations2D)
1388
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
1389
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
1390
1391
# Calculate primitive velocity variables
1392
v1_ll = rho_v1_ll / rho_ll
1393
v2_ll = rho_v2_ll / rho_ll
1394
1395
v1_rr = rho_v1_rr / rho_rr
1396
v2_rr = rho_v2_rr / rho_rr
1397
1398
# Approximate the left-most and right-most eigenvalues in the Riemann fan
1399
if orientation == 1 # x-direction
1400
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
1401
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
1402
1403
λ_min = min(v1_ll - c_f_ll, v1_rr - c_f_rr)
1404
λ_max = max(v1_ll + c_f_ll, v1_rr + c_f_rr)
1405
else # y-direction
1406
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
1407
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
1408
1409
λ_min = min(v2_ll - c_f_ll, v2_rr - c_f_rr)
1410
λ_max = max(v2_ll + c_f_ll, v1_rr + c_f_rr)
1411
end
1412
1413
return λ_min, λ_max
1414
end
1415
1416
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
1417
@inline function min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector,
1418
equations::IdealGlmMhdEquations2D)
1419
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
1420
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
1421
1422
# Calculate primitive velocity variables
1423
v1_ll = rho_v1_ll / rho_ll
1424
v2_ll = rho_v2_ll / rho_ll
1425
1426
v1_rr = rho_v1_rr / rho_rr
1427
v2_rr = rho_v2_rr / rho_rr
1428
1429
v_normal_ll = (v1_ll * normal_direction[1] + v2_ll * normal_direction[2])
1430
v_normal_rr = (v1_rr * normal_direction[1] + v2_rr * normal_direction[2])
1431
1432
c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
1433
c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
1434
1435
# Estimate the min/max eigenvalues in the normal direction
1436
λ_min = min(v_normal_ll - c_f_ll, v_normal_rr - c_f_rr)
1437
λ_max = max(v_normal_ll + c_f_ll, v_normal_rr + c_f_rr)
1438
1439
return λ_min, λ_max
1440
end
1441
1442
"""
1443
min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer, equations::IdealGlmMhdEquations2D)
1444
1445
Calculate minimum and maximum wave speeds for HLL-type fluxes as in
1446
- Li (2005)
1447
An HLLC Riemann solver for magneto-hydrodynamics
1448
[DOI: 10.1016/j.jcp.2004.08.020](https://doi.org/10.1016/j.jcp.2004.08.020).
1449
"""
1450
@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,
1451
equations::IdealGlmMhdEquations2D)
1452
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
1453
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
1454
1455
# Calculate primitive velocity variables
1456
v1_ll = rho_v1_ll / rho_ll
1457
v2_ll = rho_v2_ll / rho_ll
1458
1459
v1_rr = rho_v1_rr / rho_rr
1460
v2_rr = rho_v2_rr / rho_rr
1461
1462
# Approximate the left-most and right-most eigenvalues in the Riemann fan
1463
if orientation == 1 # x-direction
1464
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
1465
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
1466
vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)
1467
λ_min = min(v1_ll - c_f_ll, vel_roe - c_f_roe)
1468
λ_max = max(v1_rr + c_f_rr, vel_roe + c_f_roe)
1469
else # y-direction
1470
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
1471
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
1472
vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)
1473
λ_min = min(v2_ll - c_f_ll, vel_roe - c_f_roe)
1474
λ_max = max(v2_rr + c_f_rr, vel_roe + c_f_roe)
1475
end
1476
1477
return λ_min, λ_max
1478
end
1479
1480
@inline function min_max_speed_einfeldt(u_ll, u_rr, normal_direction::AbstractVector,
1481
equations::IdealGlmMhdEquations2D)
1482
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
1483
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
1484
1485
# Calculate primitive velocity variables
1486
v1_ll = rho_v1_ll / rho_ll
1487
v2_ll = rho_v2_ll / rho_ll
1488
1489
v1_rr = rho_v1_rr / rho_rr
1490
v2_rr = rho_v2_rr / rho_rr
1491
1492
v_normal_ll = (v1_ll * normal_direction[1] + v2_ll * normal_direction[2])
1493
v_normal_rr = (v1_rr * normal_direction[1] + v2_rr * normal_direction[2])
1494
1495
c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
1496
c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
1497
v_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, normal_direction, equations)
1498
1499
# Estimate the min/max eigenvalues in the normal direction
1500
λ_min = min(v_normal_ll - c_f_ll, v_roe - c_f_roe)
1501
λ_max = max(v_normal_rr + c_f_rr, v_roe + c_f_roe)
1502
1503
return λ_min, λ_max
1504
end
1505
1506
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction
1507
# has been normalized prior to this rotation of the state vector
1508
@inline function rotate_to_x(u, normal_vector, equations::IdealGlmMhdEquations2D)
1509
# cos and sin of the angle between the x-axis and the normalized normal_vector are
1510
# the normalized vector's x and y coordinates respectively (see unit circle).
1511
c = normal_vector[1]
1512
s = normal_vector[2]
1513
1514
# Apply the 2D rotation matrix with normal and tangent directions of the form
1515
# [ 1 0 0 0 0 0 0 0 0;
1516
# 0 n_1 n_2 0 0 0 0 0 0;
1517
# 0 t_1 t_2 0 0 0 0 0 0;
1518
# 0 0 0 1 0 0 0 0 0;
1519
# 0 0 0 0 1 0 0 0 0;
1520
# 0 0 0 0 0 n_1 n_2 0 0;
1521
# 0 0 0 0 0 t_1 t_2 0 0;
1522
# 0 0 0 0 0 0 0 1 0;
1523
# 0 0 0 0 0 0 0 0 1 ]
1524
# where t_1 = -n_2 and t_2 = n_1.
1525
# Note for IdealGlmMhdEquations2D only the velocities and magnetic field variables rotate
1526
1527
return SVector(u[1],
1528
c * u[2] + s * u[3],
1529
-s * u[2] + c * u[3],
1530
u[4],
1531
u[5],
1532
c * u[6] + s * u[7],
1533
-s * u[6] + c * u[7],
1534
u[8],
1535
u[9])
1536
end
1537
1538
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction
1539
# has been normalized prior to this back-rotation of the state vector
1540
@inline function rotate_from_x(u, normal_vector, equations::IdealGlmMhdEquations2D)
1541
# cos and sin of the angle between the x-axis and the normalized normal_vector are
1542
# the normalized vector's x and y coordinates respectively (see unit circle).
1543
c = normal_vector[1]
1544
s = normal_vector[2]
1545
1546
# Apply the 2D back-rotation matrix with normal and tangent directions of the form
1547
# [ 1 0 0 0 0 0 0 0 0;
1548
# 0 n_1 t_1 0 0 0 0 0 0;
1549
# 0 n_2 t_2 0 0 0 0 0 0;
1550
# 0 0 0 1 0 0 0 0 0;
1551
# 0 0 0 0 1 0 0 0 0;
1552
# 0 0 0 0 0 n_1 t_1 0 0;
1553
# 0 0 0 0 0 n_2 t_2 0 0;
1554
# 0 0 0 0 0 0 0 1 0;
1555
# 0 0 0 0 0 0 0 0 1 ]
1556
# where t_1 = -n_2 and t_2 = n_1.
1557
# Note for IdealGlmMhdEquations2D the velocities and magnetic field variables back-rotate
1558
1559
return SVector(u[1],
1560
c * u[2] - s * u[3],
1561
s * u[2] + c * u[3],
1562
u[4],
1563
u[5],
1564
c * u[6] - s * u[7],
1565
s * u[6] + c * u[7],
1566
u[8],
1567
u[9])
1568
end
1569
1570
@inline function max_abs_speeds(u, equations::IdealGlmMhdEquations2D)
1571
rho, rho_v1, rho_v2, rho_v3, _ = u
1572
v1 = rho_v1 / rho
1573
v2 = rho_v2 / rho
1574
v3 = rho_v3 / rho
1575
cf_x_direction = calc_fast_wavespeed(u, 1, equations)
1576
cf_y_direction = calc_fast_wavespeed(u, 2, equations)
1577
1578
return abs(v1) + cf_x_direction, abs(v2) + cf_y_direction
1579
end
1580
1581
# Convert conservative variables to primitive
1582
@inline function cons2prim(u, equations::IdealGlmMhdEquations2D)
1583
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u
1584
1585
v1 = rho_v1 / rho
1586
v2 = rho_v2 / rho
1587
v3 = rho_v3 / rho
1588
p = (equations.gamma - 1) * (rho_e_total -
1589
0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3
1590
+ B1 * B1 + B2 * B2 + B3 * B3
1591
+ psi * psi))
1592
1593
return SVector(rho, v1, v2, v3, p, B1, B2, B3, psi)
1594
end
1595
1596
# Convert conservative variables to entropy variables
1597
@inline function cons2entropy(u, equations::IdealGlmMhdEquations2D)
1598
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u
1599
1600
v1 = rho_v1 / rho
1601
v2 = rho_v2 / rho
1602
v3 = rho_v3 / rho
1603
v_square = v1^2 + v2^2 + v3^2
1604
p = (equations.gamma - 1) *
1605
(rho_e_total - 0.5f0 * rho * v_square - 0.5f0 * (B1^2 + B2^2 + B3^2) -
1606
0.5f0 * psi^2)
1607
s = log(p) - equations.gamma * log(rho)
1608
rho_p = rho / p
1609
1610
w1 = (equations.gamma - s) * equations.inv_gamma_minus_one -
1611
0.5f0 * rho_p * v_square
1612
w2 = rho_p * v1
1613
w3 = rho_p * v2
1614
w4 = rho_p * v3
1615
w5 = -rho_p
1616
w6 = rho_p * B1
1617
w7 = rho_p * B2
1618
w8 = rho_p * B3
1619
w9 = rho_p * psi
1620
1621
return SVector(w1, w2, w3, w4, w5, w6, w7, w8, w9)
1622
end
1623
1624
# Convert entropy variables to conservative variables
1625
@inline function entropy2cons(w, equations::IdealGlmMhdEquations2D)
1626
w1, w2, w3, w4, w5, w6, w7, w8, w9 = w
1627
1628
v1 = -w2 / w5
1629
v2 = -w3 / w5
1630
v3 = -w4 / w5
1631
1632
B1 = -w6 / w5
1633
B2 = -w7 / w5
1634
B3 = -w8 / w5
1635
psi = -w9 / w5
1636
1637
# This imitates what is done for compressible Euler 3D `entropy2cons`: we convert from
1638
# the entropy variables for `-rho * s / (gamma - 1)` to the entropy variables for the entropy
1639
# `-rho * s` used by Hughes, Franca, Mallet (1986).
1640
@unpack gamma = equations
1641
V1, V2, V3, V4, V5 = SVector(w1, w2, w3, w4, w5) * (gamma - 1)
1642
s = gamma - V1 + (V2^2 + V3^2 + V4^2) / (2 * V5)
1643
rho_iota = ((gamma - 1) / (-V5)^gamma)^(equations.inv_gamma_minus_one) *
1644
exp(-s * equations.inv_gamma_minus_one)
1645
rho = -rho_iota * V5
1646
p = -rho / w5
1647
1648
return prim2cons(SVector(rho, v1, v2, v3, p, B1, B2, B3, psi), equations)
1649
end
1650
1651
# Convert primitive to conservative variables
1652
@inline function prim2cons(prim, equations::IdealGlmMhdEquations2D)
1653
rho, v1, v2, v3, p, B1, B2, B3, psi = prim
1654
1655
rho_v1 = rho * v1
1656
rho_v2 = rho * v2
1657
rho_v3 = rho * v3
1658
rho_e_total = p * equations.inv_gamma_minus_one +
1659
0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3) +
1660
0.5f0 * (B1^2 + B2^2 + B3^2) + 0.5f0 * psi^2
1661
1662
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi)
1663
end
1664
1665
@inline function density(u, equations::IdealGlmMhdEquations2D)
1666
rho = u[1]
1667
return rho
1668
end
1669
1670
@inline function velocity(u, equations::IdealGlmMhdEquations2D)
1671
rho = u[1]
1672
v1 = u[2] / rho
1673
v2 = u[3] / rho
1674
v3 = u[4] / rho
1675
return SVector(v1, v2, v3)
1676
end
1677
1678
@inline function velocity(u, orientation::Int, equations::IdealGlmMhdEquations2D)
1679
rho = u[1]
1680
v = u[orientation + 1] / rho
1681
return v
1682
end
1683
1684
@doc raw"""
1685
pressure(u, equations::IdealGlmMhdEquations2D)
1686
1687
Computes the pressure for an ideal equation of state with
1688
isentropic exponent/adiabatic index ``\gamma`` from the conserved variables `u`.
1689
```math
1690
\begin{aligned}
1691
p &= (\gamma - 1) \left( \rho e_{\text{total}} - \rho e_{\text{kinetic}} - \rho e_{\text{magnetic}} - \frac{1}{2}\psi^2 \right) \\
1692
&= (\gamma - 1) \left( \rho e - \frac{1}{2}
1693
\left[\rho \Vert v \Vert_2^2 + \Vert B \Vert_2^2 + \psi^2 \right] \right)
1694
\end{aligned}
1695
```
1696
"""
1697
@inline function pressure(u, equations::IdealGlmMhdEquations2D)
1698
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u
1699
p = (equations.gamma - 1) * (rho_e_total -
1700
0.5f0 *
1701
((rho_v1^2 + rho_v2^2 + rho_v3^2) / rho +
1702
B1^2 + B2^2 + B3^2 + psi^2))
1703
return p
1704
end
1705
1706
# Transformation from conservative variables u to d(p)/d(u)
1707
@inline function gradient_conservative(::typeof(pressure),
1708
u, equations::IdealGlmMhdEquations2D)
1709
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u
1710
1711
v1 = rho_v1 / rho
1712
v2 = rho_v2 / rho
1713
v3 = rho_v3 / rho
1714
v_square = v1^2 + v2^2 + v3^2
1715
1716
return (equations.gamma - 1) *
1717
SVector(0.5f0 * v_square, -v1, -v2, -v3, 1, -B1, -B2, -B3, -psi)
1718
end
1719
1720
@inline function density_pressure(u, equations::IdealGlmMhdEquations2D)
1721
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = u
1722
rho_times_p = (equations.gamma - 1) * (rho * rho_e_total -
1723
0.5f0 * (rho_v1^2 + rho_v2^2 + rho_v3^2 +
1724
rho * (B1^2 + B2^2 + B3^2 + psi^2)))
1725
return rho_times_p
1726
end
1727
1728
# Compute the fastest wave speed for ideal MHD equations: c_f, the fast magnetoacoustic eigenvalue
1729
@inline function calc_fast_wavespeed(cons, orientation::Integer,
1730
equations::IdealGlmMhdEquations2D)
1731
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = cons
1732
v1 = rho_v1 / rho
1733
v2 = rho_v2 / rho
1734
v3 = rho_v3 / rho
1735
kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
1736
mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)
1737
p = (equations.gamma - 1) * (rho_e_total - kin_en - mag_en - 0.5f0 * psi^2)
1738
a_square = equations.gamma * p / rho
1739
sqrt_rho = sqrt(rho)
1740
b1 = B1 / sqrt_rho
1741
b2 = B2 / sqrt_rho
1742
b3 = B3 / sqrt_rho
1743
b_square = b1 * b1 + b2 * b2 + b3 * b3
1744
if orientation == 1 # x-direction
1745
c_f = sqrt(0.5f0 * (a_square + b_square) +
1746
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b1^2))
1747
else
1748
c_f = sqrt(0.5f0 * (a_square + b_square) +
1749
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b2^2))
1750
end
1751
return c_f
1752
end
1753
1754
@inline function calc_fast_wavespeed(cons, normal_direction::AbstractVector,
1755
equations::IdealGlmMhdEquations2D)
1756
rho, rho_v1, rho_v2, rho_v3, rho_e_total, B1, B2, B3, psi = cons
1757
v1 = rho_v1 / rho
1758
v2 = rho_v2 / rho
1759
v3 = rho_v3 / rho
1760
kin_en = 0.5f0 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
1761
mag_en = 0.5f0 * (B1 * B1 + B2 * B2 + B3 * B3)
1762
p = (equations.gamma - 1) * (rho_e_total - kin_en - mag_en - 0.5f0 * psi^2)
1763
a_square = equations.gamma * p / rho
1764
sqrt_rho = sqrt(rho)
1765
b1 = B1 / sqrt_rho
1766
b2 = B2 / sqrt_rho
1767
b3 = B3 / sqrt_rho
1768
b_square = b1 * b1 + b2 * b2 + b3 * b3
1769
norm_squared = (normal_direction[1] * normal_direction[1] +
1770
normal_direction[2] * normal_direction[2])
1771
b_dot_n_squared = (b1 * normal_direction[1] +
1772
b2 * normal_direction[2])^2 / norm_squared
1773
1774
c_f = sqrt((0.5f0 * (a_square + b_square) +
1775
0.5f0 * sqrt((a_square + b_square)^2 - 4 * a_square * b_dot_n_squared)) *
1776
norm_squared)
1777
return c_f
1778
end
1779
1780
"""
1781
calc_fast_wavespeed_roe(u_ll, u_rr, orientation_or_normal_direction, equations::IdealGlmMhdEquations2D)
1782
1783
Compute the fast magnetoacoustic wave speed using Roe averages
1784
as given by
1785
- Cargo and Gallice (1997)
1786
Roe Matrices for Ideal MHD and Systematic Construction
1787
of Roe Matrices for Systems of Conservation Laws
1788
[DOI: 10.1006/jcph.1997.5773](https://doi.org/10.1006/jcph.1997.5773)
1789
"""
1790
@inline function calc_fast_wavespeed_roe(u_ll, u_rr, orientation::Integer,
1791
equations::IdealGlmMhdEquations2D)
1792
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
1793
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
1794
1795
# Calculate primitive variables
1796
v1_ll = rho_v1_ll / rho_ll
1797
v2_ll = rho_v2_ll / rho_ll
1798
v3_ll = rho_v3_ll / rho_ll
1799
kin_en_ll = 0.5f0 * (rho_v1_ll * v1_ll + rho_v2_ll * v2_ll + rho_v3_ll * v3_ll)
1800
mag_norm_ll = B1_ll * B1_ll + B2_ll * B2_ll + B3_ll * B3_ll
1801
p_ll = (equations.gamma - 1) *
1802
(rho_e_total_ll - kin_en_ll - 0.5f0 * mag_norm_ll - 0.5f0 * psi_ll^2)
1803
1804
v1_rr = rho_v1_rr / rho_rr
1805
v2_rr = rho_v2_rr / rho_rr
1806
v3_rr = rho_v3_rr / rho_rr
1807
kin_en_rr = 0.5f0 * (rho_v1_rr * v1_rr + rho_v2_rr * v2_rr + rho_v3_rr * v3_rr)
1808
mag_norm_rr = B1_rr * B1_rr + B2_rr * B2_rr + B3_rr * B3_rr
1809
p_rr = (equations.gamma - 1) *
1810
(rho_e_total_rr - kin_en_rr - 0.5f0 * mag_norm_rr - 0.5f0 * psi_rr^2)
1811
1812
# compute total pressure which is thermal + magnetic pressures
1813
p_total_ll = p_ll + 0.5f0 * mag_norm_ll
1814
p_total_rr = p_rr + 0.5f0 * mag_norm_rr
1815
1816
# compute the Roe density averages
1817
sqrt_rho_ll = sqrt(rho_ll)
1818
sqrt_rho_rr = sqrt(rho_rr)
1819
inv_sqrt_rho_add = 1 / (sqrt_rho_ll + sqrt_rho_rr)
1820
inv_sqrt_rho_prod = 1 / (sqrt_rho_ll * sqrt_rho_rr)
1821
rho_ll_roe = sqrt_rho_ll * inv_sqrt_rho_add
1822
rho_rr_roe = sqrt_rho_rr * inv_sqrt_rho_add
1823
# Roe averages
1824
# velocities and magnetic fields
1825
v1_roe = v1_ll * rho_ll_roe + v1_rr * rho_rr_roe
1826
v2_roe = v2_ll * rho_ll_roe + v2_rr * rho_rr_roe
1827
v3_roe = v3_ll * rho_ll_roe + v3_rr * rho_rr_roe
1828
B1_roe = B1_ll * rho_ll_roe + B1_rr * rho_rr_roe
1829
B2_roe = B2_ll * rho_ll_roe + B2_rr * rho_rr_roe
1830
B3_roe = B3_ll * rho_ll_roe + B3_rr * rho_rr_roe
1831
# enthalpy
1832
H_ll = (rho_e_total_ll + p_total_ll) / rho_ll
1833
H_rr = (rho_e_total_rr + p_total_rr) / rho_rr
1834
H_roe = H_ll * rho_ll_roe + H_rr * rho_rr_roe
1835
# temporary variable see equation (4.12) in Cargo and Gallice
1836
X = 0.5f0 * ((B1_ll - B1_rr)^2 + (B2_ll - B2_rr)^2 + (B3_ll - B3_rr)^2) *
1837
inv_sqrt_rho_add^2
1838
# averaged components needed to compute c_f, the fast magnetoacoustic wave speed
1839
b_square_roe = (B1_roe^2 + B2_roe^2 + B3_roe^2) * inv_sqrt_rho_prod # scaled magnectic sum
1840
a_square_roe = ((2 - equations.gamma) * X +
1841
(equations.gamma - 1) *
1842
(H_roe - 0.5f0 * (v1_roe^2 + v2_roe^2 + v3_roe^2) -
1843
b_square_roe)) # acoustic speed
1844
# finally compute the average wave speed and set the output velocity (depends on orientation)
1845
if orientation == 1 # x-direction
1846
c_a_roe = B1_roe^2 * inv_sqrt_rho_prod # (squared) Alfvén wave speed
1847
a_star_roe = sqrt((a_square_roe + b_square_roe)^2 -
1848
4 * a_square_roe * c_a_roe)
1849
c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe))
1850
vel_out_roe = v1_roe
1851
else # y-direction
1852
c_a_roe = B2_roe^2 * inv_sqrt_rho_prod # (squared) Alfvén wave speed
1853
a_star_roe = sqrt((a_square_roe + b_square_roe)^2 -
1854
4 * a_square_roe * c_a_roe)
1855
c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe))
1856
vel_out_roe = v2_roe
1857
end
1858
1859
return vel_out_roe, c_f_roe
1860
end
1861
1862
@inline function calc_fast_wavespeed_roe(u_ll, u_rr, normal_direction::AbstractVector,
1863
equations::IdealGlmMhdEquations2D)
1864
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_total_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
1865
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_total_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
1866
1867
# Calculate primitive variables
1868
v1_ll = rho_v1_ll / rho_ll
1869
v2_ll = rho_v2_ll / rho_ll
1870
v3_ll = rho_v3_ll / rho_ll
1871
kin_en_ll = 0.5f0 * (rho_v1_ll * v1_ll + rho_v2_ll * v2_ll + rho_v3_ll * v3_ll)
1872
mag_norm_ll = B1_ll * B1_ll + B2_ll * B2_ll + B3_ll * B3_ll
1873
p_ll = (equations.gamma - 1) *
1874
(rho_e_total_ll - kin_en_ll - 0.5f0 * mag_norm_ll - 0.5f0 * psi_ll^2)
1875
1876
v1_rr = rho_v1_rr / rho_rr
1877
v2_rr = rho_v2_rr / rho_rr
1878
v3_rr = rho_v3_rr / rho_rr
1879
kin_en_rr = 0.5f0 * (rho_v1_rr * v1_rr + rho_v2_rr * v2_rr + rho_v3_rr * v3_rr)
1880
mag_norm_rr = B1_rr * B1_rr + B2_rr * B2_rr + B3_rr * B3_rr
1881
p_rr = (equations.gamma - 1) *
1882
(rho_e_total_rr - kin_en_rr - 0.5f0 * mag_norm_rr - 0.5f0 * psi_rr^2)
1883
1884
# compute total pressure which is thermal + magnetic pressures
1885
p_total_ll = p_ll + 0.5f0 * mag_norm_ll
1886
p_total_rr = p_rr + 0.5f0 * mag_norm_rr
1887
1888
# compute the Roe density averages
1889
sqrt_rho_ll = sqrt(rho_ll)
1890
sqrt_rho_rr = sqrt(rho_rr)
1891
inv_sqrt_rho_add = 1 / (sqrt_rho_ll + sqrt_rho_rr)
1892
inv_sqrt_rho_prod = 1 / (sqrt_rho_ll * sqrt_rho_rr)
1893
rho_ll_roe = sqrt_rho_ll * inv_sqrt_rho_add
1894
rho_rr_roe = sqrt_rho_rr * inv_sqrt_rho_add
1895
# Roe averages
1896
# velocities and magnetic fields
1897
v1_roe = v1_ll * rho_ll_roe + v1_rr * rho_rr_roe
1898
v2_roe = v2_ll * rho_ll_roe + v2_rr * rho_rr_roe
1899
v3_roe = v3_ll * rho_ll_roe + v3_rr * rho_rr_roe
1900
B1_roe = B1_ll * rho_ll_roe + B1_rr * rho_rr_roe
1901
B2_roe = B2_ll * rho_ll_roe + B2_rr * rho_rr_roe
1902
B3_roe = B3_ll * rho_ll_roe + B3_rr * rho_rr_roe
1903
# enthalpy
1904
H_ll = (rho_e_total_ll + p_total_ll) / rho_ll
1905
H_rr = (rho_e_total_rr + p_total_rr) / rho_rr
1906
H_roe = H_ll * rho_ll_roe + H_rr * rho_rr_roe
1907
# temporary variable see equation (4.12) in Cargo and Gallice
1908
X = 0.5f0 * ((B1_ll - B1_rr)^2 + (B2_ll - B2_rr)^2 + (B3_ll - B3_rr)^2) *
1909
inv_sqrt_rho_add^2
1910
# averaged components needed to compute c_f, the fast magnetoacoustic wave speed
1911
b_square_roe = (B1_roe^2 + B2_roe^2 + B3_roe^2) * inv_sqrt_rho_prod # scaled magnectic sum
1912
a_square_roe = ((2 - equations.gamma) * X +
1913
(equations.gamma - 1) *
1914
(H_roe - 0.5f0 * (v1_roe^2 + v2_roe^2 + v3_roe^2) -
1915
b_square_roe)) # acoustic speed
1916
1917
# finally compute the average wave speed and set the output velocity (depends on orientation)
1918
norm_squared = (normal_direction[1] * normal_direction[1] +
1919
normal_direction[2] * normal_direction[2])
1920
B_roe_dot_n_squared = (B1_roe * normal_direction[1] +
1921
B2_roe * normal_direction[2])^2 / norm_squared
1922
1923
c_a_roe = B_roe_dot_n_squared * inv_sqrt_rho_prod # (squared) Alfvén wave speed
1924
a_star_roe = sqrt((a_square_roe + b_square_roe)^2 - 4 * a_square_roe * c_a_roe)
1925
c_f_roe = sqrt(0.5f0 * (a_square_roe + b_square_roe + a_star_roe) * norm_squared)
1926
vel_out_roe = (v1_roe * normal_direction[1] +
1927
v2_roe * normal_direction[2])
1928
1929
return vel_out_roe, c_f_roe
1930
end
1931
1932
# Calculate thermodynamic entropy for a conservative state `cons`
1933
@inline function entropy_thermodynamic(cons, equations::IdealGlmMhdEquations2D)
1934
# Pressure
1935
p = (equations.gamma - 1) *
1936
(cons[5] - 0.5f0 * (cons[2]^2 + cons[3]^2 + cons[4]^2) / cons[1]
1937
-
1938
0.5f0 * (cons[6]^2 + cons[7]^2 + cons[8]^2)
1939
-
1940
0.5f0 * cons[9]^2)
1941
1942
# Thermodynamic entropy
1943
s = log(p) - equations.gamma * log(cons[1])
1944
1945
return s
1946
end
1947
1948
# Calculate mathematical entropy for a conservative state `cons`
1949
@inline function entropy_math(cons, equations::IdealGlmMhdEquations2D)
1950
S = -entropy_thermodynamic(cons, equations) * cons[1] *
1951
equations.inv_gamma_minus_one
1952
1953
return S
1954
end
1955
1956
# Default entropy is the mathematical entropy
1957
@inline entropy(cons, equations::IdealGlmMhdEquations2D) = entropy_math(cons, equations)
1958
1959
# Calculate total energy for a conservative state `cons`
1960
@inline energy_total(cons, ::IdealGlmMhdEquations2D) = cons[5]
1961
1962
# Calculate kinetic energy for a conservative state `cons`
1963
@inline function energy_kinetic(cons, equations::IdealGlmMhdEquations2D)
1964
return 0.5f0 * (cons[2]^2 + cons[3]^2 + cons[4]^2) / cons[1]
1965
end
1966
1967
# Calculate the magnetic energy for a conservative state `cons`.
1968
# OBS! For non-dinmensional form of the ideal MHD magnetic pressure ≡ magnetic energy
1969
@inline function energy_magnetic(cons, ::IdealGlmMhdEquations2D)
1970
return 0.5f0 * (cons[6]^2 + cons[7]^2 + cons[8]^2)
1971
end
1972
1973
# Calculate internal energy for a conservative state `cons`
1974
@inline function energy_internal(cons, equations::IdealGlmMhdEquations2D)
1975
return (energy_total(cons, equations)
1976
-
1977
energy_kinetic(cons, equations)
1978
-
1979
energy_magnetic(cons, equations)
1980
-
1981
cons[9]^2 / 2)
1982
end
1983
1984
# State validation for Newton-bisection method of subcell IDP limiting
1985
@inline function Base.isvalid(u, equations::IdealGlmMhdEquations2D)
1986
if u[1] <= 0 || pressure(u, equations) <= 0
1987
return false
1988
end
1989
return true
1990
end
1991
1992
# Calculate the cross helicity (\vec{v}⋅\vec{B}) for a conservative state `cons`
1993
@inline function cross_helicity(cons, ::IdealGlmMhdEquations2D)
1994
return (cons[2] * cons[6] + cons[3] * cons[7] + cons[4] * cons[8]) / cons[1]
1995
end
1996
end # @muladd
1997
1998