00001 /* 00002 * Copyright (C) 2007 Duncan Brown, Jolien Creighton, Kipp Cannon 00003 * 00004 * This program 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 * This program 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 with program; see the file COPYING. If not, write to the 00016 * Free Software Foundation, Inc., 59 Temple Place, Suite 330, Boston, 00017 * MA 02111-1307 USA 00018 */ 00019 00020 /**** <lalVerbatim file="RealFFTHV"> 00021 * $Id: RealFFT.h,v 1.18 2007/06/08 14:41:44 bema Exp $ 00022 **** </lalVerbatim> */ 00023 00024 /**** <lalLaTeX> 00025 * 00026 * \section{Header \texttt{RealFFT.h}} 00027 * \label{s:RealFFT.h} 00028 * 00029 * Performs real-to-complex and complex-to-real FFTs. 00030 * 00031 * \subsection*{Synopsis} 00032 * \begin{verbatim} 00033 * #include <lal/RealFFT.h> 00034 * \end{verbatim} 00035 * 00036 * \noindent Perform real-to-complex and complex-to-real fast Fourier 00037 * transforms of vectors, and sequences of vectors using the package 00038 * FFTW~\cite{fj:1998}. 00039 * 00040 **** </lalLaTeX> */ 00041 00042 #ifndef _REALFFT_H 00043 #define _REALFFT_H 00044 00045 #include <lal/LALDatatypes.h> 00046 00047 #ifdef __cplusplus 00048 extern "C" { 00049 #pragma } 00050 #endif 00051 00052 NRCSID( REALFFTH, "$Id: RealFFT.h,v 1.18 2007/06/08 14:41:44 bema Exp $" ); 00053 00054 00055 /**** <lalLaTeX> 00056 * \subsection*{Error conditions} 00057 **** </lalLaTeX> */ 00058 /**** <lalErrTable> */ 00059 00060 #define REALFFTH_ENULL 1 00061 #define REALFFTH_ENNUL 2 00062 #define REALFFTH_ESIZE 4 00063 #define REALFFTH_ESZMM 8 00064 #define REALFFTH_ESLEN 16 00065 #define REALFFTH_ESAME 32 00066 #define REALFFTH_ESIGN 64 00067 #define REALFFTH_EDATA 128 00068 #define REALFFTH_EALOC 256 00069 #define REALFFTH_EFFTW 512 00070 #define REALFFTH_ESNGL 1024 00071 #define REALFFTH_EINTL 2048 00072 00073 #define REALFFTH_MSGENULL "Null pointer" 00074 #define REALFFTH_MSGENNUL "Non-null pointer" 00075 #define REALFFTH_MSGESIZE "Invalid input size" 00076 #define REALFFTH_MSGESZMM "Size mismatch" 00077 #define REALFFTH_MSGESLEN "Invalid/mismatched sequence lengths" 00078 #define REALFFTH_MSGESAME "Input/Output data vectors are the same" 00079 #define REALFFTH_MSGESIGN "Incorrect plan sign" 00080 #define REALFFTH_MSGEDATA "Bad input data: DC/Nyquist should be real" 00081 #define REALFFTH_MSGEALOC "Memory allocation failed" 00082 #define REALFFTH_MSGEFFTW "Error in FFTW" 00083 #define REALFFTH_MSGESNGL "FFTW library is not single-precision" 00084 #define REALFFTH_MSGEINTL "Error in Intel FFT library" 00085 00086 00087 /**** </lalErrTable> */ 00088 /**** <lalLaTeX> 00089 * 00090 * \subsection*{Structures} 00091 * 00092 **** </lalLaTeX> */ 00093 /**** <lalVerbatim> */ 00094 typedef struct tagREAL4FFTPlan REAL4FFTPlan; 00095 typedef struct tagREAL8FFTPlan REAL8FFTPlan; 00096 #define tagRealFFTPlan tagREAL4FFTPlan 00097 #define RealFFTPlan REAL4FFTPlan 00098 /**** </lalVerbatim> */ 00099 /**** <lalLaTeX> 00100 * 00101 * This structure contains the parameters necessary for performing an FFT of a 00102 * given size and direction. The contents should not be manually adjusted. 00103 * 00104 * \newpage\input{RealFFTC} 00105 * 00106 * \newpage\subsection{XLAL Functions} 00107 * 00108 * \subsubsection*{Synopsis} 00109 * \begin{verbatim} 00110 * #include <lal/RealFFT.h> 00111 * 00112 * REAL4FFTPlan * XLALCreateREAL4FFTPlan( UINT4 size, int fwdflg, int measurelvl ); 00113 * REAL4FFTPlan * XLALCreateForwardREAL4FFTPlan( UINT4 size, int measurelvl ); 00114 * REAL4FFTPlan * XLALCreateReverseREAL4FFTPlan( UINT4 size, int measurelvl ); 00115 * void XLALDestroyREAL4FFTPlan( REAL4FFTPlan *plan ); 00116 * 00117 * int XLALREAL4ForwardFFT( COMPLEX8Vector *output, REAL4Vector *input, 00118 * REAL4FFTPlan *plan ); 00119 * int XLALREAL4ReverseFFT( REAL4Vector *output, COMPLEX8Vector *input, 00120 * REAL4FFTPlan *plan ); 00121 * int XLALREAL4VectorFFT( REAL4Vector *output, REAL4Vector *input, 00122 * REAL4FFTPlan *plan ); 00123 * int XLALREAL4PowerSpectrum( REAL4Vector *spec, REAL4Vector *data, 00124 * REAL4FFTPlan *plan ); 00125 * 00126 * REAL8FFTPlan * XLALCreateREAL8FFTPlan( UINT4 size, int fwdflg, int measurelvl ); 00127 * REAL8FFTPlan * XLALCreateForwardREAL8FFTPlan( UINT4 size, int measurelvl ); 00128 * REAL8FFTPlan * XLALCreateReverseREAL8FFTPlan( UINT4 size, int measurelvl ); 00129 * void XLALDestroyREAL8FFTPlan( REAL8FFTPlan *plan ); 00130 * 00131 * int XLALREAL8ForwardFFT( COMPLEX16Vector *output, REAL8Vector *input, 00132 * REAL8FFTPlan *plan ); 00133 * int XLALREAL8ReverseFFT( REAL8Vector *output, COMPLEX16Vector *input, 00134 * REAL8FFTPlan *plan ); 00135 * int XLALREAL8VectorFFT( REAL8Vector *output, REAL8Vector *input, 00136 * REAL8FFTPlan *plan ); 00137 * int XLALREAL8PowerSpectrum( REAL8Vector *spec, REAL8Vector *data, 00138 * REAL8FFTPlan *plan ); 00139 * \end{verbatim} 00140 * \idx{XLALCreateREAL4FFTPlan} 00141 * \idx{XLALCreateForwardREAL4FFTPlan} 00142 * \idx{XLALCreateReverseREAL4FFTPlan} 00143 * \idx{XLALDestroyREAL4FFTPlan} 00144 * \idx{XLALREAL4ForwardFFT} 00145 * \idx{XLALREAL4ReverseFFT} 00146 * \idx{XLALREAL4VectorFFT} 00147 * \idx{XLALREAL4PowerSpectrum} 00148 * \idx{XLALCreateREAL8FFTPlan} 00149 * \idx{XLALCreateForwardREAL8FFTPlan} 00150 * \idx{XLALCreateReverseREAL8FFTPlan} 00151 * \idx{XLALDestroyREAL8FFTPlan} 00152 * \idx{XLALREAL8ForwardFFT} 00153 * \idx{XLALREAL8ReverseFFT} 00154 * \idx{XLALREAL8VectorFFT} 00155 * \idx{XLALREAL8PowerSpectrum} 00156 * 00157 * \subsubsection*{Description} 00158 * 00159 * The \verb+REAL4+ routines are described below. These use single-precision 00160 * FFTs, i.e., they convert \verb+REAL4Vector+s into \verb+COMPLEX8Vector+s 00161 * and vice-versa. The \verb+REAL8+ versions of the routines are the same 00162 * but they are double-precision versions, i.e., they convert 00163 * \verb+REAL8Vector+s into \verb+COMPLEX16Vector+s. 00164 * 00165 * The routine \verb+XLALCreateREAL4FFTPlan+ creates a \verb+REAL4FFTPlan+ 00166 * structure to perform FFTs of vectors of length \verb+size+. If 00167 * \verb+fwdflg+ is non-zero then the plan is created to perform forward 00168 * (real-to-complex) FFTs with a negative exponential sign. Otherwise 00169 * the plan is created to perform reverse (complex-to-real) FFTs with a 00170 * positive exponential sign. The value of \verb+measurelvl+ determines 00171 * how much optimization of the plan FFTW will do with the most optimization 00172 * taking the most amount of time. Reasonable values for \verb+measurelvl+ 00173 * would be 0 for the fasted plan creation (FFTW does not measure the speed 00174 * of any transform with this level but rather estimates which plan will 00175 * be the fastet) or 1 to measure a few likely plans to determine the fastest. 00176 * 00177 * \verb+XLALCreateForwardREAL4FFTPlan+ is equivalent to 00178 * \verb+XLALCreateREAL4FFTPlan+ with \verb+fwdflg+ set to 1. 00179 * \verb+XLALCreateReverseREAL4FFTPlan+ is equivalent to 00180 * \verb+XLALCreateREAL4FFTPlan+ with \verb+fwdflg+ set to 0. 00181 * 00182 * \verb+XLALDestroyREAL4FFTPlan+ is used to destroy the plan, freeing all 00183 * memory that was allocated in the structure as well as the structure 00184 * itself. It can be used on either forward or reverse plans. 00185 * 00186 * \verb+XLALREAL4ForwardFFT+ and 00187 * \verb+XLALREAL4ReverseFFT+ perform forward (real to complex) and 00188 * reverse (complex to real) transforms respectively. The plan supplied 00189 * to these routines must be correctly generated for the direction of the 00190 * transform. I.e., \verb+XLALREAL4ForwardFFT+ cannot be supplied with 00191 * a plan generated by \verb+XLALCreateReverseREAL4FFTPlan+. 00192 * 00193 * \verb+XLALREAL4VectorFFT+ is a low-level routine that transforms 00194 * a real vector to a half-complex real vector (with a forward plan) or 00195 * a half-complex real vector to a real vector (with a reverse plan). 00196 * If you're not sure what this means, don't use this routine. 00197 * The input and output vectors (and their data) must be distinct pointers. 00198 * 00199 * \verb+XLALREAL4PowerSpectrum+ computes a real power spectrum of the 00200 * input real vector and a forward FFT plan. 00201 * 00202 * \subsubsection*{Return Values} 00203 * 00204 * Upon success, 00205 * \verb+XLALCreateREAL4FFTPlan+, 00206 * \verb+XLALCreateForwardREAL4FFTPlan+, and 00207 * \verb+XLALCreateReverseREAL4FFTPlan+ return a pointer to a newly-allocated 00208 * FFT plan. Upon failure, they return a \verb+NULL+ pointer and set 00209 * \verb+xlalErrno+ to one of the following values: 00210 * \verb+XLAL_EBADLEN+ if \verb+size+ is not greater than zero, 00211 * \verb+XLAL_ENOMEM+ if a memory allocation failed, or 00212 * \verb+XLAL_EFAILED+ if the FFTW plan creation routine failed. 00213 * 00214 * \verb+XLALDestroyREAL4FFTPlan+ does not return any value but, upon 00215 * failure, it will set \verb+xlalErrno+ to one of the following values: 00216 * \verb+XLAL_EFAULT+ if the routine is provided a \verb+NULL+ pointer, or 00217 * \verb+XLAL_EINVAL+ if the contents of the plan are invalid (e.g., if the 00218 * routine is provided a plan that had been previously destroyed). 00219 * 00220 * \verb+XLALREAL4ForwardFFT+, 00221 * \verb+XLALREAL4ReverseFFT+, 00222 * \verb+XLALREAL4VectorFFT+, and 00223 * \verb+XLALREAL4PowerSpectrum+ return the value 0 upon succes; upon 00224 * failure they return \verb+XLAL_FAILURE+ and set \verb+xlalErrno+ to 00225 * one of the following values: 00226 * \verb+XLAL_EFAULT+ if one of the input pointers is \verb+NULL+, 00227 * \verb+XLAL_EINVAL+ if the input, output, or plan structures appears 00228 * invalid or if the routine is passed a plan for the wrong transform 00229 * directions or if the input and output data pointers are not distinct 00230 * for \verb+XLALREAL4VectorFFT+, 00231 * \verb+XLAL_EBADLEN+ if the input and output vectors and the plan have 00232 * incompatible lengths, 00233 * \verb+XLAL_ENOMEM+ if a memory allocation of temporary internal memory 00234 * fails. 00235 * 00236 * As before, the \verb+REAL8+ versions of these routines behave the 00237 * same way but for double-precision transforms. 00238 * 00239 * \newpage\input{RealFFTTestC} 00240 **** </lalLaTeX> */ 00241 00242 /* 00243 * 00244 * XLAL REAL4 functions 00245 * 00246 */ 00247 00248 REAL4FFTPlan * XLALCreateREAL4FFTPlan( UINT4 size, int fwdflg, int measurelvl ); 00249 REAL4FFTPlan * XLALCreateForwardREAL4FFTPlan( UINT4 size, int measurelvl ); 00250 REAL4FFTPlan * XLALCreateReverseREAL4FFTPlan( UINT4 size, int measurelvl ); 00251 void XLALDestroyREAL4FFTPlan( REAL4FFTPlan *plan ); 00252 00253 int XLALREAL4ForwardFFT( COMPLEX8Vector *output, const REAL4Vector *input, 00254 const REAL4FFTPlan *plan ); 00255 int XLALREAL4ReverseFFT( REAL4Vector *output, const COMPLEX8Vector *input, 00256 const REAL4FFTPlan *plan ); 00257 int XLALREAL4VectorFFT( REAL4Vector *output, const REAL4Vector *input, 00258 const REAL4FFTPlan *plan ); 00259 int XLALREAL4PowerSpectrum( REAL4Vector *spec, const REAL4Vector *data, 00260 const REAL4FFTPlan *plan ); 00261 00262 /* 00263 * 00264 * XLAL REAL8 functions 00265 * 00266 */ 00267 00268 REAL8FFTPlan * XLALCreateREAL8FFTPlan( UINT4 size, int fwdflg, int measurelvl ); 00269 REAL8FFTPlan * XLALCreateForwardREAL8FFTPlan( UINT4 size, int measurelvl ); 00270 REAL8FFTPlan * XLALCreateReverseREAL8FFTPlan( UINT4 size, int measurelvl ); 00271 void XLALDestroyREAL8FFTPlan( REAL8FFTPlan *plan ); 00272 00273 int XLALREAL8ForwardFFT( COMPLEX16Vector *output, REAL8Vector *input, 00274 const REAL8FFTPlan *plan ); 00275 int XLALREAL8ReverseFFT( REAL8Vector *output, COMPLEX16Vector *input, 00276 const REAL8FFTPlan *plan ); 00277 int XLALREAL8VectorFFT( REAL8Vector *output, REAL8Vector *input, 00278 const REAL8FFTPlan *plan ); 00279 int XLALREAL8PowerSpectrum( REAL8Vector *spec, REAL8Vector *data, 00280 const REAL8FFTPlan *plan ); 00281 00282 /* 00283 * 00284 * LAL REAL4 functions 00285 * 00286 */ 00287 00288 void 00289 LALCreateForwardREAL4FFTPlan( 00290 LALStatus *status, 00291 REAL4FFTPlan **plan, 00292 UINT4 size, 00293 INT4 measure 00294 ); 00295 #define LALCreateForwardRealFFTPlan LALCreateForwardREAL4FFTPlan 00296 00297 void 00298 LALCreateReverseREAL4FFTPlan( 00299 LALStatus *status, 00300 REAL4FFTPlan **plan, 00301 UINT4 size, 00302 INT4 measure 00303 ); 00304 #define LALCreateReverseRealFFTPlan LALCreateReverseREAL4FFTPlan 00305 00306 void 00307 LALDestroyREAL4FFTPlan( 00308 LALStatus *status, 00309 REAL4FFTPlan **plan 00310 ); 00311 #define LALDestroyRealFFTPlan LALDestroyREAL4FFTPlan 00312 00313 void 00314 LALForwardREAL4FFT( 00315 LALStatus *status, 00316 COMPLEX8Vector *output, 00317 REAL4Vector *input, 00318 REAL4FFTPlan *plan 00319 ); 00320 #define LALForwardRealFFT LALForwardREAL4FFT 00321 00322 void 00323 LALReverseREAL4FFT( 00324 LALStatus *status, 00325 REAL4Vector *output, 00326 COMPLEX8Vector *input, 00327 REAL4FFTPlan *plan 00328 ); 00329 #define LALReverseRealFFT LALReverseREAL4FFT 00330 00331 void 00332 LALREAL4PowerSpectrum ( 00333 LALStatus *status, 00334 REAL4Vector *spec, 00335 REAL4Vector *data, 00336 REAL4FFTPlan *plan 00337 ); 00338 #define LALRealPowerSpectrum LALREAL4PowerSpectrum 00339 00340 void 00341 LALREAL4VectorFFT( 00342 LALStatus *status, 00343 REAL4Vector *output, 00344 REAL4Vector *input, 00345 REAL4FFTPlan *plan 00346 ); 00347 00348 /* 00349 * 00350 * LAL REAL8 functions 00351 * 00352 */ 00353 00354 void 00355 LALCreateForwardREAL8FFTPlan( 00356 LALStatus *status, 00357 REAL8FFTPlan **plan, 00358 UINT4 size, 00359 INT4 measure 00360 ); 00361 00362 void 00363 LALCreateReverseREAL8FFTPlan( 00364 LALStatus *status, 00365 REAL8FFTPlan **plan, 00366 UINT4 size, 00367 INT4 measure 00368 ); 00369 00370 void 00371 LALDestroyREAL8FFTPlan( 00372 LALStatus *status, 00373 REAL8FFTPlan **plan 00374 ); 00375 00376 void 00377 LALForwardREAL8FFT( 00378 LALStatus *status, 00379 COMPLEX16Vector *output, 00380 REAL8Vector *input, 00381 REAL8FFTPlan *plan 00382 ); 00383 00384 void 00385 LALReverseREAL8FFT( 00386 LALStatus *status, 00387 REAL8Vector *output, 00388 COMPLEX16Vector *input, 00389 REAL8FFTPlan *plan 00390 ); 00391 00392 void 00393 LALREAL8PowerSpectrum ( 00394 LALStatus *status, 00395 REAL8Vector *spec, 00396 REAL8Vector *data, 00397 REAL8FFTPlan *plan 00398 ); 00399 00400 void 00401 LALREAL8VectorFFT( 00402 LALStatus *status, 00403 REAL8Vector *output, 00404 REAL8Vector *input, 00405 REAL8FFTPlan *plan 00406 ); 00407 00408 #ifdef __cplusplus 00409 #pragma { 00410 } 00411 #endif 00412 00413 #endif /* _REALFFT_H */
1.5.2