GNU Octave 10.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
oct-fftw.h
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 2001-2025 The Octave Project Developers
4//
5// See the file COPYRIGHT.md in the top-level directory of this
6// distribution or <https://octave.org/copyright/>.
7//
8// This file is part of Octave.
9//
10// Octave is free software: you can redistribute it and/or modify it
11// under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// Octave is distributed in the hope that it will be useful, but
16// WITHOUT ANY WARRANTY; without even the implied warranty of
17// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18// GNU General Public License for more details.
19//
20// You should have received a copy of the GNU General Public License
21// along with Octave; see the file COPYING. If not, see
22// <https://www.gnu.org/licenses/>.
23//
24////////////////////////////////////////////////////////////////////////
25
26#if ! defined (octave_oct_fftw_h)
27#define octave_oct_fftw_h 1
28
29#include "octave-config.h"
30
31#include <cstddef>
32
33#include <string>
34
35#include "dim-vector.h"
36#include "oct-cmplx.h"
37
39
41{
42protected:
43
44 fftw_planner ();
45
46public:
47
48 OCTAVE_DISABLE_COPY_MOVE (fftw_planner)
49
51
53 {
54 UNKNOWN = -1,
59 HYBRID
60 };
61
62 static bool instance_ok ();
63
64 static void *
65 create_plan (int dir, const int rank, const dim_vector& dims,
66 octave_idx_type howmany, octave_idx_type stride,
67 octave_idx_type dist, const Complex *in,
68 Complex *out)
69 {
70 return instance_ok ()
71 ? s_instance->do_create_plan (dir, rank, dims, howmany, stride,
72 dist, in, out)
73 : nullptr;
74 }
75
76 static void *
77 create_plan (const int rank, const dim_vector& dims,
78 octave_idx_type howmany, octave_idx_type stride,
79 octave_idx_type dist, const double *in, Complex *out)
80 {
81 return instance_ok ()
82 ? s_instance->do_create_plan (rank, dims, howmany, stride, dist,
83 in, out)
84 : nullptr;
85 }
86
88 {
89 static FftwMethod dummy;
90
91 return instance_ok () ? s_instance->do_method () : dummy;
92 }
93
95 {
96 static FftwMethod dummy;
97
98 return instance_ok () ? s_instance->do_method (meth) : dummy;
99 }
100
101 static void threads (int nt);
102
103 static int threads ()
104 {
105 return instance_ok () ? s_instance->m_nthreads : 0;
106 }
107
108private:
109
110 static fftw_planner *s_instance;
111
112 static void cleanup_instance ()
113 { delete s_instance; s_instance = nullptr; }
114
115 void *
116 do_create_plan (int dir, const int rank, const dim_vector& dims,
117 octave_idx_type howmany, octave_idx_type stride,
118 octave_idx_type dist, const Complex *in,
119 Complex *out);
120
121 void *
122 do_create_plan (const int rank, const dim_vector& dims,
123 octave_idx_type howmany, octave_idx_type stride,
124 octave_idx_type dist, const double *in, Complex *out);
125
126 FftwMethod do_method ();
127
128 FftwMethod do_method (FftwMethod meth);
129
130 FftwMethod m_meth;
131
132 // FIXME: perhaps this should be split into two classes?
133
134 // Plan for fft and ifft of complex values
135 void *m_plan[2];
136
137 // dist
138 octave_idx_type m_d[2];
139
140 // stride
141 octave_idx_type m_s[2];
142
143 // rank
144 int m_r[2];
145
146 // howmany
147 octave_idx_type m_h[2];
148
149 // dims
150 dim_vector m_n[2];
151
152 bool m_simd_align[2];
153 bool m_inplace[2];
154
155 // Plan for fft of real values
156 void *m_rplan;
157
158 // dist
159 octave_idx_type m_rd;
160
161 // stride
162 octave_idx_type m_rs;
163
164 // rank
165 int m_rr;
166
167 // howmany
168 octave_idx_type m_rh;
169
170 // dims
171 dim_vector m_rn;
172
173 bool m_rsimd_align;
174
175 // number of threads. Always 1 unless compiled with multi-threading
176 // support.
177 int m_nthreads;
178};
179
181{
182protected:
183
185
186public:
187
188 OCTAVE_DISABLE_COPY_MOVE (float_fftw_planner)
189
191
193 {
194 UNKNOWN = -1,
199 HYBRID
200 };
201
202 static bool instance_ok ();
203
204 static void *
205 create_plan (int dir, const int rank, const dim_vector& dims,
206 octave_idx_type howmany, octave_idx_type stride,
207 octave_idx_type dist, const FloatComplex *in,
208 FloatComplex *out)
209 {
210 return instance_ok ()
211 ? s_instance->do_create_plan (dir, rank, dims, howmany, stride,
212 dist, in, out)
213 : nullptr;
214 }
215
216 static void *
217 create_plan (const int rank, const dim_vector& dims,
218 octave_idx_type howmany, octave_idx_type stride,
219 octave_idx_type dist, const float *in, FloatComplex *out)
220 {
221 return instance_ok ()
222 ? s_instance->do_create_plan (rank, dims, howmany, stride, dist,
223 in, out)
224 : nullptr;
225 }
226
228 {
229 static FftwMethod dummy;
230
231 return instance_ok () ? s_instance->do_method () : dummy;
232 }
233
235 {
236 static FftwMethod dummy;
237
238 return instance_ok () ? s_instance->do_method (meth) : dummy;
239 }
240
241 static void threads (int nt);
242
243 static int threads ()
244 {
245 return instance_ok () ? s_instance->m_nthreads : 0;
246 }
247
248private:
249
250 static float_fftw_planner *s_instance;
251
252 static void cleanup_instance ()
253 { delete s_instance; s_instance = nullptr; }
254
255 void *
256 do_create_plan (int dir, const int rank, const dim_vector& dims,
257 octave_idx_type howmany, octave_idx_type stride,
258 octave_idx_type dist, const FloatComplex *in,
259 FloatComplex *out);
260
261 void *
262 do_create_plan (const int rank, const dim_vector& dims,
263 octave_idx_type howmany, octave_idx_type stride,
264 octave_idx_type dist, const float *in, FloatComplex *out);
265
266 FftwMethod do_method ();
267
268 FftwMethod do_method (FftwMethod meth);
269
270 FftwMethod m_meth;
271
272 // FIXME: perhaps this should be split into two classes?
273
274 // Plan for fft and ifft of complex values
275 void *m_plan[2];
276
277 // dist
278 octave_idx_type m_d[2];
279
280 // stride
281 octave_idx_type m_s[2];
282
283 // rank
284 int m_r[2];
285
286 // howmany
287 octave_idx_type m_h[2];
288
289 // dims
290 dim_vector m_n[2];
291
292 bool m_simd_align[2];
293 bool m_inplace[2];
294
295 // Plan for fft of real values
296 void *m_rplan;
297
298 // dist
299 octave_idx_type m_rd;
300
301 // stride
302 octave_idx_type m_rs;
303
304 // rank
305 int m_rr;
306
307 // howmany
308 octave_idx_type m_rh;
309
310 // dims
311 dim_vector m_rn;
312
313 bool m_rsimd_align;
314
315 // number of threads. Always 1 unless compiled with multi-threading
316 // support.
317 int m_nthreads;
318};
319
321{
322public:
323
324 OCTAVE_DISABLE_CONSTRUCT_COPY_MOVE_DELETE (fftw)
325
326 static int fft (const double *in, Complex *out, std::size_t npts,
327 std::size_t nsamples = 1, octave_idx_type stride = 1,
328 octave_idx_type dist = -1);
329 static int fft (const Complex *in, Complex *out, std::size_t npts,
330 std::size_t nsamples = 1, octave_idx_type stride = 1,
331 octave_idx_type dist = -1);
332 static int ifft (const Complex *in, Complex *out, std::size_t npts,
333 std::size_t nsamples = 1, octave_idx_type stride = 1,
334 octave_idx_type dist = -1);
335
336 static int fftNd (const double *, Complex *, const int, const dim_vector&);
337 static int fftNd (const Complex *, Complex *, const int,
338 const dim_vector&);
339 static int ifftNd (const Complex *, Complex *, const int,
340 const dim_vector&);
341
342 static int fft (const float *in, FloatComplex *out, std::size_t npts,
343 std::size_t nsamples = 1, octave_idx_type stride = 1,
344 octave_idx_type dist = -1);
345 static int fft (const FloatComplex *in, FloatComplex *out, std::size_t npts,
346 std::size_t nsamples = 1, octave_idx_type stride = 1,
347 octave_idx_type dist = -1);
348 static int ifft (const FloatComplex *in, FloatComplex *out, std::size_t npts,
349 std::size_t nsamples = 1, octave_idx_type stride = 1,
350 octave_idx_type dist = -1);
351
352 static int fftNd (const float *, FloatComplex *, const int,
353 const dim_vector&);
354 static int fftNd (const FloatComplex *, FloatComplex *, const int,
355 const dim_vector&);
356 static int ifftNd (const FloatComplex *, FloatComplex *, const int,
357 const dim_vector&);
358};
359
360extern OCTAVE_API std::string fftw_version ();
361extern OCTAVE_API std::string fftwf_version ();
362
363OCTAVE_END_NAMESPACE(octave)
364
365#endif
Vector representing the dimensions (size) of an Array.
Definition dim-vector.h:90
static FftwMethod method(FftwMethod meth)
Definition oct-fftw.h:94
static FftwMethod method()
Definition oct-fftw.h:87
static void * create_plan(int dir, const int rank, const dim_vector &dims, octave_idx_type howmany, octave_idx_type stride, octave_idx_type dist, const Complex *in, Complex *out)
Definition oct-fftw.h:65
static int threads()
Definition oct-fftw.h:103
static void * create_plan(const int rank, const dim_vector &dims, octave_idx_type howmany, octave_idx_type stride, octave_idx_type dist, const double *in, Complex *out)
Definition oct-fftw.h:77
static void * create_plan(const int rank, const dim_vector &dims, octave_idx_type howmany, octave_idx_type stride, octave_idx_type dist, const float *in, FloatComplex *out)
Definition oct-fftw.h:217
static FftwMethod method(FftwMethod meth)
Definition oct-fftw.h:234
static int threads()
Definition oct-fftw.h:243
static FftwMethod method()
Definition oct-fftw.h:227
static void * create_plan(int dir, const int rank, const dim_vector &dims, octave_idx_type howmany, octave_idx_type stride, octave_idx_type dist, const FloatComplex *in, FloatComplex *out)
Definition oct-fftw.h:205
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
#define OCTAVE_API
Definition main.in.cc:55
std::complex< double > Complex
Definition oct-cmplx.h:33
std::complex< float > FloatComplex
Definition oct-cmplx.h:34
std::string fftwf_version()
Definition oct-fftw.cc:1144
std::string fftw_version()
Definition oct-fftw.cc:1134