public class DoubleFFT_2D extends Object
| Constructor and Description |
|---|
DoubleFFT_2D(int rows,
int columns)
Creates new instance of DoubleFFT_2D.
|
| Modifier and Type | Method and Description |
|---|---|
void |
complexForward(double[] a)
Computes 2D forward DFT of complex data leaving the result in
a. |
void |
complexForward(double[][] a)
Computes 2D forward DFT of complex data leaving the result in
a. |
void |
complexInverse(double[][] a,
boolean scale)
Computes 2D inverse DFT of complex data leaving the result in
a. |
void |
complexInverse(double[] a,
boolean scale)
Computes 2D inverse DFT of complex data leaving the result in
a. |
void |
realForward(double[] a)
Computes 2D forward DFT of real data leaving the result in
a . |
void |
realForward(double[][] a)
Computes 2D forward DFT of real data leaving the result in
a . |
void |
realForwardFull(double[] a)
Computes 2D forward DFT of real data leaving the result in
a . |
void |
realForwardFull(double[][] a)
Computes 2D forward DFT of real data leaving the result in
a . |
void |
realInverse(double[][] a,
boolean scale)
Computes 2D inverse DFT of real data leaving the result in
a . |
void |
realInverse(double[] a,
boolean scale)
Computes 2D inverse DFT of real data leaving the result in
a . |
void |
realInverseFull(double[][] a,
boolean scale)
Computes 2D inverse DFT of real data leaving the result in
a . |
void |
realInverseFull(double[] a,
boolean scale)
Computes 2D inverse DFT of real data leaving the result in
a . |
public DoubleFFT_2D(int rows,
int columns)
rows - number of rowscolumns - number of columnspublic void complexForward(double[] a)
a. The data is stored in 1D array in
row-major order. Complex number is stored as two double values in sequence: the real and imaginary part, i.e. the
input array must be of size rows*2*columns. The physical layout of the input data has to be as follows:a[k1*2*columns+2*k2] = Re[k1][k2], a[k1*2*columns+2*k2+1] = Im[k1][k2], 0<=k1<rows, 0<=k2<columns,
a - data to transformpublic void complexForward(double[][] a)
a. The data is stored in 2D 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 rows by 2*columns. The physical layout of the input data has to be as follows:a[k1][2*k2] = Re[k1][k2], a[k1][2*k2+1] = Im[k1][k2], 0<=k1<rows, 0<=k2<columns,
a - data to transformpublic void complexInverse(double[] a,
boolean scale)
a. The data is stored in 1D array in
row-major order. Complex number is stored as two double values in sequence: the real and imaginary part, i.e. the
input array must be of size rows*2*columns. The physical layout of the input data has to be as follows:a[k1*2*columns+2*k2] = Re[k1][k2], a[k1*2*columns+2*k2+1] = Im[k1][k2], 0<=k1<rows, 0<=k2<columns,
a - data to transformscale - if true then scaling is performedpublic void complexInverse(double[][] a,
boolean scale)
a. The data is stored in 2D 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 rows by 2*columns. The physical layout of the input data has to be as follows:a[k1][2*k2] = Re[k1][k2], a[k1][2*k2+1] = Im[k1][k2], 0<=k1<rows, 0<=k2<columns,
a - data to transformscale - if true then scaling is performedpublic void realForward(double[] a)
a . This method only works when the sizes
of both dimensions are power-of-two numbers. The physical layout of the output data is as follows:
a[k1*columns+2*k2] = Re[k1][k2] = Re[rows-k1][columns-k2],
a[k1*columns+2*k2+1] = Im[k1][k2] = -Im[rows-k1][columns-k2],
0<k1<rows, 0<k2<columns/2,
a[2*k2] = Re[0][k2] = Re[0][columns-k2],
a[2*k2+1] = Im[0][k2] = -Im[0][columns-k2],
0<k2<columns/2,
a[k1*columns] = Re[k1][0] = Re[rows-k1][0],
a[k1*columns+1] = Im[k1][0] = -Im[rows-k1][0],
a[(rows-k1)*columns+1] = Re[k1][columns/2] = Re[rows-k1][columns/2],
a[(rows-k1)*columns] = -Im[k1][columns/2] = Im[rows-k1][columns/2],
0<k1<rows/2,
a[0] = Re[0][0],
a[1] = Re[0][columns/2],
a[(rows/2)*columns] = Re[rows/2][0],
a[(rows/2)*columns+1] = Re[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.a - data to transformpublic void realForward(double[][] a)
a . This method only works when the sizes
of both dimensions are power-of-two numbers. The physical layout of the output data is as follows:
a[k1][2*k2] = Re[k1][k2] = Re[rows-k1][columns-k2],
a[k1][2*k2+1] = Im[k1][k2] = -Im[rows-k1][columns-k2],
0<k1<rows, 0<k2<columns/2,
a[0][2*k2] = Re[0][k2] = Re[0][columns-k2],
a[0][2*k2+1] = Im[0][k2] = -Im[0][columns-k2],
0<k2<columns/2,
a[k1][0] = Re[k1][0] = Re[rows-k1][0],
a[k1][1] = Im[k1][0] = -Im[rows-k1][0],
a[rows-k1][1] = Re[k1][columns/2] = Re[rows-k1][columns/2],
a[rows-k1][0] = -Im[k1][columns/2] = Im[rows-k1][columns/2],
0<k1<rows/2,
a[0][0] = Re[0][0],
a[0][1] = Re[0][columns/2],
a[rows/2][0] = Re[rows/2][0],
a[rows/2][1] = Re[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.a - data to transformpublic void realForwardFull(double[] a)
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
rows*2*columns, with only the first rows*columns elements filled with real data. To get back the original data,
use complexInverse on the output of this method.a - data to transformpublic void realForwardFull(double[][] a)
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 rows by
2*columns, with only the first rows by columns elements filled with real data. To get back the original data, use
complexInverse on the output of this method.a - data to transformpublic void realInverse(double[] a,
boolean scale)
a . This method only works when the sizes
of both dimensions are power-of-two numbers. The physical layout of the input data has to be as follows:
a[k1*columns+2*k2] = Re[k1][k2] = Re[rows-k1][columns-k2],
a[k1*columns+2*k2+1] = Im[k1][k2] = -Im[rows-k1][columns-k2],
0<k1<rows, 0<k2<columns/2,
a[2*k2] = Re[0][k2] = Re[0][columns-k2],
a[2*k2+1] = Im[0][k2] = -Im[0][columns-k2],
0<k2<columns/2,
a[k1*columns] = Re[k1][0] = Re[rows-k1][0],
a[k1*columns+1] = Im[k1][0] = -Im[rows-k1][0],
a[(rows-k1)*columns+1] = Re[k1][columns/2] = Re[rows-k1][columns/2],
a[(rows-k1)*columns] = -Im[k1][columns/2] = Im[rows-k1][columns/2],
0<k1<rows/2,
a[0] = Re[0][0],
a[1] = Re[0][columns/2],
a[(rows/2)*columns] = Re[rows/2][0],
a[(rows/2)*columns+1] = Re[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.a - data to transformscale - if true then scaling is performedpublic void realInverse(double[][] a,
boolean scale)
a . This method only works when the sizes
of both dimensions are power-of-two numbers. The physical layout of the input data has to be as follows:
a[k1][2*k2] = Re[k1][k2] = Re[rows-k1][columns-k2],
a[k1][2*k2+1] = Im[k1][k2] = -Im[rows-k1][columns-k2],
0<k1<rows, 0<k2<columns/2,
a[0][2*k2] = Re[0][k2] = Re[0][columns-k2],
a[0][2*k2+1] = Im[0][k2] = -Im[0][columns-k2],
0<k2<columns/2,
a[k1][0] = Re[k1][0] = Re[rows-k1][0],
a[k1][1] = Im[k1][0] = -Im[rows-k1][0],
a[rows-k1][1] = Re[k1][columns/2] = Re[rows-k1][columns/2],
a[rows-k1][0] = -Im[k1][columns/2] = Im[rows-k1][columns/2],
0<k1<rows/2,
a[0][0] = Re[0][0],
a[0][1] = Re[0][columns/2],
a[rows/2][0] = Re[rows/2][0],
a[rows/2][1] = Re[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.a - data to transformscale - if true then scaling is performedpublic void realInverseFull(double[] a,
boolean scale)
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
rows*2*columns, with only the first rows*columns elements filled with real data.a - data to transformscale - if true then scaling is performedpublic void realInverseFull(double[][] a,
boolean scale)
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 rows by
2*columns, with only the first rows by columns elements filled with real data.a - data to transformscale - if true then scaling is performedCopyright © 2019 GSI Helmholtzzentrum für Schwerionenforschung GmbH. All rights reserved.