]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG0/dNdPt/AlidNdPtEfficiency.cxx
protection agains wrong labels
[u/mrichter/AliRoot.git] / PWG0 / dNdPt / AlidNdPtEfficiency.cxx
CommitLineData
77b5a2bd 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
50using namespace std;\r
51\r
52ClassImp(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
65AlidNdPtEfficiency::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
75AlidNdPtEfficiency::~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
84void AlidNdPtEfficiency::Init(){\r
85 //\r
86 // Init histograms\r
87 //\r
98af3b87 88 const Int_t ptNbins = 63; \r
77b5a2bd 89 const Double_t ptMin = 0.; \r
90 const Double_t ptMax = 20.; \r
91\r
98af3b87 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.,25.,30.,35.,40.,50};\r
77b5a2bd 93\r
94 // \r
95 // THnSparse track histograms\r
96 //\r
97\r
98 // TPC -> ITS matching efficiency\r
55468faf 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
77b5a2bd 103\r
55468faf 104 fRecMCTrackHistTPCITS = new THnSparseF("fRecMCTrackHistTPCITS","eta:phi:pt:isPrim:charge:isMatch:isTPC",7,binsRecMCTrackHistTPCITS,minRecMCTrackHistTPCITS,maxRecMCTrackHistTPCITS);\r
77b5a2bd 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
55468faf 112 fRecMCTrackHistTPCITS->GetAxis(6)->SetTitle("isTPC");\r
77b5a2bd 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
77b5a2bd 133}\r
134\r
135//_____________________________________________________________________________\r
136void AlidNdPtEfficiency::Process(AliESDEvent *const esdEvent, AliMCEvent * const mcEvent)\r
137{\r
138 //\r
139 // Process real and/or simulated events\r
140 //\r
141 if(!esdEvent) {\r
142 AliDebug(AliLog::kError, "esdEvent not available");\r
143 return;\r
144 }\r
145\r
146 // get selection cuts\r
55468faf 147 AlidNdPtEventCuts *evtCuts = GetEventCuts(); \r
77b5a2bd 148 AlidNdPtAcceptanceCuts *accCuts = GetAcceptanceCuts(); \r
55468faf 149 AliESDtrackCuts *esdTrackCuts = GetTrackCuts(); \r
77b5a2bd 150\r
151 if(!evtCuts || !accCuts || !esdTrackCuts) {\r
152 AliDebug(AliLog::kError, "cuts not available");\r
153 return;\r
154 }\r
155\r
55468faf 156 \r
157\r
77b5a2bd 158 // trigger selection\r
159 Bool_t isEventTriggered = kTRUE;\r
02f92c0c 160 \r
77b5a2bd 161 // use MC information\r
162 AliHeader* header = 0;\r
163 AliGenEventHeader* genHeader = 0;\r
164 AliStack* stack = 0;\r
165 TArrayF vtxMC(3);\r
166\r
167 Int_t multMCTrueTracks = 0;\r
168 if(IsUseMCInfo())\r
169 {\r
170 //\r
171 if(!mcEvent) {\r
172 AliDebug(AliLog::kError, "mcEvent not available");\r
173 return;\r
174 }\r
175 // get MC event header\r
176 header = mcEvent->Header();\r
177 if (!header) {\r
178 AliDebug(AliLog::kError, "Header not available");\r
179 return;\r
180 }\r
181 // MC particle stack\r
182 stack = mcEvent->Stack();\r
183 if (!stack) {\r
184 AliDebug(AliLog::kError, "Stack not available");\r
185 return;\r
186 }\r
187\r
188 // get MC vertex\r
189 genHeader = header->GenEventHeader();\r
190 if (!genHeader) {\r
191 AliDebug(AliLog::kError, "Could not retrieve genHeader from Header");\r
192 return;\r
193 }\r
194 genHeader->PrimaryVertex(vtxMC);\r
195\r
196 // multipliticy of all MC primary tracks\r
197 // in Zv, pt and eta ranges)\r
198 multMCTrueTracks = AlidNdPtHelper::GetMCTrueTrackMult(mcEvent,evtCuts,accCuts);\r
199\r
200 } // end bUseMC\r
201\r
202 // get reconstructed vertex \r
203 const AliESDVertex* vtxESD = 0; \r
204 if(evtCuts->IsRecVertexRequired()) \r
205 {\r
206 if(GetAnalysisMode() == AlidNdPtHelper::kTPC) {\r
207 vtxESD = esdEvent->GetPrimaryVertexTPC();\r
208 }\r
209 else if(GetAnalysisMode() == AlidNdPtHelper::kTPCITS) {\r
210 vtxESD = esdEvent->GetPrimaryVertexTracks();\r
211 }\r
212 else {\r
213 return;\r
214 }\r
215 }\r
216 Bool_t isEventOK = evtCuts->AcceptEvent(esdEvent,mcEvent,vtxESD); \r
217 //printf("isEventOK %d, isEventTriggered %d \n",isEventOK, isEventTriggered);\r
218 //printf("GetAnalysisMode() %d \n",GetAnalysisMode());\r
219\r
220 TObjArray *allChargedTracks=0;\r
221 //printf("isEventOK %d, isEventTriggered %d \n",isEventOK,isEventTriggered);\r
222\r
223 // check event cuts\r
224 if(isEventOK && isEventTriggered)\r
225 {\r
226 // get all charged tracks\r
227 allChargedTracks = AlidNdPtHelper::GetAllChargedTracks(esdEvent,GetAnalysisMode());\r
228 if(!allChargedTracks) return;\r
229\r
77b5a2bd 230 Int_t entries = allChargedTracks->GetEntries();\r
55468faf 231 Bool_t isTPC = kFALSE;\r
232 Bool_t isMatch = kFALSE;\r
77b5a2bd 233\r
234 // TPC -> ITS prolongation efficiency\r
235 for(Int_t iTrack=0; iTrack<entries;++iTrack) \r
236 {\r
237 AliESDtrack *track = (AliESDtrack*)allChargedTracks->At(iTrack);\r
238 if(!track) continue;\r
239\r
55468faf 240 isTPC = kFALSE;\r
241\r
77b5a2bd 242 if(track->Charge()==0) continue;\r
55468faf 243 if(!track->GetTPCInnerParam()) continue;\r
98af3b87 244 if(!(track->GetStatus()&AliESDtrack::kTPCrefit)) continue;\r
77b5a2bd 245\r
77b5a2bd 246 // check loose cuts for TPC tracks\r
98af3b87 247 if(!esdTrackCuts->AcceptTrack(track)) { continue; } \r
77b5a2bd 248\r
55468faf 249 isTPC = kTRUE;\r
250 isMatch = kFALSE;\r
251 if( (track->GetStatus()&AliESDtrack::kITSrefit) && \r
252 (track->HasPointOnITSLayer(0) || track->HasPointOnITSLayer(1)) ) \r
77b5a2bd 253 {\r
254 isMatch = kTRUE;\r
255 }\r
256\r
257 //\r
98af3b87 258 FillHistograms(track, stack, isMatch, isTPC, kFALSE);\r
259 //if(tpcTrack) delete tpcTrack;\r
77b5a2bd 260 } \r
261\r
262 //\r
263 // ITS -> TPC prolongation efficiency\r
264 //\r
265 for(Int_t iTrack=0; iTrack<entries;++iTrack) \r
266 {\r
267 AliESDtrack *track = (AliESDtrack*)allChargedTracks->At(iTrack);\r
268 if(!track) continue;\r
269\r
270 //\r
271 // ITS stand alone tracks\r
272 //\r
273 if(!(track->GetStatus() & AliESDtrack::kITSpureSA)) continue;\r
274 if(!(track->GetStatus() & AliESDtrack::kITSrefit)) continue;\r
275 if(track->GetNcls(0)<4) continue;\r
276 if(!track->HasPointOnITSLayer(0) && !track->HasPointOnITSLayer(1)) continue;\r
277\r
278 // Check matching with TPC only track\r
77b5a2bd 279 for(Int_t jTrack=0; jTrack<entries;++jTrack) \r
280 {\r
55468faf 281 isMatch = kFALSE;\r
282\r
77b5a2bd 283 if(iTrack==jTrack) continue;\r
284\r
285 AliESDtrack *track2 = (AliESDtrack*)allChargedTracks->At(jTrack);\r
286 if(!track2) continue;\r
287 if(track2->Charge()==0) continue;\r
55468faf 288 if(!track2->GetTPCInnerParam()) continue;\r
289 if(!(track2->GetStatus() & AliESDtrack::kTPCrefit)) continue;\r
77b5a2bd 290\r
291 // Get TPC only tracks (must be deleted by user) \r
292 AliESDtrack* tpcTrack2 = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent, jTrack);\r
293 if(!tpcTrack2) continue;\r
55468faf 294 if(!tpcTrack2->RelateToVertex(vtxESD,esdEvent->GetMagneticField(),100.)) { delete tpcTrack2; continue; } \r
77b5a2bd 295\r
296 // check loose cuts for TPC tracks\r
55468faf 297 if(!esdTrackCuts->AcceptTrack(tpcTrack2)) { delete tpcTrack2; continue; }\r
77b5a2bd 298\r
299 // check matching\r
300 if (TMath::Abs(track->GetY() - tpcTrack2->GetY()) > 3) { delete tpcTrack2; continue; }\r
301 if (TMath::Abs(track->GetSnp() - tpcTrack2->GetSnp()) > 0.2) { delete tpcTrack2; continue; }\r
302 if (TMath::Abs(track->GetTgl() - tpcTrack2->GetTgl()) > 0.2) { delete tpcTrack2; continue; }\r
303\r
304 isMatch = kTRUE;\r
305 if(tpcTrack2) { \r
306 delete tpcTrack2;\r
307 }\r
308 break;\r
309 }\r
310\r
311 //\r
55468faf 312 FillHistograms(track, stack, isMatch, kFALSE, kTRUE);\r
77b5a2bd 313 } \r
314 }\r
315\r
316 if(allChargedTracks) delete allChargedTracks; allChargedTracks = 0;\r
317\r
318}\r
319\r
320//_____________________________________________________________________________\r
55468faf 321void AlidNdPtEfficiency::FillHistograms(AliESDtrack *const esdTrack, AliStack *const stack, const Bool_t isMatch, const Bool_t isTPC, const Bool_t isITSTPC) const\r
77b5a2bd 322{\r
323 //\r
324 // Fill ESD track and MC histograms \r
325 //\r
326 if(!esdTrack) return;\r
327 Int_t charge = esdTrack->Charge();\r
328 if(charge == 0.) return;\r
329\r
330 Float_t pt = esdTrack->Pt();\r
331 Float_t eta = esdTrack->Eta();\r
332 Float_t phi = esdTrack->Phi();\r
333\r
334 //\r
335 // Fill rec vs MC information\r
336 //\r
337 Bool_t isPrim = kTRUE;\r
338\r
339 if(IsUseMCInfo()) {\r
340 if(!stack) return;\r
341 Int_t label = esdTrack->GetLabel(); \r
342 if(label < 0.) return; // fake ITS track\r
343 TParticle* particle = stack->Particle(label);\r
344 if(!particle) return;\r
345 if(particle->GetPDG() && particle->GetPDG()->Charge()==0.) return;\r
346 isPrim = stack->IsPhysicalPrimary(label);\r
347 }\r
348\r
349 // fill histo\r
350 Double_t vRecMCTrackHist[6] = { eta,phi,pt,isPrim,charge,isMatch }; \r
55468faf 351 Double_t vRecMCTrackHistTPCITS[7] = { eta,phi,pt,isPrim,charge,isMatch,isTPC }; \r
77b5a2bd 352\r
353 if(isITSTPC) {\r
354 fRecMCTrackHistITSTPC->Fill(vRecMCTrackHist);\r
355 }\r
356 else {\r
55468faf 357 fRecMCTrackHistTPCITS->Fill(vRecMCTrackHistTPCITS);\r
77b5a2bd 358 }\r
359}\r
360\r
361//_____________________________________________________________________________\r
362Long64_t AlidNdPtEfficiency::Merge(TCollection* const list) \r
363{\r
364 // Merge list of objects (needed by PROOF)\r
365\r
366 if (!list)\r
367 return 0;\r
368\r
369 if (list->IsEmpty())\r
370 return 1;\r
371\r
372 TIterator* iter = list->MakeIterator();\r
373 TObject* obj = 0;\r
374\r
375 // collection of generated histograms\r
376 Int_t count=0;\r
377 while((obj = iter->Next()) != 0) {\r
378 AlidNdPtEfficiency* entry = dynamic_cast<AlidNdPtEfficiency*>(obj);\r
379 if (entry == 0) continue; \r
380 \r
381 // track histo\r
382 fRecMCTrackHistTPCITS->Add(entry->fRecMCTrackHistTPCITS);\r
383 fRecMCTrackHistITSTPC->Add(entry->fRecMCTrackHistITSTPC);\r
384\r
385 count++;\r
386 }\r
387\r
388return count;\r
389}\r
390\r
391//_____________________________________________________________________________\r
392void AlidNdPtEfficiency::Analyse() \r
393{\r
394 //\r
395 // Analyse histograms\r
396 //\r
397 TH1::AddDirectory(kFALSE);\r
398 TObjArray *aFolderObj = new TObjArray;\r
399 TH1D *h1Dall = 0; \r
400 TH1D *h1D = 0; \r
401 TH1D *h1Dc = 0; \r
402\r
403\r
404 //\r
405 // get cuts\r
406 //\r
407 AlidNdPtEventCuts *evtCuts = GetEventCuts(); \r
408 AlidNdPtAcceptanceCuts *accCuts = GetAcceptanceCuts(); \r
409 AliESDtrackCuts *esdTrackCuts = GetTrackCuts(); \r
410\r
411 if(!evtCuts || !accCuts || !esdTrackCuts) {\r
412 Error("AlidNdPtEfficiency::Analyse()", "cuts not available");\r
413 return;\r
414 }\r
415\r
416 //\r
417 // TPC->ITS efficiency\r
418 //\r
419\r
420 //eff vs eta\r
55468faf 421 fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(2,2); \r
77b5a2bd 422 h1Dall = (TH1D *)fRecMCTrackHistTPCITS->Projection(0);\r
423 fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(2,2); \r
424 h1D = (TH1D *)fRecMCTrackHistTPCITS->Projection(0);\r
425\r
426 h1Dc = (TH1D *)h1D->Clone("eff_vs_eta_TPCITS");\r
427 h1Dc->Divide(h1Dall);\r
428 aFolderObj->Add(h1Dc);\r
429 fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(1,2); \r
55468faf 430 fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(1,2); \r
77b5a2bd 431\r
432 //eff vs phi\r
55468faf 433 fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(2,2); \r
77b5a2bd 434 fRecMCTrackHistTPCITS->GetAxis(0)->SetRangeUser(-0.8, 0.799); \r
435 h1Dall = (TH1D *)fRecMCTrackHistTPCITS->Projection(1);\r
436 fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(2,2); \r
437 h1D = (TH1D *)fRecMCTrackHistTPCITS->Projection(1);\r
438\r
439 h1Dc = (TH1D *)h1D->Clone("eff_vs_phi_TPCITS");\r
440 h1Dc->Divide(h1Dall);\r
441 aFolderObj->Add(h1Dc);\r
442 fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(1,2); \r
55468faf 443 fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(1,2); \r
77b5a2bd 444\r
445 //eff vs pT\r
55468faf 446 fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(2,2); \r
77b5a2bd 447 fRecMCTrackHistTPCITS->GetAxis(0)->SetRangeUser(-0.8, 0.799); \r
448 h1Dall = (TH1D *)fRecMCTrackHistTPCITS->Projection(2);\r
449 fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(2,2); \r
450 h1D = (TH1D *)fRecMCTrackHistTPCITS->Projection(2);\r
451\r
452 h1Dc = (TH1D *)h1D->Clone("eff_vs_pT_TPCITS");\r
453 h1Dc->Divide(h1Dall);\r
454 aFolderObj->Add(h1Dc);\r
455 fRecMCTrackHistTPCITS->GetAxis(5)->SetRange(1,2); \r
55468faf 456 fRecMCTrackHistTPCITS->GetAxis(6)->SetRange(1,2); \r
77b5a2bd 457\r
458\r
459 //\r
55468faf 460 // ITS->TPC efficiency\r
77b5a2bd 461 //\r
462\r
463 fRecMCTrackHistITSTPC->GetAxis(0)->SetRangeUser(-1.5, 1.499); \r
464 fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(1,2); \r
465\r
466 //eff vs eta\r
467 h1Dall = (TH1D *)fRecMCTrackHistITSTPC->Projection(0);\r
468 fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(2,2); \r
469 h1D = (TH1D *)fRecMCTrackHistITSTPC->Projection(0);\r
470\r
471 h1Dc = (TH1D *)h1D->Clone("eff_vs_eta_ITSTPC");\r
472 h1Dc->Divide(h1Dall);\r
473 aFolderObj->Add(h1Dc);\r
474 fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(1,2); \r
475\r
476 //eff vs phi\r
477 fRecMCTrackHistITSTPC->GetAxis(0)->SetRangeUser(-0.8, 0.799); \r
478 h1Dall = (TH1D *)fRecMCTrackHistITSTPC->Projection(1);\r
479 fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(2,2); \r
480 h1D = (TH1D *)fRecMCTrackHistITSTPC->Projection(1);\r
481\r
482 h1Dc = (TH1D *)h1D->Clone("eff_vs_phi_ITSTPC");\r
483 h1Dc->Divide(h1Dall);\r
484 aFolderObj->Add(h1Dc);\r
485 fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(1,2); \r
486\r
487 //eff vs pT\r
488 fRecMCTrackHistITSTPC->GetAxis(0)->SetRangeUser(-0.8, 0.799); \r
489 h1Dall = (TH1D *)fRecMCTrackHistITSTPC->Projection(2);\r
490 fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(2,2); \r
491 h1D = (TH1D *)fRecMCTrackHistITSTPC->Projection(2);\r
492\r
493 h1Dc = (TH1D *)h1D->Clone("eff_vs_pT_ITSTPC");\r
494 h1Dc->Divide(h1Dall);\r
495 aFolderObj->Add(h1Dc);\r
496 fRecMCTrackHistITSTPC->GetAxis(5)->SetRange(1,2); \r
497 \r
498 // export objects to analysis folder\r
499 fAnalysisFolder = ExportToFolder(aFolderObj);\r
500\r
501 // delete only TObjArray\r
502 if(aFolderObj) delete aFolderObj;\r
503}\r
504\r
505//_____________________________________________________________________________\r
506TFolder* AlidNdPtEfficiency::ExportToFolder(TObjArray * const array) \r
507{\r
508 // recreate folder avery time and export objects to new one\r
509 //\r
510 AlidNdPtEfficiency * comp=this;\r
511 TFolder *folder = comp->GetAnalysisFolder();\r
512\r
513 TString name, title;\r
514 TFolder *newFolder = 0;\r
515 Int_t i = 0;\r
516 Int_t size = array->GetSize();\r
517\r
518 if(folder) { \r
519 // get name and title from old folder\r
520 name = folder->GetName(); \r
521 title = folder->GetTitle(); \r
522\r
523 // delete old one\r
524 delete folder;\r
525\r
526 // create new one\r
527 newFolder = CreateFolder(name.Data(),title.Data());\r
528 newFolder->SetOwner();\r
529\r
530 // add objects to folder\r
531 while(i < size) {\r
532 newFolder->Add(array->At(i));\r
533 i++;\r
534 }\r
535 }\r
536\r
537return newFolder;\r
538}\r
539\r
540//_____________________________________________________________________________\r
541TFolder* AlidNdPtEfficiency::CreateFolder(TString name,TString title) { \r
542// create folder for analysed histograms\r
543//\r
544TFolder *folder = 0;\r
545 folder = new TFolder(name.Data(),title.Data());\r
546\r
547 return folder;\r
548}\r