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