fftpack_FFT.h 9.91 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16
//
// IPPL FFT
//
// Copyright (c) 2008-2018
// Paul Scherrer Institut, Villigen PSI, Switzerland
// All rights reserved.
//
// OPAL is licensed under GNU GPL version 3.
//

/*
  Prototypes for accessing Fortran 1D FFT routines from
  Netlib, and definitions for templated class FFTPACK, which acts as an
  FFT engine for the FFT class, providing storage for trigonometric
  information and performing the 1D FFTs as needed.
*/
gsell's avatar
gsell committed
17 18 19 20 21 22 23 24 25

#ifndef IPPL_FFT_FFTPACK_FFT_H
#define IPPL_FFT_FFTPACK_FFT_H

#include "Utility/PAssert.h"
#include "Utility/IpplInfo.h"

// FFTPACK function prototypes for Fortran routines
extern "C" {
26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47
    // double-precision CC FFT
    void cffti (size_t n, double& wsave);
    void cfftf (size_t n, double& r, double& wsave);
    void cfftb (size_t n, double& r, double& wsave);
    // double-precision RC FFT
    void rffti (size_t n, double& wsave);
    void rfftf (size_t n, double& r, double& wsave);
    void rfftb (size_t n, double& r, double& wsave);
    // double-precision sine transform
    void sinti (size_t n, double& wsave);
    void sint  (size_t n, double& r, double& wsave);
    // single-precision CC FFT
    void fcffti (size_t n, float& wsave);
    void fcfftf (size_t n, float& r, float& wsave);
    void fcfftb (size_t n, float& r, float& wsave);
    // single-precision RC FFT
    void frffti (size_t n, float& wsave);
    void frfftf (size_t n, float& r, float& wsave);
    void frfftb (size_t n, float& r, float& wsave);
    // single-precision sine transform
    void fsinti (size_t n, float& wsave);
    void fsint (size_t n, float& r, float& wsave);
gsell's avatar
gsell committed
48 49 50 51 52 53 54 55 56 57 58 59 60
}


// FFTPACK_wrap provides static functions that wrap the Fortran functions
// in a common interface.  We specialize this class on precision type.
template <class T>
class FFTPACK_wrap {};

// Specialization for float
template <>
class FFTPACK_wrap<float> {

public:
61 62 63 64 65 66 67 68 69 70 71 72 73 74
    // interface functions used by class FFTPACK

    // initialization functions for CC FFT, RC FFT, and sine transform
    static void ccffti(size_t n, float* wsave) { fcffti (n, *wsave); }
    static void rcffti(size_t n, float* wsave) { frffti (n, *wsave); }
    static void rrffti(size_t n, float* wsave) { fsinti (n, *wsave); }
    // forward and backward CC FFT
    static void ccfftf(size_t n, float* r, float* wsave) { fcfftf (n, *r, *wsave); }
    static void ccfftb(size_t n, float* r, float* wsave) { fcfftb (n, *r, *wsave); }
    // forward and backward RC FFT
    static void rcfftf(size_t n, float* r, float* wsave) { frfftf (n, *r, *wsave); }
    static void rcfftb(size_t n, float* r, float* wsave) { frfftb (n, *r, *wsave); }
    // sine transform
    static void rrfft(size_t n, float* r, float* wsave) { fsint (n, *r, *wsave); }
gsell's avatar
gsell committed
75 76 77 78 79 80 81
};

// Specialization for double
template <>
class FFTPACK_wrap<double> {

public:
82 83 84 85 86 87 88 89 90 91 92 93 94 95
    // interface functions used by class FFTPACK

    // initialization functions for CC FFT, RC FFT, and sine transform
    static void ccffti(size_t n, double* wsave) { cffti (n, *wsave); }
    static void rcffti(size_t n, double* wsave) { rffti (n, *wsave); }
    static void rrffti(size_t n, double* wsave) { sinti (n, *wsave); }
    // forward and backward CC FFT
    static void ccfftf(size_t n, double* r, double* wsave) {cfftf (n, *r, *wsave);}
    static void ccfftb(size_t n, double* r, double* wsave) {cfftb (n, *r, *wsave);}
    // forward and backward RC FFT
    static void rcfftf(size_t n, double* r, double* wsave) {rfftf (n, *r, *wsave);}
    static void rcfftb(size_t n, double* r, double* wsave) {rfftb (n, *r, *wsave);}
    // sine transform
    static void rrfft(size_t n, double* r, double* wsave) { sint (n, *r, *wsave); }
gsell's avatar
gsell committed
96 97 98 99 100 101 102 103 104
};


