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