KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVSignal.h
Go to the documentation of this file.
1 
4 #ifndef __KVSIGNAL_H
5 #define __KVSIGNAL_H
6 
7 #include "TGraph.h"
8 #include "TArrayF.h"
9 #include "TH1F.h"
10 
11 class KVDetector;
12 class KVDBParameterList;
13 
14 class KVSignal : public TGraph {
15 public:
16  enum SignalType {
18  kI1,
20  kQ2,
21  kI2,
22  kQ3,
24  kUNKDT
25  };
26 
27 protected:
30 
33 
35 
47 
61 
66  void ResetIndexes();
67  virtual void BuildCubicSignal(); //Interpolazione mediante cubic
68  virtual void BuildCubicSplineSignal(); //Interpolazione mediante cubic spline
69  virtual void BuildSmoothingSplineSignal(); //Interpolazione mediante cubic spline
70  void init();
71  void TreateOldSignalName();
72 
73 public:
74  KVSignal();
75  KVSignal(const char* name, const char* title);
76  KVSignal(const TString& name, const TString& title);
77  virtual ~KVSignal();
78  void Copy(TObject& obj) const;
79 
80  Bool_t IsLongEnough() const;
81 
85  const Char_t* GetDetectorName() const
86  {
87  return fDetName.Data();
88  }
89  void SetDetectorName(const Char_t* name)
90  {
91  fDetName = name;
92  }
93 
94  void SetType(const Char_t* type)
95  {
96  fType = type;
97  }
98  const Char_t* GetType() const
99  {
100  return fType.Data();
101  }
103  void DeduceFromName();
104  Int_t GetIndex() const
105  {
106  return fIndex;
107  }
108 
109  virtual Bool_t IsCharge() const
110  {
111  return kFALSE;
112  }
113  virtual Bool_t IsCurrent() const
114  {
115  return kFALSE;
116  }
117 
118  void Print(Option_t* chopt = "") const;
119 
121  void SetData(Int_t nn, Double_t* xx, Double_t* yy);
122  virtual void Set(Int_t n);
123  void SetADCData();
125  {
126  return &fAdc;
127  }
128 
130  {
131  return fFPGAOutputNumbers;
132  }
133 
134  Bool_t HasFPGA() const
135  {
136  return (GetNFPGAValues() > 0);
137  }
138 
142  Double_t GetPSAParameter(const Char_t* parname);
143  virtual void LoadPSAParameters();
144  virtual void SetDefaultValues();
145  virtual void UpdatePSAParameter(KVDBParameterList* par);
146  KVSignal* ConvertTo(const Char_t* type);
147 
148  virtual bool IsOK();
149 
153  virtual void TreateSignal();
154 #define KVSIGNAL_GETPSARESULT_KVDETECTOR 1
155  virtual void GetPSAResult(KVDetector*) const
156  {
158  }
160  {
161  return fPSAIsDone;
162  }
163 
167  void SetChannelWidth(double width)
168  {
171  }
173  {
174  return fChannelWidth;
175  }
176  Bool_t TestWidth() const;
177  void ChangeChannelWidth(Double_t newwidth);
179  void SetMaxT(double t)
180  {
181  fAdc.Set((int)(t / fChannelWidth));
182  }
184  void SetNSamples(int nn)
185  {
186  fAdc.Set(nn);
187  }
189  {
190  return fAdc.GetSize();
191  }
192 
196  void SetBaseLineLength(Int_t length, Int_t first = 0)
197  {
198  fFirstBL = first;
199  fLastBL = length - first;
200  }
202  {
203  return fFirstBL;
204  }
206  {
207  return fLastBL - fFirstBL;
208  }
209 
210  virtual Double_t ComputeBaseLine();
211  virtual Double_t ComputeDuration(Double_t th = 0.2);
212  virtual Double_t ComputeEndLine();
213  virtual void RemoveBaseLine();
214 
215  void BuildReverseTimeSignal();
216 
217  Bool_t IsFired();
218 
220  {
221  return fBaseLine;
222  }
224  {
225  return fSigmaBase;
226  }
227 
229  {
230  return fEndLine;
231  }
233  {
234  return fSigmaEnd;
235  }
236 
242  {
243  return fRiseTime;
244  }
245 
251  {
252  return fAmplitude;
253  }
254 
259  {
260  fTrapRiseTime = rise;
261  fTrapFlatTop = flat;
262  }
264  {
265  fTrapRiseTime = rise;
266  }
268  {
269  fTrapFlatTop = flat;
270  }
272  {
273  return fTrapRiseTime;
274  }
276  {
277  return fTrapFlatTop;
278  }
279  Double_t ComputeCFDThreshold(Double_t threshold = 0.5);
280 
285  {
286  fSemiGaussSigma = sig;
287  }
289  {
290  return fSemiGaussSigma;
291  }
292 
296  void SetPoleZeroCorrection(Bool_t with = kTRUE)
297  {
299  }
300  void SetTauRC(Int_t taurc)
301  {
302  fTauRC = taurc;
303  }
305  {
306  return fTauRC;
307  }
308 
312  void SetInterpolation(Bool_t with = kTRUE)
313  {
314  fWithInterpolation = with;
315  }
316  void SetInterpolatedChannelWidth(double width)
317  {
319  }
321  {
323  }
324 
328  virtual void ComputeRawAmplitude(void);
330  {
331  return fYmax - fYmin;
332  }
334  {
335  return fYmin;
336  }
338  {
339  return fYmax;
340  }
341 
344  {
346  }
348  {
350  }
351 
353  Bool_t ComputeMeanAndSigma(Int_t start, Int_t stop, Double_t& mean, Double_t& sigma);
354  Bool_t ComputeMeanAndSigma(Double_t start, Double_t stop, Double_t& mean, Double_t& sigma);
355 
360 
362  void Multiply(Double_t fact);
363  void Add(Double_t fact);
364 
366  Double_t ARC_CFD(Double_t threshold = 0.3, Double_t tdelay = 10);
367  double FindTzeroLeadingEdgeCubic(double LEVEL, int Nrecurr);
368  Double_t FindTzeroCFDCubic(double level, int Nrecurr);
369  double FindTzeroCFDCubic_rev(double level, double tend, int Nrecurr);
370  Double_t CubicInterpolation(float* data, int x2, double fmax, int Nrecurr);
371 
372  virtual void BuildCubicSignal(double taufinal); //Interpolazione mediante cubic
373  virtual void BuildCubicSplineSignal(double taufinal);//Interpolazione mediante spline cubic
374  virtual void BuildSmoothingSplineSignal(double taufinal, double l = 1, int nbits = -1); //Interpolazione mediante smoothing spline
375 
376  virtual double GetDataInter(double t);
377  virtual double GetDataInterCubic(double t);
378  virtual double GetDataCubicSpline(double t);
379  virtual double GetDataSmoothingSplineLTI(double t);//metodo che serve a BuildSmoothingSpline
380  virtual double EvalCubicSpline(double X);//metodo che serve a BuildSmoothingSpline
381 
383  void FIR_ApplyTrapezoidal(double trise, double tflat); // trise=sqrt(12)*tausha di CR-RC^4 se tflat=trise/2
384  void FIR_ApplySemigaus(double tau_usec);
385  void FIR_ApplyRCLowPass(double time_usec, int reverse = 0);
386  void FIR_ApplyRCHighPass(double time_usec, int reverse = 0);
387  void FIR_ApplyRecursiveFilter(double a0, int N, double* a, double* b, int reverse);
388  void FIR_ApplyMovingAverage(int npoints);
389  void PoleZeroSuppression(Double_t tauRC);
390  int FIR_ApplySmoothingSpline(double l, int nbits = -1);
391  double ApplyNewton(double l, double x0);//metodo che serve a FIR_ApplySmoothingSpline
392 
394  void ApplyWindowing(int window_type = 3); // 0: barlett, 1:hanning, 2:hamming, 3: blackman
395  static int FFT(unsigned int p_nSamples, bool p_bInverseTransform, double* p_lpRealIn, double* p_lpImagIn, double* p_lpRealOut, double* p_lpImagOut); // nsamples: power of 2
396  int FFT(bool p_bInverseTransform, double* p_lpRealOut, double* p_lpImagOut);
397  TH1* FFT2Histo(int output, TH1* hh = 0); // 0 modulo, 1 modulo db (normalized), 2, re, 3 im
398 
400  void ApplyModifications(TGraph* newSignal = 0, Int_t nsa = -1);
401 
402 
404  void ShiftLeft(double);//shift in ns
405  void ShiftRight(double);//shift in ns
406  void TestDraw();
407 
408  static KVSignal* MakeSignal(const char* sig_type);
409 
410  ClassDef(KVSignal, 4) //Base class for FAZIA signal processing
411 
412 };
413 
414 #endif
int Int_t
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const char Option_t
#define ClassDef(name, id)
include TDocParser_001 C image html pict1_TDocParser_001 png width
int type
To store calibration parameters in a database ,.
Base class for detector geometry description, interface to energy-loss calculations.
Definition: KVDetector.h:121
Double_t fInterpolatedChannelWidth
channel width used to produced the interpolated signal
Definition: KVSignal.h:52
void SetTauRC(Int_t taurc)
Definition: KVSignal.h:300
Double_t ComputeCFDThreshold(Double_t threshold=0.5)
calculate the time during which the signal is higher than th*fAmplitude
Definition: KVSignal.cpp:539
void FIR_ApplyMovingAverage(int npoints)
Definition: KVSignal.cpp:1731
Double_t CubicInterpolation(float *data, int x2, double fmax, int Nrecurr)
Definition: KVSignal.cpp:1163
void Add(Double_t fact)
Definition: KVSignal.cpp:1790
const Char_t * GetType() const
Definition: KVSignal.h:98
Double_t GetSemiGaussSigma() const
Definition: KVSignal.h:288
Bool_t IsLongEnough() const
Definition: KVSignal.cpp:186
TString fType
string to identify the signal type : "QH1", "I2" etc ...
Definition: KVSignal.h:32
void Copy(TObject &obj) const
Definition: KVSignal.cpp:169
void ChangeChannelWidth(Double_t newwidth)
Definition: KVSignal.cpp:474
Double_t GetShaperRiseTime() const
Definition: KVSignal.h:271
virtual void UpdatePSAParameter(KVDBParameterList *par)
Definition: KVSignal.cpp:359
virtual Double_t ComputeDuration(Double_t th=0.2)
calculate the time during which the signal is higher than th*fAmplitude
Definition: KVSignal.cpp:509
SignalType
Definition: KVSignal.h:16
@ kQH1
Definition: KVSignal.h:17
@ kADC
Definition: KVSignal.h:23
@ kUNKDT
Definition: KVSignal.h:24
@ kQL1
Definition: KVSignal.h:19
Int_t fChannel
signal type (see KVSignal::SignalType enum)
Definition: KVSignal.h:29
static int FFT(unsigned int p_nSamples, bool p_bInverseTransform, double *p_lpRealIn, double *p_lpImagIn, double *p_lpRealOut, double *p_lpImagOut)
Definition: KVSignal.cpp:962
void SetInterpolation(Bool_t with=kTRUE)
Definition: KVSignal.h:312
void ApplyModifications(TGraph *newSignal=0, Int_t nsa=-1)
apply modifications of fAdc to the original signal
Definition: KVSignal.cpp:1765
virtual void RemoveBaseLine()
Definition: KVSignal.cpp:594
Double_t fEndLine
mean value of the signal line at the end
Definition: KVSignal.h:45
Double_t fYmax
raw min/max of the signal
Definition: KVSignal.h:42
void SetShaperRiseTime(Double_t rise)
Definition: KVSignal.h:263
virtual Bool_t IsCharge() const
Definition: KVSignal.h:109
Int_t GetIndex() const
Definition: KVSignal.h:104
void FIR_ApplyRCLowPass(double time_usec, int reverse=0)
Definition: KVSignal.cpp:1632
Double_t fAmplitude
results of signal treatement
Definition: KVSignal.h:38
double FindTzeroCFDCubic_rev(double level, double tend, int Nrecurr)
Definition: KVSignal.cpp:1591
virtual void TreateSignal()
Definition: KVSignal.cpp:406
virtual Double_t ComputeBaseLine()
Definition: KVSignal.cpp:491
virtual void BuildCubicSplineSignal()
Definition: KVSignal.cpp:1351
Bool_t ComputeMeanAndSigma(Int_t start, Int_t stop, Double_t &mean, Double_t &sigma)
compute mean value and rms of a subset of samples
Definition: KVSignal.cpp:693
Double_t GetPSAParameter(const Char_t *parname)
DeduceFromName has to be called before.
Definition: KVSignal.cpp:341
virtual double EvalCubicSpline(double X)
Definition: KVSignal.cpp:1577
Double_t fSigmaBase
base line rms
Definition: KVSignal.h:44
virtual double GetDataInter(double t)
Definition: KVSignal.cpp:1201
void SetData(Int_t nn, Double_t *xx, Double_t *yy)
operation on data arrays
Definition: KVSignal.cpp:207
void SetDetectorName(const Char_t *name)
Definition: KVSignal.h:89
TArrayF * GetArray()
Definition: KVSignal.h:124
void SetInterpolatedChannelWidth(double width)
Definition: KVSignal.h:316
virtual double GetDataCubicSpline(double t)
see HSIEH S.HOU IEEE Trans. Acoustic Speech, vol. ASSP-26, NO.6, DECEMBER 1978
Definition: KVSignal.cpp:1241
void ApplyWindowing(int window_type=3)
fast fourier transform and windowing of the signal (modify only fAdc)
Definition: KVSignal.cpp:926
Int_t GetNSamples() const
Definition: KVSignal.h:188
virtual void LoadPSAParameters()
Definition: KVSignal.cpp:387
TString fDetName
name of the detector, the signal is linked to, needed to find it in the KVMultiDetector
Definition: KVSignal.h:31
Int_t fFPGAOutputNumbers
ASsociated FPGA energy outputs.
Definition: KVSignal.h:34
Double_t GetAmplitude() const
Definition: KVSignal.h:250
void SetPoleZeroCorrection(Bool_t with=kTRUE)
Definition: KVSignal.h:296
Bool_t fWithInterpolation
use of interpolation or not
Definition: KVSignal.h:59
Double_t fChannelWidth
channel width in ns
Definition: KVSignal.h:51
Int_t fFirstBL
Definition: KVSignal.h:53
Double_t GetAmplitudeTriggerValue() const
routines to manage threshold for minimum charge in the detector
Definition: KVSignal.h:343
Double_t GetBLFirst() const
Definition: KVSignal.h:201
virtual bool IsOK()
Definition: KVSignal.cpp:149
Double_t fYmin
Definition: KVSignal.h:42
static KVSignal * MakeSignal(const char *sig_type)
Create new KVSignal instance corresponding to sig_type.
Definition: KVSignal.cpp:1852
Double_t fMinimumValueForAmplitude
Minimum value to say if detector has been hitted.
Definition: KVSignal.h:60
Double_t GetRawAmplitude() const
Definition: KVSignal.h:329
void ShiftLeft(double)
---------------— OPERATORI ------------------—//
Definition: KVSignal.cpp:1800
virtual void GetPSAResult(KVDetector *) const
Definition: KVSignal.h:155
Double_t fSigmaEnd
rms value of the signal line at the end
Definition: KVSignal.h:46
void SetNSamples(int nn)
Definition: KVSignal.h:184
Double_t fBaseLine
base line mean value
Definition: KVSignal.h:43
Double_t ARC_CFD(Double_t threshold=0.3, Double_t tdelay=10)
Interpolations.
Definition: KVSignal.cpp:654
const Char_t * GetDetectorName() const
Definition: KVSignal.h:85
Double_t fChannelWidthInt
internal parameter channel width of interpolated signal in ns
Definition: KVSignal.h:65
void SetBaseLineLength(Int_t length, Int_t first=0)
Definition: KVSignal.h:196
Double_t fIMax
position of the maximum in channel
Definition: KVSignal.h:40
virtual double GetDataInterCubic(double t)
Definition: KVSignal.cpp:1220
Double_t ComputeAmplitude()
Compute and return the absolute value of the signal amplitude.
Definition: KVSignal.cpp:619
void SetTrapShaperParameters(Double_t rise, Double_t flat)
Definition: KVSignal.h:258
KVSignal()
Default constructor.
Definition: KVSignal.cpp:88
void Multiply(Double_t fact)
multiply the signal (modify only fAdc)
Definition: KVSignal.cpp:1781
Double_t GetSigmaEndLine() const
Definition: KVSignal.h:232
virtual void BuildSmoothingSplineSignal()
Definition: KVSignal.cpp:1424
void Print(Option_t *chopt="") const
Definition: KVSignal.cpp:415
void init()
Definition: KVSignal.cpp:42
void FIR_ApplyRCHighPass(double time_usec, int reverse=0)
Definition: KVSignal.cpp:1654
virtual double GetDataSmoothingSplineLTI(double t)
Definition: KVSignal.cpp:1558
void SetShaperFlatTop(Double_t flat)
Definition: KVSignal.h:267
int FIR_ApplySmoothingSpline(double l, int nbits=-1)
Definition: KVSignal.cpp:1437
virtual void Set(Int_t n)
Definition: KVSignal.cpp:195
Double_t GetBaseLine() const
Definition: KVSignal.h:219
Int_t fLastBL
first and last channel number to compute the base line
Definition: KVSignal.h:53
TArrayF fAdc
Definition: KVSignal.h:36
void TestDraw()
Definition: KVSignal.cpp:1841
Int_t GetNFPGAValues() const
Definition: KVSignal.h:129
void SetMaxT(double t)
Definition: KVSignal.h:179
Bool_t IsFired()
ComputeBaseLine and ComputeEndLine methods have to be called before.
Definition: KVSignal.cpp:581
void SetADCData()
Definition: KVSignal.cpp:229
Double_t GetEndLine() const
Definition: KVSignal.h:228
void SetType(const Char_t *type)
Definition: KVSignal.h:94
void SetChannelWidth(double width)
Definition: KVSignal.h:167
Bool_t PSAHasBeenComputed() const
Definition: KVSignal.h:159
Int_t fIndex
index deduced from block, quartet and telescope numbering
Definition: KVSignal.h:28
Double_t GetBLLength() const
Definition: KVSignal.h:205
void DeduceFromName()
Definition: KVSignal.cpp:288
Double_t GetShaperFlatTop() const
Definition: KVSignal.h:275
Double_t fTauRC
tau_rc of the electronics. Used for pole zero cancellation.
Definition: KVSignal.h:54
virtual Bool_t IsCurrent() const
Definition: KVSignal.h:113
Double_t ComputeRiseTime()
Definition: KVSignal.cpp:636
Double_t fTMax
position of the maximum in ns
Definition: KVSignal.h:41
Double_t GetInterpolatedChannelWidth() const
Definition: KVSignal.h:320
void FIR_ApplyTrapezoidal(double trise, double tflat)
different shapers (modify only fAdc)
Definition: KVSignal.cpp:839
virtual void SetDefaultValues()
To be defined in child class.
Definition: KVSignal.cpp:397
void ShiftRight(double)
Definition: KVSignal.cpp:1819
Double_t GetYmax() const
Definition: KVSignal.h:337
void FIR_ApplyRecursiveFilter(double a0, int N, double *a, double *b, int reverse)
signal will be: y[n]=a0*x[n]+sum a[k] x[k] + sum b[k] y[k]
Definition: KVSignal.cpp:1675
Double_t GetSigmaBaseLine() const
Definition: KVSignal.h:223
KVSignal * ConvertTo(const Char_t *type)
Definition: KVSignal.cpp:133
double ApplyNewton(double l, double x0)
Definition: KVSignal.cpp:1515
virtual void ComputeRawAmplitude(void)
Definition: KVSignal.cpp:440
Bool_t HasFPGA() const
Definition: KVSignal.h:134
Double_t FindTzeroCFDCubic(double level, int Nrecurr)
Definition: KVSignal.cpp:885
void FIR_ApplySemigaus(double tau_usec)
Definition: KVSignal.cpp:1613
Bool_t fWithPoleZeroCorrection
use or nor pole zero correction
Definition: KVSignal.h:58
Double_t fTrapFlatTop
flat top of the trapezoidal shaper
Definition: KVSignal.h:56
double FindTzeroLeadingEdgeCubic(double LEVEL, int Nrecurr)
Definition: KVSignal.cpp:676
Double_t fRiseTime
rise time of the signal
Definition: KVSignal.h:39
void TreateOldSignalName()
Definition: KVSignal.cpp:243
Double_t fSemiGaussSigma
sigma of the semi-gaussian shaper
Definition: KVSignal.h:57
Bool_t fPSAIsDone
indicate if PSA has been done
Definition: KVSignal.h:64
void SetSemiGaussSigma(Double_t sig)
Definition: KVSignal.h:284
TH1 * FFT2Histo(int output, TH1 *hh=0)
Definition: KVSignal.cpp:1104
void PoleZeroSuppression(Double_t tauRC)
Definition: KVSignal.cpp:1748
Double_t fTrapRiseTime
rise time of the trapezoidal shaper
Definition: KVSignal.h:55
virtual void BuildCubicSignal()
Definition: KVSignal.cpp:1386
void BuildReverseTimeSignal()
Definition: KVSignal.cpp:606
Double_t GetTauRC() const
Definition: KVSignal.h:304
Bool_t TestWidth() const
Definition: KVSignal.cpp:458
virtual ~KVSignal()
Destructor.
Definition: KVSignal.cpp:124
virtual Double_t ComputeEndLine()
Definition: KVSignal.cpp:564
Double_t GetYmin() const
Definition: KVSignal.h:333
Double_t GetRiseTime() const
Definition: KVSignal.h:241
void ResetIndexes()
Definition: KVSignal.cpp:76
Double_t GetChannelWidth() const
Definition: KVSignal.h:172
void SetAmplitudeTriggerValue(Double_t val)
Definition: KVSignal.h:347
void Set(Int_t n)
Int_t GetSize() const
const char * Data() const