TBTK
Need a break? Support the development by playing Polarity Puzzles
TBTK::FourierTransform Class Reference

Fourier transform. More...

#include <FourierTransform.h>

Classes

class  ForwardPlan
 
class  InversePlan
 
class  Plan
 

Static Public Member Functions

static void transform (const CArray< std::complex< double >> &in, CArray< std::complex< double >> &out, const std::vector< unsigned int > &ranges, int sign)
 
template<typename DataType >
static void transform (Plan< DataType > &plan)
 
static void forward (const CArray< std::complex< double >> &in, CArray< std::complex< double >> &out, const std::vector< unsigned int > &ranges)
 
static void inverse (const CArray< std::complex< double >> &in, CArray< std::complex< double >> &out, const std::vector< unsigned int > &ranges)
 

Detailed Description

Fourier transform.

Example

#include "TBTK/TBTK.h"
#include "TBTK/Visualization/MatPlotLib/Plotter.h"
#include <complex>
#include <string>
using namespace std;
using namespace TBTK;
using namespace Visualization::MatPlotLib;
void example1D(){
const unsigned int SIZE = 100;
CArray<complex<double>> input(SIZE);
for(unsigned int n = 0; n < SIZE; n++){
if(n < SIZE/10)
input[n] = 1;
else
input[n] = 0;
}
CArray<complex<double>> output(100);
FourierTransform::forward(input, output, {SIZE});
Array<double> outputReal({SIZE});
Array<double> outputImaginary({SIZE});
for(unsigned int n = 0; n < SIZE; n++){
outputReal[{n}] = real(output[n]);
outputImaginary[{n}] = imag(output[n]);
}
Plotter plotter;
plotter.setTitle("1D");
plotter.plot(outputReal, {{"linewidth", "2"}, {"label", "Real"}});
plotter.plot(
outputImaginary,
{
{"linewidth", "2"},
{"label", "Imaginary"}
}
);
plotter.save("figures/1D.png");
}
void example2D(){
const unsigned int SIZE_X = 100;
const unsigned int SIZE_Y = 100;
Array<complex<double>> input({SIZE_X, SIZE_Y});
for(unsigned int x = 0; x < SIZE_X; x++){
for(unsigned int y = 0; y < SIZE_Y; y++){
if(x < SIZE_X/20 && y < SIZE_Y/20)
input[{x, y}] = 1;
else
input[{x, y}] = 0;
}
}
Array<complex<double>> output({SIZE_X, SIZE_Y});
input.getData(),
output.getData(),
{SIZE_X, SIZE_Y}
);
Array<double> outputReal({SIZE_X, SIZE_Y});
Array<double> outputImaginary({SIZE_X, SIZE_Y});
for(unsigned int x = 0; x < SIZE_X; x++){
for(unsigned int y = 0; y < SIZE_Y; y++){
outputReal[{x, y}] = real(output[{x, y}]);
outputImaginary[{x, y}] = imag(output[{x, y}]);
}
}
Plotter plotter;
plotter.setTitle("2D real");
plotter.plot(outputReal);
plotter.save("figures/2DReal.png");
plotter.clear();
plotter.setTitle("2D imaginary");
plotter.plot(outputImaginary);
plotter.save("figures/2DImaginary.png");
}
void examplePlan(){
const unsigned int SIZE = 100;
CArray<complex<double>> input(SIZE);
CArray<complex<double>> output(SIZE);
input,
output,
{SIZE}
);
Plotter plotter;
plotter.setTitle("With plan");
for(unsigned int n = 0; n < 10; n++){
for(unsigned int x = 0; x < SIZE; x++){
if(x < SIZE/(n+10))
input[x] = 1;
else
input[x] = 0;
}
Array<double> outputReal({SIZE});
Array<double> outputImaginary({SIZE});
for(unsigned int n = 0; n < SIZE; n++){
outputReal[{n}] = real(output[n]);
outputImaginary[{n}] = imag(output[n]);
}
if(n == 0){
plotter.plot(
outputReal,
{
{"linewidth", "2"},
{"color", "black"},
{"linestyle", "-"},
{"label", "Real"}
}
);
plotter.plot(
outputImaginary,
{
{"linewidth", "2"},
{"color", "orangered"},
{"linestyle", "-"},
{"label", "Imaginary"}
}
);
}
else{
plotter.plot(
outputReal,
{
{"linewidth", "2"},
{"color", "black"},
{"linestyle", "-"}
}
);
plotter.plot(
outputImaginary,
{
{"linewidth", "2"},
{"color", "orangered"},
{"linestyle", "-"}
}
);
}
}
plotter.save("figures/WithPlan.png");
}
int main(){
example1D();
example2D();
examplePlan();
}

Output

FourierTransformFourierTransform1D.png
FourierTransformFourierTransform2DReal.png
FourierTransformFourierTransform2DImaginary.png
FourierTransformFourierTransformWithPlan.png

Member Function Documentation

◆ forward()

void TBTK::FourierTransform::forward ( const CArray< std::complex< double >> &  in,
CArray< std::complex< double >> &  out,
const std::vector< unsigned int > &  ranges 
)
inlinestatic

N-dimensional complex forward Fourier transform.

Parameters
inInput data.
outOutput data.
rangesThe dimensions of the data.

◆ inverse()

void TBTK::FourierTransform::inverse ( const CArray< std::complex< double >> &  in,
CArray< std::complex< double >> &  out,
const std::vector< unsigned int > &  ranges 
)
inlinestatic

N-dimensional complex inverse Fourier transform.

Parameters
inInput data.
outOutput data.
rangesThe dimensions of the data.

◆ transform() [1/2]

static void TBTK::FourierTransform::transform ( const CArray< std::complex< double >> &  in,
CArray< std::complex< double >> &  out,
const std::vector< unsigned int > &  ranges,
int  sign 
)
static

N-dimensional complex Fourier transform.

Parameters
inInput data.
outOutput data.
rangesThe dimensions of the data.
signThe sign to use in the exponent of the Fourier transform.

◆ transform() [2/2]

template<typename DataType >
void TBTK::FourierTransform::transform ( Plan< DataType > &  plan)
inlinestatic

Execute a planned transform.

Parameters
planThe plan to execute.

The documentation for this class was generated from the following file: