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