FFTW的C++接口

头文件:

/*
 * Copyright (c) 2008-2011 Zhang Ming (M. Zhang), [email protected]
 *
 * This program is free software; you can redistribute it and/or modify it
 * under the terms of the GNU General Public License as published by the
 * Free Software Foundation, either version 2 or any later version.
 *
 * Redistribution and use in source and binary forms, with or without
 * modification, are permitted provided that the following conditions are met:
 *
 * 1. Redistributions of source code must retain the above copyright notice,
 *    this list of conditions and the following disclaimer.
 *
 * 2. Redistributions in binary form must reproduce the above copyright
 *    notice, this list of conditions and the following disclaimer in the
 *    documentation and/or other materials provided with the distribution.
 *
 * 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. A copy of the GNU General Public License is available at:
 * http://www.fsf.org/licensing/licenses
 */


/*****************************************************************************
 *                                    fftw.h
 *
 * A simple C++ interface for FFTW. This file only provides the one dimension
 * DFT and IDFT for data type as "Vector<Type>" in SPLab, where "Type" can
 * be "double", "complex<double>, "float", "complex<float>, "long double",
 * "complex<long double>.
 *
 * If you need execute FFT many times in your program, it's not a good idear
 * for using this interface because of losing efficiency.
 *
 * Zhang Ming, 2010-01, Xi'an Jiaotong University.
 *****************************************************************************/


#ifndef FFTW_H
#define FFTW_H


#include <complex>
#include <fftw3.h>
#include <vector.h>


namespace splab
{

    void fftw( Vector<double>&,      Vector< complex<double> >& );
    void fftw( Vector<float>&,       Vector< complex<float> >& );
    void fftw( Vector<long double>&, Vector< complex<long double> >& );

    void fftw( Vector< complex<double> >&, Vector< complex<double> >& );
    void fftw( Vector< complex<float> >&,  Vector< complex<float> >& );
    void fftw( Vector< complex<long double> >&,
               Vector< complex<long double> >& );

    void ifftw( Vector< complex<double> >&, Vector<double>& );
    void ifftw( Vector< complex<float> >&,  Vector<float>& );
    void ifftw( Vector< complex<long double> >&, Vector<long double>& );

    void ifftw( Vector< complex<double> >&, Vector< complex<double> >& );
    void ifftw( Vector< complex<float> >&,  Vector< complex<float> >& );
    void ifftw( Vector< complex<long double> >&,
                Vector< complex<long double> >& );


    #include <fftw-impl.h>

}
// namespace splab


#endif
// FFTW_H

实现文件:

/*
 * Copyright (c) 2008-2011 Zhang Ming (M. Zhang), [email protected]
 *
 * This program is free software; you can redistribute it and/or modify it
 * under the terms of the GNU General Public License as published by the
 * Free Software Foundation, either version 2 or any later version.
 *
 * Redistribution and use in source and binary forms, with or without
 * modification, are permitted provided that the following conditions are met:
 *
 * 1. Redistributions of source code must retain the above copyright notice,
 *    this list of conditions and the following disclaimer.
 *
 * 2. Redistributions in binary form must reproduce the above copyright
 *    notice, this list of conditions and the following disclaimer in the
 *    documentation and/or other materials provided with the distribution.
 *
 * 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. A copy of the GNU General Public License is available at:
 * http://www.fsf.org/licensing/licenses
 */


/*****************************************************************************
 *                                 fftw-impl.h
 *
 * Implementation for C++ interface for FFTW.
 *
 * Zhang Ming, 2010-01, Xi'an Jiaotong University.
 *****************************************************************************/


/**
 * Real to complex DFT of 1D signal. If "xn" has N points, "Xk"
 * should has N/2+1 points.
 */
inline void fftw( Vector<double> &xn, Vector< complex<double> > &Xk )
{
    if( Xk.size() < xn.size()/2+1 )
        Xk.resize( xn.size()/2+1 );

    fftw_plan r2cP;
    r2cP = fftw_plan_dft_r2c_1d( xn.dim(), xn.begin(),
           reinterpret_cast<fftw_complex*>(Xk.begin()),
           FFTW_ESTIMATE );
    fftw_execute( r2cP );
    fftw_destroy_plan( r2cP );
}

inline void fftw( Vector<float> &xn, Vector< complex<float> > &Xk )
{
    if( Xk.size() < xn.size()/2+1 )
        Xk.resize( xn.size()/2+1 );

    fftwf_plan r2cP;
    r2cP = fftwf_plan_dft_r2c_1d( xn.dim(), xn.begin(),
           reinterpret_cast<fftwf_complex*>(Xk.begin()),
           FFTW_ESTIMATE );
    fftwf_execute( r2cP );
    fftwf_destroy_plan( r2cP );
}

