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