Go to main content
Oracle Developer Studio 12.5 Man Pages

Exit Print View

Updated: June 2017
 
 

zfftd3 (3p)

Name

zfftd3 - pute the three-dimensional inverse Fast Fourier Transform of a three- dimensional double complex array.

Synopsis

SUBROUTINE ZFFTD3(IOPT, M, N, K, SCALE, X, LDX1, LDX2, Y, LDY1, LDY2, TRIGS, IFAC, WORK, LWORK, IERR)

INTEGER IOPT, M, N, K, LDX1, LDX2, LDY1, LDY2, IFAC(*), LWORK, IERR
DOUBLE COMPLEX X(LDX1, LDX2, *)
DOUBLE PRECISION SCALE, TRIGS(*), WORK(*), Y(LDY1, LDY2, *)

SUBROUTINE ZFFTD3_64(IOPT, M, N, K, SCALE, X, LDX1, LDX2, Y, LDY1, LDY2, TRIGS, IFAC, WORK, LWORK, IERR)

INTEGER*8 IOPT, M, N, K, LDX1, LDX2, LDY1, LDY2, IFAC(*), LWORK, IERR
DOUBLE COMPLEX X(LDX1, LDX2, *)
DOUBLE PRECISION SCALE, TRIGS(*), WORK(*), Y(LDY1, LDY2, *)




F95 INTERFACE
SUBROUTINE FFT3(IOPT, M, N, K, SCALE, X, LDX1, LDX2, Y, LDY1, LDY2, TRIGS, IFAC, WORK, LWORK, IERR)

INTEGER, INTENT(IN) :: IOPT, M, LDX2, LDY2
INTEGER, INTENT(IN), OPTIONAL :: N, K, LDX1, LDY1, 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 FFT3_64(IOPT, M, N, K, SCALE, X, LDX1, LDX2, Y, LDY1, LDY2, TRIGS, IFAC, WORK, LWORK, IERR)

INTEGER(8), INTENT(IN) :: IOPT, M, LDX2, LDY2
INTEGER(8), INTENT(IN), OPTIONAL :: N, K, LDX1, LDY1, 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 zfftd3_ (int *iopt, int *n1, int *n2, int *n3, double *scale, dou-
blecomplex *x, int *ldx1, int *ldx2, double  *y,  int  *ldy1,
int  *ldy2,  double  *trigs,  int  *ifac,  double  *work, int
*lwork, int *ierr);

void zfftd3_64_ (long *iopt, long  *n1,  long  *n2,  long  *n3,  double
*scale,  doublecomplex *x, long *ldx1, long *ldx2, double *y,
long *ldy1, long *ldy2, double  *trigs,  long  *ifac,  double
*work, long *lwork, long *ierr);

Description

Oracle Solaris Studio Performance Library                           zfftd3(3P)



NAME
       zfftd3  - initialize the trigonometric weight and factor tables or com-
       pute the three-dimensional inverse Fast Fourier Transform of  a  three-
       dimensional double complex array.

SYNOPSIS
       SUBROUTINE ZFFTD3(IOPT, M, N, K, SCALE, X, LDX1, LDX2, Y, LDY1, LDY2, TRIGS, IFAC, WORK, LWORK, IERR)

       INTEGER IOPT, M, N, K, LDX1, LDX2, LDY1, LDY2, IFAC(*), LWORK, IERR
       DOUBLE COMPLEX X(LDX1, LDX2, *)
       DOUBLE PRECISION SCALE, TRIGS(*), WORK(*), Y(LDY1, LDY2, *)

       SUBROUTINE ZFFTD3_64(IOPT, M, N, K, SCALE, X, LDX1, LDX2, Y, LDY1, LDY2, TRIGS, IFAC, WORK, LWORK, IERR)

       INTEGER*8 IOPT, M, N, K, LDX1, LDX2, LDY1, LDY2, IFAC(*), LWORK, IERR
       DOUBLE COMPLEX X(LDX1, LDX2, *)
       DOUBLE PRECISION SCALE, TRIGS(*), WORK(*), Y(LDY1, LDY2, *)




   F95 INTERFACE
       SUBROUTINE FFT3(IOPT, M, N, K, SCALE, X, LDX1, LDX2, Y, LDY1, LDY2, TRIGS, IFAC, WORK, LWORK, IERR)

       INTEGER, INTENT(IN) :: IOPT, M, LDX2, LDY2
       INTEGER, INTENT(IN), OPTIONAL :: N, K, LDX1, LDY1, 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 FFT3_64(IOPT, M, N, K, SCALE, X, LDX1, LDX2, Y, LDY1, LDY2, TRIGS, IFAC, WORK, LWORK, IERR)

       INTEGER(8), INTENT(IN) :: IOPT, M, LDX2, LDY2
       INTEGER(8), INTENT(IN), OPTIONAL :: N, K, LDX1, LDY1, 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 zfftd3_ (int *iopt, int *n1, int *n2, int *n3, double *scale, dou-
                 blecomplex *x, int *ldx1, int *ldx2, double  *y,  int  *ldy1,
                 int  *ldy2,  double  *trigs,  int  *ifac,  double  *work, int
                 *lwork, int *ierr);

       void zfftd3_64_ (long *iopt, long  *n1,  long  *n2,  long  *n3,  double
                 *scale,  doublecomplex *x, long *ldx1, long *ldx2, double *y,
                 long *ldy1, long *ldy2, double  *trigs,  long  *ifac,  double
                 *work, long *lwork, long *ierr);



