BornAgain
1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
|
Classes | |
class | Workspace |
Public Types | |
enum | EConvolutionMode { FFTW_LINEAR_FULL , FFTW_LINEAR_SAME_UNPADDED , FFTW_LINEAR_SAME , FFTW_LINEAR_VALID , FFTW_CIRCULAR_SAME , FFTW_CIRCULAR_SAME_SHIFTED , FFTW_UNDEFINED } |
typedef std::vector< double > | double1d_t |
typedef std::vector< double1d_t > | double2d_t |
Public Member Functions | |
Convolve () | |
void | fftconvolve (const double1d_t &source, const double1d_t &kernel, double1d_t &result) |
void | fftconvolve (const double2d_t &source, const double2d_t &kernel, double2d_t &result) |
void | init (int h_src, int w_src, int h_kernel, int w_kernel) |
void | setMode (EConvolutionMode mode) |
Private Member Functions | |
void | fftw_circular_convolution (const double2d_t &source, const double2d_t &kernel) |
int | find_closest_factor (int n) |
bool | is_optimal (int n) |
Private Attributes | |
Workspace | ws |
EConvolutionMode | m_mode |
std::vector< size_t > | m_implemented_factors |
Convolution of two real vectors (in 1D or 2D) using Fast Fourier Transform.
Usage: std::vector<double> signal, kernel, result; Convolve cv; cv.fftconvolve(signal, kernel, 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 33 of file Convolve.h.
typedef std::vector<double> Convolve::double1d_t |
definition of 1d vector of double
Definition at line 37 of file Convolve.h.
typedef std::vector<double1d_t> Convolve::double2d_t |
definition of 2d vector of double
Definition at line 40 of file Convolve.h.
convolution modes use LINEAR_SAME or CIRCULAR_SAME_SHIFTED for maximum performance
Enumerator | |
---|---|
FFTW_LINEAR_FULL | |
FFTW_LINEAR_SAME_UNPADDED | |
FFTW_LINEAR_SAME | |
FFTW_LINEAR_VALID | |
FFTW_CIRCULAR_SAME | |
FFTW_CIRCULAR_SAME_SHIFTED | |
FFTW_UNDEFINED |
Definition at line 46 of file Convolve.h.
Convolve::Convolve | ( | ) |
Definition at line 21 of file Convolve.cpp.
References m_implemented_factors.
void Convolve::fftconvolve | ( | const double1d_t & | source, |
const double1d_t & | kernel, | ||
double1d_t & | result | ||
) |
convolution in 1D
Definition at line 123 of file Convolve.cpp.
Referenced by ConvolutionDetectorResolution::apply1dConvolution(), and ConvolutionDetectorResolution::apply2dConvolution().
void Convolve::fftconvolve | ( | const double2d_t & | source, |
const double2d_t & | kernel, | ||
double2d_t & | result | ||
) |
convolution in 2D
Definition at line 87 of file Convolve.cpp.
References Convolve::Workspace::dst_fft, fftw_circular_convolution(), FFTW_CIRCULAR_SAME_SHIFTED, FFTW_LINEAR_SAME, FFTW_UNDEFINED, Convolve::Workspace::h_dst, Convolve::Workspace::h_fftw, Convolve::Workspace::h_kernel, Convolve::Workspace::h_offset, init(), m_mode, setMode(), Convolve::Workspace::w_dst, Convolve::Workspace::w_fftw, Convolve::Workspace::w_kernel, Convolve::Workspace::w_offset, and ws.
void Convolve::init | ( | int | h_src, |
int | w_src, | ||
int | h_kernel, | ||
int | w_kernel | ||
) |
prepare arrays for 2D convolution of given vectors
Definition at line 140 of file Convolve.cpp.
References Convolve::Workspace::clear(), Convolve::Workspace::dst_fft, FFTW_CIRCULAR_SAME, FFTW_CIRCULAR_SAME_SHIFTED, FFTW_LINEAR_FULL, FFTW_LINEAR_SAME, FFTW_LINEAR_SAME_UNPADDED, FFTW_LINEAR_VALID, find_closest_factor(), Convolve::Workspace::h_dst, Convolve::Workspace::h_fftw, Convolve::Workspace::h_kernel, Convolve::Workspace::h_offset, Convolve::Workspace::h_src, Convolve::Workspace::in_kernel, Convolve::Workspace::in_src, m_mode, Convolve::Workspace::out_kernel, Convolve::Workspace::out_src, Convolve::Workspace::p_back, Convolve::Workspace::p_forw_kernel, Convolve::Workspace::p_forw_src, Convolve::Workspace::w_dst, Convolve::Workspace::w_fftw, Convolve::Workspace::w_kernel, Convolve::Workspace::w_offset, Convolve::Workspace::w_src, and ws.
Referenced by fftconvolve().
|
inline |
Sets convolution mode.
Definition at line 66 of file Convolve.h.
References m_mode.
Referenced by fftconvolve().
|
private |
compute circual convolution of source and kernel using fast Fourier transformation
Definition at line 261 of file Convolve.cpp.
References Convolve::Workspace::dst_fft, Convolve::Workspace::h_fftw, Convolve::Workspace::h_kernel, Convolve::Workspace::h_src, Convolve::Workspace::in_kernel, Convolve::Workspace::in_src, Convolve::Workspace::out_kernel, Convolve::Workspace::out_src, Convolve::Workspace::p_back, Convolve::Workspace::p_forw_kernel, Convolve::Workspace::p_forw_src, Convolve::Workspace::w_fftw, Convolve::Workspace::w_kernel, Convolve::Workspace::w_src, and ws.
Referenced by fftconvolve().
|
private |
find closest number X>n that can be factorised according to fftw3 favorite factorisation
Definition at line 317 of file Convolve.cpp.
References is_optimal().
Referenced by init().
|
private |
if a number can be factorised using only favorite fftw3 factors
Definition at line 332 of file Convolve.cpp.
References m_implemented_factors.
Referenced by find_closest_factor().
|
private |
input and output data for fftw3
Definition at line 116 of file Convolve.h.
Referenced by fftconvolve(), fftw_circular_convolution(), and init().
|
private |
convolution mode
Definition at line 118 of file Convolve.h.
Referenced by fftconvolve(), init(), and setMode().
|
private |
Definition at line 119 of file Convolve.h.
Referenced by Convolve(), and is_optimal().