GNU Octave  6.2.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-2021 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 
38 namespace octave
39 {
40  class
41  OCTAVE_API
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  ? 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  ? 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 () ? instance->do_method () : dummy;
98  }
99 
100  static FftwMethod method (FftwMethod _meth)
101  {
102  static FftwMethod dummy;
103 
104  return instance_ok () ? instance->do_method (_meth) : dummy;
105  }
106 
107  static void threads (int nt);
108 
109  static int threads (void)
110  {
111  return instance_ok () ? instance->nthreads : 0;
112  }
113 
114  private:
115 
117 
118  static void cleanup_instance (void) { delete instance; instance = nullptr; }
119 
120  void *
121  do_create_plan (int dir, const int rank, const dim_vector& dims,
122  octave_idx_type howmany, octave_idx_type stride,
123  octave_idx_type dist, const Complex *in,
124  Complex *out);
125 
126  void *
127  do_create_plan (const int rank, const dim_vector& dims,
128  octave_idx_type howmany, octave_idx_type stride,
129  octave_idx_type dist, const double *in, Complex *out);
130 
131  FftwMethod do_method (void);
132 
133  FftwMethod do_method (FftwMethod _meth);
134 
136 
137  // FIXME: perhaps this should be split into two classes?
138 
139  // Plan for fft and ifft of complex values
140  void *plan[2];
141 
142  // dist
144 
145  // stride
147 
148  // rank
149  int r[2];
150 
151  // howmany
153 
154  // dims
156 
157  bool simd_align[2];
158  bool inplace[2];
159 
160  // Plan for fft of real values
161  void *rplan;
162 
163  // dist
165 
166  // stride
168 
169  // rank
170  int rr;
171 
172  // howmany
174 
175  // dims
177 
179 
180  // number of threads. Always 1 unless compiled with multi-threading
181  // support.
182  int nthreads;
183  };
184 
185  class
186  OCTAVE_API
188  {
189  protected:
190 
191  float_fftw_planner (void);
192 
193  public:
194 
195  // No copying!
196 
198 
200  operator = (const float_fftw_planner&) = delete;
201 
202  ~float_fftw_planner (void);
203 
205  {
206  UNKNOWN = -1,
211  HYBRID
212  };
213 
214  static bool instance_ok (void);
215 
216  static void *
217  create_plan (int dir, const int rank, const dim_vector& dims,
218  octave_idx_type howmany, octave_idx_type stride,
219  octave_idx_type dist, const FloatComplex *in,
220  FloatComplex *out)
221  {
222  return instance_ok ()
223  ? instance->do_create_plan (dir, rank, dims, howmany, stride,
224  dist, in, out)
225  : nullptr;
226  }
227 
228  static void *
229  create_plan (const int rank, const dim_vector& dims,
230  octave_idx_type howmany, octave_idx_type stride,
231  octave_idx_type dist, const float *in, FloatComplex *out)
232  {
233  return instance_ok ()
234  ? instance->do_create_plan (rank, dims, howmany, stride, dist,
235  in, out)
236  : nullptr;
237  }
238 
239  static FftwMethod method (void)
240  {
241  static FftwMethod dummy;
242 
243  return instance_ok () ? instance->do_method () : dummy;
244  }
245 
246  static FftwMethod method (FftwMethod _meth)
247  {
248  static FftwMethod dummy;
249 
250  return instance_ok () ? instance->do_method (_meth) : dummy;
251  }
252 
253  static void threads (int nt);
254 
255  static int threads (void)
256  {
257  return instance_ok () ? instance->nthreads : 0;
258  }
259 
260  private:
261 
263 
264  static void cleanup_instance (void) { delete instance; instance = nullptr; }
265 
266  void *
267  do_create_plan (int dir, const int rank, const dim_vector& dims,
268  octave_idx_type howmany, octave_idx_type stride,
269  octave_idx_type dist, const FloatComplex *in,
270  FloatComplex *out);
271 
272  void *
273  do_create_plan (const int rank, const dim_vector& dims,
274  octave_idx_type howmany, octave_idx_type stride,
275  octave_idx_type dist, const float *in, FloatComplex *out);
276 
277  FftwMethod do_method (void);
278 
279  FftwMethod do_method (FftwMethod _meth);
280 
282 
283  // FIXME: perhaps this should be split into two classes?
284 
285  // Plan for fft and ifft of complex values
286  void *plan[2];
287 
288  // dist
290 
291  // stride
293 
294  // rank
295  int r[2];
296 
297  // howmany
299 
300  // dims
302 
303  bool simd_align[2];
304  bool inplace[2];
305 
306  // Plan for fft of real values
307  void *rplan;
308 
309  // dist
311 
312  // stride
314 
315  // rank
316  int rr;
317 
318  // howmany
320 
321  // dims
323 
325 
326  // number of threads. Always 1 unless compiled with multi-threading
327  // support.
328  int nthreads;
329  };
330 
331  class
332  OCTAVE_API
333  fftw
334  {
335  public:
336 
337  fftw (void) = delete;
338 
339  // No copying.
340 
341  fftw (const fftw&) = delete;
342 
343  fftw& operator = (const fftw&) = delete;
344 
345  static int fft (const double *in, Complex *out, size_t npts,
346  size_t nsamples = 1, octave_idx_type stride = 1,
347  octave_idx_type dist = -1);
348  static int fft (const Complex *in, Complex *out, size_t npts,
349  size_t nsamples = 1, octave_idx_type stride = 1,
350  octave_idx_type dist = -1);
351  static int ifft (const Complex *in, Complex *out, size_t npts,
352  size_t nsamples = 1, octave_idx_type stride = 1,
353  octave_idx_type dist = -1);
354 
355  static int fftNd (const double*, Complex*, const int, const dim_vector&);
356  static int fftNd (const Complex*, Complex*, const int,
357  const dim_vector&);
358  static int ifftNd (const Complex*, Complex*, const int,
359  const dim_vector&);
360 
361  static int fft (const float *in, FloatComplex *out, size_t npts,
362  size_t nsamples = 1, octave_idx_type stride = 1,
363  octave_idx_type dist = -1);
364  static int fft (const FloatComplex *in, FloatComplex *out, size_t npts,
365  size_t nsamples = 1, octave_idx_type stride = 1,
366  octave_idx_type dist = -1);
367  static int ifft (const FloatComplex *in, FloatComplex *out, size_t npts,
368  size_t nsamples = 1, octave_idx_type stride = 1,
369  octave_idx_type dist = -1);
370 
371  static int fftNd (const float*, FloatComplex*, const int, const dim_vector&);
372  static int fftNd (const FloatComplex*, FloatComplex*, const int,
373  const dim_vector&);
374  static int ifftNd (const FloatComplex*, FloatComplex*, const int,
375  const dim_vector&);
376  };
377 
378  extern OCTAVE_API std::string fftw_version (void);
379  extern OCTAVE_API std::string fftwf_version (void);
380 }
381 
382 #endif
Vector representing the dimensions (size) of an Array.
Definition: dim-vector.h:95
fftw_planner(const fftw_planner &)=delete
static fftw_planner * instance
Definition: oct-fftw.h:116
static FftwMethod method(void)
Definition: oct-fftw.h:93
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
octave_idx_type rh
Definition: oct-fftw.h:173
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
static int threads(void)
Definition: oct-fftw.h:109
static FftwMethod method(FftwMethod _meth)
Definition: oct-fftw.h:100
octave_idx_type rd
Definition: oct-fftw.h:164
octave_idx_type rs
Definition: oct-fftw.h:167
static void cleanup_instance(void)
Definition: oct-fftw.h:118
FftwMethod meth
Definition: oct-fftw.h:135
fftw(const fftw &)=delete
fftw(void)=delete
static void cleanup_instance(void)
Definition: oct-fftw.h:264
float_fftw_planner(const float_fftw_planner &)=delete
octave_idx_type rd
Definition: oct-fftw.h:310
static FftwMethod method(void)
Definition: oct-fftw.h:239
static FftwMethod method(FftwMethod _meth)
Definition: oct-fftw.h:246
static int threads(void)
Definition: oct-fftw.h:255
octave_idx_type rs
Definition: oct-fftw.h:313
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:229
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:217
octave_idx_type rh
Definition: oct-fftw.h:319
static float_fftw_planner * instance
Definition: oct-fftw.h:262
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE * d
octave_idx_type n
Definition: mx-inlines.cc:753
T * r
Definition: mx-inlines.cc:773
std::string fftw_version(void)
Definition: oct-fftw.cc:1123
std::string fftwf_version(void)
Definition: oct-fftw.cc:1133
std::complex< double > Complex
Definition: oct-cmplx.h:33
std::complex< float > FloatComplex
Definition: oct-cmplx.h:34