Skip to content

Commit

Permalink
Refactor spherical_harmonics_forward
Browse files Browse the repository at this point in the history
  Refactored the spherical_harmonics_forward routine in MOM_spherical_harmonics
to work in rescaled units by making use of the unscale arguments to
reproducing_sum().  A total of 4 rescaling variables were moved into unscale
arguments, and a block of code that restores the scaling of the output variables
was eliminated.  All answers and diagnostics are bitwise identical, and no
interfaces are changed.
  • Loading branch information
Hallberg-NOAA committed Dec 12, 2024
1 parent a4d13e8 commit 2c26118
Showing 1 changed file with 7 additions and 18 deletions.
25 changes: 7 additions & 18 deletions src/parameterizations/lateral/MOM_spherical_harmonics.F90
Original file line number Diff line number Diff line change
Expand Up @@ -61,10 +61,7 @@ subroutine spherical_harmonics_forward(G, CS, var, Snm_Re, Snm_Im, Nd, tmp_scale
pmn, & ! Current associated Legendre polynomials of degree n and order m [nondim]
pmnm1, & ! Associated Legendre polynomials of degree n-1 and order m [nondim]
pmnm2 ! Associated Legendre polynomials of degree n-2 and order m [nondim]
real :: scale ! A rescaling factor to temporarily convert var to MKS units during the
! reproducing sums [a A-1 ~> 1]
real :: I_scale ! The inverse of scale [A a-1 ~> 1]
real :: sum_tot ! The total of all components output by the reproducing sum in arbitrary units [a]
real :: sum_tot ! The total of all components output by the reproducing sum in arbitrary units [A]
integer :: i, j, k
integer :: is, ie, js, je, isd, ied, jsd, jed
integer :: m, n, l
Expand All @@ -88,22 +85,21 @@ subroutine spherical_harmonics_forward(G, CS, var, Snm_Re, Snm_Im, Nd, tmp_scale
do l=1,Ltot ; Snm_Re(l) = 0.0; Snm_Im(l) = 0.0 ; enddo

if (CS%reprod_sum) then
scale = 1.0 ; if (present(tmp_scale)) scale = tmp_scale
do m=0,Nmax
l = order2index(m, Nmax)

do j=js,je ; do i=is,ie
CS%Snm_Re_raw(i,j,l) = (scale*var(i,j)) * CS%Pmm(i,j,m+1) * CS%cos_lonT_wtd(i,j,m+1)
CS%Snm_Im_raw(i,j,l) = (scale*var(i,j)) * CS%Pmm(i,j,m+1) * CS%sin_lonT_wtd(i,j,m+1)
CS%Snm_Re_raw(i,j,l) = var(i,j) * CS%Pmm(i,j,m+1) * CS%cos_lonT_wtd(i,j,m+1)
CS%Snm_Im_raw(i,j,l) = var(i,j) * CS%Pmm(i,j,m+1) * CS%sin_lonT_wtd(i,j,m+1)
pmnm2(i,j) = 0.0
pmnm1(i,j) = CS%Pmm(i,j,m+1)
enddo ; enddo

do n = m+1, Nmax ; do j=js,je ; do i=is,ie
pmn(i,j) = &
CS%a_recur(n+1,m+1) * CS%cos_clatT(i,j) * pmnm1(i,j) - CS%b_recur(n+1,m+1) * pmnm2(i,j)
CS%Snm_Re_raw(i,j,l+n-m) = (scale*var(i,j)) * pmn(i,j) * CS%cos_lonT_wtd(i,j,m+1)
CS%Snm_Im_raw(i,j,l+n-m) = (scale*var(i,j)) * pmn(i,j) * CS%sin_lonT_wtd(i,j,m+1)
CS%Snm_Re_raw(i,j,l+n-m) = var(i,j) * pmn(i,j) * CS%cos_lonT_wtd(i,j,m+1)
CS%Snm_Im_raw(i,j,l+n-m) = var(i,j) * pmn(i,j) * CS%sin_lonT_wtd(i,j,m+1)
pmnm2(i,j) = pmnm1(i,j)
pmnm1(i,j) = pmn(i,j)
enddo ; enddo ; enddo
Expand Down Expand Up @@ -133,15 +129,8 @@ subroutine spherical_harmonics_forward(G, CS, var, Snm_Re, Snm_Im, Nd, tmp_scale
if (id_clock_sht_global_sum>0) call cpu_clock_begin(id_clock_sht_global_sum)

if (CS%reprod_sum) then
sum_tot = reproducing_sum(CS%Snm_Re_raw(:,:,1:Ltot), sums=Snm_Re(1:Ltot))
sum_tot = reproducing_sum(CS%Snm_Im_raw(:,:,1:Ltot), sums=Snm_Im(1:Ltot))
if (scale /= 1.0) then
I_scale = 1.0 / scale
do l=1,Ltot
Snm_Re(l) = I_scale * Snm_Re(l)
Snm_Im(l) = I_scale * Snm_Im(l)
enddo
endif
sum_tot = reproducing_sum(CS%Snm_Re_raw(:,:,1:Ltot), sums=Snm_Re(1:Ltot), unscale=tmp_scale)
sum_tot = reproducing_sum(CS%Snm_Im_raw(:,:,1:Ltot), sums=Snm_Im(1:Ltot), unscale=tmp_scale)
else
call sum_across_PEs(Snm_Re, Ltot)
call sum_across_PEs(Snm_Im, Ltot)
Expand Down

0 comments on commit 2c26118

Please sign in to comment.