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