Skip to content

Commit bd9e9fa

Browse files
herrwang0Hallberg-NOAA
authored andcommitted
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 59cd564 commit bd9e9fa

File tree

3 files changed

+12
-12
lines changed

3 files changed

+12
-12
lines changed

src/core/MOM_barotropic.F90

Lines changed: 10 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -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: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1553,7 +1553,7 @@ subroutine initialize_dyn_split_RK2(u, v, h, tv, uh, vh, eta, Time, G, GV, US, p
15531553
! Copy eta into an output array.
15541554
do j=js,je ; do i=is,ie ; eta(i,j) = CS%eta(i,j) ; enddo ; enddo
15551555

1556-
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, &
15571557
CS%barotropic_CSp, restart_CS, calc_dtbt, CS%BT_cont, CS%SAL_CSp, HA_CSp)
15581558

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

src/core/MOM_dynamics_split_RK2b.F90

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1469,7 +1469,7 @@ subroutine initialize_dyn_split_RK2b(u, v, h, tv, uh, vh, eta, Time, G, GV, US,
14691469
! Copy eta into an output array.
14701470
do j=js,je ; do i=is,ie ; eta(i,j) = CS%eta(i,j) ; enddo ; enddo
14711471

1472-
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, &
14731473
CS%barotropic_CSp, restart_CS, calc_dtbt, CS%BT_cont, &
14741474
CS%SAL_CSp, HA_CSp)
14751475

0 commit comments

Comments
 (0)