Skip to content

Commit

Permalink
Modification of OBLIQUE OBC.
Browse files Browse the repository at this point in the history
- Use new time in interior value of tangential gradient computation.
- Not noticeably better or worse than old way but seems like it should
  be more consistent.
  • Loading branch information
kshedstrom committed Jul 17, 2017
1 parent 02e6e1d commit 73a5428
Show file tree
Hide file tree
Showing 2 changed files with 16 additions and 14 deletions.
16 changes: 8 additions & 8 deletions src/core/MOM_barotropic.F90
Original file line number Diff line number Diff line change
Expand Up @@ -2404,8 +2404,8 @@ subroutine apply_velocity_OBCs(OBC, ubt, vbt, uhbt, vhbt, ubt_trans, vbt_trans,
elseif (OBC%segment(OBC%segnum_u(I,j))%oblique) then
grad(I,J) = (ubt_old(I,j+1) - ubt_old(I,j)) * G%mask2dBu(I,J)
grad(I,J-1) = (ubt_old(I,j) - ubt_old(I,j-1)) * G%mask2dBu(I,J-1)
grad(I-1,J) = (ubt_old(I-1,j+1) - ubt_old(I-1,j)) * G%mask2dBu(I-1,J)
grad(I-1,J-1) = (ubt_old(I-1,j) - ubt_old(I-1,j-1)) * G%mask2dBu(I-1,J-1)
grad(I-1,J) = (ubt(I-1,j+1) - ubt(I-1,j)) * G%mask2dBu(I-1,J)
grad(I-1,J-1) = (ubt(I-1,j) - ubt(I-1,j-1)) * G%mask2dBu(I-1,J-1)
dhdt = ubt_old(I-1,j)-ubt(I-1,j) !old-new
dhdx = ubt(I-1,j)-ubt(I-2,j) !in new time backward sasha for I-1
! if (OBC%segment(OBC%segnum_u(I,j))%oblique) then
Expand Down Expand Up @@ -2448,8 +2448,8 @@ subroutine apply_velocity_OBCs(OBC, ubt, vbt, uhbt, vhbt, ubt_trans, vbt_trans,
elseif (OBC%segment(OBC%segnum_u(I,j))%oblique) then
grad(I,J) = (ubt_old(I,j+1) - ubt_old(I,j)) * G%mask2dBu(I,J)
grad(I,J-1) = (ubt_old(I,j) - ubt_old(I,j-1)) * G%mask2dBu(I,J-1)
grad(I+1,J) = (ubt_old(I+1,j+1) - ubt_old(I+1,j)) * G%mask2dBu(I+1,J)
grad(I+1,J-1) = (ubt_old(I+1,j) - ubt_old(I+1,j-1)) * G%mask2dBu(I+1,J-1)
grad(I+1,J) = (ubt(I+1,j+1) - ubt(I+1,j)) * G%mask2dBu(I+1,J)
grad(I+1,J-1) = (ubt(I+1,j) - ubt(I+1,j-1)) * G%mask2dBu(I+1,J-1)
dhdt = ubt_old(I+1,j)-ubt(I+1,j) !old-new
dhdx = ubt(I+1,j)-ubt(I+2,j) !in new time backward sasha for I+1
! if (OBC%segment(OBC%segnum_u(I,j))%oblique) then
Expand Down Expand Up @@ -2513,8 +2513,8 @@ subroutine apply_velocity_OBCs(OBC, ubt, vbt, uhbt, vhbt, ubt_trans, vbt_trans,
elseif (OBC%segment(OBC%segnum_v(i,J))%oblique) then
grad(I,J) = (vbt_old(i+1,J) - vbt_old(i,J)) * G%mask2dBu(I,J)
grad(I-1,J) = (vbt_old(i,J) - vbt_old(i-1,J)) * G%mask2dBu(I-1,J)
grad(I,J-1) = (vbt_old(i+1,J-1) - vbt_old(i,J-1)) * G%mask2dBu(I,J-1)
grad(I-1,J-1) = (vbt_old(i,J-1) - vbt_old(i-1,J-1)) * G%mask2dBu(I-1,J-1)
grad(I,J-1) = (vbt(i+1,J-1) - vbt(i,J-1)) * G%mask2dBu(I,J-1)
grad(I-1,J-1) = (vbt(i,J-1) - vbt(i-1,J-1)) * G%mask2dBu(I-1,J-1)
dhdt = vbt_old(i,J-1)-vbt(i,J-1) !old-new
dhdy = vbt(i,J-1)-vbt(i,J-2) !in new time backward sasha for J-1
! if (OBC%segment(OBC%segnum_v(i,J))%oblique) then
Expand Down Expand Up @@ -2558,8 +2558,8 @@ subroutine apply_velocity_OBCs(OBC, ubt, vbt, uhbt, vhbt, ubt_trans, vbt_trans,
elseif (OBC%segment(OBC%segnum_v(i,J))%oblique) then
grad(I,J) = (vbt_old(i+1,J) - vbt_old(i,J)) * G%mask2dBu(I,J)
grad(I-1,J) = (vbt_old(i,J) - vbt_old(i-1,J)) * G%mask2dBu(I-1,J)
grad(I,J+1) = (vbt_old(i+1,J+1) - vbt_old(i,J+1)) * G%mask2dBu(I,J+1)
grad(I-1,J+1) = (vbt_old(i,J+1) - vbt_old(i-1,J+1)) * G%mask2dBu(I-1,J+1)
grad(I,J+1) = (vbt(i+1,J+1) - vbt(i,J+1)) * G%mask2dBu(I,J+1)
grad(I-1,J+1) = (vbt(i,J+1) - vbt(i-1,J+1)) * G%mask2dBu(I-1,J+1)
dhdt = vbt_old(i,J+1)-vbt(i,J+1) !old-new
dhdy = vbt(i,J+1)-vbt(i,J+2) !in new time backward sasha for J+1
! if (OBC%segment(OBC%segnum_v(i,J))%oblique) then
Expand Down
14 changes: 8 additions & 6 deletions src/core/MOM_open_boundary.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1176,7 +1176,7 @@ subroutine radiation_open_bdry_conds(OBC, u_new, u_old, v_new, v_old, G, dt)
do n=1,OBC%number_of_segments
segment=>OBC%segment(n)
if (.not. segment%on_pe) cycle
if (segment%oblique) call gradient_at_q_points(G,segment,u_old,v_old)
if (segment%oblique) call gradient_at_q_points(G, segment, u_new, v_new, u_old, v_old)
if (segment%direction == OBC_DIRECTION_E) then
I=segment%HI%IscB
do k=1,nz ; do j=segment%HI%jsc,segment%HI%jec
Expand Down Expand Up @@ -1443,11 +1443,13 @@ subroutine open_boundary_zero_normal_flow(OBC, G, u, v)
end subroutine open_boundary_zero_normal_flow

!> Calculate the tangential gradient of the normal flow at the boundary q-points.
subroutine gradient_at_q_points(G,segment,uvel,vvel)
subroutine gradient_at_q_points(G, segment, uvel, vvel, uold, vold)
type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
type(OBC_segment_type), pointer :: segment !< OBC segment structure
real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(in) :: uvel !< zonal velocity
real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(in) :: vvel !< meridional velocity
real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(in) :: uold !< zonal velocity
real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(in) :: vold !< meridional velocity
integer :: i,j,k

if (.not. segment%on_pe) return
Expand All @@ -1463,15 +1465,15 @@ subroutine gradient_at_q_points(G,segment,uvel,vvel)
do k=1,G%ke
do J=segment%HI%JscB,segment%HI%JecB
segment%grad_normal(J,1,k) = (uvel(I-1,j+1,k)-uvel(I-1,j,k)) * G%mask2dBu(I-1,J)
segment%grad_normal(J,2,k) = (uvel(I,j+1,k)-uvel(I,j,k)) * G%mask2dBu(I,J)
segment%grad_normal(J,2,k) = (uold(I,j+1,k)-uold(I,j,k)) * G%mask2dBu(I,J)
enddo
enddo
else ! western segment
I=segment%HI%iscB
do k=1,G%ke
do J=segment%HI%JscB,segment%HI%JecB
segment%grad_normal(J,1,k) = (uvel(I+1,j+1,k)-uvel(I+1,j,k)) * G%mask2dBu(I+1,J)
segment%grad_normal(J,2,k) = (uvel(I,j+1,k)-uvel(I,j,k)) * G%mask2dBu(I,J)
segment%grad_normal(J,2,k) = (uold(I,j+1,k)-uold(I,j,k)) * G%mask2dBu(I,J)
enddo
enddo
endif
Expand All @@ -1486,15 +1488,15 @@ subroutine gradient_at_q_points(G,segment,uvel,vvel)
do k=1,G%ke
do I=segment%HI%IscB,segment%HI%IecB
segment%grad_normal(I,1,k) = (vvel(i+1,J-1,k)-vvel(i,J-1,k)) * G%mask2dBu(I,J-1)
segment%grad_normal(I,2,k) = (vvel(i+1,J,k)-vvel(i,J,k)) * G%mask2dBu(I,J)
segment%grad_normal(I,2,k) = (vold(i+1,J,k)-vold(i,J,k)) * G%mask2dBu(I,J)
enddo
enddo
else ! south segment
J=segment%HI%jscB
do k=1,G%ke
do I=segment%HI%IscB,segment%HI%IecB
segment%grad_normal(I,1,k) = (vvel(i+1,J+1,k)-vvel(i,J+1,k)) * G%mask2dBu(I,J+1)
segment%grad_normal(I,2,k) = (vvel(i+1,J,k)-vvel(i,J,k)) * G%mask2dBu(I,J)
segment%grad_normal(I,2,k) = (vold(i+1,J,k)-vold(i,J,k)) * G%mask2dBu(I,J)
enddo
enddo
endif
Expand Down

0 comments on commit 73a5428

Please sign in to comment.