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