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