KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVDigitalFilter.h
Go to the documentation of this file.
1 /*
2 different digital signal filters
3 classes copied from FClasses of INFN-Firenze
4 */
5 
6 
7 #include <math.h>
8 #include <stdlib.h>
9 #include <KVBase.h>
10 #include <KVSignal.h>
11 
12 #ifndef __KVDigitalFilter_H
13 #define __KVDigitalFilter_H
14 
15 class KVDigitalFilter: public KVBase {
16 public:
17  KVDigitalFilter(const double& tau = 10); //sconsigliato. usa Build*
18  KVDigitalFilter(const double& tau, int N, const double* Xcoeffs, const double* Ycoeffs);
19 
22 
23  virtual ~ KVDigitalFilter();
24  void PrintCoeffs() const;
25  void PrintCoeffsDSP() const;
26  void PrintCoeffs_AsC() const; //printed like a C array;
27 
28  virtual void Draw(Option_t* option = ""); // DRAW modulo!
31 
32 
33 
34  /******************** QUANTIZZA *************************************/
35  void Quantize(int nbits, int use_pow2 = 0, double* xgain = NULL, int* x_out = NULL, int* y_out = NULL, int* x_scale = NULL, int* y_scale = NULL);
36 
37  /************************** CREATORI ********************************/
38  static KVDigitalFilter BuildRCLowPassDeconv(const double& tau_usec, const double& tau_clk);
39  static KVDigitalFilter BuildRCLowPass(const double& tau_usec, const double& tau_clk);
40  static KVDigitalFilter BuildRCHighPass(const double& tau_usec, const double& tau_clk);
41  static KVDigitalFilter BuildRCHighPassWithPZ(const double& tau_usec, const double& preamp_decay_usec, const double& tau_clk);
42  static KVDigitalFilter BuildChebyshev(const double& freq_cutoff_mhz, int is_highpass,
43  const double& percent_ripple, int npoles,
44  const double& tau_clk);
45  static KVDigitalFilter BuildUnity(const double& tau_clk);
46  static KVDigitalFilter BuildIntegrator(const double& tau_clk);
47 
49  /************************** UTILS ***********************************/
51  static int Double2DSP(const double& val)
52  {
53  return (int)(val * 32768 * 256 + 0.5);
54  }
55  static double DSP2Double(const int val)
56  {
57  return val / (32768.*256.);
58  }
59 
60  void Compress();// shorten removing unused coeffs;
61 
63  const KVDigitalFilter& f2,
64  int parallel = 0); //default is cascade
65 
67  const KVDigitalFilter* f2,
68  int parallel = 0); //default is cascade
71  const KVDigitalFilter* f3 = NULL, const KVDigitalFilter* f4 = NULL,
72  const KVDigitalFilter* f5 = NULL, const KVDigitalFilter* f6 = NULL,
73  const KVDigitalFilter* f7 = NULL, const KVDigitalFilter* f8 = NULL,
74  const KVDigitalFilter* f9 = NULL, const KVDigitalFilter* f10 = NULL);
75 
76 
77  double ComputeGainDC(); // gain at 0 frequency
78  double ComputeGainNy(); // gain at Nyquist frequency (0.5)
79  inline void Multiply(double v)
80  {
81  for (int i = 0; i < Ncoeff; i++) a[i] *= v;
82  }
83 
84  /*------------*/
85  inline const double& GetXcoeff(int i) const
86  {
87  static const double zero = 0;
88  if (i >= 0 && i < Ncoeff) return a[i];
89  return zero;
90  }
91  inline const double& GetYcoeff(int i) const
92  {
93  static const double zero = 0;
94  if (i >= 0 && i < Ncoeff) return b[i];
95  return zero;
96  }
97  inline void SetXcoeff(int i, const double& val)
98  {
99  if (i >= 0 && i < Ncoeff)
100  a[i] = val;
101  else
102  printf("ERROR in %s: index %d out of bounds (0..%d).\n",
103  __PRETTY_FUNCTION__, i, Ncoeff - 1);
104  }
105  inline void SetYcoeff(int i, const double& val)
106  {
107  if (i >= 1 && i < Ncoeff)
108  b[i] = val;
109  else
110  printf("ERROR in %s: index %d out of bounds (0..%d).\n",
111  __PRETTY_FUNCTION__, i, Ncoeff - 1);
112  }
113  inline double* GetXarray()
114  {
115  return a;
116  }
117  inline double* GetYarray()
118  {
119  return b;
120  }
121 
122 
123  inline const double& GetTauClk()
124  {
125  return tau_clk;
126  }
127  inline int GetNCoeff()
128  {
129  return Ncoeff;
130  }
131  /************************ FILTRAGGIO ********************************/
132  void ApplyTo(double* data, const int N, int reverse = 0) const;
133  void ApplyTo(float* data, const int N, int reverse = 0) const;
134  void ApplyTo(int* data, const int N, int reverse = 0) const;
135  void FIRApplyTo(double* datax, const int NSamples, int reverse) const;
136  void FIRApplyTo(float* datax, const int NSamples, int reverse) const;
137  inline void ApplyTo(KVSignal* s, int reverse = 0) const
138  {
139  if (fabs(s->GetChannelWidth() - tau_clk) > 1e-6) {
140  printf("ERROR in %s: different tau_clk! %e != %e\n",
141  __PRETTY_FUNCTION__, s->GetChannelWidth(), tau_clk);
142  return;
143  }
144  ApplyTo(s->GetArray()->GetArray(), s->GetNSamples(), reverse);
145  }
146  inline void FIRApplyTo(KVSignal* s, int reverse = 0) const
147  {
148  if (fabs(s->GetChannelWidth() - tau_clk) > 1e-6) {
149  printf("ERROR in %s: different tau_clk! %e != %e\n",
150  __PRETTY_FUNCTION__, s->GetChannelWidth(), tau_clk);
151  return;
152  }
153  FIRApplyTo(s->GetArray()->GetArray(), s->GetNSamples(), reverse);
154  }
155  void Alloc(const int Ncoeff);
156  void Azzera();
157  int ReadMatlabFIR(char* filecoeff);
158  int WriteMatlabFIR(char* filecoeff);
159 
160 protected:
161  static void ComputeChebyshevCoeffs_serv(const double& fc, const double& lh,
162  const double& pr, const double& np,
163  int p,
164  double* a0, double* a1, double* a2,
165  double* b1, double* b2);
166 
167  static void ComputeChebyshevCoeffs(const double& freq_cutoff,
168  int is_highpass, const double& percent_ripple, int npoles,
169  double* a, double* b);
170 
171 
172 
173  /***************************** VARIABILI ********************/
174  double* a;
175  double* b;
176  int Ncoeff;
177 
178  double tau_clk;
179 
180  ClassDef(KVDigitalFilter, 1) // FIASCO: Class for digital filtering
181 };
182 
183 #endif
#define e(i)
const char Option_t
#define ClassDef(name, id)
#define N
Base class for KaliVeda framework.
Definition: KVBase.h:135
static KVDigitalFilter CombineStages(const KVDigitalFilter &f1, const KVDigitalFilter &f2, int parallel=0)
static KVDigitalFilter BuildRCHighPassWithPZ(const double &tau_usec, const double &preamp_decay_usec, const double &tau_clk)
double * GetYarray()
void Compress()
shorten filter. No deallocation of memory.
virtual void Draw(Option_t *option="")
void Multiply(double v)
void ApplyTo(KVSignal *s, int reverse=0) const
const double & GetYcoeff(int i) const
static void ComputeChebyshevCoeffs_serv(const double &fc, const double &lh, const double &pr, const double &np, int p, double *a0, double *a1, double *a2, double *b1, double *b2)
static int Double2DSP(const double &val)
– conversion to/from DSP 1.15? notation
void SetXcoeff(int i, const double &val)
static KVDigitalFilter BuildRCLowPassDeconv(const double &tau_usec, const double &tau_clk)
virtual ~ KVDigitalFilter()
void Quantize(int nbits, int use_pow2=0, double *xgain=NULL, int *x_out=NULL, int *y_out=NULL, int *x_scale=NULL, int *y_scale=NULL)
static KVDigitalFilter BuildRCLowPass(const double &tau_usec, const double &tau_clk)
static double DSP2Double(const int val)
void PrintCoeffsDSP() const
static void ComputeChebyshevCoeffs(const double &freq_cutoff, int is_highpass, const double &percent_ripple, int npoles, double *a, double *b)
static KVDigitalFilter BuildChebyshev(const double &freq_cutoff_mhz, int is_highpass, const double &percent_ripple, int npoles, const double &tau_clk)
double * GetXarray()
void ApplyTo(double *data, const int N, int reverse=0) const
double * a
coefficients.
KVDigitalFilter operator=(const KVDigitalFilter &)
int ReadMatlabFIR(char *filecoeff)
FILE *fin = fopen("notch_coeffs.txt","r");.
static KVDigitalFilter BuildUnity(const double &tau_clk)
void PrintCoeffs() const
const double & GetXcoeff(int i) const
static KVDigitalFilter BuildRCHighPass(const double &tau_usec, const double &tau_clk)
void FIRApplyTo(KVSignal *s, int reverse=0) const
void SetYcoeff(int i, const double &val)
const double & GetTauClk()
static KVDigitalFilter BuildIntegrator(const double &tau_clk)
KVDigitalFilter(const double &tau=10)
int WriteMatlabFIR(char *filecoeff)
void PrintCoeffs_AsC() const
void Alloc(const int Ncoeff)
printf("a=%p, N=%d, Ncoeff=%d\n", a, N, Ncoeff);
void FIRApplyTo(double *datax, const int NSamples, int reverse) const
static KVDigitalFilter BuildInverse(KVDigitalFilter *filter)
‍**************************************‍/
static KVDigitalFilter CombineStagesMany(const KVDigitalFilter *f1, const KVDigitalFilter *f2, const KVDigitalFilter *f3=NULL, const KVDigitalFilter *f4=NULL, const KVDigitalFilter *f5=NULL, const KVDigitalFilter *f6=NULL, const KVDigitalFilter *f7=NULL, const KVDigitalFilter *f8=NULL, const KVDigitalFilter *f9=NULL, const KVDigitalFilter *f10=NULL)
se ne devi combinare + di 1 IN CASCATA!
VecExpr< UnaryOp< Fabs< T >, SVector< T, D >, T >, T, D > fabs(const SVector< T, D > &rhs)
const long double s
Definition: KVUnits.h:94
v