Skip to content

Commit ef426b5

Browse files
author
Hyeoksu Lee
committed
add docs and test suite
1 parent 1334062 commit ef426b5

File tree

4 files changed

+38
-10
lines changed

4 files changed

+38
-10
lines changed

docs/documentation/case.md

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -461,7 +461,7 @@ It is recommended to set `weno_eps` to $10^{-6}$ for WENO-JS, and to $10^{-40}$
461461
`riemann_solver = 1`, `2`, and `3` correspond to HLL, HLLC, and Exact Riemann solver, respectively ([Toro, 2013](references.md)).
462462
`riemann_solver = 4` is only for MHD simulations. It resolves 5 of the full seven-wave structure of the MHD equations ([Miyoshi and Kusano, 2005](references.md)).
463463

464-
- `low_Mach` specifies the choice of the low Mach number correction scheme for the HLLC Riemann solver. `low_Mach = 0` is default value and does not apply any correction scheme. `low_Mach = 1` and `2` apply the anti-dissipation pressure correction method ([Chen et al., 2022](references.md)) and the improved velocity reconstruction method ([Thornber et al., 2008](references.md)). This feature requires `riemann_solver = 2` and `model_eqns = 2`.
464+
- `low_Mach` specifies the choice of the low Mach number correction scheme for the HLLC Riemann solver. `low_Mach = 0` is default value and does not apply any correction scheme. `low_Mach = 1` and `2` apply the anti-dissipation pressure correction method ([Chen et al., 2022](references.md)) and the improved velocity reconstruction method ([Thornber et al., 2008](references.md)). This feature requires `model_eqns = 2` or `3`. `low_Mach = 1` works for `riemann_solver = 1` and `2`, but `low_Mach = 2` only works for `riemann_solver = 2`.
465465

466466
- `avg_state` specifies the choice of the method to compute averaged variables at the cell-boundaries from the left and the right states in the Riemann solver by an integer of 1 or 2.
467467
`avg_state = 1` and `2` correspond to Roe- and arithmetic averages, respectively.

src/simulation/m_checker.fpp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -99,7 +99,7 @@ contains
9999
@:PROHIBIT(riemann_solver /= 3 .and. avg_state == dflt_int, "avg_state must be set if riemann_solver != 3")
100100
@:PROHIBIT(all(low_Mach /= (/0, 1, 2/)), "low_Mach must be 0, 1 or 2")
101101
@:PROHIBIT(riemann_solver /= 2 .and. low_Mach == 2, "low_Mach = 2 requires riemann_solver = 2")
102-
@:PROHIBIT(low_Mach /= 0 .and. model_eqns /= 2, "low_Mach = 1 or 2 requires model_eqns = 2")
102+
@:PROHIBIT(low_Mach /= 0 .and. all(model_eqns /= (/2, 3/)), "low_Mach = 1 or 2 requires model_eqns = 2 or 3")
103103
end subroutine s_check_inputs_riemann_solver
104104

105105
!> Checks constraints on geometry and precision

src/simulation/m_riemann_solvers.fpp

Lines changed: 30 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -369,7 +369,8 @@ contains
369369
!$acc xi_field_L, xi_field_R, &
370370
!$acc Cp_iL, Cp_iR, Xs_L, Xs_R, Gamma_iL, Gamma_iR, &
371371
!$acc Yi_avg, Phi_avg, h_iL, h_iR, h_avg_2, &
372-
!$acc c_fast, pres_mag, B, Ga, vdotB, B2, b4, cm)
372+
!$acc c_fast, pres_mag, B, Ga, vdotB, B2, b4, cm, &
373+
!$acc pcorr, zcoef, vel_L_tmp, vel_R_tmp)
373374
do l = is3%beg, is3%end
374375
do k = is2%beg, is2%end
375376
do j = is1%beg, is1%end
@@ -750,6 +751,7 @@ contains
750751
+ (5e-1_wp - sign(5e-1_wp, s_L)) &
751752
*(5e-1_wp + sign(5e-1_wp, s_R))
752753

