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