r148: *** empty log message ***
[ctsim.git] / libctsim / filter.cpp
1 /*****************************************************************************
2 ** File IDENTIFICATION
3 ** 
4 **     Name:                   filter.cpp
5 **     Purpose:                Routines for signal-procesing filters
6 **     Progammer:              Kevin Rosenberg
7 **     Date Started:           Aug 1984
8 **
9 **  This is part of the CTSim program
10 **  Copyright (C) 1983-2000 Kevin Rosenberg
11 **
12 **  $Id: filter.cpp,v 1.18 2000/07/15 08:36:13 kevin Exp $
13 **
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.
17 **
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.
22 **
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 ******************************************************************************/
27
28 #include "ct.h"
29
30
31 int SignalFilter::N_INTEGRAL=500;  //static member
32
33 /* NAME
34  *   SignalFilter::SignalFilter     Construct a signal
35  *
36  * SYNOPSIS
37  *   f = SignalFilter (filt_type, bw, filterMin, filterMax, n, param, domain, analytic)
38  *   double f           Generated filter vector
39  *   int filt_type      Type of filter wanted
40  *   double bw          Bandwidth of filter
41  *   double filterMin, filterMax        Filter limits
42  *   int nSignalPoints  Number of points in signal
43  *   double param       General input parameter to filters
44  *   int domain         FREQUENCY or SPATIAL domain wanted
45  */
46
47 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)
48 {
49   m_vecFilter = NULL;
50   m_vecFourierCosTable = NULL;
51   m_vecFourierSinTable = NULL;
52   m_idFilter = convertFilterNameToID (filterName);
53   if (m_idFilter == FILTER_INVALID) {
54     m_fail = true;
55     m_failMessage = "Invalid Filter name ";
56     m_failMessage += filterName;
57     return;
58   }
59   m_idFilterMethod = convertFilterMethodNameToID (filterMethodName);
60   if (m_idFilterMethod == FILTER_METHOD_INVALID) {
61     m_fail = true;
62     m_failMessage = "Invalid filter method name ";
63     m_failMessage += filterMethodName;
64     return;
65   }
66   m_idDomain = convertDomainNameToID (domainName);
67   if (m_idDomain == DOMAIN_INVALID) {
68     m_fail = true;
69     m_failMessage = "Invalid domain name ";
70     m_failMessage += domainName;
71     return;
72   }
73   init (m_idFilter, m_idFilterMethod, bw, signalIncrement, nSignalPoints, param, m_idDomain, zeropad, preinterpolationFactor);
74 }
75
76 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)
77 {
78   init (filterID, filterMethodID, bw, signalIncrement, nSignalPoints, param, domainID, zeropad, preinterpolationFactor);
79 }
80
81 SignalFilter::SignalFilter (const char* filterName, const char* domainName, double bw, double param)
82 {
83   m_bw = bw;
84   m_nSignalPoints = 0;
85   m_nFilterPoints = 0;
86   m_vecFilter = NULL;
87   m_vecFourierCosTable = NULL;
88   m_vecFourierSinTable = NULL;
89   m_filterParam = param;  
90   m_idFilter = convertFilterNameToID (filterName);
91   if (m_idFilter == FILTER_INVALID) {
92     m_fail = true;
93     m_failMessage = "Invalid Filter name ";
94     m_failMessage += filterName;
95     return;
96   }
97   m_idDomain = convertDomainNameToID (domainName);
98   if (m_idDomain == DOMAIN_INVALID) {
99     m_fail = true;
100     m_failMessage = "Invalid domain name ";
101     m_failMessage += domainName;
102     return;
103   }
104 }
105
106 void
107 SignalFilter::init (const FilterID filterID, const FilterMethodID filterMethodID, double bw, double signalIncrement, int nSignalPoints, double filterParam, const DomainID domainID, int zeropad, int preinterpolationFactor)
108 {
109   m_bw = bw;
110   m_idFilter = filterID;
111   m_idDomain = domainID;
112   m_idFilterMethod = filterMethodID;
113   if (m_idFilter == FILTER_INVALID || m_idDomain == DOMAIN_INVALID || m_idFilterMethod == FILTER_METHOD_INVALID) {
114     m_fail = true;
115     return;
116   }
117   m_traceLevel = TRACE_NONE;
118   m_nameFilter = convertFilterIDToName (m_idFilter);
119   m_nameDomain = convertDomainIDToName (m_idDomain);
120   m_nameFilterMethod = convertFilterMethodIDToName (m_idFilterMethod);
121   m_fail = false;
122   m_nSignalPoints = nSignalPoints;
123   m_signalInc = signalIncrement;
124   m_filterParam = filterParam;  
125   m_zeropad = zeropad;
126   m_preinterpolationFactor = preinterpolationFactor;
127
128   m_vecFourierCosTable = NULL;
129   m_vecFourierSinTable = NULL;
130   m_vecFilter = NULL;
131
132   if (m_idFilterMethod == FILTER_METHOD_FFT) {
133 #if HAVE_FFTW
134     m_idFilterMethod = FILTER_METHOD_RFFTW;
135 #else
136     m_fail = true;
137     m_failMessage = "FFT not yet implemented";
138     return;
139 #endif
140   }
141
142   if (m_idFilterMethod == FILTER_METHOD_FOURIER || m_idFilterMethod == FILTER_METHOD_FOURIER_TABLE || m_idFilterMethod == FILTER_METHOD_FFT
143 #if HAVE_FFTW
144       || m_idFilterMethod == FILTER_METHOD_FFTW || m_idFilterMethod == FILTER_METHOD_RFFTW
145 #endif
146       ) {
147     m_nFilterPoints = m_nSignalPoints;
148     if (m_zeropad > 0) {
149       double logBase2 = log(m_nSignalPoints) / log(2);
150       int nextPowerOf2 = static_cast<int>(floor(logBase2));
151       if (logBase2 != floor(logBase2))
152         nextPowerOf2++;
153       nextPowerOf2 += (m_zeropad - 1);
154       m_nFilterPoints = 1 << nextPowerOf2;
155       if (m_traceLevel >= TRACE_TEXT)
156         cout << "nFilterPoints = " << m_nFilterPoints << endl;
157     }
158     m_nOutputPoints = m_nFilterPoints * m_preinterpolationFactor;
159     m_filterMin = -1. / (2 * m_signalInc);
160     m_filterMax = 1. / (2 * m_signalInc);
161     m_filterInc = (m_filterMax - m_filterMin) / m_nFilterPoints;
162     m_vecFilter = new double [m_nFilterPoints];
163     int halfFilter = m_nFilterPoints / 2;
164     for (int i = 0; i <= halfFilter; i++) 
165         m_vecFilter[i] = static_cast<double>(i) / halfFilter/ (2. * m_signalInc);
166     for (int i = 1; i <= halfFilter; i++)
167         m_vecFilter[m_nFilterPoints - i] = static_cast<double>(i) / halfFilter / (2. * m_signalInc);
168   }
169
170   // precalculate sin and cosine tables for fourier transform
171   if (m_idFilterMethod == FILTER_METHOD_FOURIER_TABLE) {
172       int nFourier = max(m_nFilterPoints,m_nOutputPoints) * max(m_nFilterPoints, m_nOutputPoints) + 1;
173       double angleIncrement = (2. * PI) / m_nFilterPoints;
174       m_vecFourierCosTable = new double[ nFourier ];
175       m_vecFourierSinTable = new double[ nFourier ];
176       double angle = 0;
177       for (int i = 0; i < nFourier; i++) {
178           m_vecFourierCosTable[i] = cos (angle);
179           m_vecFourierSinTable[i] = sin (angle);
180           angle += angleIncrement;
181       }
182   }
183
184 #if HAVE_FFTW
185   if (m_idFilterMethod == FILTER_METHOD_FFTW || m_idFilterMethod == FILTER_METHOD_RFFTW) {
186       for (int i = 0; i < m_nFilterPoints; i++)  //fftw uses unnormalized fft
187         m_vecFilter[i] /= m_nFilterPoints;
188   }
189
190   if (m_idFilterMethod == FILTER_METHOD_RFFTW) {
191       m_realPlanForward = rfftw_create_plan (m_nFilterPoints, FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE);
192       m_realPlanBackward = rfftw_create_plan (m_nOutputPoints, FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE);
193       m_vecRealFftInput = new fftw_real [ m_nFilterPoints ];
194       m_vecRealFftSignal = new fftw_real [ m_nOutputPoints ];
195       for (int i = 0; i < m_nFilterPoints; i++) 
196           m_vecRealFftInput[i] = 0;
197   } else if (m_idFilterMethod == FILTER_METHOD_FFTW) {
198        m_complexPlanForward = fftw_create_plan (m_nFilterPoints, FFTW_FORWARD, FFTW_ESTIMATE);
199       m_complexPlanBackward = fftw_create_plan (m_nOutputPoints, FFTW_BACKWARD, FFTW_ESTIMATE);
200       m_vecComplexFftInput = new fftw_complex [ m_nFilterPoints ];
201       m_vecComplexFftSignal = new fftw_complex [ m_nOutputPoints ];
202       for (int i = 0; i < m_nFilterPoints; i++) 
203           m_vecComplexFftInput[i].re = m_vecComplexFftInput[i].im = 0;
204       for (int i = 0; i < m_nOutputPoints; i++) 
205           m_vecComplexFftSignal[i].re = m_vecComplexFftSignal[i].im = 0;
206   }
207 #endif
208
209  if (m_idFilterMethod == FILTER_METHOD_CONVOLUTION) {
210     m_nFilterPoints = 2 * m_nSignalPoints - 1;
211     m_filterMin = -m_signalInc * (m_nSignalPoints - 1);
212     m_filterMax = m_signalInc * (m_nSignalPoints - 1);
213     m_filterInc = (m_filterMax - m_filterMin) / (m_nFilterPoints - 1);
214     m_vecFilter = new double[ m_nFilterPoints ];
215     
216     if (m_idFilter == FILTER_SHEPP) {
217       double a = 2 * m_bw;
218       double c = - 4. / (a * a);
219       int center = (m_nFilterPoints - 1) / 2;
220       int sidelen = center;
221       m_vecFilter[center] = 4. / (a * a);
222       
223       for (int i = 1; i <= sidelen; i++ )
224         m_vecFilter [center + i] = m_vecFilter [center - i] = c / (4 * (i * i) - 1);
225     } else if (m_idDomain == DOMAIN_FREQUENCY) {
226       double x;
227       int i;
228       for (x = m_filterMin, i = 0; i < m_nFilterPoints; x += m_filterInc, i++)
229         m_vecFilter[i] = frequencyResponse (x, m_filterParam);
230     } else if (m_idDomain == DOMAIN_SPATIAL) {
231       double x;
232       int i;
233       for (x = m_filterMin, i = 0; i < m_nFilterPoints; x += m_filterInc, i++)
234         if (haveAnalyticSpatial(m_idFilter))
235           m_vecFilter[i] = spatialResponseAnalytic (x, m_filterParam);
236         else
237           m_vecFilter[i] = spatialResponseCalc (x, m_filterParam);
238     } else {
239       m_failMessage = "Illegal domain name ";
240       m_failMessage += m_idDomain;
241       m_fail = true;
242     }
243   }
244 }
245
246 SignalFilter::~SignalFilter (void)
247 {
248     delete [] m_vecFilter;
249     delete [] m_vecFourierSinTable;
250     delete [] m_vecFourierCosTable;
251
252 #if HAVE_FFTW
253     if (m_idFilterMethod == FILTER_METHOD_FFTW) {
254         fftw_destroy_plan(m_complexPlanForward);
255         fftw_destroy_plan(m_complexPlanBackward);
256         delete [] m_vecComplexFftInput;
257         delete [] m_vecComplexFftSignal;
258     }
259     if (m_idFilterMethod == FILTER_METHOD_RFFTW) {
260         rfftw_destroy_plan(m_realPlanForward);
261         rfftw_destroy_plan(m_realPlanBackward);
262         delete [] m_vecRealFftInput;
263         delete [] m_vecRealFftSignal;
264     }
265 #endif
266 }
267
268
269 const SignalFilter::FilterID
270 SignalFilter::convertFilterNameToID (const char *filterName)
271 {
272   FilterID filterID = FILTER_INVALID;
273
274   if (strcasecmp (filterName, FILTER_BANDLIMIT_STR) == 0)
275     filterID = FILTER_BANDLIMIT;
276   else if (strcasecmp (filterName, FILTER_HAMMING_STR) == 0)
277     filterID = FILTER_G_HAMMING;
278   else if (strcasecmp (filterName, FILTER_SINC_STR) == 0)
279     filterID = FILTER_SINC;
280   else if (strcasecmp (filterName, FILTER_COS_STR) == 0)
281     filterID = FILTER_COSINE;
282   else if (strcasecmp (filterName, FILTER_TRIANGLE_STR) == 0)
283     filterID = FILTER_TRIANGLE;
284   else if (strcasecmp (filterName, FILTER_ABS_BANDLIMIT_STR) == 0)
285     filterID = FILTER_ABS_BANDLIMIT;
286   else if (strcasecmp (filterName, FILTER_ABS_HAMMING_STR) == 0)
287     filterID = FILTER_ABS_G_HAMMING;
288   else if (strcasecmp (filterName, FILTER_ABS_SINC_STR) == 0)
289     filterID = FILTER_ABS_SINC;
290   else if (strcasecmp (filterName, FILTER_ABS_COS_STR) == 0)
291     filterID = FILTER_ABS_COSINE;
292   else if (strcasecmp (filterName, FILTER_SHEPP_STR) == 0)
293     filterID = FILTER_SHEPP;
294
295   return (filterID);
296 }
297
298 const char *
299 SignalFilter::convertFilterIDToName (const FilterID filterID)
300 {
301   const char *name = "";
302
303   if (filterID == FILTER_SHEPP)
304     name = FILTER_SHEPP_STR;
305   else if (filterID == FILTER_ABS_COSINE)
306     name = FILTER_ABS_COS_STR;
307   else if (filterID == FILTER_ABS_SINC)
308     name = FILTER_ABS_SINC_STR;
309   else if (filterID == FILTER_ABS_G_HAMMING)
310     name = FILTER_ABS_HAMMING_STR;
311   else if (filterID == FILTER_ABS_BANDLIMIT)
312     name = FILTER_ABS_BANDLIMIT_STR;
313   else if (filterID == FILTER_COSINE)
314     name = FILTER_COS_STR;
315   else if (filterID == FILTER_SINC)
316     name = FILTER_SINC_STR;
317   else if (filterID == FILTER_G_HAMMING)
318     name = FILTER_HAMMING_STR;
319   else if (filterID == FILTER_BANDLIMIT)
320     name = FILTER_BANDLIMIT_STR;
321   else if (filterID == FILTER_TRIANGLE)
322     name = FILTER_TRIANGLE_STR;
323             
324   return (name);
325 }
326       
327 const SignalFilter::FilterMethodID
328 SignalFilter::convertFilterMethodNameToID (const char* const filterMethodName)
329 {
330   FilterMethodID fmID = FILTER_METHOD_INVALID;
331
332   if (strcasecmp (filterMethodName, FILTER_METHOD_CONVOLUTION_STR) == 0)
333     fmID = FILTER_METHOD_CONVOLUTION;
334   else if (strcasecmp (filterMethodName, FILTER_METHOD_FOURIER_STR) == 0)
335     fmID = FILTER_METHOD_FOURIER;
336   else if (strcasecmp (filterMethodName, FILTER_METHOD_FOURIER_TABLE_STR) == 0)
337     fmID = FILTER_METHOD_FOURIER_TABLE;
338   else if (strcasecmp (filterMethodName, FILTER_METHOD_FFT_STR) == 0)
339     fmID = FILTER_METHOD_FFT;
340 #if HAVE_FFTW
341   else if (strcasecmp (filterMethodName, FILTER_METHOD_FFTW_STR) == 0)
342     fmID = FILTER_METHOD_FFTW;
343   else if (strcasecmp (filterMethodName, FILTER_METHOD_RFFTW_STR) == 0)
344     fmID = FILTER_METHOD_RFFTW;
345 #endif
346
347   return (fmID);
348 }
349
350 const char *
351 SignalFilter::convertFilterMethodIDToName (const FilterMethodID fmID)
352 {
353   const char *name = "";
354
355   if (fmID == FILTER_METHOD_CONVOLUTION)
356     return (FILTER_METHOD_CONVOLUTION_STR);
357   else if (fmID == FILTER_METHOD_FOURIER)
358     return (FILTER_METHOD_FOURIER_STR);
359   else if (fmID == FILTER_METHOD_FOURIER_TABLE)
360     return (FILTER_METHOD_FOURIER_TABLE_STR);
361   else if (fmID == FILTER_METHOD_FFT)
362     return (FILTER_METHOD_FFT_STR);
363 #if HAVE_FFTW
364   else if (fmID == FILTER_METHOD_FFTW)
365     return (FILTER_METHOD_FFTW_STR);
366   else if (fmID == FILTER_METHOD_RFFTW)
367     return (FILTER_METHOD_RFFTW_STR);
368 #endif
369
370   return (name);
371 }
372
373 const SignalFilter::DomainID
374 SignalFilter::convertDomainNameToID (const char* const domainName)
375 {
376   DomainID dID = DOMAIN_INVALID;
377
378   if (strcasecmp (domainName, DOMAIN_SPATIAL_STR) == 0)
379     dID = DOMAIN_SPATIAL;
380   else if (strcasecmp (domainName, DOMAIN_FREQUENCY_STR) == 0)
381     dID = DOMAIN_FREQUENCY;
382
383   return (dID);
384 }
385
386 const char *
387 SignalFilter::convertDomainIDToName (const DomainID domain)
388 {
389   const char *name = "";
390
391   if (domain == DOMAIN_SPATIAL)
392     return (DOMAIN_SPATIAL_STR);
393   else if (domain == DOMAIN_FREQUENCY)
394     return (DOMAIN_FREQUENCY_STR);
395
396   return (name);
397 }
398
399 void
400 SignalFilter::filterSignal (const float input[], double output[]) const
401 {
402   if (m_idFilterMethod == FILTER_METHOD_CONVOLUTION) {
403     for (int i = 0; i < m_nSignalPoints; i++)
404       output[i] = convolve (input, m_signalInc, i, m_nSignalPoints);
405   } else if (m_idFilterMethod == FILTER_METHOD_FOURIER) {
406     double inputSignal[m_nFilterPoints];
407     for (int i = 0; i < m_nSignalPoints; i++)
408       inputSignal[i] = input[i];
409     for (int i = m_nSignalPoints; i < m_nFilterPoints; i++)
410       inputSignal[i] = 0;  // zeropad
411     complex<double> fftSignal[m_nFilterPoints];
412     finiteFourierTransform (inputSignal, fftSignal, m_nFilterPoints, -1);
413     for (int i = 0; i < m_nFilterPoints; i++)
414       fftSignal[i] *= m_vecFilter[i];
415     double inverseFourier[m_nFilterPoints];
416     finiteFourierTransform (fftSignal, inverseFourier, m_nFilterPoints, 1);
417     for (int i = 0; i < m_nSignalPoints; i++) 
418       output[i] = inverseFourier[i];
419   } else if (m_idFilterMethod == FILTER_METHOD_FOURIER_TABLE) {
420     double inputSignal[m_nFilterPoints];
421     for (int i = 0; i < m_nSignalPoints; i++)
422       inputSignal[i] = input[i];
423     for (int i = m_nSignalPoints; i < m_nFilterPoints; i++)
424       inputSignal[i] = 0;  // zeropad
425     complex<double> fftSignal[m_nFilterPoints];
426     finiteFourierTransform (inputSignal, fftSignal, -1);
427     for (int i = 0; i < m_nFilterPoints; i++)
428       fftSignal[i] *= m_vecFilter[i];
429     double inverseFourier[m_nFilterPoints];
430     finiteFourierTransform (fftSignal, inverseFourier, 1);
431     for (int i = 0; i < m_nSignalPoints; i++) 
432       output[i] = inverseFourier[i];
433   }
434 #if HAVE_FFTW
435   else if (m_idFilterMethod == FILTER_METHOD_RFFTW) {
436       for (int i = 0; i < m_nSignalPoints; i++)
437           m_vecRealFftInput[i] = input[i];
438
439       fftw_real fftOutput [ m_nFilterPoints ];
440       rfftw_one (m_realPlanForward, m_vecRealFftInput, fftOutput);
441       for (int i = 0; i < m_nFilterPoints; i++)
442           m_vecRealFftSignal[i] = m_vecFilter[i] * fftOutput[i];
443       for (int i = m_nFilterPoints; i < m_nOutputPoints; i++)
444         m_vecRealFftSignal[i] = 0;
445
446       fftw_real ifftOutput [ m_nOutputPoints ];
447       rfftw_one(m_realPlanBackward, m_vecRealFftSignal, ifftOutput);
448       for (int i = 0; i < m_nSignalPoints * m_preinterpolationFactor; i++)
449           output[i] = ifftOutput[i];
450   } else if (m_idFilterMethod == FILTER_METHOD_FFTW) {
451       for (int i = 0; i < m_nSignalPoints; i++)
452           m_vecComplexFftInput[i].re = input[i];
453
454       fftw_complex fftOutput [ m_nFilterPoints ];
455       fftw_one(m_complexPlanForward, m_vecComplexFftInput, fftOutput);
456       for (int i = 0; i < m_nFilterPoints; i++) {
457           m_vecComplexFftSignal[i].re = m_vecFilter[i] * fftOutput[i].re;
458           m_vecComplexFftSignal[i].im = m_vecFilter[i] * fftOutput[i].im;
459       }
460       fftw_complex ifftOutput [ m_nOutputPoints ];
461       fftw_one(m_complexPlanBackward, m_vecComplexFftSignal, ifftOutput);
462       for (int i = 0; i < m_nSignalPoints * m_preinterpolationFactor; i++) 
463           output[i] = ifftOutput[i].re;
464   }
465 #endif
466 }
467
468 double
469 SignalFilter::response (double x)
470 {
471   double response = 0;
472
473   if (m_idDomain == DOMAIN_SPATIAL)
474     response = spatialResponse (m_idFilter, m_bw, x, m_filterParam);
475   else if (m_idDomain == DOMAIN_FREQUENCY)
476     response = frequencyResponse (m_idFilter, m_bw, x, m_filterParam);
477
478   return (response);
479 }
480
481
482 double 
483 SignalFilter::spatialResponse (FilterID filterID, double bw, double x, double param)
484 {
485   if (haveAnalyticSpatial(filterID))
486     return spatialResponseAnalytic (filterID, bw, x, param);
487   else
488     return spatialResponseCalc (filterID, bw, x, param, N_INTEGRAL);
489 }
490
491 /* NAME
492  *   filter_spatial_response_calc       Calculate filter by discrete inverse fourier
493  *                                      transform of filters's frequency
494  *                                      response
495  *
496  * SYNOPSIS
497  *   y = filter_spatial_response_calc (filt_type, x, m_bw, param, n)
498  *   double y                   Filter's response in spatial domain
499  *   int filt_type              Type of filter (definitions in ct.h)
500  *   double x                   Spatial position to evaluate filter
501  *   double m_bw                        Bandwidth of window
502  *   double param               General parameter for various filters
503  *   int n                      Number of points to calculate integrations
504  */
505
506 double 
507 SignalFilter::spatialResponseCalc (double x, double param) const
508 {
509   return (spatialResponseCalc (m_idFilter, m_bw, x, param, N_INTEGRAL));
510 }
511
512 double 
513 SignalFilter::spatialResponseCalc (FilterID filterID, double bw, double x, double param, int n)
514 {
515   double zmin, zmax;
516
517   if (filterID == FILTER_TRIANGLE) {
518     zmin = 0;
519     zmax = bw;
520   } else {
521     zmin = 0;
522     zmax = bw / 2;
523   }
524   double zinc = (zmax - zmin) / (n - 1);
525
526   double z = zmin;
527   double q [n];
528   for (int i = 0; i < n; i++, z += zinc)
529     q[i] = frequencyResponse (filterID, bw, z, param) * cos (TWOPI * z * x);
530   
531   double y = 2 * integrateSimpson (zmin, zmax, q, n);
532   
533   return (y);
534 }
535
536
537 /* NAME
538  *    filter_frequency_response                 Return filter frequency response
539  *
540  * SYNOPSIS
541  *    h = filter_frequency_response (filt_type, u, m_bw, param)
542  *    double h                  Filters frequency response at u
543  *    int filt_type             Type of filter
544  *    double u                  Frequency to evaluate filter at
545  *    double m_bw                       Bandwidth of filter
546  *    double param              General input parameter for various filters
547  */
548
549 double 
550 SignalFilter::frequencyResponse (double u, double param) const
551 {
552   return frequencyResponse (m_idFilter, m_bw, u, param);
553 }
554
555
556 double 
557 SignalFilter::frequencyResponse (FilterID filterID, double bw, double u, double param)
558 {
559   double q;
560   double au = fabs (u);
561
562   switch (filterID) {
563   case FILTER_BANDLIMIT:
564     if (au >= bw / 2)
565       q = 0.;
566     else
567       q = 1;
568     break;
569   case FILTER_ABS_BANDLIMIT:
570     if (au >= bw / 2)
571       q = 0.;
572     else
573       q = au;
574     break;
575   case FILTER_TRIANGLE:
576     if (au >= bw)
577       q = 0;
578     else
579       q = 1 - au / bw;
580     break;
581   case FILTER_COSINE:
582     if (au >= bw / 2)
583       q = 0;
584     else
585       q = cos(PI * u / bw);
586     break;
587   case FILTER_ABS_COSINE:
588     if (au >= bw / 2)
589       q = 0;
590     else
591       q = au * cos(PI * u / bw);
592     break;
593   case FILTER_SINC:
594     q = bw * sinc (PI * bw * u, 1.);
595     break;
596   case FILTER_ABS_SINC:
597     q = au * bw * sinc (PI * bw * u, 1.);
598     break;
599   case FILTER_G_HAMMING:
600     if (au >= bw / 2)
601       q = 0;
602     else
603       q = param + (1 - param) * cos (TWOPI * u / bw);
604     break;
605   case FILTER_ABS_G_HAMMING:
606     if (au >= bw / 2)
607       q = 0;
608     else
609       q = au * (param + (1 - param) * cos(TWOPI * u / bw));
610     break;
611   default:
612     q = 0;
613     sys_error (ERR_WARNING, "Frequency response for filter %d not implemented [filter_frequency_response]", filterID);
614     break;
615   }
616   return (q);
617 }
618
619
620
621 /* NAME
622  *   filter_spatial_response_analytic                   Calculate filter by analytic inverse fourier
623  *                              transform of filters's frequency
624  *                              response
625  *
626  * SYNOPSIS
627  *   y = filter_spatial_response_analytic (filt_type, x, m_bw, param)
628  *   double y                   Filter's response in spatial domain
629  *   int filt_type              Type of filter (definitions in ct.h)
630  *   double x                   Spatial position to evaluate filter
631  *   double m_bw                        Bandwidth of window
632  *   double param               General parameter for various filters
633  */
634
635 double 
636 SignalFilter::spatialResponseAnalytic (double x, double param) const
637 {
638   return spatialResponseAnalytic (m_idFilter, m_bw, x, param);
639 }
640
641 const bool
642 SignalFilter::haveAnalyticSpatial (FilterID filterID)
643 {
644   bool haveAnalytic = false;
645
646   switch (filterID) {
647   case FILTER_BANDLIMIT:
648   case FILTER_TRIANGLE:
649   case FILTER_COSINE:
650   case FILTER_G_HAMMING:
651   case FILTER_ABS_BANDLIMIT:
652   case FILTER_ABS_COSINE:
653   case FILTER_ABS_G_HAMMING:
654   case FILTER_SHEPP:
655   case FILTER_SINC:
656     haveAnalytic = true;
657     break;
658   default:
659     break;
660   }
661
662   return (haveAnalytic);
663 }
664
665 double 
666 SignalFilter::spatialResponseAnalytic (FilterID filterID, double bw, double x, double param)
667 {
668   double q, temp;
669   double u = TWOPI * x;
670   double w = bw / 2;
671   double b = PI / bw;
672   double b2 = TWOPI / bw;
673
674   switch (filterID) {
675   case FILTER_BANDLIMIT:
676     q = bw * sinc(u * w, 1.0);
677     break;
678   case FILTER_TRIANGLE:
679     temp = sinc (u * w, 1.0);
680     q = bw * temp * temp;
681     break;
682   case FILTER_COSINE:
683     q = sinc(b-u,w) + sinc(b+u,w);
684     break;
685   case FILTER_G_HAMMING:
686     q = 2 * param * sin(u*w)/u + (1-param) * (sinc(b2-u, w) + sinc(b2+u, w));
687     break;
688   case FILTER_ABS_BANDLIMIT:
689     q = 2 * integral_abscos (u, w);
690     break;
691   case FILTER_ABS_COSINE:
692     q = integral_abscos(b-u,w) + integral_abscos(b+u,w);
693     break;
694   case FILTER_ABS_G_HAMMING:
695     q = 2 * param * integral_abscos(u,w) +
696       (1-param)*(integral_abscos(u-b2,w)+integral_abscos(u+b2,w));
697     break;
698   case FILTER_SHEPP:
699     if (fabs (u) < 1E-7)
700       q = 4. / (PI * bw * bw);
701     else
702       q = fabs ((2 / bw) * sin (u * w)) * sinc (u * w, 1.) * sinc (u * w, 1.);
703     break;
704   case FILTER_SINC:
705     if (fabs (x) < bw / 2)
706       q = 1.;
707     else
708       q = 0.;
709     break;
710   case FILTER_ABS_SINC:
711   default:
712     sys_error (ERR_WARNING, "Analytic filter type %d not implemented [filter_spatial_response_analytic]", filterID);
713     q = 0;
714     break;
715   }
716   
717   return (q);
718 }
719
720
721 /* NAME
722  *   sinc                       Return sin(x)/x function
723  *
724  * SYNOPSIS
725  *   v = sinc (x, mult)
726  *   double v                   sinc value
727  *   double x, mult
728  *
729  * DESCRIPTION
730  *   v = sin(x * mult) / x;
731  */
732
733
734 /* NAME
735  *   integral_abscos                    Returns integral of u*cos(u)
736  *
737  * SYNOPSIS
738  *   q = integral_abscos (u, w)
739  *   double q                   Integral value
740  *   double u                   Integration variable
741  *   double w                   Upper integration boundary
742  *
743  * DESCRIPTION
744  *   Returns the value of integral of u*cos(u)*dV for V = 0 to w
745  */
746
747 double 
748 SignalFilter::integral_abscos (double u, double w)
749 {
750   return (fabs (u) > F_EPSILON 
751      ? (cos(u * w) - 1) / (u * u) + w / u * sin (u * w) 
752      : (w * w / 2));
753 }
754
755
756 /* NAME
757  *    convolve                  Discrete convolution of two functions
758  *
759  * SYNOPSIS
760  *    r = convolve (f1, f2, dx, n, np, func_type)
761  *    double r                  Convolved result
762  *    double f1[], f2[]         Functions to be convolved
763  *    double dx                 Difference between successive x values
764  *    int n                     Array index to center convolution about
765  *    int np                    Number of points in f1 array
766  *    int func_type             EVEN or ODD or EVEN_AND_ODD function f2
767  *
768  * NOTES
769  *    f1 is the projection data, its indices range from 0 to np - 1.
770  *    The index for f2, the filter, ranges from -(np-1) to (np-1).
771  *    There are 3 ways to handle the negative vertices of f2:
772  *      1. If we know f2 is an EVEN function, then f2[-n] = f2[n].
773  *         All filters used in reconstruction are even.
774  *      2. If we know f2 is an ODD function, then f2[-n] = -f2[n] 
775  *      3. If f2 is both ODD AND EVEN, then we must store the value of f2
776  *         for negative indices.  Since f2 must range from -(np-1) to (np-1),
777  *         if we add (np - 1) to f2's array index, then f2's index will
778  *         range from 0 to 2 * (np - 1), and the origin, x = 0, will be
779  *         stored at f2[np-1].
780  */
781
782 double 
783 SignalFilter::convolve (const double func[], const double dx, const int n, const int np) const
784 {
785   double sum = 0.0;
786
787 #if UNOPTIMIZED_CONVOLUTION
788   for (int i = 0; i < np; i++)
789     sum += func[i] * m_vecFilter[n - i + (np - 1)];
790 #else
791   double* f2 = m_vecFilter + n + (np - 1);
792   for (int i = 0; i < np; i++)
793     sum += *func++ * *f2--;
794 #endif
795
796   return (sum * dx);
797 }
798
799
800 double 
801 SignalFilter::convolve (const float func[], const double dx, const int n, const int np) const
802 {
803   double sum = 0.0;
804
805 #if UNOPTIMIZED_CONVOLUTION
806 for (int i = 0; i < np; i++)
807   sum += func[i] * m_vecFilter[n - i + (np - 1)];
808 #else
809 double* f2 = m_vecFilter + n + (np - 1);
810 for (int i = 0; i < np; i++)
811   sum += *func++ * *f2--;
812 #endif
813
814   return (sum * dx);
815 }
816
817
818 void
819 SignalFilter::finiteFourierTransform (const double input[], complex<double> output[], const int n, int direction)
820 {
821   if (direction < 0)
822     direction = -1;
823   else 
824     direction = 1;
825     
826   double angleIncrement = direction * 2 * PI / n;
827   for (int i = 0; i < n; i++) {
828     double sumReal = 0;
829     double sumImag = 0;
830     for (int j = 0; j < n; j++) {
831       double angle = i * j * angleIncrement;
832       sumReal += input[j] * cos(angle);
833       sumImag += input[j] * sin(angle);
834     }
835     if (direction < 0) {
836       sumReal /= n;
837       sumImag /= n;
838     }
839     output[i] = complex<double> (sumReal, sumImag);
840   }
841 }
842
843
844 void
845 SignalFilter::finiteFourierTransform (const complex<double> input[], complex<double> output[], const int n, int direction)
846 {
847   if (direction < 0)
848     direction = -1;
849   else 
850     direction = 1;
851     
852   double angleIncrement = direction * 2 * PI / n;
853   for (int i = 0; i < n; i++) {
854     complex<double> sum (0,0);
855     for (int j = 0; j < n; j++) {
856       double angle = i * j * angleIncrement;
857       complex<double> exponentTerm (cos(angle), sin(angle));
858       sum += input[j] * exponentTerm;
859     }
860     if (direction < 0) {
861       sum /= n;
862     }
863     output[i] = sum;
864   }
865 }
866
867 void
868 SignalFilter::finiteFourierTransform (const complex<double> input[], double output[], const int n, int direction)
869 {
870   if (direction < 0)
871     direction = -1;
872   else 
873     direction = 1;
874     
875   double angleIncrement = direction * 2 * PI / n;
876   for (int i = 0; i < n; i++) {
877       double sumReal = 0;
878     for (int j = 0; j < n; j++) {
879       double angle = i * j * angleIncrement;
880       sumReal += input[j].real() * cos(angle) - input[j].imag() * sin(angle);
881     }
882     if (direction < 0) {
883       sumReal /= n;
884     }
885     output[i] = sumReal;
886   }
887 }
888
889 void
890 SignalFilter::finiteFourierTransform (const double input[], complex<double> output[], int direction) const
891 {
892   if (direction < 0)
893     direction = -1;
894   else 
895     direction = 1;
896     
897   for (int i = 0; i < m_nFilterPoints; i++) {
898     double sumReal = 0, sumImag = 0;
899     for (int j = 0; j < m_nFilterPoints; j++) {
900       int tableIndex = i * j;
901       if (direction > 0) {
902         sumReal += input[j] * m_vecFourierCosTable[tableIndex];
903         sumImag += input[j] * m_vecFourierSinTable[tableIndex];
904       } else {
905         sumReal += input[j] * m_vecFourierCosTable[tableIndex];
906         sumImag -= input[j] * m_vecFourierSinTable[tableIndex];
907       }
908     }
909     if (direction < 0) {
910       sumReal /= m_nFilterPoints;
911       sumImag /= m_nFilterPoints;
912     }
913     output[i] = complex<double> (sumReal, sumImag);
914   }
915 }
916
917 // (a+bi) * (c + di) = (ac - bd) + (ad + bc)i
918 void
919 SignalFilter::finiteFourierTransform (const complex<double> input[], complex<double> output[], int direction) const
920 {
921   if (direction < 0)
922     direction = -1;
923   else 
924     direction = 1;
925     
926   for (int i = 0; i < m_nFilterPoints; i++) {
927     double sumReal = 0, sumImag = 0;
928     for (int j = 0; j < m_nFilterPoints; j++) {
929       int tableIndex = i * j;
930       if (direction > 0) {
931         sumReal += input[j].real() * m_vecFourierCosTable[tableIndex] 
932           - input[j].imag() * m_vecFourierSinTable[tableIndex];
933         sumImag += input[j].real() * m_vecFourierSinTable[tableIndex]
934           + input[j].imag() * m_vecFourierCosTable[tableIndex];
935       } else {
936         sumReal += input[j].real() * m_vecFourierCosTable[tableIndex] 
937           - input[j].imag() * -m_vecFourierSinTable[tableIndex];
938         sumImag += input[j].real() * -m_vecFourierSinTable[tableIndex]
939           + input[j].imag() * m_vecFourierCosTable[tableIndex];
940       }
941     }
942     if (direction < 0) {
943       sumReal /= m_nFilterPoints;
944       sumImag /= m_nFilterPoints;
945     }
946     output[i] = complex<double> (sumReal, sumImag);
947   }
948 }
949
950 void
951 SignalFilter::finiteFourierTransform (const complex<double> input[], double output[], int direction) const
952 {
953   if (direction < 0)
954     direction = -1;
955   else 
956     direction = 1;
957     
958   for (int i = 0; i < m_nFilterPoints; i++) {
959       double sumReal = 0;
960     for (int j = 0; j < m_nFilterPoints; j++) {
961       int tableIndex = i * j;
962       if (direction > 0) {
963         sumReal += input[j].real() * m_vecFourierCosTable[tableIndex] 
964           - input[j].imag() * m_vecFourierSinTable[tableIndex];
965       } else {
966         sumReal += input[j].real() * m_vecFourierCosTable[tableIndex] 
967           - input[j].imag() * -m_vecFourierSinTable[tableIndex];
968       }
969     }
970     if (direction < 0) {
971       sumReal /= m_nFilterPoints;
972     }
973     output[i] = sumReal;
974   }
975 }
976
977