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