updated top-level README and version_decl for V4.4.2 (#1795)
[WRF.git] / external / fftpack / fftpack5 / mradf2.F
blob7ab10c92e17648e3aa74028292c04f2e4b193d4d
1 subroutine mradf2 ( m, ido, l1, cc, im1, in1, ch, im2, in2, wa1 )
3 !*****************************************************************************80
5 !! MRADF2 is an FFTPACK5 auxiliary routine.
8 !    Copyright (C) 1995-2004, Scientific Computing Division,
9 !    University Corporation for Atmospheric Research
11 !  Modified:
13 !    27 March 2009
15 !  Author:
17 !    Paul Swarztrauber
18 !    Richard Valent
20 !  Reference:
22 !    Paul Swarztrauber,
23 !    Vectorizing the Fast Fourier Transforms,
24 !    in Parallel Computations,
25 !    edited by G. Rodrigue,
26 !    Academic Press, 1982.
28 !    Paul Swarztrauber,
29 !    Fast Fourier Transform Algorithms for Vector Computers,
30 !    Parallel Computing, pages 45-63, 1984.
32 !  Parameters:
34   implicit none
36   integer ( kind = 4 ) ido
37   integer ( kind = 4 ) in1
38   integer ( kind = 4 ) in2
39   integer ( kind = 4 ) l1
41   real ( kind = 4 ) cc(in1,ido,l1,2)
42   real ( kind = 4 ) ch(in2,ido,2,l1)
43   integer ( kind = 4 ) i
44   integer ( kind = 4 ) ic
45   integer ( kind = 4 ) idp2
46   integer ( kind = 4 ) im1
47   integer ( kind = 4 ) im2
48   integer ( kind = 4 ) k
49   integer ( kind = 4 ) m
50   integer ( kind = 4 ) m1
51   integer ( kind = 4 ) m1d
52   integer ( kind = 4 ) m2
53   integer ( kind = 4 ) m2s
54   real ( kind = 4 ) wa1(ido)
56   m1d = ( m - 1 ) * im1 + 1
57   m2s = 1 - im2
59   do k = 1, l1
60     m2 = m2s
61     do m1 = 1, m1d, im1
62       m2 = m2 + im2
63       ch(m2,1,1,k)   = cc(m1,1,k,1) + cc(m1,1,k,2)
64       ch(m2,ido,2,k) = cc(m1,1,k,1) - cc(m1,1,k,2)
65     end do
66   end do
68   if ( ido < 2 ) then
69     return
70   end if
72   if ( 2 < ido ) then
74     idp2 = ido + 2
76     do k = 1, l1
77       do i = 3, ido, 2
78         ic = idp2 - i
79         m2 = m2s
80         do m1 = 1, m1d, im1
81           m2 = m2 + im2
82           ch(m2,i,1,k) =    cc(m1,i,k,1)   + ( wa1(i-2) * cc(m1,i,k,2) &
83                                              - wa1(i-1) * cc(m1,i-1,k,2) )
84           ch(m2,ic,2,k)  = -cc(m1,i,k,1)   + ( wa1(i-2) * cc(m1,i,k,2) &
85                                              - wa1(i-1) * cc(m1,i-1,k,2) )
86           ch(m2,i-1,1,k)  = cc(m1,i-1,k,1) + ( wa1(i-2) * cc(m1,i-1,k,2) &
87                                              + wa1(i-1) * cc(m1,i,k,2))
88           ch(m2,ic-1,2,k) = cc(m1,i-1,k,1) - ( wa1(i-2) * cc(m1,i-1,k,2) &
89                                              + wa1(i-1) * cc(m1,i,k,2) )
90         end do
91       end do
92     end do
94     if ( mod ( ido, 2 ) == 1 ) then
95       return
96     end if
98   end if
100   do k = 1, l1
101     m2 = m2s
102     do m1 = 1, m1d, im1
103       m2 = m2 + im2
104       ch(m2,1,2,k) = -cc(m1,ido,k,2)
105       ch(m2,ido,1,k) = cc(m1,ido,k,1)
106     end do
107   end do
109   return