Parallel Colt 0.7.2

edu.emory.mathcs.jtransforms.fft
Class DoubleFFT_3D

java.lang.Object
  extended by edu.emory.mathcs.jtransforms.fft.DoubleFFT_3D

public class DoubleFFT_3D
extends Object

Computes 3D Discrete Fourier Transform (DFT) of complex and real, double precision data. The sizes of all three dimensions can be arbitrary numbers. This is a parallel implementation of split-radix and mixed-radix algorithms optimized for SMP systems.

Part of the code is derived from General Purpose FFT Package written by Takuya Ooura (http://www.kurims.kyoto-u.ac.jp/~ooura/fft.html)

Author:
Piotr Wendykier (piotr.wendykier@gmail.com)

Constructor Summary
DoubleFFT_3D(int slices, int rows, int columns)
          Creates new instance of DoubleFFT_3D.
 
Method Summary
 void complexForward(double[] a)
          Computes 3D forward DFT of complex data leaving the result in a.
 void complexForward(double[][][] a)
          Computes 3D forward DFT of complex data leaving the result in a.
 void complexInverse(double[][][] a, boolean scale)
          Computes 3D inverse DFT of complex data leaving the result in a.
 void complexInverse(double[] a, boolean scale)
          Computes 3D inverse DFT of complex data leaving the result in a.
 void realForward(double[] a)
          Computes 3D forward DFT of real data leaving the result in a .
 void realForward(double[][][] a)
          Computes 3D forward DFT of real data leaving the result in a .
 void realForwardFull(double[] a)
          Computes 3D forward DFT of real data leaving the result in a .
 void realForwardFull(double[][][] a)
          Computes 3D forward DFT of real data leaving the result in a .
 void realInverse(double[][][] a, boolean scale)
          Computes 3D inverse DFT of real data leaving the result in a .
 void realInverse(double[] a, boolean scale)
          Computes 3D inverse DFT of real data leaving the result in a .
 void realInverseFull(double[][][] a, boolean scale)
          Computes 3D inverse DFT of real data leaving the result in a .
 void realInverseFull(double[] a, boolean scale)
          Computes 3D inverse DFT of real data leaving the result in a .
 
Methods inherited from class java.lang.Object
equals, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
 

Constructor Detail

DoubleFFT_3D

public DoubleFFT_3D(int slices,
                    int rows,
                    int columns)
Creates new instance of DoubleFFT_3D.

Parameters:
slices - number of slices
rows - number of rows
columns - number of columns
Method Detail

complexForward

public void complexForward(double[] a)
Computes 3D forward DFT of complex data leaving the result in a. The data is stored in 1D array addressed in slice-major, then row-major, then column-major, in order of significance, i.e. element (i,j,k) of 3D array x[slices][rows][2*columns] is stored in a[i*sliceStride + j*rowStride + k], where sliceStride = rows * 2 * columns and rowStride = 2 * columns. Complex number is stored as two double values in sequence: the real and imaginary part, i.e. the input array must be of size slices*rows*2*columns. The physical layout of the input data is as follows:
 a[k1*sliceStride + k2*rowStride + 2*k3] = Re[k1][k2][k3], 
 a[k1*sliceStride + k2*rowStride + 2*k3+1] = Im[k1][k2][k3], 0<=k1<slices, 0<=k2<rows, 0<=k3<columns,
 

Parameters:
a - data to transform

complexForward

public void complexForward(double[][][] a)
Computes 3D forward DFT of complex data leaving the result in a. The data is stored in 3D array. Complex data is represented by 2 double values in sequence: the real and imaginary part, i.e. the input array must be of size slices by rows by 2*columns. The physical layout of the input data is as follows:
 a[k1][k2][2*k3] = Re[k1][k2][k3], 
 a[k1][k2][2*k3+1] = Im[k1][k2][k3], 0<=k1<slices, 0<=k2<rows, 0<=k3<columns,
 

Parameters:
a - data to transform

complexInverse

public void complexInverse(double[] a,
                           boolean scale)
Computes 3D inverse DFT of complex data leaving the result in a. The data is stored in a 1D array addressed in slice-major, then row-major, then column-major, in order of significance, i.e. element (i,j,k) of 3-d array x[slices][rows][2*columns] is stored in a[i*sliceStride + j*rowStride + k], where sliceStride = rows * 2 * columns and rowStride = 2 * columns. Complex number is stored as two double values in sequence: the real and imaginary part, i.e. the input array must be of size slices*rows*2*columns. The physical layout of the input data is as follows:
 a[k1*sliceStride + k2*rowStride + 2*k3] = Re[k1][k2][k3], 
 a[k1*sliceStride + k2*rowStride + 2*k3+1] = Im[k1][k2][k3], 0<=k1<slices, 0<=k2<rows, 0<=k3<columns,
 

Parameters:
a - data to transform
scale - if true then scaling is performed

complexInverse

public void complexInverse(double[][][] a,
                           boolean scale)
Computes 3D inverse DFT of complex data leaving the result in a. The data is stored in a 3D array. Complex data is represented by 2 double values in sequence: the real and imaginary part, i.e. the input array must be of size slices by rows by 2*columns. The physical layout of the input data is as follows:
 a[k1][k2][2*k3] = Re[k1][k2][k3], 
 a[k1][k2][2*k3+1] = Im[k1][k2][k3], 0<=k1<slices, 0<=k2<rows, 0<=k3<columns,
 

Parameters:
a - data to transform
scale - if true then scaling is performed

realForward

public void realForward(double[] a)
Computes 3D forward DFT of real data leaving the result in a . This method only works when the sizes of all three dimensions are power-of-two numbers. The data is stored in a 1D array addressed in slice-major, then row-major, then column-major, in order of significance, i.e. element (i,j,k) of 3-d array x[slices][rows][2*columns] is stored in a[i*sliceStride + j*rowStride + k], where sliceStride = rows * 2 * columns and rowStride = 2 * columns. The physical layout of the output data is as follows:
 a[k1*sliceStride + k2*rowStride + 2*k3] = Re[k1][k2][k3]
                 = Re[(slices-k1)%slices][(rows-k2)%rows][columns-k3], 
 a[k1*sliceStride + k2*rowStride + 2*k3+1] = Im[k1][k2][k3]
                   = -Im[(slices-k1)%slices][(rows-k2)%rows][columns-k3], 
     0<=k1<slices, 0<=k2<rows, 0<k3<columns/2, 
 a[k1*sliceStride + k2*rowStride] = Re[k1][k2][0]
              = Re[(slices-k1)%slices][rows-k2][0], 
 a[k1*sliceStride + k2*rowStride + 1] = Im[k1][k2][0]
              = -Im[(slices-k1)%slices][rows-k2][0], 
 a[k1*sliceStride + (rows-k2)*rowStride + 1] = Re[(slices-k1)%slices][k2][columns/2]
                 = Re[k1][rows-k2][columns/2], 
 a[k1*sliceStride + (rows-k2)*rowStride] = -Im[(slices-k1)%slices][k2][columns/2]
                 = Im[k1][rows-k2][columns/2], 
     0<=k1<slices, 0<k2<rows/2, 
 a[k1*sliceStride] = Re[k1][0][0]
             = Re[slices-k1][0][0], 
 a[k1*sliceStride + 1] = Im[k1][0][0]
             = -Im[slices-k1][0][0], 
 a[k1*sliceStride + (rows/2)*rowStride] = Re[k1][rows/2][0]
                = Re[slices-k1][rows/2][0], 
 a[k1*sliceStride + (rows/2)*rowStride + 1] = Im[k1][rows/2][0]
                = -Im[slices-k1][rows/2][0], 
 a[(slices-k1)*sliceStride + 1] = Re[k1][0][columns/2]
                = Re[slices-k1][0][columns/2], 
 a[(slices-k1)*sliceStride] = -Im[k1][0][columns/2]
                = Im[slices-k1][0][columns/2], 
 a[(slices-k1)*sliceStride + (rows/2)*rowStride + 1] = Re[k1][rows/2][columns/2]
                   = Re[slices-k1][rows/2][columns/2], 
 a[(slices-k1)*sliceStride + (rows/2) * rowStride] = -Im[k1][rows/2][columns/2]
                   = Im[slices-k1][rows/2][columns/2], 
     0<k1<slices/2, 
 a[0] = Re[0][0][0], 
 a[1] = Re[0][0][columns/2], 
 a[(rows/2)*rowStride] = Re[0][rows/2][0], 
 a[(rows/2)*rowStride + 1] = Re[0][rows/2][columns/2], 
 a[(slices/2)*sliceStride] = Re[slices/2][0][0], 
 a[(slices/2)*sliceStride + 1] = Re[slices/2][0][columns/2], 
 a[(slices/2)*sliceStride + (rows/2)*rowStride] = Re[slices/2][rows/2][0], 
 a[(slices/2)*sliceStride + (rows/2)*rowStride + 1] = Re[slices/2][rows/2][columns/2]
 
This method computes only half of the elements of the real transform. The other half satisfies the symmetry condition. If you want the full real forward transform, use realForwardFull. To get back the original data, use realInverse on the output of this method.

Parameters:
a - data to transform

realForward

public void realForward(double[][][] a)
Computes 3D forward DFT of real data leaving the result in a . This method only works when the sizes of all three dimensions are power-of-two numbers. The data is stored in a 3D array. The physical layout of the output data is as follows:
 a[k1][k2][2*k3] = Re[k1][k2][k3]
                 = Re[(slices-k1)%slices][(rows-k2)%rows][columns-k3], 
 a[k1][k2][2*k3+1] = Im[k1][k2][k3]
                   = -Im[(slices-k1)%slices][(rows-k2)%rows][columns-k3], 
     0<=k1<slices, 0<=k2<rows, 0<k3<columns/2, 
 a[k1][k2][0] = Re[k1][k2][0]
              = Re[(slices-k1)%slices][rows-k2][0], 
 a[k1][k2][1] = Im[k1][k2][0]
              = -Im[(slices-k1)%slices][rows-k2][0], 
 a[k1][rows-k2][1] = Re[(slices-k1)%slices][k2][columns/2]
                 = Re[k1][rows-k2][columns/2], 
 a[k1][rows-k2][0] = -Im[(slices-k1)%slices][k2][columns/2]
                 = Im[k1][rows-k2][columns/2], 
     0<=k1<slices, 0<k2<rows/2, 
 a[k1][0][0] = Re[k1][0][0]
             = Re[slices-k1][0][0], 
 a[k1][0][1] = Im[k1][0][0]
             = -Im[slices-k1][0][0], 
 a[k1][rows/2][0] = Re[k1][rows/2][0]
                = Re[slices-k1][rows/2][0], 
 a[k1][rows/2][1] = Im[k1][rows/2][0]
                = -Im[slices-k1][rows/2][0], 
 a[slices-k1][0][1] = Re[k1][0][columns/2]
                = Re[slices-k1][0][columns/2], 
 a[slices-k1][0][0] = -Im[k1][0][columns/2]
                = Im[slices-k1][0][columns/2], 
 a[slices-k1][rows/2][1] = Re[k1][rows/2][columns/2]
                   = Re[slices-k1][rows/2][columns/2], 
 a[slices-k1][rows/2][0] = -Im[k1][rows/2][columns/2]
                   = Im[slices-k1][rows/2][columns/2], 
     0<k1<slices/2, 
 a[0][0][0] = Re[0][0][0], 
 a[0][0][1] = Re[0][0][columns/2], 
 a[0][rows/2][0] = Re[0][rows/2][0], 
 a[0][rows/2][1] = Re[0][rows/2][columns/2], 
 a[slices/2][0][0] = Re[slices/2][0][0], 
 a[slices/2][0][1] = Re[slices/2][0][columns/2], 
 a[slices/2][rows/2][0] = Re[slices/2][rows/2][0], 
 a[slices/2][rows/2][1] = Re[slices/2][rows/2][columns/2]
 
This method computes only half of the elements of the real transform. The other half satisfies the symmetry condition. If you want the full real forward transform, use realForwardFull. To get back the original data, use realInverse on the output of this method.

Parameters:
a - data to transform

realForwardFull

public void realForwardFull(double[] a)
Computes 3D forward DFT of real data leaving the result in a . This method computes full real forward transform, i.e. you will get the same result as from complexForward called with all imaginary part equal 0. Because the result is stored in a, the input array must be of size slices*rows*2*columns, with only the first slices*rows*columns elements filled with real data. To get back the original data, use complexInverse on the output of this method.

Parameters:
a - data to transform

realForwardFull

public void realForwardFull(double[][][] a)
Computes 3D forward DFT of real data leaving the result in a . This method computes full real forward transform, i.e. you will get the same result as from complexForward called with all imaginary part equal 0. Because the result is stored in a, the input array must be of size slices by rows by 2*columns, with only the first slices by rows by columns elements filled with real data. To get back the original data, use complexInverse on the output of this method.

Parameters:
a - data to transform

realInverse

public void realInverse(double[] a,
                        boolean scale)
Computes 3D inverse DFT of real data leaving the result in a . This method only works when the sizes of all three dimensions are power-of-two numbers. The data is stored in a 1D array addressed in slice-major, then row-major, then column-major, in order of significance, i.e. element (i,j,k) of 3-d array x[slices][rows][2*columns] is stored in a[i*sliceStride + j*rowStride + k], where sliceStride = rows * 2 * columns and rowStride = 2 * columns. The physical layout of the input data has to be as follows:
 a[k1*sliceStride + k2*rowStride + 2*k3] = Re[k1][k2][k3]
                 = Re[(slices-k1)%slices][(rows-k2)%rows][columns-k3], 
 a[k1*sliceStride + k2*rowStride + 2*k3+1] = Im[k1][k2][k3]
                   = -Im[(slices-k1)%slices][(rows-k2)%rows][columns-k3], 
     0<=k1<slices, 0<=k2<rows, 0<k3<columns/2, 
 a[k1*sliceStride + k2*rowStride] = Re[k1][k2][0]
              = Re[(slices-k1)%slices][rows-k2][0], 
 a[k1*sliceStride + k2*rowStride + 1] = Im[k1][k2][0]
              = -Im[(slices-k1)%slices][rows-k2][0], 
 a[k1*sliceStride + (rows-k2)*rowStride + 1] = Re[(slices-k1)%slices][k2][columns/2]
                 = Re[k1][rows-k2][columns/2], 
 a[k1*sliceStride + (rows-k2)*rowStride] = -Im[(slices-k1)%slices][k2][columns/2]
                 = Im[k1][rows-k2][columns/2], 
     0<=k1<slices, 0<k2<rows/2, 
 a[k1*sliceStride] = Re[k1][0][0]
             = Re[slices-k1][0][0], 
 a[k1*sliceStride + 1] = Im[k1][0][0]
             = -Im[slices-k1][0][0], 
 a[k1*sliceStride + (rows/2)*rowStride] = Re[k1][rows/2][0]
                = Re[slices-k1][rows/2][0], 
 a[k1*sliceStride + (rows/2)*rowStride + 1] = Im[k1][rows/2][0]
                = -Im[slices-k1][rows/2][0], 
 a[(slices-k1)*sliceStride + 1] = Re[k1][0][columns/2]
                = Re[slices-k1][0][columns/2], 
 a[(slices-k1)*sliceStride] = -Im[k1][0][columns/2]
                = Im[slices-k1][0][columns/2], 
 a[(slices-k1)*sliceStride + (rows/2)*rowStride + 1] = Re[k1][rows/2][columns/2]
                   = Re[slices-k1][rows/2][columns/2], 
 a[(slices-k1)*sliceStride + (rows/2) * rowStride] = -Im[k1][rows/2][columns/2]
                   = Im[slices-k1][rows/2][columns/2], 
     0<k1<slices/2, 
 a[0] = Re[0][0][0], 
 a[1] = Re[0][0][columns/2], 
 a[(rows/2)*rowStride] = Re[0][rows/2][0], 
 a[(rows/2)*rowStride + 1] = Re[0][rows/2][columns/2], 
 a[(slices/2)*sliceStride] = Re[slices/2][0][0], 
 a[(slices/2)*sliceStride + 1] = Re[slices/2][0][columns/2], 
 a[(slices/2)*sliceStride + (rows/2)*rowStride] = Re[slices/2][rows/2][0], 
 a[(slices/2)*sliceStride + (rows/2)*rowStride + 1] = Re[slices/2][rows/2][columns/2]
 
This method computes only half of the elements of the real transform. The other half satisfies the symmetry condition. If you want the full real inverse transform, use realInverseFull.

Parameters:
a - data to transform
scale - if true then scaling is performed

realInverse

public void realInverse(double[][][] a,
                        boolean scale)
Computes 3D inverse DFT of real data leaving the result in a . This method only works when the sizes of all three dimensions are power-of-two numbers. The data is stored in a 3D array. The physical layout of the input data has to be as follows:
 a[k1][k2][2*k3] = Re[k1][k2][k3]
                 = Re[(slices-k1)%slices][(rows-k2)%rows][columns-k3], 
 a[k1][k2][2*k3+1] = Im[k1][k2][k3]
                   = -Im[(slices-k1)%slices][(rows-k2)%rows][columns-k3], 
     0<=k1<slices, 0<=k2<rows, 0<k3<columns/2, 
 a[k1][k2][0] = Re[k1][k2][0]
              = Re[(slices-k1)%slices][rows-k2][0], 
 a[k1][k2][1] = Im[k1][k2][0]
              = -Im[(slices-k1)%slices][rows-k2][0], 
 a[k1][rows-k2][1] = Re[(slices-k1)%slices][k2][columns/2]
                 = Re[k1][rows-k2][columns/2], 
 a[k1][rows-k2][0] = -Im[(slices-k1)%slices][k2][columns/2]
                 = Im[k1][rows-k2][columns/2], 
     0<=k1<slices, 0<k2<rows/2, 
 a[k1][0][0] = Re[k1][0][0]
             = Re[slices-k1][0][0], 
 a[k1][0][1] = Im[k1][0][0]
             = -Im[slices-k1][0][0], 
 a[k1][rows/2][0] = Re[k1][rows/2][0]
                = Re[slices-k1][rows/2][0], 
 a[k1][rows/2][1] = Im[k1][rows/2][0]
                = -Im[slices-k1][rows/2][0], 
 a[slices-k1][0][1] = Re[k1][0][columns/2]
                = Re[slices-k1][0][columns/2], 
 a[slices-k1][0][0] = -Im[k1][0][columns/2]
                = Im[slices-k1][0][columns/2], 
 a[slices-k1][rows/2][1] = Re[k1][rows/2][columns/2]
                   = Re[slices-k1][rows/2][columns/2], 
 a[slices-k1][rows/2][0] = -Im[k1][rows/2][columns/2]
                   = Im[slices-k1][rows/2][columns/2], 
     0<k1<slices/2, 
 a[0][0][0] = Re[0][0][0], 
 a[0][0][1] = Re[0][0][columns/2], 
 a[0][rows/2][0] = Re[0][rows/2][0], 
 a[0][rows/2][1] = Re[0][rows/2][columns/2], 
 a[slices/2][0][0] = Re[slices/2][0][0], 
 a[slices/2][0][1] = Re[slices/2][0][columns/2], 
 a[slices/2][rows/2][0] = Re[slices/2][rows/2][0], 
 a[slices/2][rows/2][1] = Re[slices/2][rows/2][columns/2]
 
This method computes only half of the elements of the real transform. The other half satisfies the symmetry condition. If you want the full real inverse transform, use realInverseFull.

Parameters:
a - data to transform
scale - if true then scaling is performed

realInverseFull

public void realInverseFull(double[] a,
                            boolean scale)
Computes 3D inverse DFT of real data leaving the result in a . This method computes full real inverse transform, i.e. you will get the same result as from complexInverse called with all imaginary part equal 0. Because the result is stored in a, the input array must be of size slices*rows*2*columns, with only the first slices*rows*columns elements filled with real data.

Parameters:
a - data to transform
scale - if true then scaling is performed

realInverseFull

public void realInverseFull(double[][][] a,
                            boolean scale)
Computes 3D inverse DFT of real data leaving the result in a . This method computes full real inverse transform, i.e. you will get the same result as from complexInverse called with all imaginary part equal 0. Because the result is stored in a, the input array must be of size slices by rows by 2*columns, with only the first slices by rows by columns elements filled with real data.

Parameters:
a - data to transform
scale - if true then scaling is performed

Parallel Colt 0.7.2

Jump to the Parallel Colt Homepage