2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
17 //------------------------------------------------------------------------------
18 // Implementation of AliPerfAnalyzeInvPt class. It analyzes the output of
19 // AliPerformancePtCalib.cxx and AliPerformancePtCalibMC.cxx: via
20 // AliPerformancePtCalib::Analyse() or AliPerformancePtCalibMC::Analyse()
22 // Projection of charge/pt vs theta and vs phi resp. histoprams will be fitted with either
23 // polynomial or gaussian fit function to extract minimum position of charge/pt.
24 // Fit options and theta, phi bins can be set by user.
25 // Attention: use the Set* functions of AliPerformancePtCalib.h and AliPerformancePtCalibMC.h
26 // when running AliPerformancePtCalib*::Analyse()
27 // The result of the analysis (histograms/graphs) are stored in the folder which is
28 // a data member of AliPerformancePtCalib*.
30 // Author: S.Schuchmann 11/13/2009
31 //------------------------------------------------------------------------------
35 // after running the performance task of AliPerformancePtCalib* or AliPerformancePtCalibMC*
36 // read the file, and get component:
38 TFile f("Output.root");
39 AliPerformancePtCalib *compObj = (AliPerformancePtCalib*)coutput->FindObject("AliPerformancePtCalib");
41 AliPerformancePtCalibMC *compObj = (AliPerformancePtCalibMC*)coutput->FindObject("AliPerformancePtCalibMC");
43 // set phi and theta bins for fitting and analyse comparison data
44 compObj->SetProjBinsTheta(thetaBins,nThetaBins,minPhi, maxPhi);
45 compObj->SetProjBinsPhi(phiBins,nPhiBins,minTheta,maxTheta);
46 compObj->SetMakeFitOption(kFALSE,exclRange,fitRange);
47 compObj->SetDoRebin(rebin);
49 //for details see functions of class AliPerformancePtCalib*
51 // the output histograms/graphs will be stored in the folder "folderRes"
52 compObj->GetAnalysisFolder()->ls("*");
54 // user can save whole comparison object (or only folder with anlysed histograms)
55 // in the seperate output file (e.g.)
56 TFile fout("Analysed_InvPt.root","recreate");
57 compObj->Write(); // compObj->GetAnalysisFolder()->Write();
68 #include "TGraphErrors.h"
70 #include "TObjArray.h"
72 #include "AliPerfAnalyzeInvPt.h"
74 ClassImp(AliPerfAnalyzeInvPt)
78 //____________________________________________________________________________________________________________________________________________
79 Double_t AliPerfAnalyzeInvPt::Polynomial(Double_t *x, const Double_t *par)
81 // fit function for fitting of 1/pt with a polynomial of 4th order, rejecting points between +/-par[4]
84 if (x[0] > -par[4] && x[0] < par[4]) {
89 return par[2]+par[0]*pow((x[0]-par[1]),2)+par[3]*pow((x[0]-par[1]),4);
92 //____________________________________________________________________________________________________________________________________________
93 Double_t AliPerfAnalyzeInvPt::PolynomialRejP(Double_t *x, const Double_t *par)
95 // 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)
97 Double_t pos = par[5];
98 Double_t neg = - par[5];
102 if (x[0] > neg && x[0] < pos) {
107 return par[2]+par[0]*pow((x[0]-par[1]),2)+par[3]*pow((x[0]-par[1]),4);
110 //____________________________________________________________________________________________________________________________________________
111 Double_t AliPerfAnalyzeInvPt::InvGauss(Double_t *x, const Double_t *par)
113 // fit function for fitting of 1/pt with gaussian, rejecting points between +/-par[4]
114 if (x[0] > -par[6] && x[0] < par[6]) {
119 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)) ;
121 //____________________________________________________________________________________________________________________________________________
122 Double_t AliPerfAnalyzeInvPt::InvGaussRejP(Double_t *x, const Double_t *par)
124 // 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)
125 Double_t pos = par[7];//0.12;
126 Double_t neg = - par[7];//0.12;
130 if (x[0] > neg && x[0] < pos) {
136 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) ;
140 //_____________________________________________________________________________________________________________________________________________
141 AliPerfAnalyzeInvPt::AliPerfAnalyzeInvPt(const Char_t* name, const Char_t* title):
154 //histograms,graphs and functions
155 fHistH2InvPtTheta(0),
164 // Default constructor
166 fFitGaus = kFALSE;// flag for gaussian fit
167 fNThetaBins = 0; // theta bins for projections
168 fNPhiBins = 0; //phi bins for projections
169 fRange = 0; //fit range
170 fExclRange = 0; //fit exclusion range
171 fDoRebin = kFALSE; // flag for rebin
172 fRebin = 0; // bins for rebin
174 for(Int_t i=0;i<100;i++){
176 fHistFitTheta[i] = NULL;
177 fHistFitPhi[i] = NULL;
183 //______________________________________________________________________________________________________________________________________
184 void AliPerfAnalyzeInvPt::InitGraphs(Double_t *binsXTheta,Double_t *fitParamTheta,Double_t *errFitParamTheta,Double_t *binsXPhi,Double_t *fitParamPhi,Double_t *errFitParamPhi){
186 // initialize graphs for storing fit parameter
188 fGrMinPosTheta = new TGraphErrors(fNThetaBins,binsXTheta,fitParamTheta,0, errFitParamTheta);
189 fGrMinPosTheta->SetMarkerStyle(20);
190 fGrMinPosTheta->SetMarkerColor(2);
191 fGrMinPosTheta->SetLineColor(2);
192 fGrMinPosTheta->GetYaxis()->SetTitle("min pos (GeV/c)^{-1}");
193 fGrMinPosTheta->GetXaxis()->SetTitle("#theta bin no.");
194 fGrMinPosTheta->GetYaxis()->SetTitleOffset(1.2);
195 fGrMinPosTheta->SetTitle("#theta bins");
197 fGrMinPosPhi = new TGraphErrors(fNPhiBins,binsXPhi,fitParamPhi,0,errFitParamPhi);
198 fGrMinPosPhi->SetMarkerStyle(20);
199 fGrMinPosPhi->SetMarkerColor(4);
200 fGrMinPosPhi->SetLineColor(4);
201 fGrMinPosPhi->GetYaxis()->SetTitle("min pos (GeV/c)^{-1}");
202 fGrMinPosPhi->GetXaxis()->SetTitle("#phi bin no.");
203 fGrMinPosPhi->GetYaxis()->SetTitleOffset(1.2);
204 fGrMinPosPhi->SetTitle("#phi bins");
207 //______________________________________________________________________________________________________________________________________
209 void AliPerfAnalyzeInvPt::InitFitFcn(){
211 fFitMinPos = new TF1("fFitMinPos", Polynomial,-4.0,4.0,5);
212 fFitMinPos->SetLineColor(4);
213 fFitMinPos->SetLineWidth(1);
214 fFitMinPos->SetParameter(0,1.0);
215 fFitMinPos->SetParameter(1,0.0);
216 fFitMinPos->SetParameter(2,0.0);
218 fFitMinPosRejP = new TF1("fFitMinPosRejP",PolynomialRejP,-4.0,4.0,6);
219 fFitMinPosRejP->SetLineColor(2);
220 fFitMinPosRejP->SetLineWidth(1);
221 fFitMinPosRejP->SetParameter(0,1.0);
222 fFitMinPosRejP->SetParameter(1,0.0);
223 fFitMinPosRejP->SetParameter(2,0.0);
226 fFitInvGauss = new TF1("fFitInvGauss", InvGauss,-4.0,4.0,6);
227 fFitInvGauss->SetLineColor(4);
228 fFitInvGauss->SetLineWidth(1);
229 fFitInvGauss->SetParameter(2,1.0);
231 fFitInvGaussRejP = new TF1("fFitInvGaussRejP", InvGaussRejP,-4.0,4.0,7);
232 fFitInvGaussRejP->SetLineColor(2);
233 fFitInvGaussRejP->SetLineWidth(1);
234 fFitInvGaussRejP->SetParameter(2,1.0);
237 //______________________________________________________________________________________________________________________________________
238 void AliPerfAnalyzeInvPt::StartAnalysis(const TH2F *histThetaInvPt, const TH2F *histPhiInvPt, TObjArray* aFolderObj){
239 //start analysis: fitting charge/pt spectra
243 if(!histThetaInvPt) {
244 Printf("warning: no 1/pt histogram to analyse in theta bins!");
248 Printf("warning: no 1/pt histogram to analyse in phit bins!");
252 Double_t thetaBins[9] = {0.77,0.97,1.17,1.37,1.57,1.77,1.97,2.17,2.37}; // default theta bins
253 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}; // default phi bins
255 Int_t nThetaBins = 9;
257 Double_t range = 1.0; // fit range
258 Double_t exclRange = 0.13; // range of point rejection in fit
259 Bool_t fitGaus=kFALSE; // fit with Gaussian or with polynomial (kFALSE)
262 fNThetaBins = nThetaBins;
263 for(Int_t k = 0;k<fNThetaBins;k++){
264 fThetaBins[k] = thetaBins[k];
266 Printf("warning: for theta bins no user intput! bins are set to default values.");
270 fNPhiBins = nPhiBins;
271 for(Int_t k = 0;k<fNPhiBins;k++){
272 fPhiBins[k] = phiBins[k];
275 Printf("warning: for phi bins no user intput! bins are set to default values.");
281 fExclRange = exclRange;
283 Printf("warning: no fit range is set. default fitting conditions are used.");
287 //fit parameter arrays
288 Double_t fitParamTheta[100],fitParamPhi[100];
289 Double_t errFitParamTheta[100], errFitParamPhi[100];
290 Double_t binsXTheta[100],binsXPhi[100];
293 fHistH2InvPtTheta = (TH2F*)histThetaInvPt->Clone("invPtVsTheta");
294 fHistH2InvPtPhi = (TH2F*)histPhiInvPt->Clone("invPtVsPhi");
296 // initialize fit functions
299 // canvas for 2D histograms (input)
300 TCanvas *thephiCan =new TCanvas("thephiCan","theta and phi vs invPt",800,400);
301 thephiCan->Divide(2,1);
303 fHistH2InvPtTheta->SetTitle("#theta vs charge/pt for selected #phi range");
304 fHistH2InvPtTheta->GetXaxis()->SetTitle("charge/p_{t} (GeV/c)^{-1}");
305 fHistH2InvPtTheta->GetYaxis()->SetTitle("#theta (rad)");
306 fHistH2InvPtTheta->DrawCopy("col1z");
309 fHistH2InvPtPhi->SetTitle("#phi vs charge/pt for selected #theta range");
310 fHistH2InvPtPhi->GetXaxis()->SetTitle("charge/p_{t} (GeV/c)^{-1}");
311 fHistH2InvPtPhi->GetYaxis()->SetTitle("#phi (rad)");
312 fHistH2InvPtPhi->DrawCopy("col1z");
314 // canvas for 1D histograms (projections)
315 TCanvas *theCan =new TCanvas("theCan","invPt theta bins",1200,900);
316 theCan->Divide((fNThetaBins+2)/2,2);
317 TCanvas *phiCan =new TCanvas("phiCan","invPt phi bins",1200,900);
318 phiCan->Divide((fNPhiBins+2)/2,2);
322 // analyse charge/pt in bins of theta
323 for(Int_t i=0;i<fNThetaBins;i++){
325 //set name for each projection
326 TString name = "fit_theta_";
328 fHistH2InvPtTheta->SetName(name.Data());
330 //find bins for projection
331 Int_t firstBin= fHistH2InvPtTheta->GetYaxis()->FindBin(fThetaBins[i]);
332 if(i>0) firstBin +=1;
333 Int_t lastBin = fHistH2InvPtTheta->GetYaxis()->FindBin(fThetaBins[i+1]);
334 if( i == fNThetaBins-1) {
335 firstBin= fHistH2InvPtTheta->GetYaxis()->FindBin(fThetaBins[0]);
336 lastBin = fHistH2InvPtTheta->GetYaxis()->FindBin(fThetaBins[i]);
340 fHistFitTheta[i] = (TH1D*)fHistH2InvPtTheta->ProjectionX("_px",firstBin,lastBin,"e");
342 Char_t titleTheta[50];
343 if(TMath::Abs(i-(fNThetaBins-1))<0.5) snprintf(titleTheta,50,"charge/pt (GeV/c) integrated over #theta");
344 else snprintf(titleTheta,50,"charge/pt (GeV/c) for #theta range: %1.3f - %1.3f",fThetaBins[i],fThetaBins[i+1]);
346 fHistFitTheta[i]->SetTitle(titleTheta);
348 Double_t invPtMinPos = 0;
349 Double_t invPtMinPosErr = 0;
350 Double_t invPtMinPosImpr = 0;
351 Double_t invPtMinPosErrImpr = 0;
353 if(fDoRebin) fHistFitTheta[i]->Rebin(fRebin);
357 Printf("making polynomial fit in charge/pt in theta bins");
358 theCan->cd(countPad);
359 MakeFit(fHistFitTheta[i],fFitMinPos, invPtMinPos,invPtMinPosErr, fExclRange,fRange);
360 MakeFitBetter(fHistFitTheta[i],fFitMinPosRejP, invPtMinPosImpr, invPtMinPosErrImpr, invPtMinPos,fExclRange,fRange);
362 //plot projection and fit
363 fHistFitTheta[i]->DrawCopy();
364 fFitMinPos->DrawCopy("L,same");
365 fFitMinPosRejP->DrawCopy("L,same");
368 Printf("making gauss fit in charge/pt in theta bins");
369 theCan->cd(countPad);
370 MakeFitInvGauss(fHistFitTheta[i],fFitInvGauss, invPtMinPos, invPtMinPosErr,fExclRange,fRange);
371 MakeFitInvGaussBetter(fHistFitTheta[i],fFitInvGaussRejP, invPtMinPosImpr, invPtMinPosErrImpr,invPtMinPos,fExclRange,fRange);
373 //plot projection and fit
374 fHistFitTheta[i]->DrawCopy();
375 fFitInvGauss->DrawCopy("L,same");
376 fFitInvGaussRejP->DrawCopy("L,same");
378 //add objects for analysis folder
379 aFolderObj->Add(fHistFitTheta[i]);
381 //store fit parameter
382 fitParamTheta[i] = invPtMinPosImpr;
383 errFitParamTheta[i] = invPtMinPosErrImpr;
385 //count bins and pad number for displaying
386 binsXTheta[i] = i+1.0;
394 // analyse charge/pt in bins of phi
396 for(Int_t i=0;i<fNPhiBins;i++){
398 //set name for each projection
399 TString name = "fit_phi_";
401 fHistH2InvPtPhi->SetName(name.Data());
403 //find bins for projection
404 Int_t firstBin = fHistH2InvPtPhi->GetYaxis()->FindBin(fPhiBins[i]);
405 if(i>0) firstBin +=1;
406 Int_t lastBin = fHistH2InvPtPhi->GetYaxis()->FindBin(fPhiBins[i+1]);
407 if(TMath::Abs(i-(fNPhiBins-1))<0.5){
408 firstBin = fHistH2InvPtPhi->GetYaxis()->FindBin(fPhiBins[0]);
409 lastBin = fHistH2InvPtPhi->GetYaxis()->FindBin(fPhiBins[i]);
413 fHistFitPhi[i] = (TH1D*) fHistH2InvPtPhi->ProjectionX("_px",firstBin,lastBin,"e");
416 if(TMath::Abs(i-(fNPhiBins-1))<0.5) snprintf(titlePhi,50,"charge/pt (GeV/c) integrated over #phi");
417 else snprintf(titlePhi,50,"charge/pt (GeV/c) for #phi range: %1.3f - %1.3f",fPhiBins[i],fPhiBins[i+1]);
419 fHistFitPhi[i]->SetTitle(titlePhi);
421 Double_t invPtMinPos = 0;
422 Double_t invPtMinPosErr = 0;
423 Double_t invPtMinPosImpr = 0;
424 Double_t invPtMinPosErrImpr = 0;
425 if(fDoRebin) fHistFitPhi[i]->Rebin(fRebin);
429 Printf("making polynomial fit in charge/pt in phi bins");
430 phiCan->cd(countPad);
431 MakeFit(fHistFitPhi[i],fFitMinPos, invPtMinPos, invPtMinPosErr,fExclRange,fRange);
432 MakeFitBetter(fHistFitPhi[i],fFitMinPosRejP, invPtMinPosImpr, invPtMinPosErrImpr,invPtMinPos,fExclRange,fRange);
434 //plot projection and fit
435 fHistFitPhi[i]->DrawCopy();
436 fFitMinPos->DrawCopy("L,same");
437 fFitMinPosRejP->DrawCopy("L,same");
441 Printf("making gauss fit in charge/pt in phi bins");
442 phiCan->cd(countPad);
443 MakeFitInvGauss(fHistFitPhi[i],fFitInvGauss, invPtMinPos, invPtMinPosErr, exclRange,fRange);
444 MakeFitInvGaussBetter(fHistFitPhi[i],fFitInvGaussRejP, invPtMinPosImpr, invPtMinPosErrImpr, invPtMinPos, fExclRange,fRange);
446 //plot projection and fit
447 fHistFitPhi[i]->DrawCopy();
448 fFitInvGauss->DrawCopy("L,same");
449 fFitInvGaussRejP->DrawCopy("L,same");
451 //add objects for analysis folder
452 aFolderObj->Add(fHistFitPhi[i]);
454 //store fit parameter
455 fitParamPhi[i] = invPtMinPosImpr;
456 errFitParamPhi[i] = invPtMinPosErrImpr;
458 //count bins and pad number for displaying
463 //initialize graphs for displaying the fit parameter
464 InitGraphs(binsXTheta,fitParamTheta,errFitParamTheta,binsXPhi,fitParamPhi,errFitParamPhi);
466 //plot fit values = minimum positions of charge/pt
467 TCanvas *canFitVal = new TCanvas("canFitVal","min pos histos",800,400);
468 canFitVal->Divide(2,1);
471 fGrMinPosTheta->Draw("AP");
473 fGrMinPosPhi->Draw("AP");
475 Printf("AliPerfAnalyzeInvPt::StartAnalysis: ***** NOTE: last bin is always fit result of integral over all angle BINS which have been set by user! angle ranges of 2D projections are kept.*****");
477 //add objects to folder
478 aFolderObj->Add(fGrMinPosTheta);
479 aFolderObj->Add(fGrMinPosPhi);
480 aFolderObj->Add(fHistH2InvPtTheta);
481 aFolderObj->Add(fHistH2InvPtPhi);
487 //____________________________________________________________________________________________________________________________________________
489 void AliPerfAnalyzeInvPt::MakeFit(TH1D *hproy, TF1 * fitpb, Double_t &mean, Double_t &errMean, Double_t &excl,Double_t &range)
491 // fit charge/pt and extract minimum position
493 fitpb->SetRange(-range,range);
494 fitpb->SetParLimits(1,-0.05,0.05);
495 fitpb->SetParameter(0,1.0);
496 fitpb->FixParameter(4,excl);
497 fitpb->FixParameter(2,0.0);
499 hproy->Fit(fitpb,"RM");
500 mean = fitpb->GetParameter(1);
501 errMean = fitpb->GetParError(1);
502 if(!mean) errMean = 0.0;
504 //____________________________________________________________________________________________________________________________________________
505 void AliPerfAnalyzeInvPt::MakeFitBetter(TH1D *hproy, TF1 * fitpb2, Double_t &mean, Double_t &errMean, Double_t &f, Double_t &excl, Double_t &range)
507 // adjust fit range to minimum position of AliPerfAnalyzeInvPt::MakeFit and fit 1/pt and extract new minimum position
509 fitpb2->FixParameter(5,excl);
510 fitpb2->FixParameter(4,f);
511 // fitpb2->FixParameter(2,0.0);
512 fitpb2->SetRange(-range+f,range+f);
513 fitpb2->SetParLimits(1,-0.05,0.05);
514 fitpb2->SetParameter(0,1.0);
515 hproy->Fit(fitpb2,"RM");
516 mean = fitpb2->GetParameter(1);
517 errMean = fitpb2->GetParError(1);
518 if(!mean) errMean = 0.0;
521 //____________________________________________________________________________________________________________________________________________
522 void AliPerfAnalyzeInvPt::MakeFitInvGauss(TH1D *hproy, TF1 * fitpb, Double_t &mean, Double_t &errMean, Double_t &excl,Double_t &range)
524 // fit charge/pt and extract minimum position
526 fitpb->FixParameter(6,excl);
527 fitpb->SetRange(-range,range);
528 fitpb->SetParameter(0,-1.0);
529 fitpb->SetParameter(2,1.0);
530 fitpb->SetParameter(3,25000.0);
531 fitpb->SetParameter(4,-1.0);
532 fitpb->SetParLimits(1,-0.02,0.02);
533 hproy->Fit(fitpb,"RM");
534 mean = fitpb->GetParameter(1);
535 errMean = fitpb->GetParError(1);
536 if(!mean) errMean = 0.0;
539 //____________________________________________________________________________________________________________________________________________
540 void AliPerfAnalyzeInvPt::MakeFitInvGaussBetter(TH1D *hproy, TF1 * fitpb2, Double_t &mean, Double_t &errMean, Double_t &f, Double_t &excl, Double_t &range)
542 // adjust fit range to minimum position of AliPerfAnalyzeInvPt::MakeFitInvGauss and fit 1/pt and extract new minimum position
544 fitpb2->FixParameter(7,excl);
545 fitpb2->FixParameter(6,f);
546 fitpb2->SetRange(-range+f,range+f);
547 fitpb2->SetParameter(0,-1.0);
548 fitpb2->SetParameter(2,1.0);
549 fitpb2->SetParameter(3,25000.0);
550 fitpb2->SetParameter(4,-1.0);
551 fitpb2->SetParLimits(1,-0.02,0.02);
552 hproy->Fit(fitpb2,"RM");
553 mean = fitpb2->GetParameter(1);
554 errMean = fitpb2->GetParError(1);
555 if(!mean) errMean = 0.0;
560 //____________________________________________________________________________________________________________________________________________
562 // set variables for analysis
563 void AliPerfAnalyzeInvPt::SetProjBinsPhi(const Double_t *phiBinArray,const Int_t nphBins){
565 // set theta bins for Analyse()
566 //set phi bins as array and set number of this array which is equal to number of bins analysed
567 //the last analysed bin will always be the projection from first to last bin in the array
571 for(Int_t k = 0;k<fNPhiBins;k++){
572 fPhiBins[k] = phiBinArray[k];
574 Printf("AliPerfAnalyzeInvPt::SetProjBinsPhi:number of bins in phi set to %i",fNPhiBins);
576 else Printf("AliPerfAnalyzeInvPt::SetProjBinsPhi: number of bins in theta NOT set. Default values are used.");
578 //____________________________________________________________________________________________________________________________________________
579 void AliPerfAnalyzeInvPt::SetProjBinsTheta(const Double_t *thetaBinArray, const Int_t nthBins){
580 // set theta bins for Analyse()
581 //set theta bins as array and set number of this array which is equal to number of bins analysed
582 //the last analysed bin will always be the projection from first to last bin in the array
584 fNThetaBins = nthBins;
585 for(Int_t k = 0;k<fNThetaBins;k++){
586 fThetaBins[k] = thetaBinArray[k];
588 Printf("AliPerfAnalyzeInvPt::SetProjBinsTheta: number of bins in theta set to %i",fNThetaBins);
590 else Printf("AliPerfAnalyzeInvPt::SetProjBinsTheta: number of bins in theta NOT set. Default values are used.");
592 //____________________________________________________________________________________________________________________________________________
593 void AliPerfAnalyzeInvPt::SetMakeFitOption(const Bool_t setGausFit, const Double_t exclusionR,const Double_t fitR ){
595 //set the fit options:
596 //for usage of gaussian function instead of polynomial (default) set setGausFit=kTRUE
597 //set the range of rejection of points around 0 via exclusionR
598 //set the fit range around 0 with fitR
600 fFitGaus = setGausFit;
601 fExclRange = exclusionR;
604 if(fFitGaus) Printf("AliPerfAnalyzeInvPt::SetMakeFitOption: set MakeGausFit with fit range %2.3f and exclusion range in fabs(1/pt): %2.3f",fRange,fExclRange);
605 else Printf("AliPerfAnalyzeInvPt::SetMakeFitOption: set standard polynomial fit with fit range %2.3f and exclusion range in fabs(1/pt): %2.3f",fRange,fExclRange);
607 else Printf(" AliPerfAnalyzeInvPt::SetMakeFitOption: no user input. Set standard polynomial fit with fit range 1.0 and exclusion range 0.13.");