MFFM FFTw Wrapper
real2DFFTExample.cc

This is an example of how to use the class.

/* Copyright 2001,2002 Matt Flax <flatmax@ieee.org>
This file is part of the MFFM FFTw Wrapper library.
MFFM MFFM FFTw Wrapper library is free software; you can
redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation; either version 2 of the License, or
(at your option) any later version.
MFFM FFTw Wrapper library is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You have received a copy of the GNU General Public License
along with the MFFM FFTw Wrapper library
*/
#include <string.h>
#include "real2DFFT.H"
#include <iostream>
using namespace std;
int main(void){
int x=8, y=8;
real2DFFTData *fftData = new real2DFFTData(x,y);
real2DFFT *fft= new real2DFFT(fftData);
// clear the data
fftData->clearInput();
fftData->clearOutput();
int temp=x/2, temp2=y/2;
for (int j=0;j<x;j++)
fftData->in[temp2+j*x]=10000.0;
for (int j=0;j<y;j++)
fftData->in[temp*y+j]=10000.0;
// fftData->in[temp*y+(y-1)/2]=20000.0;
for (int i=0;i<fftData->getXSize();i++){
for (int j=0;j<fftData->getYSize();j++)
cout<<fftData->in[i*x+j]<<'\t';
cout<<endl;
}
cout<<'\n'<<endl;
fft->fwdTransform();
fftData->reScale();
fftData->compPowerSpec();
fft->invTransform();
for (int i=0;i<fftData->getXSize();i++){
for (int j=0;j<fftData->getYSize();j++)
cout<<fftData->in[i*x+j]<<'\t';
cout<<endl;
}
cout<<'\n'<<endl;
/* for (int i=0;i<fftData.getXSize();i++){
for (int j=0;j<fftData.getYHalfSize();j++)
cout<<fftData.out[i][j].im<<'\t';
cout<<endl;
}*/
for (int i=0;i<x;i++){
for (int j=0;j<y/2+1;j++)
cout<<fftData->power[i*(y/2+1)+j]<<'\t';
cout<<endl;
}
delete fftData;
delete fft;
}
class real2DFFTData controls and manipulates real 2D fft data
Definition real2DFFT.H:36
int getXSize()
The row count.
Definition real2DFFT.H:68
int getYSize()
The column count.
Definition real2DFFT.H:70
void reScale(void)
Scales the output down by the number of elements.
void clearOutput(void)
Zeros the out awway.
Definition real2DFFT.H:98
void compPowerSpec()
This function computes the power spectrum and updates the totalPower, maxPower and minPower.
void clearInput(void)
Zeros the in array.
Definition real2DFFT.H:96
fftw_real * in
The input data and power spectrum.
Definition real2DFFT.H:45
class real2DFFT controls fftw plans and executes fwd/inv transforms
Definition real2DFFT.H:102
void fwdTransform()
Forward transform the data (in to out)
void invTransform()
Inverse transform the data (out to in)