GNU Octave  6.2.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
fftw.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 2006-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 (HAVE_CONFIG_H)
27 # include "config.h"
28 #endif
29 
30 #include <algorithm>
31 #include <string>
32 
33 #if defined (HAVE_FFTW3_H)
34 # include <fftw3.h>
35 #endif
36 
37 #include "oct-fftw.h"
38 
39 #include "defun-dld.h"
40 #include "error.h"
41 #include "errwarn.h"
42 #include "ov.h"
43 
44 DEFUN_DLD (fftw, args, ,
45  doc: /* -*- texinfo -*-
46 @deftypefn {} {@var{method} =} fftw ("planner")
47 @deftypefnx {} {} fftw ("planner", @var{method})
48 @deftypefnx {} {@var{wisdom} =} fftw ("dwisdom")
49 @deftypefnx {} {} fftw ("dwisdom", @var{wisdom})
50 @deftypefnx {} {} fftw ("threads", @var{nthreads})
51 @deftypefnx {} {@var{nthreads} =} fftw ("threads")
52 
53 Manage @sc{fftw} wisdom data.
54 
55 Wisdom data can be used to significantly accelerate the calculation of the
56 FFTs, but implies an initial cost in its calculation. When the @sc{fftw}
57 libraries are initialized, they read a system wide wisdom file (typically in
58 @file{/etc/fftw/wisdom}), allowing wisdom to be shared between applications
59 other than Octave. Alternatively, the @code{fftw} function can be used to
60 import wisdom. For example,
61 
62 @example
63 @var{wisdom} = fftw ("dwisdom")
64 @end example
65 
66 @noindent
67 will save the existing wisdom used by Octave to the string @var{wisdom}.
68 This string can then be saved to a file and restored using the @code{save}
69 and @code{load} commands respectively. This existing wisdom can be
70 re-imported as follows
71 
72 @example
73 fftw ("dwisdom", @var{wisdom})
74 @end example
75 
76 If @var{wisdom} is an empty string, then the wisdom used is cleared.
77 
78 During the calculation of Fourier transforms further wisdom is generated.
79 The fashion in which this wisdom is generated is also controlled by
80 the @code{fftw} function. There are five different manners in which the
81 wisdom can be treated:
82 
83 @table @asis
84 @item @qcode{"estimate"}
85 Specifies that no run-time measurement of the optimal means of
86 calculating a particular is performed, and a simple heuristic is used
87 to pick a (probably sub-optimal) plan. The advantage of this method is
88 that there is little or no overhead in the generation of the plan, which
89 is appropriate for a Fourier transform that will be calculated once.
90 
91 @item @qcode{"measure"}
92 In this case a range of algorithms to perform the transform is considered
93 and the best is selected based on their execution time.
94 
95 @item @qcode{"patient"}
96 Similar to @qcode{"measure"}, but a wider range of algorithms is
97 considered.
98 
99 @item @qcode{"exhaustive"}
100 Like @qcode{"measure"}, but all possible algorithms that may be used to
101 treat the transform are considered.
102 
103 @item @qcode{"hybrid"}
104 As run-time measurement of the algorithm can be expensive, this is a
105 compromise where @qcode{"measure"} is used for transforms up to the size
106 of 8192 and beyond that the @qcode{"estimate"} method is used.
107 @end table
108 
109 The default method is @qcode{"estimate"}. The current method can
110 be queried with
111 
112 @example
113 @var{method} = fftw ("planner")
114 @end example
115 
116 @noindent
117 or set by using
118 
119 @example
120 fftw ("planner", @var{method})
121 @end example
122 
123 Note that calculated wisdom will be lost when restarting Octave. However,
124 the wisdom data can be reloaded if it is saved to a file as described
125 above. Saved wisdom files should not be used on different platforms since
126 they will not be efficient and the point of calculating the wisdom is lost.
127 
128 The number of threads used for computing the plans and executing the
129 transforms can be set with
130 
131 @example
132 fftw ("threads", @var{NTHREADS})
133 @end example
134 
135 Note that octave must be compiled with multi-threaded @sc{fftw} support for
136 this feature. The number of processors available to the current process is
137 used per default.
138 
139 @seealso{fft, ifft, fft2, ifft2, fftn, ifftn}
140 @end deftypefn */)
141 {
142 #if defined (HAVE_FFTW)
143 
144  int nargin = args.length ();
145 
146  if (nargin < 1 || nargin > 2)
147  print_usage ();
148 
150 
151  std::string arg0 = args(0).xstring_value ("fftw: first argument must be a string");
152 
153  if (arg0 == "planner")
154  {
155  if (nargin == 2) // planner setter
156  {
157  // Use STL function to convert to lower case
158  std::transform (arg0.begin (), arg0.end (), arg0.begin (), tolower);
159 
160  std::string arg1 = args(1).xstring_value ("fftw: METHOD must be a string");
161 
162  std::transform (arg1.begin (), arg1.end (), arg1.begin (), tolower);
167 
168  if (arg1 == "estimate")
169  {
172  }
173  else if (arg1 == "measure")
174  {
177  }
178  else if (arg1 == "patient")
179  {
182  }
183  else if (arg1 == "exhaustive")
184  {
187  }
188  else if (arg1 == "hybrid")
189  {
192  }
193  else
194  error ("fftw: unrecognized planner METHOD");
195 
196  meth = octave::fftw_planner::method (meth);
198 
199  if (meth == octave::fftw_planner::MEASURE)
200  retval = octave_value ("measure");
201  else if (meth == octave::fftw_planner::PATIENT)
202  retval = octave_value ("patient");
203  else if (meth == octave::fftw_planner::EXHAUSTIVE)
204  retval = octave_value ("exhaustive");
205  else if (meth == octave::fftw_planner::HYBRID)
206  retval = octave_value ("hybrid");
207  else
208  retval = octave_value ("estimate");
209  }
210  else //planner getter
211  {
214 
215  if (meth == octave::fftw_planner::MEASURE)
216  retval = octave_value ("measure");
217  else if (meth == octave::fftw_planner::PATIENT)
218  retval = octave_value ("patient");
219  else if (meth == octave::fftw_planner::EXHAUSTIVE)
220  retval = octave_value ("exhaustive");
221  else if (meth == octave::fftw_planner::HYBRID)
222  retval = octave_value ("hybrid");
223  else
224  retval = octave_value ("estimate");
225  }
226  }
227  else if (arg0 == "dwisdom")
228  {
229  if (nargin == 2) //dwisdom setter
230  {
231  // Use STL function to convert to lower case
232  std::transform (arg0.begin (), arg0.end (), arg0.begin (),
233  tolower);
234 
235  std::string arg1 = args(1).xstring_value ("fftw: WISDOM must be a string");
236 
237  char *str = fftw_export_wisdom_to_string ();
238  if (! str)
239  error ("fftw: could not get current FFTW wisdom");
240 
241  std::string wisdom_str (str);
242  free (str);
243 
244  if (arg1.length () < 1)
245  fftw_forget_wisdom ();
246  else if (! fftw_import_wisdom_from_string (arg1.c_str ()))
247  error ("fftw: could not import supplied WISDOM");
248 
249  retval = octave_value (wisdom_str);
250  }
251  else //dwisdom getter
252  {
253  char *str = fftw_export_wisdom_to_string ();
254  if (! str)
255  error ("fftw: could not get current FFTW wisdom");
256 
257  std::string wisdom_str (str);
258  free (str);
259  retval = octave_value (wisdom_str);
260  }
261  }
262  else if (arg0 == "swisdom")
263  {
264  //swisdom uses fftwf_ functions (float), dwisdom fftw_ (real)
265  if (nargin == 2) //swisdom setter
266  {
267  // Use STL function to convert to lower case
268  std::transform (arg0.begin (), arg0.end (), arg0.begin (),
269  tolower);
270 
271  std::string arg1 = args(1).xstring_value ("fftw: WISDOM must be a string");
272 
273  char *str = fftwf_export_wisdom_to_string ();
274  if (! str)
275  error ("fftw: could not get current FFTW wisdom");
276 
277  std::string wisdom_str (str);
278  free (str);
279 
280  if (arg1.length () < 1)
281  fftwf_forget_wisdom ();
282  else if (! fftwf_import_wisdom_from_string (arg1.c_str ()))
283  error ("fftw: could not import supplied WISDOM");
284 
285  retval = octave_value (wisdom_str);
286  }
287  else //swisdom getter
288  {
289  char *str = fftwf_export_wisdom_to_string ();
290  if (! str)
291  error ("fftw: could not get current FFTW wisdom");
292 
293  std::string wisdom_str (str);
294  free (str);
295  retval = octave_value (wisdom_str);
296  }
297  }
298  else if (arg0 == "threads")
299  {
300  if (nargin == 2) //threads setter
301  {
302  if (! args(1).is_real_scalar ())
303  error ("fftw: setting threads needs one integer argument");
304 
305  int nthreads = args(1).int_value();
306  if (nthreads < 1)
307  error ("fftw: number of threads must be >=1");
308 
309 #if defined (HAVE_FFTW3_THREADS)
311 #else
312  err_disabled_feature ("fftw", "multithreaded FFTW");
313 #endif
314 #if defined (HAVE_FFTW3F_THREADS)
316 #else
317  err_disabled_feature ("fftw", "multithreaded FFTW");
318 #endif
319  }
320  else //threads getter
321 #if defined (HAVE_FFTW3_THREADS)
323 #else
324  retval = 1;
325 #endif
326  }
327  else
328  error ("fftw: unrecognized argument");
329 
330  return retval;
331 
332 #else
333 
334  octave_unused_parameter (args);
335 
336  err_disabled_feature ("fftw", "the FFTW3 planner");
337 
338 #endif
339 }
340 
341 /*
342 %!testif HAVE_FFTW
343 %! def_method = fftw ("planner");
344 %! unwind_protect
345 %! method = "estimate";
346 %! fftw ("planner", method);
347 %! assert (fftw ("planner"), method);
348 %! method = "measure";
349 %! fftw ("planner", method);
350 %! assert (fftw ("planner"), method);
351 %! method = "patient";
352 %! fftw ("planner", method);
353 %! assert (fftw ("planner"), method);
354 %! method = "exhaustive";
355 %! fftw ("planner", method);
356 %! assert (fftw ("planner"), method);
357 %! method = "hybrid";
358 %! fftw ("planner", method);
359 %! assert (fftw ("planner"), method);
360 %! unwind_protect_cleanup
361 %! fftw ("planner", def_method);
362 %! end_unwind_protect
363 
364 %!testif HAVE_FFTW
365 %! def_dwisdom = fftw ("dwisdom");
366 %! def_swisdom = fftw ("swisdom");
367 %! unwind_protect
368 %! wisdom = fftw ("dwisdom");
369 %! assert (ischar (wisdom));
370 %! fftw ("dwisdom", wisdom);
371 %! assert (fftw ("dwisdom"), wisdom);
372 %! wisdom = fftw ("swisdom");
373 %! assert (ischar (wisdom));
374 %! fftw ("swisdom", wisdom);
375 %! assert (fftw ("swisdom"), wisdom);
376 %! unwind_protect_cleanup
377 %! fftw ("dwisdom", def_dwisdom);
378 %! fftw ("swisdom", def_swisdom);
379 %! end_unwind_protect
380 
381 %!testif HAVE_FFTW3_THREADS
382 %! n = fftw ("threads");
383 %! unwind_protect
384 %! fftw ("threads", 3);
385 %! assert (fftw ("threads"), 3);
386 %! unwind_protect_cleanup
387 %! fftw ("threads", n);
388 %! end_unwind_protect
389 
390 %!error <Invalid call to fftw|was unavailable or disabled> fftw ()
391 %!error <Invalid call to fftw|was unavailable or disabled> fftw ("planner", "estimate", "measure")
392 %!error fftw (3)
393 %!error fftw ("invalid")
394 %!error fftw ("planner", "invalid")
395 %!error fftw ("planner", 2)
396 %!error fftw ("dwisdom", "invalid")
397 %!error fftw ("swisdom", "invalid")
398 %!error fftw ("threads", "invalid")
399 %!error fftw ("threads", -3)
400  */
static FftwMethod method(void)
Definition: oct-fftw.h:93
static int threads(void)
Definition: oct-fftw.h:109
static FftwMethod method(void)
Definition: oct-fftw.h:239
static int threads(void)
Definition: oct-fftw.h:255
#define DEFUN_DLD(name, args_name, nargout_name, doc)
Macro to define an at run time dynamically loadable builtin function.
Definition: defun-dld.h:61
OCTINTERP_API void print_usage(void)
Definition: defun.cc:53
void error(const char *fmt,...)
Definition: error.cc:968
void err_disabled_feature(const std::string &fcn, const std::string &feature, const std::string &pkg)
Definition: errwarn.cc:53
ColumnVector transform(const Matrix &m, double x, double y, double z)
Definition: graphics.cc:5814
void free(void *)
return octave_value(v1.char_array_value() . concat(v2.char_array_value(), ra_idx),((a1.is_sq_string()||a2.is_sq_string()) ? '\'' :'"'))
octave_value::octave_value(const Array< char > &chm, char type) return retval
Definition: ov.cc:811