]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/TPC/AliPerformanceEff.cxx
Updates in D+ histos and ntuples (Renu, Francesco, Elena)
[u/mrichter/AliRoot.git] / PWGPP / TPC / AliPerformanceEff.cxx
1 //------------------------------------------------------------------------------\r
2 // Implementation of AliPerformanceEff class. It keeps information from \r
3 // comparison of reconstructed and MC particle tracks. In addtion, \r
4 // it keeps selection cuts used during comparison. The comparison \r
5 // information is stored in the ROOT histograms. Analysis of these \r
6 // histograms can be done by using Analyse() class function. The result of \r
7 // the analysis (histograms/graphs) are stored in the folder which is \r
8 // a data member of AliPerformanceEff.\r
9 // \r
10 // Author: J.Otwinowski 04/02/2008 \r
11 //------------------------------------------------------------------------------\r
12 \r
13 /*\r
14  \r
15   // after running comparison task, read the file, and get component\r
16   gROOT->LoadMacro("$ALICE_ROOT/PWGPP/Macros/LoadMyLibs.C");\r
17   LoadMyLibs();\r
18   TFile f("Output.root");\r
19   AliPerformanceEff * compObj = (AliPerformanceEff*)coutput->FindObject("AliPerformanceEff");\r
20 \r
21   // Analyse comparison data\r
22   compObj->Analyse();\r
23 \r
24   // the output histograms/graphs will be stored in the folder "folderEff" \r
25   compObj->GetAnalysisFolder()->ls("*");\r
26 \r
27   // user can save whole comparison object (or only folder with anlysed histograms) \r
28   // in the seperate output file (e.g.)\r
29   TFile fout("Analysed_Eff.root","recreate");\r
30   compObj->Write(); // compObj->GetAnalysisFolder()->Write();\r
31   fout.Close();\r
32 \r
33 */\r
34 \r
35 #include <TAxis.h>\r
36 #include <TH1D.h>\r
37 #include "THnSparse.h"\r
38 \r
39 // \r
40 #include "AliESDtrack.h"\r
41 #include "AliRecInfoCuts.h" \r
42 #include "AliMCInfoCuts.h" \r
43 #include "AliLog.h" \r
44 #include "AliESDVertex.h" \r
45 #include "AliExternalTrackParam.h" \r
46 #include "AliTracker.h" \r
47 #include "AliESDEvent.h" \r
48 #include "AliMCEvent.h" \r
49 #include "AliMCParticle.h" \r
50 #include "AliHeader.h" \r
51 #include "AliGenEventHeader.h" \r
52 #include "AliStack.h" \r
53 #include "AliPerformanceEff.h" \r
54 \r
55 using namespace std;\r
56 \r
57 ClassImp(AliPerformanceEff)\r
58 \r
59 //_____________________________________________________________________________\r
60 AliPerformanceEff::AliPerformanceEff():\r
61   AliPerformanceObject("AliPerformanceEff"),\r
62 \r
63   // histograms\r
64   fEffHisto(0),\r
65   fEffSecHisto(0),\r
66 \r
67   // Cuts \r
68   fCutsRC(0), \r
69   fCutsMC(0),\r
70 \r
71   // histogram folder \r
72   fAnalysisFolder(0)\r
73 {\r
74   // default consttructor       \r
75   Init();\r
76 }\r
77 \r
78 //_____________________________________________________________________________\r
79 AliPerformanceEff::AliPerformanceEff(Char_t* name="AliPerformanceEff",Char_t*title="AliPerformanceEff",Int_t analysisMode=0, Bool_t hptGenerator=kFALSE):\r
80   AliPerformanceObject(name,title),\r
81 \r
82   // histograms\r
83   fEffHisto(0),\r
84   fEffSecHisto(0),\r
85 \r
86   // Cuts \r
87   fCutsRC(0), \r
88   fCutsMC(0),\r
89 \r
90   // histogram folder \r
91   fAnalysisFolder(0)\r
92 {\r
93   // named constructor\r
94   //\r
95   SetAnalysisMode(analysisMode);\r
96   SetHptGenerator(hptGenerator);\r
97 \r
98   Init();\r
99 }\r
100 \r
101 \r
102 //_____________________________________________________________________________\r
103 AliPerformanceEff::~AliPerformanceEff()\r
104 {\r
105 // destructor\r
106 \r
107   if(fEffHisto)  delete  fEffHisto; fEffHisto=0;\r
108   if(fEffSecHisto)  delete  fEffSecHisto; fEffSecHisto=0;\r
109   if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;\r
110 }\r
111 \r
112 //_____________________________________________________________________________\r
113 void AliPerformanceEff::Init()\r
114 {\r
115   // Init histograms\r
116   //\r
117 \r
118   // set pt bins\r
119   Int_t nPtBins = 50;\r
120   Double_t ptMin = 1.e-2, ptMax = 20.;\r
121 \r
122   Double_t *binsPt = 0;\r
123 \r
124   if (IsHptGenerator())  { \r
125         ptMax = 100.;\r
126   } \r
127    binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);\r
128 \r
129   /*\r
130   Int_t nPtBins = 31;\r
131   Double_t binsPt[32] = {0.,0.05,0.1,0.15,0.2,0.25,0.3,0.35,0.4,0.45,0.5,0.55,0.6,0.7,0.8,0.9,1.0,1.2,1.4,1.6,1.8,2.0,2.25,2.5,2.75,3.,3.5,4.,5.,6.,8.,10.};\r
132   Double_t ptMin = 0., ptMax = 10.;\r
133 \r
134   if(IsHptGenerator() == kTRUE) {\r
135     nPtBins = 100;\r
136     ptMin = 0.; ptMax = 100.;\r
137   }\r
138   */\r
139 \r
140   //mceta:mcphi:mcpt:pid:recStatus:findable:charge\r
141   Int_t binsEffHisto[7]={30,144,nPtBins,5,2,2,3};\r
142   Double_t minEffHisto[7]={-1.5,0.,ptMin,0.,0.,0.,-1.5};\r
143   Double_t maxEffHisto[7]={ 1.5,2.*TMath::Pi(), ptMax,5.,2.,2.,1.5};\r
144 \r
145   fEffHisto = new THnSparseF("fEffHisto","mceta:mcphi:mcpt:pid:recStatus:findable:charge",7,binsEffHisto,minEffHisto,maxEffHisto);\r
146   fEffHisto->SetBinEdges(2,binsPt);\r
147 \r
148   fEffHisto->GetAxis(0)->SetTitle("#eta_{mc}");\r
149   fEffHisto->GetAxis(1)->SetTitle("#phi_{mc} (rad)");\r
150   fEffHisto->GetAxis(2)->SetTitle("p_{Tmc} (GeV/c)");\r
151   fEffHisto->GetAxis(3)->SetTitle("pid");\r
152   fEffHisto->GetAxis(4)->SetTitle("recStatus");\r
153   fEffHisto->GetAxis(5)->SetTitle("findable");\r
154   fEffHisto->GetAxis(6)->SetTitle("charge");\r
155   fEffHisto->Sumw2();\r
156 \r
157   //mceta:mcphi:mcpt:pid:recStatus:findable:mcR:mother_phi:mother_eta:charge\r
158   Int_t binsEffSecHisto[10]={30,60,nPtBins,5,2,2,100,60,30,3};\r
159   Double_t minEffSecHisto[10]={-1.5,0.,ptMin,0.,0.,0.,0.,0.,-1.5,-1.5};\r
160   Double_t maxEffSecHisto[10]={ 1.5,2.*TMath::Pi(), ptMax,5.,2.,2.,200,2.*TMath::Pi(),1.5,1.5};\r
161 \r
162   fEffSecHisto = new THnSparseF("fEffSecHisto","mceta:mcphi:mcpt:pid:recStatus:findable:mcR:mother_phi:mother_eta:charge",10,binsEffSecHisto,minEffSecHisto,maxEffSecHisto);\r
163   fEffSecHisto->SetBinEdges(2,binsPt);\r
164 \r
165   fEffSecHisto->GetAxis(0)->SetTitle("#eta_{mc}");\r
166   fEffSecHisto->GetAxis(1)->SetTitle("#phi_{mc} (rad)");\r
167   fEffSecHisto->GetAxis(2)->SetTitle("p_{Tmc} (GeV/c)");\r
168   fEffSecHisto->GetAxis(3)->SetTitle("pid");\r
169   fEffSecHisto->GetAxis(4)->SetTitle("recStatus");\r
170   fEffSecHisto->GetAxis(5)->SetTitle("findable");\r
171   fEffSecHisto->GetAxis(6)->SetTitle("mcR (cm)");\r
172   fEffSecHisto->GetAxis(7)->SetTitle("mother_phi (rad)");\r
173   fEffSecHisto->GetAxis(8)->SetTitle("mother_eta");\r
174   fEffSecHisto->GetAxis(9)->SetTitle("charge");\r
175   fEffSecHisto->Sumw2();\r
176 \r
177   // init cuts\r
178   if(!fCutsMC) \r
179     AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");\r
180   if(!fCutsRC) \r
181     AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");\r
182 \r
183   // init folder\r
184   fAnalysisFolder = CreateFolder("folderEff","Analysis Efficiency Folder");\r
185 }\r
186 \r
187 //_____________________________________________________________________________\r
188 void AliPerformanceEff::ProcessTPC(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent)\r
189 {\r
190   // Fill TPC only efficiency comparison information \r
191   if(!esdEvent) return;\r
192   if(!mcEvent) return;\r
193 \r
194   AliStack *stack = mcEvent->Stack();\r
195   if (!stack) {\r
196     AliDebug(AliLog::kError, "Stack not available");\r
197     return;\r
198   }\r
199 \r
200   Int_t *labelsRec =  NULL;\r
201   labelsRec =  new Int_t[esdEvent->GetNumberOfTracks()];\r
202   if(!labelsRec) \r
203   {\r
204      Printf("Cannot create labelsRec");\r
205      return;\r
206   }\r
207   for(Int_t i=0;i<esdEvent->GetNumberOfTracks();i++) { labelsRec[i] = 0; }\r
208 \r
209   Int_t *labelsAllRec =  NULL;\r
210   labelsAllRec =  new Int_t[esdEvent->GetNumberOfTracks()];\r
211   if(!labelsAllRec) { \r
212      delete  [] labelsRec;\r
213      Printf("Cannot create labelsAllRec");\r
214      return;\r
215   }\r
216   for(Int_t i=0;i<esdEvent->GetNumberOfTracks();i++) { labelsAllRec[i] = 0; }\r
217 \r
218   // loop over rec. tracks\r
219   AliESDtrack *track=0;\r
220   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) \r
221   { \r
222     track = esdEvent->GetTrack(iTrack);\r
223     if(!track) continue;\r
224     if(track->Charge()==0) continue;\r
225 \r
226     // if not fUseKinkDaughters don't use tracks with kink index > 0\r
227     if(!fUseKinkDaughters && track->GetKinkIndex(0) > 0) continue;\r
228 \r
229     Int_t label = TMath::Abs(track->GetLabel()); \r
230     labelsAllRec[iTrack]=label;\r
231 \r
232     // TPC only\r
233     if(IsRecTPC(track) != 0) { \r
234       labelsRec[iTrack]=label;\r
235     }\r
236 \r
237   }\r
238 \r
239   // \r
240   // MC histograms for efficiency studies\r
241   //\r
242   Int_t nPart  = stack->GetNtrack();\r
243   //Int_t nPart  = stack->GetNprimary();\r
244   for (Int_t iMc = 0; iMc < nPart; ++iMc) \r
245   {\r
246     TParticle* particle = stack->Particle(iMc);\r
247     if (!particle) continue;\r
248     if (!particle->GetPDG()) continue; \r
249     if (particle->GetPDG()->Charge() == 0.0) continue;\r
250       \r
251     // physical primary\r
252     Bool_t prim = stack->IsPhysicalPrimary(iMc);\r
253     if(!prim) continue;\r
254 \r
255     // --- check for double filling in stack\r
256     // use only particles with no daughters in the list of primaries\r
257     Int_t nDaughters =  particle->GetNDaughters();\r
258     \r
259     for( Int_t iDaught=0; iDaught < particle->GetNDaughters(); iDaught++ ) {\r
260       if( particle->GetDaughter(iDaught) < stack->GetNprimary() )\r
261         nDaughters++;\r
262     }\r
263     if( nDaughters > 0 ) \r
264       continue;\r
265     // --- check for double filling in stack\r
266 \r
267     Bool_t findable = kFALSE;\r
268     for(Int_t iRec=0; iRec<esdEvent->GetNumberOfTracks(); ++iRec) \r
269     {\r
270       // check findable\r
271       if(iMc == labelsAllRec[iRec]) \r
272       {\r
273         findable = IsFindable(mcEvent,iMc);\r
274         break;\r
275       }\r
276     }  \r
277 \r
278     Bool_t recStatus = kFALSE;\r
279     for(Int_t iRec=0; iRec<esdEvent->GetNumberOfTracks(); ++iRec) \r
280     {\r
281       // check reconstructed\r
282       if(iMc == labelsRec[iRec]) \r
283       {\r
284         recStatus = kTRUE;\r
285         break;\r
286       }\r
287     }\r
288 \r
289     // Only 5 charged particle species (e,mu,pi,K,p)\r
290     if (fCutsMC->IsPdgParticle(TMath::Abs(particle->GetPdgCode())) == kFALSE) continue; \r
291 \r
292     // transform Pdg to Pid\r
293     Int_t pid = TransformToPID(particle);\r
294 \r
295     Float_t mceta =  particle->Eta();\r
296     Float_t mcphi =  particle->Phi();\r
297     if(mcphi<0) mcphi += 2.*TMath::Pi();\r
298     Float_t mcpt = particle->Pt();\r
299     Float_t charge = 0.;\r
300     if (particle->GetPDG()->Charge() < 0)  charge = -1.;    \r
301     else if (particle->GetPDG()->Charge() > 0)  charge = 1.;\r
302 \r
303     // Fill histograms\r
304     Double_t vEffHisto[7] = {mceta, mcphi, mcpt, pid, recStatus, findable, charge}; \r
305     fEffHisto->Fill(vEffHisto);\r
306   }\r
307 \r
308   if(labelsRec) delete [] labelsRec; labelsRec = 0;\r
309   if(labelsAllRec) delete [] labelsAllRec; labelsAllRec = 0;\r
310 }\r
311 \r
312 //_____________________________________________________________________________\r
313 void AliPerformanceEff::ProcessTPCSec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent)\r
314 {\r
315   // Fill TPC only efficiency comparison information for secondaries\r
316 \r
317   if(!esdEvent) return;\r
318 \r
319   Int_t *labelsRecSec =  NULL;\r
320   labelsRecSec =  new Int_t[esdEvent->GetNumberOfTracks()];\r
321   if(!labelsRecSec) \r
322   {\r
323      Printf("Cannot create labelsRecSec");\r
324      return;\r
325   }\r
326   for(Int_t i=0;i<esdEvent->GetNumberOfTracks();i++) { labelsRecSec[i] = 0; }\r
327 \r
328   Int_t *labelsAllRecSec =  NULL;\r
329   labelsAllRecSec =  new Int_t[esdEvent->GetNumberOfTracks()];\r
330   if(!labelsAllRecSec) { \r
331      delete [] labelsRecSec;\r
332      Printf("Cannot create labelsAllRecSec");\r
333      return;\r
334   }\r
335   for(Int_t i=0;i<esdEvent->GetNumberOfTracks();i++) { labelsAllRecSec[i] = 0; }\r
336 \r
337   // loop over rec. tracks\r
338   AliESDtrack *track=0;\r
339   Int_t multAll=0, multRec=0;\r
340   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) \r
341   { \r
342     track = esdEvent->GetTrack(iTrack);\r
343     if(!track) continue;\r
344     if(track->Charge()==0) continue;\r
345 \r
346     // if not fUseKinkDaughters don't use tracks with kink index > 0\r
347     if(!fUseKinkDaughters && track->GetKinkIndex(0) > 0) continue;\r
348 \r
349     Int_t label = TMath::Abs(track->GetLabel()); \r
350     labelsAllRecSec[multAll]=label;\r
351     multAll++;\r
352 \r
353     // TPC only\r
354     if(IsRecTPC(track) != 0) {\r
355       labelsRecSec[multRec]=label;\r
356       multRec++;\r
357     }\r
358   }\r
359 \r
360   // \r
361   // MC histograms for efficiency studies\r
362   //\r
363   if(mcEvent)  { \r
364  \r
365   AliStack *stack = mcEvent->Stack();\r
366   if (stack) {\r
367 \r
368   Int_t nPart  = stack->GetNtrack();\r
369   //Int_t nPart  = stack->GetNprimary();\r
370   for (Int_t iMc = 0; iMc < nPart; ++iMc) \r
371   {\r
372     TParticle* particle = stack->Particle(iMc);\r
373     if (!particle) continue;\r
374     if (!particle->GetPDG()) continue; \r
375     if (particle->GetPDG()->Charge() == 0.0) continue;\r
376       \r
377     // physical primary\r
378     Bool_t prim = stack->IsPhysicalPrimary(iMc);\r
379 \r
380     // only secondaries which can be reconstructed at TPC\r
381     if(prim) continue;\r
382 \r
383     //Float_t radius = TMath::Sqrt(particle->Vx()*particle->Vx()+particle->Vy()*particle->Vy()+particle->Vz()*particle->Vz());\r
384     //if(radius > fCutsMC->GetMaxR()) continue;\r
385 \r
386     // only secondary electrons from gamma conversion\r
387     //if( TMath::Abs(particle->GetPdgCode())!=fCutsMC->GetEM() ||   particle->GetUniqueID() != 5) continue;\r
388 \r
389     Bool_t findable = kFALSE;\r
390     for(Int_t iRec=0; iRec<multAll; ++iRec) \r
391     {\r
392       // check findable\r
393       if(iMc == labelsAllRecSec[iRec]) \r
394       {\r
395         findable = IsFindable(mcEvent,iMc);\r
396         break;\r
397       }\r
398     }  \r
399 \r
400     Bool_t recStatus = kFALSE;\r
401     for(Int_t iRec=0; iRec<multRec; ++iRec) \r
402     {\r
403       // check reconstructed\r
404       if(iMc == labelsRecSec[iRec]) \r
405       {\r
406         recStatus = kTRUE;\r
407         break;\r
408       }\r
409     }\r
410 \r
411     // Only 5 charged particle species (e,mu,pi,K,p)\r
412     if (fCutsMC->IsPdgParticle(TMath::Abs(particle->GetPdgCode())) == kFALSE) continue; \r
413 \r
414     // transform Pdg to Pid\r
415     Int_t pid = TransformToPID(particle);\r
416 \r
417     Float_t mceta =  particle->Eta();\r
418     Float_t mcphi =  particle->Phi();\r
419     if(mcphi<0) mcphi += 2.*TMath::Pi();\r
420     Float_t mcpt = particle->Pt();\r
421     Float_t mcR = particle->R();\r
422 \r
423     // get info about mother\r
424     Int_t motherLabel = particle->GetMother(0);\r
425     if(motherLabel < 0) continue;\r
426     TParticle *mother = stack->Particle(motherLabel);\r
427     if(!mother) continue; \r
428 \r
429     Float_t mother_eta = mother->Eta();\r
430     Float_t mother_phi = mother->Phi();\r
431     if(mother_phi<0) mother_phi += 2.*TMath::Pi();\r
432 \r
433     Float_t charge = 0.;\r
434     if (particle->GetPDG()->Charge() < 0)  charge = -1.;    \r
435     else if (particle->GetPDG()->Charge() > 0)  charge = 1.;\r
436 \r
437     // Fill histograms\r
438     Double_t vEffSecHisto[10] = { mceta, mcphi, mcpt, pid, recStatus, findable, mcR, mother_phi, mother_eta, charge }; \r
439     fEffSecHisto->Fill(vEffSecHisto);\r
440   }\r
441   }\r
442   }\r
443 \r
444   if(labelsRecSec) delete [] labelsRecSec; labelsRecSec = 0;\r
445   if(labelsAllRecSec) delete [] labelsAllRecSec; labelsAllRecSec = 0;\r
446 }\r
447 \r
448 \r
449 \r
450 \r
451 //_____________________________________________________________________________\r
452 void AliPerformanceEff::ProcessTPCITS(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent)\r
453 {\r
454   // Fill efficiency comparison information\r
455 \r
456   if(!esdEvent) return;\r
457   if(!mcEvent) return;\r
458 \r
459   AliStack *stack = mcEvent->Stack();\r
460   if (!stack) {\r
461     AliDebug(AliLog::kError, "Stack not available");\r
462     return;\r
463   }\r
464 \r
465   Int_t *labelsRecTPCITS =  NULL;\r
466   labelsRecTPCITS =  new Int_t[esdEvent->GetNumberOfTracks()];\r
467   if(!labelsRecTPCITS) \r
468   {\r
469      Printf("Cannot create labelsRecTPCITS");\r
470      return;\r
471   }\r
472   for(Int_t i=0;i<esdEvent->GetNumberOfTracks();i++) { labelsRecTPCITS[i] = 0; }\r
473 \r
474   Int_t *labelsAllRecTPCITS =  NULL;\r
475   labelsAllRecTPCITS =  new Int_t[esdEvent->GetNumberOfTracks()];\r
476   if(!labelsAllRecTPCITS) { \r
477      delete [] labelsRecTPCITS;\r
478      Printf("Cannot create labelsAllRecTPCITS");\r
479      return;\r
480   }\r
481   for(Int_t i=0;i<esdEvent->GetNumberOfTracks();i++) { labelsAllRecTPCITS[i] = 0; }\r
482 \r
483   // loop over rec. tracks\r
484   AliESDtrack *track=0;\r
485   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) \r
486   { \r
487     track = esdEvent->GetTrack(iTrack);\r
488     if(!track) continue;\r
489     if(track->Charge()==0) continue;\r
490 \r
491     // if not fUseKinkDaughters don't use tracks with kink index > 0\r
492     if(!fUseKinkDaughters && track->GetKinkIndex(0) > 0) continue;\r
493 \r
494     Int_t label = TMath::Abs(track->GetLabel()); \r
495     labelsAllRecTPCITS[iTrack]=label;\r
496 \r
497     // iTPC+ITS\r
498     if(IsRecTPCITS(track) != 0) \r
499       labelsRecTPCITS[iTrack]=label;\r
500   }\r
501 \r
502   // \r
503   // MC histograms for efficiency studies\r
504   //\r
505   //Int_t nPart  = stack->GetNtrack();\r
506   Int_t nPart  = stack->GetNprimary();\r
507   for (Int_t iMc = 0; iMc < nPart; ++iMc) \r
508   {\r
509     TParticle* particle = stack->Particle(iMc);\r
510     if (!particle) continue;\r
511     if (!particle->GetPDG()) continue; \r
512     if (particle->GetPDG()->Charge() == 0.0) continue;\r
513       \r
514     // physical primary\r
515     Bool_t prim = stack->IsPhysicalPrimary(iMc);\r
516     if(!prim) continue;\r
517 \r
518     Bool_t findable = kFALSE;\r
519     for(Int_t iRec=0; iRec<esdEvent->GetNumberOfTracks(); ++iRec) \r
520     {\r
521       // check findable\r
522       if(iMc == labelsAllRecTPCITS[iRec]) \r
523       {\r
524         findable = IsFindable(mcEvent,iMc);\r
525         break;\r
526       }\r
527     }  \r
528 \r
529     Bool_t recStatus = kFALSE;\r
530     for(Int_t iRec=0; iRec<esdEvent->GetNumberOfTracks(); ++iRec) \r
531     {\r
532       // check reconstructed\r
533       if(iMc == labelsRecTPCITS[iRec]) \r
534       {\r
535         recStatus = kTRUE;\r
536         break;\r
537       }\r
538     }\r
539 \r
540     // Only 5 charged particle species (e,mu,pi,K,p)\r
541     if (fCutsMC->IsPdgParticle(TMath::Abs(particle->GetPdgCode())) == kFALSE) continue; \r
542 \r
543     // transform Pdg to Pid\r
544     Int_t pid = TransformToPID(particle);\r
545 \r
546     Float_t mceta =  particle->Eta();\r
547     Float_t mcphi =  particle->Phi();\r
548     if(mcphi<0) mcphi += 2.*TMath::Pi();\r
549     Float_t mcpt = particle->Pt();\r
550 \r
551     Float_t charge = 0.;\r
552     if (particle->GetPDG()->Charge() < 0)  charge = -1.;    \r
553     else if (particle->GetPDG()->Charge() > 0)  charge = 1.;\r
554     \r
555     // Fill histograms\r
556     Double_t vEffHisto[7] = { mceta, mcphi, mcpt, pid, recStatus, findable, charge}; \r
557     fEffHisto->Fill(vEffHisto);\r
558   }\r
559 \r
560   if(labelsRecTPCITS) delete [] labelsRecTPCITS; labelsRecTPCITS = 0;\r
561   if(labelsAllRecTPCITS) delete [] labelsAllRecTPCITS; labelsAllRecTPCITS = 0;\r
562 }\r
563 \r
564 //_____________________________________________________________________________\r
565 void AliPerformanceEff::ProcessConstrained(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent)\r
566 {\r
567   // Process comparison information \r
568   if(!esdEvent) return;\r
569   if(!mcEvent) return;\r
570 \r
571   AliStack *stack = mcEvent->Stack();\r
572   if (!stack) {\r
573     AliDebug(AliLog::kError, "Stack not available");\r
574     return;\r
575   }\r
576 \r
577   Int_t *labelsRecConstrained =  NULL;\r
578   labelsRecConstrained =  new Int_t[esdEvent->GetNumberOfTracks()];\r
579   if(!labelsRecConstrained) \r
580   {\r
581      Printf("Cannot create labelsRecConstrained");\r
582      return;\r
583   }\r
584   for(Int_t i=0;i<esdEvent->GetNumberOfTracks();i++) { labelsRecConstrained[i] = 0; }\r
585 \r
586   Int_t *labelsAllRecConstrained =  NULL;\r
587   labelsAllRecConstrained =  new Int_t[esdEvent->GetNumberOfTracks()];\r
588   if(!labelsAllRecConstrained) { \r
589      delete [] labelsRecConstrained;\r
590      Printf("Cannot create labelsAllRecConstrained");\r
591      return;\r
592   }\r
593   for(Int_t i=0;i<esdEvent->GetNumberOfTracks();i++) { labelsAllRecConstrained[i] = 0; }\r
594 \r
595   // loop over rec. tracks\r
596   AliESDtrack *track=0;\r
597   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) \r
598   { \r
599     track = esdEvent->GetTrack(iTrack);\r
600     if(!track) continue;\r
601     if(track->Charge()==0) continue;\r
602 \r
603     // if not fUseKinkDaughters don't use tracks with kink index > 0\r
604     if(!fUseKinkDaughters && track->GetKinkIndex(0) > 0) continue;\r
605 \r
606     Int_t label = TMath::Abs(track->GetLabel()); \r
607     labelsAllRecConstrained[iTrack]=label;\r
608 \r
609     // Constrained\r
610     if(IsRecConstrained(track) != 0) \r
611       labelsRecConstrained[iTrack]=label;\r
612 \r
613   }\r
614 \r
615   // \r
616   // MC histograms for efficiency studies\r
617   //\r
618  \r
619 \r
620   //Int_t nPart  = stack->GetNtrack();\r
621   Int_t nPart  = stack->GetNprimary();\r
622   for (Int_t iMc = 0; iMc < nPart; ++iMc) \r
623   {\r
624     TParticle* particle = stack->Particle(iMc);\r
625     if (!particle) continue;\r
626     if (!particle->GetPDG()) continue; \r
627     if (particle->GetPDG()->Charge() == 0.0) continue;\r
628       \r
629     // physical primary\r
630     Bool_t prim = stack->IsPhysicalPrimary(iMc);\r
631     if(!prim) continue;\r
632 \r
633     Bool_t findable = kFALSE;\r
634     for(Int_t iRec=0; iRec<esdEvent->GetNumberOfTracks(); ++iRec) \r
635     {\r
636       // check findable\r
637       if(iMc == labelsAllRecConstrained[iRec]) \r
638       {\r
639         findable = IsFindable(mcEvent,iMc);\r
640         break;\r
641       }\r
642     }  \r
643 \r
644     Bool_t recStatus = kFALSE;\r
645     for(Int_t iRec=0; iRec<esdEvent->GetNumberOfTracks(); ++iRec) \r
646     {\r
647       // check reconstructed\r
648       if(iMc == labelsRecConstrained[iRec]) \r
649       {\r
650         recStatus = kTRUE;\r
651         break;\r
652       }\r
653     }\r
654 \r
655     // Only 5 charged particle species (e,mu,pi,K,p)\r
656     if (fCutsMC->IsPdgParticle(TMath::Abs(particle->GetPdgCode())) == kFALSE) continue; \r
657 \r
658     // transform Pdg to Pid\r
659     Int_t pid = TransformToPID(particle);\r
660 \r
661     Float_t mceta =  particle->Eta();\r
662     Float_t mcphi =  particle->Phi();\r
663     if(mcphi<0) mcphi += 2.*TMath::Pi();\r
664     Float_t mcpt = particle->Pt();\r
665 \r
666     Float_t charge = 0.;\r
667     if (particle->GetPDG()->Charge() < 0)  charge = -1.;    \r
668     else if (particle->GetPDG()->Charge() > 0)  charge = 1.;\r
669 \r
670     // Fill histograms\r
671     Double_t vEffHisto[7] = { mceta, mcphi, mcpt, pid, recStatus, findable, charge}; \r
672     fEffHisto->Fill(vEffHisto);\r
673   }\r
674 \r
675   if(labelsRecConstrained) delete [] labelsRecConstrained; labelsRecConstrained = 0;\r
676   if(labelsAllRecConstrained) delete [] labelsAllRecConstrained; labelsAllRecConstrained = 0;\r
677 }\r
678 \r
679 //_____________________________________________________________________________\r
680 void AliPerformanceEff::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend *const esdFriend, const Bool_t bUseMC, const Bool_t bUseESDfriend)\r
681 {\r
682   // Process comparison information \r
683   //\r
684   if(!esdEvent) \r
685   {\r
686     Error("Exec","esdEvent not available");\r
687     return;\r
688   }\r
689   AliHeader* header = 0;\r
690   AliGenEventHeader* genHeader = 0;\r
691   AliStack* stack = 0;\r
692   TArrayF vtxMC(3);\r
693   \r
694   if(bUseMC)\r
695   {\r
696     if(!mcEvent) {\r
697       Error("Exec","mcEvent not available");\r
698       return;\r
699     }\r
700     // get MC event header\r
701     header = mcEvent->Header();\r
702     if (!header) {\r
703       Error("Exec","Header not available");\r
704       return;\r
705     }\r
706     // MC particle stack\r
707     stack = mcEvent->Stack();\r
708     if (!stack) {\r
709       Error("Exec","Stack not available");\r
710       return;\r
711     }\r
712     // get MC vertex\r
713     genHeader = header->GenEventHeader();\r
714     if (!genHeader) {\r
715       Error("Exec","Could not retrieve genHeader from Header");\r
716       return;\r
717     }\r
718     genHeader->PrimaryVertex(vtxMC);\r
719   } \r
720   else {\r
721     Error("Exec","MC information required!");\r
722     return;\r
723   } \r
724 \r
725   // use ESD friends\r
726   if(bUseESDfriend) {\r
727     if(!esdFriend) {\r
728       Error("Exec","esdFriend not available");\r
729       return;\r
730     }\r
731   }\r
732 \r
733   //\r
734   //  Process events\r
735   //\r
736   if(GetAnalysisMode() == 0) ProcessTPC(mcEvent,esdEvent);\r
737   else if(GetAnalysisMode() == 1) ProcessTPCITS(mcEvent,esdEvent);\r
738   else if(GetAnalysisMode() == 2) ProcessConstrained(mcEvent,esdEvent);\r
739   else if(GetAnalysisMode() == 5) ProcessTPCSec(mcEvent,esdEvent);\r
740   else {\r
741     printf("ERROR: AnalysisMode %d \n",fAnalysisMode);\r
742     return;\r
743   }\r
744 }\r
745 \r
746 //_____________________________________________________________________________\r
747 Int_t AliPerformanceEff::TransformToPID(TParticle *particle) \r
748 {\r
749 // transform Pdg to Pid\r
750 // Pdg convension is different for hadrons and leptons \r
751 // (e.g. K+/K- = 321/-321; e+/e- = -11/11 ) \r
752 \r
753   Int_t pid = -1;\r
754   if( TMath::Abs(particle->GetPdgCode())==fCutsMC->GetEM() ) pid = 0; \r
755   if( TMath::Abs(particle->GetPdgCode())==fCutsMC->GetMuM() ) pid = 1; \r
756   if( TMath::Abs(particle->GetPdgCode())==fCutsMC->GetPiP() ) pid = 2; \r
757   if( TMath::Abs(particle->GetPdgCode())==fCutsMC->GetKP() ) pid = 3; \r
758   if( TMath::Abs(particle->GetPdgCode())==fCutsMC->GetProt() ) pid = 4; \r
759 \r
760 return pid;\r
761 }\r
762 \r
763 //_____________________________________________________________________________\r
764 Bool_t AliPerformanceEff::IsFindable(const AliMCEvent *mcEvent, Int_t label) \r
765 {\r
766 //\r
767 // Findfindable tracks\r
768 //\r
769 if(!mcEvent) return kFALSE;\r
770 \r
771   AliMCParticle *mcParticle = (AliMCParticle*) mcEvent->GetTrack(label);\r
772   if(!mcParticle) return kFALSE;\r
773 \r
774   Int_t counter; \r
775   Float_t tpcTrackLength = mcParticle->GetTPCTrackLength(AliTracker::GetBz(),0.05,counter,3.0); \r
776   //printf("tpcTrackLength %f \n", tpcTrackLength);\r
777 \r
778 return (tpcTrackLength>fCutsMC->GetMinTrackLength());    \r
779 }\r
780 \r
781 //_____________________________________________________________________________\r
782 Bool_t AliPerformanceEff::IsRecTPC(AliESDtrack *esdTrack) \r
783 {\r
784 //\r
785 // Check whether track is reconstructed in TPC\r
786 //\r
787 if(!esdTrack) return kFALSE;\r
788 \r
789   const AliExternalTrackParam *track = esdTrack->GetTPCInnerParam();\r
790   if(!track) return kFALSE;\r
791 \r
792   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z\r
793   esdTrack->GetImpactParametersTPC(dca,cov);\r
794 \r
795   Bool_t recStatus = kFALSE;\r
796   if(esdTrack->GetTPCNcls()>fCutsRC->GetMinNClustersTPC()) recStatus = kTRUE; \r
797 \r
798   /*\r
799   if( TMath::Abs(dca[0])<fCutsRC->GetMaxDCAToVertexXY() && \r
800       TMath::Abs(dca[1])<fCutsRC->GetMaxDCAToVertexZ())\r
801   {\r
802     recStatus = kTRUE;\r
803   }\r
804   */\r
805 \r
806 return recStatus;\r
807 }\r
808 \r
809 //_____________________________________________________________________________\r
810 Bool_t AliPerformanceEff::IsRecTPCITS(AliESDtrack *esdTrack) \r
811 {\r
812 //\r
813 // Check whether track is reconstructed in TPCITS\r
814 //\r
815 if(!esdTrack) return kFALSE;\r
816 \r
817   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z\r
818   esdTrack->GetImpactParameters(dca,cov);\r
819 \r
820   Bool_t recStatus = kFALSE;\r
821 \r
822   if ((esdTrack->GetStatus()&AliESDtrack::kTPCrefit)==0) return kFALSE; // TPC refit\r
823   if ((esdTrack->GetStatus()&AliESDtrack::kITSrefit)==0) return kFALSE; // ITS refit\r
824   if (esdTrack->GetITSclusters(0)<fCutsRC->GetMinNClustersITS()) return kFALSE;  // min. nb. ITS clusters\r
825   //if ((esdTrack->GetStatus()&AliESDtrack::kITSrefit)==0) return kFALSE; // ITS refit\r
826   //Int_t clusterITS[200];\r
827   //if(esdTrack->GetITSclusters(clusterITS)<fCutsRC->GetMinNClustersITS()) return kFALSE;  // min. nb. ITS clusters\r
828 \r
829   recStatus = kTRUE;\r
830   /*\r
831   if(TMath::Abs(dca[0])<fCutsRC->GetMaxDCAToVertexXY() && \r
832      TMath::Abs(dca[1])<fCutsRC->GetMaxDCAToVertexZ())\r
833   {\r
834     recStatus = kTRUE;\r
835   }\r
836   */\r
837 \r
838 return recStatus;\r
839 }\r
840 \r
841 //_____________________________________________________________________________\r
842 Bool_t AliPerformanceEff::IsRecConstrained(AliESDtrack *esdTrack) \r
843 {\r
844 //\r
845 // Check whether track is reconstructed in IsRecConstrained\r
846 //\r
847   if(!esdTrack) return kFALSE;\r
848 \r
849   const AliExternalTrackParam * track = esdTrack->GetConstrainedParam();\r
850   if(!track) return kFALSE;\r
851 \r
852   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z\r
853   esdTrack->GetImpactParameters(dca,cov);\r
854   //Int_t label = TMath::Abs(esdTrack->GetLabel()); \r
855 \r
856   Bool_t recStatus = kFALSE;\r
857 \r
858   if ((esdTrack->GetStatus()&AliESDtrack::kTPCrefit)==0) return kFALSE; // TPC refit\r
859   if (esdTrack->GetTPCNcls()<fCutsRC->GetMinNClustersTPC()) return kFALSE; // min. nb. TPC clusters\r
860   Int_t clusterITS[200];\r
861   if(esdTrack->GetITSclusters(clusterITS)<fCutsRC->GetMinNClustersITS()) return kFALSE;  // min. nb. ITS clusters\r
862 \r
863   recStatus = kTRUE;\r
864 \r
865   /*\r
866   if(TMath::Abs(dca[0])<fCutsRC->GetMaxDCAToVertexXY() && \r
867      TMath::Abs(dca[1])<fCutsRC->GetMaxDCAToVertexZ())\r
868   {\r
869     recStatus = kTRUE;\r
870   }\r
871   */\r
872 \r
873 return recStatus;\r
874 }\r
875 \r
876 //_____________________________________________________________________________\r
877 Long64_t AliPerformanceEff::Merge(TCollection* const list) \r
878 {\r
879   // Merge list of objects (needed by PROOF)\r
880 \r
881   if (!list)\r
882   return 0;\r
883 \r
884   if (list->IsEmpty())\r
885   return 1;\r
886 \r
887   TIterator* iter = list->MakeIterator();\r
888   TObject* obj = 0;\r
889 \r
890   // collection of generated histograms\r
891 \r
892   Int_t count=0;\r
893   while((obj = iter->Next()) != 0) \r
894   {\r
895     AliPerformanceEff* entry = dynamic_cast<AliPerformanceEff*>(obj);\r
896     if (entry == 0) continue; \r
897   \r
898      fEffHisto->Add(entry->fEffHisto);\r
899      fEffSecHisto->Add(entry->fEffSecHisto);\r
900   count++;\r
901   }\r
902 \r
903 return count;\r
904 }\r
905  \r
906 //_____________________________________________________________________________\r
907 void AliPerformanceEff::Analyse() \r
908 {\r
909   // Analyse comparison information and store output histograms\r
910   // in the folder "folderEff" \r
911   //\r
912   TH1::AddDirectory(kFALSE);\r
913   TObjArray *aFolderObj = new TObjArray;\r
914   if(!aFolderObj) return;\r
915   char title[256];\r
916 \r
917   //\r
918   // efficiency vs pt\r
919   //\r
920 \r
921   if(GetAnalysisMode() != 5) {\r
922 \r
923   fEffHisto->GetAxis(0)->SetRangeUser(-0.9,0.89); // eta range\r
924   fEffHisto->GetAxis(2)->SetRangeUser(0.1,10.);   // pt range\r
925 \r
926   // rec efficiency vs pt\r
927   fEffHisto->GetAxis(3)->SetRangeUser(0.,3.99);  // reconstructed \r
928 \r
929   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
930   aFolderObj->Add(AddHistoEff(2, "ptRecEff", "rec. efficiency"));\r
931 \r
932   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
933   aFolderObj->Add(AddHistoEff(2, "ptRecEffNeg", "rec. efficiency neg."));\r
934 \r
935   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
936   aFolderObj->Add(AddHistoEff(2, "ptRecEffPos", "rec. efficiency pos."));\r
937 \r
938   // rec efficiency vs pid vs pt\r
939   fEffHisto->GetAxis(3)->SetRangeUser(2.,2.);    // pions\r
940   \r
941   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
942   aFolderObj->Add(AddHistoEff(2, "ptRecEffPi", "rec. efficiency (pions)"));\r
943  \r
944   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
945   aFolderObj->Add(AddHistoEff(2, "ptRecEffPiNeg", "rec. efficiency (pions) neg."));\r
946 \r
947   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
948   aFolderObj->Add(AddHistoEff(2, "ptRecEffPiPos", "rec. efficiency (pions) pos."));\r
949 \r
950 \r
951   fEffHisto->GetAxis(3)->SetRangeUser(3.,3.);    // kaons\r
952   \r
953   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
954   aFolderObj->Add(AddHistoEff(2, "ptRecEffK", "rec. efficiency (kaons)"));\r
955   \r
956   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
957   aFolderObj->Add(AddHistoEff(2, "ptRecEffKNeg", "rec. efficiency (kaons) neg."));\r
958  \r
959   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
960   aFolderObj->Add(AddHistoEff(2, "ptRecEffKPos", "rec. efficiency (kaons) pos."));\r
961  \r
962   \r
963   fEffHisto->GetAxis(3)->SetRangeUser(4.,4.);    // protons\r
964 \r
965   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
966   aFolderObj->Add(AddHistoEff(2, "ptRecEffP", "rec. efficiency (protons)"));\r
967 \r
968   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
969   aFolderObj->Add(AddHistoEff(2, "ptRecEffPNeg", "rec. efficiency (protons) neg."));\r
970  \r
971   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
972   aFolderObj->Add(AddHistoEff(2, "ptRecEffPPos", "rec. efficiency (protons) pos."));\r
973 \r
974   // findable efficiency vs pt\r
975   fEffHisto->GetAxis(3)->SetRangeUser(0.,4.); \r
976   fEffHisto->GetAxis(5)->SetRangeUser(1.,1.); // findable\r
977 \r
978   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
979   aFolderObj->Add(AddHistoEff(2, "ptRecEffF", "rec. efficiency (findable)"));\r
980 \r
981   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
982   aFolderObj->Add(AddHistoEff(2, "ptRecEffFNeg", "rec. efficiency (findable) neg."));\r
983 \r
984   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
985   aFolderObj->Add(AddHistoEff(2, "ptRecEffFPos", "rec. efficiency (findable) pos."));\r
986 \r
987   //\r
988   // efficiency vs eta\r
989   //\r
990 \r
991   fEffHisto->GetAxis(0)->SetRangeUser(-1.5,1.5); // eta range\r
992   fEffHisto->GetAxis(2)->SetRangeUser(0.2,10.); // pt range\r
993   fEffHisto->GetAxis(5)->SetRangeUser(0.,1.);   // all\r
994 \r
995   // rec efficiency vs eta\r
996   fEffHisto->GetAxis(3)->SetRangeUser(0.,4.);  // reconstructed \r
997 \r
998   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
999   aFolderObj->Add(AddHistoEff(0, "etaRecEff", "rec. efficiency"));\r
1000 \r
1001   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
1002   aFolderObj->Add(AddHistoEff(0, "etaRecEffNeg", "rec. efficiency neg."));\r
1003 \r
1004   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
1005   aFolderObj->Add(AddHistoEff(0, "etaRecEffPos", "rec. efficiency pos."));\r
1006 \r
1007   // rec efficiency vs pid vs eta\r
1008   fEffHisto->GetAxis(3)->SetRangeUser(2.,2.);    // pions\r
1009   \r
1010   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
1011   aFolderObj->Add(AddHistoEff(0, "etaRecEffPi", "rec. efficiency (pions)"));\r
1012 \r
1013   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
1014   aFolderObj->Add(AddHistoEff(0, "etaRecEffPiNeg", "rec. efficiency (pions) neg."));\r
1015 \r
1016   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
1017   aFolderObj->Add(AddHistoEff(0, "etaRecEffPiPos", "rec. efficiency (pions) pos."));\r
1018 \r
1019 \r
1020   fEffHisto->GetAxis(3)->SetRangeUser(3.,3.);    // kaons\r
1021   \r
1022   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
1023   aFolderObj->Add(AddHistoEff(0, "etaRecEffK", "rec. efficiency (kaons)"));\r
1024 \r
1025   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
1026   aFolderObj->Add(AddHistoEff(0, "etaRecEffKNeg", "rec. efficiency (kaons) neg."));\r
1027   \r
1028   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
1029   aFolderObj->Add(AddHistoEff(0, "etaRecEffKPos", "rec. efficiency (kaons) pos."));\r
1030   \r
1031   \r
1032   fEffHisto->GetAxis(3)->SetRangeUser(4.,4.);    // protons\r
1033 \r
1034   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
1035   aFolderObj->Add(AddHistoEff(0, "etaRecEffP", "rec. efficiency (protons)"));\r
1036  \r
1037   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
1038   aFolderObj->Add(AddHistoEff(0, "etaRecEffPNeg", "rec. efficiency (protons) neg."));\r
1039   \r
1040   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
1041   aFolderObj->Add(AddHistoEff(0, "etaRecEffPPos", "rec. efficiency (protons) pos."));\r
1042 \r
1043   // findable efficiency vs eta\r
1044   fEffHisto->GetAxis(3)->SetRangeUser(0.,4.); \r
1045   fEffHisto->GetAxis(5)->SetRangeUser(1.,1.); // findable\r
1046 \r
1047   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
1048   aFolderObj->Add(AddHistoEff(0, "etaRecEffF", "rec. efficiency (findable)"));\r
1049 \r
1050   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
1051   aFolderObj->Add(AddHistoEff(0, "etaRecEffFNeg", "rec. efficiency (findable) neg."));\r
1052 \r
1053   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
1054   aFolderObj->Add(AddHistoEff(0, "etaRecEffFPos", "rec. efficiency (findable) pos."));\r
1055 \r
1056   //\r
1057   // efficiency vs phi\r
1058   //\r
1059 \r
1060   fEffHisto->GetAxis(0)->SetRangeUser(-0.9,0.9); // eta range\r
1061   fEffHisto->GetAxis(2)->SetRangeUser(0.2,10.); // pt range\r
1062   fEffHisto->GetAxis(5)->SetRangeUser(0.,1.);   // all\r
1063 \r
1064   // rec efficiency vs phi\r
1065   fEffHisto->GetAxis(3)->SetRangeUser(0.,4.);  // reconstructed \r
1066 \r
1067   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
1068   aFolderObj->Add(AddHistoEff(1, "phiRecEff", "rec. efficiency"));\r
1069 \r
1070   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
1071   aFolderObj->Add(AddHistoEff(1, "phiRecEffNeg", "rec. efficiency neg."));\r
1072 \r
1073   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
1074   aFolderObj->Add(AddHistoEff(1, "phiRecEffPos", "rec. efficiency pos."));\r
1075 \r
1076   // rec efficiency vs pid vs phi\r
1077   fEffHisto->GetAxis(3)->SetRangeUser(2.,2.);    // pions\r
1078   \r
1079   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
1080   aFolderObj->Add(AddHistoEff(1, "phiRecEffPi", "rec. efficiency (pions)"));\r
1081 \r
1082   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
1083   aFolderObj->Add(AddHistoEff(1, "phiRecEffPiNeg", "rec. efficiency (pions) neg."));\r
1084 \r
1085   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
1086   aFolderObj->Add(AddHistoEff(1, "phiRecEffPiPos", "rec. efficiency (pions) pos."));\r
1087 \r
1088 \r
1089   fEffHisto->GetAxis(3)->SetRangeUser(3.,3.);    // kaons\r
1090   \r
1091   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
1092   aFolderObj->Add(AddHistoEff(1, "phiRecEffK", "rec. efficiency (kaons)"));\r
1093 \r
1094   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
1095   aFolderObj->Add(AddHistoEff(1, "phiRecEffKNeg", "rec. efficiency (kaons) neg."));\r
1096 \r
1097   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
1098   aFolderObj->Add(AddHistoEff(1, "phiRecEffKPos", "rec. efficiency (kaons) pos."));\r
1099  \r
1100   \r
1101   fEffHisto->GetAxis(3)->SetRangeUser(4.,4.);    // protons\r
1102 \r
1103   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
1104   aFolderObj->Add(AddHistoEff(1, "phiRecEffP", "rec. efficiency (protons)"));\r
1105 \r
1106   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
1107   aFolderObj->Add(AddHistoEff(1, "phiRecEffPNeg", "rec. efficiency (protons) neg."));\r
1108  \r
1109   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
1110   aFolderObj->Add(AddHistoEff(1, "phiRecEffPPos", "rec. efficiency (protons) pos."));\r
1111 \r
1112   // findable efficiency vs phi\r
1113   fEffHisto->GetAxis(3)->SetRangeUser(0.,4.); \r
1114   fEffHisto->GetAxis(5)->SetRangeUser(1.,1.); // findable\r
1115 \r
1116   fEffHisto->GetAxis(6)->SetRangeUser(-2.,2.);   // charge all\r
1117   aFolderObj->Add(AddHistoEff(1, "phiRecEffF", "rec. efficiency (findable)"));\r
1118 \r
1119   fEffHisto->GetAxis(6)->SetRangeUser(-2.,0.);   // charge negativ\r
1120   aFolderObj->Add(AddHistoEff(1, "phiRecEffFNeg", "rec. efficiency (findable) neg."));\r
1121  \r
1122   fEffHisto->GetAxis(6)->SetRangeUser(0.,2.);    // charge positiv\r
1123   aFolderObj->Add(AddHistoEff(1, "phiRecEffFPos", "rec. efficiency (findable) pos."));\r
1124   }\r
1125   else {\r
1126   // \r
1127   Float_t minEta=-1.5, maxEta=1.5;\r
1128   Float_t minR=0.0, maxR=150.0;\r
1129   Float_t minPt=0.15, maxPt=100.0;\r
1130 \r
1131   // mother eta range\r
1132   fEffSecHisto->GetAxis(8)->SetRangeUser(minEta,maxEta);\r
1133 \r
1134   // particle creation radius range \r
1135   fEffSecHisto->GetAxis(6)->SetRangeUser(minR,maxR);\r
1136 \r
1137   //\r
1138   fEffSecHisto->GetAxis(0)->SetRangeUser(minEta,maxEta);\r
1139   fEffSecHisto->GetAxis(2)->SetRangeUser(minPt,maxPt);\r
1140 \r
1141   // rec efficiency vs pt\r
1142   TH1D *ptAll = fEffSecHisto->Projection(2);\r
1143 \r
1144   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.);  // reconstructed \r
1145   TH1D *ptRec = fEffSecHisto->Projection(2);\r
1146   TH1D *ptRecc = (TH1D*)ptRec->Clone();\r
1147   ptRecc->Divide(ptRec,ptAll,1,1,"B");\r
1148   ptRecc->SetName("ptRecEff");\r
1149 \r
1150   ptRecc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(2)->GetTitle());\r
1151   ptRecc->GetYaxis()->SetTitle("efficiency");\r
1152   snprintf(title,256,"%s vs %s","rec. efficiency",fEffSecHisto->GetAxis(2)->GetTitle());\r
1153   ptRecc->SetTitle(title);\r
1154 \r
1155   ptRecc->SetBit(TH1::kLogX);\r
1156   aFolderObj->Add(ptRecc);\r
1157 \r
1158   // rec efficiency vs pid vs pt\r
1159   \r
1160   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1161   fEffSecHisto->GetAxis(3)->SetRangeUser(0.,0.); // electrons\r
1162 \r
1163   TH1D *ptAllEle = fEffSecHisto->Projection(2);\r
1164 \r
1165   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1166   TH1D *ptRecEle = fEffSecHisto->Projection(2);\r
1167   TH1D *ptRecElec = (TH1D*)ptRecEle->Clone();\r
1168   ptRecElec->Divide(ptRecEle,ptAllEle,1,1,"B");\r
1169   ptRecElec->SetName("ptRecEffEle");\r
1170 \r
1171   ptRecElec->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(2)->GetTitle());\r
1172   ptRecElec->GetYaxis()->SetTitle("efficiency");\r
1173   snprintf(title,256,"%s vs %s","rec. efficiency (electrons)",fEffSecHisto->GetAxis(2)->GetTitle());\r
1174   ptRecElec->SetTitle(title);\r
1175 \r
1176   ptRecElec->SetBit(TH1::kLogX);\r
1177   aFolderObj->Add(ptRecElec);\r
1178 \r
1179   //\r
1180 \r
1181   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1182   fEffSecHisto->GetAxis(3)->SetRangeUser(2.,2.); // pions\r
1183 \r
1184   TH1D *ptAllPi = fEffSecHisto->Projection(2);\r
1185 \r
1186   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1187   TH1D *ptRecPi = fEffSecHisto->Projection(2);\r
1188   TH1D *ptRecPic = (TH1D*)ptRecPi->Clone();\r
1189   ptRecPic->Divide(ptRecPi,ptAllPi,1,1,"B");\r
1190   ptRecPic->SetName("ptRecEffPi");\r
1191 \r
1192   ptRecPic->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(2)->GetTitle());\r
1193   ptRecPic->GetYaxis()->SetTitle("efficiency");\r
1194   snprintf(title,256,"%s vs %s","rec. efficiency (pions)",fEffSecHisto->GetAxis(2)->GetTitle());\r
1195   ptRecPic->SetTitle(title);\r
1196 \r
1197   ptRecPic->SetBit(TH1::kLogX);\r
1198   aFolderObj->Add(ptRecPic);\r
1199 \r
1200   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1201   fEffSecHisto->GetAxis(3)->SetRangeUser(3.,3.); // kaons\r
1202   TH1D *ptAllK = fEffSecHisto->Projection(2);\r
1203 \r
1204   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1205   TH1D *ptRecK = fEffSecHisto->Projection(2);\r
1206 \r
1207   TH1D *ptRecKc = (TH1D*)ptRecK->Clone();\r
1208   ptRecKc->Divide(ptRecK,ptAllK,1,1,"B");\r
1209   ptRecKc->SetName("ptRecEffK");\r
1210 \r
1211   ptRecKc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(2)->GetTitle());\r
1212   ptRecKc->GetYaxis()->SetTitle("efficiency");\r
1213   snprintf(title,256,"%s vs %s","rec. efficiency (kaons)",fEffSecHisto->GetAxis(2)->GetTitle());\r
1214   ptRecKc->SetTitle(title);\r
1215 \r
1216 \r
1217   ptRecKc->SetBit(TH1::kLogX);\r
1218   aFolderObj->Add(ptRecKc);\r
1219 \r
1220   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1221   fEffSecHisto->GetAxis(3)->SetRangeUser(4.,4.); // protons\r
1222   TH1D *ptAllP = fEffSecHisto->Projection(2);\r
1223 \r
1224   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1225   TH1D *ptRecP = fEffSecHisto->Projection(2);\r
1226   TH1D *ptRecPc = (TH1D*)ptRecP->Clone();\r
1227   ptRecPc->Divide(ptRecP,ptAllP,1,1,"B");\r
1228   ptRecPc->SetName("ptRecEffP");\r
1229 \r
1230   ptRecPc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(2)->GetTitle());\r
1231   ptRecPc->GetYaxis()->SetTitle("efficiency");\r
1232   snprintf(title,256,"%s vs %s","rec. efficiency (protons)",fEffSecHisto->GetAxis(2)->GetTitle());\r
1233   ptRecPc->SetTitle(title);\r
1234 \r
1235   ptRecPc->SetBit(TH1::kLogX);\r
1236   aFolderObj->Add(ptRecPc);\r
1237 \r
1238   // findable efficiency vs pt\r
1239 \r
1240   fEffSecHisto->GetAxis(3)->SetRangeUser(0.,4.); \r
1241   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1242   fEffSecHisto->GetAxis(5)->SetRangeUser(1.,1.); // findable\r
1243   TH1D *ptAllF = fEffSecHisto->Projection(2);\r
1244 \r
1245   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.);\r
1246   fEffSecHisto->GetAxis(5)->SetRangeUser(1.,1.);\r
1247 \r
1248   TH1D *ptRecF = fEffSecHisto->Projection(2); // rec findable\r
1249   TH1D *ptRecFc = (TH1D*)ptRecF->Clone();\r
1250   ptRecFc->Divide(ptRecF,ptAllF,1,1,"B");\r
1251   ptRecFc->SetName("ptRecEffF");\r
1252 \r
1253   ptRecFc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(2)->GetTitle());\r
1254   ptRecFc->GetYaxis()->SetTitle("efficiency");\r
1255   snprintf(title,256,"%s vs %s","rec. efficiency (findable)",fEffSecHisto->GetAxis(2)->GetTitle());\r
1256   ptRecFc->SetTitle(title);\r
1257 \r
1258   ptRecFc->SetBit(TH1::kLogX);\r
1259   aFolderObj->Add(ptRecFc);\r
1260 \r
1261   //\r
1262   // efficiency vs eta\r
1263   //\r
1264   fEffSecHisto->GetAxis(2)->SetRangeUser(minPt,maxPt);\r
1265   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.);   // all\r
1266   fEffSecHisto->GetAxis(5)->SetRangeUser(0.,1.);   // all\r
1267 \r
1268   TH1D *etaAll = fEffSecHisto->Projection(0);\r
1269 \r
1270   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.);  // reconstructed \r
1271   TH1D *etaRec = fEffSecHisto->Projection(0);\r
1272   TH1D *etaRecc = (TH1D*)etaRec->Clone();\r
1273   etaRecc->Divide(etaRec,etaAll,1,1,"B");\r
1274   etaRecc->SetName("etaRecEff");\r
1275 \r
1276   etaRecc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(0)->GetTitle());\r
1277   etaRecc->GetYaxis()->SetTitle("efficiency");\r
1278   snprintf(title,256,"%s vs %s","rec. efficiency",fEffSecHisto->GetAxis(0)->GetTitle());\r
1279   etaRecc->SetTitle(title);\r
1280 \r
1281   aFolderObj->Add(etaRecc);\r
1282 \r
1283   // rec efficiency vs pid vs eta\r
1284   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1285   fEffSecHisto->GetAxis(3)->SetRangeUser(0.,0.); // electrons\r
1286 \r
1287   TH1D *etaAllEle = fEffSecHisto->Projection(0);\r
1288 \r
1289   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1290   TH1D *etaRecEle = fEffSecHisto->Projection(0);\r
1291   TH1D *etaRecElec = (TH1D*)etaRecEle->Clone();\r
1292   etaRecElec->Divide(etaRecEle,etaAllEle,1,1,"B");\r
1293   etaRecElec->SetName("etaRecEffEle");\r
1294 \r
1295   etaRecElec->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(0)->GetTitle());\r
1296   etaRecElec->GetYaxis()->SetTitle("efficiency");\r
1297   snprintf(title,256,"%s vs %s","rec. efficiency (electrons)",fEffSecHisto->GetAxis(0)->GetTitle());\r
1298   etaRecElec->SetTitle(title);\r
1299 \r
1300   aFolderObj->Add(etaRecElec);\r
1301 \r
1302   //\r
1303   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1304   fEffSecHisto->GetAxis(3)->SetRangeUser(2.,2.); // pions\r
1305 \r
1306   TH1D *etaAllPi = fEffSecHisto->Projection(0);\r
1307 \r
1308   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1309   TH1D *etaRecPi = fEffSecHisto->Projection(0);\r
1310   TH1D *etaRecPic = (TH1D*)etaRecPi->Clone();\r
1311   etaRecPic->Divide(etaRecPi,etaAllPi,1,1,"B");\r
1312   etaRecPic->SetName("etaRecEffPi");\r
1313 \r
1314   etaRecPic->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(0)->GetTitle());\r
1315   etaRecPic->GetYaxis()->SetTitle("efficiency");\r
1316   snprintf(title,256,"%s vs %s","rec. efficiency (pions)",fEffSecHisto->GetAxis(0)->GetTitle());\r
1317   etaRecPic->SetTitle(title);\r
1318 \r
1319   aFolderObj->Add(etaRecPic);\r
1320 \r
1321   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1322   fEffSecHisto->GetAxis(3)->SetRangeUser(3.,3.); // kaons\r
1323   TH1D *etaAllK = fEffSecHisto->Projection(0);\r
1324 \r
1325   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1326   TH1D *etaRecK = fEffSecHisto->Projection(0);\r
1327 \r
1328   TH1D *etaRecKc = (TH1D*)etaRecK->Clone();\r
1329   etaRecKc->Divide(etaRecK,etaAllK,1,1,"B");\r
1330   etaRecKc->SetName("etaRecEffK");\r
1331 \r
1332   etaRecKc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(0)->GetTitle());\r
1333   etaRecKc->GetYaxis()->SetTitle("efficiency");\r
1334   snprintf(title,256,"%s vs %s","rec. efficiency (kaons)",fEffSecHisto->GetAxis(0)->GetTitle());\r
1335   etaRecKc->SetTitle(title);\r
1336 \r
1337 \r
1338   aFolderObj->Add(etaRecKc);\r
1339 \r
1340   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1341   fEffSecHisto->GetAxis(3)->SetRangeUser(4.,4.); // protons\r
1342   TH1D *etaAllP = fEffSecHisto->Projection(0);\r
1343 \r
1344   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1345   TH1D *etaRecP = fEffSecHisto->Projection(0);\r
1346   TH1D *etaRecPc = (TH1D*)etaRecP->Clone();\r
1347   etaRecPc->Divide(etaRecP,etaAllP,1,1,"B");\r
1348   etaRecPc->SetName("etaRecEffP");\r
1349 \r
1350   etaRecPc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(0)->GetTitle());\r
1351   etaRecPc->GetYaxis()->SetTitle("efficiency");\r
1352   snprintf(title,256,"%s vs %s","rec. efficiency (protons)",fEffSecHisto->GetAxis(0)->GetTitle());\r
1353   etaRecPc->SetTitle(title);\r
1354 \r
1355   aFolderObj->Add(etaRecPc);\r
1356 \r
1357   // findable efficiency vs eta\r
1358 \r
1359   fEffSecHisto->GetAxis(3)->SetRangeUser(0.,4.); \r
1360   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1361   fEffSecHisto->GetAxis(5)->SetRangeUser(1.,1.); // findable\r
1362   TH1D *etaAllF = fEffSecHisto->Projection(0);\r
1363 \r
1364   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.);\r
1365   fEffSecHisto->GetAxis(5)->SetRangeUser(1.,1.);\r
1366 \r
1367   TH1D *etaRecF = fEffSecHisto->Projection(0); // rec findable\r
1368   TH1D *etaRecFc = (TH1D*)etaRecF->Clone();\r
1369   etaRecFc->Divide(etaRecF,etaAllF,1,1,"B");\r
1370   etaRecFc->SetName("etaRecEffF");\r
1371 \r
1372   etaRecFc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(0)->GetTitle());\r
1373   etaRecFc->GetYaxis()->SetTitle("efficiency");\r
1374   snprintf(title,256,"%s vs %s","rec. efficiency (findable)",fEffSecHisto->GetAxis(0)->GetTitle());\r
1375   etaRecFc->SetTitle(title);\r
1376 \r
1377   aFolderObj->Add(etaRecFc);\r
1378 \r
1379   //\r
1380   // efficiency vs phi\r
1381   //\r
1382 \r
1383   fEffSecHisto->GetAxis(0)->SetRangeUser(minEta,maxEta);\r
1384   fEffSecHisto->GetAxis(2)->SetRangeUser(minPt,maxPt);\r
1385 \r
1386   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.);   // all\r
1387   fEffSecHisto->GetAxis(5)->SetRangeUser(0.,1.);   // all\r
1388 \r
1389   TH1D *phiAll = fEffSecHisto->Projection(1);\r
1390 \r
1391   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.);  // reconstructed \r
1392   TH1D *phiRec = fEffSecHisto->Projection(1);\r
1393   TH1D *phiRecc = (TH1D*)phiRec->Clone();\r
1394   phiRecc->Divide(phiRec,phiAll,1,1,"B");\r
1395   phiRecc->SetName("phiRecEff");\r
1396 \r
1397   phiRecc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(1)->GetTitle());\r
1398   phiRecc->GetYaxis()->SetTitle("efficiency");\r
1399   snprintf(title,256,"%s vs %s","rec. efficiency",fEffSecHisto->GetAxis(1)->GetTitle());\r
1400   phiRecc->SetTitle(title);\r
1401 \r
1402   aFolderObj->Add(phiRecc);\r
1403 \r
1404   // rec efficiency vs pid vs phi\r
1405   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1406   fEffSecHisto->GetAxis(3)->SetRangeUser(2.,2.); // pions\r
1407 \r
1408   TH1D *phiAllEle = fEffSecHisto->Projection(1);\r
1409 \r
1410   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1411   TH1D *phiRecEle = fEffSecHisto->Projection(1);\r
1412   TH1D *phiRecElec = (TH1D*)phiRecEle->Clone();\r
1413   phiRecElec->Divide(phiRecEle,phiAllEle,1,1,"B");\r
1414   phiRecElec->SetName("phiRecEffEle");\r
1415 \r
1416   phiRecElec->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(1)->GetTitle());\r
1417   phiRecElec->GetYaxis()->SetTitle("efficiency");\r
1418   snprintf(title,256,"%s vs %s","rec. efficiency (electrons)",fEffSecHisto->GetAxis(1)->GetTitle());\r
1419   phiRecElec->SetTitle(title);\r
1420 \r
1421   aFolderObj->Add(phiRecElec);\r
1422 \r
1423   //\r
1424   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1425   fEffSecHisto->GetAxis(3)->SetRangeUser(2.,2.); // pions\r
1426 \r
1427   TH1D *phiAllPi = fEffSecHisto->Projection(1);\r
1428 \r
1429   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1430   TH1D *phiRecPi = fEffSecHisto->Projection(1);\r
1431   TH1D *phiRecPic = (TH1D*)phiRecPi->Clone();\r
1432   phiRecPic->Divide(phiRecPi,phiAllPi,1,1,"B");\r
1433   phiRecPic->SetName("phiRecEffPi");\r
1434 \r
1435   phiRecPic->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(1)->GetTitle());\r
1436   phiRecPic->GetYaxis()->SetTitle("efficiency");\r
1437   snprintf(title,256,"%s vs %s","rec. efficiency (pions)",fEffSecHisto->GetAxis(1)->GetTitle());\r
1438   phiRecPic->SetTitle(title);\r
1439 \r
1440   aFolderObj->Add(phiRecPic);\r
1441 \r
1442   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1443   fEffSecHisto->GetAxis(3)->SetRangeUser(3.,3.); // kaons\r
1444   TH1D *phiAllK = fEffSecHisto->Projection(1);\r
1445 \r
1446   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1447   TH1D *phiRecK = fEffSecHisto->Projection(1);\r
1448 \r
1449   TH1D *phiRecKc = (TH1D*)phiRecK->Clone();\r
1450   phiRecKc->Divide(phiRecK,phiAllK,1,1,"B");\r
1451   phiRecKc->SetName("phiRecEffK");\r
1452 \r
1453   phiRecKc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(1)->GetTitle());\r
1454   phiRecKc->GetYaxis()->SetTitle("efficiency");\r
1455   snprintf(title,256,"%s vs %s","rec. efficiency (kaons)",fEffSecHisto->GetAxis(1)->GetTitle());\r
1456   phiRecKc->SetTitle(title);\r
1457 \r
1458 \r
1459   aFolderObj->Add(phiRecKc);\r
1460 \r
1461   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1462   fEffSecHisto->GetAxis(3)->SetRangeUser(4.,4.); // protons\r
1463   TH1D *phiAllP = fEffSecHisto->Projection(1);\r
1464 \r
1465   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.); // reconstructed\r
1466   TH1D *phiRecP = fEffSecHisto->Projection(1);\r
1467   TH1D *phiRecPc = (TH1D*)phiRecP->Clone();\r
1468   phiRecPc->Divide(phiRecP,phiAllP,1,1,"B");\r
1469   phiRecPc->SetName("phiRecEffP");\r
1470 \r
1471   phiRecPc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(1)->GetTitle());\r
1472   phiRecPc->GetYaxis()->SetTitle("efficiency");\r
1473   snprintf(title,256,"%s vs %s","rec. efficiency (protons)",fEffSecHisto->GetAxis(1)->GetTitle());\r
1474   phiRecPc->SetTitle(title);\r
1475 \r
1476   aFolderObj->Add(phiRecPc);\r
1477 \r
1478   // findable efficiency vs phi\r
1479 \r
1480   fEffSecHisto->GetAxis(3)->SetRangeUser(0.,4.); \r
1481   fEffSecHisto->GetAxis(4)->SetRangeUser(0.,1.); \r
1482   fEffSecHisto->GetAxis(5)->SetRangeUser(1.,1.); // findable\r
1483   TH1D *phiAllF = fEffSecHisto->Projection(1);\r
1484 \r
1485   fEffSecHisto->GetAxis(4)->SetRangeUser(1.,1.);\r
1486   fEffSecHisto->GetAxis(5)->SetRangeUser(1.,1.);\r
1487 \r
1488   TH1D *phiRecF = fEffSecHisto->Projection(1); // rec findable\r
1489   TH1D *phiRecFc = (TH1D*)phiRecF->Clone();\r
1490   phiRecFc->Divide(phiRecF,phiAllF,1,1,"B");\r
1491   phiRecFc->SetName("phiRecEffF");\r
1492 \r
1493   phiRecFc->GetXaxis()->SetTitle(fEffSecHisto->GetAxis(1)->GetTitle());\r
1494   phiRecFc->GetYaxis()->SetTitle("efficiency");\r
1495   snprintf(title,256,"%s vs %s","rec. efficiency (findable)",fEffSecHisto->GetAxis(1)->GetTitle());\r
1496   phiRecFc->SetTitle(title);\r
1497 \r
1498   aFolderObj->Add(phiRecFc);\r
1499   }\r
1500 \r
1501   // export objects to analysis folder\r
1502   fAnalysisFolder = ExportToFolder(aFolderObj);\r
1503 \r
1504   // delete only TObjArray\r
1505   if(aFolderObj) delete aFolderObj;\r
1506 }\r
1507 \r
1508 //_____________________________________________________________________________\r
1509 TFolder* AliPerformanceEff::ExportToFolder(TObjArray * array) \r
1510 {\r
1511   // recreate folder avery time and export objects to new one\r
1512   //\r
1513   AliPerformanceEff * comp=this;\r
1514   TFolder *folder = comp->GetAnalysisFolder();\r
1515 \r
1516   TString name, title;\r
1517   TFolder *newFolder = 0;\r
1518   Int_t i = 0;\r
1519   Int_t size = array->GetSize();\r
1520 \r
1521   if(folder) { \r
1522      // get name and title from old folder\r
1523      name = folder->GetName();  \r
1524      title = folder->GetTitle();  \r
1525 \r
1526          // delete old one\r
1527      delete folder;\r
1528 \r
1529          // create new one\r
1530      newFolder = CreateFolder(name.Data(),title.Data());\r
1531      newFolder->SetOwner();\r
1532 \r
1533          // add objects to folder\r
1534      while(i < size) {\r
1535            newFolder->Add(array->At(i));\r
1536            i++;\r
1537          }\r
1538   }\r
1539 \r
1540 return newFolder;\r
1541 }\r
1542 \r
1543 \r
1544 //_____________________________________________________________________________\r
1545 TFolder* AliPerformanceEff::CreateFolder(TString name,TString title) { \r
1546 // create folder for analysed histograms\r
1547 //\r
1548 TFolder *folder = 0;\r
1549   folder = new TFolder(name.Data(),title.Data());\r
1550 \r
1551   return folder;\r
1552 }\r
1553 \r
1554 \r
1555 //_____________________________________________________________________________\r
1556 TH1D* AliPerformanceEff::AddHistoEff(Int_t axis, const Char_t *name, const Char_t* vsTitle) {\r
1557   // Create and add rec efficiency vs pt, eta, phi\r
1558   \r
1559   char title[256];\r
1560 \r
1561   fEffHisto->GetAxis(4)->SetRangeUser(0.,1.);    // all\r
1562   TH1D *all = fEffHisto->Projection(axis);\r
1563 \r
1564   fEffHisto->GetAxis(4)->SetRangeUser(1.,1.);    // reconstructed \r
1565   TH1D *rec = fEffHisto->Projection(axis);\r
1566   TH1D *recc = (TH1D*)rec->Clone();\r
1567 \r
1568   recc->Divide(rec,all,1,1,"B");\r
1569   recc->SetName(name);\r
1570 \r
1571   recc->GetXaxis()->SetTitle(fEffHisto->GetAxis(axis)->GetTitle());\r
1572   recc->GetYaxis()->SetTitle("efficiency");\r
1573 \r
1574   snprintf(title,256,"%s vs %s",vsTitle, fEffHisto->GetAxis(axis)->GetTitle());  \r
1575   recc->SetTitle(title);\r
1576 \r
1577   if (axis == 2 ) recc->SetBit(TH1::kLogX);\r
1578 \r
1579   return recc;\r
1580 }\r