Contents


NAME

     zfftdm - initialize  the  trigonometric  weight  and  factor
     tables  or  compute the one-dimensional inverse Fast Fourier
     Transform of a set of double complex data  sequences  stored
     in a two-dimensional array.

SYNOPSIS

     SUBROUTINE ZFFTDM(IOPT, N1, N2, SCALE, X, LDX, Y, LDY, TRIGS, IFAC, WORK, LWORK, IERR)

     INTEGER IOPT, N1, N2, LDX, LDY, IFAC(*), LWORK, IERR
     DOUBLE COMPLEX X(LDX, *)
     DOUBLE PRECISION SCALE, Y(LDY, *), TRIGS(*), WORK(*)

     SUBROUTINE ZFFTDM_64(IOPT, N1, N2, SCALE, X, LDX, Y, LDY, TRIGS, IFAC, WORK, LWORK, IERR)

     INTEGER*8 IOPT, N1, N2, LDX, LDY, IFAC(*), LWORK, IERR
     DOUBLE COMPLEX X(LDX, *)
     DOUBLE PRECISION SCALE, Y(LDY,*), TRIGS(*), WORK(*)

  F95 INTERFACE
     SUBROUTINE FFTM(IOPT, N1, [N2], [SCALE], X, [LDX], Y, [LDY], TRIGS,
               IFAC, WORK, [LWORK], IERR)

     INTEGER, INTENT(IN) :: IOPT, N1
     INTEGER, INTENT(IN), OPTIONAL :: N2, LDX, LDY, LWORK
     REAL(8), INTENT(IN), OPTIONAL :: SCALE
     COMPLEX(8), INTENT(IN), DIMENSION(:,:) :: X
     REAL(8), INTENT(OUT), DIMENSION(:,:) :: Y
     REAL(8), INTENT(INOUT), DIMENSION(:) :: TRIGS
     INTEGER, INTENT(INOUT), DIMENSION(:) :: IFAC
     REAL(8), INTENT(OUT), DIMENSION(:) :: WORK
     INTEGER, INTENT(OUT) :: IERR

     SUBROUTINE FFTM_64(IOPT, N1, [N2], [SCALE], X, [LDX], Y, [LDY], TRIGS, IFAC, WORK, [LWORK], IERR)

     INTEGER(8), INTENT(IN) :: IOPT, N1
     INTEGER(8), INTENT(IN), OPTIONAL :: N2, LDX, LDY, LWORK
     REAL(8), INTENT(IN), OPTIONAL :: SCALE
     COMPLEX(8), INTENT(IN), DIMENSION(:,:) :: X
     REAL(8), INTENT(OUT), DIMENSION(:,:) :: Y
     REAL(8), INTENT(INOUT), DIMENSION(:) :: TRIGS
     INTEGER(8), INTENT(INOUT), DIMENSION(:) :: IFAC
     REAL(8), INTENT(OUT), DIMENSION(:) :: WORK
     INTEGER(8), INTENT(OUT) :: IERR

  C INTERFACE
     #include <sunperf.h>
     void zfftdm_ (int *iopt, int  *m,  int  *n,  double  *scale,
               doublecomplex  *x,  int *ldx, double *y, int *ldy,
               double  *trigs,  int  *ifac,  double  *work,   int
               *lwork, int *ierr);

     void zfftdm_64_  (long  *iopt,  long  *m,  long  *n,  double
               *scale,  doublecomplex  *x,  long *ldx, double *y,
               long  *ldy,  double  *trigs,  long  *ifac,  double
               *work, long *lwork, long *ierr);

PURPOSE

     zfftdm  initializes  the  trigonometric  weight  and  factor
     tables  or computes the one-dimensional inverse Fast Fourier
     Transform of a set of double complex data  sequences  stored
     in a two-dimensional array:

                      N1-1
     Y(k,l) = scale * SUM  W*X(j,l)
                      j=0

     where
     k ranges from 0 to N1-1 and l ranges from 0 to N2-1
     i = sqrt(-1)
     isign = 1 for inverse transform
     W = exp(isign*i*j*k*2*pi/N1)
     In complex-to-real transform of length N1, the (N1/2+1) com-
     plex  input  data  points  stored are the positive-frequency
     half of the spectrum of the Discrete Fourier Transform.  The
     other  half  can be obtained through complex conjugation and
     therefore is not stored.  Furthermore, due to symmetries the
     imaginary of the component of X(0,0:N2-1) and X(N1/2,0:N2-1)
     (if N1 is even in the latter) is assumed to be zero  and  is
     not referenced.

ARGUMENTS

     IOPT (input)
               Integer specifying the operation to be performed:
               IOPT = 0 computes the trigonometric  weight  table
               and factor table
               IOPT = 1 computes inverse FFT

     N1 (input)
               Integer specifying length of the input  sequences.
               N1 is most efficient when it is a product of small
               primes.  N1 >= 0.  Unchanged on exit.

     N2 (input)
               Integer specifying number of input sequences.   N2
               >= 0.  Unchanged on exit.
     SCALE (input)
               Double precision scalar by which transform results
               are   scaled.    Unchanged   on  exit.   SCALE  is
               defaulted to 1.0D0 for F95 INTERFACE.

     X (input) X is a double complex array  of  dimensions  (LDX,
               N2)  that contains the sequences to be transformed
               stored in its columns in X(0:N1/2, 0:N2-1).

     LDX (input)
               Leading dimension of X.  LDX >= (N1/2+1) Unchanged
               on exit.

     Y (output)
               Y is a double precision array of dimensions  (LDY,
               N2)  that  contains  the  transform results of the
               input sequences in Y(0:N1-1,0:N2-1).  X and Y  can
               be  the  same  array  starting  at the same memory
               location, in which case the  input  sequences  are
               overwritten  by  their  transform results.  Other-
               wise, it is  assumed  that  there  is  no  overlap
               between X and Y in memory.

     LDY (input)
               Leading dimension of Y.  If X and Y are  the  same
               array,  LDY  =  2*LDX  Else LDY >= N1 Unchanged on
               exit.

     TRIGS (input/output)
               double precision array of length  2*N1  that  con-
               tains  the trigonometric weights.  The weights are
               computed when the routine is called with IOPT =  0
               and  they are used in subsequent calls when IOPT =
               1.  Unchanged on exit.

     IFAC (input/output)
               Integer array of dimension at least 128 that  con-
               tains the factors of N1.  The factors are computed
               when the routine is called with IOPT = 0 and  they
               are  used  in  subsequent  calls  when  IOPT  = 1.
               Unchanged on exit.

     WORK (workspace)
               double precision array of dimension at  least  N1.
               The user can also choose to have the routine allo-
               cate its own workspace (see LWORK).

     LWORK (input)
               Integer specifying workspace size.  If LWORK =  0,
               the routine will allocate its own workspace.

     IERR (output)
               On exit, integer IERR has  one  of  the  following
               values:
               0 = normal return
               -1 = IOPT is not 0 or 1
               -2 = N1 < 0
               -3 = N2 < 0
               -4 = (LDX < N1/2+1)
               -5 = (LDY < N1) or (LDY not equal 2*LDX when X and
               Y are same array)
               -6 = (LWORK not equal 0) and (LWORK < N1)
               -7 = memory allocation failed

SEE ALSO

     fft