-
Notifications
You must be signed in to change notification settings - Fork 105
/
compressible_euler_2d.jl
1492 lines (1279 loc) · 56 KB
/
compressible_euler_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"""
CompressibleEulerEquations2D(gamma)
The compressible Euler equations
```math
\frac{\partial}{\partial t}
\begin{pmatrix}
\rho \\ \rho v_1 \\ \rho v_2 \\ \rho e
\end{pmatrix}
+
\frac{\partial}{\partial x}
\begin{pmatrix}
\rho v_1 \\ \rho v_1^2 + p \\ \rho v_1 v_2 \\ (\rho e +p) v_1
\end{pmatrix}
+
\frac{\partial}{\partial y}
\begin{pmatrix}
\rho v_2 \\ \rho v_1 v_2 \\ \rho v_2^2 + p \\ (\rho e +p) v_2
\end{pmatrix}
=
\begin{pmatrix}
0 \\ 0 \\ 0 \\ 0
\end{pmatrix}
```
for an ideal gas with ratio of specific heats `gamma`
in two space dimensions.
Here, ``\rho`` is the density, ``v_1``, ``v_2`` the velocities, ``e`` the specific total energy **rather than** specific internal energy, and
```math
p = (\gamma - 1) \left( \rho e - \frac{1}{2} \rho (v_1^2+v_2^2) \right)
```
the pressure.
"""
struct CompressibleEulerEquations2D{RealT <: Real} <:
AbstractCompressibleEulerEquations{2, 4}
gamma::RealT # ratio of specific heats
inv_gamma_minus_one::RealT # = inv(gamma - 1); can be used to write slow divisions as fast multiplications
function CompressibleEulerEquations2D(gamma)
γ, inv_gamma_minus_one = promote(gamma, inv(gamma - 1))
new{typeof(γ)}(γ, inv_gamma_minus_one)
end
end
function varnames(::typeof(cons2cons), ::CompressibleEulerEquations2D)
("rho", "rho_v1", "rho_v2", "rho_e")
end
varnames(::typeof(cons2prim), ::CompressibleEulerEquations2D) = ("rho", "v1", "v2", "p")
# Set initial conditions at physical location `x` for time `t`
"""
initial_condition_constant(x, t, equations::CompressibleEulerEquations2D)
A constant initial condition to test free-stream preservation.
"""
function initial_condition_constant(x, t, equations::CompressibleEulerEquations2D)
rho = 1.0
rho_v1 = 0.1
rho_v2 = -0.2
rho_e = 10.0
return SVector(rho, rho_v1, rho_v2, rho_e)
end
"""
initial_condition_convergence_test(x, t, equations::CompressibleEulerEquations2D)
A smooth initial condition used for convergence tests in combination with
[`source_terms_convergence_test`](@ref)
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
"""
function initial_condition_convergence_test(x, t,
equations::CompressibleEulerEquations2D)
c = 2
A = 0.1
L = 2
f = 1 / L
ω = 2 * pi * f
ini = c + A * sin(ω * (x[1] + x[2] - t))
rho = ini
rho_v1 = ini
rho_v2 = ini
rho_e = ini^2
return SVector(rho, rho_v1, rho_v2, rho_e)
end
"""
source_terms_convergence_test(u, x, t, equations::CompressibleEulerEquations2D)
Source terms used for convergence tests in combination with
[`initial_condition_convergence_test`](@ref)
(and [`BoundaryConditionDirichlet(initial_condition_convergence_test)`](@ref) in non-periodic domains).
"""
@inline function source_terms_convergence_test(u, x, t,
equations::CompressibleEulerEquations2D)
# Same settings as in `initial_condition`
c = 2
A = 0.1
L = 2
f = 1 / L
ω = 2 * pi * f
γ = equations.gamma
x1, x2 = x
si, co = sincos(ω * (x1 + x2 - t))
rho = c + A * si
rho_x = ω * A * co
# Note that d/dt rho = -d/dx rho = -d/dy rho.
tmp = (2 * rho - 1) * (γ - 1)
du1 = rho_x
du2 = rho_x * (1 + tmp)
du3 = du2
du4 = 2 * rho_x * (rho + tmp)
return SVector(du1, du2, du3, du4)
end
"""
initial_condition_density_wave(x, t, equations::CompressibleEulerEquations2D)
A sine wave in the density with constant velocity and pressure; reduces the
compressible Euler equations to the linear advection equations.
This setup is the test case for stability of EC fluxes from paper
- Gregor J. Gassner, Magnus Svärd, Florian J. Hindenlang (2020)
Stability issues of entropy-stable and/or split-form high-order schemes
[arXiv: 2007.09026](https://arxiv.org/abs/2007.09026)
with the following parameters
- domain [-1, 1]
- mesh = 4x4
- polydeg = 5
"""
function initial_condition_density_wave(x, t, equations::CompressibleEulerEquations2D)
v1 = 0.1
v2 = 0.2
rho = 1 + 0.98 * sinpi(2 * (x[1] + x[2] - t * (v1 + v2)))
rho_v1 = rho * v1
rho_v2 = rho * v2
p = 20
rho_e = p / (equations.gamma - 1) + 1 / 2 * rho * (v1^2 + v2^2)
return SVector(rho, rho_v1, rho_v2, rho_e)
end
"""
initial_condition_weak_blast_wave(x, t, equations::CompressibleEulerEquations2D)
A weak blast wave taken 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::CompressibleEulerEquations2D)
# From Hennemann & Gassner JCP paper 2020 (Sec. 6.3)
# Set up polar coordinates
inicenter = SVector(0.0, 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)
sin_phi, cos_phi = sincos(phi)
# 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, p), equations)
end
"""
initial_condition_eoc_test_coupled_euler_gravity(x, t, equations::CompressibleEulerEquations2D)
Setup used for convergence tests of the Euler equations with self-gravity used in
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics
[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)
in combination with [`source_terms_eoc_test_coupled_euler_gravity`](@ref)
or [`source_terms_eoc_test_euler`](@ref).
"""
function initial_condition_eoc_test_coupled_euler_gravity(x, t,
equations::CompressibleEulerEquations2D)
# OBS! this assumes that γ = 2 other manufactured source terms are incorrect
if equations.gamma != 2.0
error("adiabatic constant must be 2 for the coupling convergence test")
end
c = 2.0
A = 0.1
ini = c + A * sin(pi * (x[1] + x[2] - t))
G = 1.0 # gravitational constant
rho = ini
v1 = 1.0
v2 = 1.0
p = ini^2 * G / pi # * 2 / ndims, but ndims==2 here
return prim2cons(SVector(rho, v1, v2, p), equations)
end
"""
source_terms_eoc_test_coupled_euler_gravity(u, x, t, equations::CompressibleEulerEquations2D)
Setup used for convergence tests of the Euler equations with self-gravity used in
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics
[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)
in combination with [`initial_condition_eoc_test_coupled_euler_gravity`](@ref).
"""
@inline function source_terms_eoc_test_coupled_euler_gravity(u, x, t,
equations::CompressibleEulerEquations2D)
# Same settings as in `initial_condition_eoc_test_coupled_euler_gravity`
c = 2.0
A = 0.1
G = 1.0 # gravitational constant, must match coupling solver
C_grav = -2 * G / pi # 2 == 4 / ndims
x1, x2 = x
si, co = sincos(pi * (x1 + x2 - t))
rhox = A * pi * co
rho = c + A * si
du1 = rhox
du2 = rhox
du3 = rhox
du4 = (1.0 - C_grav * rho) * rhox
return SVector(du1, du2, du3, du4)
end
"""
source_terms_eoc_test_euler(u, x, t, equations::CompressibleEulerEquations2D)
Setup used for convergence tests of the Euler equations with self-gravity used in
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics
[arXiv: 2008.10593](https://arxiv.org/abs/2008.10593)
in combination with [`initial_condition_eoc_test_coupled_euler_gravity`](@ref).
"""
@inline function source_terms_eoc_test_euler(u, x, t,
equations::CompressibleEulerEquations2D)
# Same settings as in `initial_condition_eoc_test_coupled_euler_gravity`
c = 2.0
A = 0.1
G = 1.0
C_grav = -2 * G / pi # 2 == 4 / ndims
x1, x2 = x
si, co = sincos(pi * (x1 + x2 - t))
rhox = A * pi * co
rho = c + A * si
du1 = rhox
du2 = rhox * (1 - C_grav * rho)
du3 = rhox * (1 - C_grav * rho)
du4 = rhox * (1 - 3 * C_grav * rho)
return SVector(du1, du2, du3, du4)
end
"""
boundary_condition_slip_wall(u_inner, normal_direction, x, t, surface_flux_function,
equations::CompressibleEulerEquations2D)
Determine the boundary numerical surface flux for a slip wall condition.
Imposes a zero normal velocity at the wall.
Density is taken from the internal solution state and pressure is computed as an
exact solution of a 1D Riemann problem. Further details about this boundary state
are available in the paper:
- J. J. W. van der Vegt and H. van der Ven (2002)
Slip flow boundary conditions in discontinuous Galerkin discretizations of
the Euler equations of gas dynamics
[PDF](https://reports.nlr.nl/bitstream/handle/10921/692/TP-2002-300.pdf?sequence=1)
Details about the 1D pressure Riemann solution can be found in Section 6.3.3 of the book
- Eleuterio F. Toro (2009)
Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction
3rd edition
[DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)
Should be used together with [`UnstructuredMesh2D`](@ref).
"""
@inline function boundary_condition_slip_wall(u_inner, normal_direction::AbstractVector,
x, t,
surface_flux_function,
equations::CompressibleEulerEquations2D)
norm_ = norm(normal_direction)
# Normalize the vector without using `normalize` since we need to multiply by the `norm_` later
normal = normal_direction / norm_
# rotate the internal solution state
u_local = rotate_to_x(u_inner, normal, equations)
# compute the primitive variables
rho_local, v_normal, v_tangent, p_local = cons2prim(u_local, equations)
# Get the solution of the pressure Riemann problem
# See Section 6.3.3 of
# Eleuterio F. Toro (2009)
# Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical Introduction
# [DOI: 10.1007/b79761](https://doi.org/10.1007/b79761)
if v_normal <= 0.0
sound_speed = sqrt(equations.gamma * p_local / rho_local) # local sound speed
p_star = p_local *
(1 + 0.5 * (equations.gamma - 1) * v_normal / sound_speed)^(2 *
equations.gamma *
equations.inv_gamma_minus_one)
else # v_normal > 0.0
A = 2 / ((equations.gamma + 1) * rho_local)
B = p_local * (equations.gamma - 1) / (equations.gamma + 1)
p_star = p_local +
0.5 * v_normal / A *
(v_normal + sqrt(v_normal^2 + 4 * A * (p_local + B)))
end
# For the slip wall we directly set the flux as the normal velocity is zero
return SVector(zero(eltype(u_inner)),
p_star * normal[1],
p_star * normal[2],
zero(eltype(u_inner))) * norm_
end
"""
boundary_condition_slip_wall(u_inner, orientation, direction, x, t,
surface_flux_function, equations::CompressibleEulerEquations2D)
Should be used together with [`TreeMesh`](@ref).
"""
@inline function boundary_condition_slip_wall(u_inner, orientation,
direction, x, t,
surface_flux_function,
equations::CompressibleEulerEquations2D)
# get the appropriate normal vector from the orientation
if orientation == 1
normal_direction = SVector(1, 0)
else # orientation == 2
normal_direction = SVector(0, 1)
end
# compute and return the flux using `boundary_condition_slip_wall` routine above
return boundary_condition_slip_wall(u_inner, normal_direction, direction,
x, t, surface_flux_function, equations)
end
"""
boundary_condition_slip_wall(u_inner, normal_direction, direction, x, t,
surface_flux_function, equations::CompressibleEulerEquations2D)
Should be used together with [`StructuredMesh`](@ref).
"""
@inline function boundary_condition_slip_wall(u_inner, normal_direction::AbstractVector,
direction, x, t,
surface_flux_function,
equations::CompressibleEulerEquations2D)
# flip sign of normal to make it outward pointing, then flip the sign of the normal flux back
# to be inward pointing on the -x and -y sides due to the orientation convention used by StructuredMesh
if isodd(direction)
boundary_flux = -boundary_condition_slip_wall(u_inner, -normal_direction,
x, t, surface_flux_function,
equations)
else
boundary_flux = boundary_condition_slip_wall(u_inner, normal_direction,
x, t, surface_flux_function,
equations)
end
return boundary_flux
end
# Calculate 2D flux for a single point
@inline function flux(u, orientation::Integer, equations::CompressibleEulerEquations2D)
rho, rho_v1, rho_v2, rho_e = u
v1 = rho_v1 / rho
v2 = rho_v2 / rho
p = (equations.gamma - 1) * (rho_e - 0.5 * (rho_v1 * v1 + rho_v2 * v2))
if orientation == 1
f1 = rho_v1
f2 = rho_v1 * v1 + p
f3 = rho_v1 * v2
f4 = (rho_e + p) * v1
else
f1 = rho_v2
f2 = rho_v2 * v1
f3 = rho_v2 * v2 + p
f4 = (rho_e + p) * v2
end
return SVector(f1, f2, f3, f4)
end
# Calculate 2D flux for a single point in the normal direction
# Note, this directional vector is not normalized
@inline function flux(u, normal_direction::AbstractVector,
equations::CompressibleEulerEquations2D)
rho_e = last(u)
rho, v1, v2, p = cons2prim(u, equations)
v_normal = v1 * normal_direction[1] + v2 * normal_direction[2]
rho_v_normal = rho * v_normal
f1 = rho_v_normal
f2 = rho_v_normal * v1 + p * normal_direction[1]
f3 = rho_v_normal * v2 + p * normal_direction[2]
f4 = (rho_e + p) * v_normal
return SVector(f1, f2, f3, f4)
end
"""
flux_shima_etal(u_ll, u_rr, orientation_or_normal_direction,
equations::CompressibleEulerEquations2D)
This flux is is a modification of the original kinetic energy preserving two-point flux by
- Yuichi Kuya, Kosuke Totani and Soshi Kawai (2018)
Kinetic energy and entropy preserving schemes for compressible flows
by split convective forms
[DOI: 10.1016/j.jcp.2018.08.058](https://doi.org/10.1016/j.jcp.2018.08.058)
The modification is in the energy flux to guarantee pressure equilibrium and was developed by
- Nao Shima, Yuichi Kuya, Yoshiharu Tamaki, Soshi Kawai (JCP 2020)
Preventing spurious pressure oscillations in split convective form discretizations for
compressible flows
[DOI: 10.1016/j.jcp.2020.110060](https://doi.org/10.1016/j.jcp.2020.110060)
"""
@inline function flux_shima_etal(u_ll, u_rr, orientation::Integer,
equations::CompressibleEulerEquations2D)
# Unpack left and right state
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
# Average each factor of products in flux
rho_avg = 1 / 2 * (rho_ll + rho_rr)
v1_avg = 1 / 2 * (v1_ll + v1_rr)
v2_avg = 1 / 2 * (v2_ll + v2_rr)
p_avg = 1 / 2 * (p_ll + p_rr)
kin_avg = 1 / 2 * (v1_ll * v1_rr + v2_ll * v2_rr)
# Calculate fluxes depending on orientation
if orientation == 1
pv1_avg = 1 / 2 * (p_ll * v1_rr + p_rr * v1_ll)
f1 = rho_avg * v1_avg
f2 = f1 * v1_avg + p_avg
f3 = f1 * v2_avg
f4 = p_avg * v1_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv1_avg
else
pv2_avg = 1 / 2 * (p_ll * v2_rr + p_rr * v2_ll)
f1 = rho_avg * v2_avg
f2 = f1 * v1_avg
f3 = f1 * v2_avg + p_avg
f4 = p_avg * v2_avg * equations.inv_gamma_minus_one + f1 * kin_avg + pv2_avg
end
return SVector(f1, f2, f3, f4)
end
@inline function flux_shima_etal(u_ll, u_rr, normal_direction::AbstractVector,
equations::CompressibleEulerEquations2D)
# Unpack left and right state
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
rho_rr, v1_rr, v2_rr, p_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]
# Average each factor of products in flux
rho_avg = 1 / 2 * (rho_ll + rho_rr)
v1_avg = 1 / 2 * (v1_ll + v1_rr)
v2_avg = 1 / 2 * (v2_ll + v2_rr)
v_dot_n_avg = 1 / 2 * (v_dot_n_ll + v_dot_n_rr)
p_avg = 1 / 2 * (p_ll + p_rr)
velocity_square_avg = 0.5 * (v1_ll * v1_rr + v2_ll * v2_rr)
# Calculate fluxes depending on normal_direction
f1 = rho_avg * v_dot_n_avg
f2 = f1 * v1_avg + p_avg * normal_direction[1]
f3 = f1 * v2_avg + p_avg * normal_direction[2]
f4 = (f1 * velocity_square_avg +
p_avg * v_dot_n_avg * equations.inv_gamma_minus_one
+ 0.5 * (p_ll * v_dot_n_rr + p_rr * v_dot_n_ll))
return SVector(f1, f2, f3, f4)
end
"""
flux_kennedy_gruber(u_ll, u_rr, orientation_or_normal_direction,
equations::CompressibleEulerEquations2D)
Kinetic energy preserving two-point flux by
- Kennedy and Gruber (2008)
Reduced aliasing formulations of the convective terms within the
Navier-Stokes equations for a compressible fluid
[DOI: 10.1016/j.jcp.2007.09.020](https://doi.org/10.1016/j.jcp.2007.09.020)
"""
@inline function flux_kennedy_gruber(u_ll, u_rr, orientation::Integer,
equations::CompressibleEulerEquations2D)
# Unpack left and right state
rho_e_ll = last(u_ll)
rho_e_rr = last(u_rr)
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
# Average each factor of products in flux
rho_avg = 1 / 2 * (rho_ll + rho_rr)
v1_avg = 1 / 2 * (v1_ll + v1_rr)
v2_avg = 1 / 2 * (v2_ll + v2_rr)
p_avg = 1 / 2 * (p_ll + p_rr)
e_avg = 1 / 2 * (rho_e_ll / rho_ll + rho_e_rr / rho_rr)
# Calculate fluxes depending on orientation
if orientation == 1
f1 = rho_avg * v1_avg
f2 = rho_avg * v1_avg * v1_avg + p_avg
f3 = rho_avg * v1_avg * v2_avg
f4 = (rho_avg * e_avg + p_avg) * v1_avg
else
f1 = rho_avg * v2_avg
f2 = rho_avg * v2_avg * v1_avg
f3 = rho_avg * v2_avg * v2_avg + p_avg
f4 = (rho_avg * e_avg + p_avg) * v2_avg
end
return SVector(f1, f2, f3, f4)
end
@inline function flux_kennedy_gruber(u_ll, u_rr, normal_direction::AbstractVector,
equations::CompressibleEulerEquations2D)
# Unpack left and right state
rho_e_ll = last(u_ll)
rho_e_rr = last(u_rr)
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
# Average each factor of products in flux
rho_avg = 0.5 * (rho_ll + rho_rr)
v1_avg = 0.5 * (v1_ll + v1_rr)
v2_avg = 0.5 * (v2_ll + v2_rr)
v_dot_n_avg = v1_avg * normal_direction[1] + v2_avg * normal_direction[2]
p_avg = 0.5 * (p_ll + p_rr)
e_avg = 0.5 * (rho_e_ll / rho_ll + rho_e_rr / rho_rr)
# Calculate fluxes depending on normal_direction
f1 = rho_avg * v_dot_n_avg
f2 = f1 * v1_avg + p_avg * normal_direction[1]
f3 = f1 * v2_avg + p_avg * normal_direction[2]
f4 = f1 * e_avg + p_avg * v_dot_n_avg
return SVector(f1, f2, f3, f4)
end
"""
flux_chandrashekar(u_ll, u_rr, orientation, equations::CompressibleEulerEquations2D)
Entropy conserving two-point flux by
- Chandrashekar (2013)
Kinetic Energy Preserving and Entropy Stable Finite Volume Schemes
for Compressible Euler and Navier-Stokes Equations
[DOI: 10.4208/cicp.170712.010313a](https://doi.org/10.4208/cicp.170712.010313a)
"""
@inline function flux_chandrashekar(u_ll, u_rr, orientation::Integer,
equations::CompressibleEulerEquations2D)
# Unpack left and right state
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
beta_ll = 0.5 * rho_ll / p_ll
beta_rr = 0.5 * rho_rr / p_rr
specific_kin_ll = 0.5 * (v1_ll^2 + v2_ll^2)
specific_kin_rr = 0.5 * (v1_rr^2 + v2_rr^2)
# Compute the necessary mean values
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)
p_mean = 0.5 * rho_avg / beta_avg
velocity_square_avg = specific_kin_ll + specific_kin_rr
# Calculate fluxes depending on orientation
if orientation == 1
f1 = rho_mean * v1_avg
f2 = f1 * v1_avg + p_mean
f3 = f1 * v2_avg
f4 = f1 * 0.5 * (1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +
f2 * v1_avg + f3 * v2_avg
else
f1 = rho_mean * v2_avg
f2 = f1 * v1_avg
f3 = f1 * v2_avg + p_mean
f4 = f1 * 0.5 * (1 / (equations.gamma - 1) / beta_mean - velocity_square_avg) +
f2 * v1_avg + f3 * v2_avg
end
return SVector(f1, f2, f3, f4)
end
"""
flux_ranocha(u_ll, u_rr, orientation_or_normal_direction,
equations::CompressibleEulerEquations2D)
Entropy conserving and kinetic energy preserving two-point flux by
- 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)
See also
- 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_ranocha(u_ll, u_rr, orientation::Integer,
equations::CompressibleEulerEquations2D)
# Unpack left and right state
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
# Compute the necessary mean values
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)
p_avg = 0.5 * (p_ll + p_rr)
velocity_square_avg = 0.5 * (v1_ll * v1_rr + v2_ll * v2_rr)
# Calculate fluxes depending on orientation
if orientation == 1
f1 = rho_mean * v1_avg
f2 = f1 * v1_avg + p_avg
f3 = f1 * v2_avg
f4 = f1 *
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +
0.5 * (p_ll * v1_rr + p_rr * v1_ll)
else
f1 = rho_mean * v2_avg
f2 = f1 * v1_avg
f3 = f1 * v2_avg + p_avg
f4 = f1 *
(velocity_square_avg + inv_rho_p_mean * equations.inv_gamma_minus_one) +
0.5 * (p_ll * v2_rr + p_rr * v2_ll)
end
return SVector(f1, f2, f3, f4)
end
@inline function flux_ranocha(u_ll, u_rr, normal_direction::AbstractVector,
equations::CompressibleEulerEquations2D)
# Unpack left and right state
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
rho_rr, v1_rr, v2_rr, p_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]
# Compute the necessary mean values
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)
p_avg = 0.5 * (p_ll + p_rr)
velocity_square_avg = 0.5 * (v1_ll * v1_rr + v2_ll * v2_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 * normal_direction[1]
f3 = f1 * v2_avg + p_avg * normal_direction[2]
f4 = (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))
return SVector(f1, f2, f3, f4)
end
"""
splitting_steger_warming(u, orientation::Integer,
equations::CompressibleEulerEquations2D)
splitting_steger_warming(u, which::Union{Val{:minus}, Val{:plus}}
orientation::Integer,
equations::CompressibleEulerEquations2D)
Splitting of the compressible Euler flux of Steger and Warming.
Returns a tuple of the fluxes "minus" (associated with waves going into the
negative axis direction) and "plus" (associated with waves going into the
positive axis direction). If only one of the fluxes is required, use the
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
!!! warning "Experimental implementation (upwind SBP)"
This is an experimental feature and may change in future releases.
## References
- Joseph L. Steger and R. F. Warming (1979)
Flux Vector Splitting of the Inviscid Gasdynamic Equations
With Application to Finite Difference Methods
[NASA Technical Memorandum](https://ntrs.nasa.gov/api/citations/19790020779/downloads/19790020779.pdf)
"""
@inline function splitting_steger_warming(u, orientation::Integer,
equations::CompressibleEulerEquations2D)
fm = splitting_steger_warming(u, Val{:minus}(), orientation, equations)
fp = splitting_steger_warming(u, Val{:plus}(), orientation, equations)
return fm, fp
end
@inline function splitting_steger_warming(u, ::Val{:plus}, orientation::Integer,
equations::CompressibleEulerEquations2D)
rho, rho_v1, rho_v2, rho_e = u
v1 = rho_v1 / rho
v2 = rho_v2 / rho
p = (equations.gamma - 1) * (rho_e - 0.5 * (rho_v1 * v1 + rho_v2 * v2))
a = sqrt(equations.gamma * p / rho)
if orientation == 1
lambda1 = v1
lambda2 = v1 + a
lambda3 = v1 - a
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
lambda2_p = positive_part(lambda2)
lambda3_p = positive_part(lambda3)
alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p
rho_2gamma = 0.5 * rho / equations.gamma
f1p = rho_2gamma * alpha_p
f2p = rho_2gamma * (alpha_p * v1 + a * (lambda2_p - lambda3_p))
f3p = rho_2gamma * alpha_p * v2
f4p = rho_2gamma *
(alpha_p * 0.5 * (v1^2 + v2^2) + a * v1 * (lambda2_p - lambda3_p)
+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)
else # orientation == 2
lambda1 = v2
lambda2 = v2 + a
lambda3 = v2 - a
lambda1_p = positive_part(lambda1) # Same as (lambda_i + abs(lambda_i)) / 2, but faster :)
lambda2_p = positive_part(lambda2)
lambda3_p = positive_part(lambda3)
alpha_p = 2 * (equations.gamma - 1) * lambda1_p + lambda2_p + lambda3_p
rho_2gamma = 0.5 * rho / equations.gamma
f1p = rho_2gamma * alpha_p
f2p = rho_2gamma * alpha_p * v1
f3p = rho_2gamma * (alpha_p * v2 + a * (lambda2_p - lambda3_p))
f4p = rho_2gamma *
(alpha_p * 0.5 * (v1^2 + v2^2) + a * v2 * (lambda2_p - lambda3_p)
+ a^2 * (lambda2_p + lambda3_p) * equations.inv_gamma_minus_one)
end
return SVector(f1p, f2p, f3p, f4p)
end
@inline function splitting_steger_warming(u, ::Val{:minus}, orientation::Integer,
equations::CompressibleEulerEquations2D)
rho, rho_v1, rho_v2, rho_e = u
v1 = rho_v1 / rho
v2 = rho_v2 / rho
p = (equations.gamma - 1) * (rho_e - 0.5 * (rho_v1 * v1 + rho_v2 * v2))
a = sqrt(equations.gamma * p / rho)
if orientation == 1
lambda1 = v1
lambda2 = v1 + a
lambda3 = v1 - a
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
lambda2_m = negative_part(lambda2)
lambda3_m = negative_part(lambda3)
alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m
rho_2gamma = 0.5 * rho / equations.gamma
f1m = rho_2gamma * alpha_m
f2m = rho_2gamma * (alpha_m * v1 + a * (lambda2_m - lambda3_m))
f3m = rho_2gamma * alpha_m * v2
f4m = rho_2gamma *
(alpha_m * 0.5 * (v1^2 + v2^2) + a * v1 * (lambda2_m - lambda3_m)
+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)
else # orientation == 2
lambda1 = v2
lambda2 = v2 + a
lambda3 = v2 - a
lambda1_m = negative_part(lambda1) # Same as (lambda_i - abs(lambda_i)) / 2, but faster :)
lambda2_m = negative_part(lambda2)
lambda3_m = negative_part(lambda3)
alpha_m = 2 * (equations.gamma - 1) * lambda1_m + lambda2_m + lambda3_m
rho_2gamma = 0.5 * rho / equations.gamma
f1m = rho_2gamma * alpha_m
f2m = rho_2gamma * alpha_m * v1
f3m = rho_2gamma * (alpha_m * v2 + a * (lambda2_m - lambda3_m))
f4m = rho_2gamma *
(alpha_m * 0.5 * (v1^2 + v2^2) + a * v2 * (lambda2_m - lambda3_m)
+ a^2 * (lambda2_m + lambda3_m) * equations.inv_gamma_minus_one)
end
return SVector(f1m, f2m, f3m, f4m)
end
"""
splitting_vanleer_haenel(u, orientation::Integer,
equations::CompressibleEulerEquations2D)
splitting_vanleer_haenel(u, which::Union{Val{:minus}, Val{:plus}}
orientation::Integer,
equations::CompressibleEulerEquations2D)
Splitting of the compressible Euler flux from van Leer. This splitting further
contains a reformulation due to Hänel et al. where the energy flux uses the
enthalpy. The pressure splitting is independent from the splitting of the
convective terms. As such there are many pressure splittings suggested across
the literature. We implement the 'p4' variant suggested by Liou and Steffen as
it proved the most robust in practice.
Returns a tuple of the fluxes "minus" (associated with waves going into the
negative axis direction) and "plus" (associated with waves going into the
positive axis direction). If only one of the fluxes is required, use the
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
!!! warning "Experimental implementation (upwind SBP)"
This is an experimental feature and may change in future releases.
## References
- Bram van Leer (1982)
Flux-Vector Splitting for the Euler Equation
[DOI: 10.1007/978-3-642-60543-7_5](https://doi.org/10.1007/978-3-642-60543-7_5)
- D. Hänel, R. Schwane and G. Seider (1987)
On the accuracy of upwind schemes for the solution of the Navier-Stokes equations
[DOI: 10.2514/6.1987-1105](https://doi.org/10.2514/6.1987-1105)
- Meng-Sing Liou and Chris J. Steffen, Jr. (1991)
High-Order Polynomial Expansions (HOPE) for Flux-Vector Splitting
[NASA Technical Memorandum](https://ntrs.nasa.gov/citations/19910016425)
"""
@inline function splitting_vanleer_haenel(u, orientation::Integer,
equations::CompressibleEulerEquations2D)
fm = splitting_vanleer_haenel(u, Val{:minus}(), orientation, equations)
fp = splitting_vanleer_haenel(u, Val{:plus}(), orientation, equations)
return fm, fp
end
@inline function splitting_vanleer_haenel(u, ::Val{:plus}, orientation::Integer,
equations::CompressibleEulerEquations2D)
rho, rho_v1, rho_v2, rho_e = u
v1 = rho_v1 / rho
v2 = rho_v2 / rho
p = (equations.gamma - 1) * (rho_e - 0.5 * (rho_v1 * v1 + rho_v2 * v2))
a = sqrt(equations.gamma * p / rho)
H = (rho_e + p) / rho
if orientation == 1
M = v1 / a
p_plus = 0.5 * (1 + equations.gamma * M) * p
f1p = 0.25 * rho * a * (M + 1)^2
f2p = f1p * v1 + p_plus
f3p = f1p * v2
f4p = f1p * H
else # orientation == 2
M = v2 / a
p_plus = 0.5 * (1 + equations.gamma * M) * p
f1p = 0.25 * rho * a * (M + 1)^2
f2p = f1p * v1
f3p = f1p * v2 + p_plus
f4p = f1p * H
end
return SVector(f1p, f2p, f3p, f4p)
end
@inline function splitting_vanleer_haenel(u, ::Val{:minus}, orientation::Integer,
equations::CompressibleEulerEquations2D)
rho, rho_v1, rho_v2, rho_e = u
v1 = rho_v1 / rho
v2 = rho_v2 / rho
p = (equations.gamma - 1) * (rho_e - 0.5 * (rho_v1 * v1 + rho_v2 * v2))
a = sqrt(equations.gamma * p / rho)
H = (rho_e + p) / rho
if orientation == 1
M = v1 / a
p_minus = 0.5 * (1 - equations.gamma * M) * p
f1m = -0.25 * rho * a * (M - 1)^2
f2m = f1m * v1 + p_minus
f3m = f1m * v2
f4m = f1m * H
else # orientation == 2
M = v2 / a
p_minus = 0.5 * (1 - equations.gamma * M) * p
f1m = -0.25 * rho * a * (M - 1)^2
f2m = f1m * v1
f3m = f1m * v2 + p_minus
f4m = f1m * H
end
return SVector(f1m, f2m, f3m, f4m)
end
"""
splitting_lax_friedrichs(u, orientation::Integer,
equations::CompressibleEulerEquations2D)
splitting_lax_friedrichs(u, which::Union{Val{:minus}, Val{:plus}}
orientation::Integer,
equations::CompressibleEulerEquations2D)
Naive local Lax-Friedrichs style flux splitting of the form `f⁺ = 0.5 (f + λ u)`
and `f⁻ = 0.5 (f - λ u)` similar to a flux splitting one would apply, e.g.,
to Burgers' equation.
Returns a tuple of the fluxes "minus" (associated with waves going into the
negative axis direction) and "plus" (associated with waves going into the
positive axis direction). If only one of the fluxes is required, use the
function signature with argument `which` set to `Val{:minus}()` or `Val{:plus}()`.
!!! warning "Experimental implementation (upwind SBP)"
This is an experimental feature and may change in future releases.
"""
@inline function splitting_lax_friedrichs(u, orientation::Integer,
equations::CompressibleEulerEquations2D)
fm = splitting_lax_friedrichs(u, Val{:minus}(), orientation, equations)
fp = splitting_lax_friedrichs(u, Val{:plus}(), orientation, equations)
return fm, fp
end
@inline function splitting_lax_friedrichs(u, ::Val{:plus}, orientation::Integer,
equations::CompressibleEulerEquations2D)
rho, rho_v1, rho_v2, rho_e = u
v1 = rho_v1 / rho
v2 = rho_v2 / rho
p = (equations.gamma - 1) * (rho_e - 0.5 * (rho_v1 * v1 + rho_v2 * v2))
a = sqrt(equations.gamma * p / rho)
H = (rho_e + p) / rho
lambda = 0.5 * (sqrt(v1^2 + v2^2) + a)
if orientation == 1
#lambda = 0.5 * (abs(v1) + a)
f1p = 0.5 * rho * v1 + lambda * u[1]
f2p = 0.5 * rho * v1 * v1 + 0.5 * p + lambda * u[2]
f3p = 0.5 * rho * v1 * v2 + lambda * u[3]
f4p = 0.5 * rho * v1 * H + lambda * u[4]
else # orientation == 2
#lambda = 0.5 * (abs(v2) + a)
f1p = 0.5 * rho * v2 + lambda * u[1]
f2p = 0.5 * rho * v2 * v1 + lambda * u[2]
f3p = 0.5 * rho * v2 * v2 + 0.5 * p + lambda * u[3]
f4p = 0.5 * rho * v2 * H + lambda * u[4]
end
return SVector(f1p, f2p, f3p, f4p)
end
@inline function splitting_lax_friedrichs(u, ::Val{:minus}, orientation::Integer,
equations::CompressibleEulerEquations2D)
rho, rho_v1, rho_v2, rho_e = u
v1 = rho_v1 / rho
v2 = rho_v2 / rho
p = (equations.gamma - 1) * (rho_e - 0.5 * (rho_v1 * v1 + rho_v2 * v2))
a = sqrt(equations.gamma * p / rho)
H = (rho_e + p) / rho
lambda = 0.5 * (sqrt(v1^2 + v2^2) + a)
if orientation == 1
#lambda = 0.5 * (abs(v1) + a)
f1m = 0.5 * rho * v1 - lambda * u[1]
f2m = 0.5 * rho * v1 * v1 + 0.5 * p - lambda * u[2]
f3m = 0.5 * rho * v1 * v2 - lambda * u[3]
f4m = 0.5 * rho * v1 * H - lambda * u[4]
else # orientation == 2
#lambda = 0.5 * (abs(v2) + a)
f1m = 0.5 * rho * v2 - lambda * u[1]
f2m = 0.5 * rho * v2 * v1 - lambda * u[2]
f3m = 0.5 * rho * v2 * v2 + 0.5 * p - lambda * u[3]
f4m = 0.5 * rho * v2 * H - lambda * u[4]
end
return SVector(f1m, f2m, f3m, f4m)
end
# Calculate maximum wave speed for local Lax-Friedrichs-type dissipation as the
# maximum velocity magnitude plus the maximum speed of sound
@inline function max_abs_speed_naive(u_ll, u_rr, orientation::Integer,
equations::CompressibleEulerEquations2D)
rho_ll, v1_ll, v2_ll, p_ll = cons2prim(u_ll, equations)
rho_rr, v1_rr, v2_rr, p_rr = cons2prim(u_rr, equations)
# Get the velocity value in the appropriate direction