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