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