fftwaffar.cc 11.7 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
/*! \file fftwaffar.cc
 * \brief engine to transfrom several signals at once (implementation)
 * 
 * ----------------------------------------------------------------------------
 * 
 * \author Thomas Forbriger
 * \date 13/05/2011
 * 
 * engine to transfrom several signals at once (implementation)
 * 
 * Copyright (c) 2011 by Thomas Forbriger (BFO Schiltach) 
 *
 * ----
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
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
 * 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
 * ----
 * 
 * REVISIONS and CHANGES 
 *  - 13/05/2011   V1.0   Thomas Forbriger
 *  - 27/05/2011   V1.1   added copy constructor (required for class member
 *                        initialization) and default constructor and
 *                        assignment operator
 * 
 * ============================================================================
 */
#define TF_FFTWAFFAR_CC_VERSION \
  "TF_FFTWAFFAR_CC   V1.1"

#include <iostream>
#include <fourier/fftwaffar.h>
#include <fourier/error.h>
#include <aff/Carray.h>
#include <aff/slice.h>
#include <aff/shaper.h>

namespace fourier {

  namespace fft {

    DRFFTWAFFArrayEngine::DRFFTWAFFArrayEngine(const int& nseis,
                                               const int& nsamp)
      : Mseriesarray(TAseries(aff::Shaper(0,nsamp-1)(0,nseis-1))),
      Mspectrumarray(TAspectrum(aff::Shaper(0,DRFFTWAFFArrayEngine::ncoeff(nsamp)-1)(0,nseis-1))),
      Mplanr2c(0), Mplanc2r(0) 
    {
      this->checkconsistency();
    }

    /*----------------------------------------------------------------------*/

    DRFFTWAFFArrayEngine::DRFFTWAFFArrayEngine(const TAseries& series,
                                               const TAspectrum& spec)
      : Mseriesarray(series), Mspectrumarray(spec),
      Mplanr2c(0), Mplanc2r(0) 
    { this->checkconsistency(); }

    /*----------------------------------------------------------------------*/

    DRFFTWAFFArrayEngine::DRFFTWAFFArrayEngine(const DRFFTWAFFArrayEngine& e)
      : Mseriesarray(e.Mseriesarray), Mspectrumarray(e.Mspectrumarray),
      Mplanr2c(0), Mplanc2r(0) 
    { }

    /*----------------------------------------------------------------------*/

    DRFFTWAFFArrayEngine::DRFFTWAFFArrayEngine()
      : Mseriesarray(TAseries(aff::Shaper(0,4-1)(0,0))),
      Mspectrumarray(TAspectrum(aff::Shaper(0,DRFFTWAFFArrayEngine::ncoeff(4)-1)(0,0))),
      Mplanr2c(0), Mplanc2r(0) 
    { }

    /*----------------------------------------------------------------------*/

    DRFFTWAFFArrayEngine& DRFFTWAFFArrayEngine::operator=(const DRFFTWAFFArrayEngine& e)
    {
      this->delete_plans();
      Mseriesarray=e.Mseriesarray;
      Mspectrumarray=e.Mspectrumarray;
      return (*this);
    } // DRFFTWAFFArrayEngine& DRFFTWAFFArrayEngine::operator=(const DRFFTWAFFArrayEngine& e)

    /*----------------------------------------------------------------------*/

    //! delete plan.
    DRFFTWAFFArrayEngine::~DRFFTWAFFArrayEngine()
    {
      this->delete_plans();
    } // DRFFTWAFFArrayEngine::~DRFFTWAFFArrayEngine()

    /*----------------------------------------------------------------------*/

    //! delete plans.
    void DRFFTWAFFArrayEngine::delete_plans() 
    {
      if (Mplanr2c != 0) { fftw_destroy_plan(Mplanr2c); }
      if (Mplanc2r != 0) { fftw_destroy_plan(Mplanc2r); }
      Mplanr2c=0;
      Mplanc2r=0;
    } // DRFFTWAFFArrayEngine::delete_plans()

    /*----------------------------------------------------------------------*/

