]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STAT/TStatToolkit.cxx
consolidate zero-length arrays (aka struct hack)
[u/mrichter/AliRoot.git] / STAT / TStatToolkit.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16
17 ///////////////////////////////////////////////////////////////////////////
18 // Class TStatToolkit
19 // 
20 // Subset of  matheamtical functions  not included in the TMath
21 //
22 //
23 /////////////////////////////////////////////////////////////////////////
24 #include "TMath.h"
25 #include "Riostream.h"
26 #include "TH1F.h"
27 #include "TH2F.h"
28 #include "TH3.h"
29 #include "TF1.h"
30 #include "TTree.h"
31 #include "TChain.h"
32 #include "TObjString.h"
33 #include "TLinearFitter.h"
34 #include "TGraph2D.h"
35 #include "TGraph.h"
36 #include "TGraphErrors.h"
37 #include "TMultiGraph.h"
38 #include "TCanvas.h"
39 #include "TLatex.h"
40 #include "TCut.h"
41 //
42 // includes neccessary for test functions
43 //
44 #include "TSystem.h"
45 #include "TRandom.h"
46 #include "TStopwatch.h"
47 #include "TTreeStream.h"
48
49 #include "TStatToolkit.h"
50
51
52 using std::cout;
53 using std::cerr;
54 using std::endl;
55
56  
57 ClassImp(TStatToolkit) // Class implementation to enable ROOT I/O
58  
59 TStatToolkit::TStatToolkit() : TObject()
60 {
61   //
62   // Default constructor
63   //
64 }
65 ///////////////////////////////////////////////////////////////////////////
66 TStatToolkit::~TStatToolkit()
67 {
68   //
69   // Destructor
70   //
71 }
72
73
74 //_____________________________________________________________________________
75 void TStatToolkit::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean
76                            , Double_t &sigma, Int_t hh)
77 {
78   //
79   // Robust estimator in 1D case MI version - (faster than ROOT version)
80   //
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
86   //
87
88   if (hh==0)
89     hh=(nvectors+2)/2;
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);
93   
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)];
96   
97   Double_t sumx  =0;
98   Double_t sumx2 =0;
99   Int_t    bestindex = -1;
100   Double_t bestmean  = 0; 
101   Double_t bestsigma = (data[index[nvectors-1]]-data[index[0]]+1.);   // maximal possible sigma
102   bestsigma *=bestsigma;
103
104   for (Int_t i=0; i<hh; i++){
105     sumx  += data[index[i]];
106     sumx2 += data[index[i]]*data[index[i]];
107   }
108   
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){
115       bestmean  = cmean;
116       bestsigma = csigma;
117       bestindex = i-hh;
118     }
119     
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]];
122   }
123   
124   Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
125   mean  = bestmean;
126   sigma = bstd;
127   delete [] index;
128
129 }
130
131
132
133 void TStatToolkit::EvaluateUniExternal(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh,  Float_t externalfactor)
134 {
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
138   // 
139
140   if (hh==0)
141     hh=(nvectors+2)/2;
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);
145   //
146   Int_t    nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
147   Double_t factor = faclts[0];
148   if (nquant>0){
149     // fix proper normalization - Anja
150     factor = faclts[nquant-1];
151   }
152
153   //
154   //
155   Double_t sumx  =0;
156   Double_t sumx2 =0;
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]];
163   }
164   //   
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){
171       bestmean  = cmean;
172       bestsigma = csigma;
173       bestindex = i-hh;
174     }
175     //
176     //
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]];
179   }
180   
181   Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
182   mean  = bestmean;
183   sigma = bstd;
184   delete [] index;
185 }
186
187
188 //_____________________________________________________________________________
189 Int_t TStatToolkit::Freq(Int_t n, const Int_t *inlist
190                         , Int_t *outlist, Bool_t down)
191 {    
192   //
193   //  Sort eleements according occurancy 
194   //  The size of output array has is 2*n 
195   //
196
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;
201   //
202   TMath::Sort(n,inlist, sindexS, down);  
203   Int_t last      = inlist[sindexS[0]];
204   Int_t val       = last;
205   sindexF[0]      = 1;
206   sindexF[0+n]    = last;
207   Int_t countPos  = 0;
208   //
209   //  find frequency
210   for(Int_t i=1;i<n; i++){
211     val = inlist[sindexS[i]];
212     if (last == val)   sindexF[countPos]++;
213     else{      
214       countPos++;
215       sindexF[countPos+n] = val;
216       sindexF[countPos]++;
217       last =val;
218     }
219   }
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]];
226   }
227   delete [] sindexS;
228   delete [] sindexF;
229   
230   return countPos;
231
232 }
233
234 //___TStatToolkit__________________________________________________________________________
235 void TStatToolkit::TruncatedMean(const TH1 * his, TVectorD *param, Float_t down, Float_t up, Bool_t verbose){
236   //
237   //
238   //
239   Int_t nbins    = his->GetNbinsX();
240   Float_t nentries = his->GetEntries();
241   Float_t sum      =0;
242   Float_t mean   = 0;
243   Float_t sigma2 = 0;
244   Float_t ncumul=0;  
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);      
252     }
253   }
254   mean/=sum;
255   sigma2= TMath::Sqrt(TMath::Abs(sigma2/sum-mean*mean));
256   if (param){
257     (*param)[0] = his->GetMaximum();
258     (*param)[1] = mean;
259     (*param)[2] = sigma2;
260     
261   }
262   if (verbose)  printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma2);
263 }
264
265 void TStatToolkit::LTM(TH1 * his, TVectorD *param , Float_t fraction,  Bool_t verbose){
266   //
267   // LTM : Trimmed mean on histogram - Modified version for binned data
268   //
269   // Robust statistic to estimate properties of the distribution
270   // See http://en.wikipedia.org/w/index.php?title=Trimmed_estimator&oldid=582847999
271   //
272   // New faster version is under preparation
273   //
274   if (!param) return;
275   (*param)[0]=0;
276   (*param)[1]=0;
277   (*param)[2]=0;  
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];
282   Int_t npoints=0;
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);
291         npoints++;
292       }
293     }
294   }
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);
299   delete [] data;
300   if (verbose)  printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma);if (param){
301     (*param)[0] = his->GetMaximum();
302     (*param)[1] = mean;
303     (*param)[2] = sigma;    
304   }
305 }
306
307
308 void TStatToolkit::MedianFilter(TH1 * his1D, Int_t nmedian){
309   //
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
314   // Input parameters:
315   //    his1D - input histogam - to be modiefied by Medianfilter
316   //    nmendian - number of bins in median filter
317   //
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);
328   }  
329 }
330
331 Bool_t TStatToolkit::LTMHisto(TH1 *his1D, TVectorD &params , Float_t fraction){
332   //
333   // LTM : Trimmed mean on histogram - Modified version for binned data
334   // 
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
339   //
340   // Function parameters:
341   //     his1D   - input histogram
342   //     params  - vector with parameters
343   //             - 0 - area
344   //             - 1 - mean
345   //             - 2 - rms 
346   //             - 3 - error estimate of mean
347   //             - 4 - error estimate of RMS
348   //             - 5 - first accepted bin position
349   //             - 6 - last accepted  bin position
350   //
351   Int_t nbins    = his1D->GetNbinsX();
352   Int_t nentries = (Int_t)his1D->GetEntries();
353   const Double_t kEpsilon=0.0000000001;
354
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);
361   Double_t sumCont=0;
362   for (Int_t ibin0=1; ibin0<=nbins; ibin0++) sumCont+=his1D->GetBinContent(ibin0);
363   //
364   Double_t minRMS=his1D->GetRMS()*10000;
365   Int_t maxBin=0;
366   //
367   for (Int_t ibin0=1; ibin0<nbins; ibin0++){
368     Double_t sum0=0, sum1=0, sum2=0;
369     Int_t ibin1=ibin0;
370     for ( ibin1=ibin0; ibin1<nbins; ibin1++){
371       Double_t cont=his1D->GetBinContent(ibin1);
372       Double_t x= his1D->GetBinCenter(ibin1);
373       sum0+=cont;
374       sum1+=cont*x;
375       sum2+=cont*x*x;
376       if ( (ibin0!=ibin1) && sum0>=fraction*sumCont) break;
377     }
378     vectorX[ibin0]=his1D->GetBinCenter(ibin0);
379     if (sum0<fraction*sumCont) continue;
380     //
381     // substract fractions of bin0 and bin1 to keep sum0=fration*sumCont
382     //
383     Double_t diff = sum0-fraction*sumCont;
384     Double_t mean = sum1/sum0;
385     //
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);
390     //
391     Double_t d = y0+y1-diff;    //enties to keep 
392     Double_t w0=0,w1=0;
393     if (y0<=kEpsilon&&y1>kEpsilon){
394       w1=d/y1;
395     } 
396     if (y1<=kEpsilon&&y0>kEpsilon){
397       w0=d/y0;
398     }
399     if (y0>kEpsilon && y1>kEpsilon && x1>x0  ){
400       w0 = (d*(x1-mean))/((x1-x0)*y0);
401       w1 = (d-y0*w0)/y1;
402       //
403       if (w0>1) {w1+=(w0-1)*y0/y1; w0=1;}
404       if (w1>1) {w0+=(w1-1)*y1/y0; w1=1;}
405     }  
406     if ( (x1>x0) &&TMath::Abs(y0*w0+y1*w1-d)>kEpsilon*sum0){
407       printf(" TStatToolkit::LTMHisto error\n");
408     }
409     sum0-=y0+y1;
410     sum1-=y0*x0;
411     sum1-=y1*x1;
412     sum2-=y0*x0*x0;
413     sum2-=y1*x1*x1;
414     //
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);
417     sum0+=y0*w0+y1*w1;
418     sum1+=y0*w0*xx0;
419     sum1+=y1*w1*xx1;
420     sum2+=y0*w0*xx0*xx0;
421     sum2+=y1*w1*xx1*xx1;
422
423     //
424     // choose the bin with smallest rms
425     //
426     if (sum0>0){
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];
431         params[0]=sum0;
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???
436         params[5]=ibin0;
437         params[6]=ibin1;
438         params[7]=his1D->GetBinCenter(ibin0);
439         params[8]=his1D->GetBinCenter(ibin1);
440         maxBin=ibin0;
441       }
442     }else{
443       break;
444     }
445   }
446   return kTRUE;
447 }
448
449
450 Double_t  TStatToolkit::FitGaus(TH1* his, TVectorD *param, TMatrixD */*matrix*/, Float_t xmin, Float_t xmax, Bool_t verbose){
451   //
452   //  Fit histogram with gaussian function
453   //  
454   //  Prameters:
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
459   //  Fitting:
460   //  1. Step - make logarithm
461   //  2. Linear  fit (parabola) - more robust - always converge
462   //  3. In case of small statistic bins are averaged
463   //  
464   static TLinearFitter fitter(3,"pol2");
465   TVectorD  par(3);
466   TVectorD  sigma(3);
467   TMatrixD mat(3,3);
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));
473
474   if (maxEstimate<1) return -1;
475   Int_t nbins    = his->GetNbinsX();
476   Int_t npoints=0;
477   //
478
479
480   if (xmin>=xmax){
481     xmin = his->GetXaxis()->GetXmin();
482     xmax = his->GetXaxis()->GetXmax();
483   }
484   for (Int_t iter=0; iter<2; iter++){
485     fitter.ClearPoints();
486     npoints=0;
487     for (Int_t ibin=1;ibin<nbins+1; ibin++){
488       Int_t countB=1;
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);
493           countB++;
494         }
495       }
496       entriesI/=countB;
497       Double_t xcenter= his->GetBinCenter(ibin);
498       if (xcenter<xmin || xcenter>xmax) continue;
499       Double_t error=1./TMath::Sqrt(countB);
500       Float_t   cont=2;
501       if (iter>0){
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));
505       }
506       if (entriesI>1&&cont>1){
507         fitter.AddPoint(&xcenter,TMath::Log(Float_t(entriesI)),error);
508         npoints++;
509       }
510     }  
511     if (npoints>3){
512       fitter.Eval();
513       fitter.GetParameters(par);
514     }else{
515       break;
516     }
517   }
518   if (npoints<=3){
519     return -1;
520   }
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]);
532   if (verbose){
533     par.Print();
534     mat.Print();
535     param->Print();
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]);    
541     f1->Draw("same");
542   }
543   return chi2;
544 }
545
546 Double_t  TStatToolkit::FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param, TMatrixD */*matrix*/, Bool_t verbose){
547   //
548   //  Fit histogram with gaussian function
549   //  
550   //  Prameters:
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!!!
555   //
556   //  Return values:
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!!
562   //
563   //  Fitting:
564   //  1. Step - make logarithm
565   //  2. Linear  fit (parabola) - more robust - always converge
566   //  
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();
572   TVectorD  par(3);
573   TVectorD  sigma(3);
574   TMatrixD matA(3,3);
575   TMatrixD b(3,1);
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;
579
580   Float_t meanCOG = 0;
581   Float_t rms2COG = 0;
582   Float_t sumCOG  = 0;
583
584   Float_t entries = 0;
585   Int_t nfilled=0;
586
587   for (Int_t i=0; i<nBins; i++){
588       entries+=arr[i];
589       if (arr[i]>0) nfilled++;
590   }
591
592   if (max<4) return -4;
593   if (entries<12) return -4;
594   if (rms<kTol) return -4;
595
596   Int_t npoints=0;
597   //
598
599   //
600   for (Int_t ibin=0;ibin<nBins; ibin++){
601       Float_t entriesI = arr[ibin];
602     if (entriesI>1){
603       Double_t xcenter = xMin+(ibin+0.5)*binWidth;
604       
605       Float_t error    = 1./TMath::Sqrt(entriesI);
606       Float_t val = TMath::Log(Float_t(entriesI));
607       fitter.AddPoint(&xcenter,val,error);
608       if (npoints<3){
609           matA(npoints,0)=1;
610           matA(npoints,1)=xcenter;
611           matA(npoints,2)=xcenter*xcenter;
612           b(npoints,0)=val;
613           meanCOG+=xcenter*entriesI;
614           rms2COG +=xcenter*entriesI*xcenter;
615           sumCOG +=entriesI;
616       }
617       npoints++;
618     }
619   }
620
621   
622   Double_t chi2 = 0;
623   if (npoints>=3){
624       if ( npoints == 3 ){
625           //analytic calculation of the parameters for three points
626           matA.Invert();
627           TMatrixD res(1,3);
628           res.Mult(matA,b);
629           par[0]=res(0,0);
630           par[1]=res(0,1);
631           par[2]=res(0,2);
632           chi2 = -3.;
633       } else {
634           // use fitter for more than three points
635           fitter.Eval();
636           fitter.GetParameters(par);
637           fitter.GetCovarianceMatrix(mat);
638           chi2 = fitter.GetChisquare()/Float_t(npoints);
639       }
640       if (TMath::Abs(par[1])<kTol) return -4;
641       if (TMath::Abs(par[2])<kTol) return -4;
642
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
645
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);
651       if (verbose){
652           par.Print();
653           mat.Print();
654           param->Print();
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]);
660           f1->Draw("same");
661       }
662       return chi2;
663   }
664
665   if (npoints == 2){
666       //use center of gravity for 2 points
667       meanCOG/=sumCOG;
668       rms2COG /=sumCOG;
669       (*param)[0] = max;
670       (*param)[1] = meanCOG;
671       (*param)[2] = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
672       chi2=-2.;
673   }
674   if ( npoints == 1 ){
675       meanCOG/=sumCOG;
676       (*param)[0] = max;
677       (*param)[1] = meanCOG;
678       (*param)[2] = binWidth/TMath::Sqrt(12);
679       chi2=-1.;
680   }
681   return chi2;
682
683 }
684
685
686 Float_t TStatToolkit::GetCOG(const Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms, Float_t *sum)
687 {
688     //
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
691     //
692     Float_t meanCOG = 0;
693     Float_t rms2COG = 0;
694     Float_t sumCOG  = 0;
695     Int_t npoints   = 0;
696
697     Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
698
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;
702         if ( entriesI>0 ){
703             meanCOG += xcenter*entriesI;
704             rms2COG += xcenter*entriesI*xcenter;
705             sumCOG  += entriesI;
706             npoints++;
707         }
708     }
709     if ( sumCOG == 0 ) return xMin;
710     meanCOG/=sumCOG;
711
712     if ( rms ){
713         rms2COG /=sumCOG;
714         (*rms) = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
715         if ( npoints == 1 ) (*rms) = binWidth/TMath::Sqrt(12);
716     }
717
718     if ( sum )
719         (*sum) = sumCOG;
720
721     return meanCOG;
722 }
723
724
725
726 ///////////////////////////////////////////////////////////////
727 //////////////         TEST functions /////////////////////////
728 ///////////////////////////////////////////////////////////////
729
730
731
732
733
734 void TStatToolkit::TestGausFit(Int_t nhistos){
735   //
736   // Test performance of the parabolic - gaussian fit - compare it with 
737   // ROOT gauss fit
738   //  nhistos - number of histograms to be used for test
739   //
740   TTreeSRedirector *pcstream = new TTreeSRedirector("fitdebug.root","recreate");
741   
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];
746   TMatrixD dummy(3,3);
747   
748   
749   TH1F **h1f = new TH1F*[nhistos];
750   TF1  *myg = new TF1("myg","gaus");
751   TF1  *fit = new TF1("fit","gaus");
752   gRandom->SetSeed(0);
753   
754   //init
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();
760     gSystem->Sleep(2);
761     sTrue[i]= .75+gRandom->Rndm()*.5;
762     myg->SetParameters(1,xTrue[i],sTrue[i]);
763     h1f[i]->FillRandom("myg");
764   }
765   
766   TStopwatch s; 
767   s.Start();
768   //standard gaus fit
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);
774   }
775   s.Stop();
776   printf("Gaussian fit\t");
777   s.Print();
778   
779   s.Start();
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);
783   }
784   
785   s.Stop();
786   printf("Parabolic fit\t");
787   s.Print();
788   //write stream
789   for (Int_t i=0;i<nhistos; i++){
790     Float_t xt  = xTrue[i];
791     Float_t st  = sTrue[i];
792     (*pcstream)<<"data"
793                <<"xTrue="<<xt
794                <<"sTrue="<<st
795                <<"pg.="<<(par1[i])
796                <<"pa.="<<(par2[i])
797                <<"\n";
798   }    
799   //delete pointers
800   for (Int_t i=0;i<nhistos; i++){
801     delete par1[i];
802     delete par2[i];
803     delete h1f[i];
804   }
805   delete pcstream;
806   delete []h1f;
807   delete []xTrue;
808   delete []sTrue;
809   //
810   delete []par1;
811   delete []par2;
812
813 }
814
815
816
817 TGraph2D * TStatToolkit::MakeStat2D(TH3 * his, Int_t delta0, Int_t delta1, Int_t type){
818   //
819   //
820   //
821   // delta - number of bins to integrate
822   // type - 0 - mean value
823
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();
829   char name[1000];
830   Int_t icount=0;
831   TGraph2D  *graph = new TGraph2D(nbinx*nbiny);
832   TF1 f1("f1","gaus");
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);
839       Float_t stat= 0;
840       if (type==0) stat = projection->GetMean();
841       if (type==1) stat = projection->GetRMS();
842       if (type==2 || type==3){
843         TVectorD vec(10);
844         TStatToolkit::LTM((TH1F*)projection,&vec,0.7);
845         if (type==2) stat= vec[1];
846         if (type==3) stat= vec[0];      
847       }
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);
852       }
853       //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
854       graph->SetPoint(icount,xcenter, ycenter, stat);
855       icount++;
856     }
857   return graph;
858 }
859
860 TGraphErrors * TStatToolkit::MakeStat1D(TH2 * his, Int_t deltaBin, Double_t fraction, Int_t returnType, Int_t markerStyle, Int_t markerColor){
861   //
862   // function to retrieve the "mean and RMS estimate" of 2D histograms
863   //     
864   // Robust statistic to estimate properties of the distribution
865   // See http://en.wikipedia.org/wiki/Trimmed_estimator
866   //
867   // deltaBin - number of bins to integrate (bin+-deltaBin)
868   // fraction - fraction of values for the LTM and for the gauss fit
869   // returnType - 
870   //        0 - mean value
871   //        1 - RMS
872   //        2 - LTM mean
873   //        3 - LTM sigma
874   //        4 - Gaus fit mean  - on LTM range
875   //        5 - Gaus fit sigma - on LTM  range
876   // 
877   TAxis * xaxis  = his->GetXaxis();
878   Int_t   nbinx  = xaxis->GetNbins();
879   char name[1000];
880   Int_t icount=0;
881   //
882   TVectorD vecX(nbinx);
883   TVectorD vecXErr(nbinx);
884   TVectorD vecY(nbinx);
885   TVectorD vecYErr(nbinx);
886   //
887   TF1 f1("f1","gaus");
888   TVectorD vecLTM(10);
889
890   for (Int_t jx=1; jx<=nbinx;jx++){
891     Int_t ix=jx-1;
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));
895     Double_t stat= 0;
896     Double_t err =0;
897     TStatToolkit::LTMHisto((TH1F*)projection,vecLTM,fraction);  
898     //
899     if (returnType==0) {
900       stat = projection->GetMean();
901       err  = projection->GetMeanError();
902     }
903     if (returnType==1) {
904       stat = projection->GetRMS();
905       err = projection->GetRMSError();
906     }
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();}
910     }
911     if (returnType==4|| returnType==5){
912       projection->Fit(&f1,"QN","QN", vecLTM[7], vecLTM[8]);
913       if (returnType==4) {
914         stat= f1.GetParameter(1);
915         err=f1.GetParError(1);
916       }
917       if (returnType==5) {
918         stat= f1.GetParameter(2);
919         err=f1.GetParError(2);
920       }
921     }
922     vecX[icount]=xcenter;
923     vecY[icount]=stat;
924     vecYErr[icount]=err;
925     icount++;
926     delete projection;
927   }
928   TGraphErrors  *graph = new TGraphErrors(icount,vecX.GetMatrixArray(), vecY.GetMatrixArray(),0, vecYErr.GetMatrixArray());
929   graph->SetMarkerStyle(markerStyle);
930   graph->SetMarkerColor(markerColor);
931   return graph;
932 }
933
934
935
936
937
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){
939    //
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
943    //
944
945    TString formulaStr(formula); 
946    TString drawStr(drawCommand);
947    TString cutStr(cuts);
948    TString ferr("1");
949
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();     
955      delete valTokens;
956    }
957
958       
959    formulaStr.ReplaceAll("++", "~");
960    TObjArray* formulaTokens = formulaStr.Tokenize("~"); 
961    Int_t dim = formulaTokens->GetEntriesFast();
962    
963    fitParam.ResizeTo(dim);
964    covMatrix.ResizeTo(dim,dim);
965    
966    TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
967    fitter->StoreData(kTRUE);   
968    fitter->ClearPoints();
969    
970    Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff",  stop-start, start);
971    if (entries == -1) {
972      delete formulaTokens;
973      return new TString("An ERROR has occured during fitting!");
974    }
975    Double_t **values = new Double_t*[dim+1] ;
976    for (Int_t i=0; i<dim+1; i++) values[i]=NULL; 
977    //
978    entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff",  stop-start, start);
979    if (entries == -1) {
980      delete formulaTokens;
981      delete []values;
982      return new TString("An ERROR has occured during fitting!");
983    }
984    Double_t *errors = new Double_t[entries];
985    memcpy(errors,  tree->GetV1(), entries*sizeof(Double_t));
986    
987    for (Int_t i = 0; i < dim + 1; i++){
988       Int_t centries = 0;
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);
991       
992       if (entries != centries) {
993         delete []errors;
994         delete []values;
995         return new TString("An ERROR has occured during fitting!");
996       }
997       values[i] = new Double_t[entries];
998       memcpy(values[i],  tree->GetV1(), entries*sizeof(Double_t)); 
999    }
1000    
1001    // add points to the fitter
1002    for (Int_t i = 0; i < entries; i++){
1003       Double_t x[1000];
1004       for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
1005       fitter->AddPoint(x, values[dim][i], errors[i]);
1006    }
1007
1008    fitter->Eval();
1009    if (frac>0.5 && frac<1){
1010      fitter->EvalRobust(frac);
1011    }else{
1012      if (fix0) {
1013        fitter->FixParameter(0,0);
1014        fitter->Eval();     
1015      }
1016    }
1017    fitter->GetParameters(fitParam);
1018    fitter->GetCovarianceMatrix(covMatrix);
1019    chi2 = fitter->GetChisquare();
1020    npoints = entries;   
1021    TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula; 
1022    
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("+");
1026    }
1027    returnFormula.Append(" )");
1028    
1029    
1030    for (Int_t j=0; j<dim+1;j++) delete [] values[j];
1031
1032
1033    delete formulaTokens;
1034    delete fitter;
1035    delete[] values;
1036    delete[] errors;
1037    return preturnFormula;
1038 }
1039
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){
1041    //
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
1045    //
1046
1047    TString formulaStr(formula); 
1048    TString drawStr(drawCommand);
1049    TString cutStr(cuts);
1050    TString ferr("1");
1051
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();     
1057      delete valTokens;
1058    }
1059
1060       
1061    formulaStr.ReplaceAll("++", "~");
1062    TObjArray* formulaTokens = formulaStr.Tokenize("~"); 
1063    Int_t dim = formulaTokens->GetEntriesFast();
1064    
1065    fitParam.ResizeTo(dim);
1066    covMatrix.ResizeTo(dim,dim);
1067    
1068    TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
1069    fitter->StoreData(kTRUE);   
1070    fitter->ClearPoints();
1071    
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!");
1076    }
1077    Double_t **values = new Double_t*[dim+1] ; 
1078    for (Int_t i=0; i<dim+1; i++) values[i]=NULL; 
1079    //
1080    entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff",  stop-start, start);
1081    if (entries == -1) {
1082      delete formulaTokens;
1083      delete [] values;
1084      return new TString("An ERROR has occured during fitting!");
1085    }
1086    Double_t *errors = new Double_t[entries];
1087    memcpy(errors,  tree->GetV1(), entries*sizeof(Double_t));
1088    
1089    for (Int_t i = 0; i < dim + 1; i++){
1090       Int_t centries = 0;
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);
1093       
1094       if (entries != centries) {
1095         delete []errors;
1096         delete []values;
1097         delete formulaTokens;
1098         return new TString("An ERROR has occured during fitting!");
1099       }
1100       values[i] = new Double_t[entries];
1101       memcpy(values[i],  tree->GetV1(), entries*sizeof(Double_t)); 
1102    }
1103    
1104    // add points to the fitter
1105    for (Int_t i = 0; i < entries; i++){
1106       Double_t x[1000];
1107       for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
1108       fitter->AddPoint(x, values[dim][i], errors[i]);
1109    }
1110    if (constrain>0){
1111      for (Int_t i = 0; i < dim; i++){
1112        Double_t x[1000];
1113        for (Int_t j=0; j<dim;j++) if (i!=j) x[j]=0;
1114        x[i]=1.;
1115        fitter->AddPoint(x, 0, constrain);
1116      }
1117    }
1118
1119
1120    fitter->Eval();
1121    if (frac>0.5 && frac<1){
1122      fitter->EvalRobust(frac);   
1123    }
1124    fitter->GetParameters(fitParam);
1125    fitter->GetCovarianceMatrix(covMatrix);
1126    chi2 = fitter->GetChisquare();
1127    npoints = entries;
1128    
1129    TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula; 
1130    
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("+");
1134    }
1135    returnFormula.Append(" )");
1136    
1137    for (Int_t j=0; j<dim+1;j++) delete [] values[j];
1138    
1139
1140
1141    delete formulaTokens;
1142    delete fitter;
1143    delete[] values;
1144    delete[] errors;
1145    return preturnFormula;
1146 }
1147
1148
1149
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){
1151    //
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
1155    //
1156
1157    TString formulaStr(formula); 
1158    TString drawStr(drawCommand);
1159    TString cutStr(cuts);
1160    TString ferr("1");
1161
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();
1167      delete valTokens;
1168    }
1169
1170       
1171    formulaStr.ReplaceAll("++", "~");
1172    TObjArray* formulaTokens = formulaStr.Tokenize("~"); 
1173    Int_t dim = formulaTokens->GetEntriesFast();
1174    
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();
1182    
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!");
1187    }
1188    Double_t **values = new Double_t*[dim+1] ; 
1189    for (Int_t i=0; i<dim+1; i++) values[i]=NULL; 
1190    //
1191    entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff",  stop-start, start);
1192    if (entries == -1) {
1193      delete []values;
1194      delete formulaTokens;
1195      return new TString("An ERROR has occured during fitting!");
1196    }
1197    Double_t *errors = new Double_t[entries];
1198    memcpy(errors,  tree->GetV1(), entries*sizeof(Double_t));
1199    
1200    for (Int_t i = 0; i < dim + 1; i++){
1201       Int_t centries = 0;
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);
1204       
1205       if (entries != centries) {
1206         delete []errors;
1207         delete []values;
1208         delete formulaTokens;
1209         return new TString("An ERROR has occured during fitting!");
1210       }
1211       values[i] = new Double_t[entries];
1212       memcpy(values[i],  tree->GetV1(), entries*sizeof(Double_t)); 
1213    }
1214    
1215    // add points to the fitter
1216    for (Int_t i = 0; i < entries; i++){
1217       Double_t x[1000];
1218       for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
1219       fitter->AddPoint(x, values[dim][i], errors[i]);
1220    }
1221
1222    fitter->Eval();
1223    if (frac>0.5 && frac<1){
1224      fitter->EvalRobust(frac);
1225    }
1226    fitter->GetParameters(fitParam);
1227    fitter->GetCovarianceMatrix(covMatrix);
1228    chi2 = fitter->GetChisquare();
1229    npoints = entries;
1230    
1231    TString *preturnFormula = new TString("("), &returnFormula = *preturnFormula; 
1232    
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("+");
1236    }
1237    returnFormula.Append(" )");
1238    
1239    
1240    for (Int_t j=0; j<dim+1;j++) delete [] values[j];
1241    
1242    delete formulaTokens;
1243    delete fitter;
1244    delete[] values;
1245    delete[] errors;
1246    return preturnFormula;
1247 }
1248
1249
1250
1251
1252
1253 Int_t TStatToolkit::GetFitIndex(const TString fString, const TString subString){
1254   //
1255   // fitString - ++ separated list of fits
1256   // substring - ++ separated list of the requiered substrings
1257   //
1258   // return the last occurance of substring in fit string
1259   // 
1260   TObjArray *arrFit = fString.Tokenize("++");
1261   TObjArray *arrSub = subString.Tokenize("++");
1262   Int_t index=-1;
1263   for (Int_t i=0; i<arrFit->GetEntries(); i++){
1264     Bool_t isOK=kTRUE;
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;
1268     }
1269     if (isOK) index=i;
1270   }
1271   delete arrFit;
1272   delete arrSub;
1273   return index;
1274 }
1275
1276
1277 TString  TStatToolkit::FilterFit(const TString &input, const TString filter, TVectorD &param, TMatrixD & covar){
1278   //
1279   // Filter fit expression make sub-fit
1280   //
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++){
1286     Bool_t isOK=kTRUE;
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;      
1290     }
1291     if (isOK) {
1292       result+="+"+str;
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());    
1295     }
1296   }
1297   result+="-0.)";
1298   delete array0;
1299   delete array1;
1300   return result;
1301 }
1302
1303 void TStatToolkit::Update1D(Double_t delta, Double_t sigma, Int_t s1, TMatrixD &vecXk, TMatrixD &covXk){
1304   //
1305   // Update parameters and covariance - with one measurement
1306   // Input:
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();
1312  
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);
1323   vecZk(0,0)=delta;
1324   measR(0,0)=sigma*sigma;
1325   //
1326   // reset matHk
1327   for (Int_t iel=0;iel<knElem;iel++) 
1328     for (Int_t ip=0;ip<knMeas;ip++) matHk(ip,iel)=0; 
1329   //mat1
1330   for (Int_t iel=0;iel<knElem;iel++) {
1331     for (Int_t jel=0;jel<knElem;jel++) mat1(iel,jel)=0;
1332     mat1(iel,iel)=1;
1333   }
1334   //
1335   matHk(0, s1)=1;
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
1339   matSk.Invert();
1340   matKk = (covXk*matHkT)*matSk;            //  Optimal Kalman gain
1341   vecXk += matKk*vecYk;                    //  updated vector 
1342   covXk2= (mat1-(matKk*matHk));
1343   covXk3 =  covXk2*covXk;          
1344   covXk = covXk3;  
1345   Int_t nrows=covXk3.GetNrows();
1346   
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; 
1351     }
1352 }
1353
1354
1355
1356 void   TStatToolkit::Constrain1D(const TString &input, const TString filter, TVectorD &param, TMatrixD & covar, Double_t mean, Double_t sigma){
1357   //
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
1363   //
1364   
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);}
1369   
1370   if (filter.Length()==0){
1371     TStatToolkit::Update1D(mean, sigma, 0, paramM, covar);//
1372   }else{  
1373     for (Int_t i=0; i<array0->GetEntries(); i++){
1374       Bool_t isOK=kTRUE;
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;      
1378       }
1379       if (isOK) {
1380         TStatToolkit::Update1D(mean, sigma, i+1, paramM, covar);//
1381       }
1382     }
1383   }
1384   for (Int_t i=0; i<=array0->GetEntries(); i++){
1385     param(i)=paramM(i,0);
1386   }
1387   delete array0;
1388   delete array1;
1389 }
1390
1391 TString  TStatToolkit::MakeFitString(const TString &input, const TVectorD &param, const TMatrixD & covar, Bool_t verbose){
1392   //
1393   //
1394   //
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());
1400     result+="+"+str;
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());    
1403   }
1404   result+="-0.)";
1405   delete array0;
1406   return result;
1407 }
1408
1409 TGraphErrors * TStatToolkit::MakeGraphErrors(TTree * tree, const char * expr, const char * cut,  Int_t mstyle, Int_t mcolor, Float_t msize, Float_t offset){
1410   //
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
1415   // expr   - examp 
1416   const Int_t entries =  tree->Draw(expr,cut,"goff");
1417   if (entries<=0) {
1418     TStatToolkit t;
1419     t.Error("TStatToolkit::MakeGraphError",Form("Empty or Not valid expression (%s) or cut *%s)", expr,cut));
1420     return 0;
1421   }
1422   if (  tree->GetV2()==0){
1423     TStatToolkit t;
1424     t.Error("TStatToolkit::MakeGraphError",Form("Not valid expression (%s) ", expr));
1425     return 0;
1426   }
1427   TGraphErrors * graph=0;
1428   if ( tree->GetV3()!=0){
1429     graph = new TGraphErrors (entries, tree->GetV2(),tree->GetV1(),0,tree->GetV3());
1430   }else{
1431     graph = new TGraphErrors (entries, tree->GetV2(),tree->GetV1(),0,0);
1432   }
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());
1441   delete charray;
1442   if (msize>0) graph->SetMarkerSize(msize);
1443   for(Int_t i=0;i<graph->GetN();i++) graph->GetX()[i]+=offset;
1444   return graph;
1445   
1446 }
1447
1448
1449 TGraph * TStatToolkit::MakeGraphSparse(TTree * tree, const char * expr, const char * cut, Int_t mstyle, Int_t mcolor, Float_t msize, Float_t offset){
1450   //
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
1454   //
1455   // Written by Weilin.Yu
1456   // updated & merged with QA-code by Patrick Reichelt
1457   //
1458   const Int_t entries = tree->Draw(expr,cut,"goff");
1459   if (entries<=0) {
1460     TStatToolkit t;
1461     t.Error("TStatToolkit::MakeGraphSparse",Form("Empty or Not valid expression (%s) or cut (%s)", expr, cut));
1462     return 0;
1463   }
1464   //  TGraph * graph = (TGraph*)gPad->GetPrimitive("Graph"); // 2D
1465
1466   Double_t *graphY, *graphX;
1467   graphY = tree->GetV1();
1468   graphX = tree->GetV2();
1469
1470   // sort according to run number
1471   Int_t *index = new Int_t[entries*4];
1472   TMath::Sort(entries,graphX,index,kFALSE);
1473
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;
1478
1479   // evaluate arrays for the new graph according to the run-number
1480   Int_t icount=0;
1481   //first entry
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++)
1486   {
1487     if(graphX[index[i]]==graphX[index[i-1]])
1488       unsortedX[index[i]] = count;
1489     else if(graphX[index[i]]!=graphX[index[i-1]]){
1490       count++;
1491       icount++;
1492       unsortedX[index[i]] = count;
1493       runNumber[icount]=graphX[index[i]];
1494     }
1495   }
1496
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++){
1501     newBins[i] = i;
1502   }
1503
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());
1509     }
1510     else { graphNew = new TGraphErrors(entries,unsortedX,graphY,0,tree->GetV3()); }
1511   }
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);
1515
1516   // set the bins for the x-axis, apply shifting of points
1517   Char_t xName[50];
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;
1522   }
1523
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;
1530   delete [] index;
1531   delete [] newBins;
1532   // 
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());
1538   delete charray;
1539   return graphNew;
1540 }
1541
1542
1543
1544 //
1545 // functions used for the trending
1546 //
1547
1548 Int_t  TStatToolkit::MakeStatAlias(TTree * tree, const char * expr, const char * cut, const char * alias) 
1549 {
1550   //
1551   // Add alias using statistical values of a given variable.
1552   // (by MI, Patrick Reichelt)
1553   //
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)
1559   // 
1560   /* Example usage:
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"
1563
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)    
1571   */
1572   // 
1573   Int_t entries = tree->Draw(expr,cut,"goff");
1574   if (entries<=1){
1575     printf("Expression or cut not valid:\t%s\t%s\n", expr, cut);
1576     return 0;
1577   }
1578   //
1579   TObjArray* oaAlias = TString(alias).Tokenize(":");
1580   if (oaAlias->GetEntries()<2) return 0;
1581   Float_t entryFraction = atof( oaAlias->At(1)->GetName() );
1582   //
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);
1588   //
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));
1594   delete oaAlias; 
1595   return entries;
1596 }
1597
1598 Int_t  TStatToolkit::SetStatusAlias(TTree * tree, const char * expr, const char * cut, const char * alias) 
1599 {
1600   //
1601   // Add alias to trending tree using statistical values of a given variable.
1602   // (by MI, Patrick Reichelt)
1603   //
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)
1610   //
1611   /* Example usage:
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)
1623   */
1624   //
1625   Int_t entries = tree->Draw(expr,cut,"goff");
1626   if (entries<1){
1627     printf("Expression or cut not valid:\t%s\t%s\n", expr, cut);
1628     return 0;
1629   }
1630   //
1631   TObjArray* oaVar = TString(expr).Tokenize(":");
1632   char varname[50];
1633   snprintf(varname,50,"%s", oaVar->At(0)->GetName());
1634   //
1635   TObjArray* oaAlias = TString(alias).Tokenize(":");
1636   if (oaAlias->GetEntries()<3) return 0;
1637   Float_t entryFraction = atof( oaAlias->At(2)->GetName() );
1638   //
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);
1644   //
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());
1657   //
1658   char query[200];
1659   char aname[200];
1660   snprintf(query,200,"%s", sQuery.Data());
1661   snprintf(aname,200,"%s", sAlias.Data());
1662   tree->SetAlias(aname, query);
1663   delete oaVar;
1664   delete oaAlias;
1665   return entries;
1666 }
1667
1668 TMultiGraph*  TStatToolkit::MakeStatusMultGr(TTree * tree, const char * expr, const char * cut, const char * alias, Int_t igr) 
1669 {
1670   //
1671   // Compute a trending multigraph that shows for which runs a variable has outliers.
1672   // (by MI, Patrick Reichelt)
1673   //
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.
1680   //
1681   TObjArray* oaVar = TString(expr).Tokenize(":");
1682   if (oaVar->GetEntries()<2) return 0;
1683   char varname[50];
1684   char var_x[50];
1685   snprintf(varname,50,"%s", oaVar->At(0)->GetName());
1686   snprintf(var_x  ,50,"%s", oaVar->At(1)->GetName());
1687   //
1688   TString sAlias(alias);
1689   sAlias.ReplaceAll("varname",varname);
1690   TObjArray* oaAlias = TString(sAlias.Data()).Tokenize(":");
1691   if (oaAlias->GetEntries()<3) return 0;
1692   //
1693   char query[200];
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]) );
1703   }
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]) );
1706   //
1707   multGr->SetName(varname);
1708   multGr->SetTitle(varname); // used for y-axis labels. // details to be included!
1709   delete oaVar;
1710   delete oaAlias;
1711   return multGr;
1712 }
1713
1714
1715 void  TStatToolkit::AddStatusPad(TCanvas* c1, Float_t padratio, Float_t bottommargin)
1716 {
1717   //
1718   // add pad to bottom of canvas for Status graphs (by Patrick Reichelt)
1719   // call function "DrawStatusGraphs(...)" afterwards
1720   //
1721   TCanvas* c1_clone = (TCanvas*) c1->Clone("c1_clone");
1722   c1->Clear();
1723   // produce new pads
1724   c1->cd();
1725   TPad* pad1 = new TPad("pad1", "pad1", 0., padratio, 1., 1.); 
1726   pad1->Draw();
1727   pad1->SetNumber(1); // so it can be called via "c1->cd(1);"
1728   c1->cd();
1729   TPad* pad2 = new TPad("pad2", "pad2", 0., 0., 1., padratio);
1730   pad2->Draw();
1731   pad2->SetNumber(2);
1732   // draw original canvas into first pad
1733   c1->cd(1);
1734   c1_clone->DrawClonePad();
1735   pad1->SetBottomMargin(0.001);
1736   pad1->SetRightMargin(0.01);
1737   // set up second pad
1738   c1->cd(2);
1739   pad2->SetGrid(3);
1740   pad2->SetTopMargin(0);
1741   pad2->SetBottomMargin(bottommargin); // for the long x-axis labels (runnumbers)
1742   pad2->SetRightMargin(0.01);
1743 }
1744
1745
1746 void  TStatToolkit::DrawStatusGraphs(TObjArray* oaMultGr)
1747 {
1748   //
1749   // draw Status graphs into active pad of canvas (by MI, Patrick Reichelt)
1750   // ...into bottom pad, if called after "AddStatusPad(...)"
1751   //
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");
1761   grAxis->Draw("ap");
1762   //
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());
1769     ylabel->Draw();
1770   }
1771 }
1772
1773
1774 TH1* TStatToolkit::DrawHistogram(TTree * tree, const char* drawCommand, const char* cuts, const char* histoname, const char* histotitle, Int_t nsigma, Float_t fraction )
1775 {
1776   //
1777   // Draw histogram from TTree with robust range
1778   // Only for 1D so far!
1779   // 
1780   // Parameters:
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
1785   //
1786
1787    TString drawStr(drawCommand);
1788    TString cutStr(cuts);
1789    Int_t dim = 1;
1790
1791    if(!tree) {
1792      cerr<<" Tree pointer is NULL!"<<endl;
1793      return 0;
1794    }
1795
1796    // get entries
1797    Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff");
1798    if (entries == -1) {
1799      cerr<<"TTree draw returns -1"<<endl;
1800      return 0;
1801    }
1802
1803    // get dimension
1804    if(tree->GetV1()) dim = 1;
1805    if(tree->GetV2()) dim = 2;
1806    if(tree->GetV3()) dim = 3;
1807    if(dim > 2){
1808      cerr<<"TTree has more than 2 dimensions (not yet supported)"<<endl;
1809      return 0;
1810    }
1811
1812    // draw robust
1813    Double_t meanX, rmsX=0;
1814    Double_t meanY, rmsY=0;
1815    TStatToolkit::EvaluateUni(entries, tree->GetV1(),meanX,rmsX, fraction*entries);
1816    if(dim==2){
1817      TStatToolkit::EvaluateUni(entries, tree->GetV1(),meanY,rmsY, fraction*entries);
1818      TStatToolkit::EvaluateUni(entries, tree->GetV2(),meanX,rmsX, fraction*entries);
1819    }
1820    TH1* hOut;
1821    if(dim==1){
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());
1825      hOut->Draw();
1826    }
1827    else if(dim==2){
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());
1832      hOut->Draw("colz");
1833    }
1834    return hOut;
1835 }