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