    //! execute c2r plan
    void DRFFTWAFFArrayEngine::c2r() 
    {
      if (Mplanc2r == 0)
      {
        this->createplanc2r();
      }
      fftw_execute(Mplanc2r);
    } // DRFFTWAFFArrayEngine::c2r()

    /*----------------------------------------------------------------------*/

    //! execute r2c plan
    void DRFFTWAFFArrayEngine::r2c() 
    {
      if (Mplanr2c == 0)
      {
        this->createplanr2c();
      }
      fftw_execute(Mplanr2c);
    } // DRFFTWAFFArrayEngine::r2c()

    /*----------------------------------------------------------------------*/

    void DRFFTWAFFArrayEngine::checkconsistency()
    {
      FOURIER_assert(Mspectrumarray.size(0)==
                     DRFFTWAFFArrayEngine::ncoeff(Mseriesarray.size(0)),
                     "ERROR: sample size inconsistent");
      FOURIER_assert(Mseriesarray.size(1)==Mspectrumarray.size(1),
                     "ERROR: inconsistent number of signals");
      FOURIER_assert(((Mseriesarray.size(2)==1) && (Mspectrumarray.size(2)==1)),
                     "ERROR: two-dimensional arrays only");
      FOURIER_assert(((Mseriesarray.size(3)==1) && (Mspectrumarray.size(3)==1)),
                     "ERROR: two-dimensional arrays only");
    } // void DRFFTWAFFArrayEngine::checkconsistency()

    /*----------------------------------------------------------------------*/

    void DRFFTWAFFArrayEngine::createplanr2c()
    {
      if (Mplanr2c==0)
      {
        // overall parameters
        // ------------------
        // tranformation only along one dimension
        const int rank=1;
        // size of each series
        int n=Mseriesarray.size(0);
        // number of signals to be transformed
        const int howmany=Mseriesarray.size(1);
        // quick design
        const unsigned flags=FFTW_ESTIMATE;

        // input array
        // -----------
        aff::CArray<TAseries::Tvalue> Cseriesarray(Mseriesarray);
        // one-dimensional transform: use default
        int* inembed=0;
        // distance to next sample
        const int istride=Cseriesarray.stride(0);
        // distance to next signal
        const int idist=Cseriesarray.stride(1);
        // casted pointer
        double* in=Cseriesarray.castedpointer<double>();

        // output array
        // ------------
        aff::CArray<TAspectrum::Tvalue> Cspectrumarray(Mspectrumarray);
        // one-dimensional transform: use default
        int* onembed=0;
        // distance to next sample
        const int ostride=Cspectrumarray.stride(0);
        // distance to next signal
        const int odist=Cspectrumarray.stride(1);
        // casted pointer
        fftw_complex* out=Cspectrumarray.castedpointer<fftw_complex>();

        // create plan
        Mplanr2c=fftw_plan_many_dft_r2c(rank, &n, howmany,
                                        in, inembed, istride, idist,
                                        out, onembed, ostride, odist,
                                        flags);
        FOURIER_assert(Mplanr2c!=0, "ERROR: creating r2c plan");
      }
    } // void DRFFTWAFFArrayEngine::createplanr2c()

    /*----------------------------------------------------------------------*/

    void DRFFTWAFFArrayEngine::createplanc2r()
    {
      if (Mplanc2r==0)
      {
        // overall parameters
        // ------------------
        // tranformation only along one dimension
        const int rank=1;
        // size of each series
        int n=Mseriesarray.size(0);
        // number of signals to be transformed
        const int howmany=Mseriesarray.size(1);
        // quick design
        const unsigned flags=FFTW_ESTIMATE;

        // input array
        // -----------
        aff::CArray<TAspectrum::Tvalue> Cspectrumarray(Mspectrumarray);
        // one-dimensional transform: use default
        int* inembed=0;
        // distance to next sample
        const int istride=Cspectrumarray.stride(0);
        // distance to next signal
        const int idist=Cspectrumarray.stride(1);
        // casted pointer
        fftw_complex* in=Cspectrumarray.castedpointer<fftw_complex>();

        // output array
        // ------------
        aff::CArray<TAseries::Tvalue> Cseriesarray(Mseriesarray);
        // one-dimensional transform: use default
        int* onembed=0;
        // distance to next sample
        const int ostride=Cseriesarray.stride(0);
        // distance to next signal
        const int odist=Cseriesarray.stride(1);
        // casted pointer
        double* out=Cseriesarray.castedpointer<double>();

        // create plan
        Mplanc2r=fftw_plan_many_dft_c2r(rank, &n, howmany,
                                        in, inembed, istride, idist,
                                        out, onembed, ostride, odist,
                                        flags);
        FOURIER_assert(Mplanc2r!=0, "ERROR: creating c2r plan");
      }
    } // void DRFFTWAFFArrayEngine::createplanc2r()
    
