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