NCEPLIBS-ip 5.2.0
Loading...
Searching...
No Matches
spfft.f
Go to the documentation of this file.
1C> @file
2C> @brief Perform multiple fast fourier transforms.
3C> @author Iredell @date 96-02-20
4
5C> This subprogram performs multiple fast fourier transforms
6C> between complex amplitudes in fourier space and real values
7C> in cyclic physical space.
8C>
9C> Subprogram spfft must be invoked first with idir=0
10C> to initialize trigonemetric data. Use subprogram spfft1
11C> to perform an fft without previous initialization.
12C> This version invokes the ibm essl fft.
13C>
14C> The restrictions on imax are that it must be a multiple
15C> of 1 to 25 factors of two, up to 2 factors of three,
16C> and up to 1 factor of five, seven and eleven.
17C>
18C> If IDIR=0, then W and G need not contain any valid data.
19C> the other parameters must be supplied and cannot change
20C> in succeeding calls until the next time it is called with IDIR=0.
21C>
22C> This subprogram is not thread-safe when IDIR=0. On the other hand,
23C> when IDIR is not zero, it can be called from a threaded region.
24C>
25C> @param IMAX number of values in the cyclic physical space
26C> (see limitations on imax in remarks below.)
27C> @param INCW first dimension of the complex amplitude array
28C> (INCW >= IMAX/2+1)
29C> @param INCG first dimension of the real value array
30C> (INCG >= IMAX)
31C> @param KMAX number of transforms to perform
32C> @param[out] W complex amplitudes if IDIR>0
33C> @param[out] G real values if IDIR<0
34C> @param IDIR direction flag
35C> - IDIR=0 to initialize internal trigonometric data
36C> - IDIR>0 TO transform from Fourier to physical space
37C> - IDIR<0 TO transform from physical to fourier space
38C>
39C> @author Iredell @date 96-02-20
40 SUBROUTINE spfft(IMAX,INCW,INCG,KMAX,W,G,IDIR)
41
42 IMPLICIT NONE
43 INTEGER,INTENT(IN):: IMAX,INCW,INCG,KMAX,IDIR
44 COMPLEX,INTENT(INOUT):: W(INCW,KMAX)
45 REAL:: WREAL(INCW,KMAX)
46 REAL,INTENT(INOUT):: G(INCG,KMAX)
47 INTEGER,SAVE:: NAUX1=0
48 REAL,SAVE,ALLOCATABLE:: AUX1CR(:),AUX1RC(:)
49 INTEGER:: NAUX2
50 REAL:: AUX2(20000+INT(0.57*IMAX))
51
52 naux2=20000+int(0.57*imax)
53
54 wreal=real(w)
55C INITIALIZATION.
56C ALLOCATE AND FILL AUXILIARY ARRAYS WITH TRIGONOMETRIC DATA
57 SELECT CASE(idir)
58 CASE(0)
59 IF(naux1.GT.0) DEALLOCATE(aux1cr,aux1rc)
60 naux1=25000+int(0.82*imax)
61 ALLOCATE(aux1cr(naux1),aux1rc(naux1))
62 CALL scrft(1,wreal,incw,g,incg,imax,kmax,-1,1.,
63 & aux1cr,naux1,aux2,naux2,0.,0)
64 CALL srcft(1,g,incg,wreal,incw,imax,kmax,+1,1./imax,
65 & aux1rc,naux1,aux2,naux2,0.,0)
66
67C FOURIER TO PHYSICAL TRANSFORM.
68 CASE(1:)
69 CALL scrft(0,wreal,incw,g,incg,imax,kmax,-1,1.,
70 & aux1cr,naux1,aux2,naux2,0.,0)
71
72C PHYSICAL TO FOURIER TRANSFORM.
73 CASE(:-1)
74 CALL srcft(0,g,incg,wreal,incw,imax,kmax,+1,1./imax,
75 & aux1rc,naux1,aux2,naux2,0.,0)
76 END SELECT
77 w=cmplx(wreal)
78 END SUBROUTINE
subroutine scrft(init, x, ldx, y, ldy, n, m, isign, scale, table, n1, wrk, n2, z, nz)
scrft
Definition fftpack.F:82
subroutine srcft(init, x, ldx, y, ldy, n, m, isign, scale, table, n1, wrk, n2, z, nz)
srcft
Definition fftpack.F:215
subroutine spfft(imax, incw, incg, kmax, w, g, idir)
This subprogram performs multiple fast fourier transforms between complex amplitudes in fourier space...
Definition spfft.f:41