KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVSignal.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Mon Dec 22 15:46:46 2014
2 //Author: ,,,
3 
4 #include "KVSignal.h"
5 #include "KVString.h"
6 #include "TMath.h"
7 #include "KVDigitalFilter.h"
8 #include "KVDataSet.h"
9 #include "KVEnv.h"
10 #include "KVDBParameterList.h"
11 
12 #include "TMatrixD.h"
13 #include "TMatrixF.h"
14 #include "TClass.h"
15 
16 #define LOG2 (double)6.93147180559945286e-01
17 # define M_PI 3.14159265358979323846 /* pi */
18 
20 
21 // BEGIN_HTML <!--
23 /* -->
24 <h2>KVSignal</h2>
25 <h4>Base class for FAZIA signal processing</h4>
26 <!-- */
27 // --> END_HTML
28 //
29 // SIGNAL NAME
30 // ===========
31 // KVSignal::GetName() returns the name of the signal as defined in FAZIA
32 // raw data, e.g. "T1-Q3-B001-QH1"
33 //
34 // SIGNAL TYPE
35 // ===========
36 // KVSignal::GetType() returns one of: "QH1", "QL1", "Q2", "Q3", "I1", "I2"
38 
39 
40 
42 void KVSignal::init()
43 {
44  fPSAIsDone = kFALSE;
45  fChannel = kUNKDT;
46  fYmin = fYmax = 0;
47  fAmplitude = 0;
48  fRiseTime = 0;
49  fIMax = 0;
50  fTMax = 0;
51  fBaseLine = 0;
52  fSigmaBase = 0;
53 
54  fChannelWidth = -1;
55  fChannelWidthInt = -1;
56  fFirstBL = -1;
57  fLastBL = -1;
58  fTauRC = -1;
59  fTrapRiseTime = -1;
60  fTrapFlatTop = -1;
61  fSemiGaussSigma = -1;
62  fWithPoleZeroCorrection = kFALSE;
63  fWithInterpolation = kFALSE;
64  fMinimumValueForAmplitude = 0;
65 
66  //DeduceFromName();
67  ResetIndexes();
68 
69  SetDefaultValues();
70 }
71 
72 
73 
75 
77 {
78  fIndex = 0;
79  fDetName = "";
81 }
82 
83 
84 
87 
89 {
90  // Default constructor
91  init();
92 }
93 
94 
95 
96 
99 
100 KVSignal::KVSignal(const char* name, const char* title) : TGraph()
101 {
102  // Write your code here
103  SetNameTitle(name, title);
104  init();
105 }
106 
107 
108 
109 
111 
112 KVSignal::KVSignal(const TString& name, const TString& title) : TGraph()
113 {
114  SetNameTitle(name, title);
115  init();
116 
117 }
118 
119 
120 
123 
125 {
126  // Destructor
127 }
128 
129 
130 
132 
134 {
135  KVSignal* sig = 0;
136  TClass* cl = TClass::GetClass(Form("KV%s", type));
137  if (cl) {
138  sig = (KVSignal*)cl->New();
139  sig->SetData(this->GetN(), this->GetX(), this->GetY());
140  sig->LoadPSAParameters();
141  }
142  return sig;
143 }
144 
145 
146 
148 
150 {
151  if (fLastBL > GetN()) return kFALSE;
152  else return kTRUE;
153 }
154 
155 
156 
157 
158 
168 
170 {
171  // This method copies the current state of 'this' object into 'obj'
172  // You should add here any member variables, for example:
173  // (supposing a member variable KVSignal::fToto)
174  // CastedObj.fToto = fToto;
175  // or
176  // CastedObj.SetToto( GetToto() );
177 
178  //TGraph::Copy((TGraph&)obj);
179  //KVSignal& CastedObj = (KVSignal&)obj;
180 }
181 
182 
183 
185 
187 {
188  return (GetN() > fLastBL + 10);
189 }
190 
191 
192 
194 
196 {
197  fPSAIsDone = kFALSE;
198  TGraph::Set(n);
199  fAdc.Set(GetN());
200 
201 }
202 
203 
204 
206 
208 {
209  Set(nn);
210  if (nn == 0) {
211  Info("SetData", "called with points number=%d", nn);
212  return;
213  }
214  Int_t np = 0;
215  fYmin = fYmax = yy[np];
216  SetPoint(np, xx[np], yy[np]);
217  for (np = 1; np < nn; np += 1) {
218  SetPoint(np, xx[np], yy[np]);
219  if (yy[np] < fYmin) fYmin = yy[np];
220  if (yy[np] > fYmax) fYmax = yy[np];
221  }
222  SetADCData();
223 }
224 
225 
226 
228 
230 {
231 
233  fAdc.Set(GetN());
234  for (int ii = 0; ii < GetN(); ii++) fAdc.AddAt(fY[ii], ii);
235 
236 }
237 
238 
239 
240 
242 
244 {
245  TString stit = GetTitle();
246  stit.ToUpper();
247 
248  KVString tmp = GetName();
249  KVString part = "";
250  if (tmp.BeginsWith("B")) {
251  tmp.Begin("-");
252  part = tmp.Next();
253  part.ReplaceAll("B", "");
254  Int_t bb = part.Atoi();
255  part = tmp.Next();
256  part.ReplaceAll("Q", "");
257  Int_t qq = part.Atoi();
258  part = tmp.Next();
259  part.ReplaceAll("T", "");
260  Int_t tt = part.Atoi();
261  fType = tmp.Next();
262 
263  fIndex = 100 * bb + 10 * qq + tt;
264  fDetName.Form("%s-%d", stit.Data(), fIndex);
265  }
266  else if (tmp.BeginsWith("RUTH")) {
267  //Old FAZIA telescope denomination
268  // Rutherford telescope case for FAZIASYM experiment
269  tmp.Begin("-");
270  part = tmp.Next();
271  fType = tmp.Next();
272 
273  TString stit = GetTitle();
274  stit.ToUpper();
275 
276  fDetName.Form("%s-RUTH", stit.Data());
277  }
278 }
279 
280 
281 
287 
289 {
290  //methods used to identify to from which detector/telescope/quartet/block
291  //the signal is associated
292  //it is assumed that the name of the signal is defined as it is in the raw data
293  //generated by fazia_reader which convert raw acquisition file in raw ROOT files
294 
295  ResetIndexes();
296  KVString tmp = GetName();
297  KVString part = "";
298  if (tmp.BeginsWith("B") || tmp.BeginsWith("RUTH")) {
299  //Old FAZIA telescope denomination
300  //Info("DeduceFromName","Old format %s",GetName());
302  }
303  else {
304 
305  if (tmp.GetNValues("-") == 2) {
306  //new general denomination
307  //in KaliVeda :
308  // signal name : [signal_type]-[index]
309  // if index is digit (number), it is the telescope number : 100*b+10*q+t
310  // if index is a string of character, we kept it as it is
311  // for example RUTH for the FAZIASYM dataset
312  TString stit = GetTitle();
313  stit.ToUpper();
314  //new format
315  //Info("DeduceFromName", "New format %s", GetName());
316  tmp.Begin("-");
317  fType = tmp.Next();
318  KVString ss = tmp.Next();
319 
320  if (ss.IsDigit()) {
321  fIndex = ss.Atoi();
322  fDetName.Form("%s-%d", stit.Data(), fIndex);
323  }
324  else {
325  fDetName.Form("%s-%s", stit.Data(), ss.Data());
326  }
327  }
328  else {
329  Warning("DeduceFromName", "Unkown format %s", GetName());
330  }
331  }
332 
333 
334 }
335 
336 
337 
340 
342 {
343 
344  //DeduceFromName has to be called before
345 
346  Double_t lval = -1;
347  KVString spar;
348  spar.Form("%s.%s", GetType(), parname);
349  if (gDataSet) lval = gDataSet->GetDataSetEnv(spar.Data(), 0.0);
350  else lval = gEnv->GetValue(spar.Data(), 0.0);
351  return lval;
352 
353 }
354 
355 
356 
358 
360 {
361  for (Int_t ii = 0; ii < par->GetParameters()->GetNpar(); ii += 1) {
362  TString nameat(par->GetParameters()->GetNameAt(ii));
363  if (nameat == "BaseLineLength") SetBaseLineLength(par->GetParameters()->GetDoubleValue(ii));
364  else if (nameat == "ChannelWidth") SetChannelWidth(par->GetParameters()->GetDoubleValue(ii));
365  else if (nameat == "ShaperRiseTime") SetShaperRiseTime(par->GetParameters()->GetDoubleValue(ii));
366  else if (nameat == "ShaperFlatTop") SetShaperFlatTop(par->GetParameters()->GetDoubleValue(ii));
367  else if (nameat == "TauRC") SetTauRC(par->GetParameters()->GetDoubleValue(ii));
368  else if (nameat == "MinimumAmplitude") SetAmplitudeTriggerValue(par->GetParameters()->GetDoubleValue(ii));
369  else if (nameat == "InterpolatedChannelWidth") SetInterpolatedChannelWidth(par->GetParameters()->GetDoubleValue(ii));
370  else if (nameat == "Interpolation") SetInterpolation((par->GetParameters()->GetDoubleValue(ii) == 1));
371  else if (nameat == "PZCorrection") SetPoleZeroCorrection((par->GetParameters()->GetDoubleValue(ii) == 1));
372  else {
373  if (nameat == "Detector" || nameat == "Signal" || nameat == "RunRange") {
374 
375  }
376  else {
377  Warning("UpdatePSAParameter", "Not supported PSA parameter : %d %s\n", ii, nameat.Data());
378  }
379  }
380  }
381 }
382 
383 
384 
386 
388 {
389  Info("LoadPSAParameters", "To be defined in child class");
390 }
391 
392 
393 
396 
398 {
399  //To be defined in child class
400 }
401 
402 
403 
405 
407 {
408  Info("TreateSignal", "To be defined in child class");
409 }
410 
411 
412 
414 
416 {
417  Info("Print", "\nName: %s - Title: %s", GetName(), GetTitle());
418  if (fDetName != "") {
419  printf("\tAssociated to the detector %s\n", fDetName.Data());
420  }
421  printf("################\nPSA parameters:\n");
422  printf("\tBaseLine: length: %lf first: %lf\n", GetBLLength(), GetBLFirst());
423  printf("\tChannelWidth: %lf\n", GetChannelWidth());
424  printf("\tTauRC: %lf\n", GetTauRC());
425  printf("\tShaperRiseTime: %lf\n", GetShaperRiseTime());
426  printf("\tShaperFlatTop: %lf\n", GetShaperFlatTop());
427  printf("\tWith Interpolation: %d", Int_t(fWithInterpolation));
428  if (fWithInterpolation)
429  printf(" %1.2lf", GetInterpolatedChannelWidth());
430  printf("\n");
431  printf("\tWith PoleZero correction: %d\n", Int_t(fWithPoleZeroCorrection));
432 
433 }
434 
435 
436 
437 
439 
441 {
442  Double_t xx, yy;
443  Int_t np = 0;
444  GetPoint(np++, xx, yy);
445  fYmin = fYmax = yy;
446 
447  for (np = 1; np < GetN(); np += 1) {
448  GetPoint(np, xx, yy);
449  if (yy < fYmin) fYmin = yy;
450  if (yy > fYmax) fYmax = yy;
451  }
452 }
453 
454 
455 
457 
459 {
460  Double_t x0, x1, y0, y1;
461 
462  GetPoint(0, x0, y0);
463  GetPoint(1, x1, y1);
464 
465  Double_t actual_width = x1 - x0;
466  return ((actual_width == GetChannelWidth()));
467 }
468 
469 
470 
471 
473 
475 {
476  Double_t xx, yy;
477  for (Int_t ii = 0; ii < GetN(); ii += 1) {
478  GetPoint(ii, xx, yy);
479  SetPoint(ii, ii * newwidth, yy);
480  }
481 }
482 
483 
484 
485 
486 
490 
492 {
493  //compute mean value of the signal and the rms between
494  // limits defined by fFirstBL and fLastBL
495  if (fLastBL >= GetN()) {
496  fBaseLine = -1;
497  fSigmaBase = -1;
498  return -1;
499  }
501  return fBaseLine;
502 }
503 
504 
505 
508 
510 {
511  // calculate the time during which the signal is higher than th*fAmplitude
512  if (fAmplitude <= 0) {
513  SetADCData();
515  }
516 
517  Multiply(-1);
518  fAmplitude *= -1;
519  Double_t qstart = FindTzeroCFDCubic(th, 3);
521  Double_t qend = FindTzeroCFDCubic(th, 3);
523  Multiply(-1);
524  fAmplitude *= -1;
525 
526  //printf("qstart=%lf - qend=%lf\n",qstart,qend);
527  if (qend < qstart || qstart <= 0 || qend <= 0)
528  return -1;
529  //Double_t deltat = GetN() * fChannelWidth - qstart - qend;
530  Double_t deltat = ((GetN() * fChannelWidth - qend) - qstart);
531  return deltat;
532 }
533 
534 
535 
538 
540 {
541  // calculate the time during which the signal is higher than th*fAmplitude
542  if (fAmplitude <= 0) {
543  SetADCData();
545  }
546 
547  Multiply(-1);
548  fAmplitude *= -1;
549 
550  Double_t qstart = FindTzeroCFDCubic(threshold, 3);
551 
552  Multiply(-1);
553  fAmplitude *= -1;
554 
555  return qstart;
556 }
557 
558 
559 
563 
565 {
566  //same as ComputeBaseLine method but made on the end of the signal
567  //in the same length as for the base line
568  if ((fLastBL - fFirstBL) >= GetN()) {
569  fEndLine = -1;
570  fSigmaEnd = -1;
571  }
573  return fEndLine;
574 }
575 
576 
577 
580 
582 {
583  //ComputeBaseLine and ComputeEndLine methods have to be called before
584 
586  return kTRUE;
587  return kFALSE;
588 }
589 
590 
591 
593 
595 {
596 
597  Add(-1.*ComputeBaseLine());
599 
600 }
601 
602 
603 
605 
607 {
608  TArrayF copy(fAdc);
609  Int_t np = fAdc.GetSize();
610 
611  for (int i = 0; i < np; i++) fAdc.AddAt(copy.GetAt(np - i - 1), i);
612 }
613 
614 
615 
618 
620 {
621  //Compute and return the absolute value of the signal amplitude
622  int i_max = 0;
623  for (int i = 0; i < fAdc.GetSize(); i++) {
624  if (fAdc.At(i_max) < fAdc.At(i)) i_max = i;
625  }
626  fIMax = i_max;
628  fAmplitude = fAdc.At(i_max);
629  return fAmplitude;
630 }
631 
632 
633 
635 
637 {
638  Multiply(-1);
639  fAmplitude *= -1;
640  Double_t qt70 = FindTzeroCFDCubic(0.7, 3);
641  Double_t qt20 = FindTzeroCFDCubic_rev(0.2, qt70, 3);
642  Double_t qtrise72 = qt70 - qt20;
643  Multiply(-1);
644  fAmplitude *= -1;
645  fRiseTime = qtrise72;
646  return fRiseTime;
647 }
648 
649 
650 
653 
655 {
656  //time of passage of the threshold
657 
658  Double_t rtime = GetRiseTime();
659  if (!(tdelay < (1 - threshold)*rtime))
660  return -1;
661 
662  rtime *= 2;
663  double time = GetAmplitude() / ((1 - threshold) * (rtime / tdelay));
664  Multiply(-1);
665  double re = FindTzeroLeadingEdgeCubic(time, 3);
666  Multiply(-1);
667 
668  return re;
669 
670 }
671 
672 
673 
675 
676 double KVSignal::FindTzeroLeadingEdgeCubic(double LEVEL, int Nrecurr)
677 {
678 
679  int i = 0;
680  float* data = fAdc.GetArray();
681  for (; i < fAdc.GetSize() - 1; i++)
682  if (-data[i] > LEVEL && -data[i + 1] > LEVEL) //to skip noise.
683  break;
684  if (i >= fAdc.GetSize() - 3) return -1;
685  i--;
686 
687  return CubicInterpolation(data, i, - LEVEL, Nrecurr);
688 }
689 
690 
691 
692 //----------------------------------
694 //----------------------------------
695 
696 
704 {
705  //compute mean value and sigma values
706  //between "start" point included and "stop" point excluded
707  //for example :
708  // ComputeMeanAndSigma(0,50,mean,sigma)
709  // compute values for the first 50 points
710  //
711  Int_t np = 0;
712  Double_t xx;
713  mean = 0;
714  Double_t mean2 = 0;
715  if (stop > GetN()) {
716  Warning("ComputeMeanAndSigma",
717  "stop position greater than number of samples %d/%d, set stop to %d", GetN(), stop, GetN());
718  stop = GetN();
719  }
720  if (start < 0) {
721  Warning("ComputeMeanAndSigma",
722  "start position unrealistic %d, set start to 0", start);
723  start = 0;
724  }
725 
726  for (Int_t ii = start; ii < stop; ii += 1) {
727  xx = fAdc.At(ii);
728  mean += xx;
729  mean2 += xx * xx;
730  np += 1;
731  }
732 
733  if (np == 0) {
734  Error("ComputeMeanAndSigma", "values cannot be computed with 0 sample");
735  return kFALSE;
736  }
737  mean /= np;
738  mean2 /= np;
739  sigma = TMath::Sqrt(mean2 - mean * mean);
740  return kTRUE;
741 }
742 
743 
744 //----------------------------------
746 //----------------------------------
747 
748 
756 {
757  //
758  //assuming that X axis is in time unit (ms)
759  //divide the X values by the fChannelWidth value which allow to set the Xaxis in time units
760  //compute mean value and sigma values
761  //between "start" and "stop" point included
762  //
763  return ComputeMeanAndSigma(Int_t(start / fChannelWidth), Int_t(stop / fChannelWidth), mean, sigma);
764 }
765 
766 
767 // double KVSignal::FindMedia(double tsta, double tsto)
768 // {
769 // Info("FindMedia(double tsta, double tsto)","Appel ...");
770 // int n1 = (int)(tsta / fChannelWidth); // Non molto preciso, ma tant'e'...
771 // int n2 = (int)(tsto / fChannelWidth);
772 //
773 // return FindMedia(n1, n2);
774 // }
775 //
776 // double KVSignal::FindMedia(int tsta, int tsto)
777 // {
778 // // Calcolo della media nel tratto tra tsta e tsto.
779 // // NOTA: questo ha senso solo se il segnale e' piatto in quella regione!!
780 //
781 // int n1 = (int)(tsta);
782 // int n2 = (int)(tsto);
783 //
784 // int N = n2 - n1 + 1;
785 // //// printf("n1=%d, n2=%d, n=%d, fChannelWidth=%e \n",n1, n2, N, fChannelWidth);
786 // if (n1 < 0 || n1 >= fAdc.GetSize() ||
787 // n2 < n1 || n2 >= fAdc.GetSize() ||
788 // N <= 0 || N >= fAdc.GetSize()) {
789 // printf("--- FSignal::FindMedia: tsta=%d, tsto=%d ?? (%d)\n", tsta, tsto, fAdc.GetSize());
790 // return -1E10;//non cambiare, serve a FindSigma2!!
791 // }
792 // double media = 0;
793 // for (int i = n1; i <= n2; i++)
794 // media += fAdc.At(i);
795 //
796 // media /= N;
797 // return media;
798 // }
799 //
800 // double KVSignal::FindSigma2(double tsta, double tsto)
801 // {
802 // // Calcolo della varianza nel tratto tra tsta e tsto.
803 // // NOTA: questo ha senso solo se il segnale e' piatto in quella regione!!
804 //
805 // int n1 = (int)(tsta / fChannelWidth); // Non molto preciso, ma tant'e'...
806 // int n2 = (int)(tsto / fChannelWidth);
807 //
808 // return FindSigma2(n1, n2);
809 // }
810 //
811 // double KVSignal::FindSigma2(int tsta, int tsto)
812 // {
813 // // Calcolo della varianza nel tratto tra tsta e tsto.
814 // // NOTA: questo ha senso solo se il segnale e' piatto in quella regione!!
815 //
816 // int n1 = (int)(tsta);
817 // int n2 = (int)(tsto);
818 //
819 // int N = n2 - n1 + 1;
820 // double sigma2 = 0;
821 // double media = FindMedia(tsta, tsto);
822 // if (media == -1E10) {
823 // printf("--- FSignal::FindSigma2(double tsta, double tsto) ---: errore nella media\n");
824 // return -1;
825 // }
826 //
827 // for (int i = n1; i <= n2; i++)
828 // sigma2 += (media - fAdc.At(i)) * (media - fAdc.At(i));
829 //
830 // sigma2 /= N-1;
831 //
832 // return sigma2;
833 // }
834 
835 
836 
838 
839 void KVSignal::FIR_ApplyTrapezoidal(double trise, double tflat) // trise=sqrt(12)*tausha di CR-RC^4 se tflat=trise/2
840 {
841  if (tflat < 0) tflat = trise / 2.;
842  int irise = (int)(1e3 * trise / fChannelWidth);
843  int iflat = (int)(1e3 * tflat / fChannelWidth);
844 
845 // Info("FIR_ApplyTrapezoidal","irise %d iflat %d chw %lf",irise,iflat, fChannelWidth);
846 
847  TArrayF sorig(fAdc);
848  float* data = fAdc.GetArray();
849  float* datao = sorig.GetArray();
850  int N = sorig.GetSize();
851 
852  // 1
853  for (int i = irise; i < N; i++) data[i] -= datao[i - irise];
854  for (int i = irise + iflat; i < N; i++) data[i] -= datao[i - irise - iflat];
855  for (int i = 2 * irise + iflat; i < N; i++) data[i] += datao[i - 2 * irise - iflat];
856 
857  // normalizzazione
858  double amp = 1e3 * trise / fChannelWidth;
859  data[0] /= amp;
860  for (int i = 1; i < N; i++) data[i] = data[i] / amp + data[i - 1];
861 }
862 
863 
864 
865 //Double_t KVSignal::GetAmplitude()
866 //{
867 // //Compute and return the absolute value of the signal amplitude
868 // int i_max=0;
869 // for(int i=0;i<fAdc.GetSize();i++)
870 // {
871 // if(fAdc.At(i_max) < fAdc.At(i) ) i_max = i;
872 // }
873 // fIMax = i_max;
874 // fTMax = fIMax*fChannelWidth;
875 // fAmplitude = fAdc.At(i_max);
876 // return fAmplitude;
877 //}
878 
879 
880 
884 
885 Double_t KVSignal::FindTzeroCFDCubic(double level, int Nrecurr)
886 {
887  // recurr=1 means: linear + 1 approx
888  // recurr=0 == FindTzeroCFD
889  float* data = fAdc.GetArray();
890  int NSamples = fAdc.GetSize();
891  double fmax = level * fAmplitude;
892  /**** 1) ricerca del punto x2 tale che x2 e' il precedente di tcfd ****/
893  int x2 = 0;
894  for (; x2 < NSamples; x2++) if (data[x2] < fmax) break;
895  x2--;
896 
897  return CubicInterpolation(data, x2, fmax, Nrecurr);
898 }
899 
900 
901 /***************************************************************************************************/
902 
904 
905 inline unsigned int ReverseBits(unsigned int p_nIndex, unsigned int p_nBits)
906 {
907 
908  unsigned int i, rev;
909 
910  for (i = rev = 0; i < p_nBits; i++) {
911  rev = (rev << 1) | (p_nIndex & 1);
912  p_nIndex >>= 1;
913  }
914 
915  return rev;
916 
917 }
918 
919 
920 
921 /***************************************************************************************************/
922 
925 
926 void KVSignal::ApplyWindowing(int window_type) // 0: barlett, 1:hanning, 2:hamming, 3: blackman
927 {
928  // vedi pag. 468 oppenheim-shafer
929  int N = fAdc.GetSize() - 1;
930  switch (window_type) {
931  case 0:/*-------------------- BARLETT ----------------------*/
932  for (int n = 0; n <= N; n++)
933  fAdc.GetArray()[n] *= (n < N / 2 ? 2.*n / (double)N : 2. - 2.*n / (double)N);
934  break;
935  case 1:/*-------------------- HANNING ----------------------*/
936  for (int n = 0; n <= N; n++)
937  fAdc.GetArray()[n] *= 0.5 - 0.5 * cos(2 * M_PI * n / (double)N);
938  break;
939  case 2:/*-------------------- HAmmING ----------------------*/
940  for (int n = 0; n <= N; n++)
941  fAdc.GetArray()[n] *= 0.54 - 0.46 * cos(2 * M_PI * n / (double)N);
942  break;
943  case 3:/*-------------------- BLACKMAN --------------------*/
944  for (int n = 0; n <= N; n++)
945  fAdc.GetArray()[n] *= 0.42 - 0.5 * cos(2 * M_PI * n / (double)N) + 0.08 * cos(4 * M_PI * n / (double)N);
946  break;
947  default:
948  printf("ERROR in %s: windowtype %d not valid!\n", __PRETTY_FUNCTION__, window_type);
949  };
950 }
951 
952 
953 
954 
955 
956 /***************************************************************************************************/
957 
961 
962 int KVSignal::FFT(unsigned int p_nSamples, bool p_bInverseTransform,
963  double* p_lpRealIn, double* p_lpImagIn,
964  double* p_lpRealOut, double* p_lpImagOut)
965 {
966  // copiata e adattata da: http://www.codeproject.com/audio/waveInFFT.asp
967  // l'unico vettore che puo' essere NULL e' p_lpImagIn
968 
969  if (!p_lpRealIn || !p_lpRealOut || !p_lpImagOut) {
970  printf("ERROR in %s: NULL vectors!\n", __PRETTY_FUNCTION__);
971  return -1;
972  }
973 
974 
975  unsigned int NumBits;
976  unsigned int i, j, k, n;
977  unsigned int BlockSize, BlockEnd;
978 
979  double angle_numerator = 2.0 * M_PI;
980  double tr, ti;
981 
982  // if( !IsPowerOfTwo(p_nSamples) )
983  // {
984  // return;
985  // }
986  if (p_nSamples < 2 || p_nSamples & (p_nSamples - 1)) {
987  printf("ERROR in %s: %d not a power of two!\n", __PRETTY_FUNCTION__, p_nSamples);
988  return -1;
989  }
990 
991  if (p_bInverseTransform) angle_numerator = -angle_numerator;
992 
993  NumBits = 0; // NumberOfBitsNeeded ( p_nSamples );
994  for (NumBits = 0; ; NumBits++) {
995  if (p_nSamples & (1 << NumBits)) break;
996  }
997 
998  for (i = 0; i < p_nSamples; i++) {
999  j = ReverseBits(i, NumBits);
1000  p_lpRealOut[j] = p_lpRealIn[i];
1001  p_lpImagOut[j] = (p_lpImagIn == NULL) ? 0.0 : p_lpImagIn[i];
1002  }
1003 
1004 
1005  BlockEnd = 1;
1006  for (BlockSize = 2; BlockSize <= p_nSamples; BlockSize <<= 1) {
1007  double delta_angle = angle_numerator / (double)BlockSize;
1008  double sm2 = sin(-2 * delta_angle);
1009  double sm1 = sin(-delta_angle);
1010  double cm2 = cos(-2 * delta_angle);
1011  double cm1 = cos(-delta_angle);
1012  double w = 2 * cm1;
1013  double ar[3], ai[3];
1014 
1015  for (i = 0; i < p_nSamples; i += BlockSize) {
1016 
1017  ar[2] = cm2;
1018  ar[1] = cm1;
1019 
1020  ai[2] = sm2;
1021  ai[1] = sm1;
1022 
1023  for (j = i, n = 0; n < BlockEnd; j++, n++) {
1024 
1025  ar[0] = w * ar[1] - ar[2];
1026  ar[2] = ar[1];
1027  ar[1] = ar[0];
1028 
1029  ai[0] = w * ai[1] - ai[2];
1030  ai[2] = ai[1];
1031  ai[1] = ai[0];
1032 
1033  k = j + BlockEnd;
1034  tr = ar[0] * p_lpRealOut[k] - ai[0] * p_lpImagOut[k];
1035  ti = ar[0] * p_lpImagOut[k] + ai[0] * p_lpRealOut[k];
1036 
1037  p_lpRealOut[k] = p_lpRealOut[j] - tr;
1038  p_lpImagOut[k] = p_lpImagOut[j] - ti;
1039 
1040  p_lpRealOut[j] += tr;
1041  p_lpImagOut[j] += ti;
1042 
1043  }
1044  }
1045 
1046  BlockEnd = BlockSize;
1047 
1048  }
1049 
1050 
1051  if (p_bInverseTransform) {
1052  double denom = (double)p_nSamples;
1053 
1054  for (i = 0; i < p_nSamples; i++) {
1055  p_lpRealOut[i] /= denom;
1056  p_lpImagOut[i] /= denom;
1057  }
1058  }
1059  return 0;
1060 }
1061 
1062 
1063 
1064 
1065 
1066 
1067 
1070 
1071 int KVSignal::FFT(bool p_bInverseTransform, double* p_lpRealOut, double* p_lpImagOut)
1072 {
1073  // returns the lenght of FFT( power of 2)
1074  static double* buffer = NULL;
1075  static int bufflen = 0;
1076  int NSA = fAdc.GetSize();
1077  int ibits = (int)ceil(log((double)NSA) / LOG2);
1078  NSA = 1 << ibits;
1079 
1080  if (buffer != NULL && bufflen < NSA) {
1081  delete [] buffer;
1082  buffer = NULL;
1083  }
1084 
1085 
1086  if (buffer == NULL) {
1087  bufflen = NSA;
1088  buffer = new double [NSA];
1089  }
1090  unsigned int N = fAdc.GetSize();
1091  float* data = fAdc.GetArray();
1092  for (unsigned int i = 0; i < N; i++)
1093  buffer[i] = data[i];
1094  // 0 padding
1095  memset(buffer + N, 0, (NSA - N)*sizeof(double));
1096  int r = FFT(NSA, p_bInverseTransform, buffer, NULL, p_lpRealOut, p_lpImagOut);
1097  if (r < 0) return r;
1098  return NSA;
1099 }
1100 
1101 
1103 
1104 TH1* KVSignal::FFT2Histo(int output, TH1* hh) // 0 modulo, 1 modulo db (normalized), 2, re, 3 im
1105 {
1106  unsigned int N = fAdc.GetSize();
1107  double* re = new double [2 * N];
1108  double* im = new double [2 * N];
1109  int NFFT = FFT(0, re, im);
1110  if (NFFT < 0) {
1111  printf("ERROR in %s: FFT returned %d!\n", __PRETTY_FUNCTION__, NFFT);
1112  delete [] re;
1113  delete [] im;
1114  return NULL;
1115  }
1116  int NF = NFFT / 2;
1117  TH1* h = 0;
1118  if (!hh) h = new TH1F("hfft", "FFT of FSignal", NF, 0, 1. / fChannelWidth * 1000 / 2);
1119  else h = hh;
1120 
1121  h->SetStats(kFALSE);
1122  for (int i = 0; i < NF; i++) {
1123  switch (output) {
1124  case 0: // modulo
1125  h->Fill(h->GetBinCenter(i + 1), sqrt(re[i]*re[i] + im[i]*im[i]));
1126  break;
1127  case 1: // modulo db
1128  h->Fill(h->GetBinCenter(i + 1), log(sqrt(re[i]*re[i] + im[i]*im[i])) * 8.68588963806503500e+00); // numero=20./log(10.)
1129  break;
1130  case 2:
1131  h->Fill(h->GetBinCenter(i + 1), re[i]);
1132  break;
1133  case 3:
1134  h->Fill(h->GetBinCenter(i + 1), im[i]);
1135  break;
1136  default:
1137  printf("ERROR in %s: output=%d not valid!!\n", __PRETTY_FUNCTION__, output);
1138  break;
1139  }
1140  }
1141 // h->GetXaxis()->SetTitle("Frequency");
1142  delete [] re;
1143  delete [] im;
1144 
1145  if (output != 1) return h;
1146  /*** normalizzazione a 0 db ****/
1147  int imax = 0;
1148  for (int i = 1; i < NF; i++)
1149  if (h->GetBinContent(i + 1) > h->GetBinContent(imax + 1))
1150  imax = i;
1151  double dbmax = h->GetBinContent(imax + 1);
1152  for (int i = 0; i < NF; i++)
1153  h->SetBinContent(i + 1, h->GetBinContent(i + 1) - dbmax);
1154  h->GetYaxis()->SetTitle("Modulo (dB)");
1155  return h;
1156 }
1157 
1158 
1162 
1163 Double_t KVSignal::CubicInterpolation(float* data, int x2, double fmax, int Nrecurr)
1164 {
1165  /*
1166  NOTA: tutti i conti fatti con i tempi in "canali". aggiungero' il *fChannelWidth
1167  solo nel return.
1168  */
1169  /**** 2) normal CFD ***************************************************/
1170  double xi0 = (fmax - data[x2]) / (data[x2 + 1] - data[x2]);
1171  if (Nrecurr <= 0) return fChannelWidth * (x2 + xi0);
1172 
1173  /**** 3) approx successive. *******************************************/
1174  // scrivo il polinomio come a3*xi**3+a2*xi**2+a1*xi+a0
1175  // dove xi=tcfd-x2 (ovvero 0<xi<1)
1176  // con maple:
1177  // 3
1178  // (1/2 y2 - 1/6 y1 + 1/6 y4 - 1/2 y3) xi
1179  //
1180  // 2
1181  // + (-y2 + 1/2 y3 + 1/2 y1) xi
1182  //
1183  // + (- 1/2 y2 - 1/6 y4 + y3 - 1/3 y1) xi + y2
1184 
1185  double a3 = 0.5 * data[x2] - (1. / 6.) * data[x2 - 1] + (1. / 6.) * data[x2 + 2] - 0.5 * data[x2 + 1];
1186  double a2 = (-data[x2] + 0.5 * data[x2 + 1] + 0.5 * data[x2 - 1]);
1187  double a1 = (- 0.5 * data[x2] - 1. / 6. *data[x2 + 2] + data[x2 + 1] - 1. / 3.* data[x2 - 1]);
1188  double a0 = data[x2];
1189  double xi = xi0;
1190  for (int rec = 0; rec < Nrecurr; rec++) {
1191  xi += (fmax - a0 - a1 * xi - a2 * xi * xi - a3 * xi * xi * xi) / (a1 + 2 * a2 * xi + 3 * a3 * xi * xi);
1192  }
1193  return fChannelWidth * (x2 + xi);
1194 }
1195 
1196 
1197 
1198 
1200 
1201 double KVSignal::GetDataInter(double t)
1202 {
1203  if (fAdc.GetSize() <= 0) return 1E10;
1204 
1205  int n = (int)(floor(t / fChannelWidth));
1206  if (n <= 0) return fAdc.At(0);
1207  if (n > fAdc.GetSize() - 2) return fAdc.At(fAdc.GetSize() - 1);
1208  if (n * fChannelWidth == t) return fAdc.At(n);
1209  double y1 = fAdc.At(n);
1210  double y2 = fAdc.At(n + 1); //quello prima e quello dopo.
1211  double x1 = fChannelWidth * n;
1212 
1213  return (t - x1) / fChannelWidth * (y2 - y1) + y1;
1214 }
1215 
1216 
1217 
1219 
1221 {
1222  int x2 = (int)(t / fChannelWidth);
1223  if (x2 < 1 || x2 > fAdc.GetSize() - 2) return GetDataInter(t);
1224  float* data = fAdc.GetArray();
1225  /***** CUT & PASTE DA CubicInterpolation *****/
1226 
1227  double a3 = 0.5 * data[x2] - (1. / 6.) * data[x2 - 1] + (1. / 6.) * data[x2 + 2] - 0.5 * data[x2 + 1];
1228  double a2 = (-data[x2] + 0.5 * data[x2 + 1] + 0.5 * data[x2 - 1]);
1229  double a1 = (- 0.5 * data[x2] - 1. / 6. *data[x2 + 2] + data[x2 + 1] - 1. / 3.* data[x2 - 1]);
1230  double a0 = data[x2];
1231  double xi = (t / fChannelWidth - x2);
1232  return a3 * xi * xi * xi + a2 * xi * xi + a1 * xi + a0;
1233 }
1234 
1235 
1236 /***********************************************/
1237 
1240 
1242 {
1243  // see HSIEH S.HOU IEEE Trans. Acoustic Speech, vol. ASSP-26, NO.6, DECEMBER 1978
1245  Double_t xk = floor(t / h) * h;
1246  int xk_index = (int)(t / h);
1247  Double_t s = (t - xk) / h;
1248 
1249  float* data = fAdc.GetArray();
1250  int N = fAdc.GetSize();
1251  int dimensione = 18; //dimensione della matrice dei campioni.!!!!deve essere pari!!!!
1252  float cm1, cNm1;
1253 
1254 
1255  TMatrixD e(dimensione, dimensione);
1256  TArrayD data_e(dimensione * dimensione);
1257  for (int k = 0, i = 0; i < dimensione; i++) {
1258  data_e[k] = 4.;
1259  if ((k + 1) < pow(dimensione, 2)) data_e[k + 1] = 1.;
1260  if ((k - 1) > 0) data_e[k - 1] = 1.;
1261  k += dimensione + 1;
1262  }
1263  e.SetMatrixArray(data_e.GetArray());
1264  e *= 1. / 6 / h;
1265  e.Invert();
1266 
1267  double dati_b[] = { -1, 3, -3, 1, 3, -6, 3, 0, -3, 0, 3, 0, 1, 4, 1, 0};
1268  TMatrixD delta(4, 4, dati_b);
1269 
1270 
1271  TMatrixF samples(dimensione, 1);
1272  TMatrixF coeff(dimensione, 1);
1273  TMatrixF b(4, 1);
1274  TMatrixF coefficienti(4, 1);
1275 
1276  if (xk_index < (dimensione / 2 - 1)) {
1277  samples.SetMatrixArray(data);//prendiamo i dati a partire dal primo
1278  }
1279  else if (xk_index > (N - dimensione / 2 - 1)) {
1280  samples.SetMatrixArray(data + N - dimensione); //prendiamo 18 dati partendo dalla fine
1281  }
1282  else {
1283  samples.SetMatrixArray(data + xk_index - (dimensione / 2 - 1)); //perche l'interp deve avvenire tra i campioni centrali della matrice
1284  }
1285  float* samples_interp = samples.GetMatrixArray();
1286 
1287  // coeff=e*samples;
1288  coeff.Mult(e, samples);
1289  float* ck = coeff.GetMatrixArray();
1290 
1291  if (xk_index < (dimensione / 2 - 1)) {
1292  if (xk_index == 0) {
1293  cm1 = (*samples_interp) * 6 * h - ((*ck) * 4 + (*(ck + 1)));
1294  float caso_zero[4] = {cm1, *ck, *(ck + 1), *(ck + 2)};
1295  coefficienti.SetMatrixArray(caso_zero);
1296  }
1297  else coefficienti.SetMatrixArray(ck + xk_index - 1);
1298  }
1299  else if (xk_index > (N - dimensione / 2 - 1)) {
1300  if (xk_index == N - 2) {
1301  cNm1 = (*(samples_interp + dimensione - 1)) * 6 * h - (*(ck + dimensione - 1) * 4 + (*(ck + dimensione - 2)));
1302  float casoN[4] = {*(ck + dimensione - 3), *(ck + dimensione - 2), *(ck + dimensione - 1), cNm1};
1303  coefficienti.SetMatrixArray(casoN);
1304  }
1305  else coefficienti.SetMatrixArray(ck + dimensione - (N - xk_index + 1));
1306  }
1307  else {
1308  coefficienti.SetMatrixArray(ck + (dimensione / 2 - 2));
1309  }
1310 
1311  // b=delta*coefficienti;
1312  b.Mult(delta, coefficienti);
1313  float* b_interp = b.GetMatrixArray();
1314  float a0 = *b_interp;
1315  float a1 = *(b_interp + 1);
1316  float a2 = *(b_interp + 2);
1317  float a3 = *(b_interp + 3);
1318 
1319  return (1. / 6 / h * (a3 + a2 * s + a1 * s * s + a0 * s * s * s));
1320 }
1321 
1322 
1323 
1324 /***********************************************/
1325 
1327 
1329 {
1330  const int Nsa = fAdc.GetSize();
1331  const double tau = fChannelWidth;
1332 
1333  fChannelWidthInt = taufinal;
1334  TArrayF interpo;
1335  interpo.Set((int)(Nsa * tau / taufinal));
1336  int nlast = interpo.GetSize() - (int)(3 * tau / taufinal);
1337  if (nlast <= 0) return;
1338 
1339  for (int i = 0; i < interpo.GetSize(); i++) interpo.AddAt(GetDataCubicSpline(i * taufinal), i);
1340  fAdc.Set(0);
1341  fAdc.Set(interpo.GetSize());
1342  for (int i = 0; i < interpo.GetSize(); i++) fAdc.AddAt(interpo.At(i), i);
1343  for (int i = nlast; i < interpo.GetSize(); i++) fAdc.AddAt(interpo.At(nlast), i);
1344 
1345 }
1346 
1347 /***********************************************/
1348 
1350 
1352 {
1354 }
1355 
1356 
1357 
1358 
1359 /***********************************************/
1360 
1362 
1363 void KVSignal::BuildCubicSignal(double taufinal)
1364 {
1365  const int Nsa = fAdc.GetSize();
1366  const double tau = fChannelWidth;
1367 
1368  fChannelWidthInt = taufinal;
1369  TArrayF interpo;
1370  interpo.Set((int)(Nsa * tau / taufinal));
1371  int nlast = interpo.GetSize() - (int)(3 * tau / taufinal);
1372  if (nlast <= 0) return;
1373 
1374  for (int i = 0; i < interpo.GetSize(); i++) interpo.AddAt(GetDataInterCubic(i * taufinal), i);
1375  fAdc.Set(0);
1376  fAdc.Set(interpo.GetSize());
1377  for (int i = 0; i < nlast; i++) fAdc.AddAt(interpo.At(i), i);
1378  for (int i = nlast; i < interpo.GetSize(); i++) fAdc.AddAt(interpo.At(nlast), i);
1379 
1380 }
1381 
1382 /***********************************************/
1383 
1385 
1387 {
1389 }
1390 
1391 
1392 /***********************************************/
1393 
1395 
1396 void KVSignal::BuildSmoothingSplineSignal(double taufinal, double l, int nbits)
1397 {
1398  const int Nsa = fAdc.GetSize();
1399  const double tau = fChannelWidth;
1400 
1401  KVSignal coeff;
1402  ApplyModifications(&coeff);
1403  coeff.SetADCData();
1404  if (coeff.FIR_ApplySmoothingSpline(l, nbits) != 0) return;
1405 
1406  fChannelWidthInt = taufinal;
1407  TArrayF interpo;
1408  interpo.Set((int)(Nsa * tau / taufinal));
1409  int nlast = interpo.GetSize() - (int)(3 * tau);
1410  if (nlast <= 0) return;
1411 
1412  for (int i = 0; i < interpo.GetSize() - (int)(53 * tau / taufinal); i++) interpo.AddAt(coeff.GetDataSmoothingSplineLTI(i * taufinal), i);
1413  fAdc.Set(0);
1414  fAdc.Set(interpo.GetSize());
1415  for (int i = 0; i < nlast; i++) fAdc.AddAt(interpo.At(i), i);
1416  for (int i = nlast; i < interpo.GetSize(); i++) fAdc.AddAt(interpo.At(nlast), i);
1417 
1418 }
1419 
1420 /***********************************************/
1421 
1423 
1425 {
1427 }
1428 
1429 /***********************************************/
1430 
1431 
1436 
1437 int KVSignal::FIR_ApplySmoothingSpline(double l, int nbits)
1438 {
1439  // This method is never called with nbits>2, therefore nmax=50 ALWAYS
1440  // and dynamic arrays xvec & yvec can safely be of fixed size
1441  // If ever we want to use nbits>2, this will have to be changed
1442 
1443  if (l < 0.1) return -1;
1444  int i;
1445  double x0 = ApplyNewton(l, 1000);
1446  double r = sqrt(x0);
1447  double a = (1 - 24 * l) / (6 * l);
1448  double cphi = -a * r / (2 * (r * r + 1));
1449  if ((cphi < -1) || (cphi > 1) || (l < 0.1)) {
1450  printf("Error Aborting on FIR_ApplySmoothingSpline\n");
1451  return -1;
1452  }
1453  double sphi = sqrt(1 - cphi * cphi);
1454  double Re = sqrt(x0) * cphi;
1455  double Im = sqrt(x0) * sphi;
1456  double a1, a2, b1, b2, t;
1457  t = 4 * Re * Re * (x0 - 1) / (x0 + 1) + 1 - x0 * x0;
1458  b2 = 1 / (l * t);
1459  a2 = -2 * Re * x0 * b2 / (1 + x0);
1460  b1 = -x0 * x0 * b2;
1461  a1 = -a2;
1462  double czr, czi;
1463  czr = a1 / 2;
1464  czi = (-a1 * Re - b1) / (2 * Im);
1465  double phib, phiz;
1466  phiz = atan2(sphi, cphi);
1467  phib = atan2(czi, czr);
1468  double roB;
1469  roB = sqrt(czr * czr + czi * czi);
1470  double roZ = sqrt(x0);
1471  int nfloat = 0;
1472  int nmax(50);
1473  double imax, fmax;
1474  if (nbits > 2) {
1475  //Determination of best notation
1476  //1 bit always for the sign
1477  //#integer bit depends on output evaluated in 0//
1478  //#fixed to 0, l>0.1
1479  nfloat = nbits - 1;
1480  //1 represented as 2^(nfloats)...//
1481  imax = ((nbits + 1) * log(2) + log(roB)) / log(roZ) - 1.;
1482  nmax = floor(imax);
1483  do {
1484  fmax = std::abs(-2 * roB * cos(phib - (nmax + 1) * phiz) / pow(roZ, nmax + 1));
1485  if (fmax * pow(2, nfloat + 1) < 1) nmax--;
1486  }
1487  while (fmax * pow(2, nfloat + 1) < 1);
1488  }
1489  double* xvec = new double[2 * nmax + 1];
1490  double* yvec = new double[2 * nmax + 1];
1491  if (nbits > 2) {
1492  for (i = 0; i <= nmax; i++) {
1493  yvec[nmax + i] = yvec[nmax - i] = 0;
1494  xvec[nmax + i] = xvec[nmax - i] = ((double)floor(-2 * roB * cos(phib - (i + 1) * phiz) / pow(roZ, i + 1) * pow(2, nfloat) + 0.5)) / pow(2, nfloat);
1495  }
1496  }
1497  else {
1498  for (i = 0; i <= nmax; i++) {
1499  yvec[nmax + i] = yvec[nmax - i] = 0;
1500  xvec[nmax + i] = xvec[nmax - i] = -2 * roB * cos(phib - (i + 1) * phiz) / pow(roZ, i + 1);
1501  }
1502  }
1503  FIR_ApplyRecursiveFilter(0, 2 * nmax + 1, xvec, yvec, 0);
1504  ShiftLeft(nmax * fChannelWidth);
1505  delete [] xvec;
1506  delete [] yvec;
1507  return 0;
1508 }
1509 
1510 
1511 /***********************************************/
1512 
1514 
1515 double KVSignal::ApplyNewton(double l, double X0)
1516 {
1517  double a = (1 - 24 * l) / (6 * l);
1518  double b = (4 + 36 * l) / (6 * l);
1519  double A = 2 - b;
1520  double B = 2 + a * a - 2 * b;
1521  Double_t x1, x0, x2, x3;
1522  double der, fx;
1523  double fxold1;
1524  x0 = X0;
1525  x1 = X0;
1526  x2 = x0;
1527  fx = pow(x1, 4) + A * pow(x1, 3) + B * pow(x1, 2) + A * x1 + 1;
1528  fxold1 = fx + 1;
1529  int fexit = 0;
1530  int loopcount = 0;
1531  do {
1532  der = 4 * pow(x1, 3) + 3 * A * pow(x1, 2) + 2 * B * x1 + A;
1533  x3 = x2;
1534  x2 = x0;
1535  x0 = x1;
1536  x1 = x0 - fx / der;
1537  if (x1 == x2) {
1538  x1 = (x1 + x0) / 2.;
1539  loopcount++;
1540  }
1541  else if (x1 == x3) {
1542  x1 = (x1 + x0 + x2) / 3.;
1543  loopcount++;
1544  }
1545  if ((x1 == x0) && (x1 == x2) && (x1 == x3)) fexit = 1;
1546  if (loopcount == 100) fexit = 1;
1547  fxold1 = fx;
1548  fx = pow(x1, 4) + A * pow(x1, 3) + B * pow(x1, 2) + A * x1 + 1;
1549  }
1550  while (((fx > 0.000000001) || (fxold1 != fx)) && (fexit == 0));
1551  return x1;
1552 }
1553 
1554 /***********************************************/
1555 
1557 
1559 {
1560  double tsamp = t / fChannelWidth;
1561  int n0 = floor(tsamp);
1562  double tres = tsamp - n0;
1563  double res = 0;
1564  if (n0 == 0)return 0;
1565  else {
1566  for (int i = -1; i < 3; i++) {
1567  res += fAdc.At(n0 + i) * EvalCubicSpline(tres - i);
1568  }
1569  }
1570  return res;
1571 }
1572 
1573 /***********************************************/
1574 
1576 
1578 {
1579  double ax = TMath::Abs(X);
1580  if (ax > 2) return 0;
1581  else if (ax > 1) return pow(2 - ax, 3) / 6.;
1582  else return pow(ax, 3) / 2. - ax * ax + 2. / 3.;
1583 }
1584 
1585 /***********************************************/
1586 
1590 
1591 double KVSignal::FindTzeroCFDCubic_rev(double level, double tend, int Nrecurr)
1592 {
1593  // recurr=1 means: linear + 1 approx
1594  // recurr=0 == FindTzeroCFD
1595  /**************************************/
1596  float* data = fAdc.GetArray();
1597  int NSamples = fAdc.GetSize();
1598  double fmax = level * fAmplitude;
1599  /**** 1) ricerca del punto x2 tale che x2 e' il precedente di tcfd ****/
1600  int x2 = (int)(tend / fChannelWidth);
1601  if (x2 <= 0 || x2 >= NSamples) return -1;
1602  for (; x2 > 0 ; x2--)
1603  if (data[x2] > fmax)
1604  break;
1605  // x2--;
1606  return CubicInterpolation(data, x2, fmax, Nrecurr);
1607 }
1608 
1609 
1610 
1612 
1613 void KVSignal::FIR_ApplySemigaus(double tau_usec)
1614 {
1615  FIR_ApplyRCLowPass(tau_usec);
1616  FIR_ApplyRCLowPass(tau_usec);
1617  FIR_ApplyRCLowPass(tau_usec);
1618  FIR_ApplyRCLowPass(tau_usec);
1619 
1620  FIR_ApplyRCHighPass(tau_usec);
1621 
1622  Multiply(1. / (32. / 3.*TMath::Exp(-4.))); // normalizzazione
1623 }
1624 
1625 
1626 
1627 
1631 
1632 void KVSignal::FIR_ApplyRCLowPass(double time_usec, int reverse)
1633 {
1634  // f=1/(2*pi*tau) se tau[ns] allora f->[GHz]
1635  // fsampling=1/fChannelWidth [GHz]
1636 #define DUEPI 6.28318530717958623
1637  // printf("fChannelWidth=%f\n",fChannelWidth);
1638  double f = 1. / (DUEPI * time_usec * 1000.) * fChannelWidth;
1639  double x = TMath::Exp(-DUEPI * f);
1640  double a0 = 1 - x;
1641  double b = x;
1642  double a = 0;
1643  // printf("f=%f, x=%f\n",x,f);
1644  FIR_ApplyRecursiveFilter(a0, 1, &a, &b, reverse);
1645 }
1646 
1647 
1648 
1653 
1654 void KVSignal::FIR_ApplyRCHighPass(double time_usec, int reverse)
1655 {
1656  // f=1/(2*pi*tau) se tau[ns] allora f->[GHz]
1657  // fsampling=1/fChannelWidth [GHz]
1658 
1659  // printf("fChannelWidth=%f\n",fChannelWidth);
1660  double f = 1. / (DUEPI * time_usec * 1000.) * fChannelWidth;
1661  double x = TMath::Exp(-DUEPI * f);
1662  double a0 = (1 + x) / 2.;
1663  double a1 = -(1 + x) / 2.;
1664  double b1 = x;
1665  // printf("f=%f, x=%f\n",x,f);
1666  FIR_ApplyRecursiveFilter(a0, 1, &a1, &b1, reverse);
1667 }
1668 
1669 #undef DUEPI
1670 
1671 
1674 
1675 void KVSignal::FIR_ApplyRecursiveFilter(double a0, int N, double* a, double* b, int reverse)
1676 {
1677  // signal will be: y[n]=a0*x[n]+sum a[k] x[k] + sum b[k] y[k]
1678  int NSamples = fAdc.GetSize();
1679  double* datay = new double[NSamples];
1680  float* datax = fAdc.GetArray();
1681  // memset(datay, 0, NSamples*sizeof(float)); //azzero l'array.
1682  /*----------------------------------------------*/
1683  int i = 0, k = 0;
1684  switch (reverse) {
1685  case 0:// direct
1686  for (i = 0; i < N; i++) { // primo loop su Npunti
1687  datay[i] = a0 * datax[i];
1688  for (k = 0; k < i; k++)
1689  datay[i] += a[k] * datax[i - k - 1] + b[k] * datay[i - k - 1];
1690  }
1691  for (i = N; i < NSamples; i++) { //secondo loop. cambia l'indice interno.
1692  datay[i] = a0 * datax[i];
1693  for (k = 0; k < N; k++)
1694  datay[i] += a[k] * datax[i - k - 1] + b[k] * datay[i - k - 1];
1695  }
1696  break; // end of direct
1697  case 1: //reverse: cut & paste from direct and NSamples-1-
1698  for (i = 0; i < N; i++) { // primo loop su Npunti
1699  datay[NSamples - 1 - i] = a0 * datax[NSamples - 1 - i];
1700  for (k = 0; k < i; k++)
1701  datay[NSamples - 1 - i] += a[k] * datax[NSamples - 1 - (i - k - 1)]
1702  + b[k] * datay[NSamples - 1 - (i - k - 1)];
1703  }
1704  for (i = N; i < NSamples; i++) { //secondo loop. cambia l'indice interno.
1705  datay[NSamples - 1 - i] = a0 * datax[NSamples - 1 - i];
1706  for (k = 0; k < N; k++)
1707  datay[NSamples - 1 - i] += a[k] * datax[NSamples - 1 - (i - k - 1)]
1708  + b[k] * datay[NSamples - 1 - (i - k - 1)];
1709  }
1710  break;
1711  case -1: // bidirectional
1712  FIR_ApplyRecursiveFilter(a0, N, a, b, 0);
1713  FIR_ApplyRecursiveFilter(a0, N, a, b, 1);
1714  delete [] datay;
1715  return;
1716  default:
1717  printf("ERROR in %s: reverse=%d not supported\n", __PRETTY_FUNCTION__, reverse);
1718  }// end of reverse switch.
1719  /*----------------------------------------------*/
1720  // void *memcpy(void *dest, const void *src, size_t n);
1722  for (int i = 0; i < NSamples; i++)
1723  datax[i] = (float)datay[i];
1724  delete [] datay;
1725 }
1726 
1727 
1728 
1730 
1731 void KVSignal::FIR_ApplyMovingAverage(int npoints) // causal moving average
1732 {
1733  TArrayF sorig = fAdc;
1734  float* data = fAdc.GetArray();
1735  float* datao = sorig.GetArray();
1736 
1737  for (int n = npoints; n < GetNSamples(); n++)
1738  data[n] = data[n - 1] + (datao[n] - datao[n - npoints]) / npoints;
1739  for (int n = 0; n < npoints; n++) data[n] = data[npoints]; // first npoints samples are put equal to first good value
1740 
1741 }
1742 
1743 
1744 
1745 
1747 
1749 {
1750 
1751  KVDigitalFilter deconv_pa1, lp1, int1;
1754  deconv_pa1 = KVDigitalFilter::CombineStagesMany(&lp1, &int1);
1755 
1756  deconv_pa1.ApplyTo(this);
1757  Multiply(fChannelWidth / tauRC / 1000.);
1758 }
1759 
1760 
1761 
1764 
1766 {
1767 // Info("ApplyModifications","called with %d",((newSignal==0)?0:1));
1768  if (!newSignal) newSignal = this;
1769 
1770  Int_t nn = fAdc.GetSize();
1771  if (nsa > 0 && nsa < nn) nn = nsa;
1772  if (newSignal->InheritsFrom("KVSignal"))((KVSignal*)newSignal)->SetChannelWidth(fChannelWidthInt);
1773  for (int ii = 0; ii < nn; ii++) newSignal->SetPoint(ii, ii * fChannelWidthInt, fAdc.At(ii));
1774 }
1775 
1776 
1777 
1778 
1780 
1782 {
1783  for (int i = 0; i < fAdc.GetSize(); i++) fAdc.AddAt(fAdc.At(i)*fact, i);
1784 }
1785 
1786 
1787 
1789 
1791 {
1792  for (int i = 0; i < fAdc.GetSize(); i++) fAdc.AddAt(fAdc.At(i) + fact, i);
1793 }
1794 
1795 
1796 /***************************************************************************/
1797 
1799 
1800 void KVSignal::ShiftLeft(double tshift)//shift is in ns
1801 {
1802  if (fAdc.GetSize() <= 0) return;
1803  if (tshift < 0) {
1804  ShiftRight(-tshift);
1805  return;
1806  }
1807 
1808  int shift = (int)floor(tshift / fChannelWidth);
1809  if (shift == 0) return;
1810  for (int i = 0; i < fAdc.GetSize() - shift; i++)
1811  fAdc.AddAt(fAdc.At(i + shift), i);
1812  fAdc.Set(fAdc.GetSize() - shift);
1813 }
1814 
1815 /***************************************************************************/
1816 
1818 
1819 void KVSignal::ShiftRight(double tshift)//shift is in ns
1820 {
1821  if (fAdc.GetSize() <= 0) return;
1822  if (tshift < 0) {
1823  ShiftLeft(-tshift);
1824  return;
1825  }
1826 
1827  int shift = (int)floor(tshift / fChannelWidth);
1828  if (shift == 0) return;
1829 
1830  fAdc.Set(fAdc.GetSize() + shift);
1831  for (int i = fAdc.GetSize() - 1; i >= shift; i--)
1832  fAdc.AddAt(fAdc.At(i - shift), i);
1833  for (int i = 0; i < shift; i++)
1834  fAdc.AddAt(fAdc.At(shift), i);
1835 }
1836 
1837 /***************************************************************************/
1838 
1840 
1842 {
1843  this->Draw();
1844  getchar();
1845 }
1846 
1847 
1848 
1851 
1852 KVSignal* KVSignal::MakeSignal(const char* sig_type)
1853 {
1854  // Create new KVSignal instance corresponding to sig_type
1855 
1856  TPluginHandler* ph = KVBase::LoadPlugin("KVSignal", sig_type);
1857  if (!ph) {
1858  ::Error("KVSignal::MakeSignal", "No plugin found for : %s", sig_type);
1859  return nullptr;
1860  }
1861  KVSignal* sig = (KVSignal*)ph->ExecPlugin(0);
1862  return sig;
1863 }
1864 
1865 
1866 
int Int_t
double
KVDataSet * gDataSet
Definition: KVDataSet.cpp:29
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
#define LOG2
Definition: KVSignal.cpp:16
unsigned int ReverseBits(unsigned int p_nIndex, unsigned int p_nBits)
Definition: KVSignal.cpp:905
#define DUEPI
#define M_PI
Definition: KVSignal.cpp:17
ROOT::R::TRInterface & r
#define b(i)
#define f(i)
#define e(i)
char Char_t
const Bool_t kFALSE
bool Bool_t
double Double_t
const Bool_t kTRUE
const char Option_t
R__EXTERN TEnv * gEnv
#define N
int type
double atan2(double, double)
double ceil(double)
double cos(double)
double pow(double, double)
double floor(double)
double sin(double)
double log(double)
char * Form(const char *fmt,...)
static TPluginHandler * LoadPlugin(const Char_t *base, const Char_t *uri="0")
Definition: KVBase.cpp:793
To store calibration parameters in a database ,.
KVNameValueList * GetParameters()
const Char_t * GetDataSetEnv(const Char_t *type, const Char_t *defval="") const
Definition: KVDataSet.cpp:758
static KVDigitalFilter BuildRCLowPassDeconv(const double &tau_usec, const double &tau_clk)
void ApplyTo(double *data, const int N, int reverse=0) const
static KVDigitalFilter BuildIntegrator(const double &tau_clk)
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!
Double_t GetDoubleValue(const Char_t *name) const
const Char_t * GetNameAt(Int_t idx) const
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
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
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
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 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 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
void ShiftLeft(double)
---------------— OPERATORI ------------------—//
Definition: KVSignal.cpp:1800
Double_t fSigmaEnd
rms value of the signal line at the end
Definition: KVSignal.h:46
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
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
KVSignal()
Default constructor.
Definition: KVSignal.cpp:88
void Multiply(Double_t fact)
multiply the signal (modify only fAdc)
Definition: KVSignal.cpp:1781
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
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
Bool_t IsFired()
ComputeBaseLine and ComputeEndLine methods have to be called before.
Definition: KVSignal.cpp:581
void SetADCData()
Definition: KVSignal.cpp:229
void SetChannelWidth(double width)
Definition: KVSignal.h:167
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 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
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
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
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 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
Bool_t fPSAIsDone
indicate if PSA has been done
Definition: KVSignal.h:64
TH1 * FFT2Histo(int output, TH1 *hh=0)
Definition: KVSignal.cpp:1104
void PoleZeroSuppression(Double_t tauRC)
Definition: KVSignal.cpp:1748
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 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
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
void Begin(TString delim) const
Definition: KVString.cpp:565
KVString Next(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:695
Int_t GetNValues(TString delim) const
Definition: KVString.cpp:886
Double_t * GetArray()
Float_t * GetArray()
Double_t GetAt(Int_t i) const
void Set(Int_t n)
Float_t At(Int_t i) const
void AddAt(Float_t c, Int_t i)
Int_t GetSize() const
static TClass * GetClass(Bool_t load=kTRUE, Bool_t silent=kFALSE)
virtual const char * GetValue(const char *name, const char *dflt) const
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
virtual void SetNameTitle(const char *name="", const char *title="")
Int_t GetN() const
Double_t * fY
Double_t * GetX() const
virtual void Draw(Option_t *chopt="")
Double_t * GetY() const
virtual void Set(Int_t n)
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y) const
virtual Double_t GetBinCenter(Int_t bin) const
TAxis * GetYaxis()
virtual Int_t Fill(const char *name, Double_t w)
virtual void SetBinContent(Int_t bin, Double_t content)
virtual Double_t GetBinContent(Int_t bin) const
virtual void SetStats(Bool_t stats=kTRUE)
virtual TMatrixTBase< Element > & SetMatrixArray(const Element *data, Option_t *option="")
virtual Element * GetMatrixArray()
void Mult(const TMatrixT< Element > &a, const TMatrixT< Element > &b)
virtual const char * GetName() const
virtual void SetTitle(const char *title="")
virtual const char * GetTitle() const
virtual void Warning(const char *method, const char *msgfmt,...) const
virtual Bool_t InheritsFrom(const char *classname) const
virtual void Error(const char *method, const char *msgfmt,...) const
virtual void Info(const char *method, const char *msgfmt,...) const
Longptr_t ExecPlugin(int nargs, const T &... params)
Int_t Atoi() const
Bool_t IsDigit() const
void ToUpper()
Bool_t BeginsWith(const char *s, ECaseCompare cmp=kExact) const
const char * Data() const
void Form(const char *fmt,...)
TString & ReplaceAll(const char *s1, const char *s2)
VecExpr< UnaryOp< Sqrt< T >, SVector< T, D >, T >, T, D > sqrt(const SVector< T, D > &rhs)
const Double_t sigma
Double_t x[n]
const Int_t n
TH1 * h
const long double s
Definition: KVUnits.h:94
const long double cl
Definition: KVUnits.h:85
Int_t Nint(T x)
Double_t Exp(Double_t x)
Double_t Sqrt(Double_t x)
Double_t Abs(Double_t d)
Double_t Max(Double_t a, Double_t b)
auto * tt
auto * l
auto * a