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