]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TPC/AliPerformancePtCalib.cxx
update by Simone Schuchmann
[u/mrichter/AliRoot.git] / PWG1 / TPC / AliPerformancePtCalib.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
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  **************************************************************************/
16
17 //------------------------------------------------------------------------------
18 // Implementation of AliPerformancePtCalib class. It fills histograms with ESD or
19 // TPC track information to study charge/pt spectra.
20 // To analyse the output of AliPerformancePtCalib use AliPerfAnalyzeInvPt class:
21 // Projection of charge/pt vs theta and vs phi resp. histoprams will be fitted
22 // with either polynomial or gaussian fit function to extract minimum position of
23 // charge/pt.
24 // Fit options and theta, phi bins can be set by user.
25 // Attention: use the Set functions of AliPerformancePtCalib  when running
26 // AliPerformancePtCalib::Analyse()
27 // The result of the analysis (histograms/graphs) are stored in the folder which is
28 // a data member of AliPerformancePtCalib.
29 //
30 // Author: S.Schuchmann 11/13/2009 
31 //         sschuchm@ikf.uni-frankfurt.de
32 //------------------------------------------------------------------------------
33
34 /*
35  
36 // after running the performance task, read the file, and get component
37
38 TFile f("Output.root");
39 AliPerformancePtCalib * compObj = (AliPerformancePtCalib*)coutput->FindObject("AliPerformancePtCalib");
40  
41 // set phi and theta bins for fitting and analyse comparison data
42 compObj->SetProjBinsTheta(thetaBins,nThetaBins,minPhi, maxPhi);
43 compObj->SetProjBinsPhi(phiBins,nPhiBins,minTheta,maxTheta);
44 compObj->SetMakeFitOption(kFALSE,exclRange,fitRange);
45 compObj->SetDoRebin(rebin);
46 compObj->Analyse();
47 //details see functions of this class
48
49 // the output histograms/graphs will be stored in the folder "folderRes" 
50 compObj->GetAnalysisFolder()->ls("*");
51
52 // user can save whole comparison object (or only folder with anlysed histograms) 
53 // in the seperate output file (e.g.)
54 TFile fout("Analysed_InvPt.root","recreate");
55 compObj->Write(); // compObj->GetAnalysisFolder()->Write();
56 fout.Close();
57
58 */
59 #include "TH1F.h"
60 #include "TH2F.h"
61 #include "THnSparse.h"
62 #include "TList.h"
63 #include "TMath.h"
64 #include "TFolder.h"
65
66 #include "AliESDEvent.h" 
67 #include "AliESDtrack.h"
68 #include "AliESDfriendTrack.h"
69 #include "AliESDfriend.h"
70 #include "AliESDtrackCuts.h"
71 #include "AliESDpid.h"
72
73 #include "AliPerfAnalyzeInvPt.h"
74 #include "AliPerformancePtCalib.h"
75
76 using namespace std;
77
78 ClassImp(AliPerformancePtCalib)
79
80 //________________________________________________________________________
81    AliPerformancePtCalib::AliPerformancePtCalib():
82       AliPerformanceObject("AliPerformancePtCalib"),
83   
84       // option parameter for AliPerformancePtCalib::Analyse()
85       fNThetaBins(0), 
86       fNPhiBins(0),
87       fMaxPhi(0),
88       fMinPhi(0),
89       fMaxTheta(0),
90       fMinTheta(0),
91       fRange(0),
92       fExclRange(0),
93       fFitGaus(0) ,
94       fDoRebin(0),
95       fRebin(0),
96       // option for user defined charge/pt shift
97       fShift(0),
98       fDeltaInvP(0),
99   
100       //options for cuts
101       fOptTPC(0),
102       fESDcuts(0),
103       fPions(0),
104       fEtaAcceptance(0),
105       fCutsRC(0),
106       fCutsMC(0),
107       fList(0),
108       // histograms
109       fHistInvPtPtThetaPhi(0),
110       
111       fHistPtShift0(0),
112       fHistPrimaryVertexPosX(0),
113       fHistPrimaryVertexPosY(0),
114       fHistPrimaryVertexPosZ(0),
115       fHistTrackMultiplicity(0),
116       fHistTrackMultiplicityCuts(0),
117
118       fHistTPCMomentaPosP(0),
119       fHistTPCMomentaNegP(0),
120       fHistTPCMomentaPosPt(0),
121       fHistTPCMomentaNegPt(0),
122       fHistUserPtShift(0),
123       fHistdedxPions(0),
124       //esd track cuts
125       fESDTrackCuts(0),
126       //pid
127       fESDpid(0),
128       // analysis folder 
129       fAnalysisFolder(0)
130 {
131    
132    // Dummy constructor
133
134    fShift = kFALSE;                       // shift in charge/pt yes/no
135    fDeltaInvP = 0.00;                     // shift value
136    //options for cuts
137    fOptTPC =  kTRUE;                      // read TPC tracks yes/no
138    fESDcuts = kFALSE;
139    fPions = kFALSE;
140    fCutsRC = NULL;
141    fCutsMC = NULL;
142    
143    //set options for esd track cuts
144    fEtaAcceptance = 0.8;
145    
146    // options for function AliPerformancePtCalib::Analyse()
147    fFitGaus = kFALSE;// use gaussian function for fitting charge/pt yes/no
148    fNThetaBins = 0; //number of theta bins
149    fNPhiBins = 0; //number of phi bins
150    fRange = 0; //fit range around 0
151    fMaxPhi = 6.5;// max phi for 2D projection on theta and charge/pt axis
152    fMinPhi = 0.0;// min phi for 2D projection on theta and charge/pt axis
153    fMaxTheta = 3.0;// max theta for 2D projection on phi and charge/pt axis
154    fMinTheta = 0.0;// min theta for 2D projection on phi and charge/pt axis
155    fExclRange =0; //range of rejection of points around 0
156    fDoRebin = kFALSE;
157    fRebin = 0;
158
159    Init();
160 }
161
162 //________________________________________________________________________
163 AliPerformancePtCalib::AliPerformancePtCalib(Char_t * name="AliPerformancePtCalib",Char_t* title ="AliPerformancePtCalib"):
164    AliPerformanceObject(name,title),
165    
166    // option parameter for AliPerformancePtCalib::Analyse()
167    fNThetaBins(0), 
168    fNPhiBins(0),
169    fMaxPhi(0),
170    fMinPhi(0),
171    fMaxTheta(0),
172    fMinTheta(0),
173    fRange(0),
174    fExclRange(0),
175    fFitGaus(0) ,
176    fDoRebin(0),
177    fRebin(0),
178    // option parameter for user defined charge/pt shift
179    fShift(0),
180    fDeltaInvP(0),
181    //options for cuts
182    fOptTPC(0),
183    fESDcuts(0),
184    fPions(0),
185    fEtaAcceptance(0),
186    fCutsRC(0),
187    fCutsMC(0),
188    fList(0),
189    // histograms
190    fHistInvPtPtThetaPhi(0),
191    fHistPtShift0(0),
192    fHistPrimaryVertexPosX(0),
193    fHistPrimaryVertexPosY(0),
194    fHistPrimaryVertexPosZ(0),
195    fHistTrackMultiplicity(0),
196    fHistTrackMultiplicityCuts(0),
197
198    fHistTPCMomentaPosP(0),
199    fHistTPCMomentaNegP(0),
200    fHistTPCMomentaPosPt(0),
201    fHistTPCMomentaNegPt(0),
202    fHistUserPtShift(0), 
203    fHistdedxPions(0),
204    //esd track cuts                                                                                                      
205    fESDTrackCuts(0),
206    //pid
207    fESDpid(0),
208    // analysis folder 
209    fAnalysisFolder(0)
210    
211   
212 {
213    // Constructor
214    fShift = kFALSE;
215    fDeltaInvP = 0.00;
216    //options for cuts
217    fOptTPC =  kTRUE;                      // read TPC tracks yes/no
218    fESDcuts = kFALSE;
219    fPions = kFALSE;
220    fCutsRC = NULL;
221    fCutsMC = NULL;
222    
223    //esd track cut options
224    fEtaAcceptance = 0.8;
225    
226    // options for function AliPerformancePtCalibMC::Analyse()
227    fFitGaus = kFALSE;// use gaussian function for fitting charge/pt yes/no
228    fNThetaBins = 0; //number of theta bins
229    fNPhiBins = 0; //number of phi bins
230    fMaxPhi = 6.5;// max phi for 2D projection on theta and charge/pt axis
231    fMinPhi = 0.0;// min phi for 2D projection on theta and charge/pt axis
232    fMaxTheta = 3.0;// max theta for 2D projection on phi and charge/pt axis
233    fMinTheta = 0.0;// min theta for 2D projection on phi and charge/pt axis
234    fRange = 0; //fit range around 0
235    fExclRange =0; //range of rejection of points around 0
236    fDoRebin = kFALSE;
237    fRebin = 0;
238   
239    Init();
240 }
241 //________________________________________________________________________
242 AliPerformancePtCalib::~AliPerformancePtCalib() { 
243    //
244    // destructor
245    //
246
247    if(fList) delete fList;
248    // histograms
249    if( fHistInvPtPtThetaPhi)  delete fHistInvPtPtThetaPhi;fHistInvPtPtThetaPhi=0;
250    if(fHistPtShift0)          delete fHistPtShift0;fHistPtShift0=0; 
251    if(fHistPrimaryVertexPosX) delete fHistPrimaryVertexPosX;fHistPrimaryVertexPosX=0; 
252    if(fHistPrimaryVertexPosY) delete fHistPrimaryVertexPosY;fHistPrimaryVertexPosY=0; 
253    if(fHistPrimaryVertexPosZ) delete fHistPrimaryVertexPosZ;fHistPrimaryVertexPosZ=0; 
254    if(fHistTrackMultiplicity) delete fHistTrackMultiplicity;fHistTrackMultiplicity=0; 
255    if(fHistTrackMultiplicityCuts) delete fHistTrackMultiplicityCuts;fHistTrackMultiplicityCuts=0; 
256
257    if(fHistTPCMomentaPosP)  delete fHistTPCMomentaPosP;fHistTPCMomentaPosP=0; 
258    if(fHistTPCMomentaNegP)  delete fHistTPCMomentaNegP;fHistTPCMomentaNegP=0; 
259    if(fHistTPCMomentaPosPt) delete fHistTPCMomentaPosPt;fHistTPCMomentaPosPt=0; 
260    if(fHistTPCMomentaNegPt) delete fHistTPCMomentaNegPt ;fHistTPCMomentaNegPt=0; 
261    if(fHistUserPtShift)     delete fHistUserPtShift;fHistUserPtShift=0; 
262    //esd track cuts
263    if(fESDTrackCuts) delete fESDTrackCuts;
264    //pid
265    if(fESDpid) delete fESDpid;fESDpid=0;
266    //analysis folder 
267    if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0; 
268 }
269
270 //________________________________________________________________________
271 void AliPerformancePtCalib::Init() 
272 {
273    // Create histograms
274    // Called once
275
276    fList = new TList();
277    // init folder
278    fAnalysisFolder = CreateFolder("folderPt_TPC","Analysis Pt Resolution Folder");
279    
280    // Primary Vertex:
281    fHistPrimaryVertexPosX       = new TH1F("fHistPrimaryVertexPosX", "Primary Vertex Position X;Primary Vertex Position X (cm);Events",100,-0.5,0.5);
282    fList->Add(fHistPrimaryVertexPosX);
283    fHistPrimaryVertexPosY       = new TH1F("fHistPrimaryVertexPosY", "Primary Vertex Position Y;Primary Vertex Position Y (cm);Events",100,-0.5,0.5);
284    fList->Add(fHistPrimaryVertexPosY);
285    fHistPrimaryVertexPosZ       = new TH1F("fHistPrimaryVertexPosZ", "Primary Vertex Position Z;Primary Vertex Position Z (cm);Events",200,-2.0,2.0);
286    fList->Add(fHistPrimaryVertexPosZ);
287    // Multiplicity:
288    fHistTrackMultiplicity     = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 250, 0, 250);
289    fList->Add(fHistTrackMultiplicity);
290    fHistTrackMultiplicityCuts = new TH1F("fHistTrackMultiplicityCuts", "Multiplicity distribution;Number of tracks after cuts;Events", 250, 0, 250);
291    fList->Add(fHistTrackMultiplicityCuts);
292  
293    // momentum histos
294    //pt shift 0 only needed if shift in 1/pt is applied
295    fHistPtShift0 = new TH1F("fHistPtShift0","1/pt dN/pt vs. pt of ESD track  ",800,-40.0,40.0);
296    fList->Add(fHistPtShift0);
297  
298    // THnSparse for 1/pt and pt spectra vs angles
299    const   Int_t invPtDims = 4;
300    fMaxPhi = 6.52;
301    fMinPhi = 0.0;
302    fMaxTheta = 3.0;
303    fMinTheta = 0.0;
304    
305    Double_t xminInvPt[invPtDims] = {-4.5,-40.0,fMinTheta,fMinPhi};
306    Double_t xmaxInvPt[invPtDims] = {4.5,40.0,fMaxTheta,fMaxPhi};
307    Int_t  binsInvPt[invPtDims] = {450,400,150,163};
308
309   
310    fHistInvPtPtThetaPhi = new THnSparseF("fHistInvPtPtThetaPhi","1/pt vs pt vs #theta vs #phi ",invPtDims,binsInvPt,xminInvPt,xmaxInvPt);
311    fList->Add(fHistInvPtPtThetaPhi);
312    
313    // momentum test histos
314    fHistTPCMomentaPosP  =  new TH2F("fHistTPCMomentaPosP","TPC p vs global esd track p pos",300,0.0,15.0,300,0.0,15.0);
315    fList->Add(fHistTPCMomentaPosP);
316    fHistTPCMomentaNegP  =  new TH2F("fHistTPCMomentaNegP","TPC p vs global esd track p neg",300,0.0,15.0,300,0.0,15.0);
317    fList->Add(fHistTPCMomentaNegP);
318    fHistTPCMomentaPosPt =  new TH2F("fHistTPCMomentaPosPt","TPC pt vs global esd track pt pos",300,0.0,15.0,300,0.0,15.0);
319    fList->Add(fHistTPCMomentaPosPt);
320    fHistTPCMomentaNegPt =  new TH2F("fHistTPCMomentaNegPt","TPC pt vs global esd track pt neg",300,0.0,15.0,300,0.0,15.0);
321    fList->Add(fHistTPCMomentaNegPt);
322
323    //user pt shift check
324    fHistUserPtShift = new TH1F("fHistUserPtShift","user defined shift in 1/pt",100,-0.5,1.5);
325    fList->Add(fHistUserPtShift);
326    //pid by dedx
327    fHistdedxPions = new TH2F ("fHistdedxPions","dEdx of pions ident. via kPID vs signed Pt",300,-15.05,15.05,200,0.0,400.0);
328    fList->Add(fHistdedxPions);
329    //pid
330    fESDpid = NULL;
331
332    // esd track cuts  
333    fESDTrackCuts =NULL;
334 }
335
336 //________________________________________________________________________
337 void AliPerformancePtCalib::SetPtShift(const Double_t shiftVal ) {
338    //set user defined shift in charge/pt
339
340    if(shiftVal) { fShift=kTRUE; fDeltaInvP = shiftVal; } 
341 }
342
343 //________________________________________________________________________
344 void AliPerformancePtCalib::Exec(AliMCEvent* const /*mcEvent*/, AliESDEvent *const esdEvent, AliESDfriend * const /*esdFriend*/, const Bool_t /*bUseMC*/, const Bool_t /*bUseESDfriend*/)
345 {
346    //exec: read esd or tpc
347
348    if(!fESDTrackCuts) Printf("no esd track cut");
349    
350    if (!esdEvent) {
351       Printf("ERROR: Event not available");
352       return;
353    }
354
355    if (!(esdEvent->GetNumberOfTracks()))  return;
356
357    
358    //vertex info for cut
359    const AliESDVertex *vtx = esdEvent->GetPrimaryVertex();
360    if (!vtx->GetStatus()) return ;
361      
362    //histo fo user defined shift in charge/pt 
363    if(fShift) fHistUserPtShift->Fill(fDeltaInvP);
364    
365    //trakc multiplicity
366    fHistTrackMultiplicity->Fill(esdEvent->GetNumberOfTracks());
367
368
369    // read primary vertex info
370    Double_t tPrimaryVtxPosition[3];
371    const AliESDVertex *primaryVtx = esdEvent->GetPrimaryVertexTPC();
372  
373    tPrimaryVtxPosition[0] = primaryVtx->GetXv();
374    tPrimaryVtxPosition[1] = primaryVtx->GetYv();
375    tPrimaryVtxPosition[2] = primaryVtx->GetZv();
376   
377    fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
378    fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
379    fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
380
381
382    //_fill histos for pt spectra and shift of transverse momentum
383    Int_t count=0;
384  
385    for(Int_t j = 0;j<esdEvent->GetNumberOfTracks();j++){// track loop
386       AliESDtrack *esdTrack = esdEvent->GetTrack(j);
387       if(!esdTrack) continue;
388     
389     
390       if(fESDcuts){
391          if(!fESDTrackCuts->AcceptTrack(esdTrack))continue;
392       }
393        
394       
395       // fill histos
396       if(fOptTPC){ //TPC tracks
397          const AliExternalTrackParam *tpcTrack = esdTrack->GetTPCInnerParam(); 
398          if(!tpcTrack) continue;
399          if(fabs(tpcTrack->Eta())>= fEtaAcceptance) continue;
400       
401          Double_t signedPt = tpcTrack->GetSignedPt();
402
403          // pid
404          if(fPions){
405            fESDpid= new AliESDpid();
406            fESDpid->GetTPCResponse().SetBetheBlochParameters(0.0283086,2.63394e+01,5.04114e-11, 2.12543e+00,4.88663e+00);
407            
408            if( TMath::Abs(fESDpid->NumberOfSigmasTPC(esdTrack,AliPID::kPion)) >1) continue;
409            fHistdedxPions->Fill(signedPt,esdTrack->GetTPCsignal());
410          }
411          
412          Double_t invPt = 0.0;
413          if(signedPt) {
414             invPt = 1.0/signedPt;
415
416             fHistPtShift0->Fill(signedPt);
417         
418             if(fShift){Printf("user shift of momentum SET to non zero value!");
419                invPt += fDeltaInvP; //shift momentum for tests
420                if(invPt) signedPt = 1.0/invPt;
421                else continue;
422             }
423             Double_t theta = tpcTrack->Theta();
424             Double_t phi = tpcTrack->Phi();
425             
426             Double_t momAng[4] = {invPt,signedPt,theta,phi};
427             fHistInvPtPtThetaPhi->Fill(momAng);
428
429             Double_t pTPC = tpcTrack->GetP();
430             Double_t pESD = esdTrack->GetP();
431             Double_t ptESD  = esdTrack->GetSignedPt();
432         
433             if(esdTrack->GetSign()>0){
434                //compare momenta ESD track and TPC track
435                fHistTPCMomentaPosP->Fill(fabs(pESD),fabs(pTPC));
436                fHistTPCMomentaPosPt->Fill(fabs(ptESD),fabs(signedPt));
437             }
438             else{
439                fHistTPCMomentaNegP->Fill(fabs(pESD),fabs(pTPC));
440                fHistTPCMomentaNegPt->Fill(fabs(ptESD),fabs(signedPt));
441             }
442             count++;
443          }
444          else continue;
445       }
446    
447       else{// ESD tracks
448          if(fabs(esdTrack->Eta())> fEtaAcceptance) continue;
449          Double_t invPt = 0.0;
450          Double_t signedPt = esdTrack->GetSignedPt();
451          if(signedPt){
452             invPt = 1.0/signedPt; 
453
454             fHistPtShift0->Fill(signedPt);
455           
456             if(fShift){Printf("user shift of momentum SET to non zero value!");
457                invPt += fDeltaInvP;//shift momentum for tests
458                if(invPt) signedPt = 1.0/invPt;
459                else continue;
460             }
461             Double_t theta = esdTrack->Theta();
462             Double_t phi = esdTrack->Phi();
463             Double_t momAng[4] = {invPt,signedPt,theta,phi};
464             fHistInvPtPtThetaPhi->Fill(momAng);
465             count++;
466          }
467       }
468    }
469     
470    fHistTrackMultiplicityCuts->Fill(count);
471 }    
472 //______________________________________________________________________________________________________________________
473
474 void AliPerformancePtCalib::Analyse()
475 {
476    // analyse charge/pt spectra in bins of theta and phi. Bins can be set by user
477    
478    THnSparseF *copyTHnSparseTheta = (THnSparseF*)fHistInvPtPtThetaPhi->Clone("copyTHnSparseTheta");
479    copyTHnSparseTheta->GetAxis(3)->SetRangeUser(fMinPhi,fMaxPhi);
480    TH2F *histInvPtTheta = (TH2F*)copyTHnSparseTheta->Projection(2,0);
481       
482    THnSparseF *copyTHnSparsePhi = (THnSparseF*)fHistInvPtPtThetaPhi->Clone("copyTHnSparsePhi");
483    copyTHnSparsePhi->GetAxis(2)->SetRangeUser(fMinTheta,fMaxTheta);
484    TH2F *histInvPtPhi   = (TH2F*)copyTHnSparsePhi->Projection(3,0);
485    
486    AliPerfAnalyzeInvPt *ana = new  AliPerfAnalyzeInvPt("AliPerfAnalyzeInvPt","AliPerfAnalyzeInvPt");
487   
488      
489    TH1::AddDirectory(kFALSE);
490  
491    ana->SetProjBinsTheta(fThetaBins,fNThetaBins);
492    ana->SetProjBinsPhi(fPhiBins,fNPhiBins);
493    ana->SetMakeFitOption(fFitGaus,fExclRange,fRange);
494    if(fDoRebin) ana->SetDoRebin(fRebin);                   
495    TObjArray *aFolderObj = new TObjArray;
496    
497    ana->StartAnalysis(histInvPtTheta,histInvPtPhi,aFolderObj);
498   
499    // export objects to analysis folder
500    fAnalysisFolder = ExportToFolder(aFolderObj);
501
502    // delete only TObjArray
503    if(aFolderObj) delete aFolderObj;
504    if(ana) delete ana;
505    
506 }
507
508 //______________________________________________________________________________________________________________________
509 TFolder* AliPerformancePtCalib::ExportToFolder(TObjArray * array) 
510 {
511    // recreate folder every time and export objects to new one
512    //
513    AliPerformancePtCalib * comp=this;
514    TFolder *folder = comp->GetAnalysisFolder();
515
516    TString name, title;
517    TFolder *newFolder = 0;
518    Int_t i = 0;
519    Int_t size = array->GetSize();
520
521    if(folder) { 
522       // get name and title from old folder
523       name = folder->GetName();  
524       title = folder->GetTitle();  
525
526       // delete old one
527       delete folder;
528
529       // create new one
530       newFolder = CreateFolder(name.Data(),title.Data());
531       newFolder->SetOwner();
532
533       // add objects to folder
534       while(i < size) {
535          newFolder->Add(array->At(i));
536          i++;
537       }
538    }
539
540    return newFolder;
541 }
542
543 //______________________________________________________________________________________________________________________
544 Long64_t AliPerformancePtCalib::Merge(TCollection* const list) 
545 {
546    // Merge list of objects (needed by PROOF)
547
548    if (!list)
549       return 0;
550
551    if (list->IsEmpty())
552       return 1;
553
554    TIterator* iter = list->MakeIterator();
555    TObject* obj = 0;
556
557    // collection of generated histograms
558    Int_t count=0;
559    while((obj = iter->Next()) != 0) 
560       {
561          AliPerformancePtCalib* entry = dynamic_cast<AliPerformancePtCalib*>(obj);
562          if (entry == 0) continue; 
563          fHistInvPtPtThetaPhi->Add(entry->fHistInvPtPtThetaPhi);
564   
565          fHistPtShift0->Add(entry->fHistPtShift0);
566          fHistPrimaryVertexPosX->Add(entry->fHistPrimaryVertexPosX);
567          fHistPrimaryVertexPosY->Add(entry->fHistPrimaryVertexPosY);
568          fHistPrimaryVertexPosZ->Add(entry->fHistPrimaryVertexPosZ);
569          fHistTrackMultiplicity->Add(entry->fHistTrackMultiplicity);
570          fHistTrackMultiplicityCuts->Add(entry->fHistTrackMultiplicityCuts);
571   
572          fHistTPCMomentaPosP->Add(entry->fHistTPCMomentaPosP);
573          fHistTPCMomentaNegP->Add(entry->fHistTPCMomentaNegP);
574          fHistTPCMomentaPosPt->Add(entry->fHistTPCMomentaPosPt);
575          fHistTPCMomentaNegPt->Add(entry->fHistTPCMomentaNegPt);
576          fHistdedxPions->Add(entry->fHistdedxPions);
577          count++;
578       }
579   
580    return count;
581 }
582
583 //______________________________________________________________________________________________________________________
584 TFolder* AliPerformancePtCalib::CreateFolder(TString name,TString title) { 
585    // create folder for analysed histograms
586    //
587    TFolder *folder = 0;
588    folder = new TFolder(name.Data(),title.Data());
589
590    return folder;
591 }
592 //______________________________________________________________________________________________________________________
593 void AliPerformancePtCalib::SetProjBinsPhi(const Double_t *phiBinArray,const Int_t nphBins, const  Double_t minTheta, const  Double_t maxTheta){
594    // set phi bins for Analyse()
595    //set phi bins as array and set number of this array which is equal to number of bins analysed
596    //the last analysed bin will always be the projection from first to last bin in the array
597    if(nphBins){
598       fNPhiBins = nphBins;
599   
600       for(Int_t k = 0;k<fNPhiBins;k++){
601          fPhiBins[k] = phiBinArray[k];
602       }
603       Printf("AliPerformancePtCalib::SetProjBinsPhi: number of bins in phi set to %i",fNPhiBins);
604    }
605    else  Printf("Warning AliPerformancePtCalib::SetProjBinsPhi: number of bins in phi NOT set!!! Default values are taken.");
606
607    if(fabs(minTheta-maxTheta)<0.001){
608       Printf("AliPerformancePtCalib::SetProjBinsPhi: theta range for projection on phi and charge/pt is too small. whole range of theta selected.");
609    }
610    else{
611       Printf("AliPerformancePtCalib::SetProjBinsPhi: theta range for projection on phi and charge/pt is selected by user: %1.3f - %1.3f rad",minTheta,maxTheta);
612       fMinTheta = minTheta;
613       fMaxTheta = maxTheta;
614    }
615 }
616
617 //____________________________________________________________________________________________________________________________________________
618 void AliPerformancePtCalib::SetProjBinsTheta(const Double_t *thetaBinArray, const Int_t nthBins, const Double_t minPhi, const Double_t maxPhi){
619    // set theta bins for Analyse()
620    //set theta bins as array and set number of this array which is equal to number of bins analysed
621    //the last analysed bin will always be the projection from first to last bin in the array
622    if(nthBins){
623       fNThetaBins = nthBins;
624       for(Int_t k = 0;k<fNThetaBins;k++){
625          fThetaBins[k] = thetaBinArray[k];
626       }
627       Printf("AliPerformancePtCalib::SetProjBinsTheta: number of bins in theta set to %i",fNThetaBins);
628    }
629    else  Printf("Warning AliPerformancePtCalib::SetProjBinsTheta: number of bins in theta NOT set!!! Default values are taken.");
630    
631    if(fabs(minPhi-maxPhi)<0.001){
632       Printf("AliPerformancePtCalib::SetProjBinsTheta: phi range for projection on theta and charge/pt is too small. whole range of phi selected.");
633    }
634    else{
635       Printf("AliPerformancePtCalib::SetProjBinsTheta: phi range for projection on theta and charge/pt is selected by user: %1.3f - %1.3f rad",minPhi,maxPhi);
636       fMinPhi = minPhi;
637       fMaxPhi = maxPhi;
638    }
639 }
640 //____________________________________________________________________________________________________________________________________________
641 void AliPerformancePtCalib::SetMakeFitOption(const Bool_t setGausFit, const Double_t exclusionR,const Double_t fitR ){
642    //set the fit options:
643    //for usage of gaussian function instead of polynomial (default) set setGausFit=kTRUE
644    //set the range of rejection of points around 0 via exclusionR
645    //set the fit range around 0 with fitR
646   
647    fRange = fitR;
648    fFitGaus = setGausFit;
649    fExclRange  = exclusionR;
650   
651    if(fFitGaus) Printf("AliPerformancePtCalib::SetMakeFitOption: set MakeGausFit with fit range %2.3f and exclusion range in fabs(1/pt): %2.3f",fRange,fExclRange);
652    else  Printf("AliPerformancePtCalib::SetMakeFitOption: set standard polynomial fit with fit range %2.3f and exclusion range in fabs(1/pt): %2.3f",fRange,fExclRange);
653  
654 }