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
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 
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
92 };
93 
94 #endif // BORNAGAIN_CORE_INSTRUMENT_FOURIERTRANSFORM_H
Workspace for Fourier Transform.
double * out_fftw
result of Fourier transformation of source
int h_src
Here, h = height (# rows), w = width (# columns)
Fourier transform of vectors (in 1D or 2D) using Fast Fourier Transform (fftw package).
Workspace ws
input and output data for fftw3
void fftw_forward_FT(const double2d_t &source)
compute FT of source using Fast Fourier transformation from fftw
void fft(const double1d_t &source, double1d_t &result)
FT in 1D.
std::vector< double > double1d_t
definition of 1D vector of double
void fftshift(double1d_t &result)
Shift low frequency to the center of 1D FT array.
std::vector< double1d_t > double2d_t
definition of 2D vector of double
void init(int h_src, int w_src)
prepare arrays for 2D Fourier Transformation (FT) of the given vector