fft.c 5.55 KB
Newer Older
Emmanuel Bertin's avatar
Emmanuel Bertin committed
1
2
3
4
5
6
7
8
9
/*
                                  fft.c

*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
*
*       Part of:        A program that uses FFTs
*
*       Author:         E.BERTIN (IAP)
*
10
*       Contents:       Routines dealing with float precision FFT.
Emmanuel Bertin's avatar
Emmanuel Bertin committed
11
*
12
*       Last modify:    26/06/2009
Emmanuel Bertin's avatar
Emmanuel Bertin committed
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
*
*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
*/

#ifdef HAVE_CONFIG_H
#include "config.h"
#endif

#include <math.h>
#include <stdio.h>
#include <stdlib.h>

#include <fftw3.h>

#include "define.h"
#include "globals.h"
#include "fft.h"
#include "prefs.h"
#ifdef USE_THREADS
#include "threads.h"
#endif

Emmanuel Bertin's avatar
Emmanuel Bertin committed
35
 fftwf_plan	fplan, bplan;
Emmanuel Bertin's avatar
Emmanuel Bertin committed
36
37
38
39
 int    firsttimeflag;
#ifdef USE_THREADS
 pthread_mutex_t	fftmutex;
#endif
Emmanuel Bertin's avatar
Emmanuel Bertin committed
40
 fftwf_complex 	*fdata1;
Emmanuel Bertin's avatar
Emmanuel Bertin committed
41
42
43
44
45
46
47
48
49
#define SWAP(a,b)       tempr=(a);(a)=(b);(b)=tempr

/****** fft_init ************************************************************
PROTO	void fft_init(void)
PURPOSE	Initialize the FFT routines
INPUT	-.
OUTPUT	-.
NOTES	Global preferences are used for multhreading.
AUTHOR	E. Bertin (IAP)
50
VERSION	26/06/2009
Emmanuel Bertin's avatar
Emmanuel Bertin committed
51
 ***/
52
void    fft_init(int nthreads)
Emmanuel Bertin's avatar
Emmanuel Bertin committed
53
54
55
56
57
 {
  if (!firsttimeflag)
    {
#ifdef USE_THREADS
    QPTHREAD_MUTEX_INIT(&fftmutex, NULL);
58
#ifdef HAVE_FFTW_MP
59
60
61
62
    if (nthreads > 1)
      {
      if (!fftw_init_threads())
        error(EXIT_FAILURE, "*Error*: thread initialization failed in ", "FFTW");
63
      fftwf_plan_with_nthreads(prefs.nthreads);
64
      }
65
#endif
Emmanuel Bertin's avatar
Emmanuel Bertin committed
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
#endif
    firsttimeflag = 1;
    }

  return;
  }


/****** fft_end ************************************************************
PROTO	void fft_init(void)
PURPOSE	Clear up stuff set by FFT routines
INPUT	-.
OUTPUT	-.
NOTES	-.
AUTHOR	E. Bertin (IAP)
81
VERSION	26/06/2009
Emmanuel Bertin's avatar
Emmanuel Bertin committed
82
83
84
85
86
87
88
89
 ***/
void    fft_end(void)
 {

  if (firsttimeflag)
    {
    firsttimeflag = 0;
#ifdef USE_THREADS
90
#ifdef HAVE_FFTW_MP
91
    fftwf_cleanup_threads();
92
#endif
Emmanuel Bertin's avatar
Emmanuel Bertin committed
93
94
    QPTHREAD_MUTEX_DESTROY(&fftmutex);
#endif
95
    fftwf_cleanup();
Emmanuel Bertin's avatar
Emmanuel Bertin committed
96
97
98
99
100
101
    }

  return;
  }


Emmanuel Bertin's avatar
Emmanuel Bertin committed
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
/****** fft_reset ***********************************************************
PROTO	void fft_reset(void)
PURPOSE	Reset the FFT plans
INPUT	-.
OUTPUT	-.
NOTES	-.
AUTHOR	E. Bertin (IAP)
VERSION	26/06/2009
 ***/
