Skip to content

Commit

Permalink
Frequency-dependent internal wave drag
Browse files Browse the repository at this point in the history
Minor performance optimization.
  • Loading branch information
c2xu committed Feb 4, 2025
1 parent 3439a85 commit 888e16e
Show file tree
Hide file tree
Showing 2 changed files with 49 additions and 21 deletions.
48 changes: 38 additions & 10 deletions src/core/MOM_barotropic.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1461,8 +1461,6 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, forces, pbce,
Drag_v(i,J) = 0.0
endif
enddo ; enddo
else
Drag_u(:,:) = 0.0 ; Drag_v(:,:) = 0.0
endif

if ((Isq > is-1) .or. (Jsq > js-1)) then
Expand Down Expand Up @@ -2090,12 +2088,19 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, forces, pbce,
!$OMP do schedule(static)
do J=jsv-1,jev ; do i=isv-1,iev+1
v_accel_bt(i,J) = v_accel_bt(i,J) + wt_accel(n) * &
((Cor_v(i,J) + PFv(i,J)) - (vbt(i,J)*Rayleigh_v(i,J) + Drag_v(i,J)))
((Cor_v(i,J) + PFv(i,J)) - vbt(i,J)*Rayleigh_v(i,J))
enddo ; enddo
else
!$OMP do schedule(static)
do J=jsv-1,jev ; do i=isv-1,iev+1
v_accel_bt(i,J) = v_accel_bt(i,J) + wt_accel(n) * (Cor_v(i,J) + PFv(i,J) - Drag_v(i,J))
v_accel_bt(i,J) = v_accel_bt(i,J) + wt_accel(n) * (Cor_v(i,J) + PFv(i,J))
enddo ; enddo
endif

if (CS%use_filter .and. CS%linear_freq_drag) then ! Apply frequency-dependent drag
!$OMP do schedule(static)
do J=jsv-1,jev ; do i=isv-1,iev+1
v_accel_bt(i,J) = v_accel_bt(i,J) - wt_accel(n) * Drag_v(i,J)
enddo ; enddo
endif

Expand Down Expand Up @@ -2168,13 +2173,21 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, forces, pbce,
!$OMP do schedule(static)
do j=jsv,jev ; do I=isv-1,iev
u_accel_bt(I,j) = u_accel_bt(I,j) + wt_accel(n) * &
((Cor_u(I,j) + PFu(I,j)) - (ubt(I,j)*Rayleigh_u(I,j) + Drag_u(I,j)))
((Cor_u(I,j) + PFu(I,j)) - ubt(I,j)*Rayleigh_u(I,j))
enddo ; enddo
!$OMP end do nowait
else
!$OMP do schedule(static)
do j=jsv,jev ; do I=isv-1,iev
u_accel_bt(I,j) = u_accel_bt(I,j) + wt_accel(n) * (Cor_u(I,j) + PFu(I,j) - Drag_u(I,j))
u_accel_bt(I,j) = u_accel_bt(I,j) + wt_accel(n) * (Cor_u(I,j) + PFu(I,j))
enddo ; enddo
!$OMP end do nowait
endif

if (CS%use_filter .and. CS%linear_freq_drag) then ! Apply frequency-dependent drag
!$OMP do schedule(static)
do j=jsv,jev ; do I=isv-1,iev
u_accel_bt(I,j) = u_accel_bt(I,j) - wt_accel(n) * Drag_u(I,j)
enddo ; enddo
!$OMP end do nowait
endif
Expand Down Expand Up @@ -2245,12 +2258,19 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, forces, pbce,
!$OMP do schedule(static)
do j=jsv-1,jev+1 ; do I=isv-1,iev
u_accel_bt(I,j) = u_accel_bt(I,j) + wt_accel(n) * &
((Cor_u(I,j) + PFu(I,j)) - (ubt(I,j)*Rayleigh_u(I,j) + Drag_u(I,j)))
((Cor_u(I,j) + PFu(I,j)) - ubt(I,j)*Rayleigh_u(I,j))
enddo ; enddo
else
!$OMP do schedule(static)
do j=jsv-1,jev+1 ; do I=isv-1,iev
u_accel_bt(I,j) = u_accel_bt(I,j) + wt_accel(n) * (Cor_u(I,j) + PFu(I,j) - Drag_u(I,j))
u_accel_bt(I,j) = u_accel_bt(I,j) + wt_accel(n) * (Cor_u(I,j) + PFu(I,j))
enddo ; enddo
endif

