00001 /* 00002 * This file is part of Healpix_cxx. 00003 * 00004 * Healpix_cxx 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 * Healpix_cxx 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 Healpix_cxx; if not, write to the Free Software 00016 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA 00017 * 00018 * For more information about HEALPix, see http://healpix.jpl.nasa.gov 00019 */ 00020 00021 /* 00022 * Healpix_cxx is being developed at the Max-Planck-Institut fuer Astrophysik 00023 * and financially supported by the Deutsches Zentrum fuer Luft- und Raumfahrt 00024 * (DLR). 00025 */ 00026 00027 /*! \file ls_fft.h 00028 * Interface for the LevelS FFT package. 00029 * 00030 * Copyright (C) 2004 Max-Planck-Society 00031 * \author Martin Reinecke 00032 */ 00033 00034 #ifndef PLANCK_LS_FFT_H 00035 #define PLANCK_LS_FFT_H 00036 00037 #ifdef __cplusplus 00038 extern "C" { 00039 #endif 00040 00041 typedef struct 00042 { 00043 double *work; 00044 int length; 00045 int bluestein; 00046 } complex_plan_i; 00047 00048 /*!\defgroup fftgroup FFT interface 00049 This package is intended to calculate one-dimensional real or complex FFTs 00050 with high accuracy and good efficiency even for lengths containing large 00051 prime factors. 00052 The code is written in C, but a Fortran wrapper exists as well. 00053 00054 Before any FFT is executed, a plan must be generated for it. Plan creation 00055 is designed to be fast, so that there is no significant overhead if the 00056 plan is only used once or a few times. 00057 00058 The main component of the code is a C port of Swarztrauber's FFTPACK 00059 (http://www.netlib.org/fftpack/), which was originally done by Pekka Janhunen 00060 and reformatted by Joerg Arndt. 00061 00062 I added a few digits to the floating-point constants to achieve higher 00063 precision, split the complex transforms into separate routines for forward 00064 and backward transforms (to increase performance a bit), and replaced 00065 the iterative twiddling factor calculation in radfg() and radbg() by an exact 00066 calculation of the factors. 00067 00068 Since FFTPACK becomes quite slow for FFT lengths with large prime factors 00069 (in the worst case of prime lengths it reaches O(n*n) complexity), I 00070 implemented Bluestein's algorithm, which computes a FFT of length n by 00071 several FFTs of length n2>=2*n and a convolution. Since n2 can be chosen 00072 to be highly composite, this algorithm is more efficient if n has large 00073 prime factors. The longer FFTs themselves are then computed using the FFTPACK 00074 routines. 00075 Bluestein's algorithm was implemented according to the description at 00076 <a href="http://en.wikipedia.org/wiki/Bluestein%27s_FFT_algorithm">Wikipedia</a>. 00077 00078 \b Thread-safety: 00079 All routines can be called concurrently; all information needed by ls_fft 00080 is stored in the plan variable. However, using the same plan variable on 00081 multiple threads simultaneously is not supported and will lead to data 00082 corruption. 00083 */ 00084 /*! \{ */ 00085 00086 /*! The opaque handle type for complex-FFT plans. */ 00087 typedef complex_plan_i * complex_plan; 00088 00089 /*! Returns a plan for a complex FFT with \a length elements. */ 00090 complex_plan make_complex_plan (int length); 00091 /*! Destroys a plan for a complex FFT. */ 00092 void kill_complex_plan (complex_plan plan); 00093 /*! Computes a complex forward FFT on \a data, using \a plan. 00094 \a Data has the form r0, i0, r1, i1, ..., r[length-1], i[length-1]. */ 00095 void complex_plan_forward (complex_plan plan, double *data); 00096 /*! Computes a complex backward FFT on \a data, using \a plan. 00097 \a Data has the form r0, i0, r1, i1, ..., r[length-1], i[length-1]. */ 00098 void complex_plan_backward (complex_plan plan, double *data); 00099 00100 typedef struct 00101 { 00102 double *work; 00103 int length; 00104 int bluestein; 00105 } real_plan_i; 00106 00107 /*! The opaque handle type for real-FFT plans. */ 00108 typedef real_plan_i * real_plan; 00109 00110 /*! Returns a plan for a real FFT with \a length elements. */ 00111 real_plan make_real_plan (int length); 00112 /*! Destroys a plan for a real FFT. */ 00113 void kill_real_plan (real_plan plan); 00114 /*! Computes a real forward FFT on \a data, using \a plan 00115 and assuming the FFTPACK storage scheme: 00116 - on entry, \a data has the form r0, r1, ..., r[length-1]; 00117 - on exit, it has the form r0, r1, i1, r2, i2, ... 00118 (a total of \a length values). */ 00119 void real_plan_forward_fftpack (real_plan plan, double *data); 00120 /*! Computes a real forward FFT on \a data, using \a plan 00121 and assuming the FFTPACK storage scheme: 00122 - on entry, \a data has the form r0, r1, i1, r2, i2, ... 00123 (a total of \a length values); 00124 - on exit, it has the form r0, r1, ..., r[length-1]. */ 00125 void real_plan_backward_fftpack (real_plan plan, double *data); 00126 /*! Computes a real forward FFT on \a data, using \a plan 00127 and assuming the FFTW halfcomplex storage scheme: 00128 - on entry, \a data has the form r0, r1, ..., r[length-1]; 00129 - on exit, it has the form r0, r1, r2, ..., i2, i1. */ 00130 void real_plan_forward_fftw (real_plan plan, double *data); 00131 /*! Computes a real backward FFT on \a data, using \a plan 00132 and assuming the FFTW halfcomplex storage scheme: 00133 - on entry, \a data has the form r0, r1, r2, ..., i2, i1. 00134 - on exit, it has the form r0, r1, ..., r[length-1]. */ 00135 void real_plan_backward_fftw (real_plan plan, double *data); 00136 /*! Computes a real forward FFT on \a data, using \a plan 00137 and assuming a full-complex storage scheme: 00138 - on entry, \a data has the form r0, [ignored], r1, [ignored], ..., 00139 r[length-1], [ignored]; 00140 - on exit, it has the form r0, i0, r1, i1, ..., r[length-1], i[length-1]. 00141 */ 00142 void real_plan_forward_c (real_plan plan, double *data); 00143 /*! Computes a real backward FFT on \a data, using \a plan 00144 and assuming a full-complex storage scheme: 00145 - on entry, \a data has the form r0, i0, r1, i1, ..., 00146 r[length-1], i[length-1]; 00147 - on exit, it has the form r0, 0, r1, 0, ..., r[length-1], 0. */ 00148 void real_plan_backward_c (real_plan plan, double *data); 00149 00150 /*! \} */ 00151 00152 #ifdef __cplusplus 00153 } 00154 #endif 00155 00156 #endif