BornAgain  1.18.0
Simulate and fit neutron and x-ray scattering at grazing incidence
Convolve Class Reference
Collaboration diagram for Convolve:

Classes

class  Workspace
 

Public Types

enum  EConvolutionMode {
  FFTW_LINEAR_FULL , FFTW_LINEAR_SAME_UNPADDED , FFTW_LINEAR_SAME , FFTW_LINEAR_VALID ,
  FFTW_CIRCULAR_SAME , FFTW_CIRCULAR_SAME_SHIFTED , FFTW_UNDEFINED
}
 
typedef std::vector< double > double1d_t
 
typedef std::vector< double1d_tdouble2d_t
 

Public Member Functions

 Convolve ()
 
void fftconvolve (const double1d_t &source, const double1d_t &kernel, double1d_t &result)
 
void fftconvolve (const double2d_t &source, const double2d_t &kernel, double2d_t &result)
 
void init (int h_src, int w_src, int h_kernel, int w_kernel)
 
void setMode (EConvolutionMode mode)
 

Private Member Functions

void fftw_circular_convolution (const double2d_t &source, const double2d_t &kernel)
 
int find_closest_factor (int n)
 
bool is_optimal (int n)
 

Private Attributes

Workspace ws
 
EConvolutionMode m_mode
 
std::vector< size_t > m_implemented_factors
 

Detailed 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 33 of file Convolve.h.

Member Typedef Documentation

◆ double1d_t

typedef std::vector<double> Convolve::double1d_t

definition of 1d vector of double

Definition at line 37 of file Convolve.h.

◆ double2d_t

typedef std::vector<double1d_t> Convolve::double2d_t

definition of 2d vector of double

Definition at line 40 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 46 of file Convolve.h.

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

Constructor & Destructor Documentation

◆ Convolve()

Convolve::Convolve ( )

Definition at line 21 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:118
std::vector< size_t > m_implemented_factors
Definition: Convolve.h:119

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 123 of file Convolve.cpp.

124 {
125  // we simply create 2d arrays with length of first dimension equal to 1, and call 2d convolution
126  double2d_t source2d, kernel2d;
127  source2d.push_back(source);
128  kernel2d.push_back(kernel);
129 
130  double2d_t result2d;
131  fftconvolve(source2d, kernel2d, result2d);
132  if (result2d.size() != 1)
133  throw Exceptions::RuntimeErrorException("Convolve::fftconvolve -> Panic in 1d");
134  result = result2d[0];
135 }
std::vector< double1d_t > double2d_t
definition of 2d vector of double
Definition: Convolve.h:40
void fftconvolve(const double1d_t &source, const double1d_t &kernel, double1d_t &result)
convolution in 1D
Definition: Convolve.cpp:123

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 87 of file Convolve.cpp.

88 {
89  // set default convolution mode, if not defined
90  if (m_mode == FFTW_UNDEFINED)
92 
93  int h_src = (int)source.size();
94  int w_src = (int)(source.size() ? source[0].size() : 0);
95  int h_kernel = (int)kernel.size();
96  int w_kernel = (kernel.size() ? (int)kernel[0].size() : 0);
97 
98  // initialisation
99  init(h_src, w_src, h_kernel, w_kernel);
100  // Compute the circular convolution
101  fftw_circular_convolution(source, kernel);
102 
103  // results
104  result.clear();
105  result.resize(ws.h_dst);
106  for (int i = 0; i < ws.h_dst; i++) {
107  result[i].resize(ws.w_dst, 0);
108  for (int j = 0; j < ws.w_dst; j++) {
109  // result[i][j]=ws.dst[i*ws.w_dst+j];
111  result[i][j] = ws.dst_fft[((i + int(ws.h_kernel / 2.0)) % ws.h_fftw) * ws.w_fftw
112  + (j + int(ws.w_kernel / 2.0)) % ws.w_fftw];
113  } else {
114  result[i][j] = ws.dst_fft[(i + ws.h_offset) * ws.w_fftw + j + ws.w_offset];
115  }
116  }
117  }
118 }
double * dst_fft
result of product of FFT(source) and FFT(kernel)
Definition: Convolve.h:107
Workspace ws
input and output data for fftw3
Definition: Convolve.h:116
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 fftw_circular_convolution(const double2d_t &source, const double2d_t &kernel)
compute circual convolution of source and kernel using fast Fourier transformation
Definition: Convolve.cpp:261

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:

◆ 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 140 of file Convolve.cpp.

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

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:

