fftw.cc

Go to the documentation of this file.
00001 /*
00002 
00003 Copyright (C) 2006-2012 David Bateman
00004 
00005 This file is part of Octave.
00006 
00007 Octave is free software; you can redistribute it and/or modify it
00008 under the terms of the GNU General Public License as published by the
00009 Free Software Foundation; either version 3 of the License, or (at your
00010 option) any later version.
00011 
00012 Octave is distributed in the hope that it will be useful, but WITHOUT
00013 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00014 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00015 for more details.
00016 
00017 You should have received a copy of the GNU General Public License
00018 along with Octave; see the file COPYING.  If not, see
00019 <http://www.gnu.org/licenses/>.
00020 
00021 */
00022 
00023 #ifdef HAVE_CONFIG_H
00024 #include <config.h>
00025 #endif
00026 
00027 #include <algorithm>
00028 
00029 #include "oct-fftw.h"
00030 
00031 #include "defun-dld.h"
00032 #include "error.h"
00033 #include "ov.h"
00034 
00035 DEFUN_DLD (fftw, args, ,
00036   "-*- texinfo -*-\n\
00037 @deftypefn  {Loadable Function} {@var{method} =} fftw ('planner')\n\
00038 @deftypefnx {Loadable Function} {} fftw ('planner', @var{method})\n\
00039 @deftypefnx {Loadable Function} {@var{wisdom} =} fftw ('dwisdom')\n\
00040 @deftypefnx {Loadable Function} {} fftw ('dwisdom', @var{wisdom})\n\
00041 \n\
00042 Manage @sc{fftw} wisdom data.  Wisdom data can be used to significantly\n\
00043 accelerate the calculation of the FFTs, but implies an initial cost\n\
00044 in its calculation.  When the @sc{fftw} libraries are initialized, they read\n\
00045 a system wide wisdom file (typically in @file{/etc/fftw/wisdom}), allowing\n\
00046 wisdom to be shared between applications other than Octave.  Alternatively,\n\
00047 the @code{fftw} function can be used to import wisdom.  For example,\n\
00048 \n\
00049 @example\n\
00050 @var{wisdom} = fftw ('dwisdom')\n\
00051 @end example\n\
00052 \n\
00053 @noindent\n\
00054 will save the existing wisdom used by Octave to the string @var{wisdom}.\n\
00055 This string can then be saved to a file and restored using the @code{save}\n\
00056 and @code{load} commands respectively.  This existing wisdom can be\n\
00057 reimported as follows\n\
00058 \n\
00059 @example\n\
00060 fftw ('dwisdom', @var{wisdom})\n\
00061 @end example\n\
00062 \n\
00063 If @var{wisdom} is an empty matrix, then the wisdom used is cleared.\n\
00064 \n\
00065 During the calculation of Fourier transforms further wisdom is generated.\n\
00066 The fashion in which this wisdom is generated is also controlled by\n\
00067 the @code{fftw} function.  There are five different manners in which the\n\
00068 wisdom can be treated:\n\
00069 \n\
00070 @table @asis\n\
00071 @item 'estimate'\n\
00072 Specifies that no run-time measurement of the optimal means of\n\
00073 calculating a particular is performed, and a simple heuristic is used\n\
00074 to pick a (probably sub-optimal) plan.  The advantage of this method is\n\
00075 that there is little or no overhead in the generation of the plan, which\n\
00076 is appropriate for a Fourier transform that will be calculated once.\n\
00077 \n\
00078 @item 'measure'\n\
00079 In this case a range of algorithms to perform the transform is considered\n\
00080 and the best is selected based on their execution time.\n\
00081 \n\
00082 @item 'patient'\n\
00083 Similar to 'measure', but a wider range of algorithms is considered.\n\
00084 \n\
00085 @item 'exhaustive'\n\
00086 Like 'measure', but all possible algorithms that may be used to\n\
00087 treat the transform are considered.\n\
00088 \n\
00089 @item 'hybrid'\n\
00090 As run-time measurement of the algorithm can be expensive, this is a\n\
00091 compromise where 'measure' is used for transforms up to the size of 8192\n\
00092 and beyond that the 'estimate' method is used.\n\
00093 @end table\n\
00094 \n\
00095 The default method is 'estimate'.  The current method can\n\
00096 be queried with\n\
00097 \n\
00098 @example\n\
00099 @var{method} = fftw ('planner')\n\
00100 @end example\n\
00101 \n\
00102 @noindent\n\
00103 or set by using\n\
00104 \n\
00105 @example\n\
00106 fftw ('planner', @var{method})\n\
00107 @end example\n\
00108 \n\
00109 Note that calculated wisdom will be lost when restarting Octave.  However,\n\
00110 the wisdom data can be reloaded if it is saved to a file as described\n\
00111 above.  Saved wisdom files should not be used on different platforms since\n\
00112 they will not be efficient and the point of calculating the wisdom is lost.\n\
00113 @seealso{fft, ifft, fft2, ifft2, fftn, ifftn}\n\
00114 @end deftypefn")
00115 {
00116   octave_value retval;
00117 
00118   int nargin = args.length();
00119 
00120   if (nargin < 1 || nargin > 2)
00121     {
00122       print_usage ();
00123       return retval;
00124     }
00125 
00126 #if defined (HAVE_FFTW)
00127   if (args(0).is_string ())
00128     {
00129       std::string arg0 = args(0).string_value ();
00130 
00131       if (!error_state)
00132         {
00133           // Use STL function to convert to lower case
00134           std::transform (arg0.begin (), arg0.end (), arg0.begin (), tolower);
00135 
00136           if (nargin == 2)
00137             {
00138               std::string arg1 = args(1).string_value ();
00139               if (!error_state)
00140                 {
00141                   if (arg0 == "planner")
00142                     {
00143                       std::transform (arg1.begin (), arg1.end (),
00144                                       arg1.begin (), tolower);
00145                       octave_fftw_planner::FftwMethod meth
00146                         = octave_fftw_planner::UNKNOWN;
00147                       octave_float_fftw_planner::FftwMethod methf
00148                         = octave_float_fftw_planner::UNKNOWN;
00149 
00150                       if (arg1 == "estimate")
00151                         {
00152                           meth = octave_fftw_planner::ESTIMATE;
00153                           methf = octave_float_fftw_planner::ESTIMATE;
00154                         }
00155                       else if (arg1 == "measure")
00156                         {
00157                           meth = octave_fftw_planner::MEASURE;
00158                           methf = octave_float_fftw_planner::MEASURE;
00159                         }
00160                       else if (arg1 == "patient")
00161                         {
00162                           meth = octave_fftw_planner::PATIENT;
00163                           methf = octave_float_fftw_planner::PATIENT;
00164                         }
00165                       else if (arg1 == "exhaustive")
00166                         {
00167                           meth = octave_fftw_planner::EXHAUSTIVE;
00168                           methf = octave_float_fftw_planner::EXHAUSTIVE;
00169                         }
00170                       else if (arg1 == "hybrid")
00171                         {
00172                           meth = octave_fftw_planner::HYBRID;
00173                           methf = octave_float_fftw_planner::HYBRID;
00174                         }
00175                       else
00176                         error ("unrecognized planner METHOD");
00177 
00178                       if (!error_state)
00179                         {
00180                           meth = octave_fftw_planner::method (meth);
00181                           octave_float_fftw_planner::method (methf);
00182 
00183                           if (meth == octave_fftw_planner::MEASURE)
00184                             retval = octave_value ("measure");
00185                           else if (meth == octave_fftw_planner::PATIENT)
00186                             retval = octave_value ("patient");
00187                           else if (meth == octave_fftw_planner::EXHAUSTIVE)
00188                             retval = octave_value ("exhaustive");
00189                           else if (meth == octave_fftw_planner::HYBRID)
00190                             retval = octave_value ("hybrid");
00191                           else
00192                             retval = octave_value ("estimate");
00193                         }
00194                     }
00195                   else if (arg0 == "dwisdom")
00196                     {
00197                       char *str = fftw_export_wisdom_to_string ();
00198 
00199                       if (arg1.length() < 1)
00200                         fftw_forget_wisdom ();
00201                       else if (! fftw_import_wisdom_from_string (arg1.c_str()))
00202                         error ("could not import supplied WISDOM");
00203 
00204                       if (!error_state)
00205                         retval = octave_value (std::string (str));
00206 
00207                       free (str);
00208                     }
00209                   else if (arg0 == "swisdom")
00210                     {
00211                       char *str = fftwf_export_wisdom_to_string ();
00212 
00213                       if (arg1.length() < 1)
00214                         fftwf_forget_wisdom ();
00215                       else if (! fftwf_import_wisdom_from_string (arg1.c_str()))
00216                         error ("could not import supplied WISDOM");
00217 
00218                       if (!error_state)
00219                         retval = octave_value (std::string (str));
00220 
00221                       free (str);
00222                     }
00223                   else
00224                     error ("unrecognized argument");
00225                 }
00226             }
00227           else
00228             {
00229               if (arg0 == "planner")
00230                 {
00231                   octave_fftw_planner::FftwMethod meth =
00232                     octave_fftw_planner::method ();
00233 
00234                   if (meth == octave_fftw_planner::MEASURE)
00235                     retval = octave_value ("measure");
00236                   else if (meth == octave_fftw_planner::PATIENT)
00237                     retval = octave_value ("patient");
00238                   else if (meth == octave_fftw_planner::EXHAUSTIVE)
00239                     retval = octave_value ("exhaustive");
00240                   else if (meth == octave_fftw_planner::HYBRID)
00241                     retval = octave_value ("hybrid");
00242                   else
00243                     retval = octave_value ("estimate");
00244                 }
00245               else if (arg0 == "dwisdom")
00246                 {
00247                   char *str = fftw_export_wisdom_to_string ();
00248                   retval = octave_value (std::string (str));
00249                   free (str);
00250                 }
00251               else if (arg0 == "swisdom")
00252                 {
00253                   char *str = fftwf_export_wisdom_to_string ();
00254                   retval = octave_value (std::string (str));
00255                   free (str);
00256                 }
00257               else
00258                 error ("unrecognized argument");
00259             }
00260         }
00261     }
00262 #else
00263 
00264   warning ("fftw: this copy of Octave was not configured to use the FFTW3 planner");
00265 
00266 #endif
00267 
00268   return retval;
00269 }
00270 
00271 /*
00272 
00273 %!testif HAVE_FFTW
00274 %! def_method = fftw ("planner");
00275 %! unwind_protect
00276 %!   method = "estimate";
00277 %!   fftw ("planner", method);
00278 %!   assert (fftw ("planner"), method);
00279 %!   method = "measure";
00280 %!   fftw ("planner", method);
00281 %!   assert (fftw ("planner"), method);
00282 %!   method = "patient";
00283 %!   fftw ("planner", method);
00284 %!   assert (fftw ("planner"), method);
00285 %!   method = "exhaustive";
00286 %!   fftw ("planner", method);
00287 %!   assert (fftw ("planner"), method);
00288 %!   method = "hybrid";
00289 %!   fftw ("planner", method);
00290 %!   assert (fftw ("planner"), method);
00291 %! unwind_protect_cleanup
00292 %!   fftw ("planner", def_method);
00293 %! end_unwind_protect
00294 
00295 %!error <Invalid call to fftw> fftw ();
00296 %!error <Invalid call to fftw> fftw ("planner", "estimate", "measure");
00297 
00298  */
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines