]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/EBYE/LRC/AliLRCAnalysis.cxx
Coverity Fix:
[u/mrichter/AliRoot.git] / PWGCF / EBYE / LRC / AliLRCAnalysis.cxx
1 /**************************************************************************
2  * Author: Andrey Ivanov.                                           *
3  * Contributors are mentioned in the code where appropriate.              *
4  *                                                                        *
5  * Permission to use, copy, modify and distribute this software and its   *
6  * documentation strictly for non-commercial purposes is hereby granted   *
7  * without fee, provided that the above copyright notice appears in all   *
8  * copies and that both the copyright notice and this permission notice   *
9  * appear in the supporting documentation. The authors make no claims     *
10  * about the suitability of this software for any purpose. It is          *
11  * provided "as is" without express or implied warranty.                  *
12  **************************************************************************/
13 //-------------------------------------------------------------------------
14 //    Description: 
15 //    This class include into LRC library for Long-Range Correlation analysis
16 //    it is base class for NN, PtN, PtPt
17 //    implements base methods for thees classes
18 //    Origin: Petr Naumenko, SPbSU-CERN, Petr.Naoumenko@cern.ch ,
19 //    Andrey Ivanov (SPbSU-CERN), Igor Altsebeev (SPbSU-CERN) 
20 //-------------------------------------------------------------------------
21
22 #include "AliLRCAnalysis.h"
23 #include "TFile.h"
24 #include "AliLRCFit.h"
25 #include "TProfile.h"
26 #include "TH1D.h"
27 #include "TH2D.h"
28 #include "TPaveText.h"
29 #include "TF1.h"
30 #include "math.h"
31 #include "TStyle.h"
32 class gStyle;
33
34
35 ClassImp(AliLRCAnalysis) 
36
37 /******************************************************
38  * AliLRCAnalysis class
39  ******************************************************/
40
41 AliLRCAnalysis::AliLRCAnalysis():TObject(), fPrAbs(new TH1D()), fPrRel(new TH1D()), fPrf(new TH1D()), fPrb(new TH1D()), fileHist(new TFile()), fdptb(.0), fEntries(0), fSx((char*)" "), fSy((char*)" "), fxFitMin(.0), fxFitMax(.0), fNsigma(2.), farel(.0), fbrel(.0), farelError(.0), fbrelError(.0), fXi2rel(.0), faabs(.0), fbabs(.0), faabsError(.0), fbabsError(.0), fXi2abs(.0),fFitMethod(0){
42 //Empty constructor
43 }
44
45 AliLRCAnalysis::AliLRCAnalysis(const AliLRCAnalysis& a):TObject(a),fPrAbs(a.fPrAbs), fPrRel(a.fPrRel), fPrf(a.fPrf), fPrb(a.fPrb), fileHist(a.fileHist), fdptb(a.fdptb), fEntries(a.fEntries), fSx(a.fSx), fSy(a.fSy), fxFitMin(a.fxFitMin), fxFitMax(a.fxFitMax), fNsigma(a.fNsigma), farel(a.farel), fbrel(a.fbrel), farelError(a.farelError), fbrelError(a.fbrelError), fXi2rel(a.fXi2rel), faabs(a.faabs), fbabs(a.fbabs),faabsError(a.faabsError), fbabsError(a.fbabsError), fXi2abs(a.fXi2abs),fFitMethod(a.fFitMethod){
46 //Constructor
47 }
48
49 AliLRCAnalysis& AliLRCAnalysis::operator= (const AliLRCAnalysis& a){
50 //Operator =
51         if(this!=&a){
52                 TObject::operator= (a);
53                 fPrAbs = a.fPrAbs;
54                 fPrRel = a.fPrRel;
55                 fPrf = a.fPrf;
56                 fPrb = a.fPrb;
57                 fileHist = a.fileHist;
58                 fSx = a.fSx;
59                 fSy = a.fSy;
60                 fdptb = a.fdptb;
61                 fEntries = a.fEntries;
62                 fxFitMin = a.fxFitMin;
63                 fxFitMax = a.fxFitMax;
64                 fNsigma = a.fNsigma;
65                 farel = a.farel;
66                 fbrel = a.fbrel;
67                 farelError = a.farelError;
68                 fbrelError = a.fbrelError;
69                 fXi2rel = a.fXi2rel;
70                 faabs = a.faabs;
71                 fbabs = a.fbabs; 
72                 faabsError = a.faabsError;
73                 fbabsError = a.fbabsError; 
74                 fXi2abs = a.fXi2abs;
75                 fFitMethod = a.fFitMethod;
76         }
77         return *this;
78 }
79
80 AliLRCAnalysis::~AliLRCAnalysis() {
81 //Destructor
82     delete fPrAbs;
83     delete fPrRel;
84     delete fileHist;
85         delete fPrf;
86         delete fPrb;
87 }
88
89
90
91 void AliLRCAnalysis::CreateHist(char *name, char *nameAbs, char *nameRel, char *atitleF, char *atitleB, char *rtitleF, char *rtitleB, TH2D* sourceHist) {
92 //Create absolute and relation var histogramm
93     TProfile* profX = (TProfile*) sourceHist->ProfileX(name, 1, sourceHist->GetNbinsY());
94     fEntries = (int) sourceHist->GetEntries();
95     fPrAbs = new TH1D(nameAbs, profX->GetTitle(), profX->GetXaxis()->GetNbins(), profX->GetXaxis()->GetXmin(), profX->GetXaxis()->GetXmax());
96     fPrAbs->SetOption("E");
97     fPrRel = new TH1D(nameRel, profX->GetTitle(), profX->GetXaxis()->GetNbins(), profX->GetXaxis()->GetXmin()/sourceHist->ProjectionX()->GetMean(), profX->GetXaxis()->GetXmax()/sourceHist->ProjectionX()->GetMean());
98     fPrRel->SetOption("E");
99     fPrAbs->GetXaxis()->SetTitle(atitleF);
100     fPrAbs->GetYaxis()->SetTitle(atitleB);
101     fPrRel->GetXaxis()->SetTitle(rtitleF);
102     fPrRel->GetYaxis()->SetTitle(rtitleB);
103     fPrf = (TH1D*) sourceHist->ProjectionX();
104     fPrb = (TH1D*) sourceHist->ProjectionY();
105     fPrf->GetXaxis()->SetTitle(atitleF);
106     fPrf->GetYaxis()->SetTitle("Tracks");
107     fPrb->GetXaxis()->SetTitle(atitleB);
108     fPrb->GetYaxis()->SetTitle("Tracks");
109     fSx = atitleF;
110     fSy = atitleB;
111     double mnf = fPrf->GetMean();
112     double rms = fPrf->GetRMS();
113     fxFitMin = mnf-fNsigma*rms;//sqrt(rms);
114     fxFitMax = mnf+fNsigma*rms;//sqrt(rms);//mnf+2*sqrt(rms);
115     //delete profX;
116
117 }
118
119 void AliLRCAnalysis::SetBinsRange(int binMin, int binMax){
120 //Set the bin range
121         TH1D* h=fPrf;
122     int n=h->GetNbinsX();
123     fxFitMin = h->GetXaxis()->GetXmin();
124     fxFitMax = h->GetXaxis()->GetXmax();
125     double df = (fxFitMax-fxFitMin)/n;
126     for(int i=1; i<=n; i++)
127         if(h->GetBinContent(i) != 0){
128                 fxFitMin +=  (i + binMin) * df;
129                 break;
130         }
131     for(int i=1; i<=n; i++)
132         if(h->GetBinContent(n-i) != 0){
133                 fxFitMax -= (i + binMax) * df;
134                 break;
135         }
136 }
137
138 bool AliLRCAnalysis::SetFitRange(double xMin, double xMax){
139 //Set the fit range
140         if(xMax < xMin){
141                 return false;
142         }
143         this->fxFitMin = xMin;
144         this->fxFitMax = xMax;
145         return true;
146 }
147 bool AliLRCAnalysis::SetFitRangeMin(double xMin){
148 //Set the fit range min
149         fxFitMin = xMin;
150         return true;
151 }
152
153 bool AliLRCAnalysis::SetFitRangeMax(double xMax){
154 //Set the fit range min
155         //this->fxFitMax = xMax;
156
157         TH1D* h=fPrf;
158         double maxBorder = 0.;
159
160         for ( int binI = h->GetNbinsX(); binI > 0; binI-- )
161         {
162                 if ( h->GetBinContent(binI) != 0 )
163                 {
164                         maxBorder = h->GetBinLowEdge(binI+1) ;
165                         break;
166                 }       
167         }
168         fxFitMax = TMath::Max( fxFitMin, TMath::Min( xMax, maxBorder ) );
169
170         return true;
171 }
172
173 void AliLRCAnalysis::SetFullFitRange(){
174 //Set fitting on full range
175         TH1D* h=fPrf;
176     fxFitMin = h->GetXaxis()->GetXmin();
177         
178         for ( int binI = h->GetNbinsX(); binI > 0; binI-- )
179         {
180                 if ( h->GetBinContent(binI) != 0 )
181                 {
182                         fxFitMax = h->GetBinLowEdge(binI+1) ;
183                         break;
184                 }       
185         }
186 }
187
188
189 void AliLRCAnalysis::SetXmin(double xMin){
190         fxFitMin = xMin;
191 }
192
193 void AliLRCAnalysis::SetXmax(double xMax){
194         fxFitMax = xMax;
195 }
196 void AliLRCAnalysis::SetNsigma(double nSigma){
197 //Sets fiting range in forward value RMSs
198         fNsigma = nSigma;
199         double mnf = fPrf->GetMean();
200         double rms = fPrf->GetRMS();
201         fxFitMin = mnf-fNsigma*rms;
202         fxFitMax = mnf+fNsigma*rms;
203
204 //      cout << "### Fit params: N sigma is " << fNsigma 
205 //              << ", mean = " << mnf << ", rms = " << rms << ", minFitX = " << fxFitMin << ", maxFitX = " << fxFitMax << endl;
206 }
207
208 void AliLRCAnalysis::SetFitMethod(int id){
209 //Choose fit method
210         fFitMethod = id;
211 }
212 double AliLRCAnalysis::GetFitXmin() const
213 {
214         return fxFitMin;
215 }
216
217 double AliLRCAnalysis::GetFitXmax() const 
218 {
219         return fxFitMax;
220 }
221 double AliLRCAnalysis::GetArel() const {
222         return farel;
223 }
224
225 double AliLRCAnalysis::GetBrel() const {
226         return fbrel;
227 }
228
229 double AliLRCAnalysis::GetArelError() const {
230         return farelError;
231 }
232
233 double AliLRCAnalysis::GetBrelError() const {
234         return fbrelError;
235 }
236
237 double AliLRCAnalysis::GetXi2rel() const {
238         return fXi2rel;
239 }
240
241 double AliLRCAnalysis::GetAabs() const {
242         return faabs;
243 }
244
245 double AliLRCAnalysis::GetBabs() const {
246         return fbabs;
247 }
248
249 double AliLRCAnalysis::GetAabsError() const {
250         return faabsError;
251 }
252
253 double AliLRCAnalysis::GetBabsError() const {
254         return fbabsError;
255 }
256 double AliLRCAnalysis::GetXi2abs() const {
257         return fXi2abs;
258 }
259
260 void AliLRCAnalysis::Calculate(){
261 //Calculate all
262         double mnf = fPrf->GetMean(); 
263
264         // new fit for abs
265         AliLRCFit *fit1 = new AliLRCFit(fPrRel, fxFitMin/mnf, fxFitMax/mnf, 1.);
266         AliLRCFit *fit2 = new AliLRCFit(fPrAbs, fxFitMin, fxFitMax, 0.);
267         farel = fit1->Geta();
268         fbrel = fit1->Getb();
269
270         faabs = fit2->Geta();
271         fbabs = fit2->Getb();
272
273         if ( fFitMethod==0 )
274         {
275                 faabsError = fit2->Getda();
276                 fbabsError = fit2->Getdb();
277
278                 farelError = fit1->Getda();
279                 fbrelError = fit1->Getdb();
280         }
281         if ( fFitMethod==1 )
282         {
283                 faabsError = fit2->Getda1();
284                 fbabsError = fit2->Getdb1();
285
286                 farelError = fit1->Getda1();
287                 fbrelError = fit1->Getdb1();
288         }
289
290         fXi2rel = fit1->Gethi2();
291         fXi2abs = fit2->Gethi2();
292 }
293
294 void AliLRCAnalysis::DrawAbs() {
295 //Draw abs var hist with ALL info
296         int * mas = new int [fgkPlotFlags];
297         for ( int i = 0; i < fgkPlotFlags; i++ )
298                 mas[i] = 1;
299         DrawAbsPure( mas, 1 );
300         delete []mas;
301 }
302
303 void AliLRCAnalysis::DrawAbs( const int * const mas ) {
304 //Draw abs var hist with REQUESTED BY ARRAY info
305         DrawAbsPure( mas, 1 );
306 }
307
308 void AliLRCAnalysis::DrawAbsPure( const int * const mDrawArray, bool drawPaveLabel ) {
309         Calculate();
310 // Draw abs var histrogram
311         DrawHist( mDrawArray, drawPaveLabel, faabs, fbabs, faabsError, fbabsError, fPrAbs, 0 );
312 }
313
314 void AliLRCAnalysis::DrawRel() {
315 //Draw rel var hist with ALL info
316         int * mas = new int [fgkPlotFlags];
317         for ( int i = 0; i < fgkPlotFlags; i++ )
318                 mas[i] = 1;
319         DrawRelPure( mas, 1 );
320         delete []mas;
321 }
322
323 void AliLRCAnalysis::DrawRel( const int * const mas ) {
324 //Draw rel var hist with REQUESTED BY ARRAY info
325         DrawRelPure( mas, 1 );
326 }
327
328 void AliLRCAnalysis::DrawRelPure( const int * const mDrawArray, bool drawPaveLabel ) {
329         Calculate();
330 // Draw rel var histogram
331         DrawHist( mDrawArray, drawPaveLabel, farel, fbrel, farelError, fbrelError, fPrRel, 1 );
332 }
333 void AliLRCAnalysis::DrawHist( const int * const mDrawArray, bool drawPaveLabel, double aCoef, double bCoef,
334                 double aCoefError, double bCoefError, TH1D* profToDraw, int histType ) {
335 // Method called by DrawRelPure or DrawAbsPure to draw corresponding LRC  histo 
336         double mnf = fPrf->GetMean();
337         double y1, y2, x1, x2;
338         Int_t i, n;
339         char str[50];
340     //mnf = fPrf->GetMean();
341         profToDraw->SetStats(0);
342     //profToDraw->Fit("pol1", "", "", fxFitMin, fxFitMax);
343     //profToDraw->Fit("pol1");
344     //TF1 *fit1 = profToDraw->GetFunction("pol1");
345
346         //draw fit line
347         TF1 *f1 = 0;
348         if ( histType == 0 )    {
349                 f1 = new TF1( "f1", "[0]+[1]*x", fxFitMin, fxFitMax);
350                 f1->SetLineColor(kRed);
351                 f1->SetParameters(aCoef,bCoef); 
352         }
353         else if ( histType == 1 )       {
354                 f1 = new TF1( "f1", "[0]+[1]*(x-1)", fxFitMin/mnf, fxFitMax/mnf);
355                 f1->SetLineColor(kGreen);
356                 f1->SetParameters(aCoef,bCoef); 
357         }
358         //cout << " set draw params: a=" << aCoef << " b=" << bCoef << endl;
359
360         y1=profToDraw->GetBinContent(1)-profToDraw->GetBinError(1);
361         y2=profToDraw->GetBinContent(1)+profToDraw->GetBinError(1);
362         n=profToDraw->GetNbinsX();
363
364         for(i=2; i<=n; i++)
365         {
366                 if(profToDraw->GetBinContent(i)-profToDraw->GetBinError(i)<y1)
367                         y1=profToDraw->GetBinContent(i)-profToDraw->GetBinError(i);
368                 if(profToDraw->GetBinContent(i)+profToDraw->GetBinError(i)>y2)
369                         y2=profToDraw->GetBinContent(i)+profToDraw->GetBinError(i);
370         }
371
372         profToDraw->DrawCopy();
373         if(f1) f1->DrawCopy("same");
374         else return;
375         
376         x1 = profToDraw->GetXaxis()->GetXmin();
377         x2 = profToDraw->GetXaxis()->GetXmax();
378
379         if ( drawPaveLabel )
380         {       
381                 int nDatas = 0;
382                 for ( int j = 0; j < 9; j++)
383                         if ( mDrawArray[j] ) nDatas++;
384                 double aXshift = (x2-x1)/7;
385                 double aYshift = (y2-y1)/20;
386
387                 TPaveText *pt1 = new TPaveText(x1+(x2-x1)/2 + aXshift, y1+aYshift, x2-(x2-x1)/6 + aXshift, y1+(y2-y1)/3*2/9*nDatas + aYshift);
388                 snprintf(str,50, "Events = %i", fEntries);
389                 //sprintf(str, "Events = %i", fEntries);
390                 if ( mDrawArray[0] ) pt1->AddText(str);
391                 snprintf(str, 50,"a = %g #pm %g", GetRoundWithError( aCoef, aCoefError ), GetRoundWithPrecision(aCoefError, 2)); //fit1->GetParameter(0), fit1->GetParError(0));
392                 //sprintf(str, "a = %g #pm %g", GetRoundWithError( aCoef, aCoefError ), GetRoundWithPrecision(aCoefError, 2)); //fit1->GetParameter(0), fit1->GetParError(0));
393                 if ( mDrawArray[1] ) pt1->AddText(str);
394                 snprintf(str, 50,"b = %g #pm %g", GetRoundWithError( bCoef, bCoefError ), GetRoundWithPrecision(bCoefError, 2)); //fit1->GetParameter(1), fit1->GetParError(1));
395                 //sprintf(str, "b = %g #pm %g", GetRoundWithError( bCoef, bCoefError ), GetRoundWithPrecision(bCoefError, 2)); //fit1->GetParameter(1), fit1->GetParError(1));
396                 if ( mDrawArray[2] ) pt1->AddText(str);
397                 snprintf(str, 50,"#hat{#chi}^{2} = #chi^{2}/(n-2) = %g", GetRoundWithPrecision(fXi2abs, 3));
398                 //sprintf(str, "#hat{#chi}^{2} = #chi^{2}/(n-2) = %g", GetRoundWithPrecision(fXi2abs, 3));
399                 if ( mDrawArray[3] ) pt1->AddText(str);
400                 snprintf(str,50, "<%s> = %g " , fSx, GetRoundWithPrecision(fPrf->GetMean(), 3));
401                 //sprintf(str, "<%s> = %g " , fSx, GetRoundWithPrecision(fPrf->GetMean(), 3));
402                 if ( mDrawArray[4] ) pt1->AddText(str);
403                 
404                 snprintf(str,50, "<%s> = %g", fSy,  GetRoundWithPrecision(fPrb->GetMean(),3));
405                 //sprintf(str, "<%s> = %g", fSy,  GetRoundWithPrecision(fPrb->GetMean(),3));
406                 if ( mDrawArray[5] ) pt1->AddText(str);
407                 
408                 snprintf(str, 50,"<<%s>> = %g " , fSx, GetRoundWithPrecision(fPrf->GetRMS(), 3));
409                 //sprintf(str, "<<%s>> = %g " , fSx, GetRoundWithPrecision(fPrf->GetRMS(), 3));
410                 if ( mDrawArray[6] ) pt1->AddText(str);
411                 snprintf(str, 50,"<<%s>> = %g", fSy,  GetRoundWithPrecision(fPrb->GetRMS(), 3));
412                 //sprintf(str, "<<%s>> = %g", fSy,  GetRoundWithPrecision(fPrb->GetRMS(), 3));
413                 if ( mDrawArray[7] ) pt1->AddText(str);
414                 
415                 if ( fdptb ) {
416                         snprintf(str,50, "d%s = %g", fSy,  GetRoundWithPrecision(fdptb, 3));
417                         //sprintf(str, "d%s = %g", fSy,  GetRoundWithPrecision(fdptb, 3));
418                         if ( mDrawArray[8] ) pt1->AddText(str);
419                 }
420                 
421                 pt1->SetTextAlign(12);
422                 pt1->SetTextFont(42);
423                 pt1->SetFillColor(4000);
424                 //pt1->SetFillStyle(4100);
425                 pt1->SetShadowColor(4000);
426                 pt1->SetBorderSize(0);
427                 
428                 pt1->DrawClone("same");
429         }
430 }
431
432 void AliLRCAnalysis::SetGraphics() const {
433 // Set root graph style
434         TStyle tempSt;
435     gStyle->SetCanvasColor(10);
436     gStyle->SetFrameFillColor(10);
437     gStyle->SetStatColor(0);
438     gStyle->SetPadColor(0);
439 }
440
441 double AliLRCAnalysis::Integral(TH2D* const source, Int_t nbin) const {
442 // calculate the integrall for x bin and y bins of 2d histogramm
443     double sum = 0;
444     for(Int_t i = 1; i<=source->GetNbinsY(); i++) {
445         sum += source->GetYaxis()->GetBinCenter(i)*source->GetBinContent(nbin, i);
446     }
447     return sum;
448 }
449
450 void AliLRCAnalysis::SetErrors(TH2D* source, const char *name){
451 // Calculate errors for NN
452         TProfile* profX = (TProfile*) source->ProfileX(name, 1, source->GetNbinsY());
453         for(int i = 0; i < profX->GetNbinsX(); i++)             
454         {
455                 fPrAbs->SetBinContent(i, profX->GetBinContent(i));
456                 if(fPrf->GetBinContent(i)!=0)
457                 {
458                    fPrAbs->SetBinError(i,sqrt(profX->GetBinContent(i)/fPrf->GetBinContent(i)));
459                 }
460                 fPrRel->SetBinContent(i, fPrAbs->GetBinContent(i)/fPrb->GetMean());
461                 fPrRel->SetBinError(i,fPrAbs->GetBinError(i)/fPrb->GetMean());  
462         }       
463 }
464
465 void AliLRCAnalysis::SetErrors(TH2D* source, const char *name, double ptd, TH2D* nb){
466 //Calculate arrors for ptn and ptpt
467         TProfile* profX = (TProfile*) source->ProfileX(name, 1, source->GetNbinsY());
468         fdptb = ptd;
469         for(int i = 0; i < profX->GetNbinsX(); i++)             
470         {
471                 fPrAbs->SetBinContent(i, profX->GetBinContent(i));
472                 if(fPrf->GetBinContent(i)!=0)
473                 {
474                           pt = profX->GetBinContent(i);
475                           fPrAbs->SetBinError(i,ptd*sqrt(Integral(nb,i))/fPrf->GetBinContent(i));
476                 }
477                 fPrRel->SetBinContent(i, fPrAbs->GetBinContent(i)/fPrb->GetMean());
478                 fPrRel->SetBinError(i,fPrAbs->GetBinError(i)/fPrb->GetMean());  
479         }       
480
481 }
482
483 void AliLRCAnalysis::SetErrors(TH2D* source, const char *name, double ptd, const TProfile* nb){
484 //Calculate arrors for ptn and ptpt
485         TProfile* profX = (TProfile*) source->ProfileX(name, 1, source->GetNbinsY());
486         fdptb = ptd;
487         for(int i = 0; i < profX->GetNbinsX(); i++)             
488         {
489                 fPrAbs->SetBinContent(i, profX->GetBinContent(i));
490                 if(fPrf->GetBinContent(i)!=0)
491                 {
492                           pt = profX->GetBinContent(i);
493                           fPrAbs->SetBinError(i,ptd*sqrt(nb->GetBinContent(i)/fPrf->GetBinContent(i)));
494                 }
495                 fPrRel->SetBinContent(i, fPrAbs->GetBinContent(i)/fPrb->GetMean());
496                 fPrRel->SetBinError(i,fPrAbs->GetBinError(i)/fPrb->GetMean());  
497         }       
498
499 }
500
501 double AliLRCAnalysis::GetRoundWithError( double value, double error ) const {
502 //Rounding error and value with DEFAULT precision
503         return GetRoundValueErrorPrecision( value, error, 2 );
504 }
505 double AliLRCAnalysis::GetRoundWithError( double value, double error, int pres ) const {
506 //Rounding error and value with REQUESTED precision
507         return GetRoundValueErrorPrecision( value, error, pres );
508 }
509 double AliLRCAnalysis::GetRoundWithPrecision( double value, int pres ) const {
510 //Rounding error and value with requested precision
511         return GetRoundValueErrorPrecision( value, 0, pres );
512 }
513
514 double AliLRCAnalysis::GetRoundValueErrorPrecision( double value, double error, int pres ) const {
515         //Rounding error and value with requested precision
516         //value == value, error == error
517         //if single argument(value=0) - calculate without errors:
518
519         int i = 0;
520         double order = 1;
521         bool noError = false;
522         pres -= 1;
523         
524         //cout << "Before rounding: " << value << " " << error << endl;
525
526         if ( !error)
527         {
528                 error = value;
529                 noError = true;
530         }
531         while ( ((int)error)%10 == 0 && i < 6 )
532         {
533                 i++;
534                 error*=10;
535                 order*=10;
536         }
537
538         for ( int j = 0; j < pres; j++ )
539         {
540                 error*=10;
541                 order*=10;
542                 i++;
543         }
544         int adding = 0;
545         if ( ((int)(error*10))%10 > 4 && ((int)(error*10))%10 != 9 ) //trouble: if we round 19 to 20 - zero disappeares!
546                 adding = 1;
547         error = (double)((int)error + adding)/order;
548         
549         if ( noError )
550         {
551                 //cout << "After rounding: " << error << endl;
552                 return error;
553         }
554         else
555         {
556                 for ( int j = 0; j < i; j++ )
557                         value*=10;
558                 
559                 adding = 0;
560                 if ( ((int)(value*10))%10 > 4 && ((int)(value*10))%10 != 9 )
561                         adding = 1;
562                 value = (double)((int)value + adding)/order;
563                 
564                 //cout << "After rounding: " << value << " " << error << endl;
565                 return value; //taking into account ERROR
566         }
567 }
568
569 TH1D* AliLRCAnalysis::GetAbsHisto() const
570 {
571 //Returns final histo in absolute variables
572 return fPrAbs;
573 }
574 TH1D* AliLRCAnalysis::GetRelHisto() const
575 {
576 //Returns final histo in relative variables
577 return fPrRel;
578 }
579
580 TH1D* AliLRCAnalysis::GetForwardValueDist() const
581 {
582 //Returns destribution of value used in forward window
583 return fPrf;
584 }
585 TH1D* AliLRCAnalysis::GetBackwardValueDist() const
586 {
587 //Returns destribution of value used in backward window
588 return fPrb;
589 }
590
591