1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
17 ///////////////////////////////////////////////////////////////////////////
20 // Subset of matheamtical functions not included in the TMath
23 /////////////////////////////////////////////////////////////////////////
25 #include "Riostream.h"
32 #include "TObjString.h"
33 #include "TLinearFitter.h"
36 #include "TGraphErrors.h"
37 #include "TMultiGraph.h"
42 // includes neccessary for test functions
46 #include "TStopwatch.h"
47 #include "TTreeStream.h"
49 #include "TStatToolkit.h"
57 ClassImp(TStatToolkit) // Class implementation to enable ROOT I/O
59 TStatToolkit::TStatToolkit() : TObject()
62 // Default constructor
65 ///////////////////////////////////////////////////////////////////////////
66 TStatToolkit::~TStatToolkit()
74 //_____________________________________________________________________________
75 void TStatToolkit::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean
76 , Double_t &sigma, Int_t hh)
79 // Robust estimator in 1D case MI version - (faster than ROOT version)
81 // For the univariate case
82 // estimates of location and scatter are returned in mean and sigma parameters
83 // the algorithm works on the same principle as in multivariate case -
84 // it finds a subset of size hh with smallest sigma, and then returns mean and
85 // sigma of this subset
90 Double_t faclts[]={2.6477,2.5092,2.3826,2.2662,2.1587,2.0589,1.9660,1.879,1.7973,1.7203,1.6473};
91 Int_t *index=new Int_t[nvectors];
92 TMath::Sort(nvectors, data, index, kFALSE);
94 Int_t nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
95 Double_t factor = faclts[TMath::Max(0,nquant-1)];
100 Double_t bestmean = 0;
101 Double_t bestsigma = (data[index[nvectors-1]]-data[index[0]]+1.); // maximal possible sigma
102 bestsigma *=bestsigma;
104 for (Int_t i=0; i<hh; i++){
105 sumx += data[index[i]];
106 sumx2 += data[index[i]]*data[index[i]];
109 Double_t norm = 1./Double_t(hh);
110 Double_t norm2 = (hh-1)>0 ? 1./Double_t(hh-1):1;
111 for (Int_t i=hh; i<nvectors; i++){
112 Double_t cmean = sumx*norm;
113 Double_t csigma = (sumx2 - hh*cmean*cmean)*norm2;
114 if (csigma<bestsigma){
120 sumx += data[index[i]]-data[index[i-hh]];
121 sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
124 Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
133 void TStatToolkit::EvaluateUniExternal(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh, Float_t externalfactor)
135 // Modified version of ROOT robust EvaluateUni
136 // robust estimator in 1D case MI version
137 // added external factor to include precision of external measurement
142 Double_t faclts[]={2.6477,2.5092,2.3826,2.2662,2.1587,2.0589,1.9660,1.879,1.7973,1.7203,1.6473};
143 Int_t *index=new Int_t[nvectors];
144 TMath::Sort(nvectors, data, index, kFALSE);
146 Int_t nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
147 Double_t factor = faclts[0];
149 // fix proper normalization - Anja
150 factor = faclts[nquant-1];
157 Int_t bestindex = -1;
158 Double_t bestmean = 0;
159 Double_t bestsigma = -1;
160 for (Int_t i=0; i<hh; i++){
161 sumx += data[index[i]];
162 sumx2 += data[index[i]]*data[index[i]];
165 Double_t kfactor = 2.*externalfactor - externalfactor*externalfactor;
166 Double_t norm = 1./Double_t(hh);
167 for (Int_t i=hh; i<nvectors; i++){
168 Double_t cmean = sumx*norm;
169 Double_t csigma = (sumx2*norm - cmean*cmean*kfactor);
170 if (csigma<bestsigma || bestsigma<0){
177 sumx += data[index[i]]-data[index[i-hh]];
178 sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
181 Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
188 //_____________________________________________________________________________
189 Int_t TStatToolkit::Freq(Int_t n, const Int_t *inlist
190 , Int_t *outlist, Bool_t down)
193 // Sort eleements according occurancy
194 // The size of output array has is 2*n
197 Int_t * sindexS = new Int_t[n]; // temp array for sorting
198 Int_t * sindexF = new Int_t[2*n];
199 for (Int_t i=0;i<n;i++) sindexS[i]=0;
200 for (Int_t i=0;i<2*n;i++) sindexF[i]=0;
202 TMath::Sort(n,inlist, sindexS, down);
203 Int_t last = inlist[sindexS[0]];
210 for(Int_t i=1;i<n; i++){
211 val = inlist[sindexS[i]];
212 if (last == val) sindexF[countPos]++;
215 sindexF[countPos+n] = val;
220 if (last==val) countPos++;
221 // sort according frequency
222 TMath::Sort(countPos, sindexF, sindexS, kTRUE);
223 for (Int_t i=0;i<countPos;i++){
224 outlist[2*i ] = sindexF[sindexS[i]+n];
225 outlist[2*i+1] = sindexF[sindexS[i]];
234 //___TStatToolkit__________________________________________________________________________
235 void TStatToolkit::TruncatedMean(const TH1 * his, TVectorD *param, Float_t down, Float_t up, Bool_t verbose){
239 Int_t nbins = his->GetNbinsX();
240 Float_t nentries = his->GetEntries();
245 for (Int_t ibin=1;ibin<nbins; ibin++){
246 ncumul+= his->GetBinContent(ibin);
247 Float_t fraction = Float_t(ncumul)/Float_t(nentries);
248 if (fraction>down && fraction<up){
249 sum+=his->GetBinContent(ibin);
250 mean+=his->GetBinCenter(ibin)*his->GetBinContent(ibin);
251 sigma2+=his->GetBinCenter(ibin)*his->GetBinCenter(ibin)*his->GetBinContent(ibin);
255 sigma2= TMath::Sqrt(TMath::Abs(sigma2/sum-mean*mean));
257 (*param)[0] = his->GetMaximum();
259 (*param)[2] = sigma2;
262 if (verbose) printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma2);
265 void TStatToolkit::LTM(TH1 * his, TVectorD *param , Float_t fraction, Bool_t verbose){
267 // LTM : Trimmed mean on histogram - Modified version for binned data
269 // Robust statistic to estimate properties of the distribution
270 // See http://en.wikipedia.org/w/index.php?title=Trimmed_estimator&oldid=582847999
272 // New faster version is under preparation
278 Int_t nbins = his->GetNbinsX();
279 Int_t nentries = (Int_t)his->GetEntries();
280 if (nentries<=0) return;
281 Double_t *data = new Double_t[nentries];
283 for (Int_t ibin=1;ibin<nbins; ibin++){
284 Double_t entriesI = his->GetBinContent(ibin);
285 //Double_t xcenter= his->GetBinCenter(ibin);
286 Double_t x0 = his->GetXaxis()->GetBinLowEdge(ibin);
287 Double_t w = his->GetXaxis()->GetBinWidth(ibin);
288 for (Int_t ic=0; ic<entriesI; ic++){
289 if (npoints<nentries){
290 data[npoints]= x0+w*Double_t((ic+0.5)/entriesI);
295 Double_t mean, sigma;
296 Int_t npoints2=TMath::Min(Int_t(fraction*Float_t(npoints)),npoints-1);
297 npoints2=TMath::Max(Int_t(0.5*Float_t(npoints)),npoints2);
298 TStatToolkit::EvaluateUni(npoints, data, mean,sigma,npoints2);
300 if (verbose) printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma);if (param){
301 (*param)[0] = his->GetMaximum();
308 void TStatToolkit::MedianFilter(TH1 * his1D, Int_t nmedian){
310 // Algorithm to filter histogram
311 // author: marian.ivanov@cern.ch
312 // Details of algorithm:
313 // http://en.wikipedia.org/w/index.php?title=Median_filter&oldid=582191524
315 // his1D - input histogam - to be modiefied by Medianfilter
316 // nmendian - number of bins in median filter
318 Int_t nbins = his1D->GetNbinsX();
319 TVectorD vectorH(nbins);
320 for (Int_t ibin=0; ibin<nbins; ibin++) vectorH[ibin]=his1D->GetBinContent(ibin+1);
321 for (Int_t ibin=0; ibin<nbins; ibin++) {
322 Int_t index0=ibin-nmedian;
323 Int_t index1=ibin+nmedian;
324 if (index0<0) {index1+=-index0; index0=0;}
325 if (index1>=nbins) {index0-=index1-nbins+1; index1=nbins-1;}
326 Double_t value= TMath::Median(index1-index0,&(vectorH.GetMatrixArray()[index0]));
327 his1D->SetBinContent(ibin+1, value);
331 Bool_t TStatToolkit::LTMHisto(TH1 *his1D, TVectorD ¶ms , Float_t fraction){
333 // LTM : Trimmed mean on histogram - Modified version for binned data
335 // Robust statistic to estimate properties of the distribution
336 // To handle binning error special treatment
337 // for definition of unbinned data see:
338 // http://en.wikipedia.org/w/index.php?title=Trimmed_estimator&oldid=582847999
340 // Function parameters:
341 // his1D - input histogram
342 // params - vector with parameters
346 // - 3 - error estimate of mean
347 // - 4 - error estimate of RMS
348 // - 5 - first accepted bin position
349 // - 6 - last accepted bin position
351 Int_t nbins = his1D->GetNbinsX();
352 Int_t nentries = (Int_t)his1D->GetEntries();
353 const Double_t kEpsilon=0.0000000001;
355 if (nentries<=0) return 0;
356 if (fraction>1) fraction=0;
357 if (fraction<0) return 0;
358 TVectorD vectorX(nbins);
359 TVectorD vectorMean(nbins);
360 TVectorD vectorRMS(nbins);
362 for (Int_t ibin0=1; ibin0<=nbins; ibin0++) sumCont+=his1D->GetBinContent(ibin0);
364 Double_t minRMS=his1D->GetRMS()*10000;
367 for (Int_t ibin0=1; ibin0<nbins; ibin0++){
368 Double_t sum0=0, sum1=0, sum2=0;
370 for ( ibin1=ibin0; ibin1<nbins; ibin1++){
371 Double_t cont=his1D->GetBinContent(ibin1);
372 Double_t x= his1D->GetBinCenter(ibin1);
376 if ( (ibin0!=ibin1) && sum0>=fraction*sumCont) break;
378 vectorX[ibin0]=his1D->GetBinCenter(ibin0);
379 if (sum0<fraction*sumCont) continue;
381 // substract fractions of bin0 and bin1 to keep sum0=fration*sumCont
383 Double_t diff = sum0-fraction*sumCont;
384 Double_t mean = sum1/sum0;
386 Double_t x0=his1D->GetBinCenter(ibin0);
387 Double_t x1=his1D->GetBinCenter(ibin1);
388 Double_t y0=his1D->GetBinContent(ibin0);
389 Double_t y1=his1D->GetBinContent(ibin1);
391 Double_t d = y0+y1-diff; //enties to keep
393 if (y0<=kEpsilon&&y1>kEpsilon){
396 if (y1<=kEpsilon&&y0>kEpsilon){
399 if (y0>kEpsilon && y1>kEpsilon && x1>x0 ){
400 w0 = (d*(x1-mean))/((x1-x0)*y0);
403 if (w0>1) {w1+=(w0-1)*y0/y1; w0=1;}
404 if (w1>1) {w0+=(w1-1)*y1/y0; w1=1;}
406 if ( (x1>x0) &&TMath::Abs(y0*w0+y1*w1-d)>kEpsilon*sum0){
407 printf(" TStatToolkit::LTMHisto error\n");
415 Double_t xx0=his1D->GetXaxis()->GetBinUpEdge(ibin0)-0.5*w0*his1D->GetBinWidth(ibin0);
416 Double_t xx1=his1D->GetXaxis()->GetBinLowEdge(ibin1)+0.5*w1*his1D->GetBinWidth(ibin1);
424 // choose the bin with smallest rms
427 vectorMean[ibin0]=sum1/sum0;
428 vectorRMS[ibin0]=TMath::Sqrt(TMath::Abs(sum2/sum0-vectorMean[ibin0]*vectorMean[ibin0]));
429 if (vectorRMS[ibin0]<minRMS){
430 minRMS=vectorRMS[ibin0];
432 params[1]=vectorMean[ibin0];
433 params[2]=vectorRMS[ibin0];
434 params[3]=vectorRMS[ibin0]/TMath::Sqrt(sumCont*fraction);
435 params[4]=0; // what is the formula for error of RMS???
438 params[7]=his1D->GetBinCenter(ibin0);
439 params[8]=his1D->GetBinCenter(ibin1);
450 Double_t TStatToolkit::FitGaus(TH1* his, TVectorD *param, TMatrixD */*matrix*/, Float_t xmin, Float_t xmax, Bool_t verbose){
452 // Fit histogram with gaussian function
455 // return value- chi2 - if negative ( not enough points)
456 // his - input histogram
457 // param - vector with parameters
458 // xmin, xmax - range to fit - if xmin=xmax=0 - the full histogram range used
460 // 1. Step - make logarithm
461 // 2. Linear fit (parabola) - more robust - always converge
462 // 3. In case of small statistic bins are averaged
464 static TLinearFitter fitter(3,"pol2");
468 if (his->GetMaximum()<4) return -1;
469 if (his->GetEntries()<12) return -1;
470 if (his->GetRMS()<mat.GetTol()) return -1;
471 Float_t maxEstimate = his->GetEntries()*his->GetBinWidth(1)/TMath::Sqrt((TMath::TwoPi()*his->GetRMS()));
472 Int_t dsmooth = TMath::Nint(6./TMath::Sqrt(maxEstimate));
474 if (maxEstimate<1) return -1;
475 Int_t nbins = his->GetNbinsX();
481 xmin = his->GetXaxis()->GetXmin();
482 xmax = his->GetXaxis()->GetXmax();
484 for (Int_t iter=0; iter<2; iter++){
485 fitter.ClearPoints();
487 for (Int_t ibin=1;ibin<nbins+1; ibin++){
489 Float_t entriesI = his->GetBinContent(ibin);
490 for (Int_t delta = -dsmooth; delta<=dsmooth; delta++){
491 if (ibin+delta>1 &&ibin+delta<nbins-1){
492 entriesI += his->GetBinContent(ibin+delta);
497 Double_t xcenter= his->GetBinCenter(ibin);
498 if (xcenter<xmin || xcenter>xmax) continue;
499 Double_t error=1./TMath::Sqrt(countB);
502 if (par[0]+par[1]*xcenter+par[2]*xcenter*xcenter>20) return 0;
503 cont = TMath::Exp(par[0]+par[1]*xcenter+par[2]*xcenter*xcenter);
504 if (cont>1.) error = 1./TMath::Sqrt(cont*Float_t(countB));
506 if (entriesI>1&&cont>1){
507 fitter.AddPoint(&xcenter,TMath::Log(Float_t(entriesI)),error);
513 fitter.GetParameters(par);
521 fitter.GetParameters(par);
522 fitter.GetCovarianceMatrix(mat);
523 if (TMath::Abs(par[1])<mat.GetTol()) return -1;
524 if (TMath::Abs(par[2])<mat.GetTol()) return -1;
525 Double_t chi2 = fitter.GetChisquare()/Float_t(npoints);
526 //fitter.GetParameters();
527 if (!param) param = new TVectorD(3);
528 // if (!matrix) matrix = new TMatrixD(3,3); // Covariance matrix to be implemented
529 (*param)[1] = par[1]/(-2.*par[2]);
530 (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
531 (*param)[0] = TMath::Exp(par[0]+ par[1]* (*param)[1] + par[2]*(*param)[1]*(*param)[1]);
536 printf("Chi2=%f\n",chi2);
537 TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",his->GetXaxis()->GetXmin(),his->GetXaxis()->GetXmax());
538 f1->SetParameter(0, (*param)[0]);
539 f1->SetParameter(1, (*param)[1]);
540 f1->SetParameter(2, (*param)[2]);
546 Double_t TStatToolkit::FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param, TMatrixD */*matrix*/, Bool_t verbose){
548 // Fit histogram with gaussian function
551 // nbins: size of the array and number of histogram bins
552 // xMin, xMax: histogram range
553 // param: paramters of the fit (0-Constant, 1-Mean, 2-Sigma)
554 // matrix: covariance matrix -- not implemented yet, pass dummy matrix!!!
557 // >0: the chi2 returned by TLinearFitter
558 // -3: only three points have been used for the calculation - no fitter was used
559 // -2: only two points have been used for the calculation - center of gravity was uesed for calculation
560 // -1: only one point has been used for the calculation - center of gravity was uesed for calculation
561 // -4: invalid result!!
564 // 1. Step - make logarithm
565 // 2. Linear fit (parabola) - more robust - always converge
567 static TLinearFitter fitter(3,"pol2");
568 static TMatrixD mat(3,3);
569 static Double_t kTol = mat.GetTol();
570 fitter.StoreData(kFALSE);
571 fitter.ClearPoints();
576 Float_t rms = TMath::RMS(nBins,arr);
577 Float_t max = TMath::MaxElement(nBins,arr);
578 Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
587 for (Int_t i=0; i<nBins; i++){
589 if (arr[i]>0) nfilled++;
592 if (max<4) return -4;
593 if (entries<12) return -4;
594 if (rms<kTol) return -4;
600 for (Int_t ibin=0;ibin<nBins; ibin++){
601 Float_t entriesI = arr[ibin];
603 Double_t xcenter = xMin+(ibin+0.5)*binWidth;
605 Float_t error = 1./TMath::Sqrt(entriesI);
606 Float_t val = TMath::Log(Float_t(entriesI));
607 fitter.AddPoint(&xcenter,val,error);
610 matA(npoints,1)=xcenter;
611 matA(npoints,2)=xcenter*xcenter;
613 meanCOG+=xcenter*entriesI;
614 rms2COG +=xcenter*entriesI*xcenter;
625 //analytic calculation of the parameters for three points
634 // use fitter for more than three points
636 fitter.GetParameters(par);
637 fitter.GetCovarianceMatrix(mat);
638 chi2 = fitter.GetChisquare()/Float_t(npoints);
640 if (TMath::Abs(par[1])<kTol) return -4;
641 if (TMath::Abs(par[2])<kTol) return -4;
643 if (!param) param = new TVectorD(3);
644 //if (!matrix) matrix = new TMatrixD(3,3); // !!!!might be a memory leek. use dummy matrix pointer to call this function! // Covariance matrix to be implemented
646 (*param)[1] = par[1]/(-2.*par[2]);
647 (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
648 Double_t lnparam0 = par[0]+ par[1]* (*param)[1] + par[2]*(*param)[1]*(*param)[1];
649 if ( lnparam0>307 ) return -4;
650 (*param)[0] = TMath::Exp(lnparam0);
655 printf("Chi2=%f\n",chi2);
656 TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",xMin,xMax);
657 f1->SetParameter(0, (*param)[0]);
658 f1->SetParameter(1, (*param)[1]);
659 f1->SetParameter(2, (*param)[2]);
666 //use center of gravity for 2 points
670 (*param)[1] = meanCOG;
671 (*param)[2] = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
677 (*param)[1] = meanCOG;
678 (*param)[2] = binWidth/TMath::Sqrt(12);
686 Float_t TStatToolkit::GetCOG(const Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms, Float_t *sum)
689 // calculate center of gravity rms and sum for array 'arr' with nBins an a x range xMin to xMax
690 // return COG; in case of failure return xMin
697 Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
699 for (Int_t ibin=0; ibin<nBins; ibin++){
700 Float_t entriesI = (Float_t)arr[ibin];
701 Double_t xcenter = xMin+(ibin+0.5)*binWidth;
703 meanCOG += xcenter*entriesI;
704 rms2COG += xcenter*entriesI*xcenter;
709 if ( sumCOG == 0 ) return xMin;
714 (*rms) = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
715 if ( npoints == 1 ) (*rms) = binWidth/TMath::Sqrt(12);
726 ///////////////////////////////////////////////////////////////
727 ////////////// TEST functions /////////////////////////
728 ///////////////////////////////////////////////////////////////
734 void TStatToolkit::TestGausFit(Int_t nhistos){
736 // Test performance of the parabolic - gaussian fit - compare it with
738 // nhistos - number of histograms to be used for test
740 TTreeSRedirector *pcstream = new TTreeSRedirector("fitdebug.root","recreate");
742 Float_t *xTrue = new Float_t[nhistos];
743 Float_t *sTrue = new Float_t[nhistos];
744 TVectorD **par1 = new TVectorD*[nhistos];
745 TVectorD **par2 = new TVectorD*[nhistos];
749 TH1F **h1f = new TH1F*[nhistos];
750 TF1 *myg = new TF1("myg","gaus");
751 TF1 *fit = new TF1("fit","gaus");
755 for (Int_t i=0;i<nhistos; i++){
756 par1[i] = new TVectorD(3);
757 par2[i] = new TVectorD(3);
758 h1f[i] = new TH1F(Form("h1f%d",i),Form("h1f%d",i),20,-10,10);
759 xTrue[i]= gRandom->Rndm();
761 sTrue[i]= .75+gRandom->Rndm()*.5;
762 myg->SetParameters(1,xTrue[i],sTrue[i]);
763 h1f[i]->FillRandom("myg");
769 for (Int_t i=0; i<nhistos; i++){
770 h1f[i]->Fit(fit,"0q");
771 (*par1[i])(0) = fit->GetParameter(0);
772 (*par1[i])(1) = fit->GetParameter(1);
773 (*par1[i])(2) = fit->GetParameter(2);
776 printf("Gaussian fit\t");
780 //TStatToolkit gaus fit
781 for (Int_t i=0; i<nhistos; i++){
782 TStatToolkit::FitGaus(h1f[i]->GetArray()+1,h1f[i]->GetNbinsX(),h1f[i]->GetXaxis()->GetXmin(),h1f[i]->GetXaxis()->GetXmax(),par2[i],&dummy);
786 printf("Parabolic fit\t");
789 for (Int_t i=0;i<nhistos; i++){
790 Float_t xt = xTrue[i];
791 Float_t st = sTrue[i];
800 for (Int_t i=0;i<nhistos; i++){
817 TGraph2D * TStatToolkit::MakeStat2D(TH3 * his, Int_t delta0, Int_t delta1, Int_t type){
821 // delta - number of bins to integrate
822 // type - 0 - mean value
824 TAxis * xaxis = his->GetXaxis();
825 TAxis * yaxis = his->GetYaxis();
826 // TAxis * zaxis = his->GetZaxis();
827 Int_t nbinx = xaxis->GetNbins();
828 Int_t nbiny = yaxis->GetNbins();
831 TGraph2D *graph = new TGraph2D(nbinx*nbiny);
833 for (Int_t ix=0; ix<nbinx;ix++)
834 for (Int_t iy=0; iy<nbiny;iy++){
835 Float_t xcenter = xaxis->GetBinCenter(ix);
836 Float_t ycenter = yaxis->GetBinCenter(iy);
837 snprintf(name,1000,"%s_%d_%d",his->GetName(), ix,iy);
838 TH1 *projection = his->ProjectionZ(name,ix-delta0,ix+delta0,iy-delta1,iy+delta1);
840 if (type==0) stat = projection->GetMean();
841 if (type==1) stat = projection->GetRMS();
842 if (type==2 || type==3){
844 TStatToolkit::LTM((TH1F*)projection,&vec,0.7);
845 if (type==2) stat= vec[1];
846 if (type==3) stat= vec[0];
848 if (type==4|| type==5){
849 projection->Fit(&f1);
850 if (type==4) stat= f1.GetParameter(1);
851 if (type==5) stat= f1.GetParameter(2);
853 //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
854 graph->SetPoint(icount,xcenter, ycenter, stat);
860 TGraphErrors * TStatToolkit::MakeStat1D(TH2 * his, Int_t deltaBin, Double_t fraction, Int_t returnType, Int_t markerStyle, Int_t markerColor){
862 // function to retrieve the "mean and RMS estimate" of 2D histograms
864 // Robust statistic to estimate properties of the distribution
865 // See http://en.wikipedia.org/wiki/Trimmed_estimator
867 // deltaBin - number of bins to integrate (bin+-deltaBin)
868 // fraction - fraction of values for the LTM and for the gauss fit
874 // 4 - Gaus fit mean - on LTM range
875 // 5 - Gaus fit sigma - on LTM range
877 TAxis * xaxis = his->GetXaxis();
878 Int_t nbinx = xaxis->GetNbins();
882 TVectorD vecX(nbinx);
883 TVectorD vecXErr(nbinx);
884 TVectorD vecY(nbinx);
885 TVectorD vecYErr(nbinx);
890 for (Int_t jx=1; jx<=nbinx;jx++){
892 Float_t xcenter = xaxis->GetBinCenter(jx);
893 snprintf(name,1000,"%s_%d",his->GetName(), ix);
894 TH1 *projection = his->ProjectionY(name,TMath::Max(jx-deltaBin,1),TMath::Min(jx+deltaBin,nbinx));
897 TStatToolkit::LTMHisto((TH1F*)projection,vecLTM,fraction);
900 stat = projection->GetMean();
901 err = projection->GetMeanError();
904 stat = projection->GetRMS();
905 err = projection->GetRMSError();
907 if (returnType==2 || returnType==3){
908 if (returnType==2) {stat= vecLTM[1]; err =projection->GetRMSError();}
909 if (returnType==3) {stat= vecLTM[2]; err =projection->GetRMSError();}
911 if (returnType==4|| returnType==5){
912 projection->Fit(&f1,"QN","QN", vecLTM[7], vecLTM[8]);
914 stat= f1.GetParameter(1);
915 err=f1.GetParError(1);
918 stat= f1.GetParameter(2);
919 err=f1.GetParError(2);
922 vecX[icount]=xcenter;
928 TGraphErrors *graph = new TGraphErrors(icount,vecX.GetMatrixArray(), vecY.GetMatrixArray(),0, vecYErr.GetMatrixArray());
929 graph->SetMarkerStyle(markerStyle);
930 graph->SetMarkerColor(markerColor);
938 TString* TStatToolkit::FitPlane(TTree *tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints, TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac, Int_t start, Int_t stop,Bool_t fix0){
940 // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
941 // returns chi2, fitParam and covMatrix
942 // returns TString with fitted formula
945 TString formulaStr(formula);
946 TString drawStr(drawCommand);
947 TString cutStr(cuts);
950 TString strVal(drawCommand);
951 if (strVal.Contains(":")){
952 TObjArray* valTokens = strVal.Tokenize(":");
953 drawStr = valTokens->At(0)->GetName();
954 ferr = valTokens->At(1)->GetName();
959 formulaStr.ReplaceAll("++", "~");
960 TObjArray* formulaTokens = formulaStr.Tokenize("~");
961 Int_t dim = formulaTokens->GetEntriesFast();
963 fitParam.ResizeTo(dim);
964 covMatrix.ResizeTo(dim,dim);
966 TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
967 fitter->StoreData(kTRUE);
968 fitter->ClearPoints();
970 Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start, start);
972 delete formulaTokens;
973 return new TString("An ERROR has occured during fitting!");
975 Double_t **values = new Double_t*[dim+1] ;
976 for (Int_t i=0; i<dim+1; i++) values[i]=NULL;
978 entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff", stop-start, start);
980 delete formulaTokens;
982 return new TString("An ERROR has occured during fitting!");
984 Double_t *errors = new Double_t[entries];
985 memcpy(errors, tree->GetV1(), entries*sizeof(Double_t));
987 for (Int_t i = 0; i < dim + 1; i++){
989 if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start);
990 else centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start);
992 if (entries != centries) {
995 return new TString("An ERROR has occured during fitting!");
997 values[i] = new Double_t[entries];
998 memcpy(values[i], tree->GetV1(), entries*sizeof(Double_t));
1001 // add points to the fitter
1002 for (Int_t i = 0; i < entries; i++){
1004 for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
1005 fitter->AddPoint(x, values[dim][i], errors[i]);
1009 if (frac>0.5 && frac<1){
1010 fitter->EvalRobust(frac);
1013 fitter->FixParameter(0,0);
1017 fitter->GetParameters(fitParam);
1018 fitter->GetCovarianceMatrix(covMatrix);
1019 chi2 = fitter->GetChisquare();
1021 TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula;
1023 for (Int_t iparam = 0; iparam < dim; iparam++) {
1024 returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
1025 if (iparam < dim-1) returnFormula.Append("+");
1027 returnFormula.Append(" )");
1030 for (Int_t j=0; j<dim+1;j++) delete [] values[j];
1033 delete formulaTokens;
1037 return preturnFormula;
1040 TString* TStatToolkit::FitPlaneConstrain(TTree *tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints, TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac, Int_t start, Int_t stop,Double_t constrain){
1042 // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
1043 // returns chi2, fitParam and covMatrix
1044 // returns TString with fitted formula
1047 TString formulaStr(formula);
1048 TString drawStr(drawCommand);
1049 TString cutStr(cuts);
1052 TString strVal(drawCommand);
1053 if (strVal.Contains(":")){
1054 TObjArray* valTokens = strVal.Tokenize(":");
1055 drawStr = valTokens->At(0)->GetName();
1056 ferr = valTokens->At(1)->GetName();
1061 formulaStr.ReplaceAll("++", "~");
1062 TObjArray* formulaTokens = formulaStr.Tokenize("~");
1063 Int_t dim = formulaTokens->GetEntriesFast();
1065 fitParam.ResizeTo(dim);
1066 covMatrix.ResizeTo(dim,dim);
1068 TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
1069 fitter->StoreData(kTRUE);
1070 fitter->ClearPoints();
1072 Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start, start);
1073 if (entries == -1) {
1074 delete formulaTokens;
1075 return new TString("An ERROR has occured during fitting!");
1077 Double_t **values = new Double_t*[dim+1] ;
1078 for (Int_t i=0; i<dim+1; i++) values[i]=NULL;
1080 entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff", stop-start, start);
1081 if (entries == -1) {
1082 delete formulaTokens;
1084 return new TString("An ERROR has occured during fitting!");
1086 Double_t *errors = new Double_t[entries];
1087 memcpy(errors, tree->GetV1(), entries*sizeof(Double_t));
1089 for (Int_t i = 0; i < dim + 1; i++){
1091 if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start);
1092 else centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start);
1094 if (entries != centries) {
1097 delete formulaTokens;
1098 return new TString("An ERROR has occured during fitting!");
1100 values[i] = new Double_t[entries];
1101 memcpy(values[i], tree->GetV1(), entries*sizeof(Double_t));
1104 // add points to the fitter
1105 for (Int_t i = 0; i < entries; i++){
1107 for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
1108 fitter->AddPoint(x, values[dim][i], errors[i]);
1111 for (Int_t i = 0; i < dim; i++){
1113 for (Int_t j=0; j<dim;j++) if (i!=j) x[j]=0;
1115 fitter->AddPoint(x, 0, constrain);
1121 if (frac>0.5 && frac<1){
1122 fitter->EvalRobust(frac);
1124 fitter->GetParameters(fitParam);
1125 fitter->GetCovarianceMatrix(covMatrix);
1126 chi2 = fitter->GetChisquare();
1129 TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula;
1131 for (Int_t iparam = 0; iparam < dim; iparam++) {
1132 returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
1133 if (iparam < dim-1) returnFormula.Append("+");
1135 returnFormula.Append(" )");
1137 for (Int_t j=0; j<dim+1;j++) delete [] values[j];
1141 delete formulaTokens;
1145 return preturnFormula;
1150 TString* TStatToolkit::FitPlaneFixed(TTree *tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints, TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac, Int_t start, Int_t stop){
1152 // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
1153 // returns chi2, fitParam and covMatrix
1154 // returns TString with fitted formula
1157 TString formulaStr(formula);
1158 TString drawStr(drawCommand);
1159 TString cutStr(cuts);
1162 TString strVal(drawCommand);
1163 if (strVal.Contains(":")){
1164 TObjArray* valTokens = strVal.Tokenize(":");
1165 drawStr = valTokens->At(0)->GetName();
1166 ferr = valTokens->At(1)->GetName();
1171 formulaStr.ReplaceAll("++", "~");
1172 TObjArray* formulaTokens = formulaStr.Tokenize("~");
1173 Int_t dim = formulaTokens->GetEntriesFast();
1175 fitParam.ResizeTo(dim);
1176 covMatrix.ResizeTo(dim,dim);
1177 TString fitString="x0";
1178 for (Int_t i=1; i<dim; i++) fitString+=Form("++x%d",i);
1179 TLinearFitter* fitter = new TLinearFitter(dim, fitString.Data());
1180 fitter->StoreData(kTRUE);
1181 fitter->ClearPoints();
1183 Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start, start);
1184 if (entries == -1) {
1185 delete formulaTokens;
1186 return new TString("An ERROR has occured during fitting!");
1188 Double_t **values = new Double_t*[dim+1] ;
1189 for (Int_t i=0; i<dim+1; i++) values[i]=NULL;
1191 entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff", stop-start, start);
1192 if (entries == -1) {
1194 delete formulaTokens;
1195 return new TString("An ERROR has occured during fitting!");
1197 Double_t *errors = new Double_t[entries];
1198 memcpy(errors, tree->GetV1(), entries*sizeof(Double_t));
1200 for (Int_t i = 0; i < dim + 1; i++){
1202 if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start);
1203 else centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start);
1205 if (entries != centries) {
1208 delete formulaTokens;
1209 return new TString("An ERROR has occured during fitting!");
1211 values[i] = new Double_t[entries];
1212 memcpy(values[i], tree->GetV1(), entries*sizeof(Double_t));
1215 // add points to the fitter
1216 for (Int_t i = 0; i < entries; i++){
1218 for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
1219 fitter->AddPoint(x, values[dim][i], errors[i]);
1223 if (frac>0.5 && frac<1){
1224 fitter->EvalRobust(frac);
1226 fitter->GetParameters(fitParam);
1227 fitter->GetCovarianceMatrix(covMatrix);
1228 chi2 = fitter->GetChisquare();
1231 TString *preturnFormula = new TString("("), &returnFormula = *preturnFormula;
1233 for (Int_t iparam = 0; iparam < dim; iparam++) {
1234 returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam]));
1235 if (iparam < dim-1) returnFormula.Append("+");
1237 returnFormula.Append(" )");
1240 for (Int_t j=0; j<dim+1;j++) delete [] values[j];
1242 delete formulaTokens;
1246 return preturnFormula;
1253 Int_t TStatToolkit::GetFitIndex(const TString fString, const TString subString){
1255 // fitString - ++ separated list of fits
1256 // substring - ++ separated list of the requiered substrings
1258 // return the last occurance of substring in fit string
1260 TObjArray *arrFit = fString.Tokenize("++");
1261 TObjArray *arrSub = subString.Tokenize("++");
1263 for (Int_t i=0; i<arrFit->GetEntries(); i++){
1265 TString str =arrFit->At(i)->GetName();
1266 for (Int_t isub=0; isub<arrSub->GetEntries(); isub++){
1267 if (str.Contains(arrSub->At(isub)->GetName())==0) isOK=kFALSE;
1277 TString TStatToolkit::FilterFit(const TString &input, const TString filter, TVectorD ¶m, TMatrixD & covar){
1279 // Filter fit expression make sub-fit
1281 TObjArray *array0= input.Tokenize("++");
1282 TObjArray *array1= filter.Tokenize("++");
1283 //TString *presult=new TString("(0");
1284 TString result="(0.0";
1285 for (Int_t i=0; i<array0->GetEntries(); i++){
1287 TString str(array0->At(i)->GetName());
1288 for (Int_t j=0; j<array1->GetEntries(); j++){
1289 if (str.Contains(array1->At(j)->GetName())==0) isOK=kFALSE;
1293 result+=Form("*(%f)",param[i+1]);
1294 printf("%f\t%f\t%s\n",param[i+1], TMath::Sqrt(covar(i+1,i+1)),str.Data());
1303 void TStatToolkit::Update1D(Double_t delta, Double_t sigma, Int_t s1, TMatrixD &vecXk, TMatrixD &covXk){
1305 // Update parameters and covariance - with one measurement
1307 // vecXk - input vector - Updated in function
1308 // covXk - covariance matrix - Updated in function
1309 // delta, sigma, s1 - new measurement, rms of new measurement and the index of measurement
1310 const Int_t knMeas=1;
1311 Int_t knElem=vecXk.GetNrows();
1313 TMatrixD mat1(knElem,knElem); // update covariance matrix
1314 TMatrixD matHk(1,knElem); // vector to mesurement
1315 TMatrixD vecYk(knMeas,1); // Innovation or measurement residual
1316 TMatrixD matHkT(knElem,knMeas); // helper matrix Hk transpose
1317 TMatrixD matSk(knMeas,knMeas); // Innovation (or residual) covariance
1318 TMatrixD matKk(knElem,knMeas); // Optimal Kalman gain
1319 TMatrixD covXk2(knElem,knElem); // helper matrix
1320 TMatrixD covXk3(knElem,knElem); // helper matrix
1321 TMatrixD vecZk(1,1);
1322 TMatrixD measR(1,1);
1324 measR(0,0)=sigma*sigma;
1327 for (Int_t iel=0;iel<knElem;iel++)
1328 for (Int_t ip=0;ip<knMeas;ip++) matHk(ip,iel)=0;
1330 for (Int_t iel=0;iel<knElem;iel++) {
1331 for (Int_t jel=0;jel<knElem;jel++) mat1(iel,jel)=0;
1336 vecYk = vecZk-matHk*vecXk; // Innovation or measurement residual
1337 matHkT=matHk.T(); matHk.T();
1338 matSk = (matHk*(covXk*matHkT))+measR; // Innovation (or residual) covariance
1340 matKk = (covXk*matHkT)*matSk; // Optimal Kalman gain
1341 vecXk += matKk*vecYk; // updated vector
1342 covXk2= (mat1-(matKk*matHk));
1343 covXk3 = covXk2*covXk;
1345 Int_t nrows=covXk3.GetNrows();
1347 for (Int_t irow=0; irow<nrows; irow++)
1348 for (Int_t icol=0; icol<nrows; icol++){
1349 // rounding problems - make matrix again symteric
1350 covXk(irow,icol)=(covXk3(irow,icol)+covXk3(icol,irow))*0.5;
1356 void TStatToolkit::Constrain1D(const TString &input, const TString filter, TVectorD ¶m, TMatrixD & covar, Double_t mean, Double_t sigma){
1358 // constrain linear fit
1359 // input - string description of fit function
1360 // filter - string filter to select sub fits
1361 // param,covar - parameters and covariance matrix of the fit
1362 // mean,sigma - new measurement uning which the fit is updated
1365 TObjArray *array0= input.Tokenize("++");
1366 TObjArray *array1= filter.Tokenize("++");
1367 TMatrixD paramM(param.GetNrows(),1);
1368 for (Int_t i=0; i<=array0->GetEntries(); i++){paramM(i,0)=param(i);}
1370 if (filter.Length()==0){
1371 TStatToolkit::Update1D(mean, sigma, 0, paramM, covar);//
1373 for (Int_t i=0; i<array0->GetEntries(); i++){
1375 TString str(array0->At(i)->GetName());
1376 for (Int_t j=0; j<array1->GetEntries(); j++){
1377 if (str.Contains(array1->At(j)->GetName())==0) isOK=kFALSE;
1380 TStatToolkit::Update1D(mean, sigma, i+1, paramM, covar);//
1384 for (Int_t i=0; i<=array0->GetEntries(); i++){
1385 param(i)=paramM(i,0);
1391 TString TStatToolkit::MakeFitString(const TString &input, const TVectorD ¶m, const TMatrixD & covar, Bool_t verbose){
1395 TObjArray *array0= input.Tokenize("++");
1396 TString result=Form("(%f",param[0]);
1397 printf("%f\t%f\t\n", param[0], TMath::Sqrt(covar(0,0)));
1398 for (Int_t i=0; i<array0->GetEntries(); i++){
1399 TString str(array0->At(i)->GetName());
1401 result+=Form("*(%f)",param[i+1]);
1402 if (verbose) printf("%f\t%f\t%s\n", param[i+1], TMath::Sqrt(covar(i+1,i+1)),str.Data());
1409 TGraphErrors * TStatToolkit::MakeGraphErrors(TTree * tree, const char * expr, const char * cut, Int_t mstyle, Int_t mcolor, Float_t msize, Float_t offset){
1411 // Query a graph errors
1412 // return TGraphErrors specified by expr and cut
1413 // Example usage TStatToolkit::MakeGraphError(tree,"Y:X:ErrY","X>0", 25,2,0.4)
1414 // tree - tree with variable
1416 const Int_t entries = tree->Draw(expr,cut,"goff");
1419 t.Error("TStatToolkit::MakeGraphError",Form("Empty or Not valid expression (%s) or cut *%s)", expr,cut));
1422 if ( tree->GetV2()==0){
1424 t.Error("TStatToolkit::MakeGraphError",Form("Not valid expression (%s) ", expr));
1427 TGraphErrors * graph=0;
1428 if ( tree->GetV3()!=0){
1429 graph = new TGraphErrors (entries, tree->GetV2(),tree->GetV1(),0,tree->GetV3());
1431 graph = new TGraphErrors (entries, tree->GetV2(),tree->GetV1(),0,0);
1433 graph->SetMarkerStyle(mstyle);
1434 graph->SetMarkerColor(mcolor);
1435 graph->SetLineColor(mcolor);
1436 graph->SetTitle(expr);
1437 TString chstring(expr);
1438 TObjArray *charray = chstring.Tokenize(":");
1439 graph->GetXaxis()->SetTitle(charray->At(1)->GetName());
1440 graph->GetYaxis()->SetTitle(charray->At(0)->GetName());
1442 if (msize>0) graph->SetMarkerSize(msize);
1443 for(Int_t i=0;i<graph->GetN();i++) graph->GetX()[i]+=offset;
1449 TGraph * TStatToolkit::MakeGraphSparse(TTree * tree, const char * expr, const char * cut, Int_t mstyle, Int_t mcolor, Float_t msize, Float_t offset){
1451 // Make a sparse draw of the variables
1452 // Format of expr : Var:Run or Var:Run:ErrorY or Var:Run:ErrorY:ErrorX
1453 // offset : points can slightly be shifted in x for better visibility with more graphs
1455 // Written by Weilin.Yu
1456 // updated & merged with QA-code by Patrick Reichelt
1458 const Int_t entries = tree->Draw(expr,cut,"goff");
1461 t.Error("TStatToolkit::MakeGraphSparse",Form("Empty or Not valid expression (%s) or cut (%s)", expr, cut));
1464 // TGraph * graph = (TGraph*)gPad->GetPrimitive("Graph"); // 2D
1466 Double_t *graphY, *graphX;
1467 graphY = tree->GetV1();
1468 graphX = tree->GetV2();
1470 // sort according to run number
1471 Int_t *index = new Int_t[entries*4];
1472 TMath::Sort(entries,graphX,index,kFALSE);
1474 // define arrays for the new graph
1475 Double_t *unsortedX = new Double_t[entries];
1476 Int_t *runNumber = new Int_t[entries];
1477 Double_t count = 0.5;
1479 // evaluate arrays for the new graph according to the run-number
1482 unsortedX[index[0]] = count;
1483 runNumber[0] = graphX[index[0]];
1484 // loop the rest of entries
1485 for(Int_t i=1;i<entries;i++)
1487 if(graphX[index[i]]==graphX[index[i-1]])
1488 unsortedX[index[i]] = count;
1489 else if(graphX[index[i]]!=graphX[index[i-1]]){
1492 unsortedX[index[i]] = count;
1493 runNumber[icount]=graphX[index[i]];
1497 // count the number of xbins (run-wise) for the new graph
1498 const Int_t newNbins = int(count+0.5);
1499 Double_t *newBins = new Double_t[newNbins+1];
1500 for(Int_t i=0; i<=count+1;i++){
1504 // define and fill the new graph
1505 TGraph *graphNew = 0;
1506 if (tree->GetV3()) {
1507 if (tree->GetV4()) {
1508 graphNew = new TGraphErrors(entries,unsortedX,graphY,tree->GetV4(),tree->GetV3());
1510 else { graphNew = new TGraphErrors(entries,unsortedX,graphY,0,tree->GetV3()); }
1512 else { graphNew = new TGraphErrors(entries,unsortedX,graphY,0,0); }
1513 // with "Set(...)", the x-axis is being sorted
1514 graphNew->GetXaxis()->Set(newNbins,newBins);
1516 // set the bins for the x-axis, apply shifting of points
1518 for(Int_t i=0;i<count;i++){
1519 snprintf(xName,50,"%d",runNumber[i]);
1520 graphNew->GetXaxis()->SetBinLabel(i+1,xName);
1521 graphNew->GetX()[i]+=offset;
1524 graphNew->GetHistogram()->SetTitle("");
1525 graphNew->SetMarkerStyle(mstyle);
1526 graphNew->SetMarkerColor(mcolor);
1527 if (msize>0) graphNew->SetMarkerSize(msize);
1528 delete [] unsortedX;
1529 delete [] runNumber;
1533 graphNew->SetTitle(expr);
1534 TString chstring(expr);
1535 TObjArray *charray = chstring.Tokenize(":");
1536 graphNew->GetXaxis()->SetTitle(charray->At(1)->GetName());
1537 graphNew->GetYaxis()->SetTitle(charray->At(0)->GetName());
1545 // functions used for the trending
1548 Int_t TStatToolkit::MakeStatAlias(TTree * tree, const char * expr, const char * cut, const char * alias)
1551 // Add alias using statistical values of a given variable.
1552 // (by MI, Patrick Reichelt)
1554 // tree - input tree
1555 // expr - variable expression
1556 // cut - selection criteria
1557 // Output - return number of entries used to define variable
1558 // In addition mean, rms, median, and robust mean and rms (choosing fraction of data with smallest RMS)
1561 1.) create the robust estimators for variable expr="QA.TPC.CPass1.meanTPCncl" and create a corresponding
1562 aliases with the prefix alias[0]="ncl", calculated using fraction alias[1]="0.90"
1564 TStatToolkit::MakeStatAlias(tree,"QA.TPC.CPass1.meanTPCncl","QA.TPC.CPass1.status>0","ncl:0.9");
1565 root [4] tree->GetListOfAliases().Print()
1566 OBJ: TNamed ncl_Median (130.964333+0)
1567 OBJ: TNamed ncl_Mean (122.120387+0)
1568 OBJ: TNamed ncl_RMS (33.509623+0)
1569 OBJ: TNamed ncl_Mean90 (131.503862+0)
1570 OBJ: TNamed ncl_RMS90 (3.738260+0)
1573 Int_t entries = tree->Draw(expr,cut,"goff");
1575 printf("Expression or cut not valid:\t%s\t%s\n", expr, cut);
1579 TObjArray* oaAlias = TString(alias).Tokenize(":");
1580 if (oaAlias->GetEntries()<2) return 0;
1581 Float_t entryFraction = atof( oaAlias->At(1)->GetName() );
1583 Double_t median = TMath::Median(entries,tree->GetV1());
1584 Double_t mean = TMath::Mean(entries,tree->GetV1());
1585 Double_t rms = TMath::RMS(entries,tree->GetV1());
1586 Double_t meanEF=0, rmsEF=0;
1587 TStatToolkit::EvaluateUni(entries, tree->GetV1(), meanEF, rmsEF, entries*entryFraction);
1589 tree->SetAlias(Form("%s_Median",oaAlias->At(0)->GetName()), Form("(%f+0)",median));
1590 tree->SetAlias(Form("%s_Mean",oaAlias->At(0)->GetName()), Form("(%f+0)",mean));
1591 tree->SetAlias(Form("%s_RMS",oaAlias->At(0)->GetName()), Form("(%f+0)",rms));
1592 tree->SetAlias(Form("%s_Mean%d",oaAlias->At(0)->GetName(),Int_t(entryFraction*100)), Form("(%f+0)",meanEF));
1593 tree->SetAlias(Form("%s_RMS%d",oaAlias->At(0)->GetName(),Int_t(entryFraction*100)), Form("(%f+0)",rmsEF));
1598 Int_t TStatToolkit::SetStatusAlias(TTree * tree, const char * expr, const char * cut, const char * alias)
1601 // Add alias to trending tree using statistical values of a given variable.
1602 // (by MI, Patrick Reichelt)
1604 // format of expr : varname (e.g. meanTPCncl)
1605 // format of cut : char like in TCut
1606 // format of alias: alias:query:entryFraction(EF) (fraction of entries used for uniformity evaluation)
1607 // e.g.: varname_Out:(abs(varname-meanEF)>6.*rmsEF):0.8
1608 // available internal variables are: 'varname, Median, Mean, MeanEF, RMS, RMSEF'
1609 // in the alias, 'varname' will be replaced by its content, and 'EF' by the percentage (e.g. MeanEF -> Mean80)
1612 1.) Define robust mean (possible, but easier done with TStatToolkit::MakeStatAlias(...))
1613 TStatToolkit::SetStatusAlias(tree, "meanTPCnclF", "meanTPCnclF>0", "meanTPCnclF_MeanEF:MeanEF:0.80") ;
1614 root [10] tree->GetListOfAliases()->Print()
1615 Collection name='TList', class='TList', size=1
1616 OBJ: TNamed meanTPCnclF_Mean80 0.899308
1617 2.) create alias outlyers - 6 sigma cut
1618 TStatToolkit::SetStatusAlias(tree, "meanTPCnclF", "meanTPCnclF>0", "meanTPCnclF_Out:(abs(meanTPCnclF-MeanEF)>6.*RMSEF):0.8")
1619 meanTPCnclF_Out ==> (abs(meanTPCnclF-0.899308)>6.*0.016590)
1620 3.) the same functionality as in 2.)
1621 TStatToolkit::SetStatusAlias(tree, "meanTPCnclF", "meanTPCnclF>0", "varname_Out2:(abs(varname-MeanEF)>6.*RMSEF):0.8")
1622 meanTPCnclF_Out2 ==> (abs(meanTPCnclF-0.899308)>6.*0.016590)
1625 Int_t entries = tree->Draw(expr,cut,"goff");
1627 printf("Expression or cut not valid:\t%s\t%s\n", expr, cut);
1631 TObjArray* oaVar = TString(expr).Tokenize(":");
1633 snprintf(varname,50,"%s", oaVar->At(0)->GetName());
1635 TObjArray* oaAlias = TString(alias).Tokenize(":");
1636 if (oaAlias->GetEntries()<3) return 0;
1637 Float_t entryFraction = atof( oaAlias->At(2)->GetName() );
1639 Double_t median = TMath::Median(entries,tree->GetV1());
1640 Double_t mean = TMath::Mean(entries,tree->GetV1());
1641 Double_t rms = TMath::RMS(entries,tree->GetV1());
1642 Double_t meanEF=0, rmsEF=0;
1643 TStatToolkit::EvaluateUni(entries, tree->GetV1(), meanEF, rmsEF, entries*entryFraction);
1645 TString sAlias( oaAlias->At(0)->GetName() );
1646 sAlias.ReplaceAll("varname",varname);
1647 sAlias.ReplaceAll("MeanEF", Form("Mean%1.0f",entryFraction*100) );
1648 sAlias.ReplaceAll("RMSEF", Form("RMS%1.0f",entryFraction*100) );
1649 TString sQuery( oaAlias->At(1)->GetName() );
1650 sQuery.ReplaceAll("varname",varname);
1651 sQuery.ReplaceAll("MeanEF", Form("%f",meanEF) );
1652 sQuery.ReplaceAll("RMSEF", Form("%f",rmsEF) ); //make sure to replace 'RMSEF' before 'RMS'...
1653 sQuery.ReplaceAll("Median", Form("%f",median) );
1654 sQuery.ReplaceAll("Mean", Form("%f",mean) );
1655 sQuery.ReplaceAll("RMS", Form("%f",rms) );
1656 printf("define alias:\t%s = %s\n", sAlias.Data(), sQuery.Data());
1660 snprintf(query,200,"%s", sQuery.Data());
1661 snprintf(aname,200,"%s", sAlias.Data());
1662 tree->SetAlias(aname, query);
1668 TMultiGraph* TStatToolkit::MakeStatusMultGr(TTree * tree, const char * expr, const char * cut, const char * alias, Int_t igr)
1671 // Compute a trending multigraph that shows for which runs a variable has outliers.
1672 // (by MI, Patrick Reichelt)
1674 // format of expr : varname:xaxis (e.g. meanTPCncl:run)
1675 // format of cut : char like in TCut
1676 // format of alias: (1):(varname_Out==0):(varname_Out)[:(varname_Warning):...]
1677 // in the alias, 'varname' will be replaced by its content (e.g. varname_Out -> meanTPCncl_Out)
1678 // note: the aliases 'varname_Out' etc have to be defined by function TStatToolkit::SetStatusAlias(...)
1679 // counter igr is used to shift the multigraph in y when filling a TObjArray.
1681 TObjArray* oaVar = TString(expr).Tokenize(":");
1682 if (oaVar->GetEntries()<2) return 0;
1685 snprintf(varname,50,"%s", oaVar->At(0)->GetName());
1686 snprintf(var_x ,50,"%s", oaVar->At(1)->GetName());
1688 TString sAlias(alias);
1689 sAlias.ReplaceAll("varname",varname);
1690 TObjArray* oaAlias = TString(sAlias.Data()).Tokenize(":");
1691 if (oaAlias->GetEntries()<3) return 0;
1694 TMultiGraph* multGr = new TMultiGraph();
1695 Int_t marArr[6] = {24+igr%2, 20+igr%2, 20+igr%2, 20+igr%2, 22, 23};
1696 Int_t colArr[6] = {kBlack, kBlack, kRed, kOrange, kMagenta, kViolet};
1697 Double_t sizArr[6] = {1.2, 1.1, 1.0, 1.0, 1, 1};
1698 const Int_t ngr = oaAlias->GetEntriesFast();
1699 for (Int_t i=0; i<ngr; i++){
1700 if (i==2) continue; // the Fatal(Out) graph will be added in the end to be plotted on top!
1701 snprintf(query,200, "%f*(%s-0.5):%s", 1.+igr, oaAlias->At(i)->GetName(), var_x);
1702 multGr->Add( (TGraphErrors*) TStatToolkit::MakeGraphSparse(tree,query,cut,marArr[i],colArr[i],sizArr[i]) );
1704 snprintf(query,200, "%f*(%s-0.5):%s", 1.+igr, oaAlias->At(2)->GetName(), var_x);
1705 multGr->Add( (TGraphErrors*) TStatToolkit::MakeGraphSparse(tree,query,cut,marArr[2],colArr[2],sizArr[2]) );
1707 multGr->SetName(varname);
1708 multGr->SetTitle(varname); // used for y-axis labels. // details to be included!
1715 void TStatToolkit::AddStatusPad(TCanvas* c1, Float_t padratio, Float_t bottommargin)
1718 // add pad to bottom of canvas for Status graphs (by Patrick Reichelt)
1719 // call function "DrawStatusGraphs(...)" afterwards
1721 TCanvas* c1_clone = (TCanvas*) c1->Clone("c1_clone");
1725 TPad* pad1 = new TPad("pad1", "pad1", 0., padratio, 1., 1.);
1727 pad1->SetNumber(1); // so it can be called via "c1->cd(1);"
1729 TPad* pad2 = new TPad("pad2", "pad2", 0., 0., 1., padratio);
1732 // draw original canvas into first pad
1734 c1_clone->DrawClonePad();
1735 pad1->SetBottomMargin(0.001);
1736 pad1->SetRightMargin(0.01);
1737 // set up second pad
1740 pad2->SetTopMargin(0);
1741 pad2->SetBottomMargin(bottommargin); // for the long x-axis labels (runnumbers)
1742 pad2->SetRightMargin(0.01);
1746 void TStatToolkit::DrawStatusGraphs(TObjArray* oaMultGr)
1749 // draw Status graphs into active pad of canvas (by MI, Patrick Reichelt)
1750 // ...into bottom pad, if called after "AddStatusPad(...)"
1752 const Int_t nvars = oaMultGr->GetEntriesFast();
1753 TGraph* grAxis = (TGraph*) ((TMultiGraph*) oaMultGr->At(0))->GetListOfGraphs()->At(0);
1754 grAxis->SetMaximum(0.5*nvars+0.5);
1755 grAxis->SetMinimum(0);
1756 grAxis->GetYaxis()->SetLabelSize(0);
1757 Int_t entries = grAxis->GetN();
1758 printf("entries (via GetN()) = %d\n",entries);
1759 grAxis->GetXaxis()->SetLabelSize(5.7*TMath::Min(TMath::Max(5./entries,0.01),0.03));
1760 grAxis->GetXaxis()->LabelsOption("v");
1763 // draw multigraphs & names of status variables on the y axis
1764 for (Int_t i=0; i<nvars; i++){
1765 ((TMultiGraph*) oaMultGr->At(i))->Draw("p");
1766 TLatex* ylabel = new TLatex(-0.1, 0.5*i+0.5, ((TMultiGraph*) oaMultGr->At(i))->GetTitle());
1767 ylabel->SetTextAlign(32); //hor:right & vert:centered
1768 ylabel->SetTextSize(0.025/gPad->GetHNDC());
1774 TH1* TStatToolkit::DrawHistogram(TTree * tree, const char* drawCommand, const char* cuts, const char* histoname, const char* histotitle, Int_t nsigma, Float_t fraction )
1777 // Draw histogram from TTree with robust range
1778 // Only for 1D so far!
1781 // - histoname: name of histogram
1782 // - histotitle: title of histgram
1783 // - fraction: fraction of data to define the robust mean
1784 // - nsigma: nsigma value for range
1787 TString drawStr(drawCommand);
1788 TString cutStr(cuts);
1792 cerr<<" Tree pointer is NULL!"<<endl;
1797 Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff");
1798 if (entries == -1) {
1799 cerr<<"TTree draw returns -1"<<endl;
1804 if(tree->GetV1()) dim = 1;
1805 if(tree->GetV2()) dim = 2;
1806 if(tree->GetV3()) dim = 3;
1808 cerr<<"TTree has more than 2 dimensions (not yet supported)"<<endl;
1813 Double_t meanX, rmsX=0;
1814 Double_t meanY, rmsY=0;
1815 TStatToolkit::EvaluateUni(entries, tree->GetV1(),meanX,rmsX, fraction*entries);
1817 TStatToolkit::EvaluateUni(entries, tree->GetV1(),meanY,rmsY, fraction*entries);
1818 TStatToolkit::EvaluateUni(entries, tree->GetV2(),meanX,rmsX, fraction*entries);
1822 hOut = new TH1F(histoname, histotitle, 200, meanX-nsigma*rmsX, meanX+nsigma*rmsX);
1823 for (Int_t i=0; i<entries; i++) hOut->Fill(tree->GetV1()[i]);
1824 hOut->GetXaxis()->SetTitle(tree->GetHistogram()->GetXaxis()->GetTitle());
1828 hOut = new TH2F(histoname, histotitle, 200, meanX-nsigma*rmsX, meanX+nsigma*rmsX,200, meanY-nsigma*rmsY, meanY+nsigma*rmsY);
1829 for (Int_t i=0; i<entries; i++) hOut->Fill(tree->GetV2()[i],tree->GetV1()[i]);
1830 hOut->GetXaxis()->SetTitle(tree->GetHistogram()->GetXaxis()->GetTitle());
1831 hOut->GetYaxis()->SetTitle(tree->GetHistogram()->GetYaxis()->GetTitle());