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

Description

Convolution of two real vectors (in 1D or 2D) using Fast Fourier Transform.

Usage: std::vector<double> signal, kernel, result; Convolve cv; cv.fftconvolve(signal, kernel, 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 Convolve.h.

Collaboration diagram for Convolve:
[legend]

Classes

class  Workspace
 Workspace for Fourier convolution. 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...
 
enum  EConvolutionMode {
  FFTW_LINEAR_FULL , FFTW_LINEAR_SAME_UNPADDED , FFTW_LINEAR_SAME , FFTW_LINEAR_VALID ,
  FFTW_CIRCULAR_SAME , FFTW_CIRCULAR_SAME_SHIFTED , FFTW_UNDEFINED
}
 convolution modes use LINEAR_SAME or CIRCULAR_SAME_SHIFTED for maximum performance More...
 

Public Member Functions

 Convolve ()
 
void fftconvolve (const double1d_t &source, const double1d_t &kernel, double1d_t &result)
 convolution in 1D More...
 
void fftconvolve (const double2d_t &source, const double2d_t &kernel, double2d_t &result)
 convolution in 2D More...
 
void init (int h_src, int w_src, int h_kernel, int w_kernel)
 prepare arrays for 2D convolution of given vectors More...
 
void setMode (EConvolutionMode mode)
 Sets convolution mode. More...
 

Private Member Functions

void fftw_circular_convolution (const double2d_t &src, const double2d_t &kernel)
 compute circual convolution of source and kernel using fast Fourier transformation More...
 
int find_closest_factor (int n)
 find closest number X>n that can be factorized according to fftw3 favorite factorization More...
 
bool is_optimal (int n)
 if a number can be factorized using only favorite fftw3 factors More...
 

Private Attributes

std::vector< size_t > m_implemented_factors
 
EConvolutionMode m_mode
 convolution mode More...
 
Workspace ws
 input and output data for fftw3 More...
 

Member Typedef Documentation

◆ double1d_t

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

definition of 1d vector of double

Definition at line 39 of file Convolve.h.

◆ double2d_t

using Convolve::double2d_t = std::vector<double1d_t>

definition of 2d vector of double

Definition at line 42 of file Convolve.h.

Member Enumeration Documentation

◆ EConvolutionMode

convolution modes use LINEAR_SAME or CIRCULAR_SAME_SHIFTED for maximum performance

Enumerator
FFTW_LINEAR_FULL 
FFTW_LINEAR_SAME_UNPADDED 
FFTW_LINEAR_SAME 
FFTW_LINEAR_VALID 
FFTW_CIRCULAR_SAME 
FFTW_CIRCULAR_SAME_SHIFTED 
FFTW_UNDEFINED 

Definition at line 48 of file Convolve.h.

48  {
56  };
@ FFTW_LINEAR_FULL
Definition: Convolve.h:49
@ FFTW_CIRCULAR_SAME_SHIFTED
Definition: Convolve.h:54
@ FFTW_LINEAR_VALID
Definition: Convolve.h:52
@ FFTW_UNDEFINED
Definition: Convolve.h:55
@ FFTW_LINEAR_SAME_UNPADDED
Definition: Convolve.h:50
@ FFTW_CIRCULAR_SAME
Definition: Convolve.h:53
@ FFTW_LINEAR_SAME
Definition: Convolve.h:51

Constructor & Destructor Documentation

◆ Convolve()

Convolve::Convolve ( )

Definition at line 20 of file Convolve.cpp.

22 {
23  // storing favorite fftw3 prime factors
24  const size_t FFTW_FACTORS[] = {13, 11, 7, 5, 3, 2};
25  m_implemented_factors.assign(FFTW_FACTORS,
26  FFTW_FACTORS + sizeof(FFTW_FACTORS) / sizeof(FFTW_FACTORS[0]));
27 }
EConvolutionMode m_mode
convolution mode
Definition: Convolve.h:119
std::vector< size_t > m_implemented_factors
Definition: Convolve.h:120

References m_implemented_factors.

Member Function Documentation

◆ fftconvolve() [1/2]

void Convolve::fftconvolve ( const double1d_t source,
const double1d_t kernel,
double1d_t result 
)

convolution in 1D

Definition at line 115 of file Convolve.cpp.

116 {
117  // we simply create 2d arrays with length of first dimension equal to 1, and call 2d convolution
118  double2d_t source2d, kernel2d;
119  source2d.push_back(source);
120  kernel2d.push_back(kernel);
121 
122  double2d_t result2d;
123  fftconvolve(source2d, kernel2d, result2d);
124  if (result2d.size() != 1)
125  throw std::runtime_error("Convolve::fftconvolve -> Panic in 1d");
126  result = result2d[0];
127 }
std::vector< double1d_t > double2d_t
definition of 2d vector of double
Definition: Convolve.h:42
void fftconvolve(const double1d_t &source, const double1d_t &kernel, double1d_t &result)
convolution in 1D
Definition: Convolve.cpp:115

Referenced by ConvolutionDetectorResolution::apply1dConvolution(), and ConvolutionDetectorResolution::apply2dConvolution().

◆ fftconvolve() [2/2]

void Convolve::fftconvolve ( const double2d_t source,
const double2d_t kernel,
double2d_t result 
)

convolution in 2D

Definition at line 79 of file Convolve.cpp.

80 {
81  // set default convolution mode, if not defined
82  if (m_mode == FFTW_UNDEFINED)
84 
85  int h_src = (int)source.size();
86  int w_src = (int)(!source.empty() ? source[0].size() : 0);
87  int h_kernel = (int)kernel.size();
88  int w_kernel = (!kernel.empty() ? (int)kernel[0].size() : 0);
89 
90  // initialization
91  init(h_src, w_src, h_kernel, w_kernel);
92  // Compute the circular convolution
93  fftw_circular_convolution(source, kernel);
94 
95  // results
96  result.clear();
97  result.resize(ws.h_dst);
98  for (int i = 0; i < ws.h_dst; i++) {
99  result[i].resize(ws.w_dst, 0);
100  for (int j = 0; j < ws.w_dst; j++) {
101  // result[i][j]=ws.dst[i*ws.w_dst+j];
103  result[i][j] = ws.dst_fft[((i + int(ws.h_kernel / 2.0)) % ws.h_fftw) * ws.w_fftw
104  + (j + int(ws.w_kernel / 2.0)) % ws.w_fftw];
105  } else {
106  result[i][j] = ws.dst_fft[(i + ws.h_offset) * ws.w_fftw + j + ws.w_offset];
107  }
108  }
109  }
110 }
double * dst_fft
result of product of FFT(source) and FFT(kernel)
Definition: Convolve.h:108
Workspace ws
input and output data for fftw3
Definition: Convolve.h:117
void fftw_circular_convolution(const double2d_t &src, const double2d_t &kernel)
compute circual convolution of source and kernel using fast Fourier transformation
Definition: Convolve.cpp:250
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:132
void setMode(EConvolutionMode mode)
Sets convolution mode.
Definition: Convolve.h:68

