Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / var / da / da_spectral / da_legtra_inv.inc
blob0c77fffa79d7f73bb49a3cecf2ca9b614baa2075
1 subroutine da_legtra_inv(jds, jde, jts, jte, max_wavenumber, alp_size, m, &
2    alp, v, r_leg)
4    !-----------------------------------------------------------------------
5    ! Purpose: TBD
6    !-----------------------------------------------------------------------
8    implicit none
10    integer, intent(in)  :: jds, jde            ! Number of latitudes.
11    integer, intent(in)  :: jts, jte            ! Number of latitudes.
12    integer, intent(in)  :: max_wavenumber      ! Maximum wavenumber.
13    integer, intent(in)  :: alp_size            ! Dimension of ALPs.
14    integer, intent(in)  :: m                   ! Zonal wavenumber.
15    real,    intent(in)  :: alp(1:alp_size)     ! Associated Legendre Polynomials
17    complex, intent(in)  :: v(m:max_wavenumber) ! Output spectral coefficient.
18    complex, intent(out) :: r_leg(jts:jte)         ! Field to transform.
20    integer              :: l, j, js, je        ! Loop counters.
21    integer              :: index_m, index_j
22    complex              :: sum_legtra          ! Summation scalars.
24    integer              :: jc, iequator
26    if (trace_use) call da_trace_entry("da_legtra_inv")
28    index_m = m * (max_wavenumber + 1 - m) + m * (m + 1) / 2 + 1 - m
30    jc = (jde-jds+1)/2
32    iequator = mod(jde-jds+1, 2)
34    je = min(jc+iequator, jte)
36    do j = jts, je
37       index_j = (j - 1) * (max_wavenumber + 1) * (max_wavenumber + 2) / 2
39       r_leg(j) = sum(v(m:max_wavenumber) * &
40          alp(index_j+index_m+m:index_j+index_m+max_wavenumber))
41    end do
43    js = max(jts, jc+iequator+1)
45    do j = js, jte
46       index_j = (jds+jde - j - 1) * (max_wavenumber + 1) * (max_wavenumber + 2) / 2
48       sum_legtra = da_zero_complex
49       do l = m, max_wavenumber
50          ! Calculate second quadrant values:
51          if(mod(l+m,2) == 1) then
52             sum_legtra = sum_legtra - v(l) * alp(index_j + index_m + l)
53          else
54             sum_legtra = sum_legtra + v(l) * alp(index_j + index_m + l)
55          end if
56       end do
57       r_leg(j) = sum_legtra
58    end do
60    if (trace_use) call da_trace_exit("da_legtra_inv")
62 end subroutine da_legtra_inv