X-Git-Url: http://git.kpe.io/?p=ctsim.git;a=blobdiff_plain;f=libctsim%2Ffilter.cpp;h=d22e99f37794bd69c261865c2ea00de0f30bd757;hp=6a344d42a6273957a8d847c264bcdebd438c0357;hb=3fba6928127cd65870bdcd96c8114ad5894247ae;hpb=6435258bbafabf9a4ce4445edfd97f771318eb6d diff --git a/libctsim/filter.cpp b/libctsim/filter.cpp index 6a344d4..d22e99f 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.14 2000/07/06 18:37:24 kevin Exp $ +** $Id: filter.cpp,v 1.17 2000/07/13 07:03:21 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 @@ -28,6 +28,8 @@ #include "ct.h" +int SignalFilter::N_INTEGRAL=500; //static member + /* NAME * SignalFilter::SignalFilter Construct a signal * @@ -40,11 +42,9 @@ * int nSignalPoints Number of points in signal * double param General input parameter to filters * int domain FREQUENCY or SPATIAL domain wanted - * int numint Number if intervals for calculating discrete inverse fourier xform - * 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 zeropad = 0, 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 preinterpolationFactor = 1) { m_vecFilter = NULL; m_vecFourierCosTable = NULL; @@ -70,15 +70,15 @@ SignalFilter::SignalFilter (const char* filterName, const char* filterMethodName m_failMessage += domainName; return; } - init (m_idFilter, m_idFilterMethod, bw, signalIncrement, nSignalPoints, param, m_idDomain, zeropad, numIntegral); + init (m_idFilter, m_idFilterMethod, bw, signalIncrement, nSignalPoints, param, m_idDomain, zeropad, preinterpolationFactor); } -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) +SignalFilter::SignalFilter (const FilterID filterID, const FilterMethodID filterMethodID, double bw, double signalIncrement, int nSignalPoints, double param, const DomainID domainID, int zeropad = 0, int preinterpolationFactor = 1) { - init (filterID, filterMethodID, bw, signalIncrement, nSignalPoints, param, domainID, zeropad, numIntegral); + init (filterID, filterMethodID, bw, signalIncrement, nSignalPoints, param, domainID, zeropad, preinterpolationFactor); } -SignalFilter::SignalFilter (const char* filterName, const char* domainName, double bw, double param, int numIntegral = 0) +SignalFilter::SignalFilter (const char* filterName, const char* domainName, double bw, double param) { m_bw = bw; m_nSignalPoints = 0; @@ -87,7 +87,6 @@ SignalFilter::SignalFilter (const char* filterName, const char* domainName, doub m_vecFourierCosTable = NULL; m_vecFourierSinTable = NULL; m_filterParam = param; - m_numIntegral = numIntegral; m_idFilter = convertFilterNameToID (filterName); if (m_idFilter == FILTER_INVALID) { m_fail = true; @@ -105,7 +104,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 zeropad, int numint) +SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID, double bw, double signalIncrement, int nSignalPoints, double filterParam, const DomainID domainID, int zeropad, int preinterpolationFactor) { m_bw = bw; m_idFilter = filterID; @@ -122,17 +121,29 @@ SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID m_fail = false; m_nSignalPoints = nSignalPoints; m_signalInc = signalIncrement; - m_filterParam = param; + m_filterParam = filterParam; m_zeropad = zeropad; + m_preinterpolationFactor = preinterpolationFactor; m_vecFourierCosTable = NULL; m_vecFourierSinTable = NULL; m_vecFilter = NULL; - if (m_idFilterMethod == FILTER_METHOD_FFT) - m_idFilterMethod = FILTER_METHOD_FFTW; + 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 || m_idFilterMethod == FILTER_METHOD_FFTW || m_idFilterMethod == FILTER_METHOD_RFFTW) { + 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_zeropad > 0) { double logBase2 = log(m_nSignalPoints) / log(2); @@ -143,6 +154,7 @@ SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID m_nFilterPoints = 1 << nextPowerOf2; cout << "nFilterPoints = " << m_nFilterPoints << endl; } + m_nOutputPoints = m_nFilterPoints * m_preinterpolationFactor; m_filterMin = -1. / (2 * m_signalInc); m_filterMax = 1. / (2 * m_signalInc); m_filterInc = (m_filterMax - m_filterMin) / m_nFilterPoints; @@ -156,16 +168,16 @@ SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID // 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; - } + int nFourier = max(m_nFilterPoints,m_nOutputPoints) * max(m_nFilterPoints, m_nOutputPoints) + 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 @@ -175,19 +187,21 @@ SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID } 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_realPlanBackward = rfftw_create_plan (m_nOutputPoints, FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE); m_vecRealFftInput = new fftw_real [ m_nFilterPoints ]; + m_vecRealFftSignal = new fftw_real [ m_nOutputPoints ]; 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_complexPlanForward = fftw_create_plan (m_nFilterPoints, FFTW_FORWARD, FFTW_ESTIMATE); + m_complexPlanBackward = fftw_create_plan (m_nOutputPoints, FFTW_BACKWARD, FFTW_ESTIMATE); m_vecComplexFftInput = new fftw_complex [ m_nFilterPoints ]; + m_vecComplexFftSignal = new fftw_complex [ m_nOutputPoints ]; for (int i = 0; i < m_nFilterPoints; i++) m_vecComplexFftInput[i].re = m_vecComplexFftInput[i].im = 0; + for (int i = 0; i < m_nOutputPoints; i++) + m_vecComplexFftSignal[i].re = m_vecComplexFftSignal[i].im = 0; } #endif @@ -196,7 +210,6 @@ SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID m_filterMin = -m_signalInc * (m_nSignalPoints - 1); m_filterMax = m_signalInc * (m_nSignalPoints - 1); m_filterInc = (m_filterMax - m_filterMin) / (m_nFilterPoints - 1); - m_numIntegral = numint; m_vecFilter = new double[ m_nFilterPoints ]; if (m_idFilter == FILTER_SHEPP) { @@ -212,15 +225,15 @@ SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID double x; int i; for (x = m_filterMin, i = 0; i < m_nFilterPoints; x += m_filterInc, i++) - m_vecFilter[i] = frequencyResponse (x, param); + m_vecFilter[i] = frequencyResponse (x, m_filterParam); } else if (m_idDomain == DOMAIN_SPATIAL) { double x; int i; for (x = m_filterMin, i = 0; i < m_nFilterPoints; x += m_filterInc, i++) - if (numint == 0) - m_vecFilter[i] = spatialResponseAnalytic (x, param); + if (haveAnalyticSpatial(m_idFilter)) + m_vecFilter[i] = spatialResponseAnalytic (x, m_filterParam); else - m_vecFilter[i] = spatialResponseCalc (x, param, numint); + m_vecFilter[i] = spatialResponseCalc (x, m_filterParam); } else { m_failMessage = "Illegal domain name "; m_failMessage += m_idDomain; @@ -240,11 +253,13 @@ SignalFilter::~SignalFilter (void) fftw_destroy_plan(m_complexPlanForward); fftw_destroy_plan(m_complexPlanBackward); delete [] m_vecComplexFftInput; + delete [] m_vecComplexFftSignal; } if (m_idFilterMethod == FILTER_METHOD_RFFTW) { rfftw_destroy_plan(m_realPlanForward); rfftw_destroy_plan(m_realPlanBackward); delete [] m_vecRealFftInput; + delete [] m_vecRealFftSignal; } #endif } @@ -321,10 +336,12 @@ SignalFilter::convertFilterMethodNameToID (const char* const filterMethodName) fmID = FILTER_METHOD_FOURIER_TABLE; else if (strcasecmp (filterMethodName, FILTER_METHOD_FFT_STR) == 0) fmID = FILTER_METHOD_FFT; +#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); } @@ -342,10 +359,12 @@ SignalFilter::convertFilterMethodIDToName (const FilterMethodID fmID) return (FILTER_METHOD_FOURIER_TABLE_STR); else if (fmID == FILTER_METHOD_FFT) return (FILTER_METHOD_FFT_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); } @@ -376,7 +395,6 @@ SignalFilter::convertDomainIDToName (const DomainID domain) return (name); } - void SignalFilter::filterSignal (const float input[], double output[]) const { @@ -391,10 +409,10 @@ SignalFilter::filterSignal (const float input[], double output[]) const inputSignal[i] = 0; // zeropad complex fftSignal[m_nFilterPoints]; finiteFourierTransform (inputSignal, fftSignal, m_nFilterPoints, -1); - complex filteredSignal[m_nFilterPoints]; - dotProduct (m_vecFilter, fftSignal, filteredSignal, m_nFilterPoints); + for (int i = 0; i < m_nFilterPoints; i++) + fftSignal[i] *= m_vecFilter[i]; double inverseFourier[m_nFilterPoints]; - finiteFourierTransform (filteredSignal, inverseFourier, m_nFilterPoints, 1); + finiteFourierTransform (fftSignal, inverseFourier, m_nFilterPoints, 1); for (int i = 0; i < m_nSignalPoints; i++) output[i] = inverseFourier[i]; } else if (m_idFilterMethod == FILTER_METHOD_FOURIER_TABLE) { @@ -405,10 +423,10 @@ SignalFilter::filterSignal (const float input[], double output[]) const inputSignal[i] = 0; // zeropad complex fftSignal[m_nFilterPoints]; finiteFourierTransform (inputSignal, fftSignal, -1); - complex filteredSignal[m_nFilterPoints]; - dotProduct (m_vecFilter, fftSignal, filteredSignal, m_nFilterPoints); + for (int i = 0; i < m_nFilterPoints; i++) + fftSignal[i] *= m_vecFilter[i]; double inverseFourier[m_nFilterPoints]; - finiteFourierTransform (filteredSignal, inverseFourier, 1); + finiteFourierTransform (fftSignal, inverseFourier, 1); for (int i = 0; i < m_nSignalPoints; i++) output[i] = inverseFourier[i]; } @@ -417,29 +435,31 @@ SignalFilter::filterSignal (const float input[], double output[]) const 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]; + fftw_real fftOutput [ m_nFilterPoints ]; + rfftw_one (m_realPlanForward, m_vecRealFftInput, fftOutput); + for (int i = 0; i < m_nFilterPoints; i++) + m_vecRealFftSignal[i] = m_vecFilter[i] * fftOutput[i]; + for (int i = m_nFilterPoints; i < m_nOutputPoints; i++) + m_vecRealFftSignal[i] = 0; + + fftw_real ifftOutput [ m_nOutputPoints ]; + rfftw_one(m_realPlanBackward, m_vecRealFftSignal, ifftOutput); + for (int i = 0; i < m_nSignalPoints * m_preinterpolationFactor; i++) + output[i] = ifftOutput[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); + fftw_complex fftOutput [ m_nFilterPoints ]; + fftw_one(m_complexPlanForward, m_vecComplexFftInput, fftOutput); for (int i = 0; i < m_nFilterPoints; i++) { - out[i].re *= m_vecFilter[i]; - out[i].im *= m_vecFilter[i]; + m_vecComplexFftSignal[i].re = m_vecFilter[i] * fftOutput[i].re; + m_vecComplexFftSignal[i].im = m_vecFilter[i] * fftOutput[i].im; } - fftw_complex outFiltered[m_nFilterPoints]; - fftw_one(m_complexPlanBackward, out, outFiltered); - for (int i = 0; i < m_nSignalPoints; i++) - output[i] = outFiltered[i].re; + fftw_complex ifftOutput [ m_nOutputPoints ]; + fftw_one(m_complexPlanBackward, m_vecComplexFftSignal, ifftOutput); + for (int i = 0; i < m_nSignalPoints * m_preinterpolationFactor; i++) + output[i] = ifftOutput[i].re; } #endif } @@ -450,7 +470,7 @@ SignalFilter::response (double x) double response = 0; if (m_idDomain == DOMAIN_SPATIAL) - response = spatialResponse (m_idFilter, m_bw, x, m_filterParam, m_numIntegral); + response = spatialResponse (m_idFilter, m_bw, x, m_filterParam); else if (m_idDomain == DOMAIN_FREQUENCY) response = frequencyResponse (m_idFilter, m_bw, x, m_filterParam); @@ -459,12 +479,12 @@ SignalFilter::response (double x) double -SignalFilter::spatialResponse (FilterID filterID, double bw, double x, double param, int nIntegral = 0) +SignalFilter::spatialResponse (FilterID filterID, double bw, double x, double param) { - if (nIntegral == 0) + if (haveAnalyticSpatial(filterID)) return spatialResponseAnalytic (filterID, bw, x, param); else - return spatialResponseCalc (filterID, bw, x, param, nIntegral); + return spatialResponseCalc (filterID, bw, x, param, N_INTEGRAL); } /* NAME @@ -483,9 +503,9 @@ SignalFilter::spatialResponse (FilterID filterID, double bw, double x, double pa */ double -SignalFilter::spatialResponseCalc (double x, double param, int nIntegral) const +SignalFilter::spatialResponseCalc (double x, double param) const { - return (spatialResponseCalc (m_idFilter, m_bw, x, param, nIntegral)); + return (spatialResponseCalc (m_idFilter, m_bw, x, param, N_INTEGRAL)); } double @@ -617,6 +637,30 @@ SignalFilter::spatialResponseAnalytic (double x, double param) const return spatialResponseAnalytic (m_idFilter, m_bw, x, param); } +const bool +SignalFilter::haveAnalyticSpatial (FilterID filterID) +{ + bool haveAnalytic = false; + + switch (filterID) { + case FILTER_BANDLIMIT: + case FILTER_TRIANGLE: + case FILTER_COSINE: + case FILTER_G_HAMMING: + case FILTER_ABS_BANDLIMIT: + case FILTER_ABS_COSINE: + case FILTER_ABS_G_HAMMING: + case FILTER_SHEPP: + case FILTER_SINC: + haveAnalytic = true; + break; + default: + break; + } + + return (haveAnalytic); +} + double SignalFilter::spatialResponseAnalytic (FilterID filterID, double bw, double x, double param) { @@ -930,9 +974,3 @@ SignalFilter::finiteFourierTransform (const complex input[], double outp } -void -SignalFilter::dotProduct (const double v1[], const complex v2[], complex output[], const int n) -{ - for (int i = 0; i < n; i++) - output[i] = v1[i] * v2[i]; -}