754+
! Low Mach correction
753755
if (low_Mach == 1) then
754756
@:compute_low_Mach_correction()
755757
else
@@ -1299,7 +1301,8 @@ contains
12991301
!$acc private(vel_L, vel_R, vel_K_Star, Re_L, Re_R, rho_avg, h_avg, gamma_avg, &
13001302
!$acc s_L, s_R, s_S, vel_avg_rms, alpha_L, alpha_R, Ys_L, Ys_R, Xs_L, Xs_R, &
13011303
!$acc Gamma_iL, Gamma_iR, Cp_iL, Cp_iR, Yi_avg, Phi_avg, h_iL, h_iR, h_avg_2, &
1302-
!$acc tau_e_L, tau_e_R, G_L, G_R, flux_ene_e, xi_field_L, xi_field_R)
1304+
!$acc tau_e_L, tau_e_R, G_L, G_R, flux_ene_e, xi_field_L, xi_field_R, pcorr, &
1305+
!$acc zcoef, vel_L_tmp, vel_R_tmp)
13031306
do l = is3%beg, is3%end
13041307
do k = is2%beg, is2%end
13051308
do j = is1%beg, is1%end
@@ -1488,6 +1491,11 @@ contains
14881491
end do
14891492
end if
14901493
1494+
! Low Mach correction
1495+
if (low_Mach == 2) then
1496+
@:compute_low_Mach_correction()
1497+
end if
1498+
14911499
! COMPUTING THE DIRECT WAVE SPEEDS
14921500
if (wave_speeds == 1) then
14931501
if (elasticity) then
@@ -1563,6 +1571,13 @@ contains
15631571
vel_K_Star = vel_L(idx1)*(1_wp - xi_MP) + xi_MP*vel_R(idx1) + &
15641572
xi_MP*xi_PP*(s_S - vel_R(idx1))
15651573
1574+
! Low Mach correction
1575+
if (low_Mach == 1) then
1576+
@:compute_low_Mach_correction()
1577+
else
1578+
pcorr = 0._wp
1579+
end if
1580+
15661581
! COMPUTING FLUXES
15671582
! MASS FLUX.
15681583
!$acc loop seq
@@ -1578,12 +1593,14 @@ contains
15781593
do i = 1, num_dims
15791594
idxi = dir_idx(i)
15801595
flux_rs${XYZ}$_vf(j, k, l, contxe + idxi) = rho_Star*vel_K_Star* &
1581-
(dir_flg(idxi)*vel_K_Star + (1_wp - dir_flg(idxi))*(xi_M*vel_L(idxi) + xi_P*vel_R(idxi))) + dir_flg(idxi)*p_Star
1596+
(dir_flg(idxi)*vel_K_Star + (1_wp - dir_flg(idxi))*(xi_M*vel_L(idxi) + xi_P*vel_R(idxi))) + dir_flg(idxi)*p_Star &
1597+
+ (s_M/s_L)*(s_P/s_R)*dir_flg(idxi)*pcorr
15821598
end do
15831599
15841600
! ENERGY FLUX.
15851601
! f = u*(E-\sigma), q = E, q_star = \xi*E+(s-u)(\rho s_star - \sigma/(s-u))
1586-
flux_rs${XYZ}$_vf(j, k, l, E_idx) = (E_star + p_Star)*vel_K_Star
1602+
flux_rs${XYZ}$_vf(j, k, l, E_idx) = (E_star + p_Star)*vel_K_Star &
1603+
+ (s_M/s_L)*(s_P/s_R)*pcorr*s_S
15871604
15881605
! ELASTICITY. Elastic shear stress additions for the momentum and energy flux
15891606
if (elasticity) then
@@ -1635,7 +1652,8 @@ contains
16351652
((xi_M*qL_prim_rs${XYZ}$_vf(j, k, l, i + advxb - 1) + xi_P*qR_prim_rs${XYZ}$_vf(j + 1, k, l, i + advxb - 1))* &
16361653
(gammas(i)*p_K_Star + pi_infs(i)) + &
16371654
(xi_M*qL_prim_rs${XYZ}$_vf(j, k, l, i + contxb - 1) + xi_P*qR_prim_rs${XYZ}$_vf(j + 1, k, l, i + contxb - 1))* &
1638-
qvs(i))*vel_K_Star
1655+
qvs(i))*vel_K_Star &
1656+
+ (s_M/s_L)*(s_P/s_R)*pcorr*s_S*(xi_M*qL_prim_rs${XYZ}$_vf(j, k, l, i + advxb - 1) + xi_P*qR_prim_rs${XYZ}$_vf(j + 1, k, l, i + advxb - 1))
16391657
end do
16401658
16411659
flux_src_rs${XYZ}$_vf(j, k, l, advxb) = vel_src_rs${XYZ}$_vf(j, k, l, idx1)
@@ -1715,6 +1733,7 @@ contains
17151733
do l = is3%beg, is3%end
17161734
do k = is2%beg, is2%end
17171735
do j = is1%beg, is1%end
1736+
17181737
!$acc loop seq
17191738
do i = 1, contxe
17201739
alpha_rho_L(i) = qL_prim_rs${XYZ}$_vf(j, k, l, i)
@@ -2212,6 +2231,7 @@ contains
22122231
end do
22132232
end if
22142233
2234+
! Low Mach correction
22152235
if (low_Mach == 2) then
22162236
@:compute_low_Mach_correction()
22172237
end if
@@ -2262,6 +2282,7 @@ contains
22622282
xi_M = (5e-1_wp + sign(5e-1_wp, s_S))
22632283
xi_P = (5e-1_wp - sign(5e-1_wp, s_S))
22642284
2285+
! Low Mach correction
22652286
if (low_Mach == 1) then
22662287
@:compute_low_Mach_correction()
22672288
else
@@ -2679,6 +2700,7 @@ contains
26792700
end do
26802701
end if
26812702
2703+
! Low Mach correction
26822704
if (low_Mach == 2) then
26832705
@:compute_low_Mach_correction()
26842706
end if
@@ -2739,14 +2761,15 @@ contains
27392761
xi_M = (5e-1_wp + sign(5e-1_wp, s_S))
27402762
xi_P = (5e-1_wp - sign(5e-1_wp, s_S))
27412763
2742-
! COMPUTING THE HLLC FLUXES
2743-
! MASS FLUX.
2764+
! Low Mach correction
27442765
if (low_Mach == 1) then
27452766
@:compute_low_Mach_correction()
27462767
else
27472768
pcorr = 0._wp
27482769
end if
27492770
2771+
! COMPUTING THE HLLC FLUXES
2772+
! MASS FLUX.
27502773
!$acc loop seq
27512774
do i = 1, contxe
27522775
flux_rs${XYZ}$_vf(j, k, l, i) = &

toolchain/mfc/test/cases.py

Lines changed: 6 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -174,12 +174,17 @@ def alter_riemann_solvers(num_fluids):
174174
stack.pop()
175175

176176
def alter_low_Mach_correction():
177-
stack.push('', {'riemann_solver': 2, 'fluid_pp(1)%gamma' : 0.16, 'fluid_pp(1)%pi_inf': 3515.0, 'dt': 1e-7})
177+
stack.push('', {'fluid_pp(1)%gamma' : 0.16, 'fluid_pp(1)%pi_inf': 3515.0, 'dt': 1e-7})
178178

179179
for low_Mach in [1, 2]:
180180
stack.push(f"low_Mach={low_Mach}", {'low_Mach': low_Mach})
181+
stack.push(f"riemann_solver=1",{'riemann_solver': 1})
181182
cases.append(define_case_d(stack, '', {}))
182183
stack.pop()
184+
stack.push(f"riemann_solver=2",{'riemann_solver': 2})
185+
cases.append(define_case_d(stack, '', {}))
186+
stack.pop()
187+
stack.pop()
183188

184189
stack.pop()
185190

0 commit comments

Comments
 (0)