Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Change to pass_var for h and correct barotropic mass source #521

Merged
merged 4 commits into from
Jun 15, 2017
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
36 changes: 12 additions & 24 deletions src/core/MOM.F90
Original file line number Diff line number Diff line change
Expand Up @@ -186,7 +186,7 @@ module MOM
logical :: adiabatic !< If true, then no diapycnal mass fluxes, with no calls
!! to routines to calculate or apply diapycnal fluxes.
logical :: use_temperature !< If true, temp and saln used as state variables.
logical :: calc_rho_for_sea_lev !< If true, calculate rho to convert pressure to sea level
logical :: calc_rho_for_sea_lev !< If true, calculate rho to convert pressure to sea level
logical :: use_frazil !< If true, liquid seawater freezes if temp below freezing,
!! with accumulated heat deficit returned to surface ocean.
logical :: bound_salinity !< If true, salt is added to keep salinity above
Expand Down Expand Up @@ -413,7 +413,6 @@ module MOM

! These are used for group halo updates.
type(group_pass_type) :: pass_tau_ustar_psurf
type(group_pass_type) :: pass_h
type(group_pass_type) :: pass_ray
type(group_pass_type) :: pass_bbl_thick_kv_bbl
type(group_pass_type) :: pass_T_S_h
Expand Down Expand Up @@ -568,8 +567,6 @@ subroutine step_MOM(fluxes, state, Time_start, time_interval, CS)
call create_group_pass(CS%pass_tau_ustar_psurf, fluxes%ustar(:,:), G%Domain)
if (ASSOCIATED(fluxes%p_surf)) &
call create_group_pass(CS%pass_tau_ustar_psurf, fluxes%p_surf(:,:), G%Domain)
if ((CS%thickness_diffuse .and. (.not.CS%thickness_diffuse_first .or. CS%t_dyn_rel_adv == 0)) .OR. CS%mixedlayer_restrat) &
call create_group_pass(CS%pass_h, h, G%Domain) !###, halo=max(2,cont_stensil))

do_pass_Ray = .FALSE.
if ((.not.G%Domain%symmetric) .and. &
Expand Down Expand Up @@ -758,7 +755,7 @@ subroutine step_MOM(fluxes, state, Time_start, time_interval, CS)
CS%MEKE, CS%VarMix, CS%CDp, CS%thickness_diffuse_CSp)
call cpu_clock_end(id_clock_thick_diff)
call cpu_clock_begin(id_clock_pass)
call do_group_pass(CS%pass_h, G%Domain)
call pass_var(h, G%Domain) !###, halo=max(2,cont_stensil))
call cpu_clock_end(id_clock_pass)
call disable_averaging(CS%diag)
if (showCallTree) call callTree_waypoint("finished thickness_diffuse_first (step_MOM)")
Expand All @@ -774,17 +771,12 @@ subroutine step_MOM(fluxes, state, Time_start, time_interval, CS)
! recalculated. This always occurs at the start of a coupling time
! step because the externally prescribed stresses may have changed.
do_calc_bbl = ((CS%t_dyn_rel_adv == 0.0) .or. (n==1))
if (do_calc_bbl) then ; if (thermo_does_span_coupling) then
bbl_time_int = dt_therm
else !### This is inconsistent with corresponding expressions above. -RWH
bbl_time_int = dt*real(1+MIN(ntstep-MOD(n,ntstep),n_max-n))
endif ; endif

if (do_calc_bbl) then
! Calculate the BBL properties and store them inside visc (u,h).
call cpu_clock_begin(id_clock_BBL_visc)
bbl_time_int = max(dt, min(dt_therm - CS%t_dyn_rel_adv, dt*(1+n_max-n)) )
call enable_averaging(bbl_time_int, &
Time_local+set_time(int(bbl_time_int-dt)), CS%diag)
Time_local+set_time(int(bbl_time_int-dt+0.5)), CS%diag)
call set_viscous_BBL(u, v, h, CS%tv, CS%visc, G, GV, CS%set_visc_CSp)
call disable_averaging(CS%diag)
call cpu_clock_end(id_clock_BBL_visc)
Expand Down Expand Up @@ -852,9 +844,7 @@ subroutine step_MOM(fluxes, state, Time_start, time_interval, CS)
dtbt_reset_time = CS%rel_time
endif

mass_src_time = CS%t_dyn_rel_adv
!### This should be mass_src_time = CS%t_dyn_rel_thermo

