GRASS GIS 8 Programmer's Manual  8.2.2dev(2023)-3d2c704037
fft.c
Go to the documentation of this file.
1 
2 /**
3  * \file fft.c
4  *
5  * \brief Fast Fourier Transformation of Two Dimensional Satellite Data 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 GIS Development Team
22  *
23  * \date 2001-2006
24  */
25 
26 #include <grass/config.h>
27 
28 #if defined(HAVE_FFTW_H) || defined(HAVE_DFFTW_H) || defined(HAVE_FFTW3_H)
29 
30 #ifdef HAVE_FFTW_H
31 #include <fftw.h>
32 #endif
33 
34 #ifdef HAVE_DFFTW_H
35 #include <dfftw.h>
36 #endif
37 
38 #ifdef HAVE_FFTW3_H
39 #include <fftw3.h>
40 #define c_re(c) ((c)[0])
41 #define c_im(c) ((c)[1])
42 #endif
43 
44 #include <stdlib.h>
45 #include <stdio.h>
46 #include <math.h>
47 #include <grass/gmath.h>
48 #include <grass/gis.h>
49 
50 
51 /**
52  * \fn int fft2(int i_sign, double (*data)[2], int NN, int dimc, int dimr)
53  *
54  * \brief Fast Fourier Transform for two-dimensional array.
55  *
56  * Fast Fourier Transform for two-dimensional array.<br>
57  * <bNote:</b> If passing real data to fft() forward transform
58  * (especially when using fft() in a loop), explicitly (re-)initialize
59  * the imaginary part to zero (DATA[1][i] = 0.0). Returns 0.
60  *
61  * \param[in] i_sign Direction of transform -1 is normal, +1 is inverse
62  * \param[in,out] data Pointer to complex linear array in row major order
63  * containing data and result
64  * \param[in] NN Value of DATA dimension (dimc * dimr)
65  * \param[in] dimc Value of image column dimension (max power of 2)
66  * \param[in] dimr Value of image row dimension (max power of 2)
67  * \return int always returns 0
68  */
69 
70 int fft2(int i_sign, double (*data)[2], int NN, int dimc, int dimr)
71 {
72 #ifdef HAVE_FFTW3_H
73  fftw_plan plan;
74 #else
75  fftwnd_plan plan;
76 #endif
77  double norm;
78  int i;
79 
80  norm = 1.0 / sqrt(NN);
81 
82 #ifdef HAVE_FFTW3_H
83  plan = fftw_plan_dft_2d(dimr, dimc, data, data,
84  (i_sign < 0) ? FFTW_FORWARD : FFTW_BACKWARD,
85  FFTW_ESTIMATE);
86 
87  fftw_execute(plan);
88 
89  fftw_destroy_plan(plan);
90 #else
91  plan = fftw2d_create_plan(dimc, dimr,
92  (i_sign < 0) ? FFTW_FORWARD : FFTW_BACKWARD,
93  FFTW_ESTIMATE | FFTW_IN_PLACE);
94 
95  fftwnd_one(plan, data, data);
96 
97  fftwnd_destroy_plan(plan);
98 #endif
99 
100  for (i = 0; i < NN; i++) {
101  data[i][0] *= norm;
102  data[i][1] *= norm;
103  }
104 
105  return 0;
106 }
107 
108 /**
109  * \fn int fft(int i_sign, double *DATA[2], int NN, int dimc, int dimr)
110  *
111  * \brief Fast Fourier Transform for two-dimensional array.
112  *
113  * Fast Fourier Transform for two-dimensional array.<br>
114  * <bNote:</b> If passing real data to fft() forward transform
115  * (especially when using fft() in a loop), explicitly (re-)initialize
116  * the imaginary part to zero (DATA[1][i] = 0.0). Returns 0.
117  *
118  * \param[in] i_sign Direction of transform -1 is normal, +1 is inverse
119  * \param[in,out] DATA Pointer to complex linear array in row major order
120  * containing data and result
121  * \param[in] NN Value of DATA dimension (dimc * dimr)
122  * \param[in] dimc Value of image column dimension (max power of 2)
123  * \param[in] dimr Value of image row dimension (max power of 2)
124  * \return int always returns 0
125  */
126 
127 int fft(int i_sign, double *DATA[2], int NN, int dimc, int dimr)
128 {
129  fftw_complex *data;
130  int i;
131 
132  data = (fftw_complex *) G_malloc(NN * sizeof(fftw_complex));
133 
134  for (i = 0; i < NN; i++) {
135  c_re(data[i]) = DATA[0][i];
136  c_im(data[i]) = DATA[1][i];
137  }
138 
139  fft2(i_sign, data, NN, dimc, dimr);
140 
141  for (i = 0; i < NN; i++) {
142  DATA[0][i] = c_re(data[i]);
143  DATA[1][i] = c_im(data[i]);
144  }
145 
146  G_free(data);
147 
148  return 0;
149 }
150 
151 #endif /* HAVE_FFT */
#define G_malloc(n)
Definition: defs/gis.h:112
void G_free(void *)
Free allocated memory.
Definition: gis/alloc.c:149
int fft(int i_sign, double *DATA[2], int NN, int dimc, int dimr)
Fast Fourier Transform for two-dimensional array.
Definition: fft.c:127
#define c_im(c)
Definition: fft.c:41
int fft2(int i_sign, double(*data)[2], int NN, int dimc, int dimr)
Fast Fourier Transform for two-dimensional array.
Definition: fft.c:70
#define c_re(c)
Definition: fft.c:40