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