da_legtra_inv.inc
References to this file elsewhere.
1 subroutine da_legtra_inv(jds, jde, jts, jte, max_wavenumber, alp_size, m, &
2 alp, v, r_leg)
3
4 !-----------------------------------------------------------------------
5 ! Purpose: TBD
6 !-----------------------------------------------------------------------
7
8 implicit none
9
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
16
17 complex, intent(in) :: v(m:max_wavenumber) ! Output spectral coefficient.
18 complex, intent(out) :: r_leg(jts:jte) ! Field to transform.
19
20 integer :: l, j, js, je ! Loop counters.
21 integer :: index_m, index_j
22 complex :: sum_legtra ! Summation scalars.
23
24 integer :: jc, iequator
25
26 index_m = m * (max_wavenumber + 1 - m) + m * (m + 1) / 2 + 1 - m
27
28 jc = (jde-jds+1)/2
29
30 iequator = mod(jde-jds+1, 2)
31
32 je = min(jc+iequator, jte)
33
34 do j = jts, je
35 index_j = (j - 1) * (max_wavenumber + 1) * (max_wavenumber + 2) / 2
36
37 r_leg(j) = sum(v(m:max_wavenumber) * &
38 alp(index_j+index_m+m:index_j+index_m+max_wavenumber))
39 end do
40
41 js = max(jts, jc+iequator+1)
42
43 do j = js, jte
44 index_j = (jds+jde - j - 1) * (max_wavenumber + 1) * (max_wavenumber + 2) / 2
45
46 sum_legtra = da_zero_complex
47 do l = m, max_wavenumber
48 ! Calculate second quadrant values:
49 if(mod(l+m,2) == 1) then
50 sum_legtra = sum_legtra - v(l) * alp(index_j + index_m + l)
51 else
52 sum_legtra = sum_legtra + v(l) * alp(index_j + index_m + l)
53 end if
54 end do
55 r_leg(j) = sum_legtra
56 end do
57
58 end subroutine da_legtra_inv
59
60