Skip to content

Commit 59cd564

Browse files
herrwang0Hallberg-NOAA
authored andcommitted
Add runtime parameter SET_DTBT_USE_BT_CONT
The runtime parameter allows to use BT_CONT optional input argument in subroutine set_dtbt if BT_CONT type is available. Originally, when BT_CONT is used, set_dtbt is estimated with zero surface height. The eta optinoal input argument has no effect since it only works with NONLINEAR_BT_CONTINUITY is true, which is by default false when BT_CONT is used. The added option allows accurate estimate of dtbt when the mean surface height of the domain is not at the same level, while maintains old answers. Also, an unused field calc_dtbt in type MOM_dyn_split_RK2_CS is removed.
1 parent ced594d commit 59cd564

File tree

2 files changed

+24
-8
lines changed

2 files changed

+24
-8
lines changed

src/core/MOM_dynamics_split_RK2.F90

Lines changed: 12 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -154,8 +154,7 @@ module MOM_dynamics_split_RK2
154154
logical :: split_bottom_stress !< If true, provide the bottom stress
155155
!! calculated by the vertical viscosity to the
156156
!! barotropic solver.
157-
logical :: calc_dtbt !< If true, calculate the barotropic time-step
158-
!! dynamically.
157+
logical :: dtbt_use_bt_cont !< If true, use BT_cont to calculate DTBT.
159158
logical :: store_CAu !< If true, store the Coriolis and advective accelerations at the
160159
!! end of the timestep for use in the next predictor step.
161160
logical :: CAu_pred_stored !< If true, the Coriolis and advective accelerations at the
@@ -648,7 +647,15 @@ subroutine step_MOM_dyn_split_RK2(u_inst, v_inst, h, tv, visc, Time_local, dt, f
648647
endif
649648

650649
call cpu_clock_begin(id_clock_btstep)
651-
if (calc_dtbt) call set_dtbt(G, GV, US, CS%barotropic_CSp, CS%pbce, eta=eta)
650+
if (calc_dtbt) then
651+
if (CS%dtbt_use_bt_cont .and. associated(CS%BT_cont)) then
652+
call set_dtbt(G, GV, US, CS%barotropic_CSp, CS%pbce, BT_cont=CS%BT_cont)
653+
else
654+
! In the following call, eta is only used when NONLINEAR_BT_CONTINUITY is True. Otherwise, dtbt is effectively
655+
! calculated with eta=0. Note that NONLINEAR_BT_CONTINUITY is False if BT_CONT is used, which is the default.
656+
call set_dtbt(G, GV, US, CS%barotropic_CSp, CS%pbce, eta=eta)
657+
endif
658+
endif
652659
if (showCallTree) call callTree_enter("btstep(), MOM_barotropic.F90")
653660
! This is the predictor step call to btstep.
654661
! The CS%ADp argument here stores the weights for certain integrated diagnostics.
@@ -1410,7 +1417,8 @@ subroutine initialize_dyn_split_RK2(u, v, h, tv, uh, vh, eta, Time, G, GV, US, p
14101417
"If SPLIT is false and USE_RK2 is true, BEGW can be "//&
14111418
"between 0 and 0.5 to damp gravity waves.", &
14121419
units="nondim", default=0.0)
1413-
1420+
call get_param(param_file, mdl, "SET_DTBT_USE_BT_CONT", CS%dtbt_use_bt_cont, &
1421+
"If true, use BT_CONT to calculate DTBT if possible.", default=.false.)
14141422
call get_param(param_file, mdl, "SPLIT_BOTTOM_STRESS", CS%split_bottom_stress, &
14151423
"If true, provide the bottom stress calculated by the "//&
14161424
"vertical viscosity to the barotropic solver.", default=.false.)

src/core/MOM_dynamics_split_RK2b.F90

Lines changed: 12 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -157,8 +157,7 @@ module MOM_dynamics_split_RK2b
157157
logical :: split_bottom_stress !< If true, provide the bottom stress
158158
!! calculated by the vertical viscosity to the
159159
!! barotropic solver.
160-
logical :: calc_dtbt !< If true, calculate the barotropic time-step
161-
!! dynamically.
160+
logical :: dtbt_use_bt_cont !< If true, use BT_cont to calculate DTBT.
162161
logical :: calculate_SAL !< If true, calculate self-attraction and loading.
163162
logical :: use_tides !< If true, tidal forcing is enabled.
164163
logical :: remap_aux !< If true, apply ALE remapping to all of the auxiliary 3-D
@@ -672,7 +671,15 @@ subroutine step_MOM_dyn_split_RK2b(u_av, v_av, h, tv, visc, Time_local, dt, forc
672671
uh_ptr => uh_in ; vh_ptr => vh_in ; u_ptr => u_inst ; v_ptr => v_inst
673672

674673
call cpu_clock_begin(id_clock_btstep)
675-
if (calc_dtbt) call set_dtbt(G, GV, US, CS%barotropic_CSp, CS%pbce, eta=eta)
674+
if (calc_dtbt) then
675+
if (CS%dtbt_use_bt_cont .and. associated(CS%BT_cont)) then
676+
call set_dtbt(G, GV, US, CS%barotropic_CSp, CS%pbce, BT_cont=CS%BT_cont)
677+
else
678+
! In the following call, eta is only used when NONLINEAR_BT_CONTINUITY is True. Otherwise, dtbt is effectively
679+
! calculated with eta=0. Note that NONLINEAR_BT_CONTINUITY is False if BT_CONT is used, which is the default.
680+
call set_dtbt(G, GV, US, CS%barotropic_CSp, CS%pbce, eta=eta)
681+
endif
682+
endif
676683
if (showCallTree) call callTree_enter("btstep(), MOM_barotropic.F90")
677684
! This is the predictor step call to btstep.
678685
! The CS%ADp argument here stores the weights for certain integrated diagnostics.
@@ -1335,7 +1342,8 @@ subroutine initialize_dyn_split_RK2b(u, v, h, tv, uh, vh, eta, Time, G, GV, US,
13351342
"If SPLIT is false and USE_RK2 is true, BEGW can be "//&
13361343
"between 0 and 0.5 to damp gravity waves.", &
13371344
units="nondim", default=0.0)
1338-
1345+
call get_param(param_file, mdl, "SET_DTBT_USE_BT_CONT", CS%dtbt_use_bt_cont, &
1346+
"If true, use BT_CONT to calculate DTBT if possible.", default=.false.)
13391347
call get_param(param_file, mdl, "SPLIT_BOTTOM_STRESS", CS%split_bottom_stress, &
13401348
"If true, provide the bottom stress calculated by the "//&
13411349
"vertical viscosity to the barotropic solver.", default=.false.)

0 commit comments

Comments
 (0)