void    fft_reset(void)
 {
  if (fplan)
    {
    QFFTWFREE(fdata1);
    fftwf_destroy_plan(fplan);
    }
  if (bplan)
    fftwf_destroy_plan(bplan);
  fplan = bplan = NULL;

  return;
  }


Emmanuel Bertin's avatar
Emmanuel Bertin committed
126
/****** fft_conv ************************************************************
127
PROTO	void fft_conv(float *data1, float *fdata2, int *size)
Emmanuel Bertin's avatar
Emmanuel Bertin committed
128
129
130
131
132
133
PURPOSE	Optimized 2-dimensional FFT convolution using the FFTW library.
INPUT	ptr to the first image,
	ptr to the Fourier transform of the second image,
	image size vector.
OUTPUT	-.
NOTES	For data1 and fdata2, memory must be allocated for
134
	size[0]* ... * 2*(size[naxis-1]/2+1) floats (padding required).
Emmanuel Bertin's avatar
Emmanuel Bertin committed
135
136
137
AUTHOR	E. Bertin (IAP)
VERSION	26/03/2007
 ***/
138
void    fft_conv(float *data1, float *fdata2, int *size)
Emmanuel Bertin's avatar
Emmanuel Bertin committed
139
  {
Emmanuel Bertin's avatar
Emmanuel Bertin committed
140
   float	*fdata1p,*fdata2p,
Emmanuel Bertin's avatar
Emmanuel Bertin committed
141
142
143
144
145
146
147
148
149
150
151
		real,imag, fac;
   int		i, npix,npix2;

/* Convert axis indexing to that of FFTW */
  npix = size[0]*size[1];
  npix2 = (((size[0]/2) + 1)*2) * size[1];

/* Forward FFT "in place" for data1 */
#ifdef USE_THREADS
  QPTHREAD_MUTEX_LOCK(&fftmutex);
#endif
Emmanuel Bertin's avatar
Emmanuel Bertin committed
152
153
154
155
156
157
  if (!fplan)
    {
    QFFTWMALLOC(fdata1, fftwf_complex, npix2);
    fplan = fftwf_plan_dft_r2c_2d(size[1], size[0], data1,
        (fftwf_complex *)fdata1, FFTW_MEASURE|FFTW_DESTROY_INPUT);
    }
Emmanuel Bertin's avatar
Emmanuel Bertin committed
158
159
160
#ifdef USE_THREADS
  QPTHREAD_MUTEX_UNLOCK(&fftmutex);
#endif
Emmanuel Bertin's avatar
Emmanuel Bertin committed
161
162
163
  fftwf_execute_dft_r2c(fplan, data1, fdata1);

//  fftwf_execute(plan);
Emmanuel Bertin's avatar
Emmanuel Bertin committed
164
165
166
167

#ifdef USE_THREADS
  QPTHREAD_MUTEX_LOCK(&fftmutex);
#endif
Emmanuel Bertin's avatar
Emmanuel Bertin committed
168
//  fftwf_destroy_plan(plan);
Emmanuel Bertin's avatar
Emmanuel Bertin committed
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
#ifdef USE_THREADS
  QPTHREAD_MUTEX_UNLOCK(&fftmutex);
#endif

/* Actual convolution (Fourier product) */
  fac = 1.0/npix;  
  fdata1p = fdata1;
  fdata2p = fdata2;
  for (i=npix2/2; i--; fdata2p+=2)
    {
    real = *fdata1p **fdata2p - *(fdata1p+1)**(fdata2p+1);
    imag = *(fdata1p+1)**fdata2p + *fdata1p**(fdata2p+1);
    *(fdata1p++) = fac*real;
    *(fdata1p++) = fac*imag;
    }

/* Reverse FFT */
#ifdef USE_THREADS
  QPTHREAD_MUTEX_LOCK(&fftmutex);
#endif
Emmanuel Bertin's avatar
Emmanuel Bertin committed
189
190
191
  if (!bplan)
    bplan = fftwf_plan_dft_c2r_2d(size[1], size[0], (fftwf_complex *)fdata1, 
        data1, FFTW_MEASURE|FFTW_DESTROY_INPUT);
Emmanuel Bertin's avatar
Emmanuel Bertin committed
192
193
194
#ifdef USE_THREADS
  QPTHREAD_MUTEX_UNLOCK(&fftmutex);
#endif
Emmanuel Bertin's avatar
Emmanuel Bertin committed
195
  fftwf_execute_dft_c2r(bplan, fdata1, data1);
Emmanuel Bertin's avatar
Emmanuel Bertin committed
196

Emmanuel Bertin's avatar
Emmanuel Bertin committed
197
//  fftwf_execute(plan);
Emmanuel Bertin's avatar
Emmanuel Bertin committed
198
199
200
201

#ifdef USE_THREADS
  QPTHREAD_MUTEX_LOCK(&fftmutex);
#endif
Emmanuel Bertin's avatar
Emmanuel Bertin committed
202
//  fftwf_destroy_plan(plan);
Emmanuel Bertin's avatar
Emmanuel Bertin committed
203
204
205
206
207
208
209
210
211
212
/* Free the fdata1 scratch array */
#ifdef USE_THREADS
  QPTHREAD_MUTEX_UNLOCK(&fftmutex);
#endif

  return;
  }


