Skip to content

Commit 075f8b3

Browse files
authored
Add option to use BT_CONT to calculate dtbt (#823)
* Reorder optional input arguments in set_dtbt Change the order of optional input arguments of subroutine set_dtbt so that they are grouped logically. The decription of the subroutine is also updated. * 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. * Minor fixes on setting dtbt in barotropic_init * Fix a bug that local variable calc_dtbt is not initialized. * Rename dtbt_tmp to dtbt_restart for clarity * Add a comment about the usage of subroutine call set_dtbt * Remove an unused input argument eta in barotropic_init
1 parent ec60bf1 commit 075f8b3

File tree

3 files changed

+52
-36
lines changed

3 files changed

+52
-36
lines changed

src/core/MOM_barotropic.F90

Lines changed: 26 additions & 26 deletions
Original file line numberDiff line numberDiff line change
@@ -2847,26 +2847,26 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, forces, pbce,
28472847

28482848
end subroutine btstep
28492849

2850-
!> This subroutine automatically determines an optimal value for dtbt based
2851-
!! on some state of the ocean.
2852-
subroutine set_dtbt(G, GV, US, CS, eta, pbce, BT_cont, gtot_est, SSH_add)
2850+
!> This subroutine automatically determines an optimal value for dtbt based on some state of the ocean. Either pbce or
2851+
!! gtot_est is required to calculate gravitational acceleration. Column thickness can be estimated using BT_cont, eta,
2852+
!! and SSH_add (default=0), with priority given in that order. The subroutine sets CS%dtbt_max and CS%dtbt.
2853+
subroutine set_dtbt(G, GV, US, CS, pbce, gtot_est, BT_cont, eta, SSH_add)
28532854
type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
28542855
type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
28552856
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
28562857
type(barotropic_CS), intent(inout) :: CS !< Barotropic control structure
2857-
real, dimension(SZI_(G),SZJ_(G)), optional, intent(in) :: eta !< The barotropic free surface
2858-
!! height anomaly or column mass anomaly [H ~> m or kg m-2].
2859-
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), optional, intent(in) :: pbce !< The baroclinic pressure
2860-
!! anomaly in each layer due to free surface
2861-
!! height anomalies [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
2862-
type(BT_cont_type), optional, pointer :: BT_cont !< A structure with elements that describe
2863-
!! the effective open face areas as a
2864-
!! function of barotropic flow.
2865-
real, optional, intent(in) :: gtot_est !< An estimate of the total gravitational
2866-
!! acceleration [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
2867-
real, optional, intent(in) :: SSH_add !< An additional contribution to SSH to
2868-
!! provide a margin of error when
2869-
!! calculating the external wave speed [Z ~> m].
2858+
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
2859+
optional, intent(in) :: pbce !< The baroclinic pressure anomaly in each layer due to free
2860+
!! surface height anomalies [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
2861+
real, optional, intent(in) :: gtot_est !< An estimate of the total gravitational acceleration
2862+
!! [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
2863+
type(BT_cont_type), optional, pointer :: BT_cont !< A structure with elements that describe the effective open
2864+
!! face areas as a function of barotropic flow.
2865+
real, dimension(SZI_(G),SZJ_(G)), &
2866+
optional, intent(in) :: eta !< The barotropic free surface height anomaly or column mass
2867+
!! anomaly [H ~> m or kg m-2].
2868+
real, optional, intent(in) :: SSH_add !< An additional contribution to SSH to provide a margin of
2869+
!! error when calculating the external wave speed [Z ~> m].
28702870

28712871
! Local variables
28722872
real, dimension(SZI_(G),SZJ_(G)) :: &
@@ -4424,7 +4424,7 @@ end subroutine bt_mass_source
44244424
!> barotropic_init initializes a number of time-invariant fields used in the
44254425
!! barotropic calculation and initializes any barotropic fields that have not
44264426
!! already been initialized.
4427-
subroutine barotropic_init(u, v, h, eta, Time, G, GV, US, param_file, diag, CS, &
4427+
subroutine barotropic_init(u, v, h, Time, G, GV, US, param_file, diag, CS, &
44284428
restart_CS, calc_dtbt, BT_cont, SAL_CSp, HA_CSp)
44294429
type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
44304430
type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
@@ -4435,9 +4435,6 @@ subroutine barotropic_init(u, v, h, eta, Time, G, GV, US, param_file, diag, CS,
44354435
intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1].
44364436
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
44374437
intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
4438-
real, dimension(SZI_(G),SZJ_(G)), &
4439-
intent(in) :: eta !< Free surface height or column mass anomaly
4440-
!! [Z ~> m] or [H ~> kg m-2].
44414438
type(time_type), target, intent(in) :: Time !< The current model time.
44424439
type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters.
44434440
type(diag_ctrl), target, intent(inout) :: diag !< A structure that is used to regulate diagnostic
@@ -4465,7 +4462,7 @@ subroutine barotropic_init(u, v, h, eta, Time, G, GV, US, param_file, diag, CS,
44654462
real :: SSH_extra ! An estimate of how much higher SSH might get, for use
44664463
! in calculating the safe external wave speed [Z ~> m].
44674464
real :: dtbt_input ! The input value of DTBT, [nondim] if negative or [s] if positive.
4468-
real :: dtbt_tmp ! A temporary copy of CS%dtbt read from a restart file [T ~> s]
4465+
real :: dtbt_restart ! A temporary copy of CS%dtbt read from a restart file [T ~> s]
44694466
real :: wave_drag_scale ! A scaling factor for the barotropic linear wave drag
44704467
! piston velocities [nondim].
44714468
character(len=200) :: inputdir ! The directory in which to find input files.
@@ -4978,9 +4975,9 @@ subroutine barotropic_init(u, v, h, eta, Time, G, GV, US, param_file, diag, CS,
49784975

49794976
CS%dtbt_fraction = 0.98 ; if (dtbt_input < 0.0) CS%dtbt_fraction = -dtbt_input
49804977

4981-
dtbt_tmp = -1.0
4978+
dtbt_restart = -1.0
49824979
if (query_initialized(CS%dtbt, "DTBT", restart_CS)) then
4983-
dtbt_tmp = CS%dtbt
4980+
dtbt_restart = CS%dtbt
49844981
endif
49854982

49864983
! Estimate the maximum stable barotropic time step.
@@ -4991,14 +4988,17 @@ subroutine barotropic_init(u, v, h, eta, Time, G, GV, US, param_file, diag, CS,
49914988
H_to_Z = GV%H_to_RZ / CS%Rho_BT_lin
49924989
do k=1,GV%ke ; gtot_estimate = gtot_estimate + H_to_Z*GV%g_prime(K) ; enddo
49934990
endif
4991+
4992+
! CS%dtbt calculated here by set_dtbt is only used when dtbt is not reset during the run, i.e. DTBT_RESET_PERIOD<0.
49944993
call set_dtbt(G, GV, US, CS, gtot_est=gtot_estimate, SSH_add=SSH_extra)
49954994

49964995
if (dtbt_input > 0.0) then
49974996
CS%dtbt = US%s_to_T * dtbt_input
4998-
elseif (dtbt_tmp > 0.0) then
4999-
CS%dtbt = dtbt_tmp
4997+
elseif (dtbt_restart > 0.0) then
4998+
CS%dtbt = dtbt_restart
50004999
endif
5001-
if ((dtbt_tmp > 0.0) .and. (dtbt_input > 0.0)) calc_dtbt = .false.
5000+
5001+
calc_dtbt = .true. ; if ((dtbt_restart > 0.0) .and. (dtbt_input > 0.0)) calc_dtbt = .false.
50025002

50035003
call log_param(param_file, mdl, "DTBT as used", CS%dtbt, units="s", unscale=US%T_to_s)
50045004
call log_param(param_file, mdl, "estimated maximum DTBT", CS%dtbt_max, units="s", unscale=US%T_to_s)

src/core/MOM_dynamics_split_RK2.F90

Lines changed: 13 additions & 5 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, eta, CS%pbce)
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.)
@@ -1545,7 +1553,7 @@ subroutine initialize_dyn_split_RK2(u, v, h, tv, uh, vh, eta, Time, G, GV, US, p
15451553
! Copy eta into an output array.
15461554
do j=js,je ; do i=is,ie ; eta(i,j) = CS%eta(i,j) ; enddo ; enddo
15471555

1548-
call barotropic_init(u, v, h, CS%eta, Time, G, GV, US, param_file, diag, &
1556+
call barotropic_init(u, v, h, Time, G, GV, US, param_file, diag, &
15491557
CS%barotropic_CSp, restart_CS, calc_dtbt, CS%BT_cont, CS%SAL_CSp, HA_CSp)
15501558

15511559
if (.not. query_initialized(CS%diffu, "diffu", restart_CS) .or. &

src/core/MOM_dynamics_split_RK2b.F90

Lines changed: 13 additions & 5 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, eta, CS%pbce)
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.)
@@ -1461,7 +1469,7 @@ subroutine initialize_dyn_split_RK2b(u, v, h, tv, uh, vh, eta, Time, G, GV, US,
14611469
! Copy eta into an output array.
14621470
do j=js,je ; do i=is,ie ; eta(i,j) = CS%eta(i,j) ; enddo ; enddo
14631471

1464-
call barotropic_init(u, v, h, CS%eta, Time, G, GV, US, param_file, diag, &
1472+
call barotropic_init(u, v, h, Time, G, GV, US, param_file, diag, &
14651473
CS%barotropic_CSp, restart_CS, calc_dtbt, CS%BT_cont, &
14661474
CS%SAL_CSp, HA_CSp)
14671475

0 commit comments

Comments
 (0)