Skip to content

Commit

Permalink
Normalize wt_[uv] in MOM_barotropic
Browse files Browse the repository at this point in the history
Add an option to normalize wt_[uv] in the barotropic solver. This is
fix step 1 of 2 to address the issue that visc_rem is applied
incorrectly to the barotropic solver.

A runtime flag BAROTROPIC_VERTICAL_WEIGHT_FIX is added to control the
behavior of this commit. The current default is false (not applying the
fix).
  • Loading branch information
herrwang0 authored and marshallward committed Jun 6, 2024
1 parent 8521222 commit a296080
Showing 1 changed file with 33 additions and 1 deletion.
34 changes: 33 additions & 1 deletion src/core/MOM_barotropic.F90
Original file line number Diff line number Diff line change
Expand Up @@ -282,6 +282,7 @@ module MOM_barotropic
logical :: tidal_sal_flather !< Apply adjustment to external gravity wave speed
!! consistent with tidal self-attraction and loading
!! used within the barotropic solver
logical :: wt_uv_fix !< If true, use a normalized wt_[uv] for vertical averages.
type(time_type), pointer :: Time => NULL() !< A pointer to the ocean models clock.
type(diag_ctrl), pointer :: diag => NULL() !< A structure that is used to regulate
!! the timing of diagnostic output.
Expand Down Expand Up @@ -506,6 +507,9 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, forces, pbce,
real :: wt_v(SZI_(G),SZJB_(G),SZK_(GV)) ! normalized weights to
! be used in calculating barotropic velocities, possibly with
! sums less than one due to viscous losses [nondim]
real :: Iwt_u_tot(SZIB_(G),SZJ_(G)) ! Iwt_u_tot and Iwt_v_tot are the
real :: Iwt_v_tot(SZI_(G),SZJB_(G)) ! inverses of wt_u and wt_v vertical integrals,
! used to normalize wt_u and wt_v [nondim]
real, dimension(SZIB_(G),SZJ_(G)) :: &
av_rem_u, & ! The weighted average of visc_rem_u [nondim]
tmp_u, & ! A temporary array at u points [L T-2 ~> m s-2] or [nondim]
Expand Down Expand Up @@ -1052,6 +1056,30 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, forces, pbce,
wt_v(i,J,k) = CS%frhatv(i,J,k) * visc_rem
enddo ; enddo ; enddo

if (CS%wt_uv_fix) then
do j=js,je ; do I=is-1,ie ; Iwt_u_tot(I,j) = wt_u(I,j,1) ; enddo ; enddo
do k=2,nz ; do j=js,je ; do I=is-1,ie
Iwt_u_tot(I,j) = Iwt_u_tot(I,j) + wt_u(I,j,k)
enddo ; enddo ; enddo
do j=js,je ; do I=is-1,ie
Iwt_u_tot(I,j) = G%mask2dCu(I,j) / (Iwt_u_tot(I,j) + h_neglect)
enddo ; enddo
do k=1,nz ; do j=js,je ; do I=is-1,ie
wt_u(I,j,k) = wt_u(I,j,k) * Iwt_u_tot(I,j)
enddo ; enddo ; enddo

do J=js-1,je ; do i=is,ie ; Iwt_v_tot(i,J) = wt_v(i,J,1) ; enddo ; enddo
do k=2,nz ; do J=js-1,je ; do i=is,ie
Iwt_v_tot(i,J) = Iwt_v_tot(i,J) + wt_v(i,J,k)
enddo ; enddo ; enddo
do J=js-1,je ; do i=is,ie
Iwt_v_tot(i,J) = G%mask2dCv(i,J) / (Iwt_v_tot(i,J) + h_neglect)
enddo ; enddo
do k=1,nz ; do J=js-1,je ; do i=is,ie
wt_v(i,J,k) = wt_v(i,J,k) * Iwt_v_tot(i,J)
enddo ; enddo ; enddo
endif

! Use u_Cor and v_Cor as the reference values for the Coriolis terms,
! including the viscous remnant.
!$OMP parallel do default(shared)
Expand Down Expand Up @@ -4580,10 +4608,14 @@ subroutine barotropic_init(u, v, h, eta, Time, G, GV, US, param_file, diag, CS,
call get_param(param_file, mdl, "BAROTROPIC_ANSWER_DATE", CS%answer_date, &
"The vintage of the expressions in the barotropic solver. "//&
"Values below 20190101 recover the answers from the end of 2018, "//&
"while higher values uuse more efficient or general expressions.", &
"while higher values use more efficient or general expressions.", &
default=default_answer_date, do_not_log=.not.GV%Boussinesq)
if (.not.GV%Boussinesq) CS%answer_date = max(CS%answer_date, 20230701)

call get_param(param_file, mdl, "BAROTROPIC_VERTICAL_WEIGHT_FIX", CS%wt_uv_fix, &
"If true, use a normalized weight function for vertical averages of "//&
"baroclinic velocity and forcing.", &
default=.false.)
call get_param(param_file, mdl, "TIDES", use_tides, &
"If true, apply tidal momentum forcing.", default=.false.)
call get_param(param_file, mdl, "CALCULATE_SAL", CS%calculate_SAL, &
Expand Down

0 comments on commit a296080

Please sign in to comment.