BornAgain
1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FourierTransform.h
Go to the documentation of this file.
1
// ************************************************************************** //
2
//
3
// BornAgain: simulate and fit scattering at grazing incidence
4
//
5
//! @file Device/Instrument/FourierTransform.h
6
//! @brief Defines class MathFunctions::FourierTransform.
7
//!
8
//! @homepage http://www.bornagainproject.org
9
//! @license GNU General Public License v3 or higher (see COPYING)
10
//! @copyright Forschungszentrum Jülich GmbH 2015
11
//! @authors Scientific Computing Group at MLZ Garching
12
//! @authors C. Durniak, M. Ganeva, G. Pospelov, W. Van Herck, J. Wuttke
13
//
14
// ************************************************************************** //
15
16
#ifndef BORNAGAIN_CORE_INSTRUMENT_FOURIERTRANSFORM_H
17
#define BORNAGAIN_CORE_INSTRUMENT_FOURIERTRANSFORM_H
18
19
#include <fftw3.h>
20
#include <vector>
21
22
//! Fourier transform of vectors (in 1D or 2D) using Fast Fourier Transform (fftw package).
23
//!
24
//! Usage:
25
//! std::vector<double> signal, result;
26
//! FourierTransform ft;
27
//! ft.fft(signal, result)
28
//!
29
//! Given code rely on code from Jeremy Fix page, http://jeremy.fix.free.fr/spip.php?article15,
30
//! see also "Efficient convolution using the Fast Fourier Transform, Application in C++"
31
//! by Jeremy Fix, May 30, 2011
32
//!
33
//! @ingroup tools_internal
34
class
FourierTransform
35
{
36
public
:
37
//! definition of 1D vector of double
38
typedef
std::vector<double>
double1d_t
;
39
40
//! definition of 2D vector of double
41
typedef
std::vector<double1d_t>
double2d_t
;
42
43
FourierTransform
();
44
45
//! FT in 1D
46
void
fft
(
const
double1d_t
& source,
double1d_t
& result);
47
48
//! Shift low frequency to the center of 1D FT array
49
void
fftshift
(
double1d_t
& result);
50
51
//! FT in 2D
52
void
fft
(
const
double2d_t
& source,
double2d_t
& result);
53
54
//! Shift low frequency to the center of 2D FT array
55
void
fftshift
(
double2d_t
& result);
56
57
//! prepare arrays for 2D Fourier Transformation (FT) of the given vector
58
void
init
(
int
h_src,
int
w_src);
59
60
private
:
61
//! compute FT of source using Fast Fourier transformation from fftw
62
void
fftw_forward_FT
(
const
double2d_t
& source);
63
64
//! Workspace for Fourier Transform.
65
66
//! Workspace contains input (src), intermediate and output (out)
67
//! arrays to run FT via fft; 'source' is our signal
68
//! Output arrays are allocated via fftw3 allocation for maximum performance.
69
class
Workspace
70
{
71
public
:
72
Workspace
();
73
~Workspace
();
74
void
clear
();
75
friend
class
FourierTransform
;
76
77
private
:
78
//! Here, h = height (# rows), w = width (# columns)
79
int
h_src
,
w_src
;
// size of input 'source' array in 2D
80
int
h_fftw
,
w_fftw
;
// size of output 'FT' array in 2D
81
82
double
*
in_src
;
// pointer to input 'source' array
83
84
//! result of Fourier transformation of source
85
double
*
out_fftw
;
// pointer to output 'FT' array
86
87
fftw_plan
p_forw_src
;
88
};
89
90
//! input and output data for fftw3
91
Workspace
ws
;
92
};
93
94
#endif
// BORNAGAIN_CORE_INSTRUMENT_FOURIERTRANSFORM_H
FourierTransform::Workspace
Workspace for Fourier Transform.
Definition:
FourierTransform.h:70
FourierTransform::Workspace::h_fftw
int h_fftw
Definition:
FourierTransform.h:80
FourierTransform::Workspace::p_forw_src
fftw_plan p_forw_src
Definition:
FourierTransform.h:87
FourierTransform::Workspace::out_fftw
double * out_fftw
result of Fourier transformation of source
Definition:
FourierTransform.h:85
FourierTransform::Workspace::Workspace
Workspace()
Definition:
FourierTransform.cpp:26
FourierTransform::Workspace::w_fftw
int w_fftw
Definition:
FourierTransform.h:80
FourierTransform::Workspace::clear
void clear()
Definition:
FourierTransform.cpp:36
FourierTransform::Workspace::in_src
double * in_src
Definition:
FourierTransform.h:82
FourierTransform::Workspace::~Workspace
~Workspace()
Definition:
FourierTransform.cpp:31
FourierTransform::Workspace::w_src
int w_src
Definition:
FourierTransform.h:79
FourierTransform::Workspace::h_src
int h_src
Here, h = height (# rows), w = width (# columns)
Definition:
FourierTransform.h:79
FourierTransform
Fourier transform of vectors (in 1D or 2D) using Fast Fourier Transform (fftw package).
Definition:
FourierTransform.h:35
FourierTransform::ws
Workspace ws
input and output data for fftw3
Definition:
FourierTransform.h:91
FourierTransform::fftw_forward_FT
void fftw_forward_FT(const double2d_t &source)
compute FT of source using Fast Fourier transformation from fftw
Definition:
FourierTransform.cpp:192
FourierTransform::fft
void fft(const double1d_t &source, double1d_t &result)
FT in 1D.
Definition:
FourierTransform.cpp:124
FourierTransform::double1d_t
std::vector< double > double1d_t
definition of 1D vector of double
Definition:
FourierTransform.h:38
FourierTransform::fftshift
void fftshift(double1d_t &result)
Shift low frequency to the center of 1D FT array.
Definition:
FourierTransform.cpp:142
FourierTransform::double2d_t
std::vector< double1d_t > double2d_t
definition of 2D vector of double
Definition:
FourierTransform.h:41
FourierTransform::FourierTransform
FourierTransform()
FourierTransform::init
void init(int h_src, int w_src)
prepare arrays for 2D Fourier Transformation (FT) of the given vector
Definition:
FourierTransform.cpp:157
Device
Instrument
FourierTransform.h
Generated by
1.9.1