◆ setMode()

void Convolve::setMode ( EConvolutionMode  mode)
inline

Sets convolution mode.

Definition at line 66 of file Convolve.h.

66 { m_mode = mode; }

References m_mode.

Referenced by fftconvolve().

◆ fftw_circular_convolution()

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

compute circual convolution of source and kernel using fast Fourier transformation

Definition at line 261 of file Convolve.cpp.

262 {
263  if (ws.h_fftw <= 0 || ws.w_fftw <= 0)
265  "Convolve::fftw_convolve() -> Panic! Initialisation is missed.");
266 
267  double *ptr, *ptr_end, *ptr2;
268 
269  // Reset the content of ws.in_src
270  for (ptr = ws.in_src, ptr_end = ws.in_src + ws.h_fftw * ws.w_fftw; ptr != ptr_end; ++ptr)
271  *ptr = 0.0;
272  for (ptr = ws.in_kernel, ptr_end = ws.in_kernel + ws.h_fftw * ws.w_fftw; ptr != ptr_end; ++ptr)
273  *ptr = 0.0;
274 
275  // Then we build our periodic signals
276  // ptr = src;
277  for (int i = 0; i < ws.h_src; ++i)
278  for (int j = 0; j < ws.w_src; ++j, ++ptr)
279  ws.in_src[(i % ws.h_fftw) * ws.w_fftw + (j % ws.w_fftw)] += src[i][j];
280 
281  // ptr = kernel;
282  for (int i = 0; i < ws.h_kernel; ++i)
283  for (int j = 0; j < ws.w_kernel; ++j, ++ptr)
284  ws.in_kernel[(i % ws.h_fftw) * ws.w_fftw + (j % ws.w_fftw)] += kernel[i][j];
285 
286  // And we compute their packed FFT
287  fftw_execute(ws.p_forw_src);
288  fftw_execute(ws.p_forw_kernel);
289 
290  // Compute the element-wise product on the packed terms
291  // Let's put the element wise products in ws.in_kernel
292  double re_s, im_s, re_k, im_k;
293  for (ptr = ws.out_src, ptr2 = ws.out_kernel,
294  ptr_end = ws.out_src + 2 * ws.h_fftw * (ws.w_fftw / 2 + 1);
295  ptr != ptr_end; ++ptr, ++ptr2) {
296  re_s = *ptr;
297  im_s = *(++ptr);
298  re_k = *ptr2;
299  im_k = *(++ptr2);
300  *(ptr2 - 1) = re_s * re_k - im_s * im_k;
301  *ptr2 = re_s * im_k + im_s * re_k;
302  }
303 
304  // Compute the backward FFT
305  // Carefull, The backward FFT does not preserve the output
306  fftw_execute(ws.p_back);
307  // Scale the transform
308  for (ptr = ws.dst_fft, ptr_end = ws.dst_fft + ws.w_fftw * ws.h_fftw; ptr != ptr_end; ++ptr)
309  *ptr /= double(ws.h_fftw * ws.w_fftw);
310 }

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 factorised according to fftw3 favorite factorisation

Definition at line 317 of file Convolve.cpp.

318 {
319  if (is_optimal(n)) {
320  return n;
321  } else {
322  int j = n + 1;
323  while (!is_optimal(j))
324  ++j;
325  return j;
326  }
327 }
bool is_optimal(int n)
if a number can be factorised using only favorite fftw3 factors
Definition: Convolve.cpp:332

References is_optimal().

Referenced by init().

Here is the call graph for this function:

◆ is_optimal()

bool Convolve::is_optimal ( int  n)
private

if a number can be factorised using only favorite fftw3 factors

Definition at line 332 of file Convolve.cpp.

333 {
334  if (n == 1)
335  return false;
336  size_t ntest = n;
337  for (size_t i = 0; i < m_implemented_factors.size(); i++) {
338  while ((ntest % m_implemented_factors[i]) == 0) {
339  ntest = ntest / m_implemented_factors[i];
340  }
341  }
342  return ntest == 1;
343 }

References m_implemented_factors.

Referenced by find_closest_factor().

Member Data Documentation

◆ ws

Workspace Convolve::ws
private

input and output data for fftw3

Definition at line 116 of file Convolve.h.

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

◆ m_mode

EConvolutionMode Convolve::m_mode
private

convolution mode

Definition at line 118 of file Convolve.h.

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

◆ m_implemented_factors

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

Definition at line 119 of file Convolve.h.

Referenced by Convolve(), and is_optimal().


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