References Convolve::Workspace::dst_fft, fftw_circular_convolution(), FFTW_CIRCULAR_SAME_SHIFTED, FFTW_LINEAR_SAME, FFTW_UNDEFINED, Convolve::Workspace::h_dst, Convolve::Workspace::h_fftw, Convolve::Workspace::h_kernel, Convolve::Workspace::h_offset, init(), m_mode, setMode(), Convolve::Workspace::w_dst, Convolve::Workspace::w_fftw, Convolve::Workspace::w_kernel, Convolve::Workspace::w_offset, and ws.

Here is the call graph for this function:

◆ fftw_circular_convolution()

void Convolve::fftw_circular_convolution ( const double2d_t src,
const double2d_t kernel 
)
private

compute circual convolution of source and kernel using fast Fourier transformation

Definition at line 250 of file Convolve.cpp.

251 {
252  if (ws.h_fftw <= 0 || ws.w_fftw <= 0)
253  throw std::runtime_error("Convolve::fftw_convolve() -> Panic! Initialization is missed.");
254 
255  double *ptr, *ptr_end, *ptr2;
256 
257  // Reset the content of ws.in_src
258  for (ptr = ws.in_src, ptr_end = ws.in_src + ws.h_fftw * ws.w_fftw; ptr != ptr_end; ++ptr)
259  *ptr = 0.0;
260  for (ptr = ws.in_kernel, ptr_end = ws.in_kernel + ws.h_fftw * ws.w_fftw; ptr != ptr_end; ++ptr)
261  *ptr = 0.0;
262 
263  // Then we build our periodic signals
264  // ptr = src;
265  for (int i = 0; i < ws.h_src; ++i)
266  for (int j = 0; j < ws.w_src; ++j, ++ptr)
267  ws.in_src[(i % ws.h_fftw) * ws.w_fftw + (j % ws.w_fftw)] += src[i][j];
268 
269  // ptr = kernel;
270  for (int i = 0; i < ws.h_kernel; ++i)
271  for (int j = 0; j < ws.w_kernel; ++j, ++ptr)
272  ws.in_kernel[(i % ws.h_fftw) * ws.w_fftw + (j % ws.w_fftw)] += kernel[i][j];
273 
274  // And we compute their packed FFT
275  fftw_execute(ws.p_forw_src);
276  fftw_execute(ws.p_forw_kernel);
277 
278  // Compute the element-wise product on the packed terms
279  // Let's put the elementwise products in ws.in_kernel
280  for (ptr = ws.out_src, ptr2 = ws.out_kernel,
281  ptr_end = ws.out_src + 2 * ws.h_fftw * (ws.w_fftw / 2 + 1);
282  ptr != ptr_end; ++ptr, ++ptr2) {
283  double re_s = *ptr;
284  double im_s = *(++ptr);
285  double re_k = *ptr2;
286  double im_k = *(++ptr2);
287  *(ptr2 - 1) = re_s * re_k - im_s * im_k;
288  *ptr2 = re_s * im_k + im_s * re_k;
289  }
290 
291  // Compute the backward FFT
292  // Carefull, The backward FFT does not preserve the output
293  fftw_execute(ws.p_back);
294  // Scale the transform
295  for (ptr = ws.dst_fft, ptr_end = ws.dst_fft + ws.w_fftw * ws.h_fftw; ptr != ptr_end; ++ptr)
296  *ptr /= double(ws.h_fftw * ws.w_fftw);
297 }
double * out_src
result of Fourier transformation of source
Definition: Convolve.h:102
fftw_plan p_back
Definition: Convolve.h:113
double * in_src
adjusted input 'source' array
Definition: Convolve.h:100
double * in_kernel
adjusted input 'kernel' array
Definition: Convolve.h:104
fftw_plan p_forw_src
Definition: Convolve.h:111
fftw_plan p_forw_kernel
Definition: Convolve.h:112
double * out_kernel
result of Fourier transformation of kernel
Definition: Convolve.h:106