inline void fftw( Vector<long double> &xn,
                  Vector< complex<long double> > &Xk )
{
    if( Xk.size() < xn.size()/2+1 )
        Xk.resize( xn.size()/2+1 );

    fftwl_plan r2cP;
    r2cP = fftwl_plan_dft_r2c_1d( xn.dim(), xn.begin(),
           reinterpret_cast<fftwl_complex*>(Xk.begin()),
           FFTW_ESTIMATE );
    fftwl_execute( r2cP );
    fftwl_destroy_plan( r2cP );
}


/**
 * Complex to complex DFT of 1D signal.
 */
inline void fftw( Vector< complex<double> > &xn,
                  Vector< complex<double> > &Xk )
{
    if( Xk.size() < xn.size() )
        Xk.resize( xn.size() );

    fftw_plan c2cP;
    c2cP = fftw_plan_dft_1d( xn.dim(),
           reinterpret_cast<fftw_complex*>(xn.begin()),
           reinterpret_cast<fftw_complex*>(Xk.begin()),
           FFTW_FORWARD, FFTW_ESTIMATE );
    fftw_execute( c2cP );
    fftw_destroy_plan( c2cP );
}

inline void fftw( Vector< complex<float> > &xn,
                  Vector< complex<float> > &Xk )
{
    if( Xk.size() < xn.size() )
        Xk.resize( xn.size() );

    fftwf_plan c2cP;
    c2cP = fftwf_plan_dft_1d( xn.dim(),
           reinterpret_cast<fftwf_complex*>(xn.begin()),
           reinterpret_cast<fftwf_complex*>(Xk.begin()),
           FFTW_FORWARD, FFTW_ESTIMATE );
    fftwf_execute( c2cP );
    fftwf_destroy_plan( c2cP );
}

inline void fftw( Vector< complex<long double> > &xn,
                  Vector< complex<long double> > &Xk )
{
    if( Xk.size() < xn.size() )
        Xk.resize( xn.size() );

    fftwl_plan c2cP;
    c2cP = fftwl_plan_dft_1d( xn.dim(),
           reinterpret_cast<fftwl_complex*>(xn.begin()),
           reinterpret_cast<fftwl_complex*>(Xk.begin()),
           FFTW_FORWARD, FFTW_ESTIMATE );
    fftwl_execute( c2cP );
    fftwl_destroy_plan( c2cP );
}


/**
 * Complex to real IDFT of 1D. If "xn" has N points, "Xk" should
 * has N/2+1 points.
 */
inline void ifftw( Vector< complex<double> > &Xk, Vector<double> &xn )
{
    assert( Xk.size() == xn.size()/2+1 );

    fftw_plan c2rP;
    c2rP = fftw_plan_dft_c2r_1d( xn.dim(),
           reinterpret_cast<fftw_complex*>(Xk.begin()),
           xn.begin(), FFTW_ESTIMATE );
    fftw_execute( c2rP );
    fftw_destroy_plan( c2rP );

    xn /= double( xn.dim() );
}

inline void ifftw( Vector< complex<float> > &Xk, Vector<float> &xn )
{
    assert( Xk.size() == xn.size()/2+1 );

    fftwf_plan c2rP;
    c2rP = fftwf_plan_dft_c2r_1d( xn.dim(),
           reinterpret_cast<fftwf_complex*>(Xk.begin()),
           xn.begin(), FFTW_ESTIMATE );
    fftwf_execute( c2rP );
    fftwf_destroy_plan( c2rP );

    xn /= float( xn.dim() );
}

inline void ifftw( Vector< complex<long double> > &Xk,
                   Vector<long double> &xn )
{
    assert( Xk.size() == xn.size()/2+1 );

    fftwl_plan c2rP;
    c2rP = fftwl_plan_dft_c2r_1d( xn.dim(),
           reinterpret_cast<fftwl_complex*>(Xk.begin()),
           xn.begin(), FFTW_ESTIMATE );
    fftwl_execute( c2rP );
    fftwl_destroy_plan( c2rP );

    xn /= (long double)( xn.dim() );
}


/**
 * Complex to complex IDFT of 1D.
 */
inline void ifftw( Vector< complex<double> > &Xk,
                   Vector< complex<double> > &xn )
{
    if( xn.size() < Xk.size() )
        xn.resize( Xk.size() );

    fftw_plan c2cP;
    c2cP = fftw_plan_dft_1d( xn.dim(),
           reinterpret_cast<fftw_complex*>(Xk.begin()),
           reinterpret_cast<fftw_complex*>(xn.begin()),
           FFTW_BACKWARD, FFTW_ESTIMATE );
    fftw_execute( c2cP );
    fftw_destroy_plan( c2cP );

    xn /= complex<double>( xn.dim(), 0.0 );
}

