NCEPLIBS-sp 2.4.0
spffte.f
Go to the documentation of this file.
1C> @file
2C> @brief Perform multiple fast Fourier transforms.
3C>
4C> ### Program History Log
5C> Date | Programmer | Comments
6C> -----|------------|---------
7C> 1998-12-18 | Iredell | Initial.
8C> 2012-11-12 | Mirvis | fixing hard-wired types problem on Intel/Linux.
9C>
10C> @author Iredell @date 96-02-20
11
12C> This subprogram performs multiple fast Fourier transforms
13C> between complex amplitudes in Fourier space and real values
14C> in cyclic physical space.
15C>
16C> This subprogram must be invoked first with IDIR=0
17C> to initialize trigonemetric data. Use subprogram spfft1()
18C> to perform an FFT without previous initialization.
19C>
20C> This version invokes the IBM ESSL FFT.
21C>
22C> @note The restrictions on IMAX are that it must be a multiple
23C> of 1 to 25 factors of two, up to 2 factors of three,
24C> and up to 1 factor of five, seven and eleven.
25C>
26C> If IDIR=0, then W and G need not contain any valid data.
27C> The other parameters must be supplied and cannot change
28C> in succeeding calls until the next time it is called with IDIR=0.
29C>
30C> This subprogram is thread-safe.
31C>
32C> @param IMAX number of values in the cyclic physical space
33C> (see limitations on imax in remarks below.)
34C> @param INCW first dimension of the complex amplitude array
35C> (INCW >= IMAX/2+1)
36C> @param INCG first dimension of the real value array
37C> (INCG >= IMAX)
38C> @param KMAX number of transforms to perform
39C> @param[out] W complex amplitudes if IDIR>0
40C> @param[out] G real values if IDIR<0
41C> @param IDIR direction flag
42C> - IDIR=0 to initialize trigonometric data
43C> - IDIR>0 to transform from Fourier to physical space
44C> - IDIR<0 to transform from physical to Fourier space
45C> @param[out] AFFT auxiliary array if IDIR<>0
46C>
47C> @author Iredell @date 96-02-20
48 SUBROUTINE spffte(IMAX,INCW,INCG,KMAX,W,G,IDIR,AFFT)
49 IMPLICIT NONE
50 INTEGER,INTENT(IN):: IMAX,INCW,INCG,KMAX,IDIR
51 REAL,INTENT(INOUT):: W(2*INCW,KMAX)
52 REAL,INTENT(INOUT):: G(INCG,KMAX)
53 REAL(8),INTENT(INOUT):: AFFT(50000+4*IMAX)
54 INTEGER:: INIT,INC2X,INC2Y,N,M,ISIGN,NAUX1,NAUX2,NAUX3
55C ==EM== ^(4)
56 REAL:: SCALE
57 REAL(8):: AUX2(20000+2*IMAX),AUX3
58 INTEGER:: IACR,IARC
59
60 naux1=25000+2*imax
61 naux2=20000+2*imax
62 naux3=1
63 iacr=1
64 iarc=1+naux1
65
66C INITIALIZATION.
67C FILL AUXILIARY ARRAYS WITH TRIGONOMETRIC DATA
68 SELECT CASE(idir)
69 CASE(0)
70 init=1
71 inc2x=incw
72 inc2y=incg
73 n=imax
74 m=kmax
75 isign=-1
76 scale=1.
77 IF(digits(1.).LT.digits(1._8)) THEN
78 CALL scrft(init,w,inc2x,g,inc2y,n,m,isign,scale,
79 & afft(iacr),naux1,aux2,naux2,aux3,naux3)
80 ELSE
81 CALL dcrft(init,w,inc2x,g,inc2y,n,m,isign,scale,
82 & afft(iacr),naux1,aux2,naux2)
83 ENDIF
84 init=1
85 inc2x=incg
86 inc2y=incw
87 n=imax
88 m=kmax
89 isign=+1
90 scale=1./imax
91 IF(digits(1.).LT.digits(1._8)) THEN
92 CALL srcft(init,g,inc2x,w,inc2y,n,m,isign,scale,
93 & afft(iarc),naux1,aux2,naux2,aux3,naux3)
94 ELSE
95 CALL drcft(init,g,inc2x,w,inc2y,n,m,isign,scale,
96 & afft(iarc),naux1,aux2,naux2)
97 ENDIF
98
99C FOURIER TO PHYSICAL TRANSFORM.
100 CASE(1:)
101 init=0
102 inc2x=incw
103 inc2y=incg
104 n=imax
105 m=kmax
106 isign=-1
107 scale=1.
108 IF(digits(1.).LT.digits(1._8)) THEN
109 CALL scrft(init,w,inc2x,g,inc2y,n,m,isign,scale,
110 & afft(iacr),naux1,aux2,naux2,aux3,naux3)
111 ELSE
112 CALL dcrft(init,w,inc2x,g,inc2y,n,m,isign,scale,
113 & afft(iacr),naux1,aux2,naux2)
114 ENDIF
115
116C PHYSICAL TO FOURIER TRANSFORM.
117 CASE(:-1)
118 init=0
119 inc2x=incg
120 inc2y=incw
121 n=imax
122 m=kmax
123 isign=+1
124 scale=1./imax
125 IF(digits(1.).LT.digits(1._8)) THEN
126 CALL srcft(init,g,inc2x,w,inc2y,n,m,isign,scale,
127 & afft(iarc),naux1,aux2,naux2,aux3,naux3)
128 ELSE
129 CALL drcft(init,g,inc2x,w,inc2y,n,m,isign,scale,
130 & afft(iarc),naux1,aux2,naux2)
131 ENDIF
132 END SELECT
133 END SUBROUTINE
subroutine drcft(init, x, ldx, y, ldy, n, m, isign, scale, table, n1, wrk, n2, z, nz)
drcft
Definition: fftpack.F:160
subroutine dcrft(init, x, ldx, y, ldy, n, m, isign, scale, table, n1, wrk, n2, z, nz)
dcrft
Definition: fftpack.F:34
subroutine scrft(init, x, ldx, y, ldy, n, m, isign, scale, table, n1, wrk, n2, z, nz)
scrft
Definition: fftpack.F:80
subroutine srcft(init, x, ldx, y, ldy, n, m, isign, scale, table, n1, wrk, n2, z, nz)
srcft
Definition: fftpack.F:209
subroutine spffte(IMAX, INCW, INCG, KMAX, W, G, IDIR, AFFT)
This subprogram performs multiple fast Fourier transforms between complex amplitudes in Fourier space...
Definition: spffte.f:49