gVirtualXRay  2.0.10
VirtualX-RayImagingLibraryonGPU
FFT.inl
Go to the documentation of this file.
1 /*
2 
3 Copyright (c) 2016-2023, Dr Franck P. Vidal, Bangor University, All rights reserved.
4 Copyright (c) 2023-present, Prof Franck P. Vidal (franck.vidal@stfc.ac.uk),
5 UK Research and Innovation, All rights reserved.
6 
7 
8 Redistribution and use in source and binary forms, with or without modification,
9 are permitted provided that the following conditions are met:
10 
11 1. Redistributions of source code must retain the above copyright notice,
12 this list of conditions and the following disclaimer.
13 
14 2. Redistributions in binary form must reproduce the above copyright notice,
15 this list of conditions and the following disclaimer in the documentation and/or
16 other materials provided with the distribution.
17 
18 3. Neither the name of Bangor University, UK Research and Innovation nor the
19 names of its contributors may be used to endorse or promote products derived
20 from this software without specific prior written permission.
21 
22 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
23 ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO,
24 THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
25 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE
26 FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
27 DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR
28 SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
29 CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
30 OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF
31 THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 
33 */
34 
35 
61 //******************************************************************************
62 // Define
63 //******************************************************************************
64 
65 
66 //******************************************************************************
67 // Include
68 //******************************************************************************
69 #include <cmath>
70 #include <complex>
71 
72 #ifdef HAS_FFTW
73 #include <fftw3.h>
74 #endif
75 
76 #ifndef __Exception_h
77 #include "Exception.h"
78 #endif
79 
80 
81 //******************************************************************************
82 // namespace
83 //******************************************************************************
84 namespace gVirtualXRay {
85 
86 
87 //---------------------------------
88 template<typename T> FFT<T>::FFT():
89 //---------------------------------
90  m_origin_image_width(0),
91  m_origin_image_height(0)
92 //---------------------------------
93 {}
94 
95 
96 //---------------------------------------------------
97 template<typename T> FFT<T>::FFT(const FFT<T>& aFFT):
98 //---------------------------------------------------
99  m_real(aFFT.m_real),
100  m_img(aFFT.m_img),
101  m_origin_image_width(aFFT.m_origin_image_width),
102  m_origin_image_height(aFFT.m_origin_image_height)
103 //---------------------------------------------------
104 {}
105 
106 
107 //----------------------------------------------------------------
108 template<typename T> FFT<T>& FFT<T>::operator=(const FFT<T>& aFFT)
109 //----------------------------------------------------------------
110 {
111  m_real = aFFT.m_real;
112  m_img = aFFT.m_img;
113  m_origin_image_width = aFFT.m_origin_image_width;
114  m_origin_image_height = aFFT.m_origin_image_height;
115 
116  return (*this);
117 }
118 
119 
120 //---------------------------------------------------------------------
121 template<typename T> FFT<T> FFT<T>::computeFFT(const Image<T>& anImage)
122 //---------------------------------------------------------------------
123 {
124 #ifdef HAS_FFTW
125  FFT fft;
126  fft.m_origin_image_width = anImage.getWidth();
127  fft.m_origin_image_height = anImage.getHeight();
128 
129  unsigned int width_power_of_two(powerOfTwo(fft.m_origin_image_width));
130  unsigned int height_power_of_two(powerOfTwo(fft.m_origin_image_height));
131 
132  if (width_power_of_two != 1)
133  {
134  width_power_of_two = (std::max(width_power_of_two, height_power_of_two));
135  }
136 
137  if (height_power_of_two != 1)
138  {
139  height_power_of_two = (std::max(width_power_of_two, height_power_of_two));
140  }
141 
142  fft.m_real = Image<T>(width_power_of_two, height_power_of_two, anImage.getDepth());
143  fft.m_img = Image<T>(width_power_of_two, height_power_of_two, anImage.getDepth());
144 
145  fftw_complex *in_src, *out_src;
146  fftw_plan p_forw;
147 
148  in_src = reinterpret_cast<fftw_complex*>(fftw_malloc(sizeof(fftw_complex) * width_power_of_two * height_power_of_two));
149  out_src = reinterpret_cast<fftw_complex*>(fftw_malloc(sizeof(fftw_complex) * width_power_of_two * height_power_of_two));
150 
151  if (width_power_of_two == 1 || height_power_of_two == 1)
152  {
153  p_forw = fftw_plan_dft_1d(std::max(width_power_of_two, height_power_of_two), in_src, out_src, FFTW_FORWARD, FFTW_ESTIMATE);
154  }
155  else
156  {
157  p_forw = fftw_plan_dft_2d(width_power_of_two, height_power_of_two, in_src, out_src, FFTW_FORWARD, FFTW_ESTIMATE);
158  }
159 
160  //float max_value(anImage.getMaxValue());
161  bool y_forward_copy(false);
162  int j2(0);
163 
164  for (unsigned int k(0); k < anImage.getDepth(); ++k)
165  {
166  for (int j1 = 0 ; j1 < height_power_of_two ; ++j1)
167  {
168  if (!(j1 % fft.m_origin_image_height))
169  {
170  y_forward_copy = !y_forward_copy;
171 
172  if (y_forward_copy)
173  {
174  j2 = 0;
175  }
176  else
177  {
178  j2 = fft.m_origin_image_height - 1;
179  }
180  }
181 
182  bool x_forward_copy(false);
183  int i2(0);
184  for (int i1 = 0 ; i1 < width_power_of_two ; ++i1)
185  {
186  if (!(i1 % fft.m_origin_image_width))
187  {
188  x_forward_copy = !x_forward_copy;
189 
190  if (x_forward_copy)
191  {
192  i2 = 0;
193  }
194  else
195  {
196  i2 = fft.m_origin_image_width - 1;
197  }
198  }
199 
200  in_src[j1 * width_power_of_two + i1][0] = anImage(i2, j2, k);
201  in_src[j1 * width_power_of_two + i1][1] = 0.0;
202 
203  if (x_forward_copy)
204  {
205  ++i2;
206  }
207  else
208  {
209  --i2;
210  }
211  }
212 
213  if (y_forward_copy)
214  {
215  ++j2;
216  }
217  else
218  {
219  --j2;
220  }
221  }
222 
223  // Compute the forward fft
224  fftw_execute(p_forw);
225 
226  //double max_value(getMaxValue());
227  double square_size(width_power_of_two * height_power_of_two);
228 
229  for (int j = 0 ; j < height_power_of_two ; ++j)
230  {
231  for (int i = 0 ; i < width_power_of_two ; ++i)
232  {
233  fft.m_real(i, j, k) = out_src[j * width_power_of_two + i][0] / square_size;
234  fft.m_img( i, j, k) = out_src[j * width_power_of_two + i][1] / square_size;
235  }
236  }
237  }
238 
239  fftw_destroy_plan(p_forw);
240  fftw_free(in_src);
241  fftw_free(out_src);
242 
243  return (fft);
244 #else // HAS_FFTW
245 throw Exception(__FILE__, __FUNCTION__, __LINE__, "FFTW not supported");
246 #endif // HAS_FFTW
247 }
248 
249 
250 //---------------------------------------------------------
251 template<typename T> Image<T> FFT<T>::getInverseFFT() const
252 //---------------------------------------------------------
253 {
254 #ifdef HAS_FFTW
255  Image<T> temp(m_real.getWidth(), m_real.getHeight(), m_real.getDepth());
256 
257  // allocate input arrays
258  fftw_complex *in_src, *out_src;
259  fftw_plan p_back;
260 
261  in_src = reinterpret_cast<fftw_complex*>(fftw_malloc(sizeof(fftw_complex) * m_real.getWidth() * m_real.getHeight()));
262  out_src = reinterpret_cast<fftw_complex*>(fftw_malloc(sizeof(fftw_complex) * m_real.getWidth() * m_real.getHeight()));
263 
264 
265  if (m_real.getWidth() == 1 || m_real.getHeight() == 1)
266  {
267  p_back = fftw_plan_dft_1d(std::max(m_real.getWidth(), m_real.getHeight()), in_src, out_src, FFTW_BACKWARD, FFTW_ESTIMATE);
268  }
269  else
270  {
271  p_back = fftw_plan_dft_2d(m_real.getWidth(), m_real.getHeight(), in_src, out_src, FFTW_BACKWARD, FFTW_ESTIMATE);
272  }
273 
274  // transform magnitude/phase to real/imaginary
275  for (unsigned int k(0); k < m_real.getDepth(); ++k)
276  {
277  for (unsigned int j(0); j < m_real.getHeight(); ++j)
278  {
279  for (unsigned int i(0); i < m_real.getWidth(); ++i)
280  {
281  in_src[j * m_real.getWidth() + i][0] = m_real(i, j, k);
282  in_src[j * m_real.getWidth() + i][1] = m_img(i, j, k);
283  }
284  }
285 
286  // perform ift
287  fftw_execute(p_back);
288 
289  // save real parts to output
290  for (unsigned int j(0); j < m_real.getHeight(); ++j)
291  {
292  for (unsigned int i(0); i < m_real.getWidth(); ++i)
293  {
294  temp(i, j, k) = out_src[j * m_real.getWidth() + i][0];
295  }
296  }
297  }
298 
299 
300  // free memory
301  fftw_destroy_plan(p_back);
302  fftw_free(in_src); fftw_free(out_src);
303 
304  return (temp.getROI(0, 0, 0, m_origin_image_width, m_origin_image_height, m_real.getDepth()));
305 #else // HAS_FFTW
306  throw Exception(__FILE__, __FUNCTION__, __LINE__, "FFTW not supported");
307 #endif // HAS_FFTW
308 }
309 
310 
311 //----------------------------------------------------
312 template<typename T> FFT<T> FFT<T>::filterRamp() const
313 //----------------------------------------------------
314 {
315  FFT filtered_fft(*this);
316 
317  if (filtered_fft.m_real.getWidth() == 1 || filtered_fft.m_real.getHeight() == 1)
318  {
319  double half_size(std::max(filtered_fft.m_real.getWidth(), filtered_fft.m_real.getHeight()) / 2.0);
320 
321  for (unsigned int i(0); i < filtered_fft.m_real.getWidth() * filtered_fft.m_real.getHeight(); ++i)
322  {
323  double distance(std::abs(half_size - i));
324  distance /= half_size;
325 
326  filtered_fft.m_real[i] *= 1.0 - distance;
327  filtered_fft.m_img[ i] *= 1.0 - distance;
328  }
329  }
330  else
331  {
332  double half_diagonal_size(std::sqrt(std::pow(filtered_fft.m_real.getWidth(), 2) + std::pow(filtered_fft.m_real.getHeight(), 2)) / 2.0);
333  double half_size(filtered_fft.m_real.getWidth() / 2.0);
334 
335  for (unsigned int j(0); j < filtered_fft.m_real.getHeight(); ++j)
336  {
337  for (unsigned int i(0); i < filtered_fft.m_real.getWidth(); ++i)
338  {
339  double distance(std::sqrt(std::pow(half_size - i, 2) + std::pow(half_size - j, 2)));
340  distance /= half_diagonal_size;
341 
342  filtered_fft.m_real(i, j) *= 1.0 - distance;
343  filtered_fft.m_img(i, j) *= 1.0 - distance;
344  }
345  }
346  }
347 
348  return (filtered_fft);
349 }
350 
351 
352 template<typename T> FFT<T> FFT<T>::rhoFilter(const std::string& aFilterType, float aScalingFactor)
353 {
354  FFT<T> filtered_fft(*this);
355 
356  unsigned int size(filtered_fft.m_real.getWidth());
357  unsigned int half_size(size / 2);
358 
359  std::vector<double> p_filter(half_size, 1.0);
360 
361  // Apply the ramp filter (Ram-Lal)
362 
363  if ((aFilterType != "none" && aFilterType != "None" && aFilterType != "NONE") ||
364  aFilterType == "ramp" || aFilterType == "Ramp" || aFilterType == "RAMP")
365  {
366  for (unsigned int i(0); i < p_filter.size(); ++i)
367  {
368  double distance(std::abs(double(half_size - i)));
369  distance /= half_size;
370  distance = aScalingFactor * (1.0 - distance);
371 
372  p_filter[i] = distance;
373  }
374  }
375 
376  if (aFilterType == "hamming" || aFilterType == "Hamming" || aFilterType == "HAMMING")
377  {
378  int m(size);
379  --m;
380 
381  for (unsigned int i(0); i < half_size; ++i)
382  {
383  //unsigned int index_i(i + half_size);
384 
385  float c1 = 1.0 - (0.54 - 0.46 * std::cos (2.0 * Pi * float(i) / float(m)));
386  //float c2 = 0.54 - 0.46 * std::cos (2.0 * Pi * float(index_i) / float(m));
387 
388  p_filter[i] *= c1;
389  //p_filter[index_i] *= c1;
390  }
391  }
392  else if (aFilterType == "hann" || aFilterType == "Hann" || aFilterType == "HANN" || aFilterType == "hanning" || aFilterType == "Hanning" || aFilterType == "HANNING")
393  {
394  throw Exception(__FILE__, __FUNCTION__, __LINE__, "rho_filter: Unknown window type");
395  }
396  else if (aFilterType == "cosine" || aFilterType == "Cosine" || aFilterType == "COSINE")
397  {
398  //f = 0.5 * (0:length (rho) - 1)' / length (rho);
399  //filt = fftshift (sin (2 * pi * f));
400  throw Exception(__FILE__, __FUNCTION__, __LINE__, "rho_filter: Unknown window type");
401  }
402  else if (aFilterType == "shepp-logan" || aFilterType == "Shepp-logan" || aFilterType == "Shepp-Logan" || aFilterType == "SHEPP-LOGAN")
403  {
404  //f = (0:length (rho) / 2)' / length (rho);
405  //filt = sin (pi * f) ./ (pi * f);
406  //filt (1) = 1;
407  //filt = [filt; filt(end - 1:-1:2)];
408  throw Exception(__FILE__, __FUNCTION__, __LINE__, "rho_filter: Unknown window type");
409  }
410  else if (!(aFilterType == "ramp" || aFilterType == "Ramp" || aFilterType == "RAMP"))
411  {
412  throw Exception(__FILE__, __FUNCTION__, __LINE__, "rho_filter: Unknown window type");
413  }
414 
415  for (unsigned int k(0); k < filtered_fft.m_real.getDepth(); ++k)
416  {
417  for (unsigned int j(0); j < filtered_fft.m_real.getHeight(); ++j)
418  {
419  for (unsigned int i(0); i < p_filter.size(); ++i)
420  {
421  filtered_fft.m_real(i, j, k) *= p_filter[i];
422  filtered_fft.m_img( i, j, k) *= p_filter[i];
423 
424  filtered_fft.m_real(size - 1 - i, j, k) *= p_filter[i];
425  filtered_fft.m_img( size - 1 - i, j, k) *= p_filter[i];
426  }
427  }
428  }
429 
430  return (filtered_fft);
431 }
432 
433 
434 template<typename T> FFT<T> FFT<T>::rhoFilter(const char* aFilterType, float aScalingFactor)
435 {
436  return (rhoFilter(std::string(aFilterType), aScalingFactor));
437 }
438 
439 
440 //-----------------------------
441 template<typename T> Image<T> FFT<T>::getMagnitude() const
442 //-----------------------------
443 {
444  Image<T> temp(m_real.getWidth(), m_real.getHeight());
445 
446  double square_size(m_real.getWidth() * m_real.getHeight());
447 
448  for(int j = 0 ; j < m_real.getHeight() ; ++j)
449  {
450  for(int i = 0 ; i < m_real.getWidth() ; ++i)
451  {
452  unsigned int index(j * m_real.getWidth() + i);
453  double real(m_real[index] / square_size);
454  double img(m_img[index] / square_size);
455 
456  temp[index] = std::sqrt(real * real + img * img);
457  }
458  }
459 
460  return (swapQuadrants(temp));
461 }
462 
463 
464 
465 //-------------------------
466 template<typename T> Image<T> FFT<T>::getPhase() const
467 //-------------------------
468 {
469  Image<T> temp(m_real.getWidth(), m_real.getHeight());
470 
471  double square_size(m_real.getWidth() * m_real.getHeight());
472 
473  for(int j = 0 ; j < m_real.getHeight() ; ++j)
474  {
475  for(int i = 0 ; i < m_real.getWidth() ; ++i)
476  {
477  unsigned int index(j * m_real.getWidth() + i);
478  double real(m_real[index] / square_size);
479  double img(m_img[index] / square_size);
480 
481  // phase
482  std::complex<double> complex_value(real, img);
483  double phase_value = arg(complex_value) + Pi;
484 
485  temp[index] = (phase_value / (2.0 * Pi));
486  }
487  }
488 
489  return (swapQuadrants(temp));
490 }
491 
492 
493 template<typename T> unsigned int FFT<T>::powerOfTwo(unsigned int i)
494 {
495  unsigned int j(1);
496 
497  while (j < i) j *= 2;
498 
499  return (j);
500 }
501 
502 
503 template<typename T> Image<T> FFT<T>::swapQuadrants(const Image<T>& anImage)
504 {
505  Image<T> temp_image(anImage);
506 
507  // It is correct as the width is a power of two (even number)
508  int half_width(temp_image.getWidth() / 2);
509 
510  // It is correct as the height is a power of two (even number)
511  int half_height(temp_image.getHeight() / 2);
512 
513  if (temp_image.getWidth() == 1)
514  {
515  for(unsigned int i = 0; i < half_width; ++i)
516  {
517  // It is correct as the width is a power of two (even number)
518  unsigned int index_i(i + half_width);
519 
520  float temp(temp_image[i]);
521  temp_image[i] = temp_image[index_i];
522  temp_image[index_i] = temp;
523  }
524  }
525  else if (temp_image.getHeight() == 1)
526  {
527  for(unsigned int j = 0; j < half_height; ++j)
528  {
529  // It is correct as the height is a power of two (even number)
530  unsigned int index_j(j + half_height);
531 
532  float temp(temp_image[j * temp_image.getWidth()]);
533  temp_image[j * temp_image.getWidth()] = temp_image[index_j * temp_image.getWidth()];
534  temp_image[index_j * temp_image.getWidth()] = temp;
535  }
536  }
537  else
538  {
539  // swap quadrants diagonally
540  for(unsigned int j = 0; j < half_height; ++j)
541  {
542  // It is correct as the height is a power of two (even number)
543  unsigned int index_j(j + half_height);
544 
545  for(unsigned int i = 0; i < half_width; ++i)
546  {
547  // It is correct as the width is a power of two (even number)
548  unsigned int index_i(i + half_width);
549 
550  float temp(temp_image[j * temp_image.getWidth() + i]);
551  temp_image[j * temp_image.getWidth() + i] = temp_image[index_j * temp_image.getWidth() + index_i];
552  temp_image[index_j * temp_image.getWidth() + index_i] = temp;
553 
554  temp = temp_image[index_j * temp_image.getWidth() + i];
555  temp_image[index_j * temp_image.getWidth() + i] = temp_image[j * temp_image.getWidth() + index_i];
556  temp_image[j * temp_image.getWidth() + index_i] = temp;
557  }
558  }
559  }
560 
561  return (temp_image);
562 }
563 
564 
565 } // namespace gVirtualXRay
FFT()
Default constructor.
Definition: FFT.inl:88
static FFT< T > computeFFT(const Image< T > &anImage)
Definition: FFT.inl:121
FFT< T > rhoFilter(const std::string &aFilterType, float aScalingFactor)
Definition: FFT.inl:352
Exception is a class to handle exceptions.
Definition: Exception.h:109
Generic class to handle exceptions.
Image is a class to manage a greyscale image.
Definition: Image.h:92
Image< T > getPhase() const
Definition: FFT.inl:466
unsigned int getHeight() const
Number of pixels along the vertical axis.
Definition: Image.inl:792
unsigned int getDepth() const
Number of pixels along the Z axis.
Definition: Image.inl:800
const double m
meter
Definition: Units.h:104
unsigned int getWidth() const
Number of pixels along the horizontal axis.
Definition: Image.inl:784
Image< T > abs(const Image< T > &anImage)
Definition: Image.inl:5177
Image< T > getInverseFFT() const
Definition: FFT.inl:251
const double Pi
Pi.
FFT is a class to compute the FFT of a greyscale image.
Definition: FFT.h:78
Image< T > getMagnitude() const
Definition: FFT.inl:441
FFT< T > & operator=(const FFT< T > &aFFT)
Assignment operator (also called copy operator).
Definition: FFT.inl:108
FFT< T > filterRamp() const
Definition: FFT.inl:312