comparison liboctave/oct-fftw.h @ 7789:82be108cc558

First attempt at single precision tyeps * * * corrections to qrupdate single precision routines * * * prefer demotion to single over promotion to double * * * Add single precision support to log2 function * * * Trivial PROJECT file update * * * Cache optimized hermitian/transpose methods * * * Add tests for tranpose/hermitian and ChangeLog entry for new transpose code
author David Bateman <dbateman@free.fr>
date Sun, 27 Apr 2008 22:34:17 +0200
parents a1dbe9d80eee
children eb63fbe60fab
comparison
equal deleted inserted replaced
7788:45f5faba05a2 7789:82be108cc558
104 dim_vector rn; 104 dim_vector rn;
105 105
106 bool rsimd_align; 106 bool rsimd_align;
107 }; 107 };
108 108
109 class
110 OCTAVE_API
111 octave_float_fftw_planner
112 {
113 public:
114
115 octave_float_fftw_planner (void);
116
117 fftwf_plan create_plan (int dir, const int rank, const dim_vector dims,
118 octave_idx_type howmany, octave_idx_type stride, octave_idx_type dist,
119 const FloatComplex *in, FloatComplex *out);
120
121 fftwf_plan create_plan (const int rank, const dim_vector dims,
122 octave_idx_type howmany, octave_idx_type stride, octave_idx_type dist,
123 const float *in, FloatComplex *out);
124
125 enum FftwMethod {
126 UNKNOWN = -1,
127 ESTIMATE,
128 MEASURE,
129 PATIENT,
130 EXHAUSTIVE,
131 HYBRID
132 };
133
134 FftwMethod method (void);
135
136 FftwMethod method (FftwMethod _meth);
137
138 private:
139
140 FftwMethod meth;
141
142 // FIXME -- perhaps this should be split into two classes?
143
144 // Plan for fft and ifft of complex values
145 fftwf_plan plan[2];
146
147 // dist
148 octave_idx_type d[2];
149
150 // stride
151 octave_idx_type s[2];
152
153 // rank
154 int r[2];
155
156 // howmany
157 octave_idx_type h[2];
158
159 // dims
160 dim_vector n[2];
161
162 bool simd_align[2];
163 bool inplace[2];
164
165 // Plan for fft of real values
166 fftwf_plan rplan;
167
168 // dist
169 octave_idx_type rd;
170
171 // stride
172 octave_idx_type rs;
173
174 // rank
175 int rr;
176
177 // howmany
178 octave_idx_type rh;
179
180 // dims
181 dim_vector rn;
182
183 bool rsimd_align;
184 };
185
109 // FIXME -- maybe octave_fftw_planner should be a singleton object? 186 // FIXME -- maybe octave_fftw_planner should be a singleton object?
110 extern OCTAVE_API octave_fftw_planner fftw_planner; 187 extern OCTAVE_API octave_fftw_planner fftw_planner;
188 extern OCTAVE_API octave_float_fftw_planner float_fftw_planner;
111 189
112 class 190 class
113 OCTAVE_API 191 OCTAVE_API
114 octave_fftw 192 octave_fftw
115 { 193 {
125 static int fftNd (const Complex*, Complex*, const int, 203 static int fftNd (const Complex*, Complex*, const int,
126 const dim_vector &); 204 const dim_vector &);
127 static int ifftNd (const Complex*, Complex*, const int, 205 static int ifftNd (const Complex*, Complex*, const int,
128 const dim_vector &); 206 const dim_vector &);
129 207
208 static int fft (const float *in, FloatComplex *out, size_t npts,
209 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
210 static int fft (const FloatComplex *in, FloatComplex *out, size_t npts,
211 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
212 static int ifft (const FloatComplex *in, FloatComplex *out, size_t npts,
213 size_t nsamples = 1, octave_idx_type stride = 1, octave_idx_type dist = -1);
214
215 static int fftNd (const float*, FloatComplex*, const int, const dim_vector &);
216 static int fftNd (const FloatComplex*, FloatComplex*, const int,
217 const dim_vector &);
218 static int ifftNd (const FloatComplex*, FloatComplex*, const int,
219 const dim_vector &);
220
130 private: 221 private:
131 octave_fftw (); 222 octave_fftw ();
132 octave_fftw (const octave_fftw&); 223 octave_fftw (const octave_fftw&);
133 octave_fftw& operator = (const octave_fftw&); 224 octave_fftw& operator = (const octave_fftw&);
134 }; 225 };