]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TPC/AliPerfAnalyzeInvPt.cxx
Change names of the output conatainer and files
[u/mrichter/AliRoot.git] / PWG1 / TPC / AliPerfAnalyzeInvPt.cxx
1 //------------------------------------------------------------------------------
2 // Implementation of AliPerfAnalyzeInvPt class. It analyzes the output of
3 // AliPerformancePtCalib.cxx  and AliPerformancePtCalibMC.cxx:
4 // Projection of 1/pt vs theta and vs phi resp. histoprams will be fitted with either
5 // polynomial or gaussian fit function to extract minimum position of 1/pt.
6 // Fit options and theta, phi bins can be set by user.
7 // Attention: use the Set* functions of AliPerformancePtCalib.h and AliPerformancePtCalibMC.h
8 // when running AliPerformancePtCalib*::Analyse()
9 // The result of the analysis (histograms/graphs) are stored in the folder which is
10 // a data member of AliPerformancePtCalib*.
11 //
12 // Author: S.Schuchmann 11/13/2009 
13 //------------------------------------------------------------------------------
14
15 /*
16  
17 // after running comparison task, read the file, and get component
18 gROOT->LoadMacro("$ALICE_ROOT/PWG1/Macros/LoadMyLibs.C");
19 LoadMyLibs();
20
21 TFile f("Output.root");
22 AliPerformancePtCalib * compObj = (AliPerformancePtCalib*)coutput->FindObject("AliPerformancePtCalib");
23  
24 // analyse comparison data
25 compObj->Analyse();
26
27 // the output histograms/graphs will be stored in the folder "folderRes" 
28 compObj->GetAnalysisFolder()->ls("*");
29
30 // user can save whole comparison object (or only folder with anlysed histograms) 
31 // in the seperate output file (e.g.)
32 TFile fout("Analysed_InvPt.root","recreate");
33 compObj->Write(); // compObj->GetAnalysisFolder()->Write();
34 fout.Close();
35
36 */
37
38 #include "TNamed.h"
39 #include "TF1.h"
40 #include "TH1F.h"
41 #include "TH2F.h"
42 #include "TMath.h"
43 #include "TH1D.h"
44 #include "TGraphErrors.h"
45 #include "TCanvas.h"
46 #include "TObjArray.h"
47
48 #include "AliPerfAnalyzeInvPt.h"
49
50 ClassImp(AliPerfAnalyzeInvPt)
51
52
53 // fit functions
54 //____________________________________________________________________________________________________________________________________________
55    Double_t AliPerfAnalyzeInvPt::Polynomial(Double_t *x, const Double_t *par)
56 {
57    // fit function for fitting of 1/pt with a polynomial of 4th order, rejecting points between  +/-par[4]
58
59   
60    if (x[0] > -par[4] && x[0] < par[4]) {
61       TF1::RejectPoint();
62       return 0;
63    }
64    
65    return  par[2]+par[0]*pow((x[0]-par[1]),2)+par[3]*pow((x[0]-par[1]),4);
66
67 }
68 //____________________________________________________________________________________________________________________________________________
69 Double_t AliPerfAnalyzeInvPt::PolynomialRejP(Double_t *x, const Double_t *par)
70 {
71    // fit function for fitting of 1/pt with a polynomial of 4th order to improve result (fit range and rejection zone is adjusted to first guess of minimum position)
72      
73    Double_t pos  = par[5];
74    Double_t neg = - par[5];
75    pos += par[4];
76    neg += par[4];
77    
78    if (x[0] > neg && x[0] < pos) {
79       TF1::RejectPoint();
80       return 0;
81    }
82    
83    return   par[2]+par[0]*pow((x[0]-par[1]),2)+par[3]*pow((x[0]-par[1]),4);
84  
85 }
86 //____________________________________________________________________________________________________________________________________________
87 Double_t AliPerfAnalyzeInvPt::InvGauss(Double_t *x, const Double_t *par)
88 {
89    // fit function for fitting of 1/pt with gaussian, rejecting points between  +/-par[4]
90    if (x[0] > -par[6] && x[0] < par[6]) {
91       TF1::RejectPoint();
92       return 0;
93    }
94
95    return par[3]+par[0]*(TMath::Exp(-0.5*(TMath::Power((x[0]-par[1])/par[2], 2.0)))+par[4]*pow((x[0]-par[1]),2)+par[5]*pow((x[0]-par[1]),4)) ;
96 }
97 //____________________________________________________________________________________________________________________________________________
98 Double_t AliPerfAnalyzeInvPt::InvGaussRejP(Double_t *x, const Double_t *par)
99 {
100    // fit function for fitting of 1/pt with gaussian to improve result (fit range and rejection zone is adjusted to first guess of minimum position)
101    Double_t pos  = par[7];//0.12;
102    Double_t neg = - par[7];//0.12;
103    pos += par[6];
104    neg += par[6];
105   
106    if (x[0] > neg && x[0] < pos) {
107       TF1::RejectPoint();
108       return 0;
109    }
110
111
112    return par[3]+par[0]*TMath::Exp(-0.5*(TMath::Power((x[0]-par[1])/par[2], 2.0)))+par[4]*pow((x[0]-par[1]),2)+par[5]*pow((x[0]-par[1]),4) ;
113 }
114
115
116 //_____________________________________________________________________________________________________________________________________________
117 AliPerfAnalyzeInvPt::AliPerfAnalyzeInvPt():
118    TNamed("AliPerfAnalyzeInvPt","AliPerfAnalyzeInvPt"),
119    fNThetaBins(0), 
120    fNPhiBins(0),
121    fRange(0),
122    fExclRange(0),
123    fFitGaus(0) ,
124    fDoRebin(0),
125    fRebin(0),
126    fHistH2InvPtTheta(0),
127    fHistH2InvPtPhi(0), 
128    fGrMinPosTheta(0),
129    fGrMinPosPhi(0),
130    fFitMinPos(0),
131    fFitMinPosRejP(0),
132    fFitInvGauss(0),
133    fFitInvGaussRejP(0)
134 {
135    // Default constructor
136   
137    fFitGaus = kFALSE;
138    fNThetaBins = 0;
139    fNPhiBins = 0;
140    fRange = 0;
141    fExclRange = 0;
142    fFitGaus = 0;
143    fDoRebin = kFALSE;
144    fRebin = 0;
145    
146    for(Int_t i=0;i<100;i++){
147       
148       fHistFitTheta[i] = NULL;
149       fHistFitPhi[i] = NULL;
150    }
151  
152    
153 }
154 //_____________________________________________________________________________________________________________________________________________
155 AliPerfAnalyzeInvPt::AliPerfAnalyzeInvPt(Char_t* name="AliAnalyzeInvPt",Char_t* title="AliAnalyzeInvPt"):
156    TNamed(name, title),
157    fNThetaBins(0), 
158    fNPhiBins(0),
159    fRange(0),
160    fExclRange(0),
161    fFitGaus(0) ,
162    fDoRebin(0),
163    fRebin(0),
164    fHistH2InvPtTheta(0),
165    fHistH2InvPtPhi(0), 
166    fGrMinPosTheta(0),
167    fGrMinPosPhi(0),
168    fFitMinPos(0),
169    fFitMinPosRejP(0),
170    fFitInvGauss(0),
171    fFitInvGaussRejP(0)
172 {
173    
174    fFitGaus = kFALSE;
175    fNThetaBins = 0;
176    fNPhiBins =0;
177    fRange = 0;
178    fExclRange = 0;
179    fFitGaus = 0;
180    fDoRebin = kFALSE;
181    fRebin = 0;
182     
183    for(Int_t i=0;i<100;i++){
184     
185       fHistFitTheta[i] = NULL;
186       fHistFitPhi[i] = NULL;
187    }
188   
189 }
190
191
192 //______________________________________________________________________________________________________________________________________
193 void AliPerfAnalyzeInvPt::InitGraphs(Double_t *binsXTheta,Double_t *fitParamTheta,Double_t *errFitParamTheta,Double_t *binsXPhi,Double_t *fitParamPhi,Double_t *errFitParamPhi){
194
195    // initialize graphs for storing fit parameter
196   
197    fGrMinPosTheta = new TGraphErrors(fNThetaBins,binsXTheta,fitParamTheta,0, errFitParamTheta);  
198    fGrMinPosTheta->SetMarkerStyle(20);
199    fGrMinPosTheta->SetMarkerColor(2);
200    fGrMinPosTheta->SetLineColor(2);
201    fGrMinPosTheta->GetYaxis()->SetTitle("min pos (Gev/c)^{-1}");
202    fGrMinPosTheta->GetXaxis()->SetTitle("#theta bin no.");
203    fGrMinPosTheta->GetYaxis()->SetTitleOffset(1.2);   
204    fGrMinPosTheta->SetTitle("#theta bins ");
205
206    fGrMinPosPhi = new TGraphErrors(fNPhiBins,binsXPhi,fitParamPhi,0,errFitParamPhi);  
207    fGrMinPosPhi->SetMarkerStyle(20);
208    fGrMinPosPhi->SetMarkerColor(4);
209    fGrMinPosPhi->SetLineColor(4);
210    fGrMinPosPhi->GetYaxis()->SetTitle("min pos (Gev/c)^{-1}");
211    fGrMinPosPhi->GetXaxis()->SetTitle("#phi bin no.");
212    fGrMinPosPhi->GetYaxis()->SetTitleOffset(1.2);   
213    fGrMinPosPhi->SetTitle("#phi bins ");
214 }
215
216 //______________________________________________________________________________________________________________________________________
217
218 void AliPerfAnalyzeInvPt::InitFitFcn(){
219    // fit functions
220    fFitMinPos = new TF1("fFitMinPos", Polynomial,-4.0,4.0,5);
221    fFitMinPos->SetLineColor(4);
222    fFitMinPos->SetLineWidth(1);
223    fFitMinPos->SetParameter(0,1.0);
224    fFitMinPos->SetParameter(1,0.0);
225    fFitMinPos->SetParameter(2,0.0);
226
227    fFitMinPosRejP = new TF1("fFitMinPosRejP",PolynomialRejP,-4.0,4.0,6);
228    fFitMinPosRejP->SetLineColor(2);
229    fFitMinPosRejP->SetLineWidth(1);
230    fFitMinPosRejP->SetParameter(0,1.0);
231    fFitMinPosRejP->SetParameter(1,0.0);
232    fFitMinPosRejP->SetParameter(2,0.0);
233  
234   
235    fFitInvGauss = new TF1("fFitInvGauss", InvGauss,-4.0,4.0,6);
236    fFitInvGauss->SetLineColor(4);
237    fFitInvGauss->SetLineWidth(1);
238    fFitInvGauss->SetParameter(2,1.0);
239   
240    fFitInvGaussRejP = new TF1("fFitInvGaussRejP", InvGaussRejP,-4.0,4.0,7);
241    fFitInvGaussRejP->SetLineColor(2);
242    fFitInvGaussRejP->SetLineWidth(1);
243    fFitInvGaussRejP->SetParameter(2,1.0);
244  
245 }
246 //______________________________________________________________________________________________________________________________________
247 void AliPerfAnalyzeInvPt::StartAnalysis(const TH2F *histThetaInvPt, const TH2F *histPhiInvPt, TObjArray* aFolderObj){
248    //start analysis: fitting 1/pt spectra
249
250   
251    
252    if(!histThetaInvPt) {
253       Printf("warning: no 1/pt histogram to analyse in theta bins!");
254    }
255    if(!histPhiInvPt) {
256       Printf("warning: no 1/pt histogram to analyse in phit bins!");
257    }
258
259    Double_t thetaBins[9] = {0.77,0.97,1.17,1.37,1.57,1.77,1.97,2.17,2.37};                 // theta bins
260    Double_t phiBins[13] = {0.0,0.5,1.0,1.5,2.0,2.5,3.0,3.5,4.0,4.5,5.0,5.5,6.5};           // phi bins
261   
262    Int_t nThetaBins = 9;
263    Int_t nPhiBins = 13;
264    Double_t range = 1.0;        // fit range
265    Double_t  exclRange  = 0.13; // range of point rejection in fit
266    Bool_t fitGaus=kFALSE;       // fit with Gaussian or with polynomial (kFALSE)
267    
268    if(!fNThetaBins){
269       fNThetaBins = nThetaBins;
270       for(Int_t k = 0;k<fNThetaBins;k++){
271          fThetaBins[k] = thetaBins[k];
272       }
273       Printf("warning: for theta bins no user intput! bins are set to default values.");
274    }
275
276    if(!fNPhiBins){
277       fNPhiBins = nPhiBins;
278       for(Int_t k = 0;k<fNPhiBins;k++){
279          fPhiBins[k] = phiBins[k];
280       }
281     
282       Printf("warning: for phi bins no user intput! bins are set to default values.");
283     
284    }
285
286    if(!fRange){
287       fRange = range;
288       fExclRange = exclRange;
289       fFitGaus = fitGaus;
290       Printf("warning: no fit range is set. default fitting conditions are used.");
291    }
292
293       
294    //param arrays
295    Double_t fitParamTheta[100],fitParamPhi[100];
296    Double_t errFitParamTheta[100], errFitParamPhi[100];
297    Double_t binsXTheta[100],binsXPhi[100];
298
299    fHistH2InvPtTheta  = (TH2F*)histThetaInvPt->Clone("invPtVsTheta");
300    fHistH2InvPtPhi    = (TH2F*)histPhiInvPt->Clone("invPtVsPhi");
301
302    InitFitFcn();
303          
304    TCanvas *theCan =new TCanvas("theCan","invPt theta bins",1200,900);
305    theCan->Divide((fNThetaBins+2)/2,2);
306    TCanvas *phiCan =new TCanvas("phiCan","invPt phi bins",1200,900);
307    phiCan->Divide((fNPhiBins+2)/2,2);
308    Int_t countPad = 1;
309          
310    // analyse 1/pt in bins of theta 
311    for(Int_t i=0;i<fNThetaBins;i++){
312       TString name = "fit_theta_";
313       name +=i;
314       Int_t firstBin= fHistH2InvPtTheta->GetYaxis()->FindBin(fThetaBins[i]);
315       if(i>0) firstBin +=1;
316       Int_t lastBin = fHistH2InvPtTheta->GetYaxis()->FindBin(fThetaBins[i+1]);
317       if( i == fNThetaBins-1) {
318          firstBin= fHistH2InvPtTheta->GetYaxis()->FindBin(fThetaBins[0]);
319          lastBin = fHistH2InvPtTheta->GetYaxis()->FindBin(fThetaBins[i]);
320       }
321     
322       fHistH2InvPtTheta->SetName(name.Data());
323       fHistFitTheta[i] =  (TH1F*)fHistH2InvPtTheta->ProjectionX("_px",firstBin,lastBin,"e");
324       
325       Char_t titleTheta[50];
326       if(TMath::Abs(i-(fNThetaBins-1))<0.5) sprintf(titleTheta,"1/pt (GeV/c) integrated over #theta");
327       else  sprintf(titleTheta,"1/pt (GeV/c) for #theta range: %1.3f - %1.3f",fThetaBins[i],fThetaBins[i+1]);
328       
329       fHistFitTheta[i]->SetTitle(titleTheta);
330    
331       Double_t invPtMinPos  = 0;
332       Double_t invPtMinPosErr = 0;
333       Double_t invPtMinPosImpr  = 0;
334       Double_t invPtMinPosErrImpr = 0;
335    
336       if(fDoRebin) fHistFitTheta[i]->Rebin(fRebin);
337       //start fitting
338       if(!fFitGaus){
339          Printf("making polynomial fit in 1/pt in theta bins");
340          theCan->cd(countPad);
341          MakeFit(fHistFitTheta[i],fFitMinPos, invPtMinPos,invPtMinPosErr, fExclRange,fRange);
342          MakeFitBetter(fHistFitTheta[i],fFitMinPosRejP, invPtMinPosImpr, invPtMinPosErrImpr, invPtMinPos,fExclRange,fRange);
343
344         
345          fHistFitTheta[i]->DrawCopy();
346          fFitMinPos->DrawCopy("L,same");
347          fFitMinPosRejP->DrawCopy("L,same");
348       }
349       else{
350          Printf("making gauss fit in 1/pt in theta bins");
351          theCan->cd(countPad);
352          MakeFitInvGauss(fHistFitTheta[i],fFitInvGauss, invPtMinPos, invPtMinPosErr,fExclRange,fRange);
353          MakeFitInvGaussBetter(fHistFitTheta[i],fFitInvGaussRejP, invPtMinPosImpr, invPtMinPosErrImpr,invPtMinPos,fExclRange,fRange);
354       
355          fHistFitTheta[i]->DrawCopy();
356          fFitInvGauss->DrawCopy("L,same");
357          fFitInvGaussRejP->DrawCopy("L,same");
358       }
359     
360       aFolderObj->Add(fHistFitTheta[i]);
361     
362       fitParamTheta[i] = invPtMinPosImpr;
363       errFitParamTheta[i] = invPtMinPosErrImpr;
364     
365       binsXTheta[i] = i+1.0;
366       countPad++;
367    }
368       
369       
370    countPad = 1;
371
372    
373    // analyse 1/pt in bins of phi 
374   
375    for(Int_t i=0;i<fNPhiBins;i++){
376       TString name = "fit_phi_";
377       name +=i;
378     
379       fHistH2InvPtPhi->SetName(name.Data());
380       Int_t  firstBin = fHistH2InvPtPhi->GetYaxis()->FindBin(fPhiBins[i]);
381       if(i>0) firstBin +=1;
382       Int_t   lastBin =  fHistH2InvPtPhi->GetYaxis()->FindBin(fPhiBins[i+1]);
383       if(TMath::Abs(i-(fNPhiBins-1))<0.5){
384          firstBin = fHistH2InvPtPhi->GetYaxis()->FindBin(fPhiBins[0]);
385          lastBin =  fHistH2InvPtPhi->GetYaxis()->FindBin(fPhiBins[i]);
386       }
387       fHistFitPhi[i] =  (TH1F*) fHistH2InvPtPhi->ProjectionX("_px",firstBin,lastBin,"e");
388       
389       Char_t titlePhi[50];
390       if(TMath::Abs(i-(fNPhiBins-1))<0.5) sprintf(titlePhi,"1/pt (GeV/c) integrated over #phi");
391       else  sprintf(titlePhi,"1/pt (GeV/c) for #phi range: %1.3f - %1.3f",fPhiBins[i],fPhiBins[i+1]);
392      
393       fHistFitPhi[i]->SetTitle(titlePhi);
394   
395       Double_t invPtMinPos  = 0;
396       Double_t invPtMinPosErr = 0;
397       Double_t invPtMinPosImpr  = 0;
398       Double_t invPtMinPosErrImpr = 0;
399      if(fDoRebin) fHistFitPhi[i]->Rebin(fRebin);
400       if(!fFitGaus){
401          Printf("making polynomial fit in 1/pt in phi bins");
402          phiCan->cd(countPad);
403          MakeFit(fHistFitPhi[i],fFitMinPos, invPtMinPos, invPtMinPosErr,fExclRange,fRange);
404          MakeFitBetter(fHistFitPhi[i],fFitMinPosRejP, invPtMinPosImpr, invPtMinPosErrImpr,invPtMinPos,fExclRange,fRange);
405          
406          fHistFitPhi[i]->DrawCopy();
407          fFitMinPos->DrawCopy("L,same");
408          fFitMinPosRejP->DrawCopy("L,same");
409
410       }
411       else {
412          Printf("making gauss fit in 1/pt in phi bins");
413          phiCan->cd(countPad);
414          MakeFitInvGauss(fHistFitPhi[i],fFitInvGauss, invPtMinPos, invPtMinPosErr, exclRange,fRange);
415          MakeFitInvGaussBetter(fHistFitPhi[i],fFitInvGaussRejP, invPtMinPosImpr, invPtMinPosErrImpr, invPtMinPos, fExclRange,fRange);
416          
417          fHistFitPhi[i]->DrawCopy();
418          fFitInvGauss->DrawCopy("L,same");
419          fFitInvGaussRejP->DrawCopy("L,same");
420       }
421     
422       aFolderObj->Add(fHistFitPhi[i]);
423     
424       fitParamPhi[i] = invPtMinPosImpr;
425       errFitParamPhi[i] = invPtMinPosErrImpr;
426     
427       binsXPhi[i] = i+1.0;
428       countPad++;
429    }
430       
431    InitGraphs(binsXTheta,fitParamTheta,errFitParamTheta,binsXPhi,fitParamPhi,errFitParamPhi);
432
433    //plot fit values = minimum positions of charge/pt
434    TCanvas *canFitVal = new TCanvas("canFitVal","min pos histos",800,400);
435    canFitVal->Divide(2,1);
436
437    canFitVal->cd(1);
438    fGrMinPosTheta->Draw("AP");
439    canFitVal->cd(2);
440    fGrMinPosPhi->Draw("AP");
441
442    Printf("AliPerfAnalyzeInvPt: NOTE: last bin is always fit result  of integral over all angle ranges which have been set by user!");
443
444    //add objects to folder
445    aFolderObj->Add(fGrMinPosTheta);
446    aFolderObj->Add(fGrMinPosPhi);
447    aFolderObj->Add(fHistH2InvPtTheta);
448    aFolderObj->Add(fHistH2InvPtPhi);
449   
450   
451 }
452
453
454 //____________________________________________________________________________________________________________________________________________
455
456 void AliPerfAnalyzeInvPt::MakeFit(TH1F *hproy, TF1 * fitpb, Double_t &mean, Double_t &errMean, Double_t &excl,Double_t &range)
457 {
458    // fit 1/pt and extract minimum position
459    
460    fitpb->SetRange(-range,range);
461    fitpb->SetParLimits(1,-0.05,0.05);
462    fitpb->SetParameter(0,1.0);
463    fitpb->FixParameter(4,excl);
464    fitpb->FixParameter(2,0.0);
465     
466    hproy->Fit(fitpb,"RM");
467    mean = fitpb->GetParameter(1);
468    errMean = fitpb->GetParError(1);
469    if(!mean)  errMean = 0.0;
470 }
471 //____________________________________________________________________________________________________________________________________________
472 void AliPerfAnalyzeInvPt::MakeFitBetter(TH1F *hproy, TF1 * fitpb2, Double_t &mean, Double_t &errMean, Double_t &f, Double_t &excl, Double_t &range)
473 {
474    // adjust fit range to minimum position of AliPerfAnalyzeInvPt::MakeFit and fit 1/pt and extract new minimum position
475    
476    fitpb2->FixParameter(5,excl);
477    fitpb2->FixParameter(4,f);
478    // fitpb2->FixParameter(2,0.0);
479    fitpb2->SetRange(-range+f,range+f);
480    fitpb2->SetParLimits(1,-0.05,0.05);
481    fitpb2->SetParameter(0,1.0);
482    hproy->Fit(fitpb2,"RM");
483    mean = fitpb2->GetParameter(1);
484    errMean = fitpb2->GetParError(1);
485    if(!mean)   errMean = 0.0;
486
487 }
488 //____________________________________________________________________________________________________________________________________________
489 void AliPerfAnalyzeInvPt::MakeFitInvGauss(TH1F *hproy, TF1 * fitpb, Double_t &mean, Double_t &errMean, Double_t &excl,Double_t &range)
490 {
491    // fit 1/pt and extract minimum position
492    
493    fitpb->FixParameter(6,excl);
494    fitpb->SetRange(-range,range);
495    fitpb->SetParameter(0,-1.0);
496    fitpb->SetParameter(2,1.0);
497    fitpb->SetParameter(3,25000.0);
498    fitpb->SetParameter(4,-1.0);
499    fitpb->SetParLimits(1,-0.02,0.02);
500    hproy->Fit(fitpb,"RM");
501    mean = fitpb->GetParameter(1);
502    errMean = fitpb->GetParError(1);
503    if(!mean)   errMean = 0.0;
504   
505 }
506 //____________________________________________________________________________________________________________________________________________
507 void AliPerfAnalyzeInvPt::MakeFitInvGaussBetter(TH1F *hproy, TF1 * fitpb2, Double_t &mean, Double_t &errMean, Double_t &f, Double_t &excl, Double_t &range)
508 {
509    // adjust fit range to minimum position of AliPerfAnalyzeInvPt::MakeFitInvGauss and fit 1/pt and extract new minimum position
510    
511    fitpb2->FixParameter(7,excl);
512    fitpb2->FixParameter(6,f);
513    fitpb2->SetRange(-range+f,range+f);
514    fitpb2->SetParameter(0,-1.0);
515    fitpb2->SetParameter(2,1.0);
516    fitpb2->SetParameter(3,25000.0);
517    fitpb2->SetParameter(4,-1.0);
518    fitpb2->SetParLimits(1,-0.02,0.02);
519    hproy->Fit(fitpb2,"RM");
520    mean = fitpb2->GetParameter(1);
521    errMean = fitpb2->GetParError(1);
522    if(!mean)   errMean = 0.0;
523    
524 }
525
526
527 //____________________________________________________________________________________________________________________________________________
528
529 // set variables for analysis
530 void AliPerfAnalyzeInvPt::SetProjBinsPhi(const Double_t *phiBinArray,const Int_t nphBins){
531
532    // set theta bins for Analyse()
533    //set phi bins as array and set number of this array which is equal to number of bins analysed
534    //the last analysed bin will always be the projection from first to last bin in the array
535    if(nphBins){
536       fNPhiBins = nphBins;
537   
538       for(Int_t k = 0;k<fNPhiBins;k++){
539          fPhiBins[k] = phiBinArray[k];
540       }
541       Printf("AliPerfAnalyzeInvPt::SetProjBinsPhi:number of bins in phi set to %i",fNPhiBins);
542    }
543    else Printf("AliPerfAnalyzeInvPt::SetProjBinsPhi: number of bins in theta NOT set. Default values are used.");
544 }
545 //____________________________________________________________________________________________________________________________________________
546 void AliPerfAnalyzeInvPt::SetProjBinsTheta(const Double_t *thetaBinArray, const Int_t nthBins){
547    // set theta bins for Analyse()
548    //set theta bins as array and set number of this array which is equal to number of bins analysed
549    //the last analysed bin will always be the projection from first to last bin in the array
550    if(nthBins){
551       fNThetaBins = nthBins;
552       for(Int_t k = 0;k<fNThetaBins;k++){
553          fThetaBins[k] = thetaBinArray[k];
554       }
555       Printf("AliPerfAnalyzeInvPt::SetProjBinsTheta: number of bins in theta set to %i",fNThetaBins);
556    }
557    else Printf("AliPerfAnalyzeInvPt::SetProjBinsTheta: number of bins in theta NOT set. Default values are used.");
558 }
559 //____________________________________________________________________________________________________________________________________________
560 void AliPerfAnalyzeInvPt::SetMakeFitOption(const Bool_t setGausFit, const Double_t exclusionR,const Double_t fitR ){
561
562    //set the fit options:
563    //for usage of gaussian function instead of polynomial (default) set setGausFit=kTRUE
564    //set the range of rejection of points around 0 via exclusionR
565    //set the fit range around 0 with fitR
566    if(fitR){
567       fFitGaus = setGausFit;
568       fExclRange  = exclusionR;
569       fRange = fitR;
570   
571       if(fFitGaus) Printf("set MakeGausFit with fit range %2.3f and exclusion range in 1/pt: %2.3f",fRange,fExclRange);
572       else  Printf("set standard polynomial fit with fit range %2.3f and exclusion range in 1/pt: %2.3f",fRange,fExclRange);
573    }
574    else Printf(" AliPerfAnalyzeInvPt::SetMakeFitOption: no user input. Set standard polynomial fit with fit range 1.0 and exclusion range 0.13.");
575 }