BornAgain  1.19.79
Open-source research software to simulate and fit neutron and x-ray reflectometry and grazing-incidence small-angle scattering
FourierTransform Class Reference

Description

Fourier transform of vectors (in 1D or 2D) using Fast Fourier Transform (fftw package).

Usage: std::vector<double> signal, result; FourierTransform ft; ft.fft(signal, result)

Given code rely on code from Jeremy Fix page, http://jeremy.fix.free.fr/spip.php?article15, see also "Efficient convolution using the Fast Fourier Transform, Application in C++" by Jeremy Fix, May 30, 2011

Definition at line 36 of file FourierTransform.h.

Collaboration diagram for FourierTransform:
[legend]

Classes

class  Workspace
 Workspace for Fourier Transform. More...
 

Public Types

using double1d_t = std::vector< double >
 definition of 1D vector of double More...
 
using double2d_t = std::vector< double1d_t >
 definition of 2D vector of double More...
 

Public Member Functions

 FourierTransform ()
 
void fft (const double1d_t &source, double1d_t &result)
 FT in 1D. More...
 
void fft (const double2d_t &source, double2d_t &result)
 FT in 2D. More...
 
void fftshift (double1d_t &result) const
 Shift low frequency to the center of 1D FT array. More...
 
void fftshift (double2d_t &result) const
 Shift low frequency to the center of 2D FT array. More...
 
void init (int h_src, int w_src)
 prepare arrays for 2D Fourier Transformation (FT) of the given vector More...
 

Private Member Functions

void fftw_forward_FT (const double2d_t &src)
 compute FT of source using Fast Fourier transformation from fftw More...
 

Private Attributes

Workspace ws
 input and output data for fftw3 More...
 

Member Typedef Documentation

◆ double1d_t

using FourierTransform::double1d_t = std::vector<double>

definition of 1D vector of double

Definition at line 39 of file FourierTransform.h.

◆ double2d_t

definition of 2D vector of double

Definition at line 42 of file FourierTransform.h.

Constructor & Destructor Documentation

◆ FourierTransform()

FourierTransform::FourierTransform ( )
default

Member Function Documentation

◆ fft() [1/2]

void FourierTransform::fft ( const double1d_t source,
double1d_t result 
)

FT in 1D.

Definition at line 114 of file FourierTransform.cpp.

115 {
116  // we simply create 2d arrays with length of first dimension equal to 1, and call 2d FT
117  double2d_t source2d;
118  source2d.push_back(source);
119 
120  double2d_t result2d;
121  fft(source2d, result2d);
122 
123  ASSERT(result2d.size() == 1);
124  result = result2d[0];
125 }
#define ASSERT(condition)
Definition: Assert.h:45
void fft(const double1d_t &source, double1d_t &result)
FT in 1D.
std::vector< double1d_t > double2d_t
definition of 2D vector of double

References ASSERT.

◆ fft() [2/2]

void FourierTransform::fft ( const double2d_t source,
double2d_t result 
)

FT in 2D.

Definition at line 49 of file FourierTransform.cpp.

50 {
51  // rows (h) and columns (w) of the input 'source'
52  int h_src = static_cast<int>(source.size());
53  int w_src = static_cast<int>((!source.empty() ? source[0].size() : 0));
54 
55  // initialization
56  init(h_src, w_src);
57 
58  // Compute the forward FT
59  fftw_forward_FT(source);
60 
61  double* ptr = ws.out_fftw;
62  result.clear();
63 
64  const size_t nh = ws.h_fftw;
65 
66  // Resize the array for holding the FT output to correct dimensions
67  result.resize(nh, std::vector<double>(static_cast<size_t>(ws.w_fftw)));
68 
69  // Storing FT output
70  for (size_t i = 0; i < nh; i++) {
71  size_t k = nh - i;
72  if (i == 0)
73  k -= nh;
74  for (size_t j = 0; j < static_cast<size_t>(ws.w_fftw / 2 + 1); j++) {
75  result[i][j] = *ptr;
76  size_t l = static_cast<size_t>(ws.w_fftw) - j;
77  if (j != 0)
78  result[k][l] = result[i][j];
79  ptr += 2; // Only interested in the magnitudes of the complex Fourier coefficients
80  }
81  }
82 }
double * out_fftw
result of Fourier transformation of source
Workspace ws
input and output data for fftw3
void fftw_forward_FT(const double2d_t &src)
compute FT of source using Fast Fourier transformation from fftw
void init(int h_src, int w_src)
prepare arrays for 2D Fourier Transformation (FT) of the given vector

References fftw_forward_FT(), FourierTransform::Workspace::h_fftw, init(), FourierTransform::Workspace::out_fftw, FourierTransform::Workspace::w_fftw, and ws.

Here is the call graph for this function:

◆ fftshift() [1/2]

void FourierTransform::fftshift ( double1d_t result) const

Shift low frequency to the center of 1D FT array.

Definition at line 130 of file FourierTransform.cpp.

131 {
132  // Centering FT around zero frequency
133  size_t col_shift;
134  if (result.size() % 2 == 0)
135  col_shift = result.size() / 2;
136  else
137  col_shift = result.size() / 2 + 1;
138 
139  std::rotate(result.begin(), result.begin() + static_cast<int>(col_shift), result.end());
140 }

