/
ideal_glm_mhd_2d.jl
1414 lines (1246 loc) · 58.1 KB
/
ideal_glm_mhd_2d.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).
# Since these FMAs can increase the performance of many numerical algorithms,
# we need to opt-in explicitly.
# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.
@muladd begin
#! format: noindent
@doc raw"""
IdealGlmMhdEquations2D(gamma)
The ideal compressible GLM-MHD equations for an ideal gas with ratio of
specific heats `gamma` in two space dimensions.
"""
mutable struct IdealGlmMhdEquations2D{RealT <: Real} <:
AbstractIdealGlmMhdEquations{2, 9}
gamma::RealT # ratio of specific heats
inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications
c_h::RealT # GLM cleaning speed
function IdealGlmMhdEquations2D(gamma, c_h)
γ, inv_gamma_minus_one, c_h = promote(gamma, inv(gamma - 1), c_h)
new{typeof(γ)}(γ, inv_gamma_minus_one, c_h)
end
end
function IdealGlmMhdEquations2D(gamma; initial_c_h = convert(typeof(gamma), NaN))
# Use `promote` to ensure that `gamma` and `initial_c_h` have the same type
IdealGlmMhdEquations2D(promote(gamma, initial_c_h)...)
end
have_nonconservative_terms(::IdealGlmMhdEquations2D) = True()
n_nonconservative_terms(::IdealGlmMhdEquations2D) = 2
function varnames(::typeof(cons2cons), ::IdealGlmMhdEquations2D)
("rho", "rho_v1", "rho_v2", "rho_v3", "rho_e", "B1", "B2", "B3", "psi")
end
function varnames(::typeof(cons2prim), ::IdealGlmMhdEquations2D)
("rho", "v1", "v2", "v3", "p", "B1", "B2", "B3", "psi")
end
function default_analysis_integrals(::IdealGlmMhdEquations2D)
(entropy_timederivative, Val(:l2_divb), Val(:linf_divb))
end
# Set initial conditions at physical location `x` for time `t`
"""
initial_condition_constant(x, t, equations::IdealGlmMhdEquations2D)
A constant initial condition to test free-stream preservation.
"""
function initial_condition_constant(x, t, equations::IdealGlmMhdEquations2D)
rho = 1.0
rho_v1 = 0.1
rho_v2 = -0.2
rho_v3 = -0.5
rho_e = 50.0
B1 = 3.0
B2 = -1.2
B3 = 0.5
psi = 0.0
return SVector(rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi)
end
"""
initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations2D)
An Alfvén wave as smooth initial condition used for convergence tests.
"""
function initial_condition_convergence_test(x, t, equations::IdealGlmMhdEquations2D)
# smooth Alfvén wave test from Derigs et al. FLASH (2016)
# domain must be set to [0, 1/cos(α)] x [0, 1/sin(α)], γ = 5/3
alpha = 0.25 * pi
x_perp = x[1] * cos(alpha) + x[2] * sin(alpha)
B_perp = 0.1 * sin(2.0 * pi * x_perp)
rho = 1.0
v1 = -B_perp * sin(alpha)
v2 = B_perp * cos(alpha)
v3 = 0.1 * cos(2.0 * pi * x_perp)
p = 0.1
B1 = cos(alpha) + v1
B2 = sin(alpha) + v2
B3 = v3
psi = 0.0
return prim2cons(SVector(rho, v1, v2, v3, p, B1, B2, B3, psi), equations)
end
"""
initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations2D)
A weak blast wave adapted from
- Sebastian Hennemann, Gregor J. Gassner (2020)
A provably entropy stable subcell shock capturing approach for high order split form DG
[arXiv: 2008.12044](https://arxiv.org/abs/2008.12044)
"""
function initial_condition_weak_blast_wave(x, t, equations::IdealGlmMhdEquations2D)
# Adapted MHD version of the weak blast wave from Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
# Same discontinuity in the velocities but with magnetic fields
# Set up polar coordinates
inicenter = (0, 0)
x_norm = x[1] - inicenter[1]
y_norm = x[2] - inicenter[2]
r = sqrt(x_norm^2 + y_norm^2)
phi = atan(y_norm, x_norm)
# Calculate primitive variables
rho = r > 0.5 ? 1.0 : 1.1691
v1 = r > 0.5 ? 0.0 : 0.1882 * cos(phi)
v2 = r > 0.5 ? 0.0 : 0.1882 * sin(phi)
p = r > 0.5 ? 1.0 : 1.245
return prim2cons(SVector(rho, v1, v2, 0.0, p, 1.0, 1.0, 1.0, 0.0), equations)
end
# Pre-defined source terms should be implemented as
# function source_terms_WHATEVER(u, x, t, equations::IdealGlmMhdEquations2D)
# Calculate 1D flux in for a single point
@inline function flux(u, orientation::Integer, equations::IdealGlmMhdEquations2D)
rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
v1 = rho_v1 / rho
v2 = rho_v2 / rho
v3 = rho_v3 / rho
kin_en = 0.5 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
mag_en = 0.5 * (B1 * B1 + B2 * B2 + B3 * B3)
p_over_gamma_minus_one = (rho_e - kin_en - mag_en - 0.5 * psi^2)
p = (equations.gamma - 1) * p_over_gamma_minus_one
if orientation == 1
f1 = rho_v1
f2 = rho_v1 * v1 + p + mag_en - B1^2
f3 = rho_v1 * v2 - B1 * B2
f4 = rho_v1 * v3 - B1 * B3
f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v1 -
B1 * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B1
f6 = equations.c_h * psi
f7 = v1 * B2 - v2 * B1
f8 = v1 * B3 - v3 * B1
f9 = equations.c_h * B1
else #if orientation == 2
f1 = rho_v2
f2 = rho_v2 * v1 - B2 * B1
f3 = rho_v2 * v2 + p + mag_en - B2^2
f4 = rho_v2 * v3 - B2 * B3
f5 = (kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v2 -
B2 * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B2
f6 = v2 * B1 - v1 * B2
f7 = equations.c_h * psi
f8 = v2 * B3 - v3 * B2
f9 = equations.c_h * B2
end
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
end
# Calculate 1D flux for a single point in the normal direction
# Note, this directional vector is not normalized
@inline function flux(u, normal_direction::AbstractVector,
equations::IdealGlmMhdEquations2D)
rho, rho_v1, rho_v2, rho_v3, rho_e, B1, B2, B3, psi = u
v1 = rho_v1 / rho
v2 = rho_v2 / rho
v3 = rho_v3 / rho
kin_en = 0.5 * (rho_v1 * v1 + rho_v2 * v2 + rho_v3 * v3)
mag_en = 0.5 * (B1 * B1 + B2 * B2 + B3 * B3)
p_over_gamma_minus_one = (rho_e - kin_en - mag_en - 0.5 * psi^2)
p = (equations.gamma - 1) * p_over_gamma_minus_one
v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]
B_normal = B1 * normal_direction[1] + B2 * normal_direction[2]
rho_v_normal = rho * v_normal
f1 = rho_v_normal
f2 = rho_v_normal * v1 - B1 * B_normal + (p + mag_en) * normal_direction[1]
f3 = rho_v_normal * v2 - B2 * B_normal + (p + mag_en) * normal_direction[2]
f4 = rho_v_normal * v3 - B3 * B_normal
f5 = ((kin_en + equations.gamma * p_over_gamma_minus_one + 2 * mag_en) * v_normal
-
B_normal * (v1 * B1 + v2 * B2 + v3 * B3) + equations.c_h * psi * B_normal)
f6 = equations.c_h * psi * normal_direction[1] +
(v2 * B1 - v1 * B2) * normal_direction[2]
f7 = equations.c_h * psi * normal_direction[2] +
(v1 * B2 - v2 * B1) * normal_direction[1]
f8 = v_normal * B3 - v3 * B_normal
f9 = equations.c_h * B_normal
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
end
"""
flux_nonconservative_powell(u_ll, u_rr, orientation::Integer,
equations::IdealGlmMhdEquations2D)
flux_nonconservative_powell(u_ll, u_rr,
normal_direction_ll ::AbstractVector,
normal_direction_average::AbstractVector,
equations::IdealGlmMhdEquations2D)
Non-symmetric two-point flux discretizing the nonconservative (source) term of
Powell and the Galilean nonconservative term associated with the GLM multiplier
of the [`IdealGlmMhdEquations2D`](@ref).
On curvilinear meshes, this nonconservative flux depends on both the
contravariant vector (normal direction) at the current node and the averaged
one. This is different from numerical fluxes used to discretize conservative
terms.
## References
- Marvin Bohm, Andrew R.Winters, Gregor J. Gassner, Dominik Derigs,
Florian Hindenlang, Joachim Saur
An entropy stable nodal discontinuous Galerkin method for the resistive MHD
equations. Part I: Theory and numerical verification
[DOI: 10.1016/j.jcp.2018.06.027](https://doi.org/10.1016/j.jcp.2018.06.027)
"""
@inline function flux_nonconservative_powell(u_ll, u_rr, orientation::Integer,
equations::IdealGlmMhdEquations2D)
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v3_ll = rho_v3_ll / rho_ll
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
if orientation == 1
f = SVector(0,
B1_ll * B1_rr,
B2_ll * B1_rr,
B3_ll * B1_rr,
v_dot_B_ll * B1_rr + v1_ll * psi_ll * psi_rr,
v1_ll * B1_rr,
v2_ll * B1_rr,
v3_ll * B1_rr,
v1_ll * psi_rr)
else # orientation == 2
f = SVector(0,
B1_ll * B2_rr,
B2_ll * B2_rr,
B3_ll * B2_rr,
v_dot_B_ll * B2_rr + v2_ll * psi_ll * psi_rr,
v1_ll * B2_rr,
v2_ll * B2_rr,
v3_ll * B2_rr,
v2_ll * psi_rr)
end
return f
end
@inline function flux_nonconservative_powell(u_ll, u_rr,
normal_direction_ll::AbstractVector,
normal_direction_average::AbstractVector,
equations::IdealGlmMhdEquations2D)
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v3_ll = rho_v3_ll / rho_ll
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
# Note that `v_dot_n_ll` uses the `normal_direction_ll` (contravariant vector
# at the same node location) while `B_dot_n_rr` uses the averaged normal
# direction. The reason for this is that `v_dot_n_ll` depends only on the left
# state and multiplies some gradient while `B_dot_n_rr` is used to compute
# the divergence of B.
v_dot_n_ll = v1_ll * normal_direction_ll[1] + v2_ll * normal_direction_ll[2]
B_dot_n_rr = B1_rr * normal_direction_average[1] +
B2_rr * normal_direction_average[2]
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
f = SVector(0,
B1_ll * B_dot_n_rr,
B2_ll * B_dot_n_rr,
B3_ll * B_dot_n_rr,
v_dot_B_ll * B_dot_n_rr + v_dot_n_ll * psi_ll * psi_rr,
v1_ll * B_dot_n_rr,
v2_ll * B_dot_n_rr,
v3_ll * B_dot_n_rr,
v_dot_n_ll * psi_rr)
return f
end
"""
flux_nonconservative_powell_local_symmetric(u_ll, u_rr,
orientation::Integer,
equations::IdealGlmMhdEquations2D)
Non-symmetric two-point flux discretizing the nonconservative (source) term of
Powell and the Galilean nonconservative term associated with the GLM multiplier
of the [`IdealGlmMhdEquations2D`](@ref).
This implementation uses a non-conservative term that can be written as the product
of local and symmetric parts. It is equivalent to the non-conservative flux of Bohm
et al. (`flux_nonconservative_powell`) for conforming meshes but it yields different
results on non-conforming meshes(!).
The two other flux functions with the same name return either the local
or symmetric portion of the non-conservative flux based on the type of the
nonconservative_type argument, employing multiple dispatch. They are used to
compute the subcell fluxes in dg_2d_subcell_limiters.jl.
## References
- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts
Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.
"""
@inline function flux_nonconservative_powell_local_symmetric(u_ll, u_rr,
orientation::Integer,
equations::IdealGlmMhdEquations2D)
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v3_ll = rho_v3_ll / rho_ll
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
psi_avg = (psi_ll + psi_rr) #* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
if orientation == 1
B1_avg = (B1_ll + B1_rr) #* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
f = SVector(0,
B1_ll * B1_avg,
B2_ll * B1_avg,
B3_ll * B1_avg,
v_dot_B_ll * B1_avg + v1_ll * psi_ll * psi_avg,
v1_ll * B1_avg,
v2_ll * B1_avg,
v3_ll * B1_avg,
v1_ll * psi_avg)
else # orientation == 2
B2_avg = (B2_ll + B2_rr) #* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
f = SVector(0,
B1_ll * B2_avg,
B2_ll * B2_avg,
B3_ll * B2_avg,
v_dot_B_ll * B2_avg + v2_ll * psi_ll * psi_avg,
v1_ll * B2_avg,
v2_ll * B2_avg,
v3_ll * B2_avg,
v2_ll * psi_avg)
end
return f
end
"""
flux_nonconservative_powell_local_symmetric(u_ll, orientation::Integer,
equations::IdealGlmMhdEquations2D,
nonconservative_type::NonConservativeLocal,
nonconservative_term::Integer)
Local part of the Powell and GLM non-conservative terms. Needed for the calculation of
the non-conservative staggered "fluxes" for subcell limiting. See, e.g.,
- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts
Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.
This function is used to compute the subcell fluxes in dg_2d_subcell_limiters.jl.
"""
@inline function flux_nonconservative_powell_local_symmetric(u_ll, orientation::Integer,
equations::IdealGlmMhdEquations2D,
nonconservative_type::NonConservativeLocal,
nonconservative_term::Integer)
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
if nonconservative_term == 1
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v3_ll = rho_v3_ll / rho_ll
v_dot_B_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
f = SVector(0,
B1_ll,
B2_ll,
B3_ll,
v_dot_B_ll,
v1_ll,
v2_ll,
v3_ll,
0)
else #nonconservative_term ==2
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
if orientation == 1
v1_ll = rho_v1_ll / rho_ll
f = SVector(0,
0,
0,
0,
v1_ll * psi_ll,
0,
0,
0,
v1_ll)
else #orientation == 2
v2_ll = rho_v2_ll / rho_ll
f = SVector(0,
0,
0,
0,
v2_ll * psi_ll,
0,
0,
0,
v2_ll)
end
end
return f
end
"""
flux_nonconservative_powell_local_symmetric(u_ll, orientation::Integer,
equations::IdealGlmMhdEquations2D,
nonconservative_type::NonConservativeSymmetric,
nonconservative_term::Integer)
Symmetric part of the Powell and GLM non-conservative terms. Needed for the calculation of
the non-conservative staggered "fluxes" for subcell limiting. See, e.g.,
- Rueda-Ramírez, Gassner (2023). A Flux-Differencing Formula for Split-Form Summation By Parts
Discretizations of Non-Conservative Systems. https://arxiv.org/pdf/2211.14009.pdf.
This function is used to compute the subcell fluxes in dg_2d_subcell_limiters.jl.
"""
@inline function flux_nonconservative_powell_local_symmetric(u_ll, u_rr,
orientation::Integer,
equations::IdealGlmMhdEquations2D,
nonconservative_type::NonConservativeSymmetric,
nonconservative_term::Integer)
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
if nonconservative_term == 1
# Powell nonconservative term: (0, B_1, B_2, B_3, v⋅B, v_1, v_2, v_3, 0)
if orientation == 1
B1_avg = (B1_ll + B1_rr)#* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
f = SVector(0,
B1_avg,
B1_avg,
B1_avg,
B1_avg,
B1_avg,
B1_avg,
B1_avg,
0)
else # orientation == 2
B2_avg = (B2_ll + B2_rr)#* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
f = SVector(0,
B2_avg,
B2_avg,
B2_avg,
B2_avg,
B2_avg,
B2_avg,
B2_avg,
0)
end
else #nonconservative_term == 2
# Galilean nonconservative term: (0, 0, 0, 0, ψ v_{1,2}, 0, 0, 0, v_{1,2})
psi_avg = (psi_ll + psi_rr)#* 0.5 # The flux is already multiplied by 0.5 wherever it is used in the code
f = SVector(0,
0,
0,
0,
psi_avg,
0,
0,
0,
psi_avg)
end
return f
end
"""
flux_derigs_etal(u_ll, u_rr, orientation, equations::IdealGlmMhdEquations2D)
Entropy conserving two-point flux by
- Derigs et al. (2018)
Ideal GLM-MHD: About the entropy consistent nine-wave magnetic field
divergence diminishing ideal magnetohydrodynamics equations
[DOI: 10.1016/j.jcp.2018.03.002](https://doi.org/10.1016/j.jcp.2018.03.002)
"""
function flux_derigs_etal(u_ll, u_rr, orientation::Integer,
equations::IdealGlmMhdEquations2D)
# Unpack left and right states to get velocities, pressure, and inverse temperature (called beta)
rho_ll, rho_v1_ll, rho_v2_ll, rho_v3_ll, rho_e_ll, B1_ll, B2_ll, B3_ll, psi_ll = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, rho_v3_rr, rho_e_rr, B1_rr, B2_rr, B3_rr, psi_rr = u_rr
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v3_ll = rho_v3_ll / rho_ll
v1_rr = rho_v1_rr / rho_rr
v2_rr = rho_v2_rr / rho_rr
v3_rr = rho_v3_rr / rho_rr
vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2
vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2
mag_norm_ll = B1_ll^2 + B2_ll^2 + B3_ll^2
mag_norm_rr = B1_rr^2 + B2_rr^2 + B3_rr^2
p_ll = (equations.gamma - 1) *
(rho_e_ll - 0.5 * rho_ll * vel_norm_ll - 0.5 * mag_norm_ll - 0.5 * psi_ll^2)
p_rr = (equations.gamma - 1) *
(rho_e_rr - 0.5 * rho_rr * vel_norm_rr - 0.5 * mag_norm_rr - 0.5 * psi_rr^2)
beta_ll = 0.5 * rho_ll / p_ll
beta_rr = 0.5 * rho_rr / p_rr
# for convenience store v⋅B
vel_dot_mag_ll = v1_ll * B1_ll + v2_ll * B2_ll + v3_ll * B3_ll
vel_dot_mag_rr = v1_rr * B1_rr + v2_rr * B2_rr + v3_rr * B3_rr
# Compute the necessary mean values needed for either direction
rho_avg = 0.5 * (rho_ll + rho_rr)
rho_mean = ln_mean(rho_ll, rho_rr)
beta_mean = ln_mean(beta_ll, beta_rr)
beta_avg = 0.5 * (beta_ll + beta_rr)
v1_avg = 0.5 * (v1_ll + v1_rr)
v2_avg = 0.5 * (v2_ll + v2_rr)
v3_avg = 0.5 * (v3_ll + v3_rr)
p_mean = 0.5 * rho_avg / beta_avg
B1_avg = 0.5 * (B1_ll + B1_rr)
B2_avg = 0.5 * (B2_ll + B2_rr)
B3_avg = 0.5 * (B3_ll + B3_rr)
psi_avg = 0.5 * (psi_ll + psi_rr)
vel_norm_avg = 0.5 * (vel_norm_ll + vel_norm_rr)
mag_norm_avg = 0.5 * (mag_norm_ll + mag_norm_rr)
vel_dot_mag_avg = 0.5 * (vel_dot_mag_ll + vel_dot_mag_rr)
# Calculate fluxes depending on orientation with specific direction averages
if orientation == 1
f1 = rho_mean * v1_avg
f2 = f1 * v1_avg + p_mean + 0.5 * mag_norm_avg - B1_avg * B1_avg
f3 = f1 * v2_avg - B1_avg * B2_avg
f4 = f1 * v3_avg - B1_avg * B3_avg
f6 = equations.c_h * psi_avg
f7 = v1_avg * B2_avg - v2_avg * B1_avg
f8 = v1_avg * B3_avg - v3_avg * B1_avg
f9 = equations.c_h * B1_avg
# total energy flux is complicated and involves the previous eight components
psi_B1_avg = 0.5 * (B1_ll * psi_ll + B1_rr * psi_rr)
v1_mag_avg = 0.5 * (v1_ll * mag_norm_ll + v1_rr * mag_norm_rr)
f5 = (f1 * 0.5 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +
f2 * v1_avg + f3 * v2_avg +
f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg + f9 * psi_avg -
0.5 * v1_mag_avg +
B1_avg * vel_dot_mag_avg - equations.c_h * psi_B1_avg)
else
f1 = rho_mean * v2_avg
f2 = f1 * v1_avg - B1_avg * B2_avg
f3 = f1 * v2_avg + p_mean + 0.5 * mag_norm_avg - B2_avg * B2_avg
f4 = f1 * v3_avg - B2_avg * B3_avg
f6 = v2_avg * B1_avg - v1_avg * B2_avg
f7 = equations.c_h * psi_avg
f8 = v2_avg * B3_avg - v3_avg * B2_avg
f9 = equations.c_h * B2_avg
# total energy flux is complicated and involves the previous eight components
psi_B2_avg = 0.5 * (B2_ll * psi_ll + B2_rr * psi_rr)
v2_mag_avg = 0.5 * (v2_ll * mag_norm_ll + v2_rr * mag_norm_rr)
f5 = (f1 * 0.5 * (1 / (equations.gamma - 1) / beta_mean - vel_norm_avg) +
f2 * v1_avg + f3 * v2_avg +
f4 * v3_avg + f6 * B1_avg + f7 * B2_avg + f8 * B3_avg + f9 * psi_avg -
0.5 * v2_mag_avg +
B2_avg * vel_dot_mag_avg - equations.c_h * psi_B2_avg)
end
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
end
"""
flux_hindenlang_gassner(u_ll, u_rr, orientation_or_normal_direction,
equations::IdealGlmMhdEquations2D)
Entropy conserving and kinetic energy preserving two-point flux of
Hindenlang and Gassner (2019), extending [`flux_ranocha`](@ref) to the MHD equations.
## References
- Florian Hindenlang, Gregor Gassner (2019)
A new entropy conservative two-point flux for ideal MHD equations derived from
first principles.
Presented at HONOM 2019: European workshop on high order numerical methods
for evolutionary PDEs, theory and applications
- Hendrik Ranocha (2018)
Generalised Summation-by-Parts Operators and Entropy Stability of Numerical Methods
for Hyperbolic Balance Laws
[PhD thesis, TU Braunschweig](https://cuvillier.de/en/shop/publications/7743)
- Hendrik Ranocha (2020)
Entropy Conserving and Kinetic Energy Preserving Numerical Methods for
the Euler Equations Using Summation-by-Parts Operators
[Proceedings of ICOSAHOM 2018](https://doi.org/10.1007/978-3-030-39647-3_42)
"""
@inline function flux_hindenlang_gassner(u_ll, u_rr, orientation::Integer,
equations::IdealGlmMhdEquations2D)
# Unpack left and right states
rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll, psi_ll = cons2prim(u_ll,
equations)
rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr, psi_rr = cons2prim(u_rr,
equations)
# Compute the necessary mean values needed for either direction
rho_mean = ln_mean(rho_ll, rho_rr)
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
# in exact arithmetic since
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
v1_avg = 0.5 * (v1_ll + v1_rr)
v2_avg = 0.5 * (v2_ll + v2_rr)
v3_avg = 0.5 * (v3_ll + v3_rr)
p_avg = 0.5 * (p_ll + p_rr)
psi_avg = 0.5 * (psi_ll + psi_rr)
velocity_square_avg = 0.5 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
magnetic_square_avg = 0.5 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)
# Calculate fluxes depending on orientation with specific direction averages
if orientation == 1
f1 = rho_mean * v1_avg
f2 = f1 * v1_avg + p_avg + magnetic_square_avg -
0.5 * (B1_ll * B1_rr + B1_rr * B1_ll)
f3 = f1 * v2_avg - 0.5 * (B1_ll * B2_rr + B1_rr * B2_ll)
f4 = f1 * v3_avg - 0.5 * (B1_ll * B3_rr + B1_rr * B3_ll)
#f5 below
f6 = equations.c_h * psi_avg
f7 = 0.5 * (v1_ll * B2_ll - v2_ll * B1_ll + v1_rr * B2_rr - v2_rr * B1_rr)
f8 = 0.5 * (v1_ll * B3_ll - v3_ll * B1_ll + v1_rr * B3_rr - v3_rr * B1_rr)
f9 = equations.c_h * 0.5 * (B1_ll + B1_rr)
# total energy flux is complicated and involves the previous components
f5 = (f1 *
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
+
0.5 * (+p_ll * v1_rr + p_rr * v1_ll
+ (v1_ll * B2_ll * B2_rr + v1_rr * B2_rr * B2_ll)
+ (v1_ll * B3_ll * B3_rr + v1_rr * B3_rr * B3_ll)
-
(v2_ll * B1_ll * B2_rr + v2_rr * B1_rr * B2_ll)
-
(v3_ll * B1_ll * B3_rr + v3_rr * B1_rr * B3_ll)
+
equations.c_h * (B1_ll * psi_rr + B1_rr * psi_ll)))
else # orientation == 2
f1 = rho_mean * v2_avg
f2 = f1 * v1_avg - 0.5 * (B2_ll * B1_rr + B2_rr * B1_ll)
f3 = f1 * v2_avg + p_avg + magnetic_square_avg -
0.5 * (B2_ll * B2_rr + B2_rr * B2_ll)
f4 = f1 * v3_avg - 0.5 * (B2_ll * B3_rr + B2_rr * B3_ll)
#f5 below
f6 = 0.5 * (v2_ll * B1_ll - v1_ll * B2_ll + v2_rr * B1_rr - v1_rr * B2_rr)
f7 = equations.c_h * psi_avg
f8 = 0.5 * (v2_ll * B3_ll - v3_ll * B2_ll + v2_rr * B3_rr - v3_rr * B2_rr)
f9 = equations.c_h * 0.5 * (B2_ll + B2_rr)
# total energy flux is complicated and involves the previous components
f5 = (f1 *
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
+
0.5 * (+p_ll * v2_rr + p_rr * v2_ll
+ (v2_ll * B1_ll * B1_rr + v2_rr * B1_rr * B1_ll)
+ (v2_ll * B3_ll * B3_rr + v2_rr * B3_rr * B3_ll)
-
(v1_ll * B2_ll * B1_rr + v1_rr * B2_rr * B1_ll)
-
(v3_ll * B2_ll * B3_rr + v3_rr * B2_rr * B3_ll)
+
equations.c_h * (B2_ll * psi_rr + B2_rr * psi_ll)))
end
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
end
@inline function flux_hindenlang_gassner(u_ll, u_rr, normal_direction::AbstractVector,
equations::IdealGlmMhdEquations2D)
# Unpack left and right states
rho_ll, v1_ll, v2_ll, v3_ll, p_ll, B1_ll, B2_ll, B3_ll, psi_ll = cons2prim(u_ll,
equations)
rho_rr, v1_rr, v2_rr, v3_rr, p_rr, B1_rr, B2_rr, B3_rr, psi_rr = cons2prim(u_rr,
equations)
v_dot_n_ll = v1_ll * normal_direction[1] + v2_ll * normal_direction[2]
v_dot_n_rr = v1_rr * normal_direction[1] + v2_rr * normal_direction[2]
B_dot_n_ll = B1_ll * normal_direction[1] + B2_ll * normal_direction[2]
B_dot_n_rr = B1_rr * normal_direction[1] + B2_rr * normal_direction[2]
# Compute the necessary mean values needed for either direction
rho_mean = ln_mean(rho_ll, rho_rr)
# Algebraically equivalent to `inv_ln_mean(rho_ll / p_ll, rho_rr / p_rr)`
# in exact arithmetic since
# log((ϱₗ/pₗ) / (ϱᵣ/pᵣ)) / (ϱₗ/pₗ - ϱᵣ/pᵣ)
# = pₗ pᵣ log((ϱₗ pᵣ) / (ϱᵣ pₗ)) / (ϱₗ pᵣ - ϱᵣ pₗ)
inv_rho_p_mean = p_ll * p_rr * inv_ln_mean(rho_ll * p_rr, rho_rr * p_ll)
v1_avg = 0.5 * (v1_ll + v1_rr)
v2_avg = 0.5 * (v2_ll + v2_rr)
v3_avg = 0.5 * (v3_ll + v3_rr)
p_avg = 0.5 * (p_ll + p_rr)
psi_avg = 0.5 * (psi_ll + psi_rr)
velocity_square_avg = 0.5 * (v1_ll * v1_rr + v2_ll * v2_rr + v3_ll * v3_rr)
magnetic_square_avg = 0.5 * (B1_ll * B1_rr + B2_ll * B2_rr + B3_ll * B3_rr)
# Calculate fluxes depending on normal_direction
f1 = rho_mean * 0.5 * (v_dot_n_ll + v_dot_n_rr)
f2 = (f1 * v1_avg + (p_avg + magnetic_square_avg) * normal_direction[1]
-
0.5 * (B_dot_n_ll * B1_rr + B_dot_n_rr * B1_ll))
f3 = (f1 * v2_avg + (p_avg + magnetic_square_avg) * normal_direction[2]
-
0.5 * (B_dot_n_ll * B2_rr + B_dot_n_rr * B2_ll))
f4 = (f1 * v3_avg
-
0.5 * (B_dot_n_ll * B3_rr + B_dot_n_rr * B3_ll))
#f5 below
f6 = (equations.c_h * psi_avg * normal_direction[1]
+
0.5 * (v_dot_n_ll * B1_ll - v1_ll * B_dot_n_ll +
v_dot_n_rr * B1_rr - v1_rr * B_dot_n_rr))
f7 = (equations.c_h * psi_avg * normal_direction[2]
+
0.5 * (v_dot_n_ll * B2_ll - v2_ll * B_dot_n_ll +
v_dot_n_rr * B2_rr - v2_rr * B_dot_n_rr))
f8 = +0.5 * (v_dot_n_ll * B3_ll - v3_ll * B_dot_n_ll +
v_dot_n_rr * B3_rr - v3_rr * B_dot_n_rr)
f9 = equations.c_h * 0.5 * (B_dot_n_ll + B_dot_n_rr)
# total energy flux is complicated and involves the previous components
f5 = (f1 * (velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one)
+
0.5 * (+p_ll * v_dot_n_rr + p_rr * v_dot_n_ll
+ (v_dot_n_ll * B1_ll * B1_rr + v_dot_n_rr * B1_rr * B1_ll)
+ (v_dot_n_ll * B2_ll * B2_rr + v_dot_n_rr * B2_rr * B2_ll)
+ (v_dot_n_ll * B3_ll * B3_rr + v_dot_n_rr * B3_rr * B3_ll)
-
(v1_ll * B_dot_n_ll * B1_rr + v1_rr * B_dot_n_rr * B1_ll)
-
(v2_ll * B_dot_n_ll * B2_rr + v2_rr * B_dot_n_rr * B2_ll)
-
(v3_ll * B_dot_n_ll * B3_rr + v3_rr * B_dot_n_rr * B3_ll)
+
equations.c_h * (B_dot_n_ll * psi_rr + B_dot_n_rr * psi_ll)))
return SVector(f1, f2, f3, f4, f5, f6, f7, f8, f9)
end
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
equations::IdealGlmMhdEquations2D)
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
# Calculate the left/right velocities and fast magnetoacoustic wave speeds
if orientation == 1
v_ll = rho_v1_ll / rho_ll
v_rr = rho_v1_rr / rho_rr
else # orientation == 2
v_ll = rho_v2_ll / rho_ll
v_rr = rho_v2_rr / rho_rr
end
cf_ll = calc_fast_wavespeed(u_ll, orientation, equations)
cf_rr = calc_fast_wavespeed(u_rr, orientation, equations)
return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)
end
@inline function max_abs_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
equations::IdealGlmMhdEquations2D)
# return max(v_mag_ll, v_mag_rr) + max(cf_ll, cf_rr)
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
# Calculate normal velocities and fast magnetoacoustic wave speeds
# left
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v_ll = (v1_ll * normal_direction[1]
+
v2_ll * normal_direction[2])
cf_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
# right
v1_rr = rho_v1_rr / rho_rr
v2_rr = rho_v2_rr / rho_rr
v_rr = (v1_rr * normal_direction[1]
+
v2_rr * normal_direction[2])
cf_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
# wave speeds already scaled by norm(normal_direction) in [`calc_fast_wavespeed`](@ref)
return max(abs(v_ll), abs(v_rr)) + max(cf_ll, cf_rr)
end
# Calculate estimate for minimum and maximum wave speeds for HLL-type fluxes
@inline function min_max_speed_naive(u_ll, u_rr, orientation::Integer,
equations::IdealGlmMhdEquations2D)
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
# Calculate primitive velocity variables
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v1_rr = rho_v1_rr / rho_rr
v2_rr = rho_v2_rr / rho_rr
# Approximate the left-most and right-most eigenvalues in the Riemann fan
if orientation == 1 # x-direction
λ_min = v1_ll - calc_fast_wavespeed(u_ll, orientation, equations)
λ_max = v1_rr + calc_fast_wavespeed(u_rr, orientation, equations)
else # y-direction
λ_min = v2_ll - calc_fast_wavespeed(u_ll, orientation, equations)
λ_max = v2_rr + calc_fast_wavespeed(u_rr, orientation, equations)
end
return λ_min, λ_max
end
@inline function min_max_speed_naive(u_ll, u_rr, normal_direction::AbstractVector,
equations::IdealGlmMhdEquations2D)
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
# Calculate primitive velocity variables
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v1_rr = rho_v1_rr / rho_rr
v2_rr = rho_v2_rr / rho_rr
v_normal_ll = (v1_ll * normal_direction[1] + v2_ll * normal_direction[2])
v_normal_rr = (v1_rr * normal_direction[1] + v2_rr * normal_direction[2])
c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
# Estimate the min/max eigenvalues in the normal direction
λ_min = min(v_normal_ll - c_f_ll, v_normal_rr - c_f_rr)
λ_max = max(v_normal_rr + c_f_rr, v_normal_rr + c_f_rr)
return λ_min, λ_max
end
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
@inline function min_max_speed_davis(u_ll, u_rr, orientation::Integer,
equations::IdealGlmMhdEquations2D)
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
# Calculate primitive velocity variables
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v1_rr = rho_v1_rr / rho_rr
v2_rr = rho_v2_rr / rho_rr
# Approximate the left-most and right-most eigenvalues in the Riemann fan
if orientation == 1 # x-direction
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
λ_min = min(v1_ll - c_f_ll, v1_rr - c_f_rr)
λ_max = max(v1_ll + c_f_ll, v1_rr + c_f_rr)
else # y-direction
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
λ_min = min(v2_ll - c_f_ll, v2_rr - c_f_rr)
λ_max = max(v2_ll + c_f_ll, v1_rr + c_f_rr)
end
return λ_min, λ_max
end
# More refined estimates for minimum and maximum wave speeds for HLL-type fluxes
@inline function min_max_speed_davis(u_ll, u_rr, normal_direction::AbstractVector,
equations::IdealGlmMhdEquations2D)
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
# Calculate primitive velocity variables
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v1_rr = rho_v1_rr / rho_rr
v2_rr = rho_v2_rr / rho_rr
v_normal_ll = (v1_ll * normal_direction[1] + v2_ll * normal_direction[2])
v_normal_rr = (v1_rr * normal_direction[1] + v2_rr * normal_direction[2])
c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
# Estimate the min/max eigenvalues in the normal direction
λ_min = min(v_normal_ll - c_f_ll, v_normal_rr - c_f_rr)
λ_max = max(v_normal_ll + c_f_ll, v_normal_rr + c_f_rr)
return λ_min, λ_max
end
"""
min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer, equations::IdealGlmMhdEquations2D)
Calculate minimum and maximum wave speeds for HLL-type fluxes as in
- Li (2005)
An HLLC Riemann solver for magneto-hydrodynamics
[DOI: 10.1016/j.jcp.2004.08.020](https://doi.org/10.1016/j.jcp.2004.08.020).
"""
@inline function min_max_speed_einfeldt(u_ll, u_rr, orientation::Integer,
equations::IdealGlmMhdEquations2D)
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
# Calculate primitive velocity variables
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v1_rr = rho_v1_rr / rho_rr
v2_rr = rho_v2_rr / rho_rr
# Approximate the left-most and right-most eigenvalues in the Riemann fan
if orientation == 1 # x-direction
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)
λ_min = min(v1_ll - c_f_ll, vel_roe - c_f_roe)
λ_max = max(v1_rr + c_f_rr, vel_roe + c_f_roe)
else # y-direction
c_f_ll = calc_fast_wavespeed(u_ll, orientation, equations)
c_f_rr = calc_fast_wavespeed(u_rr, orientation, equations)
vel_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, orientation, equations)
λ_min = min(v2_ll - c_f_ll, vel_roe - c_f_roe)
λ_max = max(v2_rr + c_f_rr, vel_roe + c_f_roe)
end
return λ_min, λ_max
end
@inline function min_max_speed_einfeldt(u_ll, u_rr, normal_direction::AbstractVector,
equations::IdealGlmMhdEquations2D)
rho_ll, rho_v1_ll, rho_v2_ll, _ = u_ll
rho_rr, rho_v1_rr, rho_v2_rr, _ = u_rr
# Calculate primitive velocity variables
v1_ll = rho_v1_ll / rho_ll
v2_ll = rho_v2_ll / rho_ll
v1_rr = rho_v1_rr / rho_rr
v2_rr = rho_v2_rr / rho_rr
v_normal_ll = (v1_ll * normal_direction[1] + v2_ll * normal_direction[2])
v_normal_rr = (v1_rr * normal_direction[1] + v2_rr * normal_direction[2])
c_f_ll = calc_fast_wavespeed(u_ll, normal_direction, equations)
c_f_rr = calc_fast_wavespeed(u_rr, normal_direction, equations)
v_roe, c_f_roe = calc_fast_wavespeed_roe(u_ll, u_rr, normal_direction, equations)
# Estimate the min/max eigenvalues in the normal direction
λ_min = min(v_normal_ll - c_f_ll, v_roe - c_f_roe)
λ_max = max(v_normal_rr + c_f_rr, v_roe + c_f_roe)
return λ_min, λ_max
end
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction
# has been normalized prior to this rotation of the state vector
@inline function rotate_to_x(u, normal_vector, equations::IdealGlmMhdEquations2D)
# cos and sin of the angle between the x-axis and the normalized normal_vector are
# the normalized vector's x and y coordinates respectively (see unit circle).
c = normal_vector[1]
s = normal_vector[2]
# Apply the 2D rotation matrix with normal and tangent directions of the form
# [ 1 0 0 0 0 0 0 0 0;
# 0 n_1 n_2 0 0 0 0 0 0;
# 0 t_1 t_2 0 0 0 0 0 0;
# 0 0 0 1 0 0 0 0 0;
# 0 0 0 0 1 0 0 0 0;
# 0 0 0 0 0 n_1 n_2 0 0;
# 0 0 0 0 0 t_1 t_2 0 0;
# 0 0 0 0 0 0 0 1 0;
# 0 0 0 0 0 0 0 0 1 ]
# where t_1 = -n_2 and t_2 = n_1.
# Note for IdealGlmMhdEquations2D only the velocities and magnetic field variables rotate
return SVector(u[1],
c * u[2] + s * u[3],
-s * u[2] + c * u[3],
u[4],
u[5],
c * u[6] + s * u[7],
-s * u[6] + c * u[7],
u[8],
u[9])
end
# Called inside `FluxRotated` in `numerical_fluxes.jl` so the direction
# has been normalized prior to this back-rotation of the state vector
@inline function rotate_from_x(u, normal_vector, equations::IdealGlmMhdEquations2D)
# cos and sin of the angle between the x-axis and the normalized normal_vector are
# the normalized vector's x and y coordinates respectively (see unit circle).
c = normal_vector[1]
s = normal_vector[2]
# Apply the 2D back-rotation matrix with normal and tangent directions of the form
# [ 1 0 0 0 0 0 0 0 0;
# 0 n_1 t_1 0 0 0 0 0 0;
# 0 n_2 t_2 0 0 0 0 0 0;
# 0 0 0 1 0 0 0 0 0;
# 0 0 0 0 1 0 0 0 0;
# 0 0 0 0 0 n_1 t_1 0 0;
# 0 0 0 0 0 n_2 t_2 0 0;
# 0 0 0 0 0 0 0 1 0;
# 0 0 0 0 0 0 0 0 1 ]
# where t_1 = -n_2 and t_2 = n_1.
# Note for IdealGlmMhdEquations2D the velocities and magnetic field variables back-rotate