X-Git-Url: http://git.kpe.io/?a=blobdiff_plain;f=libctsim%2Ffilter.cpp;h=65163f3185b0be2de8beb7342ad643f0701ff875;hb=dfa390de2efc04d85b03718a6480f735516df0e8;hp=640bfb2d6b84783904d5e42034052ee4cabc7845;hpb=b847a34057751df3a3f4e2bf795227a13ef4788c;p=ctsim.git diff --git a/libctsim/filter.cpp b/libctsim/filter.cpp index 640bfb2..65163f3 100644 --- a/libctsim/filter.cpp +++ b/libctsim/filter.cpp @@ -9,7 +9,7 @@ ** This is part of the CTSim program ** Copyright (C) 1983-2000 Kevin Rosenberg ** -** $Id: filter.cpp,v 1.10 2000/07/05 01:34:46 kevin Exp $ +** $Id: filter.cpp,v 1.15 2000/07/07 15:30:59 kevin Exp $ ** ** This program is free software; you can redistribute it and/or modify ** it under the terms of the GNU General Public License (version 2) as @@ -44,7 +44,7 @@ * for spatial domain filters. For ANALYTIC solutions, use numint = 0 */ -SignalFilter::SignalFilter (const char* filterName, const char* filterMethodName, double bw, double signalIncrement, int nSignalPoints, double param, const char* domainName, int numIntegral = 0) +SignalFilter::SignalFilter (const char* filterName, const char* filterMethodName, double bw, double signalIncrement, int nSignalPoints, double param, const char* domainName, int zeropad = 0, int numIntegral = 0) { m_vecFilter = NULL; m_vecFourierCosTable = NULL; @@ -70,12 +70,12 @@ SignalFilter::SignalFilter (const char* filterName, const char* filterMethodName m_failMessage += domainName; return; } - init (m_idFilter, m_idFilterMethod, bw, signalIncrement, nSignalPoints, param, m_idDomain, numIntegral); + init (m_idFilter, m_idFilterMethod, bw, signalIncrement, nSignalPoints, param, m_idDomain, zeropad, numIntegral); } -SignalFilter::SignalFilter (const FilterID filterID, const FilterMethodID filterMethodID, double bw, double signalIncrement, int nSignalPoints, double param, const DomainID domainID, int numIntegral = 0) +SignalFilter::SignalFilter (const FilterID filterID, const FilterMethodID filterMethodID, double bw, double signalIncrement, int nSignalPoints, double param, const DomainID domainID, int zeropad = 0, int numIntegral = 0) { - init (filterID, filterMethodID, bw, signalIncrement, nSignalPoints, param, domainID, numIntegral); + init (filterID, filterMethodID, bw, signalIncrement, nSignalPoints, param, domainID, zeropad, numIntegral); } SignalFilter::SignalFilter (const char* filterName, const char* domainName, double bw, double param, int numIntegral = 0) @@ -105,7 +105,7 @@ SignalFilter::SignalFilter (const char* filterName, const char* domainName, doub } void -SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID, double bw, double signalIncrement, int nSignalPoints, double param, const DomainID domainID, int numint) +SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID, double bw, double signalIncrement, int nSignalPoints, double param, const DomainID domainID, int zeropad, int numint) { m_bw = bw; m_idFilter = filterID; @@ -123,37 +123,34 @@ SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID m_nSignalPoints = nSignalPoints; m_signalInc = signalIncrement; m_filterParam = param; - - if (m_idFilterMethod == FILTER_METHOD_FOURIER) { - int nFourier = m_nSignalPoints * m_nSignalPoints + 1; - double angleIncrement = (2. * PI) / m_nSignalPoints; - m_vecFourierCosTable = new double[ nFourier ]; - m_vecFourierSinTable = new double[ nFourier ]; - for (int i = 0; i < nFourier; i++) { - m_vecFourierCosTable[i] = cos (angleIncrement * i); - m_vecFourierSinTable[i] = sin (angleIncrement * i); - } - m_nFilterPoints = m_nSignalPoints; - m_filterMin = -1. / (2 * m_signalInc); - m_filterMax = 1. / (2 * m_signalInc); - m_filterInc = (m_filterMax - m_filterMin) / m_nFilterPoints; - m_vecFilter = new double [m_nFilterPoints]; - int halfFilter = m_nFilterPoints / 2; - for (int i = 0; i < halfFilter; i++) - m_vecFilter[i] = static_cast(i) / (halfFilter - 1) / (2 * m_signalInc); - for (int i = 0; i < halfFilter; i++) - m_vecFilter[m_nFilterPoints - i - 1] = static_cast(i) / (halfFilter - 1) / (2 * m_signalInc); - if (halfFilter % 2) // odd - m_vecFilter[halfFilter] = 1 / (2 * m_signalInc); - } else if (m_idFilterMethod == FILTER_METHOD_FFT || m_idFilterMethod == FILTER_METHOD_FFT_ZEROPAD_2 || m_idFilterMethod == FILTER_METHOD_FFT_ZEROPAD_4) { + m_zeropad = zeropad; + + m_vecFourierCosTable = NULL; + m_vecFourierSinTable = NULL; + m_vecFilter = NULL; + + if (m_idFilterMethod == FILTER_METHOD_FFT) { +#if HAVE_FFTW + m_idFilterMethod = FILTER_METHOD_RFFTW; +#else + m_fail = true; + m_failMessage = "FFT not yet implemented"; + return; +#endif + } + + if (m_idFilterMethod == FILTER_METHOD_FOURIER || FILTER_METHOD_FOURIER_TABLE || m_idFilterMethod == FILTER_METHOD_FFT +#if HAVE_FFTW + || m_idFilterMethod == FILTER_METHOD_FFTW || m_idFilterMethod == FILTER_METHOD_RFFTW +#endif + ) { m_nFilterPoints = m_nSignalPoints; - if (m_idFilterMethod == FILTER_METHOD_FFT_ZEROPAD_2 || m_idFilterMethod == FILTER_METHOD_FFT_ZEROPAD_4) { + if (m_zeropad > 0) { double logBase2 = log(m_nSignalPoints) / log(2); - int nextPowerOf2 = static_cast(floor(logBase2)) + 1; - if (m_idFilterMethod == FILTER_METHOD_FFT_ZEROPAD_4) - nextPowerOf2++; + int nextPowerOf2 = static_cast(floor(logBase2)); if (logBase2 != floor(logBase2)) nextPowerOf2++; + nextPowerOf2 += (m_zeropad - 1); m_nFilterPoints = 1 << nextPowerOf2; cout << "nFilterPoints = " << m_nFilterPoints << endl; } @@ -162,19 +159,49 @@ SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID m_filterInc = (m_filterMax - m_filterMin) / m_nFilterPoints; m_vecFilter = new double [m_nFilterPoints]; int halfFilter = m_nFilterPoints / 2; - for (int i = 0; i < halfFilter; i++) - m_vecFilter[i] = static_cast(i) / (halfFilter - 1) / (2 * m_signalInc); - for (int i = 0; i < halfFilter; i++) - m_vecFilter[m_nFilterPoints - i - 1] = static_cast(i) / (halfFilter - 1) / (2 * m_signalInc); - if (halfFilter % 2) // odd - m_vecFilter[halfFilter] = 1 / (2 * m_signalInc); + for (int i = 0; i <= halfFilter; i++) + m_vecFilter[i] = static_cast(i) / halfFilter/ (2. * m_signalInc); + for (int i = 1; i <= halfFilter; i++) + m_vecFilter[m_nFilterPoints - i] = static_cast(i) / halfFilter / (2. * m_signalInc); + } + + // precalculate sin and cosine tables for fourier transform + if (m_idFilterMethod == FILTER_METHOD_FOURIER_TABLE) { + int nFourier = m_nFilterPoints * m_nFilterPoints + 1; + double angleIncrement = (2. * PI) / m_nFilterPoints; + m_vecFourierCosTable = new double[ nFourier ]; + m_vecFourierSinTable = new double[ nFourier ]; + double angle = 0; + for (int i = 0; i < nFourier; i++) { + m_vecFourierCosTable[i] = cos (angle); + m_vecFourierSinTable[i] = sin (angle); + angle += angleIncrement; + } + } #if HAVE_FFTW - m_planForward = fftw_create_plan (m_nFilterPoints, FFTW_FORWARD, FFTW_ESTIMATE); - m_planBackward = fftw_create_plan (m_nFilterPoints, FFTW_BACKWARD, FFTW_ESTIMATE); -#endif + if (m_idFilterMethod == FILTER_METHOD_FFTW || m_idFilterMethod == FILTER_METHOD_RFFTW) { + for (int i = 0; i < m_nFilterPoints; i++) //fftw uses unnormalized fft + m_vecFilter[i] /= m_nFilterPoints; } + if (m_idFilterMethod == FILTER_METHOD_RFFTW) { + m_complexPlanForward = m_complexPlanBackward = NULL; + m_realPlanForward = rfftw_create_plan (m_nFilterPoints, FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE); + m_realPlanBackward = rfftw_create_plan (m_nFilterPoints, FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE); + m_vecRealFftInput = new fftw_real [ m_nFilterPoints ]; + for (int i = 0; i < m_nFilterPoints; i++) + m_vecRealFftInput[i] = 0; + } else if (m_idFilterMethod == FILTER_METHOD_FFTW) { + m_realPlanForward = m_realPlanBackward = NULL; + m_complexPlanForward = fftw_create_plan (m_nFilterPoints, FFTW_FORWARD, FFTW_ESTIMATE); + m_complexPlanBackward = fftw_create_plan (m_nFilterPoints, FFTW_BACKWARD, FFTW_ESTIMATE); + m_vecComplexFftInput = new fftw_complex [ m_nFilterPoints ]; + for (int i = 0; i < m_nFilterPoints; i++) + m_vecComplexFftInput[i].re = m_vecComplexFftInput[i].im = 0; + } +#endif + if (m_idFilterMethod == FILTER_METHOD_CONVOLUTION) { m_nFilterPoints = 2 * m_nSignalPoints - 1; m_filterMin = -m_signalInc * (m_nSignalPoints - 1); @@ -215,13 +242,20 @@ SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID SignalFilter::~SignalFilter (void) { - delete m_vecFilter; - delete m_vecFourierSinTable; - delete m_vecFourierCosTable; + delete [] m_vecFilter; + delete [] m_vecFourierSinTable; + delete [] m_vecFourierCosTable; + #if HAVE_FFTW - if (m_idFilterMethod == FILTER_METHOD_FFT) { - fftw_destroy_plan(m_planForward); - fftw_destroy_plan(m_planBackward); + if (m_idFilterMethod == FILTER_METHOD_FFTW) { + fftw_destroy_plan(m_complexPlanForward); + fftw_destroy_plan(m_complexPlanBackward); + delete [] m_vecComplexFftInput; + } + if (m_idFilterMethod == FILTER_METHOD_RFFTW) { + rfftw_destroy_plan(m_realPlanForward); + rfftw_destroy_plan(m_realPlanBackward); + delete [] m_vecRealFftInput; } #endif } @@ -294,14 +328,16 @@ SignalFilter::convertFilterMethodNameToID (const char* const filterMethodName) fmID = FILTER_METHOD_CONVOLUTION; else if (strcasecmp (filterMethodName, FILTER_METHOD_FOURIER_STR) == 0) fmID = FILTER_METHOD_FOURIER; + else if (strcasecmp (filterMethodName, FILTER_METHOD_FOURIER_TABLE_STR) == 0) + fmID = FILTER_METHOD_FOURIER_TABLE; else if (strcasecmp (filterMethodName, FILTER_METHOD_FFT_STR) == 0) fmID = FILTER_METHOD_FFT; - else if (strcasecmp (filterMethodName, FILTER_METHOD_FFT_ZEROPAD_2_STR) == 0) - fmID = FILTER_METHOD_FFT_ZEROPAD_2; - else if (strcasecmp (filterMethodName, FILTER_METHOD_FFT_ZEROPAD_4_STR) == 0) - fmID = FILTER_METHOD_FFT_ZEROPAD_4; - else if (strcasecmp (filterMethodName, FILTER_METHOD_FFT_ZEROPAD_6_STR) == 0) - fmID = FILTER_METHOD_FFT_ZEROPAD_6; +#if HAVE_FFTW + else if (strcasecmp (filterMethodName, FILTER_METHOD_FFTW_STR) == 0) + fmID = FILTER_METHOD_FFTW; + else if (strcasecmp (filterMethodName, FILTER_METHOD_RFFTW_STR) == 0) + fmID = FILTER_METHOD_RFFTW; +#endif return (fmID); } @@ -315,14 +351,16 @@ SignalFilter::convertFilterMethodIDToName (const FilterMethodID fmID) return (FILTER_METHOD_CONVOLUTION_STR); else if (fmID == FILTER_METHOD_FOURIER) return (FILTER_METHOD_FOURIER_STR); + else if (fmID == FILTER_METHOD_FOURIER_TABLE) + return (FILTER_METHOD_FOURIER_TABLE_STR); else if (fmID == FILTER_METHOD_FFT) return (FILTER_METHOD_FFT_STR); - else if (fmID == FILTER_METHOD_FFT_ZEROPAD_2) - return (FILTER_METHOD_FFT_ZEROPAD_2_STR); - else if (fmID == FILTER_METHOD_FFT_ZEROPAD_4) - return (FILTER_METHOD_FFT_ZEROPAD_4_STR); - else if (fmID == FILTER_METHOD_FFT_ZEROPAD_6) - return (FILTER_METHOD_FFT_ZEROPAD_6_STR); +#if HAVE_FFTW + else if (fmID == FILTER_METHOD_FFTW) + return (FILTER_METHOD_FFTW_STR); + else if (fmID == FILTER_METHOD_RFFTW) + return (FILTER_METHOD_RFFTW_STR); +#endif return (name); } @@ -361,67 +399,64 @@ SignalFilter::filterSignal (const float input[], double output[]) const for (int i = 0; i < m_nSignalPoints; i++) output[i] = convolve (input, m_signalInc, i, m_nSignalPoints); } else if (m_idFilterMethod == FILTER_METHOD_FOURIER) { - complex fftSignal[m_nSignalPoints]; - complex complexOutput[m_nSignalPoints]; - complex filteredSignal[m_nSignalPoints]; - finiteFourierTransform (input, fftSignal, m_nSignalPoints, -1); - if (m_traceLevel >= TRACE_PLOT) { - double test[m_nSignalPoints]; - for (int i = 0; i < m_nSignalPoints; i++) - test[i] = abs(fftSignal[i]); - ezplot_1d(test, m_nSignalPoints); - cio_kb_getc(); - } - dotProduct (m_vecFilter, fftSignal, filteredSignal, m_nSignalPoints); - if (m_traceLevel >= TRACE_PLOT) { - double test[m_nSignalPoints]; - for (int i = 0; i < m_nSignalPoints; i++) - test[i] = abs(filteredSignal[i]); - ezplot_1d(test, m_nSignalPoints); - cio_kb_getc(); - } - finiteFourierTransform (filteredSignal, complexOutput, m_nSignalPoints, 1); + double inputSignal[m_nFilterPoints]; + for (int i = 0; i < m_nSignalPoints; i++) + inputSignal[i] = input[i]; + for (int i = m_nSignalPoints; i < m_nFilterPoints; i++) + inputSignal[i] = 0; // zeropad + complex fftSignal[m_nFilterPoints]; + finiteFourierTransform (inputSignal, fftSignal, m_nFilterPoints, -1); + for (int i = 0; i < m_nFilterPoints; i++) + fftSignal[i] *= m_vecFilter[i]; + double inverseFourier[m_nFilterPoints]; + finiteFourierTransform (fftSignal, inverseFourier, m_nFilterPoints, 1); for (int i = 0; i < m_nSignalPoints; i++) - output[i] = abs( complexOutput[i] ); - } else if (m_idFilterMethod == FILTER_METHOD_FFT || FILTER_METHOD_FFT_ZEROPAD_2 || FILTER_METHOD_FFT_ZEROPAD_4) { - fftw_complex in[m_nFilterPoints], out[m_nFilterPoints]; - for (int i = 0; i < m_nSignalPoints; i++) { - in[i].re = input[i]; - in[i].im = 0; - } - for (int i = m_nSignalPoints; i < m_nFilterPoints; i++) { - in[i].re = in[i].im = 0; // ZeroPad - } - fftw_one(m_planForward, in, out); - if (m_traceLevel >= TRACE_PLOT) { - double test[m_nFilterPoints]; - for (int i = 0; i < m_nFilterPoints; i++) - test[i] = sqrt(out[i].re * out[i].re + out[i].im * out[i].im); - ezplot_1d(test, m_nFilterPoints); - cio_kb_getc(); - } - for (int i = 0; i < m_nFilterPoints; i++) { - out[i].re = m_vecFilter[i] * out[i].re / m_nSignalPoints; - out[i].im = m_vecFilter[i] * out[i].im / m_nSignalPoints; - } - if (m_traceLevel >= TRACE_PLOT) { - double test[m_nFilterPoints]; - for (int i = 0; i < m_nFilterPoints; i++) - test[i] = sqrt(out[i].re * out[i].re + out[i].im * out[i].im); - ezplot_1d(test, m_nFilterPoints); - cio_kb_getc(); - } - fftw_one(m_planBackward, out, in); - if (m_traceLevel >= TRACE_PLOT) { - double test[m_nFilterPoints]; - for (int i = 0; i < m_nFilterPoints; i++) - test[i] = sqrt(in[i].re * in[i].re + in[i].im * in[i].im); - ezplot_1d(test, m_nFilterPoints); - cio_kb_getc(); - } + output[i] = inverseFourier[i]; + } else if (m_idFilterMethod == FILTER_METHOD_FOURIER_TABLE) { + double inputSignal[m_nFilterPoints]; + for (int i = 0; i < m_nSignalPoints; i++) + inputSignal[i] = input[i]; + for (int i = m_nSignalPoints; i < m_nFilterPoints; i++) + inputSignal[i] = 0; // zeropad + complex fftSignal[m_nFilterPoints]; + finiteFourierTransform (inputSignal, fftSignal, -1); + for (int i = 0; i < m_nFilterPoints; i++) + fftSignal[i] *= m_vecFilter[i]; + double inverseFourier[m_nFilterPoints]; + finiteFourierTransform (fftSignal, inverseFourier, 1); for (int i = 0; i < m_nSignalPoints; i++) - output[i] = sqrt (in[i].re * in[i].re + in[i].im * in[i].im); + output[i] = inverseFourier[i]; + } +#if HAVE_FFTW + else if (m_idFilterMethod == FILTER_METHOD_RFFTW) { + for (int i = 0; i < m_nSignalPoints; i++) + m_vecRealFftInput[i] = input[i]; + + fftw_real out[m_nFilterPoints]; + rfftw_one (m_realPlanForward, m_vecRealFftInput, out); + for (int i = 0; i < m_nFilterPoints; i++) { + out[i] *= m_vecFilter[i]; + } + fftw_real outFiltered[m_nFilterPoints]; + rfftw_one(m_realPlanBackward, out, outFiltered); + for (int i = 0; i < m_nSignalPoints; i++) + output[i] = outFiltered[i]; + } else if (m_idFilterMethod == FILTER_METHOD_FFTW) { + for (int i = 0; i < m_nSignalPoints; i++) + m_vecComplexFftInput[i].re = input[i]; + + fftw_complex out[m_nFilterPoints]; + fftw_one(m_complexPlanForward, m_vecComplexFftInput, out); + for (int i = 0; i < m_nFilterPoints; i++) { + out[i].re *= m_vecFilter[i]; + out[i].im *= m_vecFilter[i]; + } + fftw_complex outFiltered[m_nFilterPoints]; + fftw_one(m_complexPlanBackward, out, outFiltered); + for (int i = 0; i < m_nSignalPoints; i++) + output[i] = outFiltered[i].re; } +#endif } double @@ -751,19 +786,19 @@ for (int i = 0; i < np; i++) void -SignalFilter::finiteFourierTransform (const float input[], complex output[], const int n, int direction) +SignalFilter::finiteFourierTransform (const double input[], complex output[], const int n, int direction) { if (direction < 0) direction = -1; else direction = 1; - double angleIncrement = 2 * PI / n; + double angleIncrement = direction * 2 * PI / n; for (int i = 0; i < n; i++) { double sumReal = 0; double sumImag = 0; for (int j = 0; j < n; j++) { - double angle = i * j * angleIncrement * direction; + double angle = i * j * angleIncrement; sumReal += input[j] * cos(angle); sumImag += input[j] * sin(angle); } @@ -784,11 +819,11 @@ SignalFilter::finiteFourierTransform (const complex input[], complex sum (0,0); for (int j = 0; j < n; j++) { - double angle = i * j * angleIncrement * direction; + double angle = i * j * angleIncrement; complex exponentTerm (cos(angle), sin(angle)); sum += input[j] * exponentTerm; } @@ -800,35 +835,56 @@ SignalFilter::finiteFourierTransform (const complex input[], complex output[], int direction) const +SignalFilter::finiteFourierTransform (const complex input[], double output[], const int n, int direction) +{ + if (direction < 0) + direction = -1; + else + direction = 1; + + double angleIncrement = direction * 2 * PI / n; + for (int i = 0; i < n; i++) { + double sumReal = 0; + for (int j = 0; j < n; j++) { + double angle = i * j * angleIncrement; + sumReal += input[j].real() * cos(angle) - input[j].imag() * sin(angle); + } + if (direction < 0) { + sumReal /= n; + } + output[i] = sumReal; + } +} + +void +SignalFilter::finiteFourierTransform (const double input[], complex output[], int direction) const { if (direction < 0) direction = -1; else direction = 1; - for (int i = 0; i < m_nSignalPoints; i++) { + for (int i = 0; i < m_nFilterPoints; i++) { double sumReal = 0, sumImag = 0; - for (int j = 0; j < m_nSignalPoints; j++) { + for (int j = 0; j < m_nFilterPoints; j++) { int tableIndex = i * j; if (direction > 0) { - sumReal += input[i] * m_vecFourierCosTable[tableIndex]; - sumImag += input[i] * m_vecFourierSinTable[tableIndex]; + sumReal += input[j] * m_vecFourierCosTable[tableIndex]; + sumImag += input[j] * m_vecFourierSinTable[tableIndex]; } else { - sumReal += input[i] * m_vecFourierCosTable[tableIndex]; - sumImag -= input[i] * m_vecFourierSinTable[tableIndex]; + sumReal += input[j] * m_vecFourierCosTable[tableIndex]; + sumImag -= input[j] * m_vecFourierSinTable[tableIndex]; } } if (direction < 0) { - sumReal /= m_nSignalPoints; - sumImag /= m_nSignalPoints; + sumReal /= m_nFilterPoints; + sumImag /= m_nFilterPoints; } output[i] = complex (sumReal, sumImag); } } -// (a+bi) * (c + di) = (ac - db) + (bc + da)i -#if 0 +// (a+bi) * (c + di) = (ac - bd) + (ad + bc)i void SignalFilter::finiteFourierTransform (const complex input[], complex output[], int direction) const { @@ -837,30 +893,55 @@ SignalFilter::finiteFourierTransform (const complex input[], complex 0) { - sumReal += input[i] * m_vecFourierCosTable[tableIndex]; - sumImag += input[i] * m_vecFourierSinTable[tableIndex]; + sumReal += input[j].real() * m_vecFourierCosTable[tableIndex] + - input[j].imag() * m_vecFourierSinTable[tableIndex]; + sumImag += input[j].real() * m_vecFourierSinTable[tableIndex] + + input[j].imag() * m_vecFourierCosTable[tableIndex]; } else { - sumReal += input[i] * m_vecFourierCosTable[tableIndex]; - sumImag -= input[i] * m_vecFourierSinTable[tableIndex]; + sumReal += input[j].real() * m_vecFourierCosTable[tableIndex] + - input[j].imag() * -m_vecFourierSinTable[tableIndex]; + sumImag += input[j].real() * -m_vecFourierSinTable[tableIndex] + + input[j].imag() * m_vecFourierCosTable[tableIndex]; } } - if (direction > 0) { - sumReal /= m_nSignalPoints; - sumImag /= m_nSignalPoints; + if (direction < 0) { + sumReal /= m_nFilterPoints; + sumImag /= m_nFilterPoints; } output[i] = complex (sumReal, sumImag); } } -#endif -void -SignalFilter::dotProduct (const double v1[], const complex v2[], complex output[], const int n) +void +SignalFilter::finiteFourierTransform (const complex input[], double output[], int direction) const { - for (int i = 0; i < n; i++) - output[i] = v1[i] * v2[i]; + if (direction < 0) + direction = -1; + else + direction = 1; + + for (int i = 0; i < m_nFilterPoints; i++) { + double sumReal = 0; + for (int j = 0; j < m_nFilterPoints; j++) { + int tableIndex = i * j; + if (direction > 0) { + sumReal += input[j].real() * m_vecFourierCosTable[tableIndex] + - input[j].imag() * m_vecFourierSinTable[tableIndex]; + } else { + sumReal += input[j].real() * m_vecFourierCosTable[tableIndex] + - input[j].imag() * -m_vecFourierSinTable[tableIndex]; + } + } + if (direction < 0) { + sumReal /= m_nFilterPoints; + } + output[i] = sumReal; + } } + +