ls_fft.h

Go to the documentation of this file.
00001 /*
00002  *  This file is part of libfftpack.
00003  *
00004  *  libfftpack is free software; you can redistribute it and/or modify
00005  *  it under the terms of the GNU General Public License as published by
00006  *  the Free Software Foundation; either version 2 of the License, or
00007  *  (at your option) any later version.
00008  *
00009  *  libfftpack is distributed in the hope that it will be useful,
00010  *  but WITHOUT ANY WARRANTY; without even the implied warranty of
00011  *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00012  *  GNU General Public License for more details.
00013  *
00014  *  You should have received a copy of the GNU General Public License
00015  *  along with libfftpack; if not, write to the Free Software
00016  *  Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00017  */
00018 
00019 /*
00020  *  libfftpack is being developed at the Max-Planck-Institut fuer Astrophysik
00021  *  and financially supported by the Deutsches Zentrum fuer Luft- und Raumfahrt
00022  *  (DLR).
00023  */
00024 
00025 /*! \file ls_fft.h
00026  *  Interface for the LevelS FFT package.
00027  *
00028  *  Copyright (C) 2004 Max-Planck-Society
00029  *  \author Martin Reinecke
00030  */
00031 
00032 #ifndef PLANCK_LS_FFT_H
00033 #define PLANCK_LS_FFT_H
00034 
00035 #include "c_utils.h"
00036 
00037 #ifdef __cplusplus
00038 extern "C" {
00039 #endif
00040 
00041 /*!\defgroup fftgroup FFT interface
00042 This package is intended to calculate one-dimensional real or complex FFTs
00043 with high accuracy and good efficiency even for lengths containing large
00044 prime factors.
00045 The code is written in C, but a Fortran wrapper exists as well.
00046 
00047 Before any FFT is executed, a plan must be generated for it. Plan creation
00048 is designed to be fast, so that there is no significant overhead if the
00049 plan is only used once or a few times.
00050 
00051 The main component of the code is based on Paul N. Swarztrauber's FFTPACK in the
00052 double precision incarnation by Hugh C. Pumphrey
00053 (http://www.netlib.org/fftpack/dp.tgz).
00054 
00055 I replaced the iterative sine and cosine calculations in radfg() and radbg()
00056 by an exact calculation, which slightly improves the transform accuracy for
00057 real FFTs with lengths containing large prime factors.
00058 
00059 Since FFTPACK becomes quite slow for FFT lengths with large prime factors
00060 (in the worst case of prime lengths it reaches \f$\mathcal{O}(n^2)\f$
00061 complexity), I implemented Bluestein's algorithm, which computes a FFT of length
00062 \f$n\f$ by several FFTs of length \f$n_2\ge 2n-1\f$ and a convolution. Since
00063 \f$n_2\f$ can be chosen to be highly composite, this algorithm is more efficient
00064 if \f$n\f$ has large prime factors. The longer FFTs themselves are then computed
00065 using the FFTPACK routines.
00066 Bluestein's algorithm was implemented according to the description on Wikipedia
00067 (<a href="http://en.wikipedia.org/wiki/Bluestein%27s_FFT_algorithm">
00068 http://en.wikipedia.org/wiki/Bluestein%27s_FFT_algorithm</a>).
00069 
00070 \b Thread-safety:
00071 All routines can be called concurrently; all information needed by
00072 <tt>ls_fft</tt> is stored in the plan variable. However, using the same plan
00073 variable on multiple threads simultaneously is not supported and will lead to
00074 data corruption.
00075 */
00076 /*! \{ */
00077 
00078 typedef struct
00079   {
00080   double *work;
00081   size_t length, worksize;
00082   int bluestein;
00083   } complex_plan_i;
00084 
00085 /*! The opaque handle type for complex-FFT plans. */
00086 typedef complex_plan_i * complex_plan;
00087 
00088 /*! Returns a plan for a complex FFT with \a length elements. */
00089 complex_plan make_complex_plan (size_t length);
00090 /*! Constructs a copy of \a plan. */
00091 complex_plan copy_complex_plan (complex_plan plan);
00092 /*! Destroys a plan for a complex FFT. */
00093 void kill_complex_plan (complex_plan plan);
00094 /*! Computes a complex forward FFT on \a data, using \a plan.
00095     \a Data has the form <tt>r0, i0, r1, i1, ...,
00096     r[length-1], i[length-1]</tt>. */
00097 void complex_plan_forward (complex_plan plan, double *data);
00098 /*! Computes a complex backward FFT on \a data, using \a plan.
00099     \a Data has the form <tt>r0, i0, r1, i1, ...,
00100     r[length-1], i[length-1]</tt>. */
00101 void complex_plan_backward (complex_plan plan, double *data);
00102 
00103 typedef struct
00104   {
00105   double *work;
00106   size_t length, worksize;
00107   int bluestein;
00108   } real_plan_i;
00109 
00110 /*! The opaque handle type for real-FFT plans. */
00111 typedef real_plan_i * real_plan;
00112 
00113 /*! Returns a plan for a real FFT with \a length elements. */
00114 real_plan make_real_plan (size_t length);
00115 /*! Constructs a copy of \a plan. */
00116 real_plan copy_real_plan (real_plan plan);
00117 /*! Destroys a plan for a real FFT. */
00118 void kill_real_plan (real_plan plan);
00119 /*! Computes a real forward FFT on \a data, using \a plan
00120     and assuming the FFTPACK storage scheme:
00121     - on entry, \a data has the form <tt>r0, r1, ..., r[length-1]</tt>;
00122     - on exit, it has the form <tt>r0, r1, i1, r2, i2, ...</tt>
00123       (a total of \a length values). */
00124 void real_plan_forward_fftpack (real_plan plan, double *data);
00125 /*! Computes a real backward FFT on \a data, using \a plan
00126     and assuming the FFTPACK storage scheme:
00127     - on entry, \a data has the form <tt>r0, r1, i1, r2, i2, ...</tt>
00128     (a total of \a length values);
00129     - on exit, it has the form <tt>r0, r1, ..., r[length-1]</tt>. */
00130 void real_plan_backward_fftpack (real_plan plan, double *data);
00131 /*! Computes a real forward FFT on \a data, using \a plan
00132     and assuming the FFTW halfcomplex storage scheme:
00133     - on entry, \a data has the form <tt>r0, r1, ..., r[length-1]</tt>;
00134     - on exit, it has the form <tt>r0, r1, r2, ..., i2, i1</tt>. */
00135 void real_plan_forward_fftw (real_plan plan, double *data);
00136 /*! Computes a real backward FFT on \a data, using \a plan
00137     and assuming the FFTW halfcomplex storage scheme:
00138     - on entry, \a data has the form <tt>r0, r1, r2, ..., i2, i1</tt>.
00139     - on exit, it has the form <tt>r0, r1, ..., r[length-1]</tt>. */
00140 void real_plan_backward_fftw (real_plan plan, double *data);
00141 /*! Computes a real forward FFT on \a data, using \a plan
00142     and assuming a full-complex storage scheme:
00143     - on entry, \a data has the form <tt>r0, [ignored], r1, [ignored], ...,
00144       r[length-1], [ignored]</tt>;
00145     - on exit, it has the form <tt>r0, i0, r1, i1, ...,
00146       r[length-1], i[length-1]</tt>. */
00147 void real_plan_forward_c (real_plan plan, double *data);
00148 /*! Computes a real backward FFT on \a data, using \a plan
00149     and assuming a full-complex storage scheme:
00150     - on entry, \a data has the form <tt>r0, i0, r1, i1, ...,
00151       r[length-1], i[length-1]</tt>;
00152     - on exit, it has the form <tt>r0, 0, r1, 0, ..., r[length-1], 0</tt>. */
00153 void real_plan_backward_c (real_plan plan, double *data);
00154 
00155 /*! \} */
00156 
00157 #ifdef __cplusplus
00158 }
00159 #endif
00160 
00161 #endif

Generated on Thu Oct 8 14:48:49 2015 for LevelS FFT library