Go to main content
Oracle Developer Studio 12.5 Man Pages

Exit Print View

Updated: June 2017
 
 

dfftz2 (3p)

Name

dfftz2 - pute the two-dimensional forward Fast Fourier Transform of a two-dimen- sional double precision array.

Synopsis

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

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

SUBROUTINE DFFTZ2_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 Y(LDY, *)
DOUBLE PRECISION X(LDX, *), SCALE, TRIGS(*), WORK(*)




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

INTEGER, INTENT(IN) :: IOPT
INTEGER, INTENT(IN), OPTIONAL :: M, N, LDX, LDY, LWORK
REAL(8), INTENT(IN), OPTIONAL :: SCALE
REAL(8), INTENT(IN), DIMENSION(:,:) :: X
COMPLEX(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
INTEGER(8), INTENT(IN), OPTIONAL :: M, N, LDX, LDY, LWORK
REAL(8), INTENT(IN), OPTIONAL :: SCALE
REAL(8), INTENT(IN), DIMENSION(:,:) :: X
COMPLEX(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  dfftz2_  (int  *iopt, int *n1, int *n2, double *scale, double *x,
int *ldx, doublecomplex *y,  int  *ldy,  double  *trigs,  int
*ifac, double *work, int *lwork, int *ierr);

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

Description

Oracle Solaris Studio Performance Library                           dfftz2(3P)



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

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

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

       SUBROUTINE DFFTZ2_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 Y(LDY, *)
       DOUBLE PRECISION X(LDX, *), SCALE, TRIGS(*), WORK(*)




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

       INTEGER, INTENT(IN) :: IOPT
       INTEGER, INTENT(IN), OPTIONAL :: M, N, LDX, LDY, LWORK
       REAL(8), INTENT(IN), OPTIONAL :: SCALE
       REAL(8), INTENT(IN), DIMENSION(:,:) :: X
       COMPLEX(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
       INTEGER(8), INTENT(IN), OPTIONAL :: M, N, LDX, LDY, LWORK
       REAL(8), INTENT(IN), OPTIONAL :: SCALE
       REAL(8), INTENT(IN), DIMENSION(:,:) :: X
       COMPLEX(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  dfftz2_  (int  *iopt, int *n1, int *n2, double *scale, double *x,
                 int *ldx, doublecomplex *y,  int  *ldy,  double  *trigs,  int
                 *ifac, double *work, int *lwork, int *ierr);

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



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

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

       where
       k1 ranges from 0 to M-1 and k2 ranges from 0 to N-1
       i = sqrt(-1)
       isign = -1 for forward transform
       W1 = exp(isign*i*j1*k1*2*pi/M)
       W2 = exp(isign*i*j2*k2*2*pi/N)
       In real-to-complex transform of length M, the  (M/2+1)  complex  output
       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 forward 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  precision array of dimensions (LDX, N) that
                 contains input data to be transformed.  X and Y  can  be  the
                 same array.

       LDX (input)
                 Leading  dimension of X.  LDX >= M if X is not the same array
                 as Y.  Else, LDX = 2*LDY.  Unchanged on exit.

       Y (output)
                 Y is a double complex array of dimensions (LDY, N) that  con-
                 tains  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.  Otherwise,
                 it is assumed that there is no overlap between  X  and  Y  in
                 memory.

       LDY (input)
                 Leading dimension of Y.  LDY >= M/2+1 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 rou-
                 tine.   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  = (LDX < M) or (LDX not equal 2*LDY when X and Y are same
                 array)
                 -5 = (LDY < M/2+1)
                 -6 = (LWORK not equal 0) and (LWORK < MAX(M,2*N))
                 -7 = memory allocation failed

SEE ALSO
       fft

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



                                  7 Nov 2015                        dfftz2(3P)