fftwaff.h 7.85 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
/*! \file fftwaff.h
 * \brief use fftw together with aff containers (prototypes)
 * 
 * ----------------------------------------------------------------------------
 * 
 * \author Thomas Forbriger
 * \date 11/07/2006
 * 
 * use fftw together with aff containers (prototypes)
 *
 * link with -lrfftw -lfftw -lm -laff
 *
 * ----
14
 * libfourier is free software; you can redistribute it and/or modify
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation; either version 2 of the License, or
 * (at your option) any later version. 
 * 
 * This program is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 * 
 * You should have received a copy of the GNU General Public License
 * along with this program; if not, write to the Free Software
 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
 * ----
 *
 * 
 * Copyright (c) 2006 by Thomas Forbriger (BFO Schiltach) 
 * 
 * REVISIONS and CHANGES 
 *  - 11/07/2006   V1.0   Thomas Forbriger
 *  - 12/09/2007   V1.1   first running version
 *  - 07/10/2010   V1.2   
 *                        - migrate to FFTW3:
 *                        - use different fftw header file (fftw3.h)
 *                        - type of FFTW plan has changed
 *  - 02/10/2012   V1.3  
 *                        - provide size calculation functions
 *
 * Migration to FFTW3:
 * - The location of the arrays in memory are part of the plan. 
 *   Consequently we have to either create a new plan for each transform or to
 *   allocate the working array together with the plan and keep it available in
 *   the background. The solution of choice is to keep arrays. We will
 *   introduce a control parameter to the arguments of the constructor which
 *   optionally allows to get rid of the arrays after each transformation.
 * - FFTW3 uses wisdom by default.
 *   \code
 *   void fftw_forget_wisdom(void);
 *   \endcode
 *   should be called by the destructor of the class. Alternatively
 *   \code
 *   void fftw_cleanup(void);
 *   \endcode 
 *   could be called which apparently is even more rigorous. Plans, however,
 *   must be destroyed prior to calling \c fftw_cleanup.
 * - The FFTW documentation states:
 *   To the extent that this is true, if you have a variable 
 *   \code
 *   complex<double> *x,
 *   \endcode
 *   you can pass it directly to FFTW via
 *   \code
 *   reinterpret_cast<fftw_complex*>(x). 
 *   \endcode
 * 
 * ============================================================================
 */

// include guard
#ifndef TF_FFTWAFF_H_VERSION

#define TF_FFTWAFF_H_VERSION \
  "TF_FFTWAFF_H   V1.3"

#include<complex>
#ifdef FFTWFALLBACK
#include<drfftw.h>
#else
#include<fftw3.h>
#endif
#include<aff/series.h>

namespace fourier {

  /*! All Fourier transform stuff is collected here.
   */
  namespace fft {

    /*! A rigid class to do simple transforms using libdrfftw.a.
     *
     * uses real double arrays
     *
     * How to use this class:
     *
     * You may create one instance of this class and use it to transform
     * several signals in both directions in turn. The class itself takes care
     * of the transform size and creates a new plan if necessary. FFTs are
     * invoked by the bracket operators. You should use the class object like
     * a function. The scaling operators (taking the sampling interval as one
     * of their arguments) return a series or spectrum scaled appropriately to
     * match the values of samples from the corresponding Fourier integral
     * transform (usual convention with \f$dt\f$ and \f$df\f$
     * integrals - not \f$d \omega\f$).
     *
     * \note
     * The appropriate number of samples for the time series obtained from a
     * given set of Fourier coefficients is non-unique, if only the Fourier
     * coefficients for positive frequencies are given, as is the case here.
     * A time series with \f$2n\f$ samples and a time series with \f$2n+1\f$
     * samples both will results in a set of \f$n+1\f$ Fourier coefficients.
     * The method to determine the required number of time samples from the
     * imaginary part of the Nyquist Fourier coefficients only works if the
     * Nyquist coefficient (at least real part) is finite - which is not the
     * case for most of our signals. This way we will obtain \f$n+1\f$ Fourier
     * coefficients from \f$2n\f$ time series samples. When reconstructing the
     * time series, we will obtain \f$2n+1\f$ samples with a sampling interval
     * now being \f$2n/(2n+1)\f$ of the original interval. For this reason it
     * is appropriate to set the size of the expected time series explicitely
     * either by using the constructor 
     * DRFFTWAFF::DRFFTWAFF(const unsigned int& n) or by
     * calling DRFFTWAFF::size(const unsigned int& s) prior to
     * DRFFTWAFF::operator()(const Tspectrum::Tcoc& s).
     *
     * \sa \ref page_fftw3
     */
    class DRFFTWAFF {
      public:
        typedef double Tsample;
        typedef std::complex<Tsample> Tcoeff;
        typedef aff::Series<Tsample> Tseries;
        typedef aff::Series<Tcoeff> Tspectrum;
#ifdef FFTWFALLBACK
        DRFFTWAFF(const unsigned int& n):
          Msize(n), Mplan_forward(0), Mplan_backward(0) { }
        DRFFTWAFF():
          Msize(0), Mplan_forward(0), Mplan_backward(0) { }
#else
        DRFFTWAFF(const unsigned int& n, const bool& deletearrays=false):
          Msize(n), Mplan_forward(0), Mplan_backward(0),
          Mseriesarray(0), Mspectrumarray(0),
          Mdeletearrays(deletearrays) { }
        DRFFTWAFF(const bool& deletearrays=false):
          Msize(0), Mplan_forward(0), Mplan_backward(0),
          Mseriesarray(0), Mspectrumarray(0),
          Mdeletearrays(deletearrays) { }
#endif
        ~DRFFTWAFF();
        Tspectrum operator()(const Tseries::Tcoc& s,
                             const bool& debug=false) const;
        Tseries operator()(const Tspectrum::Tcoc& s,
                           const bool& debug=false) const;
        Tspectrum operator()(const Tseries::Tcoc& s,
                             const Tsample& dt,
                             const bool& debug=false) const;
        Tseries operator()(const Tspectrum::Tcoc& s,
                           const Tsample& dt,
                           const bool& debug=false) const;
        Tsample scale_series(const Tsample& dt) const;
        Tsample scale_spectrum(const Tsample& dt) const;
        unsigned int size() const { return(Msize); }
        void size(const unsigned int& s) const { this->set_size(s); }

        //! return number of coefficients for given number of samples
        inline
        static unsigned int spectrumsize(const unsigned int& n) 
        { return(n/2+1); }

        //! return number of samples for given number of coefficients
        inline
        static unsigned int seriessize(const unsigned int& n) 
        { return(n*2-1); }

      private:
        void create_plan_forward() const;
        void create_plan_backward() const;
        void delete_plans() const;
#ifndef FFTWFALLBACK
        void create_arrays() const;
        void delete_arrays() const;
        unsigned int ssize() const 
        { return(DRFFTWAFF::spectrumsize(this->size())); }
#endif
        void set_size(const unsigned int& n) const;
        mutable unsigned int Msize;
#ifdef FFTWFALLBACK
        mutable rfftw_plan Mplan_forward;
        mutable rfftw_plan Mplan_backward;
#else
        mutable fftw_plan Mplan_forward;
        mutable fftw_plan Mplan_backward;
        mutable double *Mseriesarray;
        mutable fftw_complex *Mspectrumarray;
        mutable Tspectrum Mspectrum;
        mutable Tseries Mseries;
        bool Mdeletearrays;
#endif
    }; // class DRFFTWAFF

  } // namespace ftt

} // namespace fourier

#endif // TF_FFTWAFF_H_VERSION (includeguard)

/* ----- END OF fftwaff.h ----- */