Skip to content

Commit

Permalink
case-opt for wenoz_q and teno_CT
Browse files Browse the repository at this point in the history
  • Loading branch information
ChrisZYJ committed Oct 28, 2024
1 parent 07a8456 commit c7a6a2d
Show file tree
Hide file tree
Showing 5 changed files with 17 additions and 11 deletions.
2 changes: 1 addition & 1 deletion src/common/m_checker_common.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -194,7 +194,7 @@ contains
!> Checks constraints regarding WENO order.
!! Called by s_check_inputs_common for all three stages
subroutine s_check_inputs_weno
@:PROHIBIT(all(weno_order /= (/1, 3, 5, 7/)), "weno_order must be 1, 3, or 5")
@:PROHIBIT(all(weno_order /= (/1, 3, 5, 7/)), "weno_order must be 1, 3, 5, or 7")
@:PROHIBIT(m + 1 < weno_order, "m must be at least weno_order - 1")
@:PROHIBIT(n > 0 .and. n + 1 < weno_order, "n must be at least weno_order - 1")
@:PROHIBIT(p > 0 .and. p + 1 < weno_order, "p must be at least weno_order - 1")
Expand Down
17 changes: 10 additions & 7 deletions src/simulation/m_global_parameters.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -135,6 +135,8 @@ module m_global_parameters
logical, parameter :: mapped_weno = (${mapped_weno}$ /= 0) !< WENO-M (WENO with mapping of nonlinear weights)
logical, parameter :: wenoz = (${wenoz}$ /= 0) !< WENO-Z
logical, parameter :: teno = (${teno}$ /= 0) !< TENO (Targeted ENO)
real(kind(0d0)), parameter :: wenoz_q = ${wenoz_q}$ !< Power constant for WENO-Z
real(kind(0d0)), parameter :: teno_CT = ${teno_CT}$ !< Smoothness threshold for TENO
#:else
integer :: weno_polyn !< Degree of the WENO polynomials (polyn)
integer :: weno_order !< Order of the WENO reconstruction
Expand All @@ -144,11 +146,11 @@ module m_global_parameters
logical :: mapped_weno !< WENO-M (WENO with mapping of nonlinear weights)
logical :: wenoz !< WENO-Z
logical :: teno !< TENO (Targeted ENO)
real(kind(0d0)) :: wenoz_q !< Power constant for WENO-Z
real(kind(0d0)) :: teno_CT !< Smoothness threshold for TENO
#:endif

real(kind(0d0)) :: weno_eps !< Binding for the WENO nonlinear weights
real(kind(0d0)) :: teno_CT !< Smoothness threshold for TENO
real(kind(0d0)) :: wenoz_q !< Power constant for WENO-Z
logical :: mp_weno !< Monotonicity preserving (MP) WENO
logical :: weno_avg ! Average left/right cell-boundary states
logical :: weno_Re_flux !< WENO reconstruct velocity gradients for viscous stress tensor
Expand Down Expand Up @@ -177,10 +179,10 @@ module m_global_parameters
integer :: cpu_start, cpu_end, cpu_rate

#:if not MFC_CASE_OPTIMIZATION
!$acc declare create(num_dims, weno_polyn, weno_order, weno_num_stencils, num_fluids, wenojs, mapped_weno, wenoz, teno)
!$acc declare create(num_dims, weno_polyn, weno_order, weno_num_stencils, num_fluids, wenojs, mapped_weno, wenoz, teno, wenoz_q, teno_CT)
#:endif

!$acc declare create(mpp_lim, model_eqns, mixture_err, alt_soundspeed, avg_state, mp_weno, weno_eps, teno_CT, wenoz_q, hypoelasticity, low_Mach)
!$acc declare create(mpp_lim, model_eqns, mixture_err, alt_soundspeed, avg_state, mp_weno, weno_eps, hypoelasticity, low_Mach)

logical :: relax !< activate phase change
integer :: relax_model !< Relaxation model
Expand Down Expand Up @@ -527,8 +529,6 @@ contains
mpp_lim = .false.
time_stepper = dflt_int
weno_eps = dflt_real
teno_CT = dflt_real
wenoz_q = dflt_real
mp_weno = .false.
weno_avg = .false.
weno_Re_flux = .false.
Expand All @@ -555,6 +555,8 @@ contains
mapped_weno = .false.
wenoz = .false.
teno = .false.
teno_CT = dflt_real
wenoz_q = dflt_real
#:endif

