changelog shortlog tags changeset files revisions annotate raw

src/DLD-FUNCTIONS/fftw.cc

changeset 10289: 4b124317dc38
parent:0ce82753dd72
author: John W. Eaton <jwe@octave.org>
date: Tue Feb 09 20:58:55 2010 -0500 (47 minutes ago)
permissions: -rw-r--r--
description: base_properties::set_children: account for hidden children
1/*
2
3Copyright (C) 2006, 2007, 2008, 2009 David Bateman
4
5This file is part of Octave.
6
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.
11
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
15for more details.
16
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/>.
20
21*/
22
23#ifdef HAVE_CONFIG_H
24#include <config.h>
25#endif
26
27#include <algorithm>
28
29#include "oct-fftw.h"
30
31#include "defun-dld.h"
32#include "error.h"
33#include "ov.h"
34
35DEFUN_DLD (fftw, args, ,
36 "-*- texinfo -*-\n\
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\
41\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\
48\n\
49@example\n\
50@var{wisdom} = fftw ('dwisdom')\n\
51@end example\n\
52\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\
56as follows\n\
57\n\
58@example\n\
59fftw ('dwisdom', @var{wisdom})\n\
60@end example \n\
61\n\
62If @var{wisdom} is an empty matrix, then the wisdom used is cleared.\n\
63\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\
68\n\
69@table @asis\n\
70@item 'estimate'\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\
76\n\
77@item 'measure'\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\
80\n\
81@item 'patient'\n\
82This is like 'measure', but a wider range of algorithms is considered.\n\
83\n\
84@item 'exhaustive'\n\
85This is like 'measure', but all possible algorithms that may be used to\n\
86treat the transform are considered.\n\
87\n\
88@item 'hybrid'\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\
92@end table\n\
93\n\
94The default method is 'estimate', and the method currently being used can\n\
95be probed with\n\
96\n\
97@example\n\
98@var{method} = fftw ('planner')\n\
99@end example\n\
100\n\
101and the method used can be set using\n\
102\n\
103@example\n\
104fftw ('planner', @var{method})\n\
105@end example\n\
106\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\
112@end deftypefn")
113{
114 octave_value retval;
115
116 int nargin = args.length();
117
118 if (nargin < 1 || nargin > 2)
119 {
120 print_usage ();
121 return retval;
122 }
123
124#if defined (HAVE_FFTW)
125 if (args(0).is_string ())
126 {
127 std::string arg0 = args(0).string_value ();
128
129 if (!error_state)
130 {
131 // Use STL function to convert to lower case
132 std::transform (arg0.begin (), arg0.end (), arg0.begin (), tolower);
133
134 if (nargin == 2)
135 {
136 std::string arg1 = args(1).string_value ();
137 if (!error_state)
138 {
139 if (arg0 == "planner")
140 {
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;
147
148 if (arg1 == "estimate")
149 {
150 meth = octave_fftw_planner::ESTIMATE;
151 methf = octave_float_fftw_planner::ESTIMATE;
152 }
153 else if (arg1 == "measure")
154 {
155 meth = octave_fftw_planner::MEASURE;
156 methf = octave_float_fftw_planner::MEASURE;
157 }
158 else if (arg1 == "patient")
159 {
160 meth = octave_fftw_planner::PATIENT;
161 methf = octave_float_fftw_planner::PATIENT;
162 }
163 else if (arg1 == "exhaustive")
164 {
165 meth = octave_fftw_planner::EXHAUSTIVE;
166 methf = octave_float_fftw_planner::EXHAUSTIVE;
167 }
168 else if (arg1 == "hybrid")
169 {
170 meth = octave_fftw_planner::HYBRID;
171 methf = octave_float_fftw_planner::HYBRID;
172 }
173 else
174 error ("unrecognized planner method");
175
176 if (!error_state)
177 {
178 meth = octave_fftw_planner::method (meth);
179 octave_float_fftw_planner::method (methf);
180
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");
189 else
190 retval = octave_value ("estimate");
191 }
192 }
193 else if (arg0 == "dwisdom")
194 {
195 char *str = fftw_export_wisdom_to_string ();
196
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");
201
202 if (!error_state)
203 retval = octave_value (std::string (str));
204
205 free (str);
206 }
207 else if (arg0 == "swisdom")
208 {
209 char *str = fftwf_export_wisdom_to_string ();
210
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");
215
216 if (!error_state)
217 retval = octave_value (std::string (str));
218
219 free (str);
220 }
221 else
222 error ("unrecognized argument");
223 }
224 }
225 else
226 {
227 if (arg0 == "planner")
228 {
229 octave_fftw_planner::FftwMethod meth =
230 octave_fftw_planner::method ();
231
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");
240 else
241 retval = octave_value ("estimate");
242 }
243 else if (arg0 == "dwisdom")
244 {
245 char *str = fftw_export_wisdom_to_string ();
246 retval = octave_value (std::string (str));
247 free (str);
248 }
249 else if (arg0 == "swisdom")
250 {
251 char *str = fftwf_export_wisdom_to_string ();
252 retval = octave_value (std::string (str));
253 free (str);
254 }
255 else
256 error ("unrecognized argument");
257 }
258 }
259 }
260#else
261
262 warning ("fftw: this copy of Octave was not configured to use the FFTW3 planner");
263
264#endif
265
266 return retval;
267}