inline void ifftw( Vector< complex<float> > &Xk,
                   Vector< complex<float> > &xn )
{
    if( xn.size() < Xk.size() )
        xn.resize( Xk.size() );

    fftwf_plan c2cP;
    c2cP = fftwf_plan_dft_1d( xn.dim(),
           reinterpret_cast<fftwf_complex*>(Xk.begin()),
           reinterpret_cast<fftwf_complex*>(xn.begin()),
           FFTW_BACKWARD, FFTW_ESTIMATE );
    fftwf_execute( c2cP );
    fftwf_destroy_plan( c2cP );

    xn /= complex<float>( xn.dim(), 0.0 );
}

inline void ifftw( Vector< complex<long double> > &Xk,
                   Vector< complex<long double> > &xn )
{
    if( xn.size() < Xk.size() )
        xn.resize( Xk.size() );

    fftwl_plan c2cP;
    c2cP = fftwl_plan_dft_1d( xn.dim(),
           reinterpret_cast<fftwl_complex*>(Xk.begin()),
           reinterpret_cast<fftwl_complex*>(xn.begin()),
           FFTW_BACKWARD, FFTW_ESTIMATE );
    fftwl_execute( c2cP );
    fftwl_destroy_plan( c2cP );

    xn /= complex<long double>( xn.dim(), 0.0 );
}

测试代码:

/*****************************************************************************
 *                               fftw_test.cpp
 *
 * FFTW interface testing.
 *
 * Zhang Ming, 2010-01, Xi'an Jiaotong University.
 *****************************************************************************/


#define BOUNDS_CHECK

#include <iostream>
#include <iomanip>
#include <fftw.h>


using namespace std;
using namespace splab;


typedef float   Type;
const   int     N = 7;


int main()
{
    // complex to complex dft test...
    Vector< complex<Type> > xn( N );
    Vector< complex<Type> > yn( N );
    Vector< complex<Type> > Xk( N );

    for( int i=0; i<N; ++i )
    {
        Type theta = Type( 2*PI * i / N );
        xn[i] = complex<Type>( cos(theta), sin(theta) );
    }

    cout << setiosflags(ios::fixed) << setiosflags(ios::showpos) << setprecision(8);
    cout << "xn:   " << xn << endl << endl;

    fftw( xn, Xk );
    cout << "Xk=fft(xn):   " << Xk << endl << endl;
    ifftw( Xk, yn );
    cout << "xn-ifft(Xk):   " << xn-yn << endl << endl;

    // real to complex and complex to real dft test...
    Vector<Type> sn(N), tn(N);
    Vector< complex<Type> > Sk;
    for( int i=0; i<N; ++i )
    {
        Type theta = Type( 2*PI * i / N );
        sn[i] = sin(theta);
    }
    cout << "sn:   " << sn << endl;

    fftw( sn, Sk );
    cout << "Sk=fft(sn):   " << Sk << endl << endl;
    ifftw( Sk, tn );
    cout << "sn-ifft(Sk):   " << sn-tn << endl;

    return 0;
}

运行结果:

xn:   size: +7 by 1
(+1.00000000,+0.00000000)
(+0.62348980,+0.78183150)
(-0.22252095,+0.97492790)
(-0.90096885,+0.43388382)
(-0.90096885,-0.43388376)
(-0.22252101,-0.97492790)
(+0.62348968,-0.78183162)


Xk=fft(xn):   size: +7 by 1
(-0.00000018,-0.00000006)
(+7.00000000,-0.00000028)
(+0.00000021,-0.00000003)
(+0.00000018,+0.00000009)
(-0.00000003,+0.00000010)
(+0.00000005,+0.00000015)
(-0.00000028,+0.00000002)


xn-ifft(Xk):   size: +7 by 1
(+0.00000006,-0.00000000)
(+0.00000000,-0.00000006)
(-0.00000001,+0.00000000)
(+0.00000000,+0.00000006)
(+0.00000000,-0.00000003)
(-0.00000001,+0.00000000)
(+0.00000000,+0.00000000)


sn:   size: +7 by 1
+0.00000000
+0.78183150
+0.97492790
+0.43388382
-0.43388376
-0.97492790
-0.78183162

Sk=fft(sn):   size: +4 by 1
(-0.00000006,+0.00000000)
(-0.00000013,-3.50000024)
(+0.00000006,-0.00000008)
(+0.00000009,-0.00000011)


sn-ifft(Sk):   size: +7 by 1
+0.00000000
-0.00000006
+0.00000000
+0.00000000
-0.00000003
+0.00000000
+0.00000006


Process returned 0 (0x0)   execution time : 0.109 s
Press any key to continue.

你可能感兴趣的:(FFTW的C++接口)