Subversion Repositories shark

Rev

Details | Last modification | View Log | RSS feed

Rev Author Line No. Line
2 pj 1
/* -*- C -*- */
2
/*
3
 * Copyright (c) 1997-1999 Massachusetts Institute of Technology
4
 *
5
 * This program is free software; you can redistribute it and/or modify
6
 * it under the terms of the GNU General Public License as published by
7
 * the Free Software Foundation; either version 2 of the License, or
8
 * (at your option) any later version.
9
 *
10
 * This program is distributed in the hope that it will be useful,
11
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
13
 * GNU General Public License for more details.
14
 *
15
 * You should have received a copy of the GNU General Public License
16
 * along with this program; if not, write to the Free Software
17
 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
18
 *
19
 */
20
 
21
/* fftw.h -- system-wide definitions */
22
/* $Id: fftw.h,v 1.1.1.1 2002-03-29 14:12:51 pj Exp $ */
23
 
24
#ifndef FFTW_H
25
#define FFTW_H
26
 
27
#include <stdlib.h>
28
//#include <stdio.h>
29
//#include "halsa.h"
30
 
31
#ifdef __cplusplus
32
extern "C" {
33
#endif                          /* __cplusplus */
34
 
35
/* Define for using single precision */
36
/*
37
 * If you can, use configure --enable-float instead of changing this
38
 * flag directly
39
 */
40
/* #undef FFTW_ENABLE_FLOAT */
41
 
42
/* our real numbers */
43
#ifdef FFTW_ENABLE_FLOAT
44
typedef float fftw_real;
45
#else
46
typedef double fftw_real;
47
#endif
48
 
49
/*********************************************
50
 * Complex numbers and operations
51
 *********************************************/
52
typedef struct {
53
     fftw_real re, im;
54
} fftw_complex;
55
 
56
#define c_re(c)  ((c).re)
57
#define c_im(c)  ((c).im)
58
 
59
typedef enum {
60
     FFTW_FORWARD = -1, FFTW_BACKWARD = 1
61
} fftw_direction;
62
 
63
/* backward compatibility with FFTW-1.3 */
64
typedef fftw_complex FFTW_COMPLEX;
65
typedef fftw_real FFTW_REAL;
66
 
67
#ifndef FFTW_1_0_COMPATIBILITY
68
#define FFTW_1_0_COMPATIBILITY 0
69
#endif
70
 
71
#if FFTW_1_0_COMPATIBILITY
72
/* backward compatibility with FFTW-1.0 */
73
#define REAL fftw_real
74
#define COMPLEX fftw_complex
75
#endif
76
 
77
/*********************************************
78
 * Success or failure status
79
 *********************************************/
80
 
81
typedef enum {
82
     FFTW_SUCCESS = 0, FFTW_FAILURE = -1
83
} fftw_status;
84
 
85
/*********************************************
86
 *              Codelets
87
 *********************************************/
88
typedef void (fftw_notw_codelet)
89
     (const fftw_complex *, fftw_complex *, int, int);
90
typedef void (fftw_twiddle_codelet)
91
     (fftw_complex *, const fftw_complex *, int,
92
      int, int);
93
typedef void (fftw_generic_codelet)
94
     (fftw_complex *, const fftw_complex *, int,
95
      int, int, int);
96
typedef void (fftw_real2hc_codelet)
97
     (const fftw_real *, fftw_real *, fftw_real *,
98
      int, int, int);
99
typedef void (fftw_hc2real_codelet)
100
     (const fftw_real *, const fftw_real *,
101
      fftw_real *, int, int, int);
102
typedef void (fftw_hc2hc_codelet)
103
     (fftw_real *, const fftw_complex *,
104
      int, int, int);
105
typedef void (fftw_rgeneric_codelet)
106
     (fftw_real *, const fftw_complex *, int,
107
      int, int, int);
108
 
109
/*********************************************
110
 *     Configurations
111
 *********************************************/
112
/*
113
 * A configuration is a database of all known codelets
114
 */
115
 
116
enum fftw_node_type {
117
     FFTW_NOTW, FFTW_TWIDDLE, FFTW_GENERIC, FFTW_RADER,
118
     FFTW_REAL2HC, FFTW_HC2REAL, FFTW_HC2HC, FFTW_RGENERIC
119
};
120
 
121
/* description of a codelet */
122
typedef struct {
123
     const char *name;          /* name of the codelet */
124
     void (*codelet) ();        /* pointer to the codelet itself */
125
     int size;                  /* size of the codelet */
126
     fftw_direction dir;        /* direction */
127
     enum fftw_node_type type;  /* TWIDDLE or NO_TWIDDLE */
128
     int signature;             /* unique id */
129
     int ntwiddle;              /* number of twiddle factors */
130
     const int *twiddle_order;  /*
131
                                 * array that determines the order
132
                                 * in which the codelet expects
133
                                 * the twiddle factors
134
                                 */
135
} fftw_codelet_desc;
136
 
137
/* On Win32, you need to do funny things to access global variables
138
   in shared libraries.  Thanks to Andrew Sterian for this hack. */
139
#if defined(__WIN32__) || defined(WIN32) || defined(_WINDOWS)
140
#  if defined(BUILD_FFTW_DLL)
141
#    define DL_IMPORT(type) __declspec(dllexport) type
142
#  elif defined(USE_FFTW_DLL)
143
#    define DL_IMPORT(type) __declspec(dllimport) type
144
#  else
145
#    define DL_IMPORT(type) type
146
#  endif
147
#else
148
#  define DL_IMPORT(type) type
149
#endif
150
 
151
extern DL_IMPORT(const char *) fftw_version;
152
 
153
/*****************************
154
 *        Plans
155
 *****************************/
156
/*
157
 * A plan is a sequence of reductions to compute a FFT of
158
 * a given size.  At each step, the FFT algorithm can:
159
 *
160
 * 1) apply a notw codelet, or
161
 * 2) recurse and apply a twiddle codelet, or
162
 * 3) apply the generic codelet.
163
 */
164
 
165
/* structure that contains twiddle factors */
166
typedef struct fftw_twiddle_struct {
167
     int n;
168
     const fftw_codelet_desc *cdesc;
169
     fftw_complex *twarray;
170
     struct fftw_twiddle_struct *next;
171
     int refcnt;
172
} fftw_twiddle;
173
 
174
typedef struct fftw_rader_data_struct {
175
     struct fftw_plan_struct *plan;
176
     fftw_complex *omega;
177
     int g, ginv;
178
     int p, flags, refcount;
179
     struct fftw_rader_data_struct *next;
180
     fftw_codelet_desc *cdesc;
181
} fftw_rader_data;
182
 
183
typedef void (fftw_rader_codelet)
184
     (fftw_complex *, const fftw_complex *, int,
185
      int, int, fftw_rader_data *);
186
 
187
/* structure that holds all the data needed for a given step */
188
typedef struct fftw_plan_node_struct {
189
     enum fftw_node_type type;
190
 
191
     union {
192
          /* nodes of type FFTW_NOTW */
193
          struct {
194
               int size;
195
               fftw_notw_codelet *codelet;
196
               const fftw_codelet_desc *codelet_desc;
197
          } notw;
198
 
199
          /* nodes of type FFTW_TWIDDLE */
200
          struct {
201
               int size;
202
               fftw_twiddle_codelet *codelet;
203
               fftw_twiddle *tw;
204
               struct fftw_plan_node_struct *recurse;
205
               const fftw_codelet_desc *codelet_desc;
206
          } twiddle;
207
 
208
          /* nodes of type FFTW_GENERIC */
209
          struct {
210
               int size;
211
               fftw_generic_codelet *codelet;
212
               fftw_twiddle *tw;
213
               struct fftw_plan_node_struct *recurse;
214
          } generic;
215
 
216
          /* nodes of type FFTW_RADER */
217
          struct {
218
               int size;
219
               fftw_rader_codelet *codelet;
220
               fftw_rader_data *rader_data;
221
               fftw_twiddle *tw;
222
               struct fftw_plan_node_struct *recurse;
223
          } rader;
224
 
225
          /* nodes of type FFTW_REAL2HC */
226
          struct {
227
               int size;
228
               fftw_real2hc_codelet *codelet;
229
               const fftw_codelet_desc *codelet_desc;
230
          } real2hc;
231
 
232
          /* nodes of type FFTW_HC2REAL */
233
          struct {
234
               int size;
235
               fftw_hc2real_codelet *codelet;
236
               const fftw_codelet_desc *codelet_desc;
237
          } hc2real;
238
 
239
          /* nodes of type FFTW_HC2HC */
240
          struct {
241
               int size;
242
               fftw_direction dir;
243
               fftw_hc2hc_codelet *codelet;
244
               fftw_twiddle *tw;
245
               struct fftw_plan_node_struct *recurse;
246
               const fftw_codelet_desc *codelet_desc;
247
          } hc2hc;
248
 
249
          /* nodes of type FFTW_RGENERIC */
250
          struct {
251
               int size;
252
               fftw_direction dir;
253
               fftw_rgeneric_codelet *codelet;
254
               fftw_twiddle *tw;
255
               struct fftw_plan_node_struct *recurse;
256
          } rgeneric;
257
     } nodeu;
258
 
259
     int refcnt;
260
} fftw_plan_node;
261
 
262
struct fftw_plan_struct {
263
     int n;
264
     int refcnt;
265
     fftw_direction dir;
266
     int flags;
267
     int wisdom_signature;
268
     enum fftw_node_type wisdom_type;
269
     struct fftw_plan_struct *next;
270
     fftw_plan_node *root;
271
     double cost;
272
};
273
 
274
/* a plan is just an array of instructions */
275
typedef struct fftw_plan_struct *fftw_plan;
276
 
277
/* flags for the planner */
278
#define  FFTW_ESTIMATE (0)
279
#define  FFTW_MEASURE  (1)
280
 
281
#define FFTW_OUT_OF_PLACE (0)
282
#define FFTW_IN_PLACE (8)
283
#define FFTW_USE_WISDOM (16)
284
 
285
#define FFTW_THREADSAFE (128)  /* guarantee plan is read-only so that the
286
                                  same plan can be used in parallel by
287
                                  multiple threads */
288
 
289
#define FFTWND_FORCE_BUFFERED (256)     /* internal, undocumented flag */
290
 
291
extern fftw_plan fftw_create_plan_specific(int n, fftw_direction dir,
292
                                           int flags,
293
                                           fftw_complex *in, int istride,
294
                                         fftw_complex *out, int ostride);
295
#define FFTW_HAS_PLAN_SPECIFIC
296
extern fftw_plan fftw_create_plan(int n, fftw_direction dir, int flags);
297
extern void fftw_print_plan(fftw_plan plan);
298
extern void fftw_destroy_plan(fftw_plan plan);
299
extern void fftw(fftw_plan plan, int howmany, fftw_complex *in, int istride,
300
                 int idist, fftw_complex *out, int ostride, int odist);
301
extern void fftw_one(fftw_plan plan, fftw_complex *in, fftw_complex *out);
302
extern void fftw_die(const char *s);
303
extern void *fftw_malloc(size_t n);
304
extern void fftw_free(void *p);
305
extern void fftw_check_memory_leaks(void);
306
extern void fftw_print_max_memory_usage(void);
307
 
308
typedef void *(*fftw_malloc_type_function) (size_t n);
309
typedef void  (*fftw_free_type_function) (void *p);
310
typedef void  (*fftw_die_type_function) (const char *errString);
311
extern DL_IMPORT(fftw_malloc_type_function) fftw_malloc_hook;
312
extern DL_IMPORT(fftw_free_type_function) fftw_free_hook;
313
extern DL_IMPORT(fftw_die_type_function) fftw_die_hook;
314
 
315
extern size_t fftw_sizeof_fftw_real(void);
316
 
317
/* Wisdom: */
318
/*
319
 * define this symbol so that users know we are using a version of FFTW
320
 * with wisdom
321
 */
322
//#define FFTW_HAS_WISDOM
323
extern void fftw_forget_wisdom(void);
324
extern void fftw_export_wisdom(void (*emitter) (char c, void *), void *data);
325
extern fftw_status fftw_import_wisdom(int (*g) (void *), void *data);
326
//extern void fftw_export_wisdom_to_file(FILE *output_file);
327
//extern fftw_status fftw_import_wisdom_from_file(FILE *input_file);
328
extern char *fftw_export_wisdom_to_string(void);
329
extern fftw_status fftw_import_wisdom_from_string(const char *input_string);
330
 
331
/*
332
 * define symbol so we know this function is available (it is not in
333
 * older FFTWs)
334
 */
335
//#define FFTW_HAS_FPRINT_PLAN
336
//extern void fftw_fprint_plan(FILE *f, fftw_plan plan);
337
 
338
/*****************************
339
 *    N-dimensional code
340
 *****************************/
341
typedef struct {
342
     int is_in_place;           /* 1 if for in-place FFTs, 0 otherwise */
343
 
344
     int rank;                  /*
345
                                 * the rank (number of dimensions) of the
346
                                 * array to be FFTed
347
                                 */
348
     int *n;                    /*
349
                                 * the dimensions of the array to the
350
                                 * FFTed
351
                                 */
352
     fftw_direction dir;
353
 
354
     int *n_before;             /*
355
                                 * n_before[i] = product of n[j] for j < i
356
                                 */
357
     int *n_after;              /* n_after[i] = product of n[j] for j > i */
358
 
359
     fftw_plan *plans;          /* 1d fftw plans for each dimension */
360
 
361
     int nbuffers, nwork;
362
     fftw_complex *work;        /*
363
                                 * work array big enough to hold
364
                                 * nbuffers+1 of the largest dimension
365
                                 * (has nwork elements)
366
                                 */
367
} fftwnd_data;
368
 
369
typedef fftwnd_data *fftwnd_plan;
370
 
371
/* Initializing the FFTWND plan: */
372
extern fftwnd_plan fftw2d_create_plan(int nx, int ny, fftw_direction dir,
373
                                      int flags);
374
extern fftwnd_plan fftw3d_create_plan(int nx, int ny, int nz,
375
                                      fftw_direction dir, int flags);
376
extern fftwnd_plan fftwnd_create_plan(int rank, const int *n,
377
                                      fftw_direction dir,
378
                                      int flags);
379
 
380
extern fftwnd_plan fftw2d_create_plan_specific(int nx, int ny,
381
                                               fftw_direction dir,
382
                                               int flags,
383
                                           fftw_complex *in, int istride,
384
                                         fftw_complex *out, int ostride);
385
extern fftwnd_plan fftw3d_create_plan_specific(int nx, int ny, int nz,
386
                                           fftw_direction dir, int flags,
387
                                           fftw_complex *in, int istride,
388
                                         fftw_complex *out, int ostride);
389
extern fftwnd_plan fftwnd_create_plan_specific(int rank, const int *n,
390
                                               fftw_direction dir,
391
                                               int flags,
392
                                           fftw_complex *in, int istride,
393
                                         fftw_complex *out, int ostride);
394
 
395
/* Freeing the FFTWND plan: */
396
extern void fftwnd_destroy_plan(fftwnd_plan plan);
397
 
398
/* Printing the plan: */
399
//extern void fftwnd_fprint_plan(FILE *f, fftwnd_plan p);
400
extern void fftwnd_print_plan(fftwnd_plan p);
401
#define FFTWND_HAS_PRINT_PLAN
402
 
403
/* Computing the N-Dimensional FFT */
404
extern void fftwnd(fftwnd_plan plan, int howmany,
405
                   fftw_complex *in, int istride, int idist,
406
                   fftw_complex *out, int ostride, int odist);
407
extern void fftwnd_one(fftwnd_plan p, fftw_complex *in, fftw_complex *out);
408
 
409
#ifdef __cplusplus
410
}                               /* extern "C" */
411
 
412
#endif                          /* __cplusplus */
413
#endif                          /* FFTW_H */