00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #if !defined (octave_oct_fftw_h)
00024 #define octave_oct_fftw_h 1
00025
00026 #include <cstddef>
00027
00028 #if defined (HAVE_FFTW3_H)
00029 #include <fftw3.h>
00030 #endif
00031
00032 #include "oct-cmplx.h"
00033 #include "dim-vector.h"
00034
00035 #if defined (HAVE_FFTW)
00036
00037 class
00038 OCTAVE_API
00039 octave_fftw_planner
00040 {
00041 protected:
00042
00043 octave_fftw_planner (void);
00044
00045 public:
00046
00047 ~octave_fftw_planner (void);
00048
00049 enum FftwMethod
00050 {
00051 UNKNOWN = -1,
00052 ESTIMATE,
00053 MEASURE,
00054 PATIENT,
00055 EXHAUSTIVE,
00056 HYBRID
00057 };
00058
00059 static bool instance_ok (void);
00060
00061 static fftw_plan
00062 create_plan (int dir, const int rank, const dim_vector dims,
00063 octave_idx_type howmany, octave_idx_type stride,
00064 octave_idx_type dist, const Complex *in,
00065 Complex *out)
00066 {
00067 static fftw_plan dummy;
00068
00069 return instance_ok ()
00070 ? instance->do_create_plan (dir, rank, dims, howmany, stride,
00071 dist, in, out)
00072 : dummy;
00073 }
00074
00075 static fftw_plan
00076 create_plan (const int rank, const dim_vector dims,
00077 octave_idx_type howmany, octave_idx_type stride,
00078 octave_idx_type dist, const double *in, Complex *out)
00079 {
00080 static fftw_plan dummy;
00081
00082 return instance_ok ()
00083 ? instance->do_create_plan (rank, dims, howmany, stride, dist, in, out)
00084 : dummy;
00085 }
00086
00087 static FftwMethod method (void)
00088 {
00089 static FftwMethod dummy;
00090
00091 return instance_ok () ? instance->do_method () : dummy;
00092 }
00093
00094 static FftwMethod method (FftwMethod _meth)
00095 {
00096 static FftwMethod dummy;
00097
00098 return instance_ok () ? instance->do_method (_meth) : dummy;
00099 }
00100
00101 private:
00102
00103
00104
00105 octave_fftw_planner (const octave_fftw_planner&);
00106
00107 octave_fftw_planner& operator = (const octave_fftw_planner&);
00108
00109 static octave_fftw_planner *instance;
00110
00111 static void cleanup_instance (void) { delete instance; instance = 0; }
00112
00113 fftw_plan
00114 do_create_plan (int dir, const int rank, const dim_vector dims,
00115 octave_idx_type howmany, octave_idx_type stride,
00116 octave_idx_type dist, const Complex *in,
00117 Complex *out);
00118
00119 fftw_plan
00120 do_create_plan (const int rank, const dim_vector dims,
00121 octave_idx_type howmany, octave_idx_type stride,
00122 octave_idx_type dist, const double *in, Complex *out);
00123
00124 FftwMethod do_method (void);
00125
00126 FftwMethod do_method (FftwMethod _meth);
00127
00128 FftwMethod meth;
00129
00130
00131
00132
00133 fftw_plan plan[2];
00134
00135
00136 octave_idx_type d[2];
00137
00138
00139 octave_idx_type s[2];
00140
00141
00142 int r[2];
00143
00144
00145 octave_idx_type h[2];
00146
00147
00148 dim_vector n[2];
00149
00150 bool simd_align[2];
00151 bool inplace[2];
00152
00153
00154 fftw_plan rplan;
00155
00156
00157 octave_idx_type rd;
00158
00159
00160 octave_idx_type rs;
00161
00162
00163 int rr;
00164
00165
00166 octave_idx_type rh;
00167
00168
00169 dim_vector rn;
00170
00171 bool rsimd_align;
00172 };
00173
00174 class
00175 OCTAVE_API
00176 octave_float_fftw_planner
00177 {
00178 protected:
00179
00180 octave_float_fftw_planner (void);
00181
00182 public:
00183
00184 ~octave_float_fftw_planner (void);
00185
00186 enum FftwMethod
00187 {
00188 UNKNOWN = -1,
00189 ESTIMATE,
00190 MEASURE,
00191 PATIENT,
00192 EXHAUSTIVE,
00193 HYBRID
00194 };
00195
00196 static bool instance_ok (void);
00197
00198 static fftwf_plan
00199 create_plan (int dir, const int rank, const dim_vector dims,
00200 octave_idx_type howmany, octave_idx_type stride,
00201 octave_idx_type dist, const FloatComplex *in,
00202 FloatComplex *out)
00203 {
00204 static fftwf_plan dummy;
00205
00206 return instance_ok ()
00207 ? instance->do_create_plan (dir, rank, dims, howmany, stride,
00208 dist, in, out)
00209 : dummy;
00210 }
00211
00212 static fftwf_plan
00213 create_plan (const int rank, const dim_vector dims,
00214 octave_idx_type howmany, octave_idx_type stride,
00215 octave_idx_type dist, const float *in, FloatComplex *out)
00216 {
00217 static fftwf_plan dummy;
00218
00219 return instance_ok ()
00220 ? instance->do_create_plan (rank, dims, howmany, stride, dist, in, out)
00221 : dummy;
00222 }
00223
00224 static FftwMethod method (void)
00225 {
00226 static FftwMethod dummy;
00227
00228 return instance_ok () ? instance->do_method () : dummy;
00229 }
00230
00231 static FftwMethod method (FftwMethod _meth)
00232 {
00233 static FftwMethod dummy;
00234
00235 return instance_ok () ? instance->do_method (_meth) : dummy;
00236 }
00237
00238 private:
00239
00240
00241
00242 octave_float_fftw_planner (const octave_float_fftw_planner&);
00243
00244 octave_float_fftw_planner& operator = (const octave_float_fftw_planner&);
00245
00246 static octave_float_fftw_planner *instance;
00247
00248 static void cleanup_instance (void) { delete instance; instance = 0; }
00249
00250 fftwf_plan
00251 do_create_plan (int dir, const int rank, const dim_vector dims,
00252 octave_idx_type howmany, octave_idx_type stride,
00253 octave_idx_type dist, const FloatComplex *in,
00254 FloatComplex *out);
00255
00256 fftwf_plan
00257 do_create_plan (const int rank, const dim_vector dims,
00258 octave_idx_type howmany, octave_idx_type stride,
00259 octave_idx_type dist, const float *in, FloatComplex *out);
00260
00261 FftwMethod do_method (void);
00262
00263 FftwMethod do_method (FftwMethod _meth);
00264
00265 FftwMethod meth;
00266
00267
00268
00269
00270 fftwf_plan plan[2];
00271
00272
00273 octave_idx_type d[2];
00274
00275
00276 octave_idx_type s[2];
00277
00278
00279 int r[2];
00280
00281
00282 octave_idx_type h[2];
00283
00284
00285 dim_vector n[2];
00286
00287 bool simd_align[2];
00288 bool inplace[2];
00289
00290
00291 fftwf_plan rplan;
00292
00293
00294 octave_idx_type rd;
00295
00296
00297 octave_idx_type rs;
00298
00299
00300 int rr;
00301
00302
00303 octave_idx_type rh;
00304
00305
00306 dim_vector rn;
00307
00308 bool rsimd_align;
00309 };
00310
00311 class
00312 OCTAVE_API
00313 octave_fftw
00314 {
00315 public:
00316
00317 static int fft (const double *in, Complex *out, size_t npts,
00318 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00319 static int fft (const Complex *in, Complex *out, size_t npts,
00320 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00321 static int ifft (const Complex *in, Complex *out, size_t npts,
00322 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00323
00324 static int fftNd (const double*, Complex*, const int, const dim_vector &);
00325 static int fftNd (const Complex*, Complex*, const int,
00326 const dim_vector &);
00327 static int ifftNd (const Complex*, Complex*, const int,
00328 const dim_vector &);
00329
00330 static int fft (const float *in, FloatComplex *out, size_t npts,
00331 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00332 static int fft (const FloatComplex *in, FloatComplex *out, size_t npts,
00333 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00334 static int ifft (const FloatComplex *in, FloatComplex *out, size_t npts,
00335 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00336
00337 static int fftNd (const float*, FloatComplex*, const int, const dim_vector &);
00338 static int fftNd (const FloatComplex*, FloatComplex*, const int,
00339 const dim_vector &);
00340 static int ifftNd (const FloatComplex*, FloatComplex*, const int,
00341 const dim_vector &);
00342
00343 private:
00344 octave_fftw (void);
00345 octave_fftw (const octave_fftw&);
00346 octave_fftw& operator = (const octave_fftw&);
00347 };
00348
00349 #endif
00350
00351 #endif