KaliVeda
Toolkit for HIC analysis
KVSignal.h
1 
4 #ifndef __KVSIGNAL_H
5 #define __KVSIGNAL_H
6 
7 #include "TGraph.h"
8 #include "TArrayF.h"
9 #include "TH1F.h"
10 #include "KVNameValueList.h"
11 
12 class KVSignal : public TGraph {
13 
14 protected:
26 
40 
42 
47 
48  virtual void BuildCubicSignal(); //Interpolazione mediante cubic
49  virtual void BuildCubicSplineSignal(); //Interpolazione mediante cubic spline
50  virtual void BuildSmoothingSplineSignal(); //Interpolazione mediante cubic spline
51  void init();
52 
53  void print_psa_parameters() const;
54 
55  public:
56  KVSignal();
57  KVSignal(const char* name, const char* title);
58  KVSignal(const TString& name, const TString& title);
59 
60  Bool_t IsLongEnough() const;
61 
62  virtual Bool_t IsCharge() const
63  {
64  return kFALSE;
65  }
66  virtual Bool_t IsCurrent() const
67  {
68  return kFALSE;
69  }
73  virtual const Char_t* GetDetectorName() const
74  {
75  return "";
76  }
77  void SetType(const Char_t* type)
78  {
79  fType = type;
80  SetTitle(type);
81  }
82  const Char_t* GetType() const
83  {
84  return fType.Data();
85  }
86  virtual Int_t GetIndex() const
87  {
88  return -1;
89  }
90  virtual Int_t GetNFPGAValues() const
91  {
92  return -1;
93  }
94  virtual Bool_t HasFPGA() const
95  {
96  return false;
97  }
98 
99  void Print(Option_t* chopt = "") const override;
100 
102  void SetData(Int_t nn, Double_t* xx, Double_t* yy);
103  void Set(Int_t n) override;
104  void SetADCData();
106  {
107  return &fAdc;
108  }
109 
113  virtual Double_t GetPSAParameter(const Char_t* parname)
114  {
115  return 0;
116  }
117  virtual void LoadPSAParameters();
118  virtual void SetDefaultValues();
119  virtual void UpdatePSAParameter(KVNameValueList* par);
120  KVSignal* ConvertTo(const Char_t* type);
121 
122  bool IsOK();
123 
127  virtual void TreateSignal();
128  virtual void GetPSAResults(KVNameValueList&) const {}
130  {
131  return fPSAIsDone;
132  }
133 
137  void SetChannelWidth(double width)
138  {
141  }
143  {
144  return fChannelWidth;
145  }
146  Bool_t TestWidth() const;
147  void ChangeChannelWidth(Double_t newwidth);
149  void SetMaxT(double t)
150  {
151  fAdc.Set((int)(t / fChannelWidth));
152  }
154  void SetNSamples(int nn)
155  {
156  fAdc.Set(nn);
157  }
159  {
160  return fAdc.GetSize();
161  }
162 
166  void SetBaseLineLength(Int_t length, Int_t first = 0)
167  {
168  fFirstBL = first;
169  fLastBL = length - first;
170  }
172  {
173  return fFirstBL;
174  }
176  {
177  return fLastBL - fFirstBL;
178  }
179 
180  virtual Double_t ComputeBaseLine();
181  virtual Double_t ComputeDuration(Double_t th = 0.2);
182  virtual Double_t ComputeEndLine();
183  virtual void RemoveBaseLine();
184 
185  void BuildReverseTimeSignal();
186 
187  Bool_t IsFired();
188 
190  {
191  return fBaseLine;
192  }
194  {
195  return fSigmaBase;
196  }
197 
199  {
200  return fEndLine;
201  }
203  {
204  return fSigmaEnd;
205  }
206 
212  {
213  return fRiseTime;
214  }
215 
221  {
222  return fAmplitude;
223  }
224 
229  {
230  fTrapRiseTime = rise;
231  fTrapFlatTop = flat;
232  }
234  {
235  fTrapRiseTime = rise;
236  }
238  {
239  fTrapFlatTop = flat;
240  }
242  {
243  return fTrapRiseTime;
244  }
246  {
247  return fTrapFlatTop;
248  }
249  Double_t ComputeCFDThreshold(Double_t threshold = 0.5);
250 
255  {
256  fSemiGaussSigma = sig;
257  }
259  {
260  return fSemiGaussSigma;
261  }
262 
266  void SetPoleZeroCorrection(Bool_t with = kTRUE)
267  {
269  }
270  void SetTauRC(Int_t taurc)
271  {
272  fTauRC = taurc;
273  }
275  {
276  return fTauRC;
277  }
278 
282  void SetInterpolation(Bool_t with = kTRUE)
283  {
284  fWithInterpolation = with;
285  }
286  void SetInterpolatedChannelWidth(double width)
287  {
289  }
291  {
293  }
294 
298  virtual void ComputeRawAmplitude(void);
300  {
301  return fYmax - fYmin;
302  }
304  {
305  return fYmin;
306  }
308  {
309  return fYmax;
310  }
311 
314  {
316  }
318  {
320  }
321 
323  Bool_t ComputeMeanAndSigma(Int_t start, Int_t stop, Double_t& mean, Double_t& sigma);
324  Bool_t ComputeMeanAndSigma(Double_t start, Double_t stop, Double_t& mean, Double_t& sigma);
325 
327  void Multiply(Double_t fact);
328  void Add(Double_t fact);
329 
331  Double_t ARC_CFD(Double_t threshold = 0.3, Double_t tdelay = 10);
332  double FindTzeroLeadingEdgeCubic(double LEVEL, int Nrecurr);
333  Double_t FindTzeroCFDCubic(double level, int Nrecurr);
334  double FindTzeroCFDCubic_rev(double level, double tend, int Nrecurr);
335  Double_t CubicInterpolation(float* data, int x2, double fmax, int Nrecurr);
336 
337  virtual void BuildCubicSignal(double taufinal); //Interpolazione mediante cubic
338  virtual void BuildCubicSplineSignal(double taufinal);//Interpolazione mediante spline cubic
339  virtual void BuildSmoothingSplineSignal(double taufinal, double l = 1, int nbits = -1); //Interpolazione mediante smoothing spline
340 
341  virtual double GetDataInter(double t);
342  virtual double GetDataInterCubic(double t);
343  virtual double GetDataCubicSpline(double t);
344  virtual double GetDataSmoothingSplineLTI(double t);//metodo che serve a BuildSmoothingSpline
345  virtual double EvalCubicSpline(double X);//metodo che serve a BuildSmoothingSpline
346 
348  void FIR_ApplyTrapezoidal(double trise, double tflat); // trise=sqrt(12)*tausha di CR-RC^4 se tflat=trise/2
349  void FIR_ApplySemigaus(double tau_usec);
350  void FIR_ApplyRCLowPass(double time_usec, int reverse = 0);
351  void FIR_ApplyRCHighPass(double time_usec, int reverse = 0);
352  void FIR_ApplyRecursiveFilter(double a0, int N, double* a, double* b, int reverse);
353  void FIR_ApplyMovingAverage(int npoints);
354  void PoleZeroSuppression(Double_t tauRC);
355  int FIR_ApplySmoothingSpline(double l, int nbits = -1);
356  double ApplyNewton(double l, double x0);//metodo che serve a FIR_ApplySmoothingSpline
357 
359  void ApplyWindowing(int window_type = 3); // 0: barlett, 1:hanning, 2:hamming, 3: blackman
360  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
361  int FFT(bool p_bInverseTransform, double* p_lpRealOut, double* p_lpImagOut);
362  TH1* FFT2Histo(int output, TH1* hh = 0); // 0 modulo, 1 modulo db (normalized), 2, re, 3 im
363 
365  void ApplyModifications(TGraph* newSignal = 0, Int_t nsa = -1);
366 
367 
369  void ShiftLeft(double);//shift in ns
370  void ShiftRight(double);//shift in ns
371  void TestDraw();
372 
373  static KVSignal* MakeSignal(const char* sig_type);
374 
375  ClassDefOverride(KVSignal, 5) //Base class for signal processing
376 
377 };
378 
379 #endif
int Int_t
bool Bool_t
char Char_t
constexpr Bool_t kFALSE
double Double_t
const char Option_t
#define ClassDefOverride(name, id)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h length
Option_t Option_t width
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Double_t fInterpolatedChannelWidth
channel width used to produced the interpolated signal
Definition: KVSignal.h:31
void SetTauRC(Int_t taurc)
Definition: KVSignal.h:270
Double_t ComputeCFDThreshold(Double_t threshold=0.5)
calculate the time during which the signal is higher than th*fAmplitude
Definition: KVSignal.cpp:351
void FIR_ApplyMovingAverage(int npoints)
Definition: KVSignal.cpp:1543
Double_t CubicInterpolation(float *data, int x2, double fmax, int Nrecurr)
Definition: KVSignal.cpp:975
void Add(Double_t fact)
Definition: KVSignal.cpp:1602
const Char_t * GetType() const
Definition: KVSignal.h:82
Double_t GetSemiGaussSigma() const
Definition: KVSignal.h:258
Bool_t IsLongEnough() const
Definition: KVSignal.cpp:109
TString fType
string to identify the signal type
Definition: KVSignal.h:41
void ChangeChannelWidth(Double_t newwidth)
Definition: KVSignal.cpp:286
Double_t GetShaperRiseTime() const
Definition: KVSignal.h:241
virtual Double_t ComputeDuration(Double_t th=0.2)
calculate the time during which the signal is higher than th*fAmplitude
Definition: KVSignal.cpp:321
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:774
void SetInterpolation(Bool_t with=kTRUE)
Definition: KVSignal.h:282
void ApplyModifications(TGraph *newSignal=0, Int_t nsa=-1)
apply modifications of fAdc to the original signal
Definition: KVSignal.cpp:1577
virtual void RemoveBaseLine()
Definition: KVSignal.cpp:406
virtual const Char_t * GetDetectorName() const
Definition: KVSignal.h:73
Double_t fEndLine
mean value of the signal line at the end
Definition: KVSignal.h:24
Double_t fYmax
raw min/max of the signal
Definition: KVSignal.h:21
void SetShaperRiseTime(Double_t rise)
Definition: KVSignal.h:233
virtual void GetPSAResults(KVNameValueList &) const
Definition: KVSignal.h:128
virtual Bool_t IsCharge() const
Definition: KVSignal.h:62
void FIR_ApplyRCLowPass(double time_usec, int reverse=0)
Definition: KVSignal.cpp:1444
Double_t fAmplitude
results of signal treatement
Definition: KVSignal.h:17
double FindTzeroCFDCubic_rev(double level, double tend, int Nrecurr)
Definition: KVSignal.cpp:1403
virtual void TreateSignal()
Definition: KVSignal.cpp:212
virtual Double_t ComputeBaseLine()
Definition: KVSignal.cpp:303
virtual void BuildCubicSplineSignal()
Definition: KVSignal.cpp:1163
void Set(Int_t n) override
Definition: KVSignal.cpp:118
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:505
virtual double EvalCubicSpline(double X)
Definition: KVSignal.cpp:1389
Double_t fSigmaBase
base line rms
Definition: KVSignal.h:23
virtual double GetDataInter(double t)
Definition: KVSignal.cpp:1013
void SetData(Int_t nn, Double_t *xx, Double_t *yy)
operation on data arrays
Definition: KVSignal.cpp:130
virtual Bool_t HasFPGA() const
Definition: KVSignal.h:94
TArrayF * GetArray()
Definition: KVSignal.h:105
void SetInterpolatedChannelWidth(double width)
Definition: KVSignal.h:286
virtual Int_t GetNFPGAValues() const
Definition: KVSignal.h:90
virtual double GetDataCubicSpline(double t)
see HSIEH S.HOU IEEE Trans. Acoustic Speech, vol. ASSP-26, NO.6, DECEMBER 1978
Definition: KVSignal.cpp:1053
void ApplyWindowing(int window_type=3)
fast fourier transform and windowing of the signal (modify only fAdc)
Definition: KVSignal.cpp:738
Int_t GetNSamples() const
Definition: KVSignal.h:158
void print_psa_parameters() const
Definition: KVSignal.cpp:221
virtual void LoadPSAParameters()
Definition: KVSignal.cpp:193
Double_t GetAmplitude() const
Definition: KVSignal.h:220
void SetPoleZeroCorrection(Bool_t with=kTRUE)
Definition: KVSignal.h:266
Bool_t fWithInterpolation
use of interpolation or not
Definition: KVSignal.h:38
Double_t fChannelWidth
channel width in ns
Definition: KVSignal.h:30
Int_t fFirstBL
Definition: KVSignal.h:32
virtual Double_t GetPSAParameter(const Char_t *parname)
Definition: KVSignal.h:113
Double_t GetAmplitudeTriggerValue() const
routines to manage threshold for minimum charge in the detector
Definition: KVSignal.h:313
Double_t GetBLFirst() const
Definition: KVSignal.h:171
bool IsOK()
Definition: KVSignal.cpp:99
Double_t fYmin
Definition: KVSignal.h:21
static KVSignal * MakeSignal(const char *sig_type)
Create new KVSignal instance corresponding to sig_type.
Definition: KVSignal.cpp:1664
Double_t fMinimumValueForAmplitude
Minimum value to say if detector has been hitted.
Definition: KVSignal.h:39
Double_t GetRawAmplitude() const
Definition: KVSignal.h:299
void ShiftLeft(double)
---------------— OPERATORI ------------------—//
Definition: KVSignal.cpp:1612
Double_t fSigmaEnd
rms value of the signal line at the end
Definition: KVSignal.h:25
void SetNSamples(int nn)
Definition: KVSignal.h:154
Double_t fBaseLine
base line mean value
Definition: KVSignal.h:22
Double_t ARC_CFD(Double_t threshold=0.3, Double_t tdelay=10)
Interpolations.
Definition: KVSignal.cpp:466
Double_t fChannelWidthInt
internal parameter channel width of interpolated signal in ns
Definition: KVSignal.h:46
void SetBaseLineLength(Int_t length, Int_t first=0)
Definition: KVSignal.h:166
Double_t fIMax
position of the maximum in channel
Definition: KVSignal.h:19
virtual double GetDataInterCubic(double t)
Definition: KVSignal.cpp:1032
Double_t ComputeAmplitude()
Compute and return the absolute value of the signal amplitude.
Definition: KVSignal.cpp:431
void SetTrapShaperParameters(Double_t rise, Double_t flat)
Definition: KVSignal.h:228
KVSignal()
Default constructor.
Definition: KVSignal.cpp:50
void Multiply(Double_t fact)
multiply the signal (modify only fAdc)
Definition: KVSignal.cpp:1593
Double_t GetSigmaEndLine() const
Definition: KVSignal.h:202
virtual void BuildSmoothingSplineSignal()
Definition: KVSignal.cpp:1236
void init()
Definition: KVSignal.cpp:21
void FIR_ApplyRCHighPass(double time_usec, int reverse=0)
Definition: KVSignal.cpp:1466
virtual double GetDataSmoothingSplineLTI(double t)
Definition: KVSignal.cpp:1370
void SetShaperFlatTop(Double_t flat)
Definition: KVSignal.h:237
int FIR_ApplySmoothingSpline(double l, int nbits=-1)
Definition: KVSignal.cpp:1249
Double_t GetBaseLine() const
Definition: KVSignal.h:189
Int_t fLastBL
first and last channel number to compute the base line
Definition: KVSignal.h:32
TArrayF fAdc
Definition: KVSignal.h:15
void TestDraw()
Definition: KVSignal.cpp:1653
void Print(Option_t *chopt="") const override
Definition: KVSignal.cpp:240
void SetMaxT(double t)
Definition: KVSignal.h:149
Bool_t IsFired()
ComputeBaseLine and ComputeEndLine methods have to be called before.
Definition: KVSignal.cpp:393
void SetADCData()
Definition: KVSignal.cpp:152
Double_t GetEndLine() const
Definition: KVSignal.h:198
void SetType(const Char_t *type)
Definition: KVSignal.h:77
void SetChannelWidth(double width)
Definition: KVSignal.h:137
Bool_t PSAHasBeenComputed() const
Definition: KVSignal.h:129
Double_t GetBLLength() const
Definition: KVSignal.h:175
Double_t GetShaperFlatTop() const
Definition: KVSignal.h:245
Double_t fTauRC
tau_rc of the electronics. Used for pole zero cancellation.
Definition: KVSignal.h:33
virtual Bool_t IsCurrent() const
Definition: KVSignal.h:66
Double_t ComputeRiseTime()
Definition: KVSignal.cpp:448
Double_t fTMax
position of the maximum in ns
Definition: KVSignal.h:20
Double_t GetInterpolatedChannelWidth() const
Definition: KVSignal.h:290
void FIR_ApplyTrapezoidal(double trise, double tflat)
different shapers (modify only fAdc)
Definition: KVSignal.cpp:651
virtual void SetDefaultValues()
To be defined in child class.
Definition: KVSignal.cpp:203
void ShiftRight(double)
Definition: KVSignal.cpp:1631
Double_t GetYmax() const
Definition: KVSignal.h:307
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:1487
Double_t GetSigmaBaseLine() const
Definition: KVSignal.h:193
KVSignal * ConvertTo(const Char_t *type)
Definition: KVSignal.cpp:85
double ApplyNewton(double l, double x0)
Definition: KVSignal.cpp:1327
virtual void ComputeRawAmplitude(void)
Definition: KVSignal.cpp:252
Double_t FindTzeroCFDCubic(double level, int Nrecurr)
Definition: KVSignal.cpp:697
void FIR_ApplySemigaus(double tau_usec)
Definition: KVSignal.cpp:1425
Bool_t fWithPoleZeroCorrection
use or nor pole zero correction
Definition: KVSignal.h:37
Double_t fTrapFlatTop
flat top of the trapezoidal shaper
Definition: KVSignal.h:35
double FindTzeroLeadingEdgeCubic(double LEVEL, int Nrecurr)
Definition: KVSignal.cpp:488
Double_t fRiseTime
rise time of the signal
Definition: KVSignal.h:18
virtual void UpdatePSAParameter(KVNameValueList *par)
Definition: KVSignal.cpp:165
Double_t fSemiGaussSigma
sigma of the semi-gaussian shaper
Definition: KVSignal.h:36
Bool_t fPSAIsDone
indicate if PSA has been done
Definition: KVSignal.h:45
void SetSemiGaussSigma(Double_t sig)
Definition: KVSignal.h:254
TH1 * FFT2Histo(int output, TH1 *hh=0)
Definition: KVSignal.cpp:916
void PoleZeroSuppression(Double_t tauRC)
Definition: KVSignal.cpp:1560
Double_t fTrapRiseTime
rise time of the trapezoidal shaper
Definition: KVSignal.h:34
virtual void BuildCubicSignal()
Definition: KVSignal.cpp:1198
void BuildReverseTimeSignal()
Definition: KVSignal.cpp:418
Double_t GetTauRC() const
Definition: KVSignal.h:274
Bool_t TestWidth() const
Definition: KVSignal.cpp:270
virtual Double_t ComputeEndLine()
Definition: KVSignal.cpp:376
Double_t GetYmin() const
Definition: KVSignal.h:303
Double_t GetRiseTime() const
Definition: KVSignal.h:211
virtual Int_t GetIndex() const
Definition: KVSignal.h:86
Double_t GetChannelWidth() const
Definition: KVSignal.h:142
void SetAmplitudeTriggerValue(Double_t val)
Definition: KVSignal.h:317
void Set(Int_t n) override
Int_t GetSize() const
void SetTitle(const char *title="") override
const char * Data() const