]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/STEERBase/AliMathBase.cxx
Fix in the last caall to CleanOwnPrimaryVertex
[u/mrichter/AliRoot.git] / STEER / STEERBase / AliMathBase.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 AliMathBase
19 // 
20 // Subset of  matheamtical functions  not included in the TMath
21 //
22
23 ///////////////////////////////////////////////////////////////////////////
24 #include "TMath.h"
25 #include "AliMathBase.h"
26 #include "Riostream.h"
27 #include "TH1F.h"
28 #include "TH3.h"
29 #include "TF1.h"
30 #include "TLinearFitter.h"
31
32 #include "AliExternalTrackParam.h"
33
34 //
35 // includes neccessary for test functions
36 //
37
38 #include "TSystem.h"
39 #include "TRandom.h"
40 #include "TStopwatch.h"
41 #include "TTreeStream.h"
42  
43 ClassImp(AliMathBase) // Class implementation to enable ROOT I/O
44  
45 AliMathBase::AliMathBase() : TObject()
46 {
47   //
48   // Default constructor
49   //
50 }
51 ///////////////////////////////////////////////////////////////////////////
52 AliMathBase::~AliMathBase()
53 {
54   //
55   // Destructor
56   //
57 }
58
59
60 //_____________________________________________________________________________
61 void AliMathBase::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean
62                            , Double_t &sigma, Int_t hh)
63 {
64   //
65   // Robust estimator in 1D case MI version - (faster than ROOT version)
66   //
67   // For the univariate case
68   // estimates of location and scatter are returned in mean and sigma parameters
69   // the algorithm works on the same principle as in multivariate case -
70   // it finds a subset of size hh with smallest sigma, and then returns mean and
71   // sigma of this subset
72   //
73
74   if (hh==0)
75     hh=(nvectors+2)/2;
76   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};
77   Int_t *index=new Int_t[nvectors];
78   TMath::Sort(nvectors, data, index, kFALSE);
79   
80   Int_t    nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
81   Double_t factor = faclts[TMath::Max(0,nquant-1)];
82   
83   Double_t sumx  =0;
84   Double_t sumx2 =0;
85   Int_t    bestindex = -1;
86   Double_t bestmean  = 0; 
87   Double_t bestsigma = (data[index[nvectors-1]]-data[index[0]]+1.);   // maximal possible sigma
88   bestsigma *=bestsigma;
89
90   for (Int_t i=0; i<hh; i++){
91     sumx  += data[index[i]];
92     sumx2 += data[index[i]]*data[index[i]];
93   }
94   
95   Double_t norm = 1./Double_t(hh);
96   Double_t norm2 = 1./Double_t(hh-1);
97   for (Int_t i=hh; i<nvectors; i++){
98     Double_t cmean  = sumx*norm;
99     Double_t csigma = (sumx2 - hh*cmean*cmean)*norm2;
100     if (csigma<bestsigma){
101       bestmean  = cmean;
102       bestsigma = csigma;
103       bestindex = i-hh;
104     }
105     
106     sumx  += data[index[i]]-data[index[i-hh]];
107     sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
108   }
109   
110   Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
111   mean  = bestmean;
112   sigma = bstd;
113   delete [] index;
114
115 }
116
117
118
119 void AliMathBase::EvaluateUniExternal(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh,  Float_t externalfactor)
120 {
121   // Modified version of ROOT robust EvaluateUni
122   // robust estimator in 1D case MI version
123   // added external factor to include precision of external measurement
124   // 
125
126   if (hh==0)
127     hh=(nvectors+2)/2;
128   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};
129   Int_t *index=new Int_t[nvectors];
130   TMath::Sort(nvectors, data, index, kFALSE);
131   //
132   Int_t    nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
133   Double_t factor = faclts[0];
134   if (nquant>0){
135     // fix proper normalization - Anja
136     factor = faclts[nquant-1];
137   }
138
139   //
140   //
141   Double_t sumx  =0;
142   Double_t sumx2 =0;
143   Int_t    bestindex = -1;
144   Double_t bestmean  = 0; 
145   Double_t bestsigma = -1;
146   for (Int_t i=0; i<hh; i++){
147     sumx  += data[index[i]];
148     sumx2 += data[index[i]]*data[index[i]];
149   }
150   //   
151   Double_t kfactor = 2.*externalfactor - externalfactor*externalfactor;
152   Double_t norm = 1./Double_t(hh);
153   for (Int_t i=hh; i<nvectors; i++){
154     Double_t cmean  = sumx*norm;
155     Double_t csigma = (sumx2*norm - cmean*cmean*kfactor);
156     if (csigma<bestsigma ||  bestsigma<0){
157       bestmean  = cmean;
158       bestsigma = csigma;
159       bestindex = i-hh;
160     }
161     //
162     //
163     sumx  += data[index[i]]-data[index[i-hh]];
164     sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]];
165   }
166   
167   Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma));
168   mean  = bestmean;
169   sigma = bstd;
170   delete [] index;
171 }
172
173
174 //_____________________________________________________________________________
175 Int_t AliMathBase::Freq(Int_t n, const Int_t *inlist
176                         , Int_t *outlist, Bool_t down)
177 {    
178   //
179   //  Sort eleements according occurancy 
180   //  The size of output array has is 2*n 
181   //
182
183   Int_t * sindexS = new Int_t[n];     // temp array for sorting
184   Int_t * sindexF = new Int_t[2*n];   
185   for (Int_t i=0;i<n;i++) sindexF[i]=0;
186   //
187   TMath::Sort(n,inlist, sindexS, down);  
188   Int_t last      = inlist[sindexS[0]];
189   Int_t val       = last;
190   sindexF[0]      = 1;
191   sindexF[0+n]    = last;
192   Int_t countPos  = 0;
193   //
194   //  find frequency
195   for(Int_t i=1;i<n; i++){
196     val = inlist[sindexS[i]];
197     if (last == val)   sindexF[countPos]++;
198     else{      
199       countPos++;
200       sindexF[countPos+n] = val;
201       sindexF[countPos]++;
202       last =val;
203     }
204   }
205   if (last==val) countPos++;
206   // sort according frequency
207   TMath::Sort(countPos, sindexF, sindexS, kTRUE);
208   for (Int_t i=0;i<countPos;i++){
209     outlist[2*i  ] = sindexF[sindexS[i]+n];
210     outlist[2*i+1] = sindexF[sindexS[i]];
211   }
212   delete [] sindexS;
213   delete [] sindexF;
214   
215   return countPos;
216
217 }
218
219 //___AliMathBase__________________________________________________________________________
220 void AliMathBase::TruncatedMean(TH1F * his, TVectorD *param, Float_t down, Float_t up, Bool_t verbose){
221   //
222   //
223   //
224   Int_t nbins    = his->GetNbinsX();
225   Float_t nentries = his->GetEntries();
226   Float_t sum      =0;
227   Float_t mean   = 0;
228   Float_t sigma2 = 0;
229   Float_t ncumul=0;  
230   for (Int_t ibin=1;ibin<nbins; ibin++){
231     ncumul+= his->GetBinContent(ibin);
232     Float_t fraction = Float_t(ncumul)/Float_t(nentries);
233     if (fraction>down && fraction<up){
234       sum+=his->GetBinContent(ibin);
235       mean+=his->GetBinCenter(ibin)*his->GetBinContent(ibin);
236       sigma2+=his->GetBinCenter(ibin)*his->GetBinCenter(ibin)*his->GetBinContent(ibin);      
237     }
238   }
239   mean/=sum;
240   sigma2= TMath::Sqrt(TMath::Abs(sigma2/sum-mean*mean));
241   if (param){
242     (*param)[0] = his->GetMaximum();
243     (*param)[1] = mean;
244     (*param)[2] = sigma2;
245     
246   }
247   if (verbose)  printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma2);
248 }
249
250 void AliMathBase::LTM(TH1F * his, TVectorD *param , Float_t fraction,  Bool_t verbose){
251   //
252   // LTM
253   //
254   Int_t nbins    = his->GetNbinsX();
255   Int_t nentries = (Int_t)his->GetEntries();
256   Double_t *data  = new Double_t[nentries];
257   Int_t npoints=0;
258   for (Int_t ibin=1;ibin<nbins; ibin++){
259     Float_t entriesI = his->GetBinContent(ibin);
260     Float_t xcenter= his->GetBinCenter(ibin);
261     for (Int_t ic=0; ic<entriesI; ic++){
262       if (npoints<nentries){
263         data[npoints]= xcenter;
264         npoints++;
265       }
266     }
267   }
268   Double_t mean, sigma;
269   Int_t npoints2=TMath::Min(Int_t(fraction*Float_t(npoints)),npoints-1);
270   npoints2=TMath::Max(Int_t(0.5*Float_t(npoints)),npoints2);
271   AliMathBase::EvaluateUni(npoints, data, mean,sigma,npoints2);
272   delete [] data;
273   if (verbose)  printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma);if (param){
274     (*param)[0] = his->GetMaximum();
275     (*param)[1] = mean;
276     (*param)[2] = sigma;    
277   }
278 }
279
280 Double_t  AliMathBase::FitGaus(TH1F* his, TVectorD *param, TMatrixD */*matrix*/, Float_t xmin, Float_t xmax, Bool_t verbose){
281   //
282   //  Fit histogram with gaussian function
283   //  
284   //  Prameters:
285   //       return value- chi2 - if negative ( not enough points)
286   //       his        -  input histogram
287   //       param      -  vector with parameters 
288   //       xmin, xmax -  range to fit - if xmin=xmax=0 - the full histogram range used
289   //  Fitting:
290   //  1. Step - make logarithm
291   //  2. Linear  fit (parabola) - more robust - always converge
292   //  3. In case of small statistic bins are averaged
293   //  
294   static TLinearFitter fitter(3,"pol2");
295   TVectorD  par(3);
296   TVectorD  sigma(3);
297   TMatrixD mat(3,3);
298   if (his->GetMaximum()<4) return -1;  
299   if (his->GetEntries()<12) return -1;  
300   if (his->GetRMS()<mat.GetTol()) return -1;
301   Float_t maxEstimate   = his->GetEntries()*his->GetBinWidth(1)/TMath::Sqrt((TMath::TwoPi()*his->GetRMS()));
302   Int_t dsmooth = TMath::Nint(6./TMath::Sqrt(maxEstimate));
303
304   if (maxEstimate<1) return -1;
305   Int_t nbins    = his->GetNbinsX();
306   Int_t npoints=0;
307   //
308
309
310   if (xmin>=xmax){
311     xmin = his->GetXaxis()->GetXmin();
312     xmax = his->GetXaxis()->GetXmax();
313   }
314   for (Int_t iter=0; iter<2; iter++){
315     fitter.ClearPoints();
316     npoints=0;
317     for (Int_t ibin=1;ibin<nbins+1; ibin++){
318       Int_t countB=1;
319       Float_t entriesI =  his->GetBinContent(ibin);
320       for (Int_t delta = -dsmooth; delta<=dsmooth; delta++){
321         if (ibin+delta>1 &&ibin+delta<nbins-1){
322           entriesI +=  his->GetBinContent(ibin+delta);
323           countB++;
324         }
325       }
326       entriesI/=countB;
327       Double_t xcenter= his->GetBinCenter(ibin);
328       if (xcenter<xmin || xcenter>xmax) continue;
329       Double_t error=1./TMath::Sqrt(countB);
330       Float_t   cont=2;
331       if (iter>0){
332         if (par[0]+par[1]*xcenter+par[2]*xcenter*xcenter>20) return 0;
333         cont = TMath::Exp(par[0]+par[1]*xcenter+par[2]*xcenter*xcenter);
334         if (cont>1.) error = 1./TMath::Sqrt(cont*Float_t(countB));
335       }
336       if (entriesI>1&&cont>1){
337         fitter.AddPoint(&xcenter,TMath::Log(Float_t(entriesI)),error);
338         npoints++;
339       }
340     }  
341     if (npoints>3){
342       fitter.Eval();
343       fitter.GetParameters(par);
344     }else{
345       break;
346     }
347   }
348   if (npoints<=3){
349     return -1;
350   }
351   fitter.GetParameters(par);
352   fitter.GetCovarianceMatrix(mat);
353   if (TMath::Abs(par[1])<mat.GetTol()) return -1;
354   if (TMath::Abs(par[2])<mat.GetTol()) return -1;
355   Double_t chi2 = fitter.GetChisquare()/Float_t(npoints);
356   //fitter.GetParameters();
357   if (!param)  param  = new TVectorD(3);
358   //if (!matrix) matrix = new TMatrixD(3,3);
359   (*param)[1] = par[1]/(-2.*par[2]);
360   (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
361   (*param)[0] = TMath::Exp(par[0]+ par[1]* (*param)[1] +  par[2]*(*param)[1]*(*param)[1]);
362   if (verbose){
363     par.Print();
364     mat.Print();
365     param->Print();
366     printf("Chi2=%f\n",chi2);
367     TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",his->GetXaxis()->GetXmin(),his->GetXaxis()->GetXmax());
368     f1->SetParameter(0, (*param)[0]);
369     f1->SetParameter(1, (*param)[1]);
370     f1->SetParameter(2, (*param)[2]);    
371     f1->Draw("same");
372   }
373   return chi2;
374 }
375
376 Double_t  AliMathBase::FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param, TMatrixD */*matrix*/, Bool_t verbose){
377   //
378   //  Fit histogram with gaussian function
379   //  
380   //  Prameters:
381   //     nbins: size of the array and number of histogram bins
382   //     xMin, xMax: histogram range
383   //     param: paramters of the fit (0-Constant, 1-Mean, 2-Sigma, 3-Sum)
384   //     matrix: covariance matrix -- not implemented yet, pass dummy matrix!!!
385   //
386   //  Return values:
387   //    >0: the chi2 returned by TLinearFitter
388   //    -3: only three points have been used for the calculation - no fitter was used
389   //    -2: only two points have been used for the calculation - center of gravity was uesed for calculation
390   //    -1: only one point has been used for the calculation - center of gravity was uesed for calculation
391   //    -4: invalid result!!
392   //
393   //  Fitting:
394   //  1. Step - make logarithm
395   //  2. Linear  fit (parabola) - more robust - always converge
396   //  
397   static TLinearFitter fitter(3,"pol2");
398   static TMatrixD mat(3,3);
399   static Double_t kTol = mat.GetTol();
400   fitter.StoreData(kFALSE);
401   fitter.ClearPoints();
402   TVectorD  par(3);
403   TVectorD  sigma(3);
404   TMatrixD A(3,3);
405   TMatrixD b(3,1);
406   Float_t rms = TMath::RMS(nBins,arr);
407   Float_t max = TMath::MaxElement(nBins,arr);
408   Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
409
410   Float_t meanCOG = 0;
411   Float_t rms2COG = 0;
412   Float_t sumCOG  = 0;
413
414   Float_t entries = 0;
415   Int_t nfilled=0;
416
417   if (!param)  param  = new TVectorD(4);
418
419   for (Int_t i=0; i<nBins; i++){
420       entries+=arr[i];
421       if (arr[i]>0) nfilled++;
422   }
423   (*param)[0] = 0;
424   (*param)[1] = 0;
425   (*param)[2] = 0;
426   (*param)[3] = 0;
427
428   if (max<4) return -4;
429   if (entries<12) return -4;
430   if (rms<kTol) return -4;
431
432   (*param)[3] = entries;
433
434   Int_t npoints=0;
435   for (Int_t ibin=0;ibin<nBins; ibin++){
436       Float_t entriesI = arr[ibin];
437     if (entriesI>1){
438       Double_t xcenter = xMin+(ibin+0.5)*binWidth;
439       Float_t error    = 1./TMath::Sqrt(entriesI);
440       Float_t val = TMath::Log(Float_t(entriesI));
441       fitter.AddPoint(&xcenter,val,error);
442       if (npoints<3){
443           A(npoints,0)=1;
444           A(npoints,1)=xcenter;
445           A(npoints,2)=xcenter*xcenter;
446           b(npoints,0)=val;
447           meanCOG+=xcenter*entriesI;
448           rms2COG +=xcenter*entriesI*xcenter;
449           sumCOG +=entriesI;
450       }
451       npoints++;
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(4);
477       if ( param->GetNrows()<4 ) param->ResizeTo(4);
478       //if (!matrix) matrix = new TMatrixD(3,3);  // !!!!might be a memory leek. use dummy matrix pointer to call this function!
479
480       (*param)[1] = par[1]/(-2.*par[2]);
481       (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2]));
482       Double_t lnparam0 = par[0]+ par[1]* (*param)[1] +  par[2]*(*param)[1]*(*param)[1];
483       if ( lnparam0>307 ) return -4;
484       (*param)[0] = TMath::Exp(lnparam0);
485       if (verbose){
486           par.Print();
487           mat.Print();
488           param->Print();
489           printf("Chi2=%f\n",chi2);
490           TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",xMin,xMax);
491           f1->SetParameter(0, (*param)[0]);
492           f1->SetParameter(1, (*param)[1]);
493           f1->SetParameter(2, (*param)[2]);
494           f1->Draw("same");
495       }
496       return chi2;
497   }
498
499   if (npoints == 2){
500       //use center of gravity for 2 points
501       meanCOG/=sumCOG;
502       rms2COG /=sumCOG;
503       (*param)[0] = max;
504       (*param)[1] = meanCOG;
505       (*param)[2] = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
506       chi2=-2.;
507   }
508   if ( npoints == 1 ){
509       meanCOG/=sumCOG;
510       (*param)[0] = max;
511       (*param)[1] = meanCOG;
512       (*param)[2] = binWidth/TMath::Sqrt(12);
513       chi2=-1.;
514   }
515   return chi2;
516
517 }
518
519
520 Float_t AliMathBase::GetCOG(Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms, Float_t *sum)
521 {
522     //
523     //  calculate center of gravity rms and sum for array 'arr' with nBins an a x range xMin to xMax
524     //  return COG; in case of failure return xMin
525     //
526     Float_t meanCOG = 0;
527     Float_t rms2COG = 0;
528     Float_t sumCOG  = 0;
529     Int_t npoints   = 0;
530
531     Float_t binWidth = (xMax-xMin)/(Float_t)nBins;
532
533     for (Int_t ibin=0; ibin<nBins; ibin++){
534         Float_t entriesI = (Float_t)arr[ibin];
535         Double_t xcenter = xMin+(ibin+0.5)*binWidth;
536         if ( entriesI>0 ){
537             meanCOG += xcenter*entriesI;
538             rms2COG += xcenter*entriesI*xcenter;
539             sumCOG  += entriesI;
540             npoints++;
541         }
542     }
543     if ( sumCOG == 0 ) return xMin;
544     meanCOG/=sumCOG;
545
546     if ( rms ){
547         rms2COG /=sumCOG;
548         (*rms) = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG));
549         if ( npoints == 1 ) (*rms) = binWidth/TMath::Sqrt(12);
550     }
551
552     if ( sum )
553         (*sum) = sumCOG;
554
555     return meanCOG;
556 }
557
558
559 Double_t AliMathBase::ErfcFast(Double_t x){
560   // Fast implementation of the complementary error function
561   // The error of the approximation is |eps(x)| < 5E-4
562   // See Abramowitz and Stegun, p.299, 7.1.27
563
564   Double_t z = TMath::Abs(x);
565   Double_t ans = 1+z*(0.278393+z*(0.230389+z*(0.000972+z*0.078108)));
566   ans = 1.0/ans;
567   ans *= ans;
568   ans *= ans;
569
570   return (x>=0.0 ? ans : 2.0 - ans);
571 }
572
573 ///////////////////////////////////////////////////////////////
574 //////////////         TEST functions /////////////////////////
575 ///////////////////////////////////////////////////////////////
576
577
578
579
580
581 void AliMathBase::TestGausFit(Int_t nhistos){
582   //
583   // Test performance of the parabolic - gaussian fit - compare it with 
584   // ROOT gauss fit
585   //  nhistos - number of histograms to be used for test
586   //
587   TTreeSRedirector *pcstream = new TTreeSRedirector("fitdebug.root");
588   
589   Float_t  *xTrue = new Float_t[nhistos];
590   Float_t  *sTrue = new Float_t[nhistos];
591   TVectorD **par1  = new TVectorD*[nhistos];
592   TVectorD **par2  = new TVectorD*[nhistos];
593   TMatrixD dummy(3,3);
594   
595   
596   TH1F **h1f = new TH1F*[nhistos];
597   TF1  *myg = new TF1("myg","gaus");
598   TF1  *fit = new TF1("fit","gaus");
599   gRandom->SetSeed(0);
600   
601   //init
602   for (Int_t i=0;i<nhistos; i++){
603     par1[i] = new TVectorD(3);
604     par2[i] = new TVectorD(3);
605     h1f[i]  = new TH1F(Form("h1f%d",i),Form("h1f%d",i),20,-10,10);
606     xTrue[i]= gRandom->Rndm();
607     gSystem->Sleep(2);
608     sTrue[i]= .75+gRandom->Rndm()*.5;
609     myg->SetParameters(1,xTrue[i],sTrue[i]);
610     h1f[i]->FillRandom("myg");
611   }
612   
613   TStopwatch s;
614   s.Start();
615   //standard gaus fit
616   for (Int_t i=0; i<nhistos; i++){
617     h1f[i]->Fit(fit,"0q");
618     (*par1[i])(0) = fit->GetParameter(0);
619     (*par1[i])(1) = fit->GetParameter(1);
620     (*par1[i])(2) = fit->GetParameter(2);
621   }
622   s.Stop();
623   printf("Gaussian fit\t");
624   s.Print();
625   
626   s.Start();
627   //AliMathBase gaus fit
628   for (Int_t i=0; i<nhistos; i++){
629     AliMathBase::FitGaus(h1f[i]->GetArray()+1,h1f[i]->GetNbinsX(),h1f[i]->GetXaxis()->GetXmin(),h1f[i]->GetXaxis()->GetXmax(),par2[i],&dummy);
630   }
631   
632   s.Stop();
633   printf("Parabolic fit\t");
634   s.Print();
635   //write stream
636   for (Int_t i=0;i<nhistos; i++){
637     Float_t xt  = xTrue[i];
638     Float_t st  = sTrue[i];
639     (*pcstream)<<"data"
640                <<"xTrue="<<xt
641                <<"sTrue="<<st
642                <<"pg.="<<(par1[i])
643                <<"pa.="<<(par2[i])
644                <<"\n";
645   }    
646   //delete pointers
647   for (Int_t i=0;i<nhistos; i++){
648     delete par1[i];
649     delete par2[i];
650     delete h1f[i];
651   }
652   delete pcstream;
653   delete []h1f;
654   delete []xTrue;
655   delete []sTrue;
656   //
657   delete []par1;
658   delete []par2;
659
660 }
661
662
663
664 TGraph2D * AliMathBase::MakeStat2D(TH3 * his, Int_t delta0, Int_t delta1, Int_t type){
665   //
666   //
667   //
668   // delta - number of bins to integrate
669   // type - 0 - mean value
670
671   TAxis * xaxis  = his->GetXaxis();
672   TAxis * yaxis  = his->GetYaxis();
673   //  TAxis * zaxis  = his->GetZaxis();
674   Int_t   nbinx  = xaxis->GetNbins();
675   Int_t   nbiny  = yaxis->GetNbins();
676   const Int_t nc=1000;
677   char name[nc];
678   Int_t icount=0;
679   TGraph2D  *graph = new TGraph2D(nbinx*nbiny);
680   TF1 f1("f1","gaus");
681   for (Int_t ix=0; ix<nbinx;ix++)
682     for (Int_t iy=0; iy<nbiny;iy++){
683       Float_t xcenter = xaxis->GetBinCenter(ix); 
684       Float_t ycenter = yaxis->GetBinCenter(iy); 
685       snprintf(name,nc,"%s_%d_%d",his->GetName(), ix,iy);
686       TH1 *projection = his->ProjectionZ(name,ix-delta0,ix+delta0,iy-delta1,iy+delta1);
687       Float_t stat= 0;
688       if (type==0) stat = projection->GetMean();
689       if (type==1) stat = projection->GetRMS();
690       if (type==2 || type==3){
691         TVectorD vec(3);
692         AliMathBase::LTM((TH1F*)projection,&vec,0.7);
693         if (type==2) stat= vec[1];
694         if (type==3) stat= vec[0];      
695       }
696       if (type==4|| type==5){
697         projection->Fit(&f1);
698         if (type==4) stat= f1.GetParameter(1);
699         if (type==5) stat= f1.GetParameter(2);
700       }
701       //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
702       graph->SetPoint(icount,xcenter, ycenter, stat);
703       icount++;
704     }
705   return graph;
706 }
707
708 TGraph * AliMathBase::MakeStat1D(TH3 * his, Int_t delta1, Int_t type){
709   //
710   //
711   //
712   // delta - number of bins to integrate
713   // type - 0 - mean value
714
715   TAxis * xaxis  = his->GetXaxis();
716   TAxis * yaxis  = his->GetYaxis();
717   //  TAxis * zaxis  = his->GetZaxis();
718   Int_t   nbinx  = xaxis->GetNbins();
719   Int_t   nbiny  = yaxis->GetNbins();
720   const Int_t nc=1000;
721   char name[nc];
722   Int_t icount=0;
723   TGraph  *graph = new TGraph(nbinx);
724   TF1 f1("f1","gaus");
725   for (Int_t ix=0; ix<nbinx;ix++){
726     Float_t xcenter = xaxis->GetBinCenter(ix); 
727     //    Float_t ycenter = yaxis->GetBinCenter(iy); 
728     snprintf(name,nc,"%s_%d",his->GetName(), ix);
729     TH1 *projection = his->ProjectionZ(name,ix-delta1,ix+delta1,0,nbiny);
730     Float_t stat= 0;
731     if (type==0) stat = projection->GetMean();
732     if (type==1) stat = projection->GetRMS();
733     if (type==2 || type==3){
734       TVectorD vec(3);
735         AliMathBase::LTM((TH1F*)projection,&vec,0.7);
736         if (type==2) stat= vec[1];
737         if (type==3) stat= vec[0];      
738     }
739     if (type==4|| type==5){
740       projection->Fit(&f1);
741       if (type==4) stat= f1.GetParameter(1);
742       if (type==5) stat= f1.GetParameter(2);
743     }
744       //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat);
745     graph->SetPoint(icount,xcenter, stat);
746     icount++;
747   }
748   return graph;
749 }
750
751 Double_t AliMathBase::TruncatedGaus(Double_t mean, Double_t sigma, Double_t cutat)
752 {
753   // return number generated according to a gaussian distribution N(mean,sigma) truncated at cutat
754   //
755   Double_t value;
756   do{
757     value=gRandom->Gaus(mean,sigma);
758   }while(TMath::Abs(value-mean)>cutat);
759   return value;
760 }
761
762 Double_t AliMathBase::TruncatedGaus(Double_t mean, Double_t sigma, Double_t leftCut, Double_t rightCut)
763 {
764   // return number generated according to a gaussian distribution N(mean,sigma)
765   // truncated at leftCut and rightCut
766   //
767   Double_t value;
768   do{
769     value=gRandom->Gaus(mean,sigma);
770   }while((value-mean)<-leftCut || (value-mean)>rightCut);
771   return value;
772 }
773
774 Double_t AliMathBase::BetheBlochAleph(Double_t bg,
775          Double_t kp1,
776          Double_t kp2,
777          Double_t kp3,
778          Double_t kp4,
779          Double_t kp5) {
780   //
781   // This is the empirical ALEPH parameterization of the Bethe-Bloch formula.
782   // It is normalized to 1 at the minimum.
783   //
784   // bg - beta*gamma
785   // 
786   // The default values for the kp* parameters are for ALICE TPC.
787   // The returned value is in MIP units
788   //
789
790   return AliExternalTrackParam::BetheBlochAleph(bg,kp1,kp2,kp3,kp4,kp5);
791 }