References Convolve::Workspace::dst_fft, Convolve::Workspace::h_fftw, Convolve::Workspace::h_kernel, Convolve::Workspace::h_src, Convolve::Workspace::in_kernel, Convolve::Workspace::in_src, Convolve::Workspace::out_kernel, Convolve::Workspace::out_src, Convolve::Workspace::p_back, Convolve::Workspace::p_forw_kernel, Convolve::Workspace::p_forw_src, Convolve::Workspace::w_fftw, Convolve::Workspace::w_kernel, Convolve::Workspace::w_src, and ws.

Referenced by fftconvolve().

◆ find_closest_factor()

int Convolve::find_closest_factor ( int  n)
private

find closest number X>n that can be factorized according to fftw3 favorite factorization

Definition at line 304 of file Convolve.cpp.

305 {
306  if (is_optimal(n))
307  return n;
308  int j = n + 1;
309  while (!is_optimal(j))
310  ++j;
311  return j;
312 }
bool is_optimal(int n)
if a number can be factorized using only favorite fftw3 factors
Definition: Convolve.cpp:317

References is_optimal().

Referenced by init().

Here is the call graph for this function:

◆ init()

void Convolve::init ( int  h_src,
int  w_src,
int  h_kernel,
int  w_kernel 
)

prepare arrays for 2D convolution of given vectors

Definition at line 132 of file Convolve.cpp.

