NCEPLIBS-sp  2.3.3
spfft1.f
Go to the documentation of this file.
1 C> @file
2 C>
3 C> Perform multiple fast fourier transforms
4 C> @author IREDELL @date 96-02-20
5 
6 C> This subprogram performs multiple fast fourier transforms
7 C> between complex amplitudes in fourier space and real values
8 C> in cyclic physical space.
9 C> Subprogram SPFFT1 initializes trigonometric data each call.
10 C> Use subprogram SPFFT to save time and initialize once.
11 C> this version invokes the ibm essl fft.
12 C>
13 C> @param IMAX - INTEGER NUMBER OF VALUES IN THE CYCLIC> PHYSICAL SPACE
14 C> (SEE LIMITATIONS ON IMAX IN REMARKS BELOW.)
15 C> @param INCW - INTEGER FIRST DIMENSION OF THE COMPLEX AMPLITUDE ARRAY
16 C> (INCW >= IMAX/2+1)
17 C> @param INCG - INTEGER FIRST DIMENSION OF THE REAL VALUE ARRAY
18 C> (INCG >= IMAX)
19 C> @param KMAX - INTEGER NUMBER OF TRANSFORMS TO PERFORM
20 C> @param[out] W - COMPLEX(INCW,KMAX) COMPLEX AMPLITUDES IF IDIR>0
21 C> @param[out] G - REAL(INCG,KMAX) REAL VALUES IF IDIR<0
22 C> @param IDIR - INTEGER DIRECTION FLAG
23 C> IDIR>0 TO TRANSFORM FROM FOURIER TO PHYSICAL SPACE
24 C> IDIR<0 TO TRANSFORM FROM PHYSICAL TO FOURIER SPACE
25 C>
26 C> SUBPROGRAMS CALLED:
27 C> - SCRFT IBM ESSL COMPLEX TO REAL FOURIER TRANSFORM
28 C> - SRCFT IBM ESSL REAL TO COMPLEX FOURIER TRANSFORM
29 C>
30 C> @note The restrictions on IMAX are that it must be a multiple of 1
31 C> to 25 factors of two, up to 2 factors of three, and up to 1 factor of
32 C> five, seven and eleven.
33 C>
34 C> @note This subprogram is thread-safe.
35  SUBROUTINE spfft1(IMAX,INCW,INCG,KMAX,W,G,IDIR)
36  IMPLICIT NONE
37  INTEGER,INTENT(IN):: IMAX,INCW,INCG,KMAX,IDIR
38  COMPLEX,INTENT(INOUT):: W(INCW,KMAX)
39  REAL,INTENT(INOUT):: G(INCG,KMAX)
40  REAL:: AUX1(25000+INT(0.82*IMAX))
41  REAL:: AUX2(20000+INT(0.57*IMAX))
42  INTEGER:: NAUX1,NAUX2
43 C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44  naux1=25000+int(0.82*imax)
45  naux2=20000+int(0.57*imax)
46 C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47 C FOURIER TO PHYSICAL TRANSFORM.
48  SELECT CASE(idir)
49  CASE(1:)
50  CALL scrft(1,w,incw,g,incg,imax,kmax,-1,1.,
51  & aux1,naux1,aux2,naux2,0.,0)
52  CALL scrft(0,w,incw,g,incg,imax,kmax,-1,1.,
53  & aux1,naux1,aux2,naux2,0.,0)
54 C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55 C PHYSICAL TO FOURIER TRANSFORM.
56  CASE(:-1)
57  CALL srcft(1,g,incg,w,incw,imax,kmax,+1,1./imax,
58  & aux1,naux1,aux2,naux2,0.,0)
59  CALL srcft(0,g,incg,w,incw,imax,kmax,+1,1./imax,
60  & aux1,naux1,aux2,naux2,0.,0)
61  END SELECT
62  END SUBROUTINE
spfft1
subroutine spfft1(IMAX, INCW, INCG, KMAX, W, G, IDIR)
This subprogram performs multiple fast fourier transforms between complex amplitudes in fourier space...
Definition: spfft1.f:36