    /*----------------------------------------------------------------------*/

253
254
255
256
257
258
259
260
261
    /*! \brief Return appropriate scaling factor for sampling interval dt.
     *
     * Factor to be applied when transforming to time domain.
     *
     * \param[in] dt sampling interval
     * \return scalar factor to be applied to all samples
     *
     * \sa \ref sec_fftw3_integral_transform
     */
262
263
264
265
266
267
268
269
    DRFFTWAFFArrayEngine::Tsample
      DRFFTWAFFArrayEngine::scale_series(const Tsample& dt) const
    {
      return(1./(Mseriesarray.size(0)*dt));
    } // Tsample DRFFTWAFFArrayEngine::scale_series(const Tsample& dt) const

    /*----------------------------------------------------------------------*/
    
270
271
272
273
274
275
276
277
278
    /*! \brief Return appropriate scaling factor for sampling interval dt.
     *
     * Factor to be applied when transforming to Fourier domain.
     *
     * \param[in] dt sampling interval
     * \return scalar factor to be applied to all samples
     *
     * \sa \ref sec_fftw3_integral_transform
     */
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
    DRFFTWAFFArrayEngine::Tsample
      DRFFTWAFFArrayEngine::scale_spectrum(const Tsample& dt) const
    {
      return(dt);
    } // Tsample DRFFTWAFFArrayEngine::scale_spectrum(const Tsample& dt) const

    /*----------------------------------------------------------------------*/

    unsigned int DRFFTWAFFArrayEngine::nsamples() const
    {
      return Mseriesarray.size(0);
    } // unsigned int DRFFTWAFFArrayEngine::nsamples() const

    /*----------------------------------------------------------------------*/

    unsigned int DRFFTWAFFArrayEngine::nfrequencies() const
    {
      return Mspectrumarray.size(0);
    } // unsigned int DRFFTWAFFArrayEngine::nfrequencies() const

    /*----------------------------------------------------------------------*/

    unsigned int DRFFTWAFFArrayEngine::nseries() const
    {
      return Mseriesarray.size(1);
    } // unsigned int DRFFTWAFFArrayEngine::nseries() const

    /*----------------------------------------------------------------------*/

    //! \brief return a reference to the time series i
    DRFFTWAFFArrayEngine::TAseries 
      DRFFTWAFFArrayEngine::series(const unsigned int& i) const
    {
      FOURIER_assert(((i>=0) && (i<this->nseries())),
                     "ERROR: signal index is out of range");
      return(aff::slice(Mseriesarray)()(i));
    } // DRFFTWAFFArrayEngine::TAseries DRFFTWAFFArrayEngine::series(i) const

    /*----------------------------------------------------------------------*/

    //! \brief return a reference to the Fourier coefficients of series i
    DRFFTWAFFArrayEngine::TAspectrum 
      DRFFTWAFFArrayEngine::spectrum(const unsigned int& i) const
    {
      FOURIER_assert(((i>=0) && (i<this->nseries())),
                     "ERROR: signal index is out of range");
      return(aff::slice(Mspectrumarray)()(i));
    } // DRFFTWAFFArrayEngine::TAspectrum DRFFTWAFFArrayEngine::spectrum(i) const

  } // namespace fft

} // namespace fourier

/* ----- END OF fftwaffar.cc ----- */