BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
Convolve.h
Go to the documentation of this file.
1 // ************************************************************************** //
2 //
3 // BornAgain: simulate and fit scattering at grazing incidence
4 //
5 //! @file Device/Resolution/Convolve.h
6 //! @brief Defines class MathFunctions::Convolve.
7 //!
8 //! @homepage http://www.bornagainproject.org
9 //! @license GNU General Public License v3 or higher (see COPYING)
10 //! @copyright Forschungszentrum Jülich GmbH 2018
11 //! @authors Scientific Computing Group at MLZ (see CITATION, AUTHORS)
12 //
13 // ************************************************************************** //
14 
15 #ifndef BORNAGAIN_CORE_DETECTOR_CONVOLVE_H
16 #define BORNAGAIN_CORE_DETECTOR_CONVOLVE_H
17 
18 #include <fftw3.h>
19 #include <vector>
20 
21 //! Convolution of two real vectors (in 1D or 2D) using Fast Fourier Transform.
22 //!
23 //! Usage:
24 //! std::vector<double> signal, kernel, result;
25 //! Convolve cv;
26 //! cv.fftconvolve(signal, kernel, result)
27 //!
28 //! Given code rely on code from Jeremy Fix page, http://jeremy.fix.free.fr/spip.php?article15,
29 //! see also "Efficient convolution using the Fast Fourier Transform, Application in C++"
30 //! by Jeremy Fix, May 30, 2011
31 //!
32 //! @ingroup tools_internal
33 class Convolve
34 {
35 public:
36  //! definition of 1d vector of double
37  typedef std::vector<double> double1d_t;
38 
39  //! definition of 2d vector of double
40  typedef std::vector<double1d_t> double2d_t;
41 
42  Convolve();
43 
44  //! convolution modes
45  //! use LINEAR_SAME or CIRCULAR_SAME_SHIFTED for maximum performance
47  FFTW_LINEAR_FULL,
48  FFTW_LINEAR_SAME_UNPADDED,
49  FFTW_LINEAR_SAME,
50  FFTW_LINEAR_VALID,
51  FFTW_CIRCULAR_SAME,
52  FFTW_CIRCULAR_SAME_SHIFTED,
53  FFTW_UNDEFINED
54  };
55 
56  //! convolution in 1D
57  void fftconvolve(const double1d_t& source, const double1d_t& kernel, double1d_t& result);
58 
59  //! convolution in 2D
60  void fftconvolve(const double2d_t& source, const double2d_t& kernel, double2d_t& result);
61 
62  //! prepare arrays for 2D convolution of given vectors
63  void init(int h_src, int w_src, int h_kernel, int w_kernel);
64 
65  //! Sets convolution mode
66  void setMode(EConvolutionMode mode) { m_mode = mode; }
67 
68 private:
69  //! compute circual convolution of source and kernel using fast Fourier transformation
70  void fftw_circular_convolution(const double2d_t& source, const double2d_t& kernel);
71 
72  //! find closest number X>n that can be factorised according to fftw3 favorite factorisation
73  int find_closest_factor(int n);
74 
75  //! if a number can be factorised using only favorite fftw3 factors
76  bool is_optimal(int n);
77 
78  //! Workspace for Fourier convolution.
79 
80  //! Workspace contains input (source and kernel), intermediate and output
81  //! arrays to run convolution via fft; 'source' it is our signal, 'kernel'
82  //! it is our resolution function.
83  //! Sizes of input arrays are adjusted; output arrays are alocated via
84  //! fftw3 allocation for maximum performance.
85  class Workspace
86  {
87  public:
88  Workspace();
89  ~Workspace();
90  void clear();
91  friend class Convolve;
92 
93  private:
94  int h_src, w_src; // size of original 'source' array in 2 dimensions
95  int h_kernel, w_kernel; // size of original 'kernel' array in 2 dimensions
96  // size of adjusted source and kernel arrays (in_src, out_src, in_kernel, out_kernel)
97  int w_fftw, h_fftw;
98  //! adjusted input 'source' array
99  double* in_src;
100  //! result of Fourier transformation of source
101  double* out_src;
102  //! adjusted input 'kernel' array
103  double* in_kernel;
104  //! result of Fourier transformation of kernel
105  double* out_kernel;
106  //! result of product of FFT(source) and FFT(kernel)
107  double* dst_fft;
108  int h_dst, w_dst; // size of resulting array
109  int h_offset, w_offset; // offsets to copy result into output arrays
110  fftw_plan p_forw_src;
111  fftw_plan p_forw_kernel;
112  fftw_plan p_back;
113  };
114 
115  //! input and output data for fftw3
116  Workspace ws;
117  //! convolution mode
118  EConvolutionMode m_mode;
119  std::vector<size_t> m_implemented_factors; // favorite factorization terms of fftw3
120 };
121 
122 #endif // BORNAGAIN_CORE_DETECTOR_CONVOLVE_H
Convolution of two real vectors (in 1D or 2D) using Fast Fourier Transform.
Definition: Convolve.h:34
std::vector< double1d_t > double2d_t
definition of 2d vector of double
Definition: Convolve.h:40
std::vector< double > double1d_t
definition of 1d vector of double
Definition: Convolve.h:37
EConvolutionMode
convolution modes use LINEAR_SAME or CIRCULAR_SAME_SHIFTED for maximum performance
Definition: Convolve.h:46
void init(int h_src, int w_src, int h_kernel, int w_kernel)
prepare arrays for 2D convolution of given vectors
Definition: Convolve.cpp:140
void setMode(EConvolutionMode mode)
Sets convolution mode.
Definition: Convolve.h:66
void fftconvolve(const double1d_t &source, const double1d_t &kernel, double1d_t &result)
convolution in 1D
Definition: Convolve.cpp:123