if (CS%use_filter .and. CS%linear_freq_drag) then ! Apply frequency-dependent drag
!$OMP do schedule(static)
do j=jsv-1,jev+1 ; do I=isv-1,iev
u_accel_bt(I,j) = u_accel_bt(I,j) - wt_accel(n) * Drag_u(I,j)
enddo ; enddo
endif

Expand Down Expand Up @@ -2334,13 +2354,21 @@ subroutine btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, forces, pbce,
!$OMP do schedule(static)
do J=jsv-1,jev ; do i=isv,iev
v_accel_bt(i,J) = v_accel_bt(i,J) + wt_accel(n) * &
((Cor_v(i,J) + PFv(i,J)) - (vbt(i,J)*Rayleigh_v(i,J) + Drag_v(i,J)))
((Cor_v(i,J) + PFv(i,J)) - vbt(i,J)*Rayleigh_v(i,J))
enddo ; enddo
!$OMP end do nowait
else
!$OMP do schedule(static)
do J=jsv-1,jev ; do i=isv,iev
v_accel_bt(i,J) = v_accel_bt(i,J) + wt_accel(n) * (Cor_v(i,J) + PFv(i,J) - Drag_v(i,J))
v_accel_bt(i,J) = v_accel_bt(i,J) + wt_accel(n) * (Cor_v(i,J) + PFv(i,J))
enddo ; enddo
!$OMP end do nowait
endif

if (CS%use_filter .and. CS%linear_freq_drag) then ! Apply frequency-dependent drag
!$OMP do schedule(static)
do J=jsv-1,jev ; do i=isv,iev
v_accel_bt(i,J) = v_accel_bt(i,J) - wt_accel(n) * Drag_v(i,J)
enddo ; enddo
!$OMP end do nowait
endif
Expand Down
22 changes: 11 additions & 11 deletions src/parameterizations/lateral/MOM_wave_drag.F90
Original file line number Diff line number Diff line change
Expand Up @@ -18,10 +18,9 @@ module MOM_wave_drag

!> Control structure for the MOM_wave_drag module
type, public :: wave_drag_CS ; private
integer :: nf !< Number of filters to be used in the simulation
!>@{ Spatially varying, frequency-dependent drag coefficients [H T-1 ~> m s-1]
real, allocatable, dimension(:,:,:) :: coef_u, coef_v
!>@}
integer :: nf !< Number of filters to be used in the simulation
real, allocatable, dimension(:,:,:) :: coef_u !< frequency-dependent drag coefficients [H T-1 ~> m s-1]
real, allocatable, dimension(:,:,:) :: coef_v !< frequency-dependent drag coefficients [H T-1 ~> m s-1]
end type wave_drag_CS

contains
Expand Down Expand Up @@ -92,13 +91,14 @@ end subroutine wave_drag_init
subroutine wave_drag_calc(u, v, drag_u, drag_v, G, CS)
type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
type(wave_drag_CS), intent(in) :: CS !< Control structure of MOM_wave_drag
!>@{ Tidal velocities from the output of streaming band-pass filters [L T-1 ~> m s-1]
real, dimension(:,:,:), pointer, intent(in) :: u, v
!>@}
!>@{ Sum of products of tidal velocities and scaled frequency-dependent drag [L2 T-2 ~> m2 s-2]
real, dimension(G%IsdB:G%IedB,G%jsd:G%jed), intent(out) :: drag_u
real, dimension(G%isd:G%ied,G%JsdB:G%JedB), intent(out) :: drag_v
!>@}
real, dimension(:,:,:), pointer, intent(in) :: u !< Zonal velocity from the output of
!! streaming band-pass filters [L T-1 ~> m s-1]
real, dimension(:,:,:), pointer, intent(in) :: v !< Meridional velocity from the output of
!! streaming band-pass filters [L T-1 ~> m s-1]
real, dimension(G%IsdB:G%IedB,G%jsd:G%jed), intent(out) :: drag_u !< Sum of products of filtered velocities
!! and scaled frequency-dependent drag [L2 T-2 ~> m2 s-2]
real, dimension(G%isd:G%ied,G%JsdB:G%JedB), intent(out) :: drag_v !< Sum of products of filtered velocities
!! and scaled frequency-dependent drag [L2 T-2 ~> m2 s-2]

! Local variables
integer :: is, ie, js, je, i, j, k
Expand Down

0 comments on commit 888e16e

Please sign in to comment.