GNU Octave 7.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
oct-fftw.h
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 2001-2022 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
38namespace octave
39{
40 class
43 {
44 protected:
45
46 fftw_planner (void);
47
48 public:
49
50 // No copying!
51
52 fftw_planner (const fftw_planner&) = delete;
53
54 fftw_planner& operator = (const fftw_planner&) = delete;
55
56 ~fftw_planner (void);
57
59 {
60 UNKNOWN = -1,
65 HYBRID
66 };
67
68 static bool instance_ok (void);
69
70 static void *
71 create_plan (int dir, const int rank, const dim_vector& dims,
72 octave_idx_type howmany, octave_idx_type stride,
73 octave_idx_type dist, const Complex *in,
74 Complex *out)
75 {
76 return instance_ok ()
77 ? s_instance->do_create_plan (dir, rank, dims, howmany, stride,
78 dist, in, out)
79 : nullptr;
80 }
81
82 static void *
83 create_plan (const int rank, const dim_vector& dims,
84 octave_idx_type howmany, octave_idx_type stride,
85 octave_idx_type dist, const double *in, Complex *out)
86 {
87 return instance_ok ()
88 ? s_instance->do_create_plan (rank, dims, howmany, stride, dist,
89 in, out)
90 : nullptr;
91 }
92
93 static FftwMethod method (void)
94 {
95 static FftwMethod dummy;
96
97 return instance_ok () ? s_instance->do_method () : dummy;
98 }
99
101 {
102 static FftwMethod dummy;
103
104 return instance_ok () ? s_instance->do_method (meth) : dummy;
105 }
106
107 static void threads (int nt);
108
109 static int threads (void)
110 {
111 return instance_ok () ? s_instance->m_nthreads : 0;
112 }
113
114 private:
115
117
118 static void cleanup_instance (void)
119 { delete s_instance; s_instance = nullptr; }
120
121 void *
122 do_create_plan (int dir, const int rank, const dim_vector& dims,
123 octave_idx_type howmany, octave_idx_type stride,
124 octave_idx_type dist, const Complex *in,
125 Complex *out);
126
127 void *
128 do_create_plan (const int rank, const dim_vector& dims,
129 octave_idx_type howmany, octave_idx_type stride,
130 octave_idx_type dist, const double *in, Complex *out);
131
132 FftwMethod do_method (void);
133
134 FftwMethod do_method (FftwMethod meth);
135
137
138 // FIXME: perhaps this should be split into two classes?
139
140 // Plan for fft and ifft of complex values
141 void *m_plan[2];
142
143 // dist
145
146 // stride
148
149 // rank
150 int m_r[2];
151
152 // howmany
154
155 // dims
157
158 bool m_simd_align[2];
159 bool m_inplace[2];
160
161 // Plan for fft of real values
162 void *m_rplan;
163
164 // dist
166
167 // stride
169
170 // rank
171 int m_rr;
172
173 // howmany
175
176 // dims
178
180
181 // number of threads. Always 1 unless compiled with multi-threading
182 // support.
184 };
185
186 class
189 {
190 protected:
191
192 float_fftw_planner (void);
193
194 public:
195
196 // No copying!
197
199
201 operator = (const float_fftw_planner&) = delete;
202
203 ~float_fftw_planner (void);
204
206 {
207 UNKNOWN = -1,
212 HYBRID
213 };
214
215 static bool instance_ok (void);
216
217 static void *
218 create_plan (int dir, const int rank, const dim_vector& dims,
219 octave_idx_type howmany, octave_idx_type stride,
220 octave_idx_type dist, const FloatComplex *in,
221 FloatComplex *out)
222 {
223 return instance_ok ()
224 ? s_instance->do_create_plan (dir, rank, dims, howmany, stride,
225 dist, in, out)
226 : nullptr;
227 }
228
229 static void *
230 create_plan (const int rank, const dim_vector& dims,
231 octave_idx_type howmany, octave_idx_type stride,
232 octave_idx_type dist, const float *in, FloatComplex *out)
233 {
234 return instance_ok ()
235 ? s_instance->do_create_plan (rank, dims, howmany, stride, dist,
236 in, out)
237 : nullptr;
238 }
239
240 static FftwMethod method (void)
241 {
242 static FftwMethod dummy;
243
244 return instance_ok () ? s_instance->do_method () : dummy;
245 }
246
248 {
249 static FftwMethod dummy;
250
251 return instance_ok () ? s_instance->do_method (meth) : dummy;
252 }
253
254 static void threads (int nt);
255
256 static int threads (void)
257 {
258 return instance_ok () ? s_instance->m_nthreads : 0;
259 }
260
261 private:
262
264
265 static void cleanup_instance (void)
266 { delete s_instance; s_instance = nullptr; }
267
268 void *
269 do_create_plan (int dir, const int rank, const dim_vector& dims,
270 octave_idx_type howmany, octave_idx_type stride,
271 octave_idx_type dist, const FloatComplex *in,
272 FloatComplex *out);
273
274 void *
275 do_create_plan (const int rank, const dim_vector& dims,
276 octave_idx_type howmany, octave_idx_type stride,
277 octave_idx_type dist, const float *in, FloatComplex *out);
278
279 FftwMethod do_method (void);
280
281 FftwMethod do_method (FftwMethod meth);
282
284
285 // FIXME: perhaps this should be split into two classes?
286
287 // Plan for fft and ifft of complex values
288 void *m_plan[2];
289
290 // dist
292
293 // stride
295
296 // rank
297 int m_r[2];
298
299 // howmany
301
302 // dims
304
305 bool m_simd_align[2];
306 bool m_inplace[2];
307
308 // Plan for fft of real values
309 void *m_rplan;
310
311 // dist
313
314 // stride
316
317 // rank
318 int m_rr;
319
320 // howmany
322
323 // dims
325
327
328 // number of threads. Always 1 unless compiled with multi-threading
329 // support.
331 };
332
333 class
335 fftw
336 {
337 public:
338
339 fftw (void) = delete;
340
341 // No copying.
342
343 fftw (const fftw&) = delete;
344
345 fftw& operator = (const fftw&) = delete;
346
347 static int fft (const double *in, Complex *out, std::size_t npts,
348 std::size_t nsamples = 1, octave_idx_type stride = 1,
349 octave_idx_type dist = -1);
350 static int fft (const Complex *in, Complex *out, std::size_t npts,
351 std::size_t nsamples = 1, octave_idx_type stride = 1,
352 octave_idx_type dist = -1);
353 static int ifft (const Complex *in, Complex *out, std::size_t npts,
354 std::size_t nsamples = 1, octave_idx_type stride = 1,
355 octave_idx_type dist = -1);
356
357 static int fftNd (const double *, Complex *, const int, const dim_vector&);
358 static int fftNd (const Complex *, Complex *, const int,
359 const dim_vector&);
360 static int ifftNd (const Complex *, Complex *, const int,
361 const dim_vector&);
362
363 static int fft (const float *in, FloatComplex *out, std::size_t npts,
364 std::size_t nsamples = 1, octave_idx_type stride = 1,
365 octave_idx_type dist = -1);
366 static int fft (const FloatComplex *in, FloatComplex *out, std::size_t npts,
367 std::size_t nsamples = 1, octave_idx_type stride = 1,
368 octave_idx_type dist = -1);
369 static int ifft (const FloatComplex *in, FloatComplex *out, std::size_t npts,
370 std::size_t nsamples = 1, octave_idx_type stride = 1,
371 octave_idx_type dist = -1);
372
373 static int fftNd (const float *, FloatComplex *, const int,
374 const dim_vector&);
375 static int fftNd (const FloatComplex *, FloatComplex *, const int,
376 const dim_vector&);
377 static int ifftNd (const FloatComplex *, FloatComplex *, const int,
378 const dim_vector&);
379 };
380
381 extern OCTAVE_API std::string fftw_version (void);
382 extern OCTAVE_API std::string fftwf_version (void);
383}
384
385#endif
Vector representing the dimensions (size) of an Array.
Definition: dim-vector.h:94
fftw_planner(const fftw_planner &)=delete
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:83
static FftwMethod method(FftwMethod meth)
Definition: oct-fftw.h:100
octave_idx_type m_rs
Definition: oct-fftw.h:168
FftwMethod m_meth
Definition: oct-fftw.h:136
static fftw_planner * s_instance
Definition: oct-fftw.h:116
static FftwMethod method(void)
Definition: oct-fftw.h:93
static int threads(void)
Definition: oct-fftw.h:109
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:71
octave_idx_type m_rd
Definition: oct-fftw.h:165
dim_vector m_rn
Definition: oct-fftw.h:177
octave_idx_type m_rh
Definition: oct-fftw.h:174
static void cleanup_instance(void)
Definition: oct-fftw.h:118
fftw(const fftw &)=delete
fftw(void)=delete
static void cleanup_instance(void)
Definition: oct-fftw.h:265
octave_idx_type m_rd
Definition: oct-fftw.h:312
float_fftw_planner(const float_fftw_planner &)=delete
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:218
static FftwMethod method(void)
Definition: oct-fftw.h:240
static float_fftw_planner * s_instance
Definition: oct-fftw.h:263
octave_idx_type m_rs
Definition: oct-fftw.h:315
static FftwMethod method(FftwMethod meth)
Definition: oct-fftw.h:247
static int threads(void)
Definition: oct-fftw.h:256
octave_idx_type m_rh
Definition: oct-fftw.h:321
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:230
#define OCTAVE_API
Definition: main.in.cc:55
std::string fftw_version(void)
Definition: oct-fftw.cc:1133
std::string fftwf_version(void)
Definition: oct-fftw.cc:1143
std::complex< double > Complex
Definition: oct-cmplx.h:33
std::complex< float > FloatComplex
Definition: oct-cmplx.h:34