1 /*****************************************************************************
5 ** Purpose: Signal filter header file
6 ** Programmer: Kevin Rosenberg
7 ** Date Started: June 2000
9 ** This is part of the CTSim program
10 ** Copyright (C) 1983-2000 Kevin Rosenberg
12 ** $Id: filter.h,v 1.13 2000/07/11 10:32:44 kevin Exp $
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.
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.
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 ******************************************************************************/
42 FILTER_ABS_BANDLIMIT, // filter times |x|
50 FILTER_METHOD_INVALID,
51 FILTER_METHOD_CONVOLUTION,
52 FILTER_METHOD_FOURIER,
53 FILTER_METHOD_FOURIER_TABLE,
67 static const char FILTER_ABS_BANDLIMIT_STR[]= "abs_bandlimit";
68 static const char FILTER_ABS_SINC_STR[]= "abs_sinc";
69 static const char FILTER_ABS_COS_STR[]= "abs_cos";
70 static const char FILTER_ABS_HAMMING_STR[]= "abs_hamming";
71 static const char FILTER_SHEPP_STR[]= "shepp";
72 static const char FILTER_BANDLIMIT_STR[]= "bandlimit";
73 static const char FILTER_SINC_STR[]= "sinc";
74 static const char FILTER_COS_STR[]= "cos";
75 static const char FILTER_HAMMING_STR[]= "hamming";
76 static const char FILTER_TRIANGLE_STR[]= "triangle";
78 static const char FILTER_METHOD_CONVOLUTION_STR[]= "convolution";
79 static const char FILTER_METHOD_FOURIER_STR[]= "fourier";
80 static const char FILTER_METHOD_FOURIER_TABLE_STR[]="fourier_table";
81 static const char FILTER_METHOD_FFT_STR[]= "fft";
83 static const char FILTER_METHOD_FFTW_STR[]= "fftw";
84 static const char FILTER_METHOD_RFFTW_STR[]= "rfftw";
87 static const char DOMAIN_FREQUENCY_STR[]="frequency";
88 static const char DOMAIN_SPATIAL_STR[]="spatial";
91 SignalFilter (const char* filterName, const char* filterMethodName,double bw, double signalIncrement, int n, double param, const char* domainName, const int zeropad = 0, const int preinterpolationFactor = 1);
93 SignalFilter (const FilterID filt_type, FilterMethodID filterMethodID, double bw, double signalIncrement, int n, double param, const DomainID domain, const int zeropad = 0, const int preinterpolationFactor = 1);
95 SignalFilter (const char* filterName, const char* domainName, double bw, double param);
99 double* getFilter (void) const
100 { return m_vecFilter; }
102 int getNFilterPoints (void) const
103 { return m_nFilterPoints; }
105 double convolve (const double f[], const double dx, const int n, const int np) const;
107 double convolve (const float f[], const double dx, const int n, const int np) const;
109 void filterSignal (const double input[], double output[]) const;
110 void filterSignal (const float input[], double output[]) const;
112 static void finiteFourierTransform (const double input[], complex<double> output[], const int n, const int direction);
113 static void finiteFourierTransform (const complex<double> input[], complex<double> output[], const int n, const int direction);
114 static void finiteFourierTransform (const complex<double> input[], double output[], const int n, const int direction);
116 void finiteFourierTransform (const double input[], complex<double> output[], const int direction) const;
117 void finiteFourierTransform (const complex<double> input[], complex<double> output[], const int direction) const;
118 void finiteFourierTransform (const complex<double> input[], double output[], const int direction) const;
120 void setTraceLevel (int traceLevel) {m_traceLevel = traceLevel; }
122 bool fail(void) const {return m_fail;}
123 const string& failMessage(void) const {return m_failMessage;}
125 const string& nameFilter(void) const { return m_nameFilter;}
126 const string& nameDomain(void) const { return m_nameDomain;}
127 const FilterID idFilter(void) const { return m_idFilter;}
128 const DomainID idDomain(void) const { return m_idDomain;}
129 const double getFilterMin(void) const {return m_filterMin;}
130 const double getFilterMax(void) const {return m_filterMax;}
131 const double getFilterIncrement(void) const {return m_filterInc;}
133 double response (double x);
135 static double spatialResponse (FilterID fType, double bw, double x, double param);
137 static double frequencyResponse (FilterID fType, double bw, double u, double param);
139 static double spatialResponseAnalytic (FilterID fType, double bw, double x, double param);
141 static double spatialResponseCalc (FilterID fType, double bw, double x, double param, int nIntegral);
143 static void setNumIntegral(int nIntegral) {N_INTEGRAL = nIntegral;}
154 double* m_vecFourierCosTable;
155 double* m_vecFourierSinTable;
156 complex<double>* m_complexVecFilter;
158 fftw_real* m_vecRealFftInput, *m_vecRealFftSignal;
159 rfftw_plan m_realPlanForward, m_realPlanBackward;
160 fftw_complex* m_vecComplexFftInput, *m_vecComplexFftSignal;
161 fftw_plan m_complexPlanForward, m_complexPlanBackward;
165 string m_failMessage;
167 string m_nameFilterMethod;
170 FilterMethodID m_idFilterMethod;
172 double m_filterParam;
176 int m_preinterpolationFactor;
178 static int N_INTEGRAL;
180 static const bool haveAnalyticSpatial (const FilterID filterID);
181 static const FilterID convertFilterNameToID (const char* filterName);
182 static const char* convertFilterIDToName (const FilterID filterID);
183 static const FilterMethodID convertFilterMethodNameToID (const char* filterMethodName);
184 static const char* convertFilterMethodIDToName (const FilterMethodID filterMethodID);
185 static const DomainID convertDomainNameToID (const char* domainName);
186 static const char* convertDomainIDToName (const DomainID domainID);
188 void init (const FilterID filt_type, const FilterMethodID filterMethod, double bw, double signalIncrement, int n, double param, const DomainID domain, const int zeropad, const int preInterpScale);
190 double spatialResponseCalc (double x, double param) const;
192 double spatialResponseAnalytic (double x, double param) const;
194 double frequencyResponse (double u, double param) const;
196 static double sinc (double x, double mult)
197 { return (fabs(x) > F_EPSILON ? (sin (x * mult) / x) : 1.0); }
199 static double integral_abscos (double u, double w);