chem_params%advection = .false.
Expand Down Expand Up @@ -1109,10 +1111,11 @@ contains
!$acc update device(cfl_target, m, n, p)

!$acc update device(alt_soundspeed, acoustic_source, num_source)
!$acc update device(dt, sys_size, buff_size, pref, rhoref, gamma_idx, pi_inf_idx, E_idx, alf_idx, stress_idx, mpp_lim, bubbles, hypoelasticity, alt_soundspeed, avg_state, num_fluids, model_eqns, num_dims, mixture_err, grid_geometry, cyl_coord, mp_weno, weno_eps, teno_CT, wenoz_q, low_Mach)
!$acc update device(dt, sys_size, buff_size, pref, rhoref, gamma_idx, pi_inf_idx, E_idx, alf_idx, stress_idx, mpp_lim, bubbles, hypoelasticity, alt_soundspeed, avg_state, num_fluids, model_eqns, num_dims, mixture_err, grid_geometry, cyl_coord, mp_weno, weno_eps, low_Mach)

#:if not MFC_CASE_OPTIMIZATION
!$acc update device(wenojs, mapped_weno, wenoz, teno)
!$acc update device(teno_CT, wenoz_q)
#:endif

!$acc enter data copyin(nb, R0ref, Ca, Web, Re_inv, weight, R0, V0, bubbles, polytropic, polydisperse, qbmm, R0_type, ptil, bubble_model, thermal, poly_sigma)
Expand Down
4 changes: 2 additions & 2 deletions src/simulation/m_start_up.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -136,7 +136,7 @@ contains
t_step_start, t_step_stop, t_step_save, t_step_print, &
model_eqns, mpp_lim, time_stepper, weno_eps, weno_flat, &
riemann_flat, rdma_mpi, cu_tensor, &
teno_CT, wenoz_q, mp_weno, weno_avg, &
mp_weno, weno_avg, &
riemann_solver, low_Mach, wave_speeds, avg_state, &
bc_x, bc_y, bc_z, &
x_domain, y_domain, z_domain, &
Expand All @@ -149,7 +149,7 @@ contains
rhoref, pref, bubbles, bubble_model, &
R0ref, chem_params, &
#:if not MFC_CASE_OPTIMIZATION
nb, mapped_weno, wenoz, teno, weno_order, num_fluids, &
nb, mapped_weno, wenoz, teno, teno_CT, wenoz_q, weno_order, num_fluids, &
#:endif
Ca, Web, Re_inv, &
acoustic_source, acoustic, num_source, &
Expand Down
3 changes: 3 additions & 0 deletions toolchain/mfc/case.py
Original file line number Diff line number Diff line change
Expand Up @@ -199,6 +199,7 @@ def __get_sim_fpp(self, print: bool) -> str:
else:
weno_num_stencils = weno_polyn

# Throw error if wenoz_q or teno_CT are required but not set
return f"""\
#:set MFC_CASE_OPTIMIZATION = {ARG("case_optimization")}
#:set weno_order = {weno_order}
Expand All @@ -212,6 +213,8 @@ def __get_sim_fpp(self, print: bool) -> str:
#:set mapped_weno = {mapped_weno}
#:set wenoz = {wenoz}
#:set teno = {teno}
#:set wenoz_q = {int(self.params.get("wenoz_q", -1))}
#:set teno_CT = {self.params.get("teno_CT", -1.0)}
"""

return """\
Expand Down
2 changes: 1 addition & 1 deletion toolchain/mfc/run/case_dicts.py
Original file line number Diff line number Diff line change
Expand Up @@ -361,7 +361,7 @@ def analytic(self):
ALL.update(SIMULATION)
ALL.update(POST_PROCESS)

CASE_OPTIMIZATION = [ "mapped_weno", "wenoz", "teno", "nb", "weno_order", "num_fluids" ]
CASE_OPTIMIZATION = [ "mapped_weno", "wenoz", "teno", "wenoz_q", "teno_CT", "nb", "weno_order", "num_fluids" ]

_properties = { k: v.value for k, v in ALL.items() }

Expand Down

0 comments on commit c7a6a2d

Please sign in to comment.