BornAgain
1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
|
Classes | |
class | Workspace |
Public Types | |
typedef std::vector< double > | double1d_t |
typedef std::vector< double1d_t > | double2d_t |
Public Member Functions | |
FourierTransform () | |
void | fft (const double1d_t &source, double1d_t &result) |
void | fftshift (double1d_t &result) |
void | fft (const double2d_t &source, double2d_t &result) |
void | fftshift (double2d_t &result) |
void | init (int h_src, int w_src) |
Private Member Functions | |
void | fftw_forward_FT (const double2d_t &source) |
Private Attributes | |
Workspace | ws |
Fourier transform of vectors (in 1D or 2D) using Fast Fourier Transform (fftw package).
Usage: std::vector<double> signal, result; FourierTransform ft; ft.fft(signal, result)
Given code rely on code from Jeremy Fix page, http://jeremy.fix.free.fr/spip.php?article15, see also "Efficient convolution using the Fast Fourier Transform, Application in C++" by Jeremy Fix, May 30, 2011
Definition at line 34 of file FourierTransform.h.
typedef std::vector<double> FourierTransform::double1d_t |
definition of 1D vector of double
Definition at line 38 of file FourierTransform.h.
typedef std::vector<double1d_t> FourierTransform::double2d_t |
definition of 2D vector of double
Definition at line 41 of file FourierTransform.h.
|
default |
void FourierTransform::fft | ( | const double1d_t & | source, |
double1d_t & | result | ||
) |
FT in 1D.
Definition at line 124 of file FourierTransform.cpp.
Referenced by IntensityDataFunctions::FT2DArray().
void FourierTransform::fftshift | ( | double1d_t & | result | ) |
Shift low frequency to the center of 1D FT array.
Definition at line 142 of file FourierTransform.cpp.
Referenced by IntensityDataFunctions::FT2DArray().
void FourierTransform::fft | ( | const double2d_t & | source, |
double2d_t & | result | ||
) |
FT in 2D.
Definition at line 59 of file FourierTransform.cpp.
References fftw_forward_FT(), FourierTransform::Workspace::h_fftw, init(), FourierTransform::Workspace::out_fftw, FourierTransform::Workspace::w_fftw, and ws.
void FourierTransform::fftshift | ( | double2d_t & | result | ) |
Shift low frequency to the center of 2D FT array.
Definition at line 96 of file FourierTransform.cpp.
References FourierTransform::Workspace::h_fftw, and ws.
void FourierTransform::init | ( | int | h_src, |
int | w_src | ||
) |
prepare arrays for 2D Fourier Transformation (FT) of the given vector
Definition at line 157 of file FourierTransform.cpp.
References FourierTransform::Workspace::clear(), FourierTransform::Workspace::h_fftw, FourierTransform::Workspace::h_src, FourierTransform::Workspace::in_src, FourierTransform::Workspace::out_fftw, FourierTransform::Workspace::p_forw_src, FourierTransform::Workspace::w_fftw, FourierTransform::Workspace::w_src, and ws.
Referenced by fft().
|
private |
compute FT of source using Fast Fourier transformation from fftw
Definition at line 192 of file FourierTransform.cpp.
References FourierTransform::Workspace::h_fftw, FourierTransform::Workspace::h_src, FourierTransform::Workspace::in_src, FourierTransform::Workspace::out_fftw, FourierTransform::Workspace::p_forw_src, FourierTransform::Workspace::w_fftw, FourierTransform::Workspace::w_src, and ws.
Referenced by fft().
|
private |
input and output data for fftw3
Definition at line 91 of file FourierTransform.h.
Referenced by fft(), fftshift(), fftw_forward_FT(), and init().