PURPOSE
       zfftd3  initializes  the trigonometric weight and factor tables or com-
       putes  the  three-dimensional  inverse  Fast  Fourier  Transform  of  a
       three-dimensional double complex array.

                             K-1  N-1  M-1
       Y(k1,k2,k3) = scale * SUM   SUM   SUM   W3*W2*W1*X(j1,j2,j3)
                             j3=0  j2=0  j1=0

       where
       k1  ranges  from 0 to M-1; k2 ranges from 0 to N-1 and k3 ranges from 0
       to K-1
       i = sqrt(-1)
       isign = 1 for inverse transform
       W1 = exp(isign*i*j1*k1*2*pi/M)
       W2 = exp(isign*i*j2*k2*2*pi/N)
       W3 = exp(isign*i*j3*k3*2*pi/K)


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

       M (input)
                 Integer  specifying  length  of  the  transform  in the first
                 dimension.  M is most efficient when it is a product of small
                 primes.  M >= 0.  Unchanged on exit.

       N (input)
                 Integer  specifying  length  of  the  transform in the second
                 dimension.  N is most efficient when it is a product of small
                 primes.  N >= 0.  Unchanged on exit.

       K (input)
                 Integer  specifying  length  of  the  transform  in the third
                 dimension.  K is most efficient when it is a product of small
                 primes.  K >= 0.  Unchanged on exit.

       SCALE (input)
                 Double  precision  scalar  by  which  transform  results  are
                 scaled.  Unchanged on exit.

       X (input) X is a double complex array of  dimensions  (LDX1,  LDX2,  K)
                 that contains input data to be transformed.

       LDX1 (input)
                 first dimension of X.  LDX1 >= M/2+1 Unchanged on exit.

       LDX2 (input)
                 second dimension of X.  LDX2 >= N Unchanged on exit.

       Y (output)
                 Y  is  a double precision array of dimensions (LDY1, LDY2, K)
                 that contains the transform results.  X and Y can be the same
                 array starting at the same memory location, in which case the
                 input data are overwritten by their transform results.   Oth-
                 erwise,  it is assumed that there is no overlap between X and
                 Y in memory.

       LDY1 (input)
                 first dimension of Y.  If X and Y are the same array, LDY1  =
                 2*LDX1  Else  LDY1  >=  2*LDX1  and LDY1 is even Unchanged on
                 exit.

       LDY2 (input)
                 second dimension of Y.  If X and Y are the same array, LDY2 =
                 LDX2 Else LDY2 >= N Unchanged on exit.

       TRIGS (input/output)
                 Double  precision array of length 2*(M+N+K) that contains the
                 trigonometric weights.  The weights  are  computed  when  the
                 routine  is  called with IOPT = 0 and they are used in subse-
                 quent calls when IOPT = 1.  Unchanged on exit.

       IFAC (input/output)
                 Integer array of dimension at least 3*128 that  contains  the
                 factors  of  M,  N  and K.  The factors are computed when the
                 routine is called with IOPT = 0 and they are used  in  subse-
                 quent calls when IOPT = 1.  Unchanged on exit.

       WORK (workspace)
                 Double  precision array of dimension at least (MAX(N,2*N,2*K)
                 + 16*K) * NCPUS where NCPUS is the number of threads used  to
                 execute  the  routine.   The user can also choose to have the
                 routine allocate 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 = M < 0
                 -3 = N < 0
                 -4 = K < 0
                 -5 = (LDX1 < M/2+1)
                 -6 = (LDX2 < N)
                 -7 = LDY1 not equal 2*LDX1 when X and Y are same array
                 -8  =  (LDY1  < 2*LDX1) or (LDY1 is odd) when X and Y are not
                 same array
                 -9 = (LDY2 < N) or (LDY2 not equal LDX2) when  X  and  Y  are
                 same array
                 -10  =  (LWORK  not  equal  0) and ((LWORK < MAX(N,2*N,2*K) +
                 16*K)*NCPUS)
                 -11 = memory allocation failed

SEE ALSO
       fft

CAUTIONS
       On exit, output subarray Y(1:LDY1, 1:N, 1:K) is overwritten.



                                  7 Nov 2015                        zfftd3(3P)