133 {
134  if (!h_src || !w_src || !h_kernel || !w_kernel) {
135  std::ostringstream os;
136  os << "Convolve::init() -> Panic! Wrong dimensions " << h_src << " " << w_src << " "
137  << h_kernel << " " << w_kernel << std::endl;
138  throw std::runtime_error(os.str());
139  }
140 
141  ws.clear();
142  ws.h_src = h_src;
143  ws.w_src = w_src;
144  ws.h_kernel = h_kernel;
145  ws.w_kernel = w_kernel;
146  switch (m_mode) {
147  case FFTW_LINEAR_FULL:
148  // Full Linear convolution
149  ws.h_fftw = find_closest_factor(h_src + h_kernel - 1);
150  ws.w_fftw = find_closest_factor(w_src + w_kernel - 1);
151  ws.h_dst = h_src + h_kernel - 1;
152  ws.w_dst = w_src + w_kernel - 1;
153  // Here we just keep the first [0:h_dst-1 ; 0:w_dst-1] real part elements of out_src
154  ws.h_offset = 0;
155  ws.w_offset = 0;
156  break;
158  // Same Linear convolution
159  ws.h_fftw = h_src + int(h_kernel / 2.0);
160  ws.w_fftw = w_src + int(w_kernel / 2.0);
161  ws.h_dst = h_src;
162  ws.w_dst = w_src;
163  // Here we just keep the first [h_filt/2:h_filt/2+h_dst-1 ; w_filt/2:w_filt/2+w_dst-1]
164  // real part elements of out_src
165  ws.h_offset = int(ws.h_kernel / 2.0);
166  ws.w_offset = int(ws.w_kernel / 2.0);
167  break;
168  case FFTW_LINEAR_SAME:
169  // Same Linear convolution
170  ws.h_fftw = find_closest_factor(h_src + int(h_kernel / 2.0));
171  ws.w_fftw = find_closest_factor(w_src + int(w_kernel / 2.0));
172  ws.h_dst = h_src;
173  ws.w_dst = w_src;
174  // Here we just keep the first [h_filt/2:h_filt/2+h_dst-1 ; w_filt/2:w_filt/2+w_dst-1]
175  // real part elements of out_src
176  ws.h_offset = int(ws.h_kernel / 2.0);
177  ws.w_offset = int(ws.w_kernel / 2.0);
178  break;
179  case FFTW_LINEAR_VALID:
180  // Valid Linear convolution
181  if (ws.h_kernel > ws.h_src || ws.w_kernel > ws.w_src) {
182  ws.h_fftw = 0;
183  ws.w_fftw = 0;
184  ws.h_dst = 0;
185  ws.w_dst = 0;
186  std::cout << "The 'valid' convolution results in an empty matrix" << std::endl;
187  throw std::runtime_error("The 'valid' convolution results in an empty matrix");
188  } else {
189  ws.h_fftw = find_closest_factor(h_src);
190  ws.w_fftw = find_closest_factor(w_src);
191  ws.h_dst = h_src - h_kernel + 1;
192  ws.w_dst = w_src - w_kernel + 1;
193  }
194  // Here we just take [h_dst x w_dst] elements starting at [h_kernel-1;w_kernel-1]
195  ws.h_offset = ws.h_kernel - 1;
196  ws.w_offset = ws.w_kernel - 1;
197  break;
198  case FFTW_CIRCULAR_SAME:
200  // Circular convolution, modulo N
201  // We don't padd with zeros because if we do, we need to padd with at least h_kernel/2;
202  // w_kernel/2 elements plus the wrapp around, which in facts leads to too much
203  // computations compared to the gain obtained with the optimal size
204  ws.h_fftw = h_src;
205  ws.w_fftw = w_src;
206  ws.h_dst = h_src;
207  ws.w_dst = w_src;
208  // We copy the first [0:h_dst-1 ; 0:w_dst-1] real part elements of out_src
209  ws.h_offset = 0;
210  ws.w_offset = 0;
211  break;
212  default:
213  std::cout
214  << "Unrecognized convolution mode, possible modes are "
215  << "FFTW_LINEAR_FULL, FFTW_LINEAR_SAME, FFTW_LINEAR_SAME_UNPADDED, FFTW_LINEAR_VALID, "
216  << "FFTW_CIRCULAR_SAME, FFTW_CIRCULAR_SHIFTED " << std::endl;
217  break;
218  }
219 
220  ws.in_src = new double[ws.h_fftw * ws.w_fftw];
221  ws.out_src = (double*)fftw_malloc(sizeof(fftw_complex) * ws.h_fftw * (ws.w_fftw / 2 + 1));
222  ws.in_kernel = new double[ws.h_fftw * ws.w_fftw];
223  ws.out_kernel = (double*)fftw_malloc(sizeof(fftw_complex) * ws.h_fftw * (ws.w_fftw / 2 + 1));
224 
225  ws.dst_fft = new double[ws.h_fftw * ws.w_fftw];
226  // ws.dst = new double[ws.h_dst * ws.w_dst];
227 
228  // Initialization of the plans
229  ws.p_forw_src = fftw_plan_dft_r2c_2d(ws.h_fftw, ws.w_fftw, ws.in_src, (fftw_complex*)ws.out_src,
230  FFTW_ESTIMATE);
231  if (ws.p_forw_src == nullptr)
232  throw std::runtime_error("Convolve::init() -> Error! Can't initialise p_forw_src plan.");
233 
234  ws.p_forw_kernel = fftw_plan_dft_r2c_2d(ws.h_fftw, ws.w_fftw, ws.in_kernel,
235  (fftw_complex*)ws.out_kernel, FFTW_ESTIMATE);
236  if (ws.p_forw_kernel == nullptr)
237  throw std::runtime_error("Convolve::init() -> Error! Can't initialise p_forw_kernel plan.");
238 
239  // The backward FFT takes ws.out_kernel as input
240  ws.p_back = fftw_plan_dft_c2r_2d(ws.h_fftw, ws.w_fftw, (fftw_complex*)ws.out_kernel, ws.dst_fft,
241  FFTW_ESTIMATE);
242  if (ws.p_back == nullptr)
243  throw std::runtime_error("Convolve::init() -> Error! Can't initialise p_back plan.");
244 }
int find_closest_factor(int n)
find closest number X>n that can be factorized according to fftw3 favorite factorization
Definition: Convolve.cpp:304