/****** fft_rtf ************************************************************
213
PROTO	float *fft_rtf(float *data, int *size)
Emmanuel Bertin's avatar
Emmanuel Bertin committed
214
215
216
217
218
219
220
221
PURPOSE	Optimized 2-dimensional FFT "in place" using the FFTW library.
INPUT	ptr to the image,
	ptr to image size vector.
OUTPUT	Pointer to the compressed, memory-allocated Fourier transform.
NOTES	Input data may end up corrupted.
AUTHOR	E. Bertin (IAP)
VERSION	26/03/2007
 ***/
222
float	*fft_rtf(float *data, int *size)
Emmanuel Bertin's avatar
Emmanuel Bertin committed
223
  {
224
225
226
   fftwf_plan   	plan;
   fftwf_complex	*fdata;
   int			npix2;
Emmanuel Bertin's avatar
Emmanuel Bertin committed
227
228
229
230
231
232
233
234

/* Convert axis indexing to that of FFTW */
  npix2 = ((size[0]/2) + 1) * size[1];

/* Forward FFT "in place" for data1 */
#ifdef USE_THREADS
  QPTHREAD_MUTEX_LOCK(&fftmutex);
#endif
235
236
  QFFTWMALLOC(fdata, fftwf_complex, npix2);
  plan = fftwf_plan_dft_r2c_2d(size[1], size[0], data,
Emmanuel Bertin's avatar
Emmanuel Bertin committed
237
238
239
240
        fdata, FFTW_ESTIMATE|FFTW_DESTROY_INPUT);
#ifdef USE_THREADS
  QPTHREAD_MUTEX_UNLOCK(&fftmutex);
#endif
241
  fftwf_execute(plan);
Emmanuel Bertin's avatar
Emmanuel Bertin committed
242
243
244
#ifdef USE_THREADS
  QPTHREAD_MUTEX_LOCK(&fftmutex);
#endif
245
  fftwf_destroy_plan(plan);
Emmanuel Bertin's avatar
Emmanuel Bertin committed
246
247
248
249
#ifdef USE_THREADS
  QPTHREAD_MUTEX_UNLOCK(&fftmutex);
#endif

250
  return (float *)fdata;
Emmanuel Bertin's avatar
Emmanuel Bertin committed
251
252
253
  }