1 /*****************************************************************************
5 ** Purpose: Routines for signal-procesing filters
6 ** Progammer: Kevin Rosenberg
7 ** Date Started: Aug 1984
9 ** This is part of the CTSim program
10 ** Copyright (C) 1983-2000 Kevin Rosenberg
12 ** $Id: filter.cpp,v 1.20 2000/07/22 15:45:33 kevin Exp $
14 ** This program is free software; you can redistribute it and/or modify
15 ** it under the terms of the GNU General Public License (version 2) as
16 ** published by the Free Software Foundation.
18 ** This program is distributed in the hope that it will be useful,
19 ** but WITHOUT ANY WARRANTY; without even the implied warranty of
20 ** MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 ** GNU General Public License for more details.
23 ** You should have received a copy of the GNU General Public License
24 ** along with this program; if not, write to the Free Software
25 ** Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26 ******************************************************************************/
30 int SignalFilter::N_INTEGRAL=500; //static member
32 const int SignalFilter::FILTER_INVALID = -1 ;
33 const int SignalFilter::FILTER_BANDLIMIT = 0;
34 const int SignalFilter::FILTER_SINC = 1;
35 const int SignalFilter::FILTER_G_HAMMING = 2;
36 const int SignalFilter::FILTER_COSINE = 3;
37 const int SignalFilter::FILTER_TRIANGLE = 4;
38 const int SignalFilter::FILTER_ABS_BANDLIMIT = 5; // filter times |x = |
39 const int SignalFilter::FILTER_ABS_SINC = 6;
40 const int SignalFilter::FILTER_ABS_G_HAMMING = 7;
41 const int SignalFilter::FILTER_ABS_COSINE = 8;
42 const int SignalFilter::FILTER_SHEPP = 9;
44 const char* SignalFilter::s_aszFilterName[] = {
57 const char* SignalFilter::s_aszFilterTitle[] = {
63 {"Abs(w) * Bandlimit"},
70 const int SignalFilter::s_iFilterCount = sizeof(s_aszFilterName) / sizeof(const char*);
72 const int SignalFilter::FILTER_METHOD_INVALID = -1;
73 const int SignalFilter::FILTER_METHOD_CONVOLUTION = 0;
74 const int SignalFilter::FILTER_METHOD_FOURIER = 1;
75 const int SignalFilter::FILTER_METHOD_FOURIER_TABLE = 2;
76 const int SignalFilter::FILTER_METHOD_FFT = 3;
78 const int SignalFilter::FILTER_METHOD_FFTW = 4;
79 const int SignalFilter::FILTER_METHOD_RFFTW =5 ;
82 const char* SignalFilter::s_aszFilterMethodName[] = {
93 const char* SignalFilter::s_aszFilterMethodTitle[] = {
96 {"Fouier Trigometric Table Lookout"},
100 {"Real/Half-Complex FFTW"},
104 const int SignalFilter::s_iFilterMethodCount = sizeof(s_aszFilterMethodName) / sizeof(const char*);
107 const int SignalFilter::DOMAIN_INVALID = -1;
108 const int SignalFilter::DOMAIN_FREQUENCY = 0;
109 const int SignalFilter::DOMAIN_SPATIAL = 1;
111 const char* SignalFilter::s_aszDomainName[] = {
116 const char* SignalFilter::s_aszDomainTitle[] = {
121 const int SignalFilter::s_iDomainCount = sizeof(s_aszDomainName) / sizeof(const char*);
125 * SignalFilter::SignalFilter Construct a signal
128 * f = SignalFilter (filt_type, bw, filterMin, filterMax, n, param, domain, analytic)
129 * double f Generated filter vector
130 * int filt_type Type of filter wanted
131 * double bw Bandwidth of filter
132 * double filterMin, filterMax Filter limits
133 * int nSignalPoints Number of points in signal
134 * double param General input parameter to filters
135 * int domain FREQUENCY or SPATIAL domain wanted
138 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)
141 m_vecFourierCosTable = NULL;
142 m_vecFourierSinTable = NULL;
143 m_idFilter = convertFilterNameToID (filterName);
144 if (m_idFilter == FILTER_INVALID) {
146 m_failMessage = "Invalid Filter name ";
147 m_failMessage += filterName;
150 m_idFilterMethod = convertFilterMethodNameToID (filterMethodName);
151 if (m_idFilterMethod == FILTER_METHOD_INVALID) {
153 m_failMessage = "Invalid filter method name ";
154 m_failMessage += filterMethodName;
157 m_idDomain = convertDomainNameToID (domainName);
158 if (m_idDomain == DOMAIN_INVALID) {
160 m_failMessage = "Invalid domain name ";
161 m_failMessage += domainName;
164 init (m_idFilter, m_idFilterMethod, bw, signalIncrement, nSignalPoints, param, m_idDomain, zeropad, preinterpolationFactor);
167 SignalFilter::SignalFilter (const int filterID, const int filterMethodID, double bw, double signalIncrement, int nSignalPoints, double param, const int domainID, int zeropad = 0, int preinterpolationFactor = 1)
169 init (filterID, filterMethodID, bw, signalIncrement, nSignalPoints, param, domainID, zeropad, preinterpolationFactor);
172 SignalFilter::SignalFilter (const char* filterName, const char* domainName, double bw, double param)
178 m_vecFourierCosTable = NULL;
179 m_vecFourierSinTable = NULL;
180 m_filterParam = param;
181 m_idFilter = convertFilterNameToID (filterName);
182 if (m_idFilter == FILTER_INVALID) {
184 m_failMessage = "Invalid Filter name ";
185 m_failMessage += filterName;
188 m_idDomain = convertDomainNameToID (domainName);
189 if (m_idDomain == DOMAIN_INVALID) {
191 m_failMessage = "Invalid domain name ";
192 m_failMessage += domainName;
198 SignalFilter::init (const int filterID, const int filterMethodID, double bw, double signalIncrement, int nSignalPoints, double filterParam, const int domainID, int zeropad, int preinterpolationFactor)
201 m_idFilter = filterID;
202 m_idDomain = domainID;
203 m_idFilterMethod = filterMethodID;
204 if (m_idFilter == FILTER_INVALID || m_idDomain == DOMAIN_INVALID || m_idFilterMethod == FILTER_METHOD_INVALID) {
208 m_traceLevel = TRACE_NONE;
209 m_nameFilter = convertFilterIDToName (m_idFilter);
210 m_nameDomain = convertDomainIDToName (m_idDomain);
211 m_nameFilterMethod = convertFilterMethodIDToName (m_idFilterMethod);
213 m_nSignalPoints = nSignalPoints;
214 m_signalInc = signalIncrement;
215 m_filterParam = filterParam;
217 m_preinterpolationFactor = preinterpolationFactor;
219 m_vecFourierCosTable = NULL;
220 m_vecFourierSinTable = NULL;
223 if (m_idFilterMethod == FILTER_METHOD_FFT) {
225 m_idFilterMethod = FILTER_METHOD_RFFTW;
228 m_failMessage = "FFT not yet implemented";
233 if (m_idFilterMethod == FILTER_METHOD_FOURIER || m_idFilterMethod == FILTER_METHOD_FOURIER_TABLE || m_idFilterMethod == FILTER_METHOD_FFT
235 || m_idFilterMethod == FILTER_METHOD_FFTW || m_idFilterMethod == FILTER_METHOD_RFFTW
238 m_nFilterPoints = m_nSignalPoints;
240 double logBase2 = log(m_nSignalPoints) / log(2);
241 int nextPowerOf2 = static_cast<int>(floor(logBase2));
242 if (logBase2 != floor(logBase2))
244 nextPowerOf2 += (m_zeropad - 1);
245 m_nFilterPoints = 1 << nextPowerOf2;
246 if (m_traceLevel >= TRACE_TEXT)
247 cout << "nFilterPoints = " << m_nFilterPoints << endl;
249 m_nOutputPoints = m_nFilterPoints * m_preinterpolationFactor;
250 m_filterMin = -1. / (2 * m_signalInc);
251 m_filterMax = 1. / (2 * m_signalInc);
252 m_filterInc = (m_filterMax - m_filterMin) / m_nFilterPoints;
253 m_vecFilter = new double [m_nFilterPoints];
254 int halfFilter = m_nFilterPoints / 2;
255 for (int i = 0; i <= halfFilter; i++)
256 m_vecFilter[i] = static_cast<double>(i) / halfFilter/ (2. * m_signalInc);
257 for (int i = 1; i <= halfFilter; i++)
258 m_vecFilter[m_nFilterPoints - i] = static_cast<double>(i) / halfFilter / (2. * m_signalInc);
261 // precalculate sin and cosine tables for fourier transform
262 if (m_idFilterMethod == FILTER_METHOD_FOURIER_TABLE) {
263 int nFourier = max(m_nFilterPoints,m_nOutputPoints) * max(m_nFilterPoints, m_nOutputPoints) + 1;
264 double angleIncrement = (2. * PI) / m_nFilterPoints;
265 m_vecFourierCosTable = new double[ nFourier ];
266 m_vecFourierSinTable = new double[ nFourier ];
268 for (int i = 0; i < nFourier; i++) {
269 m_vecFourierCosTable[i] = cos (angle);
270 m_vecFourierSinTable[i] = sin (angle);
271 angle += angleIncrement;
276 if (m_idFilterMethod == FILTER_METHOD_FFTW || m_idFilterMethod == FILTER_METHOD_RFFTW) {
277 for (int i = 0; i < m_nFilterPoints; i++) //fftw uses unnormalized fft
278 m_vecFilter[i] /= m_nFilterPoints;
281 if (m_idFilterMethod == FILTER_METHOD_RFFTW) {
282 m_realPlanForward = rfftw_create_plan (m_nFilterPoints, FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE);
283 m_realPlanBackward = rfftw_create_plan (m_nOutputPoints, FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE);
284 m_vecRealFftInput = new fftw_real [ m_nFilterPoints ];
285 m_vecRealFftSignal = new fftw_real [ m_nOutputPoints ];
286 for (int i = 0; i < m_nFilterPoints; i++)
287 m_vecRealFftInput[i] = 0;
288 } else if (m_idFilterMethod == FILTER_METHOD_FFTW) {
289 m_complexPlanForward = fftw_create_plan (m_nFilterPoints, FFTW_FORWARD, FFTW_ESTIMATE);
290 m_complexPlanBackward = fftw_create_plan (m_nOutputPoints, FFTW_BACKWARD, FFTW_ESTIMATE);
291 m_vecComplexFftInput = new fftw_complex [ m_nFilterPoints ];
292 m_vecComplexFftSignal = new fftw_complex [ m_nOutputPoints ];
293 for (int i = 0; i < m_nFilterPoints; i++)
294 m_vecComplexFftInput[i].re = m_vecComplexFftInput[i].im = 0;
295 for (int i = 0; i < m_nOutputPoints; i++)
296 m_vecComplexFftSignal[i].re = m_vecComplexFftSignal[i].im = 0;
300 if (m_idFilterMethod == FILTER_METHOD_CONVOLUTION) {
301 m_nFilterPoints = 2 * m_nSignalPoints - 1;
302 m_filterMin = -m_signalInc * (m_nSignalPoints - 1);
303 m_filterMax = m_signalInc * (m_nSignalPoints - 1);
304 m_filterInc = (m_filterMax - m_filterMin) / (m_nFilterPoints - 1);
305 m_vecFilter = new double[ m_nFilterPoints ];
307 if (m_idFilter == FILTER_SHEPP) {
309 double c = - 4. / (a * a);
310 int center = (m_nFilterPoints - 1) / 2;
311 int sidelen = center;
312 m_vecFilter[center] = 4. / (a * a);
314 for (int i = 1; i <= sidelen; i++ )
315 m_vecFilter [center + i] = m_vecFilter [center - i] = c / (4 * (i * i) - 1);
316 } else if (m_idDomain == DOMAIN_FREQUENCY) {
319 for (x = m_filterMin, i = 0; i < m_nFilterPoints; x += m_filterInc, i++)
320 m_vecFilter[i] = frequencyResponse (x, m_filterParam);
321 } else if (m_idDomain == DOMAIN_SPATIAL) {
324 for (x = m_filterMin, i = 0; i < m_nFilterPoints; x += m_filterInc, i++)
325 if (haveAnalyticSpatial(m_idFilter))
326 m_vecFilter[i] = spatialResponseAnalytic (x, m_filterParam);
328 m_vecFilter[i] = spatialResponseCalc (x, m_filterParam);
330 m_failMessage = "Illegal domain name ";
331 m_failMessage += m_idDomain;
337 SignalFilter::~SignalFilter (void)
339 delete [] m_vecFilter;
340 delete [] m_vecFourierSinTable;
341 delete [] m_vecFourierCosTable;
344 if (m_idFilterMethod == FILTER_METHOD_FFTW) {
345 fftw_destroy_plan(m_complexPlanForward);
346 fftw_destroy_plan(m_complexPlanBackward);
347 delete [] m_vecComplexFftInput;
348 delete [] m_vecComplexFftSignal;
350 if (m_idFilterMethod == FILTER_METHOD_RFFTW) {
351 rfftw_destroy_plan(m_realPlanForward);
352 rfftw_destroy_plan(m_realPlanBackward);
353 delete [] m_vecRealFftInput;
354 delete [] m_vecRealFftSignal;
361 SignalFilter::convertFilterNameToID (const char *filterName)
363 int filterID = FILTER_INVALID;
365 for (int i = 0; i < s_iFilterCount; i++)
366 if (strcasecmp (filterName, s_aszFilterName[i]) == 0) {
375 SignalFilter::convertFilterIDToName (const int filterID)
377 static const char *name = "";
379 if (filterID >= 0 && filterID < s_iFilterCount)
380 return (s_aszFilterName [filterID]);
386 SignalFilter::convertFilterIDToTitle (const int filterID)
388 static const char *title = "";
390 if (filterID >= 0 && filterID < s_iFilterCount)
391 return (s_aszFilterTitle [filterID]);
397 SignalFilter::convertFilterMethodNameToID (const char* const filterMethodName)
399 int fmID = FILTER_METHOD_INVALID;
401 for (int i = 0; i < s_iFilterMethodCount; i++)
402 if (strcasecmp (filterMethodName, s_aszFilterMethodName[i]) == 0) {
411 SignalFilter::convertFilterMethodIDToName (const int fmID)
413 static const char *name = "";
415 if (fmID >= 0 && fmID < s_iFilterMethodCount)
416 return (s_aszFilterName [fmID]);
422 SignalFilter::convertFilterMethodIDToTitle (const int fmID)
424 static const char *title = "";
426 if (fmID >= 0 && fmID < s_iFilterMethodCount)
427 return (s_aszFilterTitle [fmID]);
433 SignalFilter::convertDomainNameToID (const char* const domainName)
435 int dID = DOMAIN_INVALID;
437 for (int i = 0; i < s_iDomainCount; i++)
438 if (strcasecmp (domainName, s_aszDomainName[i]) == 0) {
447 SignalFilter::convertDomainIDToName (const int domainID)
449 static const char *name = "";
451 if (domainID >= 0 && domainID < s_iDomainCount)
452 return (s_aszDomainName [domainID]);
458 SignalFilter::convertDomainIDToTitle (const int domainID)
460 static const char *title = "";
462 if (domainID >= 0 && domainID < s_iDomainCount)
463 return (s_aszDomainTitle [domainID]);
469 SignalFilter::filterSignal (const float input[], double output[]) const
471 if (m_idFilterMethod == FILTER_METHOD_CONVOLUTION) {
472 for (int i = 0; i < m_nSignalPoints; i++)
473 output[i] = convolve (input, m_signalInc, i, m_nSignalPoints);
474 } else if (m_idFilterMethod == FILTER_METHOD_FOURIER) {
475 double inputSignal[m_nFilterPoints];
476 for (int i = 0; i < m_nSignalPoints; i++)
477 inputSignal[i] = input[i];
478 for (int i = m_nSignalPoints; i < m_nFilterPoints; i++)
479 inputSignal[i] = 0; // zeropad
480 complex<double> fftSignal[m_nFilterPoints];
481 finiteFourierTransform (inputSignal, fftSignal, m_nFilterPoints, -1);
482 for (int i = 0; i < m_nFilterPoints; i++)
483 fftSignal[i] *= m_vecFilter[i];
484 double inverseFourier[m_nFilterPoints];
485 finiteFourierTransform (fftSignal, inverseFourier, m_nFilterPoints, 1);
486 for (int i = 0; i < m_nSignalPoints; i++)
487 output[i] = inverseFourier[i];
488 } else if (m_idFilterMethod == FILTER_METHOD_FOURIER_TABLE) {
489 double inputSignal[m_nFilterPoints];
490 for (int i = 0; i < m_nSignalPoints; i++)
491 inputSignal[i] = input[i];
492 for (int i = m_nSignalPoints; i < m_nFilterPoints; i++)
493 inputSignal[i] = 0; // zeropad
494 complex<double> fftSignal[m_nFilterPoints];
495 finiteFourierTransform (inputSignal, fftSignal, -1);
496 for (int i = 0; i < m_nFilterPoints; i++)
497 fftSignal[i] *= m_vecFilter[i];
498 double inverseFourier[m_nFilterPoints];
499 finiteFourierTransform (fftSignal, inverseFourier, 1);
500 for (int i = 0; i < m_nSignalPoints; i++)
501 output[i] = inverseFourier[i];
504 else if (m_idFilterMethod == FILTER_METHOD_RFFTW) {
505 for (int i = 0; i < m_nSignalPoints; i++)
506 m_vecRealFftInput[i] = input[i];
508 fftw_real fftOutput [ m_nFilterPoints ];
509 rfftw_one (m_realPlanForward, m_vecRealFftInput, fftOutput);
510 for (int i = 0; i < m_nFilterPoints; i++)
511 m_vecRealFftSignal[i] = m_vecFilter[i] * fftOutput[i];
512 for (int i = m_nFilterPoints; i < m_nOutputPoints; i++)
513 m_vecRealFftSignal[i] = 0;
515 fftw_real ifftOutput [ m_nOutputPoints ];
516 rfftw_one(m_realPlanBackward, m_vecRealFftSignal, ifftOutput);
517 for (int i = 0; i < m_nSignalPoints * m_preinterpolationFactor; i++)
518 output[i] = ifftOutput[i];
519 } else if (m_idFilterMethod == FILTER_METHOD_FFTW) {
520 for (int i = 0; i < m_nSignalPoints; i++)
521 m_vecComplexFftInput[i].re = input[i];
523 fftw_complex fftOutput [ m_nFilterPoints ];
524 fftw_one(m_complexPlanForward, m_vecComplexFftInput, fftOutput);
525 for (int i = 0; i < m_nFilterPoints; i++) {
526 m_vecComplexFftSignal[i].re = m_vecFilter[i] * fftOutput[i].re;
527 m_vecComplexFftSignal[i].im = m_vecFilter[i] * fftOutput[i].im;
529 fftw_complex ifftOutput [ m_nOutputPoints ];
530 fftw_one(m_complexPlanBackward, m_vecComplexFftSignal, ifftOutput);
531 for (int i = 0; i < m_nSignalPoints * m_preinterpolationFactor; i++)
532 output[i] = ifftOutput[i].re;
538 SignalFilter::response (double x)
542 if (m_idDomain == DOMAIN_SPATIAL)
543 response = spatialResponse (m_idFilter, m_bw, x, m_filterParam);
544 else if (m_idDomain == DOMAIN_FREQUENCY)
545 response = frequencyResponse (m_idFilter, m_bw, x, m_filterParam);
552 SignalFilter::spatialResponse (int filterID, double bw, double x, double param)
554 if (haveAnalyticSpatial(filterID))
555 return spatialResponseAnalytic (filterID, bw, x, param);
557 return spatialResponseCalc (filterID, bw, x, param, N_INTEGRAL);
561 * filter_spatial_response_calc Calculate filter by discrete inverse fourier
562 * transform of filters's frequency
566 * y = filter_spatial_response_calc (filt_type, x, m_bw, param, n)
567 * double y Filter's response in spatial domain
568 * int filt_type Type of filter (definitions in ct.h)
569 * double x Spatial position to evaluate filter
570 * double m_bw Bandwidth of window
571 * double param General parameter for various filters
572 * int n Number of points to calculate integrations
576 SignalFilter::spatialResponseCalc (double x, double param) const
578 return (spatialResponseCalc (m_idFilter, m_bw, x, param, N_INTEGRAL));
582 SignalFilter::spatialResponseCalc (int filterID, double bw, double x, double param, int n)
586 if (filterID == FILTER_TRIANGLE) {
593 double zinc = (zmax - zmin) / (n - 1);
597 for (int i = 0; i < n; i++, z += zinc)
598 q[i] = frequencyResponse (filterID, bw, z, param) * cos (TWOPI * z * x);
600 double y = 2 * integrateSimpson (zmin, zmax, q, n);
607 * filter_frequency_response Return filter frequency response
610 * h = filter_frequency_response (filt_type, u, m_bw, param)
611 * double h Filters frequency response at u
612 * int filt_type Type of filter
613 * double u Frequency to evaluate filter at
614 * double m_bw Bandwidth of filter
615 * double param General input parameter for various filters
619 SignalFilter::frequencyResponse (double u, double param) const
621 return frequencyResponse (m_idFilter, m_bw, u, param);
626 SignalFilter::frequencyResponse (int filterID, double bw, double u, double param)
629 double au = fabs (u);
632 case FILTER_BANDLIMIT:
638 case FILTER_ABS_BANDLIMIT:
644 case FILTER_TRIANGLE:
654 q = cos(PI * u / bw);
656 case FILTER_ABS_COSINE:
660 q = au * cos(PI * u / bw);
663 q = bw * sinc (PI * bw * u, 1.);
665 case FILTER_ABS_SINC:
666 q = au * bw * sinc (PI * bw * u, 1.);
668 case FILTER_G_HAMMING:
672 q = param + (1 - param) * cos (TWOPI * u / bw);
674 case FILTER_ABS_G_HAMMING:
678 q = au * (param + (1 - param) * cos(TWOPI * u / bw));
682 sys_error (ERR_WARNING, "Frequency response for filter %d not implemented [filter_frequency_response]", filterID);
691 * filter_spatial_response_analytic Calculate filter by analytic inverse fourier
692 * transform of filters's frequency
696 * y = filter_spatial_response_analytic (filt_type, x, m_bw, param)
697 * double y Filter's response in spatial domain
698 * int filt_type Type of filter (definitions in ct.h)
699 * double x Spatial position to evaluate filter
700 * double m_bw Bandwidth of window
701 * double param General parameter for various filters
705 SignalFilter::spatialResponseAnalytic (double x, double param) const
707 return spatialResponseAnalytic (m_idFilter, m_bw, x, param);
711 SignalFilter::haveAnalyticSpatial (int filterID)
713 bool haveAnalytic = false;
716 case FILTER_BANDLIMIT:
717 case FILTER_TRIANGLE:
719 case FILTER_G_HAMMING:
720 case FILTER_ABS_BANDLIMIT:
721 case FILTER_ABS_COSINE:
722 case FILTER_ABS_G_HAMMING:
731 return (haveAnalytic);
735 SignalFilter::spatialResponseAnalytic (int filterID, double bw, double x, double param)
738 double u = TWOPI * x;
741 double b2 = TWOPI / bw;
744 case FILTER_BANDLIMIT:
745 q = bw * sinc(u * w, 1.0);
747 case FILTER_TRIANGLE:
748 temp = sinc (u * w, 1.0);
749 q = bw * temp * temp;
752 q = sinc(b-u,w) + sinc(b+u,w);
754 case FILTER_G_HAMMING:
755 q = 2 * param * sin(u*w)/u + (1-param) * (sinc(b2-u, w) + sinc(b2+u, w));
757 case FILTER_ABS_BANDLIMIT:
758 q = 2 * integral_abscos (u, w);
760 case FILTER_ABS_COSINE:
761 q = integral_abscos(b-u,w) + integral_abscos(b+u,w);
763 case FILTER_ABS_G_HAMMING:
764 q = 2 * param * integral_abscos(u,w) +
765 (1-param)*(integral_abscos(u-b2,w)+integral_abscos(u+b2,w));
769 q = 4. / (PI * bw * bw);
771 q = fabs ((2 / bw) * sin (u * w)) * sinc (u * w, 1.) * sinc (u * w, 1.);
774 if (fabs (x) < bw / 2)
779 case FILTER_ABS_SINC:
781 sys_error (ERR_WARNING, "Analytic filter type %d not implemented [filter_spatial_response_analytic]", filterID);
791 * sinc Return sin(x)/x function
795 * double v sinc value
799 * v = sin(x * mult) / x;
804 * integral_abscos Returns integral of u*cos(u)
807 * q = integral_abscos (u, w)
808 * double q Integral value
809 * double u Integration variable
810 * double w Upper integration boundary
813 * Returns the value of integral of u*cos(u)*dV for V = 0 to w
817 SignalFilter::integral_abscos (double u, double w)
819 return (fabs (u) > F_EPSILON
820 ? (cos(u * w) - 1) / (u * u) + w / u * sin (u * w)
826 * convolve Discrete convolution of two functions
829 * r = convolve (f1, f2, dx, n, np, func_type)
830 * double r Convolved result
831 * double f1[], f2[] Functions to be convolved
832 * double dx Difference between successive x values
833 * int n Array index to center convolution about
834 * int np Number of points in f1 array
835 * int func_type EVEN or ODD or EVEN_AND_ODD function f2
838 * f1 is the projection data, its indices range from 0 to np - 1.
839 * The index for f2, the filter, ranges from -(np-1) to (np-1).
840 * There are 3 ways to handle the negative vertices of f2:
841 * 1. If we know f2 is an EVEN function, then f2[-n] = f2[n].
842 * All filters used in reconstruction are even.
843 * 2. If we know f2 is an ODD function, then f2[-n] = -f2[n]
844 * 3. If f2 is both ODD AND EVEN, then we must store the value of f2
845 * for negative indices. Since f2 must range from -(np-1) to (np-1),
846 * if we add (np - 1) to f2's array index, then f2's index will
847 * range from 0 to 2 * (np - 1), and the origin, x = 0, will be
848 * stored at f2[np-1].
852 SignalFilter::convolve (const double func[], const double dx, const int n, const int np) const
856 #if UNOPTIMIZED_CONVOLUTION
857 for (int i = 0; i < np; i++)
858 sum += func[i] * m_vecFilter[n - i + (np - 1)];
860 double* f2 = m_vecFilter + n + (np - 1);
861 for (int i = 0; i < np; i++)
862 sum += *func++ * *f2--;
870 SignalFilter::convolve (const float func[], const double dx, const int n, const int np) const
874 #if UNOPTIMIZED_CONVOLUTION
875 for (int i = 0; i < np; i++)
876 sum += func[i] * m_vecFilter[n - i + (np - 1)];
878 double* f2 = m_vecFilter + n + (np - 1);
879 for (int i = 0; i < np; i++)
880 sum += *func++ * *f2--;
888 SignalFilter::finiteFourierTransform (const double input[], complex<double> output[], const int n, int direction)
895 double angleIncrement = direction * 2 * PI / n;
896 for (int i = 0; i < n; i++) {
899 for (int j = 0; j < n; j++) {
900 double angle = i * j * angleIncrement;
901 sumReal += input[j] * cos(angle);
902 sumImag += input[j] * sin(angle);
908 output[i] = complex<double> (sumReal, sumImag);
914 SignalFilter::finiteFourierTransform (const complex<double> input[], complex<double> output[], const int n, int direction)
921 double angleIncrement = direction * 2 * PI / n;
922 for (int i = 0; i < n; i++) {
923 complex<double> sum (0,0);
924 for (int j = 0; j < n; j++) {
925 double angle = i * j * angleIncrement;
926 complex<double> exponentTerm (cos(angle), sin(angle));
927 sum += input[j] * exponentTerm;
937 SignalFilter::finiteFourierTransform (const complex<double> input[], double output[], const int n, int direction)
944 double angleIncrement = direction * 2 * PI / n;
945 for (int i = 0; i < n; i++) {
947 for (int j = 0; j < n; j++) {
948 double angle = i * j * angleIncrement;
949 sumReal += input[j].real() * cos(angle) - input[j].imag() * sin(angle);
959 SignalFilter::finiteFourierTransform (const double input[], complex<double> output[], int direction) const
966 for (int i = 0; i < m_nFilterPoints; i++) {
967 double sumReal = 0, sumImag = 0;
968 for (int j = 0; j < m_nFilterPoints; j++) {
969 int tableIndex = i * j;
971 sumReal += input[j] * m_vecFourierCosTable[tableIndex];
972 sumImag += input[j] * m_vecFourierSinTable[tableIndex];
974 sumReal += input[j] * m_vecFourierCosTable[tableIndex];
975 sumImag -= input[j] * m_vecFourierSinTable[tableIndex];
979 sumReal /= m_nFilterPoints;
980 sumImag /= m_nFilterPoints;
982 output[i] = complex<double> (sumReal, sumImag);
986 // (a+bi) * (c + di) = (ac - bd) + (ad + bc)i
988 SignalFilter::finiteFourierTransform (const complex<double> input[], complex<double> output[], int direction) const
995 for (int i = 0; i < m_nFilterPoints; i++) {
996 double sumReal = 0, sumImag = 0;
997 for (int j = 0; j < m_nFilterPoints; j++) {
998 int tableIndex = i * j;
1000 sumReal += input[j].real() * m_vecFourierCosTable[tableIndex]
1001 - input[j].imag() * m_vecFourierSinTable[tableIndex];
1002 sumImag += input[j].real() * m_vecFourierSinTable[tableIndex]
1003 + input[j].imag() * m_vecFourierCosTable[tableIndex];
1005 sumReal += input[j].real() * m_vecFourierCosTable[tableIndex]
1006 - input[j].imag() * -m_vecFourierSinTable[tableIndex];
1007 sumImag += input[j].real() * -m_vecFourierSinTable[tableIndex]
1008 + input[j].imag() * m_vecFourierCosTable[tableIndex];
1011 if (direction < 0) {
1012 sumReal /= m_nFilterPoints;
1013 sumImag /= m_nFilterPoints;
1015 output[i] = complex<double> (sumReal, sumImag);
1020 SignalFilter::finiteFourierTransform (const complex<double> input[], double output[], int direction) const
1027 for (int i = 0; i < m_nFilterPoints; i++) {
1029 for (int j = 0; j < m_nFilterPoints; j++) {
1030 int tableIndex = i * j;
1031 if (direction > 0) {
1032 sumReal += input[j].real() * m_vecFourierCosTable[tableIndex]
1033 - input[j].imag() * m_vecFourierSinTable[tableIndex];
1035 sumReal += input[j].real() * m_vecFourierCosTable[tableIndex]
1036 - input[j].imag() * -m_vecFourierSinTable[tableIndex];
1039 if (direction < 0) {
1040 sumReal /= m_nFilterPoints;
1042 output[i] = sumReal;