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 static octave_fftw_planner *instance;
00104
00105 fftw_plan
00106 do_create_plan (int dir, const int rank, const dim_vector dims,
00107 octave_idx_type howmany, octave_idx_type stride,
00108 octave_idx_type dist, const Complex *in,
00109 Complex *out);
00110
00111 fftw_plan
00112 do_create_plan (const int rank, const dim_vector dims,
00113 octave_idx_type howmany, octave_idx_type stride,
00114 octave_idx_type dist, const double *in, Complex *out);
00115
00116 FftwMethod do_method (void);
00117
00118 FftwMethod do_method (FftwMethod _meth);
00119
00120 FftwMethod meth;
00121
00122
00123
00124
00125 fftw_plan plan[2];
00126
00127
00128 octave_idx_type d[2];
00129
00130
00131 octave_idx_type s[2];
00132
00133
00134 int r[2];
00135
00136
00137 octave_idx_type h[2];
00138
00139
00140 dim_vector n[2];
00141
00142 bool simd_align[2];
00143 bool inplace[2];
00144
00145
00146 fftw_plan rplan;
00147
00148
00149 octave_idx_type rd;
00150
00151
00152 octave_idx_type rs;
00153
00154
00155 int rr;
00156
00157
00158 octave_idx_type rh;
00159
00160
00161 dim_vector rn;
00162
00163 bool rsimd_align;
00164 };
00165
00166 class
00167 OCTAVE_API
00168 octave_float_fftw_planner
00169 {
00170 protected:
00171
00172 octave_float_fftw_planner (void);
00173
00174 public:
00175
00176 ~octave_float_fftw_planner (void) { }
00177
00178 enum FftwMethod
00179 {
00180 UNKNOWN = -1,
00181 ESTIMATE,
00182 MEASURE,
00183 PATIENT,
00184 EXHAUSTIVE,
00185 HYBRID
00186 };
00187
00188 static bool instance_ok (void);
00189
00190 static fftwf_plan
00191 create_plan (int dir, const int rank, const dim_vector dims,
00192 octave_idx_type howmany, octave_idx_type stride,
00193 octave_idx_type dist, const FloatComplex *in,
00194 FloatComplex *out)
00195 {
00196 static fftwf_plan dummy;
00197
00198 return instance_ok ()
00199 ? instance->do_create_plan (dir, rank, dims, howmany, stride,
00200 dist, in, out)
00201 : dummy;
00202 }
00203
00204 static fftwf_plan
00205 create_plan (const int rank, const dim_vector dims,
00206 octave_idx_type howmany, octave_idx_type stride,
00207 octave_idx_type dist, const float *in, FloatComplex *out)
00208 {
00209 static fftwf_plan dummy;
00210
00211 return instance_ok ()
00212 ? instance->do_create_plan (rank, dims, howmany, stride, dist, in, out)
00213 : dummy;
00214 }
00215
00216 static FftwMethod method (void)
00217 {
00218 static FftwMethod dummy;
00219
00220 return instance_ok () ? instance->method () : dummy;
00221 }
00222
00223 static FftwMethod method (FftwMethod _meth)
00224 {
00225 static FftwMethod dummy;
00226
00227 return instance_ok () ? instance->method (_meth) : dummy;
00228 }
00229
00230 private:
00231
00232 static octave_float_fftw_planner *instance;
00233
00234 fftwf_plan
00235 do_create_plan (int dir, const int rank, const dim_vector dims,
00236 octave_idx_type howmany, octave_idx_type stride,
00237 octave_idx_type dist, const FloatComplex *in,
00238 FloatComplex *out);
00239
00240 fftwf_plan
00241 do_create_plan (const int rank, const dim_vector dims,
00242 octave_idx_type howmany, octave_idx_type stride,
00243 octave_idx_type dist, const float *in, FloatComplex *out);
00244
00245 FftwMethod do_method (void);
00246
00247 FftwMethod do_method (FftwMethod _meth);
00248
00249 FftwMethod meth;
00250
00251
00252
00253
00254 fftwf_plan plan[2];
00255
00256
00257 octave_idx_type d[2];
00258
00259
00260 octave_idx_type s[2];
00261
00262
00263 int r[2];
00264
00265
00266 octave_idx_type h[2];
00267
00268
00269 dim_vector n[2];
00270
00271 bool simd_align[2];
00272 bool inplace[2];
00273
00274
00275 fftwf_plan rplan;
00276
00277
00278 octave_idx_type rd;
00279
00280
00281 octave_idx_type rs;
00282
00283
00284 int rr;
00285
00286
00287 octave_idx_type rh;
00288
00289
00290 dim_vector rn;
00291
00292 bool rsimd_align;
00293 };
00294
00295 class
00296 OCTAVE_API
00297 octave_fftw
00298 {
00299 public:
00300
00301 static int fft (const double *in, Complex *out, size_t npts,
00302 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00303 static int fft (const Complex *in, Complex *out, size_t npts,
00304 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00305 static int ifft (const Complex *in, Complex *out, size_t npts,
00306 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00307
00308 static int fftNd (const double*, Complex*, const int, const dim_vector &);
00309 static int fftNd (const Complex*, Complex*, const int,
00310 const dim_vector &);
00311 static int ifftNd (const Complex*, Complex*, const int,
00312 const dim_vector &);
00313
00314 static int fft (const float *in, FloatComplex *out, size_t npts,
00315 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00316 static int fft (const FloatComplex *in, FloatComplex *out, size_t npts,
00317 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00318 static int ifft (const FloatComplex *in, FloatComplex *out, size_t npts,
00319 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
00320
00321 static int fftNd (const float*, FloatComplex*, const int, const dim_vector &);
00322 static int fftNd (const FloatComplex*, FloatComplex*, const int,
00323 const dim_vector &);
00324 static int ifftNd (const FloatComplex*, FloatComplex*, const int,
00325 const dim_vector &);
00326
00327 private:
00328 octave_fftw (void);
00329 octave_fftw (const octave_fftw&);
00330 octave_fftw& operator = (const octave_fftw&);
00331 };
00332
00333 #endif
00334
00335 #endif
00336
00337
00338
00339
00340
00341
00342