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