fft.c 4.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152
  1. /**
  2. * \file fft.c
  3. *
  4. * \brief Fast Fourier Transformation of Two Dimensional Satellite Data functions.
  5. *
  6. * This program is free software; you can redistribute it and/or modify
  7. * it under the terms of the GNU General Public License as published by
  8. * the Free Software Foundation; either version 2 of the License, or (at
  9. * your option) any later version.
  10. *
  11. * This program is distributed in the hope that it will be useful, but
  12. * WITHOUT ANY WARRANTY; without even the implied warranty of
  13. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  14. * General Public License for more details.
  15. *
  16. * You should have received a copy of the GNU General Public License
  17. * along with this program; if not, write to the Free Software
  18. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
  19. *
  20. * \author GRASS GIS Development Team
  21. *
  22. * \date 2001-2006
  23. */
  24. #include <grass/config.h>
  25. #if defined(HAVE_FFTW_H) || defined(HAVE_DFFTW_H) || defined(HAVE_FFTW3_H)
  26. #ifdef HAVE_FFTW_H
  27. #include <fftw.h>
  28. #endif
  29. #ifdef HAVE_DFFTW_H
  30. #include <dfftw.h>
  31. #endif
  32. #ifdef HAVE_FFTW3_H
  33. #include <fftw3.h>
  34. #define c_re(c) ((c)[0])
  35. #define c_im(c) ((c)[1])
  36. #endif
  37. #include <stdlib.h>
  38. #include <stdio.h>
  39. #include <math.h>
  40. #include <grass/gmath.h>
  41. #include <grass/gis.h>
  42. /**
  43. * \fn int fft2(int i_sign, double (*data)[2], int NN, int dimc, int dimr)
  44. *
  45. * \brief Fast Fourier Transform for two-dimensional array.
  46. *
  47. * Fast Fourier Transform for two-dimensional array.<br>
  48. * <bNote:</b> If passing real data to fft() forward transform
  49. * (especially when using fft() in a loop), explicitly (re-)initialize
  50. * the imaginary part to zero (DATA[1][i] = 0.0). Returns 0.
  51. *
  52. * \param[in] i_sign Direction of transform -1 is normal, +1 is inverse
  53. * \param[in,out] data Pointer to complex linear array in row major order
  54. * containing data and result
  55. * \param[in] NN Value of DATA dimension (dimc * dimr)
  56. * \param[in] dimc Value of image column dimension (max power of 2)
  57. * \param[in] dimr Value of image row dimension (max power of 2)
  58. * \return int always returns 0
  59. */
  60. int fft2(int i_sign, double (*data)[2], int NN, int dimc, int dimr)
  61. {
  62. #ifdef HAVE_FFTW3_H
  63. fftw_plan plan;
  64. #else
  65. fftwnd_plan plan;
  66. #endif
  67. double norm;
  68. int i;
  69. norm = 1.0 / sqrt(NN);
  70. #ifdef HAVE_FFTW3_H
  71. plan = fftw_plan_dft_2d(dimr, dimc, data, data,
  72. (i_sign < 0) ? FFTW_FORWARD : FFTW_BACKWARD,
  73. FFTW_ESTIMATE);
  74. fftw_execute(plan);
  75. fftw_destroy_plan(plan);
  76. #else
  77. plan = fftw2d_create_plan(dimc, dimr,
  78. (i_sign < 0) ? FFTW_FORWARD : FFTW_BACKWARD,
  79. FFTW_ESTIMATE | FFTW_IN_PLACE);
  80. fftwnd_one(plan, data, data);
  81. fftwnd_destroy_plan(plan);
  82. #endif
  83. for (i = 0; i < NN; i++) {
  84. data[i][0] *= norm;
  85. data[i][1] *= norm;
  86. }
  87. return 0;
  88. }
  89. /**
  90. * \fn int fft(int i_sign, double *DATA[2], int NN, int dimc, int dimr)
  91. *
  92. * \brief Fast Fourier Transform for two-dimensional array.
  93. *
  94. * Fast Fourier Transform for two-dimensional array.<br>
  95. * <bNote:</b> If passing real data to fft() forward transform
  96. * (especially when using fft() in a loop), explicitly (re-)initialize
  97. * the imaginary part to zero (DATA[1][i] = 0.0). Returns 0.
  98. *
  99. * \param[in] i_sign Direction of transform -1 is normal, +1 is inverse
  100. * \param[in,out] DATA Pointer to complex linear array in row major order
  101. * containing data and result
  102. * \param[in] NN Value of DATA dimension (dimc * dimr)
  103. * \param[in] dimc Value of image column dimension (max power of 2)
  104. * \param[in] dimr Value of image row dimension (max power of 2)
  105. * \return int always returns 0
  106. */
  107. int fft(int i_sign, double *DATA[2], int NN, int dimc, int dimr)
  108. {
  109. fftw_complex *data;
  110. int i;
  111. data = (fftw_complex *) G_malloc(NN * sizeof(fftw_complex));
  112. for (i = 0; i < NN; i++) {
  113. c_re(data[i]) = DATA[0][i];
  114. c_im(data[i]) = DATA[1][i];
  115. }
  116. fft2(i_sign, data, NN, dimc, dimr);
  117. for (i = 0; i < NN; i++) {
  118. DATA[0][i] = c_re(data[i]);
  119. DATA[1][i] = c_im(data[i]);
  120. }
  121. G_free(data);
  122. return 0;
  123. }
  124. #endif /* HAVE_FFT */