◆ fftshift() [2/2]

void FourierTransform::fftshift ( double2d_t result) const

Shift low frequency to the center of 2D FT array.

Definition at line 87 of file FourierTransform.cpp.

88 {
89  // Centering FT around zero frequency
90  size_t row_shift;
91  if (result.size() % 2 == 0)
92  row_shift = result.size() / 2;
93  else
94  row_shift = result.size() / 2 + 1;
95 
96  size_t col_shift;
97  if (result[0].size() % 2 == 0)
98  col_shift = result[0].size() / 2;
99  else
100  col_shift = result[0].size() / 2 + 1;
101 
102  // First, shifting the rows
103  std::rotate(result.begin(), result.begin() + static_cast<int>(row_shift), result.end());
104 
105  // Second, shifting the cols
106  for (size_t i = 0; i < static_cast<size_t>(ws.h_fftw); i++)
107  std::rotate(result[i].begin(), result[i].begin() + static_cast<int>(col_shift),
108  result[i].end());
109 }

References FourierTransform::Workspace::h_fftw, and ws.

◆ fftw_forward_FT()

void FourierTransform::fftw_forward_FT ( const double2d_t src)
private

compute FT of source using Fast Fourier transformation from fftw

Definition at line 174 of file FourierTransform.cpp.

175 {
176  ASSERT(ws.h_fftw > 0);
177  ASSERT(ws.w_fftw > 0);
178 
179  double *ptr, *ptr_end;
180 
181  // Initializing the content of ws.in_src to zero for all elements
182  for (ptr = ws.in_src, ptr_end = ws.in_src + ws.h_fftw * ws.w_fftw; ptr != ptr_end; ++ptr)
183  *ptr = 0.0;
184 
185  // Building the input signal for fftw algorithm
186  for (size_t row = 0; row < static_cast<size_t>(ws.h_src); ++row)
187  for (size_t col = 0; col < static_cast<size_t>(ws.w_src); ++col)
188  ws.in_src[(static_cast<int>(row) % ws.h_fftw) * ws.w_fftw
189  + (static_cast<int>(col) % ws.w_fftw)] += src[row][col];
190 
191  // Computing the FFT with fftw plan
192  fftw_execute(ws.p_forw_src);
193 
194  for (ptr = ws.out_fftw, ptr_end = ws.out_fftw + 2 * ws.h_fftw * (ws.w_fftw / 2 + 1);
195  ptr != ptr_end; ++ptr) {
196  double re_out = *ptr;
197  double im_out = *(++ptr);
198 
199  // magnitude
200  *(ptr - 1) = sqrt(re_out * re_out + im_out * im_out);
201  // phase
202  *ptr = atan2(im_out, re_out);
203  }
204 }
int h_src
Here, h = height (# rows), w = width (# columns)

References ASSERT, FourierTransform::Workspace::h_fftw, FourierTransform::Workspace::h_src, FourierTransform::Workspace::in_src, FourierTransform::Workspace::out_fftw, FourierTransform::Workspace::p_forw_src, FourierTransform::Workspace::w_fftw, FourierTransform::Workspace::w_src, and ws.

Referenced by fft().

◆ init()

void FourierTransform::init ( int  h_src,
int  w_src 
)

prepare arrays for 2D Fourier Transformation (FT) of the given vector

Definition at line 145 of file FourierTransform.cpp.

146 {
147  ASSERT(h_src);
148  ASSERT(w_src);
149 
150  ws.clear();
151  ws.h_src = h_src;
152  ws.w_src = w_src;
153 
154  ws.h_fftw = h_src;
155  ws.w_fftw = w_src;
156 
157  ws.in_src = new double[static_cast<size_t>(ws.h_fftw * ws.w_fftw)];
158  ws.out_fftw = static_cast<double*>(
159  fftw_malloc(sizeof(fftw_complex) * static_cast<size_t>(ws.h_fftw * (ws.w_fftw / 2 + 1))));
160 
161  // Initialization of the plans
162  ws.p_forw_src = fftw_plan_dft_r2c_2d(ws.h_fftw, ws.w_fftw, ws.in_src,
163  (fftw_complex*)ws.out_fftw, FFTW_ESTIMATE);
164  // ws.p_forw_src = fftw_plan_dft_r2c_2d(ws.h_fftw, ws.w_fftw, ws.in_src,
165  // static_cast<double*>(ws.out_src), FFTW_ESTIMATE);
166 
168 }

References ASSERT, FourierTransform::Workspace::clear(), FourierTransform::Workspace::h_fftw, FourierTransform::Workspace::h_src, FourierTransform::Workspace::in_src, FourierTransform::Workspace::out_fftw, FourierTransform::Workspace::p_forw_src, FourierTransform::Workspace::w_fftw, FourierTransform::Workspace::w_src, and ws.

Referenced by fft().

Here is the call graph for this function:

Member Data Documentation

◆ ws

Workspace FourierTransform::ws
private

input and output data for fftw3

Definition at line 91 of file FourierTransform.h.

Referenced by fft(), fftshift(), fftw_forward_FT(), and init().


The documentation for this class was generated from the following files: