Go to main content
Oracle Developer Studio 12.5 Man Pages

Exit Print View

Updated: June 2017
 
 

zfftd2 (3p)

Name

zfftd2 - pute the two-dimensional inverse Fast Fourier Transform of a two-dimen- sional double complex array.

Synopsis

SUBROUTINE ZFFTD2(IOPT, M, N, SCALE, X, LDX, Y, LDY, TRIGS, IFAC, WORK, LWORK, IERR)

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

SUBROUTINE ZFFTD2_64(IOPT, M, N, SCALE, X, LDX, Y, LDY, TRIGS, IFAC, WORK, LWORK, IERR)

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




F95 INTERFACE
SUBROUTINE FFT2(IOPT, M, N, SCALE, X, LDX, Y, LDY, TRIGS,
IFAC, WORK, LWORK, IERR)

INTEGER, INTENT(IN) :: IOPT, M
INTEGER, INTENT(IN), OPTIONAL :: N, 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 FFT2_64(IOPT, M, N, SCALE, X, LDX, Y, LDY, TRIGS, IFAC, WORK, LWORK, IERR)

INTEGER(8), INTENT(IN) :: IOPT, M
INTEGER(8), INTENT(IN), OPTIONAL :: N, 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 zfftd2_ (int *iopt, int *n1, int *n2, double *scale, doublecomplex
*x, int *ldx, double *y, int *ldy, double *trigs, int  *ifac,
double *work, int *lwork, int *ierr);

void zfftd2_64_ (long *iopt, long *n1, long *n2, double *scale, double-
complex *x, long *ldx, double *y, long *ldy,  double  *trigs,
long *ifac, double *work, long *lwork, long *ierr);

Description

Oracle Solaris Studio Performance Library                           zfftd2(3P)



NAME
       zfftd2  - initialize the trigonometric weight and factor tables or com-
       pute the two-dimensional inverse Fast Fourier Transform of a two-dimen-
       sional double complex array.

SYNOPSIS
       SUBROUTINE ZFFTD2(IOPT, M, N, SCALE, X, LDX, Y, LDY, TRIGS, IFAC, WORK, LWORK, IERR)

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

       SUBROUTINE ZFFTD2_64(IOPT, M, N, SCALE, X, LDX, Y, LDY, TRIGS, IFAC, WORK, LWORK, IERR)

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




   F95 INTERFACE
       SUBROUTINE FFT2(IOPT, M, N, SCALE, X, LDX, Y, LDY, TRIGS,
                 IFAC, WORK, LWORK, IERR)

       INTEGER, INTENT(IN) :: IOPT, M
       INTEGER, INTENT(IN), OPTIONAL :: N, 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 FFT2_64(IOPT, M, N, SCALE, X, LDX, Y, LDY, TRIGS, IFAC, WORK, LWORK, IERR)

       INTEGER(8), INTENT(IN) :: IOPT, M
       INTEGER(8), INTENT(IN), OPTIONAL :: N, 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 zfftd2_ (int *iopt, int *n1, int *n2, double *scale, doublecomplex
                 *x, int *ldx, double *y, int *ldy, double *trigs, int  *ifac,
                 double *work, int *lwork, int *ierr);

       void zfftd2_64_ (long *iopt, long *n1, long *n2, double *scale, double-
                 complex *x, long *ldx, double *y, long *ldy,  double  *trigs,
                 long *ifac, double *work, long *lwork, long *ierr);



PURPOSE
       zfftd2  initializes  the trigonometric weight and factor tables or com-
       putes  the  two-dimensional  inverse  Fast  Fourier  Transform   of   a
       two-dimensional double complex array.  In computing the two-dimensional
       FFT, one-dimensional FFTs are computed along  the  rows  of  the  input
       array.  One-dimensional FFTs are then computed along the columns of the
       intermediate results.

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

       where
       k1 ranges from 0 to M-1 and k2 ranges from 0 to N-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)
       In complex-to-real transform of length M,  the  (M/2+1)  complex  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.


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.

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

       X (input) X  is a double complex array of dimensions (LDX, N) that con-
                 tains input data to be transformed.

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

       Y (output)
                 Y  is  a  double  precision array of dimensions (LDY, N) 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.

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

       TRIGS (input/output)
                 Double precision array of length 2*(M+N)  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 2*128 that contains the
                 factors of M and N.  The factors are computed when  the  rou-
                 tine  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 MAX(M,2*N) 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, 1
                 -2 = M < 0
                 -3 = N < 0
                 -4 = (LDX < M/2+1)
                 -5 = LDY not equal 2*LDX when X and Y are same array
                 -6 = (LDY < 2*LDX or LDY odd) when X and Y are same array
                 -7 = (LWORK not equal 0) and (LWORK < MAX(M,2*N))
                 -8 = memory allocation failed

SEE ALSO
       fft

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



                                  7 Nov 2015                        zfftd2(3P)