]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STAT/TStatToolkit.cxx
Fixing ompilation warnings
[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++) sindexF[i]=0;
188   //
189   TMath::Sort(n,inlist, sindexS, down);  
190   Int_t last      = inlist[sindexS[0]];
191   Int_t val       = last;
192   sindexF[0]      = 1;
193   sindexF[0+n]    = last;
194   Int_t countPos  = 0;
195   //
196   //  find frequency
197   for(Int_t i=1;i<n; i++){
198     val = inlist[sindexS[i]];
199     if (last == val)   sindexF[countPos]++;
200     else{      
201       countPos++;
202       sindexF[countPos+n] = val;
203       sindexF[countPos]++;
204       last =val;
205     }
206   }
207   if (last==val) countPos++;
208   // sort according frequency
209   TMath::Sort(countPos, sindexF, sindexS, kTRUE);
210   for (Int_t i=0;i<countPos;i++){
211     outlist[2*i  ] = sindexF[sindexS[i]+n];
212     outlist[2*i+1] = sindexF[sindexS[i]];
213   }
214   delete [] sindexS;
215   delete [] sindexF;
216   
217   return countPos;
218
219 }
220
221 //___TStatToolkit__________________________________________________________________________
222 void TStatToolkit::TruncatedMean(TH1F * his, TVectorD *param, Float_t down, Float_t up, Bool_t verbose){
223   //
224   //
225   //
226   Int_t nbins    = his->GetNbinsX();
227   Float_t nentries = his->GetEntries();
228   Float_t sum      =0;
229   Float_t mean   = 0;
230   Float_t sigma2 = 0;
231   Float_t ncumul=0;  
232   for (Int_t ibin=1;ibin<nbins; ibin++){
233     ncumul+= his->GetBinContent(ibin);
234     Float_t fraction = Float_t(ncumul)/Float_t(nentries);
235     if (fraction>down && fraction<up){
236       sum+=his->GetBinContent(ibin);
237       mean+=his->GetBinCenter(ibin)*his->GetBinContent(ibin);
238       sigma2+=his->GetBinCenter(ibin)*his->GetBinCenter(ibin)*his->GetBinContent(ibin);      
239     }
240   }
241   mean/=sum;
242   sigma2= TMath::Sqrt(TMath::Abs(sigma2/sum-mean*mean));
243   if (param){
244     (*param)[0] = his->GetMaximum();
245     (*param)[1] = mean;
246     (*param)[2] = sigma2;
247     
248   }
249   if (verbose)  printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma2);
250 }
251
252 void TStatToolkit::LTM(TH1F * his, TVectorD *param , Float_t fraction,  Bool_t verbose){
253   //
254   // LTM
255   //
256   Int_t nbins    = his->GetNbinsX();
257   Int_t nentries = (Int_t)his->GetEntries();
258   Double_t *data  = new Double_t[nentries];
259   Int_t npoints=0;
260   for (Int_t ibin=1;ibin<nbins; ibin++){
261     Float_t entriesI = his->GetBinContent(ibin);
262     Float_t xcenter= his->GetBinCenter(ibin);
263     for (Int_t ic=0; ic<entriesI; ic++){
264       if (npoints<nentries){
265         data[npoints]= xcenter;
266         npoints++;
267       }
268     }
269   }
270   Double_t mean, sigma;
271   Int_t npoints2=TMath::Min(Int_t(fraction*Float_t(npoints)),npoints-1);
272   npoints2=TMath::Max(Int_t(0.5*Float_t(npoints)),npoints2);
273   TStatToolkit::EvaluateUni(npoints, data, mean,sigma,npoints2);
274   delete [] data;
275   if (verbose)  printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma);if (param){
276     (*param)[0] = his->GetMaximum();
277     (*param)[1] = mean;
278     (*param)[2] = sigma;    
279   }
280 }
281
282 Double_t  TStatToolkit::FitGaus(TH1F* his, TVectorD *param, TMatrixD *matrix, Float_t xmin, Float_t xmax, Bool_t verbose){
283   //
284   //  Fit histogram with gaussian function
285   //  
286   //  Prameters:
287   //       return value- chi2 - if negative ( not enough points)
288   //       his        -  input histogram
289   //       param      -  vector with parameters 
290   //       xmin, xmax -  range to fit - if xmin=xmax=0 - the full histogram range used
291   //  Fitting:
292   //  1. Step - make logarithm
293   //  2. Linear  fit (parabola) - more robust - always converge
294   //  3. In case of small statistic bins are averaged
295   //  
296   static TLinearFitter fitter(3,"pol2");
297   TVectorD  par(3);
298   TVectorD  sigma(3);
299   TMatrixD mat(3,3);
300   if (his->GetMaximum()<4) return -1;  
301   if (his->GetEntries()<12) return -1;  
302   if (his->GetRMS()<mat.GetTol()) return -1;
303   Float_t maxEstimate   = his->GetEntries()*his->GetBinWidth(1)/TMath::Sqrt((TMath::TwoPi()*his->GetRMS()));
304   Int_t dsmooth = TMath::Nint(6./TMath::Sqrt(maxEstimate));
305
306   if (maxEstimate<1) return -1;
307   Int_t nbins    = his->GetNbinsX();
308   Int_t npoints=0;
309   //
310
311
312   if (xmin>=xmax){
313     xmin = his->GetXaxis()->GetXmin();
314     xmax = his->GetXaxis()->GetXmax();
315   }
316   for (Int_t iter=0; iter<2; iter++){
317     fitter.ClearPoints();
318     npoints=0;
319     for (Int_t ibin=1;ibin<nbins+1; ibin++){
320       Int_t countB=1;
321       Float_t entriesI =  his->GetBinContent(ibin);
322       for (Int_t delta = -dsmooth; delta<=dsmooth; delta++){
323         if (ibin+delta>1 &&ibin+delta<nbins-1){
324           entriesI +=  his->GetBinContent(ibin+delta);
325           countB++;
326         }
327       }
328       entriesI/=countB;
329       Double_t xcenter= his->GetBinCenter(ibin);
330       if (xcenter<xmin || xcenter>xmax) continue;
331       Double_t error=1./TMath::Sqrt(countB);
332       Float_t   cont=2;
333       if (iter>0){
334         if (par[0]+par[1]*xcenter+par[2]*xcenter*xcenter>20) return 0;
335         cont = TMath::Exp(par[0]+par[1]*xcenter+par[2]*xcenter*xcenter);
336         if (cont>1.) error = 1./TMath::Sqrt(cont*Float_t(countB));
337       }
338       if (entriesI>1&&cont>1){
339         fitter.AddPoint(&xcenter,TMath::Log(Float_t(entriesI)),error);
340         npoints++;
341       }
342     }  
343     if (npoints>3){
344       fitter.Eval();
345       fitter.GetParameters(par);
346     }else{
347       break;
348     }
349   }
350   if (npoints<=3){
351     return -1;
352   }
353   fitter.GetParameters(par);
354   fitter.GetCovarianceMatrix(mat);
355   if (TMath::Abs(par[1])<mat.GetTol()) return -1;
356   if (TMath::Abs(par[2])<mat.GetTol()) return -1;
357   Double_t chi2 = fitter.GetChisquare()/Float_t(npoints);
358   //fitter.GetParameters();
359   if (!param)  param  = new TVectorD(3);
360   if (!matrix) matrix = new TMatrixD(3,3);
361   (*param)[1] = par[1]/(-2.*par[2]);
362   (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
363   (*param)[0] = TMath::Exp(par[0]+ par[1]* (*param)[1] +  par[2]*(*param)[1]*(*param)[1]);
364   if (verbose){
365     par.Print();
366     mat.Print();
367     param->Print();
368     printf("Chi2=%f\n",chi2);
369     TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",his->GetXaxis()->GetXmin(),his->GetXaxis()->GetXmax());
370     f1->SetParameter(0, (*param)[0]);
371     f1->SetParameter(1, (*param)[1]);
372     f1->SetParameter(2, (*param)[2]);    
373     f1->Draw("same");
374   }
375   return chi2;
376 }
377
378 Double_t  TStatToolkit::FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param, TMatrixD *matrix, Bool_t verbose){
379   //
380   //  Fit histogram with gaussian function
381   //  
382   //  Prameters:
383   //     nbins: size of the array and number of histogram bins
384   //     xMin, xMax: histogram range
385   //     param: paramters of the fit (0-Constant, 1-Mean, 2-Sigma)
386   //     matrix: covariance matrix -- not implemented yet, pass dummy matrix!!!
387   //
388   //  Return values:
389   //    >0: the chi2 returned by TLinearFitter
390   //    -3: only three points have been used for the calculation - no fitter was used
391   //    -2: only two points have been used for the calculation - center of gravity was uesed for calculation
392   //    -1: only one point has been used for the calculation - center of gravity was uesed for calculation
393   //    -4: invalid result!!
394   //
395   //  Fitting:
396   //  1. Step - make logarithm
397   //  2. Linear  fit (parabola) - more robust - always converge
398   //  
399   static TLinearFitter fitter(3,"pol2");
400   static TMatrixD mat(3,3);
401   static Double_t kTol = mat.GetTol();
402   fitter.StoreData(kFALSE);
403   fitter.ClearPoints();
404   TVectorD  par(3);
405   TVectorD  sigma(3);
406   TMatrixD A(3,3);
407   TMatrixD b(3,1);
408   Float_t rms = TMath::RMS(nBins,arr);
409   Float_t max = TMath::MaxElement(nBins,arr);
410   Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
411
412   Float_t meanCOG = 0;
413   Float_t rms2COG = 0;
414   Float_t sumCOG  = 0;
415
416   Float_t entries = 0;
417   Int_t nfilled=0;
418
419   for (Int_t i=0; i<nBins; i++){
420       entries+=arr[i];
421       if (arr[i]>0) nfilled++;
422   }
423
424   if (max<4) return -4;
425   if (entries<12) return -4;
426   if (rms<kTol) return -4;
427
428   Int_t npoints=0;
429   //
430
431   //
432   for (Int_t ibin=0;ibin<nBins; ibin++){
433       Float_t entriesI = arr[ibin];
434     if (entriesI>1){
435       Double_t xcenter = xMin+(ibin+0.5)*binWidth;
436       
437       Float_t error    = 1./TMath::Sqrt(entriesI);
438       Float_t val = TMath::Log(Float_t(entriesI));
439       fitter.AddPoint(&xcenter,val,error);
440       if (npoints<3){
441           A(npoints,0)=1;
442           A(npoints,1)=xcenter;
443           A(npoints,2)=xcenter*xcenter;
444           b(npoints,0)=val;
445           meanCOG+=xcenter*entriesI;
446           rms2COG +=xcenter*entriesI*xcenter;
447           sumCOG +=entriesI;
448       }
449       npoints++;
450     }
451   }
452
453   
454   Double_t chi2 = 0;
455   if (npoints>=3){
456       if ( npoints == 3 ){
457           //analytic calculation of the parameters for three points
458           A.Invert();
459           TMatrixD res(1,3);
460           res.Mult(A,b);
461           par[0]=res(0,0);
462           par[1]=res(0,1);
463           par[2]=res(0,2);
464           chi2 = -3.;
465       } else {
466           // use fitter for more than three points
467           fitter.Eval();
468           fitter.GetParameters(par);
469           fitter.GetCovarianceMatrix(mat);
470           chi2 = fitter.GetChisquare()/Float_t(npoints);
471       }
472       if (TMath::Abs(par[1])<kTol) return -4;
473       if (TMath::Abs(par[2])<kTol) return -4;
474
475       if (!param)  param  = new TVectorD(3);
476       if (!matrix) matrix = new TMatrixD(3,3);  // !!!!might be a memory leek. use dummy matrix pointer to call this function!
477
478       (*param)[1] = par[1]/(-2.*par[2]);
479       (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
480       Double_t lnparam0 = par[0]+ par[1]* (*param)[1] +  par[2]*(*param)[1]*(*param)[1];
481       if ( lnparam0>307 ) return -4;
482       (*param)[0] = TMath::Exp(lnparam0);
483       if (verbose){
484           par.Print();
485           mat.Print();
486           param->Print();
487           printf("Chi2=%f\n",chi2);
488           TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",xMin,xMax);
489           f1->SetParameter(0, (*param)[0]);
490           f1->SetParameter(1, (*param)[1]);
491           f1->SetParameter(2, (*param)[2]);
492           f1->Draw("same");
493       }
494       return chi2;
495   }
496
497   if (npoints == 2){
498       //use center of gravity for 2 points
499       meanCOG/=sumCOG;
500       rms2COG /=sumCOG;
501       (*param)[0] = max;
502       (*param)[1] = meanCOG;
503       (*param)[2] = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
504       chi2=-2.;
505   }
506   if ( npoints == 1 ){
507       meanCOG/=sumCOG;
508       (*param)[0] = max;
509       (*param)[1] = meanCOG;
510       (*param)[2] = binWidth/TMath::Sqrt(12);
511       chi2=-1.;
512   }
513   return chi2;
514
515 }
516
517
518 Float_t TStatToolkit::GetCOG(Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms, Float_t *sum)
519 {
520     //
521     //  calculate center of gravity rms and sum for array 'arr' with nBins an a x range xMin to xMax
522     //  return COG; in case of failure return xMin
523     //
524     Float_t meanCOG = 0;
525     Float_t rms2COG = 0;
526     Float_t sumCOG  = 0;
527     Int_t npoints   = 0;
528
529     Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
530
531     for (Int_t ibin=0; ibin<nBins; ibin++){
532         Float_t entriesI = (Float_t)arr[ibin];
533         Double_t xcenter = xMin+(ibin+0.5)*binWidth;
534         if ( entriesI>0 ){
535             meanCOG += xcenter*entriesI;
536             rms2COG += xcenter*entriesI*xcenter;
537             sumCOG  += entriesI;
538             npoints++;
539         }
540     }
541     if ( sumCOG == 0 ) return xMin;
542     meanCOG/=sumCOG;
543
544     if ( rms ){
545         rms2COG /=sumCOG;
546         (*rms) = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
547         if ( npoints == 1 ) (*rms) = binWidth/TMath::Sqrt(12);
548     }
549
550     if ( sum )
551         (*sum) = sumCOG;
552
553     return meanCOG;
554 }
555
556
557
558 ///////////////////////////////////////////////////////////////
559 //////////////         TEST functions /////////////////////////
560 ///////////////////////////////////////////////////////////////
561
562
563
564
565
566 void TStatToolkit::TestGausFit(Int_t nhistos){
567   //
568   // Test performance of the parabolic - gaussian fit - compare it with 
569   // ROOT gauss fit
570   //  nhistos - number of histograms to be used for test
571   //
572   TTreeSRedirector *pcstream = new TTreeSRedirector("fitdebug.root");
573   
574   Float_t  *xTrue = new Float_t[nhistos];
575   Float_t  *sTrue = new Float_t[nhistos];
576   TVectorD **par1  = new TVectorD*[nhistos];
577   TVectorD **par2  = new TVectorD*[nhistos];
578   TMatrixD dummy(3,3);
579   
580   
581   TH1F **h1f = new TH1F*[nhistos];
582   TF1  *myg = new TF1("myg","gaus");
583   TF1  *fit = new TF1("fit","gaus");
584   gRandom->SetSeed(0);
585   
586   //init
587   for (Int_t i=0;i<nhistos; i++){
588     par1[i] = new TVectorD(3);
589     par2[i] = new TVectorD(3);
590     h1f[i]  = new TH1F(Form("h1f%d",i),Form("h1f%d",i),20,-10,10);
591     xTrue[i]= gRandom->Rndm();
592     gSystem->Sleep(2);
593     sTrue[i]= .75+gRandom->Rndm()*.5;
594     myg->SetParameters(1,xTrue[i],sTrue[i]);
595     h1f[i]->FillRandom("myg");
596   }
597   
598   TStopwatch s;
599   s.Start();
600   //standard gaus fit
601   for (Int_t i=0; i<nhistos; i++){
602     h1f[i]->Fit(fit,"0q");
603     (*par1[i])(0) = fit->GetParameter(0);
604     (*par1[i])(1) = fit->GetParameter(1);
605     (*par1[i])(2) = fit->GetParameter(2);
606   }
607   s.Stop();
608   printf("Gaussian fit\t");
609   s.Print();
610   
611   s.Start();
612   //TStatToolkit gaus fit
613   for (Int_t i=0; i<nhistos; i++){
614     TStatToolkit::FitGaus(h1f[i]->GetArray()+1,h1f[i]->GetNbinsX(),h1f[i]->GetXaxis()->GetXmin(),h1f[i]->GetXaxis()->GetXmax(),par2[i],&dummy);
615   }
616   
617   s.Stop();
618   printf("Parabolic fit\t");
619   s.Print();
620   //write stream
621   for (Int_t i=0;i<nhistos; i++){
622     Float_t xt  = xTrue[i];
623     Float_t st  = sTrue[i];
624     (*pcstream)<<"data"
625                <<"xTrue="<<xt
626                <<"sTrue="<<st
627                <<"pg.="<<(par1[i])
628                <<"pa.="<<(par2[i])
629                <<"\n";
630   }    
631   //delete pointers
632   for (Int_t i=0;i<nhistos; i++){
633     delete par1[i];
634     delete par2[i];
635     delete h1f[i];
636   }
637   delete pcstream;
638   delete []h1f;
639   delete []xTrue;
640   delete []sTrue;
641   //
642   delete []par1;
643   delete []par2;
644
645 }
646
647
648
649 TGraph2D * TStatToolkit::MakeStat2D(TH3 * his, Int_t delta0, Int_t delta1, Int_t type){
650   //
651   //
652   //
653   // delta - number of bins to integrate
654   // type - 0 - mean value
655
656   TAxis * xaxis  = his->GetXaxis();
657   TAxis * yaxis  = his->GetYaxis();
658   //  TAxis * zaxis  = his->GetZaxis();
659   Int_t   nbinx  = xaxis->GetNbins();
660   Int_t   nbiny  = yaxis->GetNbins();
661   char name[1000];
662   Int_t icount=0;
663   TGraph2D  *graph = new TGraph2D(nbinx*nbiny);
664   TF1 f1("f1","gaus");
665   for (Int_t ix=0; ix<nbinx;ix++)
666     for (Int_t iy=0; iy<nbiny;iy++){
667       Float_t xcenter = xaxis->GetBinCenter(ix); 
668       Float_t ycenter = yaxis->GetBinCenter(iy); 
669       sprintf(name,"%s_%d_%d",his->GetName(), ix,iy);
670       TH1 *projection = his->ProjectionZ(name,ix-delta0,ix+delta0,iy-delta1,iy+delta1);
671       Float_t stat= 0;
672       if (type==0) stat = projection->GetMean();
673       if (type==1) stat = projection->GetRMS();
674       if (type==2 || type==3){
675         TVectorD vec(3);
676         TStatToolkit::LTM((TH1F*)projection,&vec,0.7);
677         if (type==2) stat= vec[1];
678         if (type==3) stat= vec[0];      
679       }
680       if (type==4|| type==5){
681         projection->Fit(&f1);
682         if (type==4) stat= f1.GetParameter(1);
683         if (type==5) stat= f1.GetParameter(2);
684       }
685       //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
686       graph->SetPoint(icount,xcenter, ycenter, stat);
687       icount++;
688     }
689   return graph;
690 }
691
692 TGraph * TStatToolkit::MakeStat1D(TH3 * his, Int_t delta1, Int_t type){
693   //
694   //
695   //
696   // delta - number of bins to integrate
697   // type - 0 - mean value
698
699   TAxis * xaxis  = his->GetXaxis();
700   TAxis * yaxis  = his->GetYaxis();
701   //  TAxis * zaxis  = his->GetZaxis();
702   Int_t   nbinx  = xaxis->GetNbins();
703   Int_t   nbiny  = yaxis->GetNbins();
704   char name[1000];
705   Int_t icount=0;
706   TGraph  *graph = new TGraph(nbinx);
707   TF1 f1("f1","gaus");
708   for (Int_t ix=0; ix<nbinx;ix++){
709     Float_t xcenter = xaxis->GetBinCenter(ix); 
710     //    Float_t ycenter = yaxis->GetBinCenter(iy); 
711     sprintf(name,"%s_%d",his->GetName(), ix);
712     TH1 *projection = his->ProjectionZ(name,ix-delta1,ix+delta1,0,nbiny);
713     Float_t stat= 0;
714     if (type==0) stat = projection->GetMean();
715     if (type==1) stat = projection->GetRMS();
716     if (type==2 || type==3){
717       TVectorD vec(3);
718         TStatToolkit::LTM((TH1F*)projection,&vec,0.7);
719         if (type==2) stat= vec[1];
720         if (type==3) stat= vec[0];      
721     }
722     if (type==4|| type==5){
723       projection->Fit(&f1);
724       if (type==4) stat= f1.GetParameter(1);
725       if (type==5) stat= f1.GetParameter(2);
726     }
727       //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
728     graph->SetPoint(icount,xcenter, stat);
729     icount++;
730   }
731   return graph;
732 }
733
734
735
736
737
738 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){
739    //
740    // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
741    // returns chi2, fitParam and covMatrix
742    // returns TString with fitted formula
743    //
744
745    TString formulaStr(formula); 
746    TString drawStr(drawCommand);
747    TString cutStr(cuts);
748    TString ferr("1");
749
750    TString strVal(drawCommand);
751    if (strVal.Contains(":")){
752      TObjArray* valTokens = strVal.Tokenize(":");
753      drawStr = valTokens->At(0)->GetName();
754      ferr       = valTokens->At(1)->GetName();     
755    }
756
757       
758    formulaStr.ReplaceAll("++", "~");
759    TObjArray* formulaTokens = formulaStr.Tokenize("~"); 
760    Int_t dim = formulaTokens->GetEntriesFast();
761    
762    fitParam.ResizeTo(dim);
763    covMatrix.ResizeTo(dim,dim);
764    
765    TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
766    fitter->StoreData(kTRUE);   
767    fitter->ClearPoints();
768    
769    Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff",  stop-start, start);
770    if (entries == -1) return new TString("An ERROR has occured during fitting!");
771    Double_t **values = new Double_t*[dim+1] ; 
772    //
773    entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff",  stop-start, start);
774    if (entries == -1) return new TString("An ERROR has occured during fitting!");
775    Double_t *errors = new Double_t[entries];
776    memcpy(errors,  tree->GetV1(), entries*sizeof(Double_t));
777    
778    for (Int_t i = 0; i < dim + 1; i++){
779       Int_t centries = 0;
780       if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start);
781       else  centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start);
782       
783       if (entries != centries) return new TString("An ERROR has occured during fitting!");
784       values[i] = new Double_t[entries];
785       memcpy(values[i],  tree->GetV1(), entries*sizeof(Double_t)); 
786    }
787    
788    // add points to the fitter
789    for (Int_t i = 0; i < entries; i++){
790       Double_t x[1000];
791       for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
792       fitter->AddPoint(x, values[dim][i], errors[i]);
793    }
794
795    fitter->Eval();
796    if (frac>0.5 && frac<1){
797      fitter->EvalRobust(frac);
798    }
799    fitter->GetParameters(fitParam);
800    fitter->GetCovarianceMatrix(covMatrix);
801    chi2 = fitter->GetChisquare();
802    chi2 = chi2;
803    npoints = entries;
804 //    TString *preturnFormula = new TString(Form("%f*(",fitParam[0])), &returnFormula = *preturnFormula; 
805    
806 //    for (Int_t iparam = 0; iparam < dim; iparam++) {
807 //      returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]/fitParam[0]));
808 //      if (iparam < dim-1) returnFormula.Append("+");
809 //    }
810 //    returnFormula.Append(" )");
811    
812    TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula; 
813    
814    for (Int_t iparam = 0; iparam < dim; iparam++) {
815      returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
816      if (iparam < dim-1) returnFormula.Append("+");
817    }
818    returnFormula.Append(" )");
819    
820    
821
822
823    delete formulaTokens;
824    delete fitter;
825    delete[] values;
826    return preturnFormula;
827 }