Skip to content

Commit a33ad14

Browse files
committed
Added new ES numerical fluxes
1 parent 31d27ad commit a33ad14

File tree

2 files changed

+269
-0
lines changed

2 files changed

+269
-0
lines changed

src/Trixi.jl

+1
Original file line numberDiff line numberDiff line change
@@ -178,6 +178,7 @@ export flux, flux_central, flux_lax_friedrichs, flux_hll, flux_hllc, flux_hlle,
178178
hydrostatic_reconstruction_chen_noelle, flux_nonconservative_chen_noelle,
179179
flux_hll_chen_noelle,
180180
FluxPlusDissipation, DissipationGlobalLaxFriedrichs, DissipationLocalLaxFriedrichs,
181+
DissipationEntropyStableLump, DissipationEntropyStable,
181182
FluxLaxFriedrichs, max_abs_speed_naive,
182183
FluxHLL, min_max_speed_naive, min_max_speed_davis, min_max_speed_einfeldt,
183184
min_max_speed_chen_noelle,

src/equations/ideal_mhd_multiion_2d.jl

+268
Original file line numberDiff line numberDiff line change
@@ -1043,4 +1043,272 @@ Computes the sum of the densities times the sum of the pressures
10431043
end
10441044
return rho_total * p_total
10451045
end
1046+
1047+
"""
1048+
DissipationEntropyStableLump(max_abs_speed=max_abs_speed_naive)
1049+
1050+
Create a local Lax-Friedrichs dissipation operator where the maximum absolute wave speed
1051+
is estimated as
1052+
`max_abs_speed(u_ll, u_rr, orientation_or_normal_direction, equations)`,
1053+
defaulting to [`max_abs_speed_naive`](@ref).
1054+
"""
1055+
struct DissipationEntropyStableLump{MaxAbsSpeed}
1056+
max_abs_speed::MaxAbsSpeed
1057+
end
1058+
1059+
DissipationEntropyStableLump() = DissipationEntropyStableLump(max_abs_speed_naive)
1060+
1061+
@inline function (dissipation::DissipationEntropyStableLump)(u_ll, u_rr,
1062+
orientation_or_normal_direction,
1063+
equations::IdealMhdMultiIonEquations2D)
1064+
@unpack gammas = equations
1065+
λ = dissipation.max_abs_speed(u_ll, u_rr, orientation_or_normal_direction,
1066+
equations)
1067+
1068+
w_ll = cons2entropy(u_ll, equations)
1069+
w_rr = cons2entropy(u_rr, equations)
1070+
prim_ll = cons2prim(u_ll, equations)
1071+
prim_rr = cons2prim(u_rr, equations)
1072+
B1_ll, B2_ll, B3_ll = magnetic_field(u_ll, equations)
1073+
B1_rr, B2_rr, B3_rr = magnetic_field(u_rr, equations)
1074+
psi_ll = divergence_cleaning_field(u_ll, equations)
1075+
psi_rr = divergence_cleaning_field(u_rr, equations)
1076+
1077+
# Some global averages
1078+
B1_avg = 0.5 * (B1_ll + B1_rr)
1079+
B2_avg = 0.5 * (B2_ll + B2_rr)
1080+
B3_avg = 0.5 * (B3_ll + B3_rr)
1081+
psi_avg = 0.5 * (psi_ll + psi_rr)
1082+
1083+
dissipation = zero(MVector{nvariables(equations), eltype(u_ll)})
1084+
1085+
beta_plus_ll = 0
1086+
beta_plus_rr = 0
1087+
# Get the lumped dissipation for all components
1088+
for k in eachcomponent(equations)
1089+
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = get_component(k, prim_ll, equations)
1090+
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = get_component(k, prim_rr, equations)
1091+
1092+
w1_ll, w2_ll, w3_ll, w4_ll, w5_ll = get_component(k, w_ll, equations)
1093+
w1_rr, w2_rr, w3_rr, w4_rr, w5_rr = get_component(k, w_rr, equations)
1094+
1095+
# Auxiliary variables
1096+
beta_ll = 0.5 * rho_ll / p_ll
1097+
beta_rr = 0.5 * rho_rr / p_rr
1098+
vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2
1099+
vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2
1100+
1101+
# Mean variables
1102+
rho_ln = ln_mean(rho_ll, rho_rr)
1103+
beta_ln = ln_mean(beta_ll, beta_rr)
1104+
rho_avg = 0.5 * (rho_ll + rho_rr)
1105+
v1_avg = 0.5 * (v1_ll + v1_rr)
1106+
v2_avg = 0.5 * (v2_ll + v2_rr)
1107+
v3_avg = 0.5 * (v3_ll + v3_rr)
1108+
beta_avg = 0.5 * (beta_ll + beta_rr)
1109+
tau = 1 / (beta_ll + beta_rr)
1110+
p_mean = 0.5 * rho_avg / beta_avg
1111+
p_star = 0.5 * rho_ln / beta_ln
1112+
vel_norm_avg = 0.5 * (vel_norm_ll + vel_norm_rr)
1113+
vel_avg_norm = v1_avg^2 + v2_avg^2 + v3_avg^2
1114+
E_bar = p_star / (gammas[k] - 1) + 0.5 * rho_ln * (2 * vel_avg_norm - vel_norm_avg)
1115+
1116+
h11 = rho_ln
1117+
h22 = rho_ln * v1_avg^2 + p_mean
1118+
h33 = rho_ln * v2_avg^2 + p_mean
1119+
h44 = rho_ln * v3_avg^2 + p_mean
1120+
h55 = ((p_star^2 / (gammas[k] - 1) + E_bar * E_bar) / rho_ln
1121+
+ vel_avg_norm * p_mean
1122+
+ tau * ( B1_avg^2 + B2_avg^2 + B3_avg^2 + psi_avg^2 ))
1123+
1124+
d1 = -0.5 * λ * h11 * (w1_rr - w1_ll)
1125+
d2 = -0.5 * λ * h22 * (w2_rr - w2_ll)
1126+
d3 = -0.5 * λ * h33 * (w3_rr - w3_ll)
1127+
d4 = -0.5 * λ * h44 * (w4_rr - w4_ll)
1128+
d5 = -0.5 * λ * h55 * (w5_rr - w5_ll)
1129+
1130+
beta_plus_ll += beta_ll
1131+
beta_plus_rr += beta_rr
1132+
1133+
set_component!(dissipation, k, d1, d2, d3, d4, d5, equations)
1134+
end
1135+
1136+
# Set the magnetic field and psi terms
1137+
h_B_psi = 1 / (beta_plus_ll + beta_plus_rr)
1138+
dissipation[1] = -0.5 * λ * h_B_psi * (w_rr[1] - w_ll[1])
1139+
dissipation[2] = -0.5 * λ * h_B_psi * (w_rr[2] - w_ll[2])
1140+
dissipation[3] = -0.5 * λ * h_B_psi * (w_rr[3] - w_ll[3])
1141+
dissipation[end] = -0.5 * λ * h_B_psi * (w_rr[end] - w_ll[end])
1142+
1143+
return dissipation
1144+
end
1145+
1146+
function Base.show(io::IO, d::DissipationEntropyStableLump)
1147+
print(io, "DissipationEntropyStableLump(", d.max_abs_speed, ")")
1148+
end
1149+
1150+
"""
1151+
DissipationEntropyStable(max_abs_speed=max_abs_speed_naive)
1152+
1153+
Create a local Lax-Friedrichs dissipation operator where the maximum absolute wave speed
1154+
is estimated as
1155+
`max_abs_speed(u_ll, u_rr, orientation_or_normal_direction, equations)`,
1156+
defaulting to [`max_abs_speed_naive`](@ref).
1157+
"""
1158+
struct DissipationEntropyStable{MaxAbsSpeed}
1159+
max_abs_speed::MaxAbsSpeed
1160+
end
1161+
1162+
DissipationEntropyStable() = DissipationEntropyStable(max_abs_speed_naive)
1163+
1164+
@inline function (dissipation::DissipationEntropyStable)(u_ll, u_rr,
1165+
orientation_or_normal_direction,
1166+
equations::IdealMhdMultiIonEquations2D)
1167+
@unpack gammas = equations
1168+
λ = dissipation.max_abs_speed(u_ll, u_rr, orientation_or_normal_direction,
1169+
equations)
1170+
1171+
w_ll = cons2entropy(u_ll, equations)
1172+
w_rr = cons2entropy(u_rr, equations)
1173+
prim_ll = cons2prim(u_ll, equations)
1174+
prim_rr = cons2prim(u_rr, equations)
1175+
B1_ll, B2_ll, B3_ll = magnetic_field(u_ll, equations)
1176+
B1_rr, B2_rr, B3_rr = magnetic_field(u_rr, equations)
1177+
psi_ll = divergence_cleaning_field(u_ll, equations)
1178+
psi_rr = divergence_cleaning_field(u_rr, equations)
1179+
1180+
# Some global averages
1181+
B1_avg = 0.5 * (B1_ll + B1_rr)
1182+
B2_avg = 0.5 * (B2_ll + B2_rr)
1183+
B3_avg = 0.5 * (B3_ll + B3_rr)
1184+
psi_avg = 0.5 * (psi_ll + psi_rr)
1185+
1186+
dissipation = zero(MVector{nvariables(equations), eltype(u_ll)})
1187+
1188+
beta_plus_ll = 0
1189+
beta_plus_rr = 0
1190+
# Get the lumped dissipation for all components
1191+
for k in eachcomponent(equations)
1192+
rho_ll, v1_ll, v2_ll, v3_ll, p_ll = get_component(k, prim_ll, equations)
1193+
rho_rr, v1_rr, v2_rr, v3_rr, p_rr = get_component(k, prim_rr, equations)
1194+
1195+
w1_ll, w2_ll, w3_ll, w4_ll, w5_ll = get_component(k, w_ll, equations)
1196+
w1_rr, w2_rr, w3_rr, w4_rr, w5_rr = get_component(k, w_rr, equations)
1197+
1198+
# Auxiliary variables
1199+
beta_ll = 0.5 * rho_ll / p_ll
1200+
beta_rr = 0.5 * rho_rr / p_rr
1201+
vel_norm_ll = v1_ll^2 + v2_ll^2 + v3_ll^2
1202+
vel_norm_rr = v1_rr^2 + v2_rr^2 + v3_rr^2
1203+
1204+
# Mean variables
1205+
rho_ln = ln_mean(rho_ll, rho_rr)
1206+
beta_ln = ln_mean(beta_ll, beta_rr)
1207+
rho_avg = 0.5 * (rho_ll + rho_rr)
1208+
v1_avg = 0.5 * (v1_ll + v1_rr)
1209+
v2_avg = 0.5 * (v2_ll + v2_rr)
1210+
v3_avg = 0.5 * (v3_ll + v3_rr)
1211+
beta_avg = 0.5 * (beta_ll + beta_rr)
1212+
tau = 1 / (beta_ll + beta_rr)
1213+
p_mean = 0.5 * rho_avg / beta_avg
1214+
p_star = 0.5 * rho_ln / beta_ln
1215+
vel_norm_avg = 0.5 * (vel_norm_ll + vel_norm_rr)
1216+
vel_avg_norm = v1_avg^2 + v2_avg^2 + v3_avg^2
1217+
E_bar = p_star / (gammas[k] - 1) + 0.5 * rho_ln * (2 * vel_avg_norm - vel_norm_avg)
1218+
1219+
h11 = rho_ln
1220+
h12 = rho_ln * v1_avg
1221+
h13 = rho_ln * v2_avg
1222+
h14 = rho_ln * v3_avg
1223+
h15 = E_bar
1224+
d1 = -0.5 * λ * (h11 * (w1_rr - w1_ll) +
1225+
h12 * (w2_rr - w2_ll) +
1226+
h13 * (w3_rr - w3_ll) +
1227+
h14 * (w4_rr - w4_ll) +
1228+
h15 * (w5_rr - w5_ll) )
1229+
1230+
h21 = h12
1231+
h22 = rho_ln * v1_avg^2 + p_mean
1232+
h23 = h21 * v2_avg
1233+
h24 = h21 * v3_avg
1234+
h25 = (E_bar + p_mean) * v1_avg
1235+
d2 = -0.5 * λ * (h21 * (w1_rr - w1_ll) +
1236+
h22 * (w2_rr - w2_ll) +
1237+
h23 * (w3_rr - w3_ll) +
1238+
h24 * (w4_rr - w4_ll) +
1239+
h25 * (w5_rr - w5_ll) )
1240+
1241+
h31 = h13
1242+
h32 = h23
1243+
h33 = rho_ln * v2_avg^2 + p_mean
1244+
h34 = h31 * v3_avg
1245+
h35 = (E_bar + p_mean) * v2_avg
1246+
d3 = -0.5 * λ * (h31 * (w1_rr - w1_ll) +
1247+
h32 * (w2_rr - w2_ll) +
1248+
h33 * (w3_rr - w3_ll) +
1249+
h34 * (w4_rr - w4_ll) +
1250+
h35 * (w5_rr - w5_ll) )
1251+
1252+
h41 = h14
1253+
h42 = h24
1254+
h43 = h34
1255+
h44 = rho_ln * v3_avg^2 + p_mean
1256+
h45 = (E_bar + p_mean) * v3_avg
1257+
d4 = -0.5 * λ * (h41 * (w1_rr - w1_ll) +
1258+
h42 * (w2_rr - w2_ll) +
1259+
h43 * (w3_rr - w3_ll) +
1260+
h44 * (w4_rr - w4_ll) +
1261+
h45 * (w5_rr - w5_ll) )
1262+
1263+
h51 = h15
1264+
h52 = h25
1265+
h53 = h35
1266+
h54 = h45
1267+
h55 = ((p_star^2 / (gammas[k] - 1) + E_bar * E_bar) / rho_ln
1268+
+ vel_avg_norm * p_mean
1269+
+ tau * ( B1_avg^2 + B2_avg^2 + B3_avg^2 + psi_avg^2 ))
1270+
d5 = -0.5 * λ * (h51 * (w1_rr - w1_ll) +
1271+
h52 * (w2_rr - w2_ll) +
1272+
h53 * (w3_rr - w3_ll) +
1273+
h54 * (w4_rr - w4_ll) +
1274+
h55 * (w5_rr - w5_ll) )
1275+
1276+
beta_plus_ll += beta_ll
1277+
beta_plus_rr += beta_rr
1278+
1279+
set_component!(dissipation, k, d1, d2, d3, d4, d5, equations)
1280+
end
1281+
1282+
# Set the magnetic field and psi terms
1283+
h_B_psi = 1 / (beta_plus_ll + beta_plus_rr)
1284+
1285+
# diagonal entries
1286+
dissipation[1] = -0.5 * λ * h_B_psi * (w_rr[1] - w_ll[1])
1287+
dissipation[2] = -0.5 * λ * h_B_psi * (w_rr[2] - w_ll[2])
1288+
dissipation[3] = -0.5 * λ * h_B_psi * (w_rr[3] - w_ll[3])
1289+
dissipation[end] = -0.5 * λ * h_B_psi * (w_rr[end] - w_ll[end])
1290+
# Off-diagonal entries
1291+
for k in eachcomponent(equations)
1292+
_, _, _, _, w5_ll = get_component(k, w_ll, equations)
1293+
_, _, _, _, w5_rr = get_component(k, w_rr, equations)
1294+
1295+
dissipation[1] -= 0.5 * λ * h_B_psi * B1_avg * (w5_rr - w5_ll)
1296+
dissipation[2] -= 0.5 * λ * h_B_psi * B2_avg * (w5_rr - w5_ll)
1297+
dissipation[3] -= 0.5 * λ * h_B_psi * B3_avg * (w5_rr - w5_ll)
1298+
dissipation[end] -= 0.5 * λ * h_B_psi * psi_avg * (w5_rr - w5_ll)
1299+
1300+
ind_E = 3 + (k - 1) * 5 + 5
1301+
dissipation[ind_E] -= 0.5 * λ * h_B_psi * B1_avg * (w_rr[1] - w_ll[1])
1302+
dissipation[ind_E] -= 0.5 * λ * h_B_psi * B2_avg * (w_rr[2] - w_ll[2])
1303+
dissipation[ind_E] -= 0.5 * λ * h_B_psi * B3_avg * (w_rr[3] - w_ll[3])
1304+
dissipation[ind_E] -= 0.5 * λ * h_B_psi * psi_avg * (w_rr[end] - w_ll[end])
1305+
end
1306+
1307+
return dissipation
1308+
end
1309+
1310+
function Base.show(io::IO, d::DissipationEntropyStable)
1311+
print(io, "DissipationEntropyStable(", d.max_abs_speed, ")")
1312+
end
1313+
10461314
end # @muladd

0 commit comments

Comments
 (0)