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 #include "KVDBParameterList.h"
12 
13 class KVSignal : public TGraph {
14 public:
15  enum SignalType {
17  kI1,
19  kQ2,
20  kI2,
21  kQ3,
23  kUNKDT
24  };
25 
26 protected:
29 
32 
34 
46 
60 
65  void ResetIndexes();
66  virtual void BuildCubicSignal(); //Interpolazione mediante cubic
67  virtual void BuildCubicSplineSignal(); //Interpolazione mediante cubic spline
68  virtual void BuildSmoothingSplineSignal(); //Interpolazione mediante cubic spline
69  void init();
70  void TreateOldSignalName();
71 
72 public:
73  KVSignal();
74  KVSignal(const char* name, const char* title);
75  KVSignal(const TString& name, const TString& title);
76  virtual ~KVSignal();
77  void Copy(TObject& obj) const override;
78 
79  Bool_t IsLongEnough() const;
80 
84  const Char_t* GetDetectorName() const
85  {
86  return fDetName.Data();
87  }
88  void SetDetectorName(const Char_t* name)
89  {
90  fDetName = name;
91  }
92 
93  void SetType(const Char_t* type)
94  {
95  fType = type;
96  SetTitle(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 override;
119 
121  void SetData(Int_t nn, Double_t* xx, Double_t* yy);
122  void Set(Int_t n) override;
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  virtual void GetPSAResults(KVNameValueList&) const {}
156  {
157  return fPSAIsDone;
158  }
159 
163  void SetChannelWidth(double width)
164  {
167  }
169  {
170  return fChannelWidth;
171  }
172  Bool_t TestWidth() const;
173  void ChangeChannelWidth(Double_t newwidth);
175  void SetMaxT(double t)
176  {
177  fAdc.Set((int)(t / fChannelWidth));
178  }
180  void SetNSamples(int nn)
181  {
182  fAdc.Set(nn);
183  }
185  {
186  return fAdc.GetSize();
187  }
188 
192  void SetBaseLineLength(Int_t length, Int_t first = 0)
193  {
194  fFirstBL = first;
195  fLastBL = length - first;
196  }
198  {
199  return fFirstBL;
200  }
202  {
203  return fLastBL - fFirstBL;
204  }
205 
206  virtual Double_t ComputeBaseLine();
207  virtual Double_t ComputeDuration(Double_t th = 0.2);
208  virtual Double_t ComputeEndLine();
209  virtual void RemoveBaseLine();
210 
211  void BuildReverseTimeSignal();
212 
213  Bool_t IsFired();
214 
216  {
217  return fBaseLine;
218  }
220  {
221  return fSigmaBase;
222  }
223 
225  {
226  return fEndLine;
227  }
229  {
230  return fSigmaEnd;
231  }
232 
238  {
239  return fRiseTime;
240  }
241 
247  {
248  return fAmplitude;
249  }
250 
255  {
256  fTrapRiseTime = rise;
257  fTrapFlatTop = flat;
258  }
260  {
261  fTrapRiseTime = rise;
262  }
264  {
265  fTrapFlatTop = flat;
266  }
268  {
269  return fTrapRiseTime;
270  }
272  {
273  return fTrapFlatTop;
274  }
275  Double_t ComputeCFDThreshold(Double_t threshold = 0.5);
276 
281  {
282  fSemiGaussSigma = sig;
283  }
285  {
286  return fSemiGaussSigma;
287  }
288 
292  void SetPoleZeroCorrection(Bool_t with = kTRUE)
293  {
295  }
296  void SetTauRC(Int_t taurc)
297  {
298  fTauRC = taurc;
299  }
301  {
302  return fTauRC;
303  }
304 
308  void SetInterpolation(Bool_t with = kTRUE)
309  {
310  fWithInterpolation = with;
311  }
312  void SetInterpolatedChannelWidth(double width)
313  {
315  }
317  {
319  }
320 
324  virtual void ComputeRawAmplitude(void);
326  {
327  return fYmax - fYmin;
328  }
330  {
331  return fYmin;
332  }
334  {
335  return fYmax;
336  }
337 
340  {
342  }
344  {
346  }
347 
349  Bool_t ComputeMeanAndSigma(Int_t start, Int_t stop, Double_t& mean, Double_t& sigma);
350  Bool_t ComputeMeanAndSigma(Double_t start, Double_t stop, Double_t& mean, Double_t& sigma);
351 
356 
358  void Multiply(Double_t fact);
359  void Add(Double_t fact);
360 
362  Double_t ARC_CFD(Double_t threshold = 0.3, Double_t tdelay = 10);
363  double FindTzeroLeadingEdgeCubic(double LEVEL, int Nrecurr);
364  Double_t FindTzeroCFDCubic(double level, int Nrecurr);
365  double FindTzeroCFDCubic_rev(double level, double tend, int Nrecurr);
366  Double_t CubicInterpolation(float* data, int x2, double fmax, int Nrecurr);
367 
368  virtual void BuildCubicSignal(double taufinal); //Interpolazione mediante cubic
369  virtual void BuildCubicSplineSignal(double taufinal);//Interpolazione mediante spline cubic
370  virtual void BuildSmoothingSplineSignal(double taufinal, double l = 1, int nbits = -1); //Interpolazione mediante smoothing spline
371 
372  virtual double GetDataInter(double t);
373  virtual double GetDataInterCubic(double t);
374  virtual double GetDataCubicSpline(double t);
375  virtual double GetDataSmoothingSplineLTI(double t);//metodo che serve a BuildSmoothingSpline
376  virtual double EvalCubicSpline(double X);//metodo che serve a BuildSmoothingSpline
377 
379  void FIR_ApplyTrapezoidal(double trise, double tflat); // trise=sqrt(12)*tausha di CR-RC^4 se tflat=trise/2
380  void FIR_ApplySemigaus(double tau_usec);
381  void FIR_ApplyRCLowPass(double time_usec, int reverse = 0);
382  void FIR_ApplyRCHighPass(double time_usec, int reverse = 0);
383  void FIR_ApplyRecursiveFilter(double a0, int N, double* a, double* b, int reverse);
384  void FIR_ApplyMovingAverage(int npoints);
385  void PoleZeroSuppression(Double_t tauRC);
386  int FIR_ApplySmoothingSpline(double l, int nbits = -1);
387  double ApplyNewton(double l, double x0);//metodo che serve a FIR_ApplySmoothingSpline
388 
390  void ApplyWindowing(int window_type = 3); // 0: barlett, 1:hanning, 2:hamming, 3: blackman
391  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
392  int FFT(bool p_bInverseTransform, double* p_lpRealOut, double* p_lpImagOut);
393  TH1* FFT2Histo(int output, TH1* hh = 0); // 0 modulo, 1 modulo db (normalized), 2, re, 3 im
394 
396  void ApplyModifications(TGraph* newSignal = 0, Int_t nsa = -1);
397 
398 
400  void ShiftLeft(double);//shift in ns
401  void ShiftRight(double);//shift in ns
402  void TestDraw();
403 
404  static KVSignal* MakeSignal(const char* sig_type);
405 
406  ClassDefOverride(KVSignal, 4) //Base class for FAZIA signal processing
407 
408 };
409 
410 #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
char name[80]
To store calibration parameters in a database ,.
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:51
void SetTauRC(Int_t taurc)
Definition: KVSignal.h:296
Double_t ComputeCFDThreshold(Double_t threshold=0.5)
calculate the time during which the signal is higher than th*fAmplitude
Definition: KVSignal.cpp:538
void FIR_ApplyMovingAverage(int npoints)
Definition: KVSignal.cpp:1730
Double_t CubicInterpolation(float *data, int x2, double fmax, int Nrecurr)
Definition: KVSignal.cpp:1162
void Add(Double_t fact)
Definition: KVSignal.cpp:1789
const Char_t * GetType() const
Definition: KVSignal.h:98
Double_t GetSemiGaussSigma() const
Definition: KVSignal.h:284
Bool_t IsLongEnough() const
Definition: KVSignal.cpp:185
TString fType
string to identify the signal type : "QH1", "I2" etc ...
Definition: KVSignal.h:31
void ChangeChannelWidth(Double_t newwidth)
Definition: KVSignal.cpp:473
Double_t GetShaperRiseTime() const
Definition: KVSignal.h:267
virtual void UpdatePSAParameter(KVDBParameterList *par)
Definition: KVSignal.cpp:358
virtual Double_t ComputeDuration(Double_t th=0.2)
calculate the time during which the signal is higher than th*fAmplitude
Definition: KVSignal.cpp:508
SignalType
Definition: KVSignal.h:15
@ kQH1
Definition: KVSignal.h:16
@ kADC
Definition: KVSignal.h:22
@ kUNKDT
Definition: KVSignal.h:23
@ kQL1
Definition: KVSignal.h:18
Int_t fChannel
signal type (see KVSignal::SignalType enum)
Definition: KVSignal.h:28
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:961
void SetInterpolation(Bool_t with=kTRUE)
Definition: KVSignal.h:308
void ApplyModifications(TGraph *newSignal=0, Int_t nsa=-1)
apply modifications of fAdc to the original signal
Definition: KVSignal.cpp:1764
virtual void RemoveBaseLine()
Definition: KVSignal.cpp:593
Double_t fEndLine
mean value of the signal line at the end
Definition: KVSignal.h:44
Double_t fYmax
raw min/max of the signal
Definition: KVSignal.h:41
void SetShaperRiseTime(Double_t rise)
Definition: KVSignal.h:259
virtual void GetPSAResults(KVNameValueList &) const
Definition: KVSignal.h:154
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:1631
Double_t fAmplitude
results of signal treatement
Definition: KVSignal.h:37
double FindTzeroCFDCubic_rev(double level, double tend, int Nrecurr)
Definition: KVSignal.cpp:1590
virtual void TreateSignal()
Definition: KVSignal.cpp:405
virtual Double_t ComputeBaseLine()
Definition: KVSignal.cpp:490
virtual void BuildCubicSplineSignal()
Definition: KVSignal.cpp:1350
void Set(Int_t n) override
Definition: KVSignal.cpp:194
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:692
Double_t GetPSAParameter(const Char_t *parname)
DeduceFromName has to be called before.
Definition: KVSignal.cpp:340
virtual double EvalCubicSpline(double X)
Definition: KVSignal.cpp:1576
Double_t fSigmaBase
base line rms
Definition: KVSignal.h:43
virtual double GetDataInter(double t)
Definition: KVSignal.cpp:1200
void SetData(Int_t nn, Double_t *xx, Double_t *yy)
operation on data arrays
Definition: KVSignal.cpp:206
void SetDetectorName(const Char_t *name)
Definition: KVSignal.h:88
TArrayF * GetArray()
Definition: KVSignal.h:124
void SetInterpolatedChannelWidth(double width)
Definition: KVSignal.h:312
virtual double GetDataCubicSpline(double t)
see HSIEH S.HOU IEEE Trans. Acoustic Speech, vol. ASSP-26, NO.6, DECEMBER 1978
Definition: KVSignal.cpp:1240
void ApplyWindowing(int window_type=3)
fast fourier transform and windowing of the signal (modify only fAdc)
Definition: KVSignal.cpp:925
Int_t GetNSamples() const
Definition: KVSignal.h:184
virtual void LoadPSAParameters()
Definition: KVSignal.cpp:386
TString fDetName
name of the detector, the signal is linked to, needed to find it in the KVMultiDetector
Definition: KVSignal.h:30
Int_t fFPGAOutputNumbers
ASsociated FPGA energy outputs.
Definition: KVSignal.h:33
Double_t GetAmplitude() const
Definition: KVSignal.h:246
void SetPoleZeroCorrection(Bool_t with=kTRUE)
Definition: KVSignal.h:292
Bool_t fWithInterpolation
use of interpolation or not
Definition: KVSignal.h:58
Double_t fChannelWidth
channel width in ns
Definition: KVSignal.h:50
Int_t fFirstBL
Definition: KVSignal.h:52
Double_t GetAmplitudeTriggerValue() const
routines to manage threshold for minimum charge in the detector
Definition: KVSignal.h:339
Double_t GetBLFirst() const
Definition: KVSignal.h:197
virtual bool IsOK()
Definition: KVSignal.cpp:148
Double_t fYmin
Definition: KVSignal.h:41
static KVSignal * MakeSignal(const char *sig_type)
Create new KVSignal instance corresponding to sig_type.
Definition: KVSignal.cpp:1851
Double_t fMinimumValueForAmplitude
Minimum value to say if detector has been hitted.
Definition: KVSignal.h:59
Double_t GetRawAmplitude() const
Definition: KVSignal.h:325
void ShiftLeft(double)
---------------— OPERATORI ------------------—//
Definition: KVSignal.cpp:1799
Double_t fSigmaEnd
rms value of the signal line at the end
Definition: KVSignal.h:45
void SetNSamples(int nn)
Definition: KVSignal.h:180
Double_t fBaseLine
base line mean value
Definition: KVSignal.h:42
Double_t ARC_CFD(Double_t threshold=0.3, Double_t tdelay=10)
Interpolations.
Definition: KVSignal.cpp:653
const Char_t * GetDetectorName() const
Definition: KVSignal.h:84
Double_t fChannelWidthInt
internal parameter channel width of interpolated signal in ns
Definition: KVSignal.h:64
void SetBaseLineLength(Int_t length, Int_t first=0)
Definition: KVSignal.h:192
Double_t fIMax
position of the maximum in channel
Definition: KVSignal.h:39
virtual double GetDataInterCubic(double t)
Definition: KVSignal.cpp:1219
Double_t ComputeAmplitude()
Compute and return the absolute value of the signal amplitude.
Definition: KVSignal.cpp:618
void SetTrapShaperParameters(Double_t rise, Double_t flat)
Definition: KVSignal.h:254
KVSignal()
Default constructor.
Definition: KVSignal.cpp:87
void Multiply(Double_t fact)
multiply the signal (modify only fAdc)
Definition: KVSignal.cpp:1780
Double_t GetSigmaEndLine() const
Definition: KVSignal.h:228
virtual void BuildSmoothingSplineSignal()
Definition: KVSignal.cpp:1423
void init()
Definition: KVSignal.cpp:41
void FIR_ApplyRCHighPass(double time_usec, int reverse=0)
Definition: KVSignal.cpp:1653
void Copy(TObject &obj) const override
Definition: KVSignal.cpp:168
virtual double GetDataSmoothingSplineLTI(double t)
Definition: KVSignal.cpp:1557
void SetShaperFlatTop(Double_t flat)
Definition: KVSignal.h:263
int FIR_ApplySmoothingSpline(double l, int nbits=-1)
Definition: KVSignal.cpp:1436
Double_t GetBaseLine() const
Definition: KVSignal.h:215
Int_t fLastBL
first and last channel number to compute the base line
Definition: KVSignal.h:52
TArrayF fAdc
Definition: KVSignal.h:35
void TestDraw()
Definition: KVSignal.cpp:1840
Int_t GetNFPGAValues() const
Definition: KVSignal.h:129
void Print(Option_t *chopt="") const override
Definition: KVSignal.cpp:414
void SetMaxT(double t)
Definition: KVSignal.h:175
Bool_t IsFired()
ComputeBaseLine and ComputeEndLine methods have to be called before.
Definition: KVSignal.cpp:580
void SetADCData()
Definition: KVSignal.cpp:228
Double_t GetEndLine() const
Definition: KVSignal.h:224
void SetType(const Char_t *type)
Definition: KVSignal.h:93
void SetChannelWidth(double width)
Definition: KVSignal.h:163
Bool_t PSAHasBeenComputed() const
Definition: KVSignal.h:155
Int_t fIndex
index deduced from block, quartet and telescope numbering
Definition: KVSignal.h:27
Double_t GetBLLength() const
Definition: KVSignal.h:201
void DeduceFromName()
Definition: KVSignal.cpp:287
Double_t GetShaperFlatTop() const
Definition: KVSignal.h:271
Double_t fTauRC
tau_rc of the electronics. Used for pole zero cancellation.
Definition: KVSignal.h:53
virtual Bool_t IsCurrent() const
Definition: KVSignal.h:113
Double_t ComputeRiseTime()
Definition: KVSignal.cpp:635
Double_t fTMax
position of the maximum in ns
Definition: KVSignal.h:40
Double_t GetInterpolatedChannelWidth() const
Definition: KVSignal.h:316
void FIR_ApplyTrapezoidal(double trise, double tflat)
different shapers (modify only fAdc)
Definition: KVSignal.cpp:838
virtual void SetDefaultValues()
To be defined in child class.
Definition: KVSignal.cpp:396
void ShiftRight(double)
Definition: KVSignal.cpp:1818
Double_t GetYmax() const
Definition: KVSignal.h:333
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:1674
Double_t GetSigmaBaseLine() const
Definition: KVSignal.h:219
KVSignal * ConvertTo(const Char_t *type)
Definition: KVSignal.cpp:132
double ApplyNewton(double l, double x0)
Definition: KVSignal.cpp:1514
virtual void ComputeRawAmplitude(void)
Definition: KVSignal.cpp:439
Bool_t HasFPGA() const
Definition: KVSignal.h:134
Double_t FindTzeroCFDCubic(double level, int Nrecurr)
Definition: KVSignal.cpp:884
void FIR_ApplySemigaus(double tau_usec)
Definition: KVSignal.cpp:1612
Bool_t fWithPoleZeroCorrection
use or nor pole zero correction
Definition: KVSignal.h:57
Double_t fTrapFlatTop
flat top of the trapezoidal shaper
Definition: KVSignal.h:55
double FindTzeroLeadingEdgeCubic(double LEVEL, int Nrecurr)
Definition: KVSignal.cpp:675
Double_t fRiseTime
rise time of the signal
Definition: KVSignal.h:38
void TreateOldSignalName()
Definition: KVSignal.cpp:242
Double_t fSemiGaussSigma
sigma of the semi-gaussian shaper
Definition: KVSignal.h:56
Bool_t fPSAIsDone
indicate if PSA has been done
Definition: KVSignal.h:63
void SetSemiGaussSigma(Double_t sig)
Definition: KVSignal.h:280
TH1 * FFT2Histo(int output, TH1 *hh=0)
Definition: KVSignal.cpp:1103
void PoleZeroSuppression(Double_t tauRC)
Definition: KVSignal.cpp:1747
Double_t fTrapRiseTime
rise time of the trapezoidal shaper
Definition: KVSignal.h:54
virtual void BuildCubicSignal()
Definition: KVSignal.cpp:1385
void BuildReverseTimeSignal()
Definition: KVSignal.cpp:605
Double_t GetTauRC() const
Definition: KVSignal.h:300
Bool_t TestWidth() const
Definition: KVSignal.cpp:457
virtual ~KVSignal()
Destructor.
Definition: KVSignal.cpp:123
virtual Double_t ComputeEndLine()
Definition: KVSignal.cpp:563
Double_t GetYmin() const
Definition: KVSignal.h:329
Double_t GetRiseTime() const
Definition: KVSignal.h:237
void ResetIndexes()
Definition: KVSignal.cpp:75
Double_t GetChannelWidth() const
Definition: KVSignal.h:168
void SetAmplitudeTriggerValue(Double_t val)
Definition: KVSignal.h:343
void Set(Int_t n) override
Int_t GetSize() const
void SetTitle(const char *title="") override
const char * Data() const