// Definition of FFT engine class FFTPACK
template <class T>
class FFTPACK {

public:

105 106
    // definition of complex type
    typedef std::complex<T> Complex_t;
gsell's avatar
gsell committed
107

108 109
    // Trivial constructor.  Do the real work in setup function.
    FFTPACK(void) {}
gsell's avatar
gsell committed
110

111 112
    // destructor
    ~FFTPACK(void);
gsell's avatar
gsell committed
113

114 115 116 117 118
    // setup internal storage and prepare to perform FFTs
    // inputs are number of dimensions to transform, the transform types,
    // and the lengths of these dimensions.
    void setup(unsigned numTransformDims, const int* transformTypes,
               const int* axisLengths);
gsell's avatar
gsell committed
119

120 121
    // invoke FFT on complex data for given dimension and direction
    void callFFT(unsigned transformDim, int direction, Complex_t* data);
gsell's avatar
gsell committed
122

123 124
    // invoke FFT on real data for given dimension and direction
    void callFFT(unsigned transformDim, int direction, T* data);
gsell's avatar
gsell committed
125 126 127

private:

128 129 130 131
    unsigned numTransformDims_m;  // number of dimensions to transform
    int* transformType_m;         // transform type for each dimension
    int* axisLength_m;            // length of each transform dimension
    T** trig_m;                   // trigonometric tables
gsell's avatar
gsell committed
132 133 134 135 136 137 138 139 140 141

};


// Inline member function definitions

// destructor
template <class T>
inline
FFTPACK<T>::~FFTPACK(void) {
142 143 144 145 146 147
    // delete storage
    for (unsigned d=0; d<numTransformDims_m; ++d)
        delete [] trig_m[d];
    delete [] trig_m;
    delete [] transformType_m;
    delete [] axisLength_m;
gsell's avatar
gsell committed
148 149 150 151 152 153 154 155
}

// setup internal storage to prepare for FFTs
template <class T>
inline void
FFTPACK<T>::setup(unsigned numTransformDims, const int* transformTypes,
                  const int* axisLengths) {

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
    // store transform types and lengths for each transform dim
    numTransformDims_m = numTransformDims;
    transformType_m = new int[numTransformDims_m];
    axisLength_m = new int[numTransformDims_m];
    unsigned d;
    for (d=0; d<numTransformDims_m; ++d) {
        transformType_m[d] = transformTypes[d];
        axisLength_m[d] = axisLengths[d];
    }

    // allocate and initialize trig table
    trig_m = new T*[numTransformDims_m];
    for (d=0; d<numTransformDims_m; ++d) {
        switch (transformType_m[d]) {
        case 0:  // CC FFT
            trig_m[d] = new T[4 * axisLength_m[d] + 15];
            FFTPACK_wrap<T>::ccffti(axisLength_m[d], trig_m[d]);
            break;
        case 1:  // RC FFT
            trig_m[d] = new T[2 * axisLength_m[d] + 15];
            FFTPACK_wrap<T>::rcffti(axisLength_m[d], trig_m[d]);
            break;
        case 2:  // Sine transform
            trig_m[d] = new T[static_cast<int>(2.5 * axisLength_m[d] + 0.5) + 15];
            FFTPACK_wrap<T>::rrffti(axisLength_m[d], trig_m[d]);
            break;
        default:
            ERRORMSG("Unknown transform type requested!!" << endl);
            break;
        }
gsell's avatar
gsell committed
186 187
    }

188
    return;
gsell's avatar
gsell committed
189 190 191 192 193 194 195 196
}

