]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/dNdPt/AlidNdPtEfficiency.cxx
2dfb3e77fd65d816c9b0a303f31633a0068b2c45
[u/mrichter/AliRoot.git] / PWG0 / dNdPt / AlidNdPtEfficiency.cxx
1 /**************************************************************************\r
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *\r
3  *                                                                        *\r
4  * Author: The ALICE Off-line Project.                                    *\r
5  * Contributors are mentioned in the code where appropriate.              *\r
6  *                                                                        *\r
7  * Permission to use, copy, modify and distribute this software and its   *\r
8  * documentation strictly for non-commercial purposes is hereby granted   *\r
9  * without fee, provided that the above copyright notice appears in all   *\r
10  * copies and that both the copyright notice and this permission notice   *\r
11  * appear in the supporting documentation. The authors make no claims     *\r
12  * about the suitability of this software for any purpose. It is          *\r
13  * provided "as is" without express or implied warranty.                  *\r
14  **************************************************************************/\r
15 //------------------------------------------------------------------------------\r
16 // AlidNdPtEfficiency class. \r
17 //\r
18 // a. functionality:\r
19 // - fills generic cut histograms\r
20 // - generates cuts (selection criteria)\r
21 //\r
22 // b. data members:\r
23 // - generic cut histograms\r
24 // - control histograms\r
25 //\r
26 // Author: J.Otwinowski 18/11/2010 \r
27 //------------------------------------------------------------------------------\r
28 #include "TH1.h"\r
29 #include "TH2.h"\r
30 \r
31 #include "AliHeader.h"  \r
32 #include "AliGenEventHeader.h"  \r
33 #include "AliStack.h"  \r
34 #include "AliESDEvent.h"  \r
35 #include "AliMCEvent.h"  \r
36 #include "AliESDtrackCuts.h"  \r
37 #include "AliLog.h" \r
38 #include "AliTracker.h" \r
39 \r
40 #include "AlidNdPtEventCuts.h"\r
41 #include "AlidNdPtAcceptanceCuts.h"\r
42 #include "AlidNdPtBackgroundCuts.h"\r
43 #include "AlidNdPtAnalysis.h"\r
44 #include "AliPhysicsSelection.h"\r
45 \r
46 #include "AliPWG0Helper.h"\r
47 #include "AlidNdPtHelper.h"\r
48 #include "AlidNdPtEfficiency.h"\r
49 \r
50 using namespace std;\r
51 \r
52 ClassImp(AlidNdPtEfficiency)\r
53 \r
54 //_____________________________________________________________________________\r
55   AlidNdPtEfficiency::AlidNdPtEfficiency(): AlidNdPt(),\r
56   fAnalysisFolder(0),\r
57   fRecMCTrackHistTPCITS(0),\r
58   fRecMCTrackHistITSTPC(0)\r
59 {\r
60   // default constructor\r
61   Init();\r
62 }\r
63 \r
64 //_____________________________________________________________________________\r
65 AlidNdPtEfficiency::AlidNdPtEfficiency(Char_t* name, Char_t* title): AlidNdPt(name,title),\r
66   fAnalysisFolder(0),\r
67   fRecMCTrackHistTPCITS(0),\r
68   fRecMCTrackHistITSTPC(0)\r
69 {\r
70   // constructor\r
71   Init();\r
72 }\r
73 \r
74 //_____________________________________________________________________________\r
75 AlidNdPtEfficiency::~AlidNdPtEfficiency() {\r
76   // \r
77   if(fRecMCTrackHistTPCITS) delete fRecMCTrackHistTPCITS; fRecMCTrackHistTPCITS=0;\r
78   if(fRecMCTrackHistITSTPC) delete fRecMCTrackHistITSTPC; fRecMCTrackHistITSTPC=0;\r
79 \r
80   if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;\r
81 }\r
82 \r
83 //_____________________________________________________________________________\r
84 void AlidNdPtEfficiency::Init(){\r
85   //\r
86   // Init histograms\r
87   //\r
88   const Int_t ptNbins = 58; \r
89   const Double_t ptMin = 0.; \r
90   const Double_t ptMax = 20.; \r
91 \r
92   Double_t binsPt[ptNbins+1] = {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.65,0.7,0.75,0.8,0.85,0.9,0.95,1.0,1.1,1.2,1.3,1.4,1.5,1.6,1.7,1.8,1.9,2.0,2.2,2.4,2.6,2.8,3.0,3.2,3.4,3.6,3.8,4.0,4.5,5.0,5.5,6.0,6.5,7.0,7.5,8.0,9.0,10.0,11.0,12.0,13.0,14.0,15.0,16.0,18.0, 20.};\r
93 \r
94   // \r
95   // THnSparse track histograms\r
96   //\r
97 \r
98   // TPC -> ITS matching efficiency\r
99   // eta:phi:pt:isPrim:charge:isMatch:isTPC\r
100   Int_t binsRecMCTrackHistTPCITS[7]=  { 30,  90,             ptNbins, 2,   3,  2,  2  };\r
101   Double_t minRecMCTrackHistTPCITS[7]={-1.5, 0.,             ptMin,   0., -1., 0., 0. };\r
102   Double_t maxRecMCTrackHistTPCITS[7]={ 1.5, 2.*TMath::Pi(), ptMax,   2.,  2., 2., 2. };\r
103 \r
104   fRecMCTrackHistTPCITS = new THnSparseF("fRecMCTrackHistTPCITS","eta:phi:pt:isPrim:charge:isMatch:isTPC",7,binsRecMCTrackHistTPCITS,minRecMCTrackHistTPCITS,maxRecMCTrackHistTPCITS);\r
105   fRecMCTrackHistTPCITS->SetBinEdges(2,binsPt);\r
106   fRecMCTrackHistTPCITS->GetAxis(0)->SetTitle("#eta");\r
107   fRecMCTrackHistTPCITS->GetAxis(1)->SetTitle("#phi (rad)");\r
108   fRecMCTrackHistTPCITS->GetAxis(2)->SetTitle("p_{T} (GeV/c)");\r
109   fRecMCTrackHistTPCITS->GetAxis(3)->SetTitle("isPrim");\r
110   fRecMCTrackHistTPCITS->GetAxis(4)->SetTitle("charge");\r
111   fRecMCTrackHistTPCITS->GetAxis(5)->SetTitle("isMatch");\r
112   fRecMCTrackHistTPCITS->GetAxis(6)->SetTitle("isTPC");\r
113   fRecMCTrackHistTPCITS->Sumw2();\r
114 \r
115   // ITS -> TPC matching efficiency\r
116   // eta:phi:pt:isPrim:charge:isMatch\r
117   Int_t binsRecMCTrackHistITSTPC[6]=  { 30,  90,             ptNbins,   2,   3,  2 };\r
118   Double_t minRecMCTrackHistITSTPC[6]={-1.5, 0.,             ptMin,     0., -1., 0 };\r
119   Double_t maxRecMCTrackHistITSTPC[6]={ 1.5, 2.*TMath::Pi(), ptMax,     2.,  2., 2.};\r
120 \r
121   fRecMCTrackHistITSTPC = new THnSparseF("fRecMCTrackHistITSTPC","eta:phi:pt:isPrim:charge:isMatch",6,binsRecMCTrackHistITSTPC,minRecMCTrackHistITSTPC,maxRecMCTrackHistITSTPC);\r
122   fRecMCTrackHistITSTPC->SetBinEdges(2,binsPt);\r
123   fRecMCTrackHistITSTPC->GetAxis(0)->SetTitle("#eta");\r
124   fRecMCTrackHistITSTPC->GetAxis(1)->SetTitle("#phi (rad)");\r
125   fRecMCTrackHistITSTPC->GetAxis(2)->SetTitle("p_{T} (GeV/c)");\r
126   fRecMCTrackHistITSTPC->GetAxis(3)->SetTitle("isPrim");\r
127   fRecMCTrackHistITSTPC->GetAxis(4)->SetTitle("charge");\r
128   fRecMCTrackHistITSTPC->GetAxis(5)->SetTitle("isMatch");\r
129   fRecMCTrackHistITSTPC->Sumw2();\r
130 \r
131   // init output folder\r
132   fAnalysisFolder = CreateFolder("folderdNdPt","Analysis dNdPt Folder");\r
133 \r
134 }\r
135 \r
136 //_____________________________________________________________________________\r
137 void AlidNdPtEfficiency::Process(AliESDEvent *const esdEvent, AliMCEvent * const mcEvent)\r
138 {\r
139   //\r
140   // Process real and/or simulated events\r
141   //\r
142   if(!esdEvent) {\r
143     AliDebug(AliLog::kError, "esdEvent not available");\r
144     return;\r
145   }\r
146 \r
147   // get selection cuts\r
148   AlidNdPtEventCuts *evtCuts      = GetEventCuts(); \r
149   AlidNdPtAcceptanceCuts *accCuts = GetAcceptanceCuts(); \r
150   AliESDtrackCuts *esdTrackCuts   = GetTrackCuts(); \r
151 \r
152   if(!evtCuts || !accCuts  || !esdTrackCuts) {\r
153     AliDebug(AliLog::kError, "cuts not available");\r
154     return;\r
155   }\r
156 \r
157    \r
158 \r
159   // trigger selection\r
160   Bool_t isEventTriggered = kTRUE;\r
161    \r
162   // use MC information\r
163   AliHeader* header = 0;\r
164   AliGenEventHeader* genHeader = 0;\r
165   AliStack* stack = 0;\r
166   TArrayF vtxMC(3);\r
167 \r
168   Int_t multMCTrueTracks = 0;\r
169   if(IsUseMCInfo())\r
170   {\r
171     //\r
172     if(!mcEvent) {\r
173       AliDebug(AliLog::kError, "mcEvent not available");\r
174       return;\r
175     }\r
176     // get MC event header\r
177     header = mcEvent->Header();\r
178     if (!header) {\r
179       AliDebug(AliLog::kError, "Header not available");\r
180       return;\r
181     }\r
182     // MC particle stack\r
183     stack = mcEvent->Stack();\r
184     if (!stack) {\r
185       AliDebug(AliLog::kError, "Stack not available");\r
186       return;\r
187     }\r
188 \r
189     // get MC vertex\r
190     genHeader = header->GenEventHeader();\r
191     if (!genHeader) {\r
192       AliDebug(AliLog::kError, "Could not retrieve genHeader from Header");\r
193       return;\r
194     }\r
195     genHeader->PrimaryVertex(vtxMC);\r
196 \r
197     // multipliticy of all MC primary tracks\r
198     // in Zv, pt and eta ranges)\r
199     multMCTrueTracks = AlidNdPtHelper::GetMCTrueTrackMult(mcEvent,evtCuts,accCuts);\r
200 \r
201   } // end bUseMC\r
202 \r
203   // get reconstructed vertex  \r
204   const AliESDVertex* vtxESD = 0; \r
205   if(evtCuts->IsRecVertexRequired()) \r
206   {\r
207      if(GetAnalysisMode() == AlidNdPtHelper::kTPC) {\r
208         vtxESD = esdEvent->GetPrimaryVertexTPC();\r
209     }\r
210     else if(GetAnalysisMode() == AlidNdPtHelper::kTPCITS) {\r
211       vtxESD = esdEvent->GetPrimaryVertexTracks();\r
212     }\r
213     else {\r
214         return;\r
215     }\r
216   }\r
217   Bool_t isEventOK = evtCuts->AcceptEvent(esdEvent,mcEvent,vtxESD); \r
218   //printf("isEventOK %d, isEventTriggered %d \n",isEventOK, isEventTriggered);\r
219   //printf("GetAnalysisMode() %d \n",GetAnalysisMode());\r
220 \r
221   TObjArray *allChargedTracks=0;\r
222   //printf("isEventOK %d, isEventTriggered %d \n",isEventOK,isEventTriggered);\r
223 \r
224   // check event cuts\r
225   if(isEventOK && isEventTriggered)\r
226   {\r
227     // get all charged tracks\r
228     allChargedTracks = AlidNdPtHelper::GetAllChargedTracks(esdEvent,GetAnalysisMode());\r
229     if(!allChargedTracks) return;\r
230 \r
231     Int_t entries = allChargedTracks->GetEntries();\r
232     Bool_t isTPC = kFALSE;\r
233     Bool_t isMatch = kFALSE;\r
234 \r
235     // TPC -> ITS prolongation efficiency\r
236     for(Int_t iTrack=0; iTrack<entries;++iTrack) \r
237     {\r
238       AliESDtrack *track = (AliESDtrack*)allChargedTracks->At(iTrack);\r
239       if(!track) continue;\r
240 \r
241       isTPC = kFALSE;\r
242 \r
243       if(track->Charge()==0) continue;\r
244       if(!track->GetTPCInnerParam()) continue;\r
245       if(!(track->GetStatus() & AliESDtrack::kTPCrefit)) continue;\r
246 \r
247       // Get TPC only tracks (must be deleted by user) \r
248       AliESDtrack* tpcTrack = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent,iTrack);\r
249       if(!tpcTrack) continue;\r
250       if(!tpcTrack->RelateToVertex(vtxESD,esdEvent->GetMagneticField(),100.)) { delete tpcTrack; continue; } \r
251 \r
252       // check loose cuts for TPC tracks\r
253       if(!esdTrackCuts->AcceptTrack(tpcTrack))  { delete tpcTrack; continue; } \r
254 \r
255       isTPC = kTRUE;\r
256       isMatch = kFALSE;\r
257       if( (track->GetStatus()&AliESDtrack::kITSrefit) && \r
258           (track->HasPointOnITSLayer(0) || track->HasPointOnITSLayer(1)) ) \r
259       {\r
260         isMatch = kTRUE;\r
261       }\r
262 \r
263       //\r
264       FillHistograms(tpcTrack, stack, isMatch, isTPC, kFALSE);\r
265       if(tpcTrack) delete tpcTrack;\r
266     } \r
267 \r
268     //\r
269     // ITS -> TPC prolongation efficiency\r
270     //\r
271     for(Int_t iTrack=0; iTrack<entries;++iTrack) \r
272     {\r
273       AliESDtrack *track = (AliESDtrack*)allChargedTracks->At(iTrack);\r
274       if(!track) continue;\r
275 \r
276       //\r
277       // ITS stand alone tracks\r
278       //\r
279       if(!(track->GetStatus() & AliESDtrack::kITSpureSA)) continue;\r
280       if(!(track->GetStatus() & AliESDtrack::kITSrefit)) continue;\r
281       if(track->GetNcls(0)<4) continue;\r
282       if(!track->HasPointOnITSLayer(0) && !track->HasPointOnITSLayer(1)) continue;\r
283 \r
284       // Check matching with TPC only track\r
285       for(Int_t jTrack=0; jTrack<entries;++jTrack) \r
286       {\r
287         isMatch = kFALSE;\r
288 \r
289         if(iTrack==jTrack) continue;\r
290 \r
291         AliESDtrack *track2 = (AliESDtrack*)allChargedTracks->At(jTrack);\r
292         if(!track2) continue;\r
293         if(track2->Charge()==0) continue;\r
294         if(!track2->GetTPCInnerParam()) continue;\r
295         if(!(track2->GetStatus() & AliESDtrack::kTPCrefit)) continue;\r
296 \r
297         // Get TPC only tracks (must be deleted by user) \r
298         AliESDtrack* tpcTrack2 = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent, jTrack);\r
299         if(!tpcTrack2) continue;\r
300         if(!tpcTrack2->RelateToVertex(vtxESD,esdEvent->GetMagneticField(),100.)) { delete tpcTrack2; continue; } \r
301 \r
302         // check loose cuts for TPC tracks\r
303         if(!esdTrackCuts->AcceptTrack(tpcTrack2)) { delete tpcTrack2; continue; }\r
304 \r
305         // check matching\r
306         if (TMath::Abs(track->GetY() - tpcTrack2->GetY()) > 3) { delete tpcTrack2; continue; }\r
307         if (TMath::Abs(track->GetSnp() - tpcTrack2->GetSnp()) > 0.2) { delete tpcTrack2; continue; }\r
308         if (TMath::Abs(track->GetTgl() - tpcTrack2->GetTgl()) > 0.2) { delete tpcTrack2; continue; }\r
309 \r
310         isMatch = kTRUE;\r
311         if(tpcTrack2) { \r
312           delete tpcTrack2;\r
313         }\r
314           break;\r
315       }\r
316 \r
317        //\r
318        FillHistograms(track, stack, isMatch, kFALSE, kTRUE);\r
319     } \r
320   }\r
321 \r
322   if(allChargedTracks) delete allChargedTracks; allChargedTracks = 0;\r
323 \r
324 }\r
325 \r
326 //_____________________________________________________________________________\r
327 void AlidNdPtEfficiency::FillHistograms(AliESDtrack *const esdTrack, AliStack *const stack, const Bool_t isMatch, const Bool_t isTPC,  const Bool_t isITSTPC) const\r
328 {\r
329   //\r
330   // Fill ESD track and MC histograms \r
331   //\r
332   if(!esdTrack) return;\r
333   Int_t charge = esdTrack->Charge();\r
334   if(charge == 0.) return;\r
335 \r
336   Float_t pt = esdTrack->Pt();\r
337   Float_t eta = esdTrack->Eta();\r
338   Float_t phi = esdTrack->Phi();\r
339 \r
340   //\r
341   // Fill rec vs MC information\r
342   //\r
343   Bool_t isPrim = kTRUE;\r
344 \r
345   if(IsUseMCInfo()) {\r
346     if(!stack) return;\r
347     Int_t label = esdTrack->GetLabel(); \r
348     if(label < 0.) return; // fake ITS track\r
349     TParticle* particle = stack->Particle(label);\r
350     if(!particle) return;\r
351     if(particle->GetPDG() && particle->GetPDG()->Charge()==0.) return;\r
352     isPrim = stack->IsPhysicalPrimary(label);\r
353   }\r
354 \r
355   // fill histo\r
356   Double_t vRecMCTrackHist[6] = { eta,phi,pt,isPrim,charge,isMatch }; \r
357   Double_t vRecMCTrackHistTPCITS[7] = { eta,phi,pt,isPrim,charge,isMatch,isTPC }; \r
358 \r
359   if(isITSTPC) {\r
360     fRecMCTrackHistITSTPC->Fill(vRecMCTrackHist);\r
361   }\r
362   else {\r
363     fRecMCTrackHistTPCITS->Fill(vRecMCTrackHistTPCITS);\r
364   }\r
365 }\r
366 \r
367 //_____________________________________________________________________________\r
368 Long64_t AlidNdPtEfficiency::Merge(TCollection* const list) \r
369 {\r
370   // Merge list of objects (needed by PROOF)\r
371 \r
372   if (!list)\r
373   return 0;\r
374 \r
375   if (list->IsEmpty())\r
376   return 1;\r
377 \r
378   TIterator* iter = list->MakeIterator();\r
379   TObject* obj = 0;\r
380 \r
381   // collection of generated histograms\r
382   Int_t count=0;\r
383   while((obj = iter->Next()) != 0) {\r
384     AlidNdPtEfficiency* entry = dynamic_cast<AlidNdPtEfficiency*>(obj);\r
385     if (entry == 0) continue; \r
386   \r
387     // track histo\r
388     fRecMCTrackHistTPCITS->Add(entry->fRecMCTrackHistTPCITS);\r
389     fRecMCTrackHistITSTPC->Add(entry->fRecMCTrackHistITSTPC);\r
390 \r
391   count++;\r
392   }\r
393 \r
394 return count;\r
395 }\r
396 \r
397 //_____________________________________________________________________________\r
398 void AlidNdPtEfficiency::Analyse() \r
399 {\r
400   //\r
401   // Analyse histograms\r
402   //\r
403   TH1::AddDirectory(kFALSE);\r
404   TObjArray *aFolderObj = new TObjArray;\r
405   TH1D *h1Dall = 0; \r
406   TH1D *h1D = 0; \r
407   TH1D *h1Dc = 0; \r
408 \r
409 \r
410   //\r
411   // get cuts\r
412   //\r
413   AlidNdPtEventCuts *evtCuts = GetEventCuts(); \r
414   AlidNdPtAcceptanceCuts *accCuts = GetAcceptanceCuts(); \r
415   AliESDtrackCuts *esdTrackCuts = GetTrackCuts(); \r
416 \r
417   if(!evtCuts || !accCuts || !esdTrackCuts) {\r
418     Error("AlidNdPtEfficiency::Analyse()", "cuts not available");\r
419     return;\r
420   }\r
421 \r
422   //\r
423   // TPC->ITS efficiency\r
424   //\r
425 \r
426   //eff vs eta\r
427   fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(2,2);  \r
428   h1Dall = (TH1D *)fRecMCTrackHistTPCITS->Projection(0);\r
429   fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(2,2);  \r
430   h1D = (TH1D *)fRecMCTrackHistTPCITS->Projection(0);\r
431 \r
432   h1Dc = (TH1D *)h1D->Clone("eff_vs_eta_TPCITS");\r
433   h1Dc->Divide(h1Dall);\r
434   aFolderObj->Add(h1Dc);\r
435   fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(1,2);  \r
436   fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(1,2);  \r
437 \r
438   //eff vs phi\r
439   fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(2,2);  \r
440   fRecMCTrackHistTPCITS->GetAxis(0)->SetRangeUser(-0.8, 0.799);  \r
441   h1Dall = (TH1D *)fRecMCTrackHistTPCITS->Projection(1);\r
442   fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(2,2);  \r
443   h1D = (TH1D *)fRecMCTrackHistTPCITS->Projection(1);\r
444 \r
445   h1Dc = (TH1D *)h1D->Clone("eff_vs_phi_TPCITS");\r
446   h1Dc->Divide(h1Dall);\r
447   aFolderObj->Add(h1Dc);\r
448   fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(1,2);  \r
449   fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(1,2);  \r
450 \r
451   //eff vs pT\r
452   fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(2,2);  \r
453   fRecMCTrackHistTPCITS->GetAxis(0)->SetRangeUser(-0.8, 0.799);  \r
454   h1Dall = (TH1D *)fRecMCTrackHistTPCITS->Projection(2);\r
455   fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(2,2);  \r
456   h1D = (TH1D *)fRecMCTrackHistTPCITS->Projection(2);\r
457 \r
458   h1Dc = (TH1D *)h1D->Clone("eff_vs_pT_TPCITS");\r
459   h1Dc->Divide(h1Dall);\r
460   aFolderObj->Add(h1Dc);\r
461   fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(1,2);  \r
462   fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(1,2);  \r
463 \r
464 \r
465   //\r
466   // ITS->TPC efficiency\r
467   //\r
468 \r
469   fRecMCTrackHistITSTPC->GetAxis(0)->SetRangeUser(-1.5, 1.499);  \r
470   fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(1,2);  \r
471 \r
472   //eff vs eta\r
473   h1Dall = (TH1D *)fRecMCTrackHistITSTPC->Projection(0);\r
474   fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(2,2);  \r
475   h1D = (TH1D *)fRecMCTrackHistITSTPC->Projection(0);\r
476 \r
477   h1Dc = (TH1D *)h1D->Clone("eff_vs_eta_ITSTPC");\r
478   h1Dc->Divide(h1Dall);\r
479   aFolderObj->Add(h1Dc);\r
480   fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(1,2);  \r
481 \r
482   //eff vs phi\r
483   fRecMCTrackHistITSTPC->GetAxis(0)->SetRangeUser(-0.8, 0.799);  \r
484   h1Dall = (TH1D *)fRecMCTrackHistITSTPC->Projection(1);\r
485   fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(2,2);  \r
486   h1D = (TH1D *)fRecMCTrackHistITSTPC->Projection(1);\r
487 \r
488   h1Dc = (TH1D *)h1D->Clone("eff_vs_phi_ITSTPC");\r
489   h1Dc->Divide(h1Dall);\r
490   aFolderObj->Add(h1Dc);\r
491   fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(1,2);  \r
492 \r
493   //eff vs pT\r
494   fRecMCTrackHistITSTPC->GetAxis(0)->SetRangeUser(-0.8, 0.799);  \r
495   h1Dall = (TH1D *)fRecMCTrackHistITSTPC->Projection(2);\r
496   fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(2,2);  \r
497   h1D = (TH1D *)fRecMCTrackHistITSTPC->Projection(2);\r
498 \r
499   h1Dc = (TH1D *)h1D->Clone("eff_vs_pT_ITSTPC");\r
500   h1Dc->Divide(h1Dall);\r
501   aFolderObj->Add(h1Dc);\r
502   fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(1,2);  \r
503   \r
504   // export objects to analysis folder\r
505   fAnalysisFolder = ExportToFolder(aFolderObj);\r
506 \r
507   // delete only TObjArray\r
508   if(aFolderObj) delete aFolderObj;\r
509 }\r
510 \r
511 //_____________________________________________________________________________\r
512 TFolder* AlidNdPtEfficiency::ExportToFolder(TObjArray * const array) \r
513 {\r
514   // recreate folder avery time and export objects to new one\r
515   //\r
516   AlidNdPtEfficiency * comp=this;\r
517   TFolder *folder = comp->GetAnalysisFolder();\r
518 \r
519   TString name, title;\r
520   TFolder *newFolder = 0;\r
521   Int_t i = 0;\r
522   Int_t size = array->GetSize();\r
523 \r
524   if(folder) { \r
525      // get name and title from old folder\r
526      name = folder->GetName();  \r
527      title = folder->GetTitle();  \r
528 \r
529          // delete old one\r
530      delete folder;\r
531 \r
532          // create new one\r
533      newFolder = CreateFolder(name.Data(),title.Data());\r
534      newFolder->SetOwner();\r
535 \r
536          // add objects to folder\r
537      while(i < size) {\r
538            newFolder->Add(array->At(i));\r
539            i++;\r
540          }\r
541   }\r
542 \r
543 return newFolder;\r
544 }\r
545 \r
546 //_____________________________________________________________________________\r
547 TFolder* AlidNdPtEfficiency::CreateFolder(TString name,TString title) { \r
548 // create folder for analysed histograms\r
549 //\r
550 TFolder *folder = 0;\r
551   folder = new TFolder(name.Data(),title.Data());\r
552 \r
553   return folder;\r
554 }\r