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