Skip to content

Commit

Permalink
Reverted MOM_barotropic.F90 merge with z1l branch.
Browse files Browse the repository at this point in the history
This update give reproducing results for the
following examples:

SOLO TEST CASES
===============
double_gyre
DOME
adjustment2d/layer
adjustment2d/rho
adjustment2d/z
benchmark
lock_exchange
global
global_ALE/layer

OCEAN-SIS TEST CASES
===============
GOLD_SIS
  • Loading branch information
MJHarrison-GFDL committed Sep 25, 2014
1 parent dee81f1 commit 1aa03f4
Showing 1 changed file with 3 additions and 7 deletions.
10 changes: 3 additions & 7 deletions src/core/MOM_barotropic.F90
Original file line number Diff line number Diff line change
Expand Up @@ -707,7 +707,7 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, &
"implemented for a non-Boussinesq model.")
endif ; endif

nstep = max(CEILING(dt/CS%dtbt - 0.0001),4)
nstep = CEILING(dt/CS%dtbt - 0.0001)
if (is_root_PE() .and. (nstep /= CS%nstep_last)) then
write(mesg,'("btstep is using a dynamic barotropic timestep of ", ES12.6, &
& " seconds, max ", ES12.6, ".")') (dt/nstep), CS%dtbt_max
Expand Down Expand Up @@ -776,6 +776,8 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, &
if (.not. use_BT_cont) then
call create_group_pass(pass_Dat_uv, Datu, Datv, CS%BT_Domain, To_All+Scalar_Pair)
endif
call create_group_pass(pass_eta_ubt, eta, CS%BT_Domain)
call create_group_pass(pass_eta_ubt, ubt, vbt, CS%BT_Domain)
if (find_etaav) then
call create_group_pass(pass_etaav, etaav, G%Domain)
endif
Expand Down Expand Up @@ -1440,10 +1442,6 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, &
endif
nfilter = ceiling(dt_filt / dtbt)

call create_group_pass(pass_eta_ubt, eta, CS%BT_Domain)
call create_group_pass(pass_eta_ubt, ubt, vbt, CS%BT_Domain)


! Set up the normalized weights for the filtered velocity.
sum_wt_vel = 0.0 ; sum_wt_eta = 0.0 ; sum_wt_accel = 0.0 ; sum_wt_trans = 0.0
allocate(wt_vel(nstep+nfilter)) ; allocate(wt_eta(nstep+nfilter))
Expand Down Expand Up @@ -2206,8 +2204,6 @@ subroutine set_dtbt(G, CS, eta, pbce, BT_cont, gtot_est, SSH_add)

CS%dtbt = CS%dtbt_fraction * dtbt_max
CS%dtbt_max = dtbt_max


end subroutine set_dtbt

! The following 4 subroutines apply the open boundary conditions.
Expand Down

0 comments on commit 1aa03f4

Please sign in to comment.