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