3Copyright (C) 2006, 2007, 2008, 2009 David Bateman
5This file is part of Octave.
7Octave is free software; you can redistribute it and/or modify it
8under the terms of the GNU General Public License as published by the
9Free Software Foundation; either version 3 of the License, or (at your
10option) any later version.
12Octave is distributed in the hope that it will be useful, but WITHOUT
13ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
14FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
17You should have received a copy of the GNU General Public License
18along with Octave; see the file COPYING. If not, see
19<http://www.gnu.org/licenses/>.
35DEFUN_DLD (fftw, args, ,
37@deftypefn {Loadable Function} {@var{method} =} fftw ('planner')\n\
38@deftypefnx {Loadable Function} {} fftw ('planner', @var{method})\n\
39@deftypefnx {Loadable Function} {@var{wisdom} =} fftw ('dwisdom')\n\
40@deftypefnx {Loadable Function} {@var{wisdom} =} fftw ('dwisdom', @var{wisdom})\n\
42Manage @sc{fftw} wisdom data. Wisdom data can be used to significantly\n\
43accelerate the calculation of the FFTs but implies an initial cost\n\
44in its calculation. When the @sc{fftw} libraries are initialized, they read\n\
45a system wide wisdom file (typically in @file{/etc/fftw/wisdom}), allowing wisdom\n\
46to be shared between applications other than Octave. Alternatively, the\n\
47@code{fftw} function can be used to import wisdom. For example\n\
50@var{wisdom} = fftw ('dwisdom')\n\
53will save the existing wisdom used by Octave to the string @var{wisdom}.\n\
54This string can then be saved to a file and restored using the @code{save}\n\
55and @code{load} commands respectively. This existing wisdom can be reimported\n\
59fftw ('dwisdom', @var{wisdom})\n\
62If @var{wisdom} is an empty matrix, then the wisdom used is cleared.\n\
64During the calculation of Fourier transforms further wisdom is generated.\n\
65The fashion in which this wisdom is generated is equally controlled by\n\
66the @code{fftw} function. There are five different manners in which the\n\
67wisdom can be treated, these being\n\
71This specifies that no run-time measurement of the optimal means of\n\
72calculating a particular is performed, and a simple heuristic is used\n\
73to pick a (probably sub-optimal) plan. The advantage of this method is\n\
74that there is little or no overhead in the generation of the plan, which\n\
75is appropriate for a Fourier transform that will be calculated once.\n\
78In this case a range of algorithms to perform the transform is considered\n\
79and the best is selected based on their execution time.\n\
82This is like 'measure', but a wider range of algorithms is considered.\n\
85This is like 'measure', but all possible algorithms that may be used to\n\
86treat the transform are considered.\n\
89As run-time measurement of the algorithm can be expensive, this is a\n\
90compromise where 'measure' is used for transforms up to the size of 8192\n\
91and beyond that the 'estimate' method is used.\n\
94The default method is 'estimate', and the method currently being used can\n\
98@var{method} = fftw ('planner')\n\
101and the method used can be set using\n\
104fftw ('planner', @var{method})\n\
107Note that calculated wisdom will be lost when restarting Octave. However,\n\
108the wisdom data can be reloaded if it is saved to a file as described\n\
109above. Saved wisdom files should not be used on different platforms since\n\
110they will not be efficient and the point of calculating the wisdom is lost.\n\
111@seealso{fft, ifft, fft2, ifft2, fftn, ifftn}\n\
116 int nargin = args.length();
118 if (nargin < 1 || nargin > 2)
124#if defined (HAVE_FFTW)
125 if (args(0).is_string ())
127 std::string arg0 = args(0).string_value ();
131 // Use STL function to convert to lower case
132 std::transform (arg0.begin (), arg0.end (), arg0.begin (), tolower);
136 std::string arg1 = args(1).string_value ();
139 if (arg0 == "planner")
141 std::transform (arg1.begin (), arg1.end (),
142 arg1.begin (), tolower);
143 octave_fftw_planner::FftwMethod meth
144 = octave_fftw_planner::UNKNOWN;
145 octave_float_fftw_planner::FftwMethod methf
146 = octave_float_fftw_planner::UNKNOWN;
148 if (arg1 == "estimate")
150 meth = octave_fftw_planner::ESTIMATE;
151 methf = octave_float_fftw_planner::ESTIMATE;
153 else if (arg1 == "measure")
155 meth = octave_fftw_planner::MEASURE;
156 methf = octave_float_fftw_planner::MEASURE;
158 else if (arg1 == "patient")
160 meth = octave_fftw_planner::PATIENT;
161 methf = octave_float_fftw_planner::PATIENT;
163 else if (arg1 == "exhaustive")
165 meth = octave_fftw_planner::EXHAUSTIVE;
166 methf = octave_float_fftw_planner::EXHAUSTIVE;
168 else if (arg1 == "hybrid")
170 meth = octave_fftw_planner::HYBRID;
171 methf = octave_float_fftw_planner::HYBRID;
174 error ("unrecognized planner method");
178 meth = octave_fftw_planner::method (meth);
179 octave_float_fftw_planner::method (methf);
181 if (meth == octave_fftw_planner::MEASURE)
182 retval = octave_value ("measure");
183 else if (meth == octave_fftw_planner::PATIENT)
184 retval = octave_value ("patient");
185 else if (meth == octave_fftw_planner::EXHAUSTIVE)
186 retval = octave_value ("exhaustive");
187 else if (meth == octave_fftw_planner::HYBRID)
188 retval = octave_value ("hybrid");
190 retval = octave_value ("estimate");
193 else if (arg0 == "dwisdom")
195 char *str = fftw_export_wisdom_to_string ();
197 if (arg1.length() < 1)
198 fftw_forget_wisdom ();
199 else if (! fftw_import_wisdom_from_string (arg1.c_str()))
200 error ("could not import supplied wisdom");
203 retval = octave_value (std::string (str));
207 else if (arg0 == "swisdom")
209 char *str = fftwf_export_wisdom_to_string ();
211 if (arg1.length() < 1)
212 fftwf_forget_wisdom ();
213 else if (! fftwf_import_wisdom_from_string (arg1.c_str()))
214 error ("could not import supplied wisdom");
217 retval = octave_value (std::string (str));
222 error ("unrecognized argument");
227 if (arg0 == "planner")
229 octave_fftw_planner::FftwMethod meth =
230 octave_fftw_planner::method ();
232 if (meth == octave_fftw_planner::MEASURE)
233 retval = octave_value ("measure");
234 else if (meth == octave_fftw_planner::PATIENT)
235 retval = octave_value ("patient");
236 else if (meth == octave_fftw_planner::EXHAUSTIVE)
237 retval = octave_value ("exhaustive");
238 else if (meth == octave_fftw_planner::HYBRID)
239 retval = octave_value ("hybrid");
241 retval = octave_value ("estimate");
243 else if (arg0 == "dwisdom")
245 char *str = fftw_export_wisdom_to_string ();
246 retval = octave_value (std::string (str));
249 else if (arg0 == "swisdom")
251 char *str = fftwf_export_wisdom_to_string ();
252 retval = octave_value (std::string (str));
256 error ("unrecognized argument");
262 warning ("fftw: this copy of Octave was not configured to use the FFTW3 planner");