BornAgain
1.19.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 reflection and scattering
4
//
5
//! @file Base/Math/FourierTransform.h
6
//! @brief Defines class Math::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
#ifdef SWIG
17
#error no need to expose this header to Swig
18
#endif
19
20
#ifndef USER_API
21
#ifndef BORNAGAIN_BASE_MATH_FOURIERTRANSFORM_H
22
#define BORNAGAIN_BASE_MATH_FOURIERTRANSFORM_H
23
24
#include <fftw3.h>
25
#include <vector>
26
27
//! Fourier transform of vectors (in 1D or 2D) using Fast Fourier Transform (fftw package).
28
//!
29
//! Usage:
30
//! std::vector<double> signal, result;
31
//! FourierTransform ft;
32
//! ft.fft(signal, result)
33
//!
34
//! Given code rely on code from Jeremy Fix page, http://jeremy.fix.free.fr/spip.php?article15,
35
//! see also "Efficient convolution using the Fast Fourier Transform, Application in C++"
36
//! by Jeremy Fix, May 30, 2011
37
//!
38
//! @ingroup tools_internal
39
class
FourierTransform
{
40
public
:
41
//! definition of 1D vector of double
42
typedef
std::vector<double>
double1d_t
;
43
44
//! definition of 2D vector of double
45
typedef
std::vector<double1d_t>
double2d_t
;
46
47
FourierTransform
();
48
49
//! FT in 1D
50
void
fft
(
const
double1d_t
& source,
double1d_t
& result);
51
52
//! Shift low frequency to the center of 1D FT array
53
void
fftshift
(
double1d_t
& result);
54
55
//! FT in 2D
56
void
fft
(
const
double2d_t
& source,
double2d_t
& result);
57
58
//! Shift low frequency to the center of 2D FT array
59
void
fftshift
(
double2d_t
& result);
60
61
//! prepare arrays for 2D Fourier Transformation (FT) of the given vector
62
void
init
(
int
h_src,
int
w_src);
63
64
private
:
65
//! compute FT of source using Fast Fourier transformation from fftw
66
void
fftw_forward_FT
(
const
double2d_t
& source);
67
68
//! Workspace for Fourier Transform.
69
70
//! Workspace contains input (src), intermediate and output (out)
71
//! arrays to run FT via fft; 'source' is our signal
72
//! Output arrays are allocated via fftw3 allocation for maximum performance.
73
class
Workspace
{
74
public
:
75
Workspace
();
76
~Workspace
();
77
void
clear
();
78
friend
class
FourierTransform
;
79
80
private
:
81
//! Here, h = height (# rows), w = width (# columns)
82
int
h_src
,
w_src
;
// size of input 'source' array in 2D
83
int
h_fftw
,
w_fftw
;
// size of output 'FT' array in 2D
84
85
double
*
in_src
;
// pointer to input 'source' array
86
87
//! result of Fourier transformation of source
88
double
*
out_fftw
;
// pointer to output 'FT' array
89
90
fftw_plan
p_forw_src
;
91
};
92
93
//! input and output data for fftw3
94
Workspace
ws
;
95
};
96
97
#endif
// BORNAGAIN_BASE_MATH_FOURIERTRANSFORM_H
98
#endif
// USER_API
FourierTransform::Workspace
Workspace for Fourier Transform.
Definition:
FourierTransform.h:73
FourierTransform::Workspace::h_fftw
int h_fftw
Definition:
FourierTransform.h:83
FourierTransform::Workspace::p_forw_src
fftw_plan p_forw_src
Definition:
FourierTransform.h:90
FourierTransform::Workspace::out_fftw
double * out_fftw
result of Fourier transformation of source
Definition:
FourierTransform.h:88
FourierTransform::Workspace::Workspace
Workspace()
Definition:
FourierTransform.cpp:25
FourierTransform::Workspace::w_fftw
int w_fftw
Definition:
FourierTransform.h:83
FourierTransform::Workspace::clear
void clear()
Definition:
FourierTransform.cpp:35
FourierTransform::Workspace::in_src
double * in_src
Definition:
FourierTransform.h:85
FourierTransform::Workspace::~Workspace
~Workspace()
Definition:
FourierTransform.cpp:30
FourierTransform::Workspace::w_src
int w_src
Definition:
FourierTransform.h:82
FourierTransform::Workspace::h_src
int h_src
Here, h = height (# rows), w = width (# columns)
Definition:
FourierTransform.h:82
FourierTransform
Fourier transform of vectors (in 1D or 2D) using Fast Fourier Transform (fftw package).
Definition:
FourierTransform.h:39
FourierTransform::ws
Workspace ws
input and output data for fftw3
Definition:
FourierTransform.h:94
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:191
FourierTransform::fft
void fft(const double1d_t &source, double1d_t &result)
FT in 1D.
Definition:
FourierTransform.cpp:123
FourierTransform::double1d_t
std::vector< double > double1d_t
definition of 1D vector of double
Definition:
FourierTransform.h:42
FourierTransform::fftshift
void fftshift(double1d_t &result)
Shift low frequency to the center of 1D FT array.
Definition:
FourierTransform.cpp:141
FourierTransform::double2d_t
std::vector< double1d_t > double2d_t
definition of 2D vector of double
Definition:
FourierTransform.h:45
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:156
Base
Math
FourierTransform.h
Generated by
1.9.1