References Convolve::Workspace::clear(), Convolve::Workspace::dst_fft, FFTW_CIRCULAR_SAME, FFTW_CIRCULAR_SAME_SHIFTED, FFTW_LINEAR_FULL, FFTW_LINEAR_SAME, FFTW_LINEAR_SAME_UNPADDED, FFTW_LINEAR_VALID, find_closest_factor(), Convolve::Workspace::h_dst, Convolve::Workspace::h_fftw, Convolve::Workspace::h_kernel, Convolve::Workspace::h_offset, Convolve::Workspace::h_src, Convolve::Workspace::in_kernel, Convolve::Workspace::in_src, m_mode, Convolve::Workspace::out_kernel, Convolve::Workspace::out_src, Convolve::Workspace::p_back, Convolve::Workspace::p_forw_kernel, Convolve::Workspace::p_forw_src, Convolve::Workspace::w_dst, Convolve::Workspace::w_fftw, Convolve::Workspace::w_kernel, Convolve::Workspace::w_offset, Convolve::Workspace::w_src, and ws.

Referenced by fftconvolve().

Here is the call graph for this function:

◆ is_optimal()

bool Convolve::is_optimal ( int  n)
private

if a number can be factorized using only favorite fftw3 factors

Definition at line 317 of file Convolve.cpp.

318 {
319  if (n == 1)
320  return false;
321  size_t ntest = n;
322  for (size_t i = 0; i < m_implemented_factors.size(); i++) {
323  while ((ntest % m_implemented_factors[i]) == 0)
324  ntest = ntest / m_implemented_factors[i];
325  }
326  return ntest == 1;
327 }

References m_implemented_factors.

Referenced by find_closest_factor().

◆ setMode()

void Convolve::setMode ( EConvolutionMode  mode)
inline

Sets convolution mode.

Definition at line 68 of file Convolve.h.

68 { m_mode = mode; }

References m_mode.

Referenced by fftconvolve().

Member Data Documentation

◆ m_implemented_factors

std::vector<size_t> Convolve::m_implemented_factors
private

Definition at line 120 of file Convolve.h.

Referenced by Convolve(), and is_optimal().

◆ m_mode

EConvolutionMode Convolve::m_mode
private

convolution mode

Definition at line 119 of file Convolve.h.

Referenced by fftconvolve(), init(), and setMode().

◆ ws

Workspace Convolve::ws
private

input and output data for fftw3

Definition at line 117 of file Convolve.h.

Referenced by fftconvolve(), fftw_circular_convolution(), and init().


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