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