BornAgain  1.19.0
Simulate and fit neutron and x-ray scattering at grazing incidence
FourierTransform.cpp
Go to the documentation of this file.
1 // ************************************************************************************************
2 //
3 // BornAgain: simulate and fit reflection and scattering
4 //
5 //! @file Base/Math/FourierTransform.cpp
6 //! @brief Implements class 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 
17 #include <algorithm>
18 #include <cmath>
19 #include <iostream>
20 #include <sstream>
21 #include <stdexcept>
22 
24 
26  : h_src(0), w_src(0), h_fftw(0), w_fftw(0), in_src(0), out_fftw(0), p_forw_src(nullptr)
27 {
28 }
29 
31 {
32  clear();
33 }
34 
36 {
37  // rows (h) and columns (w) of the input 'source'
38  h_src = 0;
39  w_src = 0;
40 
41  if (in_src)
42  delete[] in_src;
43  in_src = 0;
44 
45  if (out_fftw)
46  fftw_free((fftw_complex*)out_fftw);
47  out_fftw = 0;
48 
49  if (p_forw_src != nullptr)
50  fftw_destroy_plan(p_forw_src);
51 
52  fftw_cleanup();
53 }
54 
55 /* ************************************************************************* */
56 // Fourier Transform in 2D
57 /* ************************************************************************* */
58 void FourierTransform::fft(const double2d_t& source, double2d_t& result)
59 {
60  // rows (h) and columns (w) of the input 'source'
61  int h_src = static_cast<int>(source.size());
62  int w_src = static_cast<int>((source.size() ? source[0].size() : 0));
63 
64  // initialisation
65  init(h_src, w_src);
66 
67  // Compute the forward FT
68  fftw_forward_FT(source);
69 
70  double* ptr = ws.out_fftw;
71  result.clear();
72 
73  // Resize the array for holding the FT output to correct dimensions
74  result.resize(static_cast<size_t>(ws.h_fftw),
75  std::vector<double>(static_cast<size_t>(ws.w_fftw)));
76 
77  // Storing FT output
78  for (size_t i = 0; i < static_cast<size_t>(ws.h_fftw); i++) {
79  size_t k = static_cast<size_t>(ws.h_fftw) - i;
80  if (i == 0)
81  k -= static_cast<size_t>(ws.h_fftw);
82  for (size_t j = 0; j < static_cast<size_t>(ws.w_fftw / 2 + 1); j++) {
83  result[i][j] = *ptr;
84  size_t l = static_cast<size_t>(ws.w_fftw) - j;
85  if (j != 0)
86  result[k][l] = result[i][j];
87  ptr += 2; // Only interested in the magnitudes of the complex Fourier coefficients
88  }
89  }
90 }
91 
92 /* ************************************************************************* */
93 // Fourier Transform 2D shift - center array around zero frequency
94 /* ************************************************************************* */
96 {
97  // Centering FT around zero frequency
98  size_t row_shift;
99  if (result.size() % 2 == 0)
100  row_shift = result.size() / 2;
101  else
102  row_shift = result.size() / 2 + 1;
103 
104  size_t col_shift;
105  if (result[0].size() % 2 == 0)
106  col_shift = result[0].size() / 2;
107  else
108  col_shift = result[0].size() / 2 + 1;
109 
110  // First, shifting the rows
111  std::rotate(result.begin(), result.begin() + static_cast<int>(row_shift), result.end());
112 
113  // Second, shifting the cols
114  for (size_t i = 0; i < static_cast<size_t>(ws.h_fftw); i++) {
115  std::rotate(result[i].begin(), result[i].begin() + static_cast<int>(col_shift),
116  result[i].end());
117  }
118 }
119 
120 /* ************************************************************************* */
121 // Fourier Transform in 1D
122 /* ************************************************************************* */
123 void FourierTransform::fft(const double1d_t& source, double1d_t& result)
124 {
125  // we simply create 2d arrays with length of first dimension equal to 1, and call 2d FT
126  double2d_t source2d;
127  source2d.push_back(source);
128 
129  double2d_t result2d;
130  fft(source2d, result2d);
131 
132  if (result2d.size() != 1)
133  throw std::runtime_error("FourierTransform::fft -> Panic in 1d");
134 
135  result = result2d[0];
136 }
137 
138 /* ************************************************************************* */
139 // Fourier Transform 1D shift - center around zero frequency
140 /* ************************************************************************* */
142 {
143  // Centering FT around zero frequency
144  size_t col_shift;
145  if (result.size() % 2 == 0)
146  col_shift = result.size() / 2;
147  else
148  col_shift = result.size() / 2 + 1;
149 
150  std::rotate(result.begin(), result.begin() + static_cast<int>(col_shift), result.end());
151 }
152 
153 /* ************************************************************************************ */
154 // initialise input and output arrays in workspace for fast Fourier transformation (fftw)
155 /* ************************************************************************************ */
156 void FourierTransform::init(int h_src, int w_src)
157 {
158  if (!h_src || !w_src) {
159  std::ostringstream os;
160  os << "FourierTransform::init() -> Panic! Wrong dimensions " << h_src << " " << w_src
161  << std::endl;
162  throw std::runtime_error(os.str());
163  }
164 
165  ws.clear();
166  ws.h_src = h_src;
167  ws.w_src = w_src;
168 
169  ws.h_fftw = h_src;
170  ws.w_fftw = w_src;
171 
172  ws.in_src = new double[static_cast<size_t>(ws.h_fftw * ws.w_fftw)];
173  ws.out_fftw = static_cast<double*>(
174  fftw_malloc(sizeof(fftw_complex) * static_cast<size_t>(ws.h_fftw * (ws.w_fftw / 2 + 1))));
175 
176  // Initialization of the plans
177  ws.p_forw_src = fftw_plan_dft_r2c_2d(ws.h_fftw, ws.w_fftw, ws.in_src,
178  (fftw_complex*)ws.out_fftw, FFTW_ESTIMATE);
179  // ws.p_forw_src = fftw_plan_dft_r2c_2d(ws.h_fftw, ws.w_fftw, ws.in_src,
180  // static_cast<double*>(ws.out_src), FFTW_ESTIMATE);
181 
182  if (ws.p_forw_src == nullptr)
183  throw std::runtime_error(
184  "FourierTransform::init() -> Error! Can't initialise p_forw_src plan.");
185 }
186 
187 /* ************************************************************************* */
188 // initialise input and output arrays for fast Fourier transformation
189 /* ************************************************************************* */
190 
192 {
193  if (ws.h_fftw <= 0 || ws.w_fftw <= 0)
194  throw std::runtime_error(
195  "FourierTransform::fftw_forward_FT() -> Panic! Initialisation is missed.");
196 
197  double *ptr, *ptr_end;
198 
199  // Initializing the content of ws.in_src to zero for all elements
200  for (ptr = ws.in_src, ptr_end = ws.in_src + ws.h_fftw * ws.w_fftw; ptr != ptr_end; ++ptr)
201  *ptr = 0.0;
202 
203  // Building the input signal for fftw algorithm
204  for (size_t row = 0; row < static_cast<size_t>(ws.h_src); ++row)
205  for (size_t col = 0; col < static_cast<size_t>(ws.w_src); ++col)
206  ws.in_src[(static_cast<int>(row) % ws.h_fftw) * ws.w_fftw
207  + (static_cast<int>(col) % ws.w_fftw)] += src[row][col];
208 
209  // Computing the FFT with fftw plan
210  fftw_execute(ws.p_forw_src);
211 
212  double re_out, im_out;
213  for (ptr = ws.out_fftw, ptr_end = ws.out_fftw + 2 * ws.h_fftw * (ws.w_fftw / 2 + 1);
214  ptr != ptr_end; ++ptr) {
215  re_out = *ptr;
216  im_out = *(++ptr);
217 
218  // magnitude
219  *(ptr - 1) = sqrt(re_out * re_out + im_out * im_out);
220  // phase
221  *ptr = atan2(im_out, re_out);
222  }
223 }
Defines class Math::FourierTransform.
double * out_fftw
result of Fourier transformation of source
int h_src
Here, h = height (# rows), w = width (# columns)
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