// invoke FFT on complex data for given dimension and direction
template <class T>
inline void
FFTPACK<T>::callFFT(unsigned transformDim, int direction,
                    FFTPACK<T>::Complex_t* data) {

197 198 199
    // check transform dimension and direction arguments
    PAssert_LT(transformDim, numTransformDims_m);
    PAssert_EQ(std::abs(direction), 1);
gsell's avatar
gsell committed
200

201 202
    // cast complex number pointer to T* for calling Fortran routines
    T* rdata = reinterpret_cast<T*>(data);
gsell's avatar
gsell committed
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
    // branch on transform type for this dimension
    switch (transformType_m[transformDim]) {
    case 0:  // CC FFT
        if (direction == +1) {
            // call forward complex-to-complex FFT
            FFTPACK_wrap<T>::ccfftf(axisLength_m[transformDim], rdata,
                                    trig_m[transformDim]);
        }
        else {
            // call backward complex-to-complex FFT
            FFTPACK_wrap<T>::ccfftb(axisLength_m[transformDim], rdata,
                                    trig_m[transformDim]);
        }
        break;
    case 1:  // RC FFT
        if (direction == +1) {
            // call forward real-to-complex FFT
            FFTPACK_wrap<T>::rcfftf(axisLength_m[transformDim], rdata,
                                    trig_m[transformDim]);
            // rearrange output to conform with SGI format for complex result
            int clen = axisLength_m[transformDim]/2 + 1;
            data[clen-1] = Complex_t(imag(data[clen-2]),0.0);
            for (int i = clen-2; i > 0; --i)
                data[i] = Complex_t(imag(data[i-1]),real(data[i]));
            data[0] = Complex_t(real(data[0]),0.0);
        }
        else {                
            // rearrange input to conform with Netlib format for complex modes
            int clen = axisLength_m[transformDim]/2 + 1;
            data[0] = Complex_t(real(data[0]),real(data[1]));
            for (int i = 1; i < clen-1; ++i)
                data[i] = Complex_t(imag(data[i]),real(data[i+1]));
            // call backward complex-to-real FFT
            FFTPACK_wrap<T>::rcfftb(axisLength_m[transformDim], rdata,
                                    trig_m[transformDim]);
        }
        break;
    case 2:  // Sine transform
        ERRORMSG("Input for real-to-real FFT should be real!!" << endl);
        break;
    default:
        ERRORMSG("Unknown transform type requested!!" << endl);
        break;
gsell's avatar
gsell committed
247
    }
248 249

    return;
gsell's avatar
gsell committed
250 251 252 253 254 255 256
}

// invoke FFT on real data for given dimension and direction
template <class T>
inline void
FFTPACK<T>::callFFT(unsigned transformDim, int direction, T* data) {

257 258
    // check transform dimension and direction arguments
    PAssert_LT(transformDim, numTransformDims_m);
259 260 261
    // avoid unused variable warning if we compile with Release
    // :FIXME: remove direction
    (void)direction;
262
    PAssert_EQ(std::abs(direction), 1);
gsell's avatar
gsell committed
263

264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280
    // branch on transform type for this dimension
    switch (transformType_m[transformDim]) {
    case 0:  // CC FFT
        ERRORMSG("Input for complex-to-complex FFT should be complex!!" << endl);
        break;
    case 1:  // RC FFT
        ERRORMSG("real-to-complex FFT uses complex input!!" << endl);
        break;
    case 2:  // Sine transform
        // invoke the real-to-real transform on the data
        FFTPACK_wrap<T>::rrfft(axisLength_m[transformDim], data,
                               trig_m[transformDim]);
        break;
    default:
        ERRORMSG("Unknown transform type requested!!" << endl);
        break;
    }
gsell's avatar
gsell committed
281

282 283
    return;
}
gsell's avatar
gsell committed
284 285
#endif // IPPL_FFT_FFTPACK_FFT_H

286 287 288 289
// vi: set et ts=4 sw=4 sts=4:
// Local Variables:
// mode:c
// c-basic-offset: 4
gsell's avatar
gsell committed
290 291
// indent-tabs-mode: nil
// require-final-newline: nil
292
// End:
gsell's avatar
gsell committed
293