mass_src_time = CS%t_dyn_rel_thermo
if (CS%legacy_split) then
call step_MOM_dyn_legacy_split(u, v, h, CS%tv, CS%visc, &
Time_local, dt, fluxes, CS%p_surf_begin, CS%p_surf_end, &
Expand Down Expand Up @@ -905,7 +895,7 @@ subroutine step_MOM(fluxes, state, Time_start, time_interval, CS)
if (CS%debug) call hchksum(h,"Post-thickness_diffuse h", G%HI, haloshift=1, scale=GV%H_to_m)
call cpu_clock_end(id_clock_thick_diff)
call cpu_clock_begin(id_clock_pass)
call do_group_pass(CS%pass_h, G%Domain)
call pass_var(h, G%Domain) !###, halo=max(2,cont_stensil))
call cpu_clock_end(id_clock_pass)
if (showCallTree) call callTree_waypoint("finished thickness_diffuse (step_MOM)")
endif
Expand All @@ -922,7 +912,7 @@ subroutine step_MOM(fluxes, state, Time_start, time_interval, CS)
G, GV, CS%mixedlayer_restrat_CSp)
call cpu_clock_end(id_clock_ml_restrat)
call cpu_clock_begin(id_clock_pass)
call do_group_pass(CS%pass_h, G%Domain)
call pass_var(h, G%Domain) !###, halo=max(2,cont_stensil))
call cpu_clock_end(id_clock_pass)
if (CS%debug) then
call hchksum(h,"Post-mixedlayer_restrat h", G%HI, haloshift=1, scale=GV%H_to_m)
Expand Down Expand Up @@ -1403,7 +1393,7 @@ subroutine step_tracers(fluxes, state, Time_start, time_interval, CS)
fluxes%netMassOut(:,:) = CS%offline_CSp%netMassOut(:,:)
call offline_diabatic_ale(fluxes, Time_start, Time_end, time_interval, CS%offline_CSp, &
CS%h, eatr, ebtr)
call pass_var(CS%h,G%Domain)
call pass_var(CS%h, G%Domain)

! Do the transport, the final ALE remappings, horizontal diffusion if it is
! the last iteration
Expand All @@ -1422,9 +1412,7 @@ subroutine step_tracers(fluxes, state, Time_start, time_interval, CS)
call cpu_clock_begin(id_clock_ALE)
call ALE_offline_tracer_final( G, GV, CS%h, h_end, CS%tracer_Reg, CS%ALE_CSp)
call cpu_clock_end(id_clock_ALE)
call pass_var(CS%h,G%Domain)


call pass_var(CS%h, G%Domain)
endif

else ! NON-ALE MODE...NOT WELL TESTED
Expand Down Expand Up @@ -1452,9 +1440,9 @@ subroutine step_tracers(fluxes, state, Time_start, time_interval, CS)
CS%tv%S = salt_mean
CS%h = h_end

call pass_var(CS%tv%T,G%Domain)
call pass_var(CS%tv%S,G%Domain)
call pass_var(CS%h,G%Domain)
call pass_var(CS%tv%T, G%Domain)
call pass_var(CS%tv%S, G%Domain)
call pass_var(CS%h, G%Domain)

endif

Expand Down
7 changes: 3 additions & 4 deletions src/core/MOM_barotropic.F90
Original file line number Diff line number Diff line change
Expand Up @@ -2342,20 +2342,19 @@ subroutine apply_velocity_OBCs(OBC, ubt, vbt, uhbt, vhbt, ubt_trans, vbt_trans,
!! column mass anomaly, in m or kg m-2.
real, dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: ubt_old !< The starting value of ubt in a barotropic step,
!! m s-1.
real, dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: vbt_old !< The starting value of ubt in a barotropic step,
real, dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: vbt_old !< The starting value of vbt in a barotropic step,
!! m s-1.
type(BT_OBC_type), intent(in) :: BT_OBC !< A structure with the private barotropic arrays
!! related to the open boundary conditions,
!! set by set_up_BT_OBC.
integer, intent(in) :: halo !< The extra halo size to use here.
real, intent(in) :: dtbt !< The time step, in s.
real, intent(in) :: bebt !< The fractional weighting of the future velocity
!! in
!! determining the transport.
!! in determining the transport.
logical, intent(in) :: use_BT_cont !< If true, use the BT_cont_types to calculate
!! transports.
real, dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: Datu !< A fixed estimate of the face areas at u points.
real, dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: Datv !< A fixed estimate of the face areas at u points.
real, dimension(SZIW_(MS),SZJBW_(MS)), intent(in) :: Datv !< A fixed estimate of the face areas at v points.
type(local_BT_cont_u_type), dimension(SZIBW_(MS),SZJW_(MS)), intent(in) :: BTCL_u !< Structure of information used
!! for a dynamic estimate of the face areas at
!! u-points.
Expand Down
2 changes: 1 addition & 1 deletion src/core/MOM_dynamics_split_RK2.F90
Original file line number Diff line number Diff line change
Expand Up @@ -214,7 +214,7 @@ subroutine step_MOM_dyn_split_RK2(u, v, h, tv, visc, &
type(forcing), intent(in) :: fluxes !< forcing fields
real, dimension(:,:), pointer :: p_surf_begin !< surf pressure at start of this dynamic time step (Pa)
real, dimension(:,:), pointer :: p_surf_end !< surf pressure at end of this dynamic time step (Pa)
real, intent(in) :: dt_since_flux !< elapesed time since fluxes were applied (sec)
real, intent(in) :: dt_since_flux !< elapsed time since fluxes were applied (sec)
real, intent(in) :: dt_therm !< thermodynamic time step (sec)
real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), target, intent(inout) :: uh !< zonal volume/mass transport (m3/s or kg/s)
real, dimension(SZI_(G),SZJB_(G),SZK_(G)), target, intent(inout) :: vh !< merid volume/mass transport (m3/s or kg/s)
Expand Down