]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STAT/TStatToolkit.cxx
bug fix - intorduced with coverity changes
[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 "TH3.h"
28 #include "TF1.h"
29 #include "TTree.h"
30 #include "TChain.h"
31 #include "TObjString.h"
32 #include "TLinearFitter.h"
33
34 //
35 // includes neccessary for test functions
36 //
37 #include "TSystem.h"
38 #include "TRandom.h"
39 #include "TStopwatch.h"
40 #include "TTreeStream.h"
41
42 #include "TStatToolkit.h"
43
44  
45 ClassImp(TStatToolkit) // Class implementation to enable ROOT I/O
46  
47 TStatToolkit::TStatToolkit() : TObject()
48 {
49   //
50   // Default constructor
51   //
52 }
53 ///////////////////////////////////////////////////////////////////////////
54 TStatToolkit::~TStatToolkit()
55 {
56   //
57   // Destructor
58   //
59 }
60
61
62 //_____________________________________________________________________________
63 void TStatToolkit::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean
64                            , Double_t &sigma, Int_t hh)
65 {
66   //
67   // Robust estimator in 1D case MI version - (faster than ROOT version)
68   //
69   // For the univariate case
70   // estimates of location and scatter are returned in mean and sigma parameters
71   // the algorithm works on the same principle as in multivariate case -
72   // it finds a subset of size hh with smallest sigma, and then returns mean and
73   // sigma of this subset
74   //
75
76   if (hh==0)
77     hh=(nvectors+2)/2;
78   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};
79   Int_t *index=new Int_t[nvectors];
80   TMath::Sort(nvectors, data, index, kFALSE);
81   
82   Int_t    nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
83   Double_t factor = faclts[TMath::Max(0,nquant-1)];
84   
85   Double_t sumx  =0;
86   Double_t sumx2 =0;
87   Int_t    bestindex = -1;
88   Double_t bestmean  = 0; 
89   Double_t bestsigma = (data[index[nvectors-1]]-data[index[0]]+1.);   // maximal possible sigma
90   bestsigma *=bestsigma;
91
92   for (Int_t i=0; i<hh; i++){
93     sumx  += data[index[i]];
94     sumx2 += data[index[i]]*data[index[i]];
95   }
96   
97   Double_t norm = 1./Double_t(hh);
98   Double_t norm2 = 1./Double_t(hh-1);
99   for (Int_t i=hh; i<nvectors; i++){
100     Double_t cmean  = sumx*norm;
101     Double_t csigma = (sumx2 - hh*cmean*cmean)*norm2;
102     if (csigma<bestsigma){
103       bestmean  = cmean;
104       bestsigma = csigma;
105       bestindex = i-hh;
106     }
107     
108     sumx  += data[index[i]]-data[index[i-hh]];
109     sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
110   }
111   
112   Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
113   mean  = bestmean;
114   sigma = bstd;
115   delete [] index;
116
117 }
118
119
120
121 void TStatToolkit::EvaluateUniExternal(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh,  Float_t externalfactor)
122 {
123   // Modified version of ROOT robust EvaluateUni
124   // robust estimator in 1D case MI version
125   // added external factor to include precision of external measurement
126   // 
127
128   if (hh==0)
129     hh=(nvectors+2)/2;
130   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};
131   Int_t *index=new Int_t[nvectors];
132   TMath::Sort(nvectors, data, index, kFALSE);
133   //
134   Int_t    nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
135   Double_t factor = faclts[0];
136   if (nquant>0){
137     // fix proper normalization - Anja
138     factor = faclts[nquant-1];
139   }
140
141   //
142   //
143   Double_t sumx  =0;
144   Double_t sumx2 =0;
145   Int_t    bestindex = -1;
146   Double_t bestmean  = 0; 
147   Double_t bestsigma = -1;
148   for (Int_t i=0; i<hh; i++){
149     sumx  += data[index[i]];
150     sumx2 += data[index[i]]*data[index[i]];
151   }
152   //   
153   Double_t kfactor = 2.*externalfactor - externalfactor*externalfactor;
154   Double_t norm = 1./Double_t(hh);
155   for (Int_t i=hh; i<nvectors; i++){
156     Double_t cmean  = sumx*norm;
157     Double_t csigma = (sumx2*norm - cmean*cmean*kfactor);
158     if (csigma<bestsigma ||  bestsigma<0){
159       bestmean  = cmean;
160       bestsigma = csigma;
161       bestindex = i-hh;
162     }
163     //
164     //
165     sumx  += data[index[i]]-data[index[i-hh]];
166     sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
167   }
168   
169   Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
170   mean  = bestmean;
171   sigma = bstd;
172   delete [] index;
173 }
174
175
176 //_____________________________________________________________________________
177 Int_t TStatToolkit::Freq(Int_t n, const Int_t *inlist
178                         , Int_t *outlist, Bool_t down)
179 {    
180   //
181   //  Sort eleements according occurancy 
182   //  The size of output array has is 2*n 
183   //
184
185   Int_t * sindexS = new Int_t[n];     // temp array for sorting
186   Int_t * sindexF = new Int_t[2*n];   
187   for (Int_t i=0;i<n;i++) sindexS[i]=0;
188   for (Int_t i=0;i<2*n;i++) sindexF[i]=0;
189   //
190   TMath::Sort(n,inlist, sindexS, down);  
191   Int_t last      = inlist[sindexS[0]];
192   Int_t val       = last;
193   sindexF[0]      = 1;
194   sindexF[0+n]    = last;
195   Int_t countPos  = 0;
196   //
197   //  find frequency
198   for(Int_t i=1;i<n; i++){
199     val = inlist[sindexS[i]];
200     if (last == val)   sindexF[countPos]++;
201     else{      
202       countPos++;
203       sindexF[countPos+n] = val;
204       sindexF[countPos]++;
205       last =val;
206     }
207   }
208   if (last==val) countPos++;
209   // sort according frequency
210   TMath::Sort(countPos, sindexF, sindexS, kTRUE);
211   for (Int_t i=0;i<countPos;i++){
212     outlist[2*i  ] = sindexF[sindexS[i]+n];
213     outlist[2*i+1] = sindexF[sindexS[i]];
214   }
215   delete [] sindexS;
216   delete [] sindexF;
217   
218   return countPos;
219
220 }
221
222 //___TStatToolkit__________________________________________________________________________
223 void TStatToolkit::TruncatedMean(TH1F * his, TVectorD *param, Float_t down, Float_t up, Bool_t verbose){
224   //
225   //
226   //
227   Int_t nbins    = his->GetNbinsX();
228   Float_t nentries = his->GetEntries();
229   Float_t sum      =0;
230   Float_t mean   = 0;
231   Float_t sigma2 = 0;
232   Float_t ncumul=0;  
233   for (Int_t ibin=1;ibin<nbins; ibin++){
234     ncumul+= his->GetBinContent(ibin);
235     Float_t fraction = Float_t(ncumul)/Float_t(nentries);
236     if (fraction>down && fraction<up){
237       sum+=his->GetBinContent(ibin);
238       mean+=his->GetBinCenter(ibin)*his->GetBinContent(ibin);
239       sigma2+=his->GetBinCenter(ibin)*his->GetBinCenter(ibin)*his->GetBinContent(ibin);      
240     }
241   }
242   mean/=sum;
243   sigma2= TMath::Sqrt(TMath::Abs(sigma2/sum-mean*mean));
244   if (param){
245     (*param)[0] = his->GetMaximum();
246     (*param)[1] = mean;
247     (*param)[2] = sigma2;
248     
249   }
250   if (verbose)  printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma2);
251 }
252
253 void TStatToolkit::LTM(TH1F * his, TVectorD *param , Float_t fraction,  Bool_t verbose){
254   //
255   // LTM
256   //
257   Int_t nbins    = his->GetNbinsX();
258   Int_t nentries = (Int_t)his->GetEntries();
259   Double_t *data  = new Double_t[nentries];
260   Int_t npoints=0;
261   for (Int_t ibin=1;ibin<nbins; ibin++){
262     Float_t entriesI = his->GetBinContent(ibin);
263     Float_t xcenter= his->GetBinCenter(ibin);
264     for (Int_t ic=0; ic<entriesI; ic++){
265       if (npoints<nentries){
266         data[npoints]= xcenter;
267         npoints++;
268       }
269     }
270   }
271   Double_t mean, sigma;
272   Int_t npoints2=TMath::Min(Int_t(fraction*Float_t(npoints)),npoints-1);
273   npoints2=TMath::Max(Int_t(0.5*Float_t(npoints)),npoints2);
274   TStatToolkit::EvaluateUni(npoints, data, mean,sigma,npoints2);
275   delete [] data;
276   if (verbose)  printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma);if (param){
277     (*param)[0] = his->GetMaximum();
278     (*param)[1] = mean;
279     (*param)[2] = sigma;    
280   }
281 }
282
283 Double_t  TStatToolkit::FitGaus(TH1F* his, TVectorD *param, TMatrixD */*matrix*/, Float_t xmin, Float_t xmax, Bool_t verbose){
284   //
285   //  Fit histogram with gaussian function
286   //  
287   //  Prameters:
288   //       return value- chi2 - if negative ( not enough points)
289   //       his        -  input histogram
290   //       param      -  vector with parameters 
291   //       xmin, xmax -  range to fit - if xmin=xmax=0 - the full histogram range used
292   //  Fitting:
293   //  1. Step - make logarithm
294   //  2. Linear  fit (parabola) - more robust - always converge
295   //  3. In case of small statistic bins are averaged
296   //  
297   static TLinearFitter fitter(3,"pol2");
298   TVectorD  par(3);
299   TVectorD  sigma(3);
300   TMatrixD mat(3,3);
301   if (his->GetMaximum()<4) return -1;  
302   if (his->GetEntries()<12) return -1;  
303   if (his->GetRMS()<mat.GetTol()) return -1;
304   Float_t maxEstimate   = his->GetEntries()*his->GetBinWidth(1)/TMath::Sqrt((TMath::TwoPi()*his->GetRMS()));
305   Int_t dsmooth = TMath::Nint(6./TMath::Sqrt(maxEstimate));
306
307   if (maxEstimate<1) return -1;
308   Int_t nbins    = his->GetNbinsX();
309   Int_t npoints=0;
310   //
311
312
313   if (xmin>=xmax){
314     xmin = his->GetXaxis()->GetXmin();
315     xmax = his->GetXaxis()->GetXmax();
316   }
317   for (Int_t iter=0; iter<2; iter++){
318     fitter.ClearPoints();
319     npoints=0;
320     for (Int_t ibin=1;ibin<nbins+1; ibin++){
321       Int_t countB=1;
322       Float_t entriesI =  his->GetBinContent(ibin);
323       for (Int_t delta = -dsmooth; delta<=dsmooth; delta++){
324         if (ibin+delta>1 &&ibin+delta<nbins-1){
325           entriesI +=  his->GetBinContent(ibin+delta);
326           countB++;
327         }
328       }
329       entriesI/=countB;
330       Double_t xcenter= his->GetBinCenter(ibin);
331       if (xcenter<xmin || xcenter>xmax) continue;
332       Double_t error=1./TMath::Sqrt(countB);
333       Float_t   cont=2;
334       if (iter>0){
335         if (par[0]+par[1]*xcenter+par[2]*xcenter*xcenter>20) return 0;
336         cont = TMath::Exp(par[0]+par[1]*xcenter+par[2]*xcenter*xcenter);
337         if (cont>1.) error = 1./TMath::Sqrt(cont*Float_t(countB));
338       }
339       if (entriesI>1&&cont>1){
340         fitter.AddPoint(&xcenter,TMath::Log(Float_t(entriesI)),error);
341         npoints++;
342       }
343     }  
344     if (npoints>3){
345       fitter.Eval();
346       fitter.GetParameters(par);
347     }else{
348       break;
349     }
350   }
351   if (npoints<=3){
352     return -1;
353   }
354   fitter.GetParameters(par);
355   fitter.GetCovarianceMatrix(mat);
356   if (TMath::Abs(par[1])<mat.GetTol()) return -1;
357   if (TMath::Abs(par[2])<mat.GetTol()) return -1;
358   Double_t chi2 = fitter.GetChisquare()/Float_t(npoints);
359   //fitter.GetParameters();
360   if (!param)  param  = new TVectorD(3);
361   // if (!matrix) matrix = new TMatrixD(3,3); // Covariance matrix to be implemented
362   (*param)[1] = par[1]/(-2.*par[2]);
363   (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
364   (*param)[0] = TMath::Exp(par[0]+ par[1]* (*param)[1] +  par[2]*(*param)[1]*(*param)[1]);
365   if (verbose){
366     par.Print();
367     mat.Print();
368     param->Print();
369     printf("Chi2=%f\n",chi2);
370     TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",his->GetXaxis()->GetXmin(),his->GetXaxis()->GetXmax());
371     f1->SetParameter(0, (*param)[0]);
372     f1->SetParameter(1, (*param)[1]);
373     f1->SetParameter(2, (*param)[2]);    
374     f1->Draw("same");
375   }
376   return chi2;
377 }
378
379 Double_t  TStatToolkit::FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param, TMatrixD */*matrix*/, Bool_t verbose){
380   //
381   //  Fit histogram with gaussian function
382   //  
383   //  Prameters:
384   //     nbins: size of the array and number of histogram bins
385   //     xMin, xMax: histogram range
386   //     param: paramters of the fit (0-Constant, 1-Mean, 2-Sigma)
387   //     matrix: covariance matrix -- not implemented yet, pass dummy matrix!!!
388   //
389   //  Return values:
390   //    >0: the chi2 returned by TLinearFitter
391   //    -3: only three points have been used for the calculation - no fitter was used
392   //    -2: only two points have been used for the calculation - center of gravity was uesed for calculation
393   //    -1: only one point has been used for the calculation - center of gravity was uesed for calculation
394   //    -4: invalid result!!
395   //
396   //  Fitting:
397   //  1. Step - make logarithm
398   //  2. Linear  fit (parabola) - more robust - always converge
399   //  
400   static TLinearFitter fitter(3,"pol2");
401   static TMatrixD mat(3,3);
402   static Double_t kTol = mat.GetTol();
403   fitter.StoreData(kFALSE);
404   fitter.ClearPoints();
405   TVectorD  par(3);
406   TVectorD  sigma(3);
407   TMatrixD A(3,3);
408   TMatrixD b(3,1);
409   Float_t rms = TMath::RMS(nBins,arr);
410   Float_t max = TMath::MaxElement(nBins,arr);
411   Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
412
413   Float_t meanCOG = 0;
414   Float_t rms2COG = 0;
415   Float_t sumCOG  = 0;
416
417   Float_t entries = 0;
418   Int_t nfilled=0;
419
420   for (Int_t i=0; i<nBins; i++){
421       entries+=arr[i];
422       if (arr[i]>0) nfilled++;
423   }
424
425   if (max<4) return -4;
426   if (entries<12) return -4;
427   if (rms<kTol) return -4;
428
429   Int_t npoints=0;
430   //
431
432   //
433   for (Int_t ibin=0;ibin<nBins; ibin++){
434       Float_t entriesI = arr[ibin];
435     if (entriesI>1){
436       Double_t xcenter = xMin+(ibin+0.5)*binWidth;
437       
438       Float_t error    = 1./TMath::Sqrt(entriesI);
439       Float_t val = TMath::Log(Float_t(entriesI));
440       fitter.AddPoint(&xcenter,val,error);
441       if (npoints<3){
442           A(npoints,0)=1;
443           A(npoints,1)=xcenter;
444           A(npoints,2)=xcenter*xcenter;
445           b(npoints,0)=val;
446           meanCOG+=xcenter*entriesI;
447           rms2COG +=xcenter*entriesI*xcenter;
448           sumCOG +=entriesI;
449       }
450       npoints++;
451     }
452   }
453
454   
455   Double_t chi2 = 0;
456   if (npoints>=3){
457       if ( npoints == 3 ){
458           //analytic calculation of the parameters for three points
459           A.Invert();
460           TMatrixD res(1,3);
461           res.Mult(A,b);
462           par[0]=res(0,0);
463           par[1]=res(0,1);
464           par[2]=res(0,2);
465           chi2 = -3.;
466       } else {
467           // use fitter for more than three points
468           fitter.Eval();
469           fitter.GetParameters(par);
470           fitter.GetCovarianceMatrix(mat);
471           chi2 = fitter.GetChisquare()/Float_t(npoints);
472       }
473       if (TMath::Abs(par[1])<kTol) return -4;
474       if (TMath::Abs(par[2])<kTol) return -4;
475
476       if (!param)  param  = new TVectorD(3);
477       //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
478
479       (*param)[1] = par[1]/(-2.*par[2]);
480       (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
481       Double_t lnparam0 = par[0]+ par[1]* (*param)[1] +  par[2]*(*param)[1]*(*param)[1];
482       if ( lnparam0>307 ) return -4;
483       (*param)[0] = TMath::Exp(lnparam0);
484       if (verbose){
485           par.Print();
486           mat.Print();
487           param->Print();
488           printf("Chi2=%f\n",chi2);
489           TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",xMin,xMax);
490           f1->SetParameter(0, (*param)[0]);
491           f1->SetParameter(1, (*param)[1]);
492           f1->SetParameter(2, (*param)[2]);
493           f1->Draw("same");
494       }
495       return chi2;
496   }
497
498   if (npoints == 2){
499       //use center of gravity for 2 points
500       meanCOG/=sumCOG;
501       rms2COG /=sumCOG;
502       (*param)[0] = max;
503       (*param)[1] = meanCOG;
504       (*param)[2] = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
505       chi2=-2.;
506   }
507   if ( npoints == 1 ){
508       meanCOG/=sumCOG;
509       (*param)[0] = max;
510       (*param)[1] = meanCOG;
511       (*param)[2] = binWidth/TMath::Sqrt(12);
512       chi2=-1.;
513   }
514   return chi2;
515
516 }
517
518
519 Float_t TStatToolkit::GetCOG(Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms, Float_t *sum)
520 {
521     //
522     //  calculate center of gravity rms and sum for array 'arr' with nBins an a x range xMin to xMax
523     //  return COG; in case of failure return xMin
524     //
525     Float_t meanCOG = 0;
526     Float_t rms2COG = 0;
527     Float_t sumCOG  = 0;
528     Int_t npoints   = 0;
529
530     Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
531
532     for (Int_t ibin=0; ibin<nBins; ibin++){
533         Float_t entriesI = (Float_t)arr[ibin];
534         Double_t xcenter = xMin+(ibin+0.5)*binWidth;
535         if ( entriesI>0 ){
536             meanCOG += xcenter*entriesI;
537             rms2COG += xcenter*entriesI*xcenter;
538             sumCOG  += entriesI;
539             npoints++;
540         }
541     }
542     if ( sumCOG == 0 ) return xMin;
543     meanCOG/=sumCOG;
544
545     if ( rms ){
546         rms2COG /=sumCOG;
547         (*rms) = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
548         if ( npoints == 1 ) (*rms) = binWidth/TMath::Sqrt(12);
549     }
550
551     if ( sum )
552         (*sum) = sumCOG;
553
554     return meanCOG;
555 }
556
557
558
559 ///////////////////////////////////////////////////////////////
560 //////////////         TEST functions /////////////////////////
561 ///////////////////////////////////////////////////////////////
562
563
564
565
566
567 void TStatToolkit::TestGausFit(Int_t nhistos){
568   //
569   // Test performance of the parabolic - gaussian fit - compare it with 
570   // ROOT gauss fit
571   //  nhistos - number of histograms to be used for test
572   //
573   TTreeSRedirector *pcstream = new TTreeSRedirector("fitdebug.root");
574   
575   Float_t  *xTrue = new Float_t[nhistos];
576   Float_t  *sTrue = new Float_t[nhistos];
577   TVectorD **par1  = new TVectorD*[nhistos];
578   TVectorD **par2  = new TVectorD*[nhistos];
579   TMatrixD dummy(3,3);
580   
581   
582   TH1F **h1f = new TH1F*[nhistos];
583   TF1  *myg = new TF1("myg","gaus");
584   TF1  *fit = new TF1("fit","gaus");
585   gRandom->SetSeed(0);
586   
587   //init
588   for (Int_t i=0;i<nhistos; i++){
589     par1[i] = new TVectorD(3);
590     par2[i] = new TVectorD(3);
591     h1f[i]  = new TH1F(Form("h1f%d",i),Form("h1f%d",i),20,-10,10);
592     xTrue[i]= gRandom->Rndm();
593     gSystem->Sleep(2);
594     sTrue[i]= .75+gRandom->Rndm()*.5;
595     myg->SetParameters(1,xTrue[i],sTrue[i]);
596     h1f[i]->FillRandom("myg");
597   }
598   
599   TStopwatch s;
600   s.Start();
601   //standard gaus fit
602   for (Int_t i=0; i<nhistos; i++){
603     h1f[i]->Fit(fit,"0q");
604     (*par1[i])(0) = fit->GetParameter(0);
605     (*par1[i])(1) = fit->GetParameter(1);
606     (*par1[i])(2) = fit->GetParameter(2);
607   }
608   s.Stop();
609   printf("Gaussian fit\t");
610   s.Print();
611   
612   s.Start();
613   //TStatToolkit gaus fit
614   for (Int_t i=0; i<nhistos; i++){
615     TStatToolkit::FitGaus(h1f[i]->GetArray()+1,h1f[i]->GetNbinsX(),h1f[i]->GetXaxis()->GetXmin(),h1f[i]->GetXaxis()->GetXmax(),par2[i],&dummy);
616   }
617   
618   s.Stop();
619   printf("Parabolic fit\t");
620   s.Print();
621   //write stream
622   for (Int_t i=0;i<nhistos; i++){
623     Float_t xt  = xTrue[i];
624     Float_t st  = sTrue[i];
625     (*pcstream)<<"data"
626                <<"xTrue="<<xt
627                <<"sTrue="<<st
628                <<"pg.="<<(par1[i])
629                <<"pa.="<<(par2[i])
630                <<"\n";
631   }    
632   //delete pointers
633   for (Int_t i=0;i<nhistos; i++){
634     delete par1[i];
635     delete par2[i];
636     delete h1f[i];
637   }
638   delete pcstream;
639   delete []h1f;
640   delete []xTrue;
641   delete []sTrue;
642   //
643   delete []par1;
644   delete []par2;
645
646 }
647
648
649
650 TGraph2D * TStatToolkit::MakeStat2D(TH3 * his, Int_t delta0, Int_t delta1, Int_t type){
651   //
652   //
653   //
654   // delta - number of bins to integrate
655   // type - 0 - mean value
656
657   TAxis * xaxis  = his->GetXaxis();
658   TAxis * yaxis  = his->GetYaxis();
659   //  TAxis * zaxis  = his->GetZaxis();
660   Int_t   nbinx  = xaxis->GetNbins();
661   Int_t   nbiny  = yaxis->GetNbins();
662   char name[1000];
663   Int_t icount=0;
664   TGraph2D  *graph = new TGraph2D(nbinx*nbiny);
665   TF1 f1("f1","gaus");
666   for (Int_t ix=0; ix<nbinx;ix++)
667     for (Int_t iy=0; iy<nbiny;iy++){
668       Float_t xcenter = xaxis->GetBinCenter(ix); 
669       Float_t ycenter = yaxis->GetBinCenter(iy); 
670       snprintf(name,1000,"%s_%d_%d",his->GetName(), ix,iy);
671       TH1 *projection = his->ProjectionZ(name,ix-delta0,ix+delta0,iy-delta1,iy+delta1);
672       Float_t stat= 0;
673       if (type==0) stat = projection->GetMean();
674       if (type==1) stat = projection->GetRMS();
675       if (type==2 || type==3){
676         TVectorD vec(3);
677         TStatToolkit::LTM((TH1F*)projection,&vec,0.7);
678         if (type==2) stat= vec[1];
679         if (type==3) stat= vec[0];      
680       }
681       if (type==4|| type==5){
682         projection->Fit(&f1);
683         if (type==4) stat= f1.GetParameter(1);
684         if (type==5) stat= f1.GetParameter(2);
685       }
686       //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
687       graph->SetPoint(icount,xcenter, ycenter, stat);
688       icount++;
689     }
690   return graph;
691 }
692
693 TGraph * TStatToolkit::MakeStat1D(TH3 * his, Int_t delta1, Int_t type){
694   //
695   //
696   //
697   // delta - number of bins to integrate
698   // type - 0 - mean value
699
700   TAxis * xaxis  = his->GetXaxis();
701   TAxis * yaxis  = his->GetYaxis();
702   //  TAxis * zaxis  = his->GetZaxis();
703   Int_t   nbinx  = xaxis->GetNbins();
704   Int_t   nbiny  = yaxis->GetNbins();
705   char name[1000];
706   Int_t icount=0;
707   TGraph  *graph = new TGraph(nbinx);
708   TF1 f1("f1","gaus");
709   for (Int_t ix=0; ix<nbinx;ix++){
710     Float_t xcenter = xaxis->GetBinCenter(ix); 
711     //    Float_t ycenter = yaxis->GetBinCenter(iy); 
712     snprintf(name,1000,"%s_%d",his->GetName(), ix);
713     TH1 *projection = his->ProjectionZ(name,ix-delta1,ix+delta1,0,nbiny);
714     Float_t stat= 0;
715     if (type==0) stat = projection->GetMean();
716     if (type==1) stat = projection->GetRMS();
717     if (type==2 || type==3){
718       TVectorD vec(3);
719         TStatToolkit::LTM((TH1F*)projection,&vec,0.7);
720         if (type==2) stat= vec[1];
721         if (type==3) stat= vec[0];      
722     }
723     if (type==4|| type==5){
724       projection->Fit(&f1);
725       if (type==4) stat= f1.GetParameter(1);
726       if (type==5) stat= f1.GetParameter(2);
727     }
728       //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
729     graph->SetPoint(icount,xcenter, stat);
730     icount++;
731   }
732   return graph;
733 }
734
735
736
737
738
739 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){
740    //
741    // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
742    // returns chi2, fitParam and covMatrix
743    // returns TString with fitted formula
744    //
745
746    TString formulaStr(formula); 
747    TString drawStr(drawCommand);
748    TString cutStr(cuts);
749    TString ferr("1");
750
751    TString strVal(drawCommand);
752    if (strVal.Contains(":")){
753      TObjArray* valTokens = strVal.Tokenize(":");
754      drawStr = valTokens->At(0)->GetName();
755      ferr       = valTokens->At(1)->GetName();     
756    }
757
758       
759    formulaStr.ReplaceAll("++", "~");
760    TObjArray* formulaTokens = formulaStr.Tokenize("~"); 
761    Int_t dim = formulaTokens->GetEntriesFast();
762    
763    fitParam.ResizeTo(dim);
764    covMatrix.ResizeTo(dim,dim);
765    
766    TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
767    fitter->StoreData(kTRUE);   
768    fitter->ClearPoints();
769    
770    Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff",  stop-start, start);
771    if (entries == -1) return new TString("An ERROR has occured during fitting!");
772    Double_t **values = new Double_t*[dim+1] ; 
773    //
774    entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff",  stop-start, start);
775    if (entries == -1) {
776      delete []values;
777      return new TString("An ERROR has occured during fitting!");
778    }
779    Double_t *errors = new Double_t[entries];
780    memcpy(errors,  tree->GetV1(), entries*sizeof(Double_t));
781    
782    for (Int_t i = 0; i < dim + 1; i++){
783       Int_t centries = 0;
784       if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start);
785       else  centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start);
786       
787       if (entries != centries) {
788         delete []errors;
789         delete []values;
790         return new TString("An ERROR has occured during fitting!");
791       }
792       values[i] = new Double_t[entries];
793       memcpy(values[i],  tree->GetV1(), entries*sizeof(Double_t)); 
794    }
795    
796    // add points to the fitter
797    for (Int_t i = 0; i < entries; i++){
798       Double_t x[1000];
799       for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
800       fitter->AddPoint(x, values[dim][i], errors[i]);
801    }
802
803    fitter->Eval();
804    if (frac>0.5 && frac<1){
805      fitter->EvalRobust(frac);
806    }else{
807      if (fix0) {
808        fitter->FixParameter(0,0);
809        fitter->Eval();     
810      }
811    }
812    fitter->GetParameters(fitParam);
813    fitter->GetCovarianceMatrix(covMatrix);
814    chi2 = fitter->GetChisquare();
815    npoints = entries;   
816    TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula; 
817    
818    for (Int_t iparam = 0; iparam < dim; iparam++) {
819      returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
820      if (iparam < dim-1) returnFormula.Append("+");
821    }
822    returnFormula.Append(" )");
823    
824    
825    for (Int_t j=0; j<dim+1;j++) delete [] values[j];
826
827
828    delete formulaTokens;
829    delete fitter;
830    delete[] values;
831    delete[] errors;
832    return preturnFormula;
833 }
834
835 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){
836    //
837    // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
838    // returns chi2, fitParam and covMatrix
839    // returns TString with fitted formula
840    //
841
842    TString formulaStr(formula); 
843    TString drawStr(drawCommand);
844    TString cutStr(cuts);
845    TString ferr("1");
846
847    TString strVal(drawCommand);
848    if (strVal.Contains(":")){
849      TObjArray* valTokens = strVal.Tokenize(":");
850      drawStr = valTokens->At(0)->GetName();
851      ferr       = valTokens->At(1)->GetName();     
852    }
853
854       
855    formulaStr.ReplaceAll("++", "~");
856    TObjArray* formulaTokens = formulaStr.Tokenize("~"); 
857    Int_t dim = formulaTokens->GetEntriesFast();
858    
859    fitParam.ResizeTo(dim);
860    covMatrix.ResizeTo(dim,dim);
861    
862    TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
863    fitter->StoreData(kTRUE);   
864    fitter->ClearPoints();
865    
866    Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff",  stop-start, start);
867    if (entries == -1) return new TString("An ERROR has occured during fitting!");
868    Double_t **values = new Double_t*[dim+1] ; 
869    //
870    entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff",  stop-start, start);
871    if (entries == -1) {
872      delete [] values;
873      return new TString("An ERROR has occured during fitting!");
874    }
875    Double_t *errors = new Double_t[entries];
876    memcpy(errors,  tree->GetV1(), entries*sizeof(Double_t));
877    
878    for (Int_t i = 0; i < dim + 1; i++){
879       Int_t centries = 0;
880       if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start);
881       else  centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start);
882       
883       if (entries != centries) {
884         delete []errors;
885         delete []values;
886         return new TString("An ERROR has occured during fitting!");
887       }
888       values[i] = new Double_t[entries];
889       memcpy(values[i],  tree->GetV1(), entries*sizeof(Double_t)); 
890    }
891    
892    // add points to the fitter
893    for (Int_t i = 0; i < entries; i++){
894       Double_t x[1000];
895       for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
896       fitter->AddPoint(x, values[dim][i], errors[i]);
897    }
898    if (constrain>0){
899      for (Int_t i = 0; i < dim; i++){
900        Double_t x[1000];
901        for (Int_t j=0; j<dim;j++) if (i!=j) x[j]=0;
902        x[i]=1.;
903        fitter->AddPoint(x, 0, constrain);
904      }
905    }
906
907
908    fitter->Eval();
909    if (frac>0.5 && frac<1){
910      fitter->EvalRobust(frac);   
911    }
912    fitter->GetParameters(fitParam);
913    fitter->GetCovarianceMatrix(covMatrix);
914    chi2 = fitter->GetChisquare();
915    npoints = entries;
916    
917    TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula; 
918    
919    for (Int_t iparam = 0; iparam < dim; iparam++) {
920      returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
921      if (iparam < dim-1) returnFormula.Append("+");
922    }
923    returnFormula.Append(" )");
924    
925    for (Int_t j=0; j<dim+1;j++) delete [] values[j];
926    
927
928
929    delete formulaTokens;
930    delete fitter;
931    delete[] values;
932    delete[] errors;
933    return preturnFormula;
934 }
935
936
937
938 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){
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    }
956
957       
958    formulaStr.ReplaceAll("++", "~");
959    TObjArray* formulaTokens = formulaStr.Tokenize("~"); 
960    Int_t dim = formulaTokens->GetEntriesFast();
961    
962    fitParam.ResizeTo(dim);
963    covMatrix.ResizeTo(dim,dim);
964    TString fitString="x0";
965    for (Int_t i=1; i<dim; i++) fitString+=Form("++x%d",i);     
966    TLinearFitter* fitter = new TLinearFitter(dim, fitString.Data());
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) return new TString("An ERROR has occured during fitting!");
972    Double_t **values = new Double_t*[dim+1] ; 
973    //
974    entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff",  stop-start, start);
975    if (entries == -1) {
976      delete []values;
977      return new TString("An ERROR has occured during fitting!");
978    }
979    Double_t *errors = new Double_t[entries];
980    memcpy(errors,  tree->GetV1(), entries*sizeof(Double_t));
981    
982    for (Int_t i = 0; i < dim + 1; i++){
983       Int_t centries = 0;
984       if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start);
985       else  centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start);
986       
987       if (entries != centries) {
988         delete []errors;
989         delete []values;
990         return new TString("An ERROR has occured during fitting!");
991       }
992       values[i] = new Double_t[entries];
993       memcpy(values[i],  tree->GetV1(), entries*sizeof(Double_t)); 
994    }
995    
996    // add points to the fitter
997    for (Int_t i = 0; i < entries; i++){
998       Double_t x[1000];
999       for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
1000       fitter->AddPoint(x, values[dim][i], errors[i]);
1001    }
1002
1003    fitter->Eval();
1004    if (frac>0.5 && frac<1){
1005      fitter->EvalRobust(frac);
1006    }
1007    fitter->GetParameters(fitParam);
1008    fitter->GetCovarianceMatrix(covMatrix);
1009    chi2 = fitter->GetChisquare();
1010    npoints = entries;
1011    
1012    TString *preturnFormula = new TString("("), &returnFormula = *preturnFormula; 
1013    
1014    for (Int_t iparam = 0; iparam < dim; iparam++) {
1015      returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam]));
1016      if (iparam < dim-1) returnFormula.Append("+");
1017    }
1018    returnFormula.Append(" )");
1019    
1020    
1021    for (Int_t j=0; j<dim+1;j++) delete [] values[j];
1022    
1023    delete formulaTokens;
1024    delete fitter;
1025    delete[] values;
1026    delete[] errors;
1027    return preturnFormula;
1028 }
1029
1030
1031
1032
1033
1034 Int_t TStatToolkit::GetFitIndex(TString fString, TString subString){
1035   //
1036   // fitString - ++ separated list of fits
1037   // substring - ++ separated list of the requiered substrings
1038   //
1039   // return the last occurance of substring in fit string
1040   // 
1041   TObjArray *arrFit = fString.Tokenize("++");
1042   TObjArray *arrSub = subString.Tokenize("++");
1043   Int_t index=-1;
1044   for (Int_t i=0; i<arrFit->GetEntries(); i++){
1045     Bool_t isOK=kTRUE;
1046     TString str =arrFit->At(i)->GetName();
1047     for (Int_t isub=0; isub<arrSub->GetEntries(); isub++){
1048       if (str.Contains(arrSub->At(isub)->GetName())==0) isOK=kFALSE;
1049     }
1050     if (isOK) index=i;
1051   }
1052   return index;
1053 }
1054
1055
1056 TString  TStatToolkit::FilterFit(TString &input, TString filter, TVectorD &param, TMatrixD & covar){
1057   //
1058   // Filter fit expression make sub-fit
1059   //
1060   TObjArray *array0= input.Tokenize("++");
1061   TObjArray *array1= filter.Tokenize("++");
1062   //TString *presult=new TString("(0");
1063   TString result="(0.0";
1064   for (Int_t i=0; i<array0->GetEntries(); i++){
1065     Bool_t isOK=kTRUE;
1066     TString str(array0->At(i)->GetName());
1067     for (Int_t j=0; j<array1->GetEntries(); j++){
1068       if (str.Contains(array1->At(j)->GetName())==0) isOK=kFALSE;      
1069     }
1070     if (isOK) {
1071       result+="+"+str;
1072       result+=Form("*(%f)",param[i+1]);
1073       printf("%f\t%f\t%s\n",param[i+1], TMath::Sqrt(covar(i+1,i+1)),str.Data());    
1074     }
1075   }
1076   result+="-0.)";
1077   return result;
1078 }
1079
1080 void TStatToolkit::Update1D(Double_t delta, Double_t sigma, Int_t s1, TMatrixD &vecXk, TMatrixD &covXk){
1081   //
1082   // Update parameters and covariance - with one measurement
1083   // Input:
1084   // vecXk - input vector - Updated in function 
1085   // covXk - covariance matrix - Updated in function
1086   // delta, sigma, s1 - new measurement, rms of new measurement and the index of measurement
1087   const Int_t knMeas=1;
1088   Int_t knElem=vecXk.GetNrows();
1089  
1090   TMatrixD mat1(knElem,knElem);            // update covariance matrix
1091   TMatrixD matHk(1,knElem);        // vector to mesurement
1092   TMatrixD vecYk(knMeas,1);        // Innovation or measurement residual
1093   TMatrixD matHkT(knElem,knMeas);  // helper matrix Hk transpose
1094   TMatrixD matSk(knMeas,knMeas);   // Innovation (or residual) covariance
1095   TMatrixD matKk(knElem,knMeas);   // Optimal Kalman gain
1096   TMatrixD covXk2(knElem,knElem);  // helper matrix
1097   TMatrixD covXk3(knElem,knElem);  // helper matrix
1098   TMatrixD vecZk(1,1);
1099   TMatrixD measR(1,1);
1100   vecZk(0,0)=delta;
1101   measR(0,0)=sigma*sigma;
1102   //
1103   // reset matHk
1104   for (Int_t iel=0;iel<knElem;iel++) 
1105     for (Int_t ip=0;ip<knMeas;ip++) matHk(ip,iel)=0; 
1106   //mat1
1107   for (Int_t iel=0;iel<knElem;iel++) {
1108     for (Int_t jel=0;jel<knElem;jel++) mat1(iel,jel)=0;
1109     mat1(iel,iel)=1;
1110   }
1111   //
1112   matHk(0, s1)=1;
1113   vecYk = vecZk-matHk*vecXk;               // Innovation or measurement residual
1114   matHkT=matHk.T(); matHk.T();
1115   matSk = (matHk*(covXk*matHkT))+measR;    // Innovation (or residual) covariance
1116   matSk.Invert();
1117   matKk = (covXk*matHkT)*matSk;            //  Optimal Kalman gain
1118   vecXk += matKk*vecYk;                    //  updated vector 
1119   covXk2= (mat1-(matKk*matHk));
1120   covXk3 =  covXk2*covXk;          
1121   covXk = covXk3;  
1122   Int_t nrows=covXk3.GetNrows();
1123   
1124   for (Int_t irow=0; irow<nrows; irow++)
1125     for (Int_t icol=0; icol<nrows; icol++){
1126       // rounding problems - make matrix again symteric
1127       covXk(irow,icol)=(covXk3(irow,icol)+covXk3(icol,irow))*0.5; 
1128     }
1129 }
1130
1131
1132
1133 void   TStatToolkit::Constrain1D(TString &input, TString filter, TVectorD &param, TMatrixD & covar, Double_t mean, Double_t sigma){
1134   //
1135   // constrain linear fit
1136   // input  - string description of fit function
1137   // filter - string filter to select sub fits
1138   // param,covar - parameters and covariance matrix of the fit
1139   // mean,sigma  - new measurement uning which the fit is updated
1140   //
1141   TObjArray *array0= input.Tokenize("++");
1142   TObjArray *array1= filter.Tokenize("++");
1143   TMatrixD paramM(param.GetNrows(),1);
1144   for (Int_t i=0; i<=array0->GetEntries(); i++){paramM(i,0)=param(i);}
1145   
1146   for (Int_t i=0; i<array0->GetEntries(); i++){
1147     Bool_t isOK=kTRUE;
1148     TString str(array0->At(i)->GetName());
1149     for (Int_t j=0; j<array1->GetEntries(); j++){
1150       if (str.Contains(array1->At(j)->GetName())==0) isOK=kFALSE;      
1151     }
1152     if (isOK) {
1153       TStatToolkit::Update1D(mean, sigma, i+1, paramM, covar);//
1154     }
1155   }
1156   for (Int_t i=0; i<=array0->GetEntries(); i++){
1157     param(i)=paramM(i,0);
1158   }
1159 }
1160
1161 TString  TStatToolkit::MakeFitString(TString &input, TVectorD &param, TMatrixD & covar){
1162   //
1163   //
1164   //
1165   TObjArray *array0= input.Tokenize("++");
1166   TString result="(0.0";
1167   for (Int_t i=0; i<array0->GetEntries(); i++){
1168     TString str(array0->At(i)->GetName());
1169     result+="+"+str;
1170     result+=Form("*(%f)",param[i+1]);
1171     printf("%f\t%f\t%s\n", param[i+1], TMath::Sqrt(covar(i+1,i+1)),str.Data());    
1172   }
1173   result+="-0.)";
1174   return result;
1175 }
1176
1177
1178 TGraph * TStatToolkit::MakeGraphSparse(TTree * tree, const char * expr, const char * cut){
1179   //
1180   // Make a sparse draw of the variables
1181   //
1182   const Int_t entries =  tree->Draw(expr,cut,"goff");
1183   //  TGraph * graph = (TGraph*)gPad->GetPrimitive("Graph"); // 2D
1184   TGraph * graph = new TGraph (entries, tree->GetV2(),tree->GetV1());
1185   //
1186   Int_t *index = new Int_t[entries];
1187   TMath::Sort(entries,graph->GetX(),index,kFALSE);
1188   
1189   Double_t *VV = new Double_t[entries];
1190
1191   Double_t count = 0.5;
1192   vector<Int_t> vrun;
1193   VV[index[0]] = count;
1194   vrun.push_back(graph->GetX()[index[0]]);
1195   for(Int_t i=1;i<entries;i++){
1196     if(graph->GetX()[index[i]]==graph->GetX()[index[i-1]])
1197       VV[index[i]] = count; 
1198     else if(graph->GetX()[index[i]]!=graph->GetX()[index[i-1]]){
1199       count++;
1200       VV[index[i]] = count;
1201       vrun.push_back(graph->GetX()[index[i]]);
1202     }
1203   }
1204   
1205   const Int_t newNbins = int(count+0.5);
1206   Double_t *newBins = new Double_t[newNbins+1];
1207   for(Int_t i=0; i<=count+1;i++){
1208     newBins[i] = i;
1209   }
1210   
1211   TGraph *graphNew = new TGraph(entries,VV,graph->GetY());
1212   graphNew->GetXaxis()->Set(newNbins,newBins);
1213   
1214   Char_t xName[50];
1215   Double_t bin_unit = graphNew->GetXaxis()->GetNbins()/count;
1216   for(Int_t i=0;i<count;i++){
1217     snprintf(xName,50,"%d",vrun.at(i));
1218     graphNew->GetXaxis()->SetBinLabel(i+1,xName);
1219   }
1220   graphNew->GetHistogram()->SetTitle("");
1221   
1222   delete [] VV;
1223   delete [] index;
1224   delete [] newBins;
1225   return graphNew;
1226 }
1227