da_vtovv_spectral.inc
References to this file elsewhere.
1 subroutine da_vtovv_spectral(max_wave, sizec, lenr, lenwrk, lensav, inc, &
2 alp_size, alp, wsave, power, rcv, field)
3
4 !-------------------------------------------------------------------------
5 ! Purpose: Performs spectral to gridpoint transformation on a sphere.
6 !----------------------------------------------------------------------
7
8 implicit none
9
10 integer, intent(in) :: max_wave ! Max total wavenumber.
11 integer, intent(in) :: sizec ! Size of packed spectral array.
12 integer, intent(in) :: lenr ! FFT info.
13 integer, intent(in) :: lenwrk ! FFT info.
14 integer, intent(in) :: lensav ! FFT info.
15 integer, intent(in) :: inc ! FFT info.
16 integer, intent(in) :: alp_size ! Size of alp array.
17 real, intent(in) :: alp(1:alp_size) ! Associated Legendre Polynomials.
18 real, intent(in) :: wsave(1:lensav) ! Primes for FFT.
19 real, intent(in) :: power(0:max_wave) ! Power spectrum
20 real, intent(in) :: rcv(1:2*sizec) ! Spectral modes.
21 real, intent(out) :: field(its:ite,jts:jte) ! Gridpoint field.
22
23 integer :: j, l,m, n ! Loop counters.
24 integer :: index_start ! Position markers in cv.
25 integer :: index_r, index_c ! Array index for complex v_fou.
26
27 real :: r_fou(1:lenr) ! FFT array.
28 real :: work(1:lenwrk) ! FFT work array.
29 complex :: v_fou(its:ite,0:max_wave)! Intermediate Fourier state.
30 complex :: ccv(1:sizec) ! Spectral modes.
31 integer :: index_m, index_j
32 integer :: jc, js, je, iequator
33 complex :: sum_legtra ! Summation scalars.
34
35
36 if (trace_use) call da_trace_entry("da_vtovv_spectral")
37
38 !----------------------------------------------------------------------------
39 ! [1] Create complex array from read array:
40 !----------------------------------------------------------------------------
41
42 v_fou = 0.0
43 do n = 1, sizec
44 ccv(n) = CMPLX(rcv(2*n-1), rcv(2*n))
45 end do
46
47 !----------------------------------------------------------------------------
48 ! [2] Apply power spectrum
49 !----------------------------------------------------------------------------
50
51 if (.not. test_transforms) call da_apply_power(power, max_wave, ccv, sizec)
52
53 !----------------------------------------------------------------------------
54 ! [3] Perform inverse Legendre decomposition in N-S direction:
55 !----------------------------------------------------------------------------
56
57 do m = 0, max_wave
58 index_start = m * (max_wave + 1 - m) + m * (m + 1) / 2 + 1
59
60 index_m = m * (max_wave + 1 - m) + m * (m + 1) / 2 + 1 - m
61
62 jc = (jde-jds+1)/2
63
64 iequator = mod(jde-jds+1, 2)
65
66 je = min(jc+iequator, jte)
67
68 do j = jts, je
69 index_j = (j - 1) * (max_wave + 1) * (max_wave + 2) / 2
70
71 v_fou(j,m) = sum(ccv(index_start:index_start-m+max_wave) * &
72 alp(index_j+index_m+m:index_j+index_m+max_wave))
73 end do
74
75 js = max(jts, jc+iequator+1)
76
77 do j = js, jte
78 index_j = (jds+jde - j - 1) * (max_wave + 1) * (max_wave + 2) / 2
79
80 sum_legtra = da_zero_complex
81 do l = m, max_wave
82 ! Calculate second quadrant values:
83 if(mod(l+m,2) == 1) then
84 sum_legtra = sum_legtra - ccv(index_start-m+l) * alp(index_j + index_m + l)
85 else
86 sum_legtra = sum_legtra + ccv(index_start-m+l) * alp(index_j + index_m + l)
87 end if
88 end do
89 v_fou(j,m) = sum_legtra
90 end do
91 end do
92
93 !----------------------------------------------------------------------------
94 ! [4] Perform inverse Fourier decomposition in E-W direction:
95 !----------------------------------------------------------------------------
96
97 do j = jts, jte
98 r_fou(its) = real(v_fou(j,0)) ! R(m=0) is real.
99 ! r_fou(ite) = aimag(v_fou(j,0)) ! R(m=NI/2) is real, but packed in imag m = 0)
100 ! make r_fou(ide) zero as there is no power computed corresponding to this wavenumber
101 r_fou(ite) = 0.0
102
103 do m = 1, max_wave
104 index_r = 2 * m
105 index_c = 2 * m + 1
106 r_fou(index_r) = real(v_fou(j,m))
107 r_fou(index_c) = aimag(v_fou(j,m))
108 end do
109
110 #ifdef FFTPACK
111 call rfft1b(ide, inc, r_fou, lenr, wsave, lensav, work, lenwrk, ierr)
112 #else
113 call da_error(__FILE__,__LINE__,(/"Must compile with FFTPACK"/))
114 #endif
115 field(its:ite,j) = r_fou(its:ite)
116 end do
117
118 if (trace_use) call da_trace_exit("da_vtovv_spectral")
119
120 end subroutine da_vtovv_spectral
121
122