]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/EvTrkSelection/AliCFSingleTrackEfficiencyTask.cxx
Fix for ESD analysis
[u/mrichter/AliRoot.git] / PWGPP / EvTrkSelection / AliCFSingleTrackEfficiencyTask.cxx
1 #include "TCanvas.h"
2 #include "TParticle.h"
3 #include "TH1I.h"
4 #include <TDatabasePDG.h>
5
6 #include "AliAnalysisDataSlot.h"
7 #include "AliAnalysisDataContainer.h"
8 #include "AliAnalysisManager.h"
9 #include "AliAODEvent.h"
10 #include "AliAODMCParticle.h"
11 #include "AliAODTrack.h"
12 #include "AliAODTracklets.h"
13 #include "AliMultiplicity.h"
14 #include "AliCFManager.h"
15 #include "AliCFCutBase.h"
16 #include "AliCFContainer.h"
17 #include "AliESDEvent.h"
18 #include "AliESDtrack.h"
19 #include "AliESDtrackCuts.h"
20 #include "AliESDVertex.h"
21 #include "AliInputEventHandler.h"
22 #include "AliMCEvent.h"
23 #include "AliVEvent.h"
24 #include "AliAODEvent.h"
25 #include "AliVVertex.h"
26 #include "AliLog.h"
27 #include "AliStack.h"
28 #include "AliGenEventHeader.h"
29 #include "AliAODMCHeader.h"
30
31 #include "AliSingleTrackEffCuts.h"
32 #include "AliCFSingleTrackEfficiencyTask.h"
33
34
35 ClassImp(AliCFSingleTrackEfficiencyTask)
36
37 //__________________________________________________________________________
38 AliCFSingleTrackEfficiencyTask::AliCFSingleTrackEfficiencyTask() :
39   fReadAODData(0),
40   fCFManager(0x0),
41   fQAHistList(0x0),
42   fTrackCuts(0x0),
43   fMCCuts(0x0),
44   fTriggerMask(AliVEvent::kAny),
45   fSetFilterBit(kFALSE),
46   fbit(0),
47   fHistEventsProcessed(0x0)
48 {
49   //
50   //Default constructor
51   //
52 }
53
54 //___________________________________________________________________________
55 AliCFSingleTrackEfficiencyTask::AliCFSingleTrackEfficiencyTask(const Char_t* name,AliESDtrackCuts *trackcuts, AliSingleTrackEffCuts * mccuts) :
56   AliAnalysisTaskSE(name),
57   fReadAODData(0),
58   fCFManager(0x0),
59   fQAHistList(0x0),
60   fTrackCuts(trackcuts),
61   fMCCuts(mccuts),
62   fTriggerMask(AliVEvent::kAny),
63   fSetFilterBit(kFALSE),
64   fbit(0),
65   fHistEventsProcessed(0x0)
66 {
67   //
68   // Constructor. Initialization of Inputs and Outputs
69   //
70   Info("AliCFSingleTrackEfficiencyTask","Calling Constructor");
71
72   // Output slot #1 writes into a TList container (nevents histogran)
73   DefineOutput(1,TH1I::Class());
74   // Output slot #2 writes into a TList container (distributions)
75   DefineOutput(2,AliCFContainer::Class());
76   // Output slot #3 writes the QA list
77   DefineOutput(3,TList::Class());
78   // Output slot #3 writes the ESD track cuts used
79   DefineOutput(4,AliESDtrackCuts::Class());
80   // Output slot #4 writes the particle and event selection object
81   DefineOutput(5,AliSingleTrackEffCuts::Class());
82 }
83
84 //_________________________________________________________________________________________________________________
85 AliCFSingleTrackEfficiencyTask& AliCFSingleTrackEfficiencyTask::operator=(const AliCFSingleTrackEfficiencyTask& c) 
86 {
87   //
88   // Assignment operator
89   //
90   if (this!=&c) {
91     AliAnalysisTaskSE::operator=(c) ;
92
93     fReadAODData = c.fReadAODData ;
94     fCFManager  = c.fCFManager;
95     fQAHistList = c.fQAHistList ;
96
97     if(c.fTrackCuts) { delete fTrackCuts; fTrackCuts = new AliESDtrackCuts(*(c.fTrackCuts)); }
98     if(c.fMCCuts) { delete fMCCuts; fMCCuts = new AliSingleTrackEffCuts(*(c.fMCCuts)); }
99     fTriggerMask = c.fTriggerMask;
100
101     fSetFilterBit  = c.fSetFilterBit;
102     fbit = c.fbit;
103     fHistEventsProcessed = c.fHistEventsProcessed;
104   }
105   return *this;
106 }
107
108 //________________________________________________________________________________________________________
109 AliCFSingleTrackEfficiencyTask::AliCFSingleTrackEfficiencyTask(const AliCFSingleTrackEfficiencyTask& c) :
110   AliAnalysisTaskSE(c),
111   fReadAODData(c.fReadAODData),
112   fCFManager(c.fCFManager),
113   fQAHistList(c.fQAHistList),
114   fTrackCuts(c.fTrackCuts),
115   fMCCuts(c.fMCCuts),
116   fTriggerMask(c.fTriggerMask),
117   fSetFilterBit(c.fSetFilterBit),
118   fbit(c.fbit),
119   fHistEventsProcessed(c.fHistEventsProcessed)
120 {
121   //
122   // Copy Constructor
123   //
124 }
125
126 //___________________________________________________________________________
127 AliCFSingleTrackEfficiencyTask::~AliCFSingleTrackEfficiencyTask()
128 {
129   //
130   // Destructor
131   //
132   Info("~AliCFSingleTrackEfficiencyTask","Calling Destructor");
133
134   if (fCFManager)           delete fCFManager;
135   if (fHistEventsProcessed) delete fHistEventsProcessed;
136   if (fQAHistList) { fQAHistList->Clear(); delete fQAHistList; }
137   if (fTrackCuts)           delete fTrackCuts;
138   if (fMCCuts)              delete fMCCuts;
139 }
140
141 //_________________________________________________
142 void AliCFSingleTrackEfficiencyTask::Init() 
143 {
144   //
145   // Initialization, checks + copy cuts
146   //
147   if(!fMCCuts) {
148     AliFatal(" MC Cuts not defined");
149     return;
150   }
151   if(!fTrackCuts) {
152     AliFatal(" Track Cuts not defined");
153     return;
154   }
155
156   AliESDtrackCuts* copyfTrackCuts = new AliESDtrackCuts(*fTrackCuts);
157   const char* nameoutput = GetOutputSlot(4)->GetContainer()->GetName();
158   copyfTrackCuts->SetName(nameoutput);
159   // Post the data
160   PostData(4,copyfTrackCuts);
161
162   AliSingleTrackEffCuts* copyfMCCuts = new AliSingleTrackEffCuts(*fMCCuts);
163   nameoutput = GetOutputSlot(5)->GetContainer()->GetName();
164   copyfMCCuts->SetName(nameoutput);
165   // Post the data
166   PostData(5,copyfMCCuts);
167
168   return;
169 }
170
171 //_________________________________________________________________
172 void AliCFSingleTrackEfficiencyTask::UserExec(Option_t *)
173 {
174   //
175   // User Exec
176   //
177
178   Info("UserExec","") ;
179
180   AliVEvent* event = fInputEvent;
181
182   if(!fInputEvent) {
183     AliFatal("NO EVENT FOUND!");
184     return;
185   }
186   if(!fMCEvent) {
187     AliFatal("NO MC INFO FOUND");
188     return;
189   }
190
191   fHistEventsProcessed->Fill(0.5); // # of Event proceed        
192   Bool_t IsEventMCSelected = kFALSE;
193   Bool_t isAOD = fInputEvent->IsA()->InheritsFrom("AliAODEvent");
194
195   //
196   // Step 0: MC Gen Event Selection
197   //
198   if(isAOD) {
199     if(!event && AODEvent() && IsStandardAOD()) {
200       // In case there is an AOD handler writing a standard AOD, use the AOD 
201       // event in memory rather than the input (ESD) event.    
202       event = dynamic_cast<AliAODEvent*> (AODEvent());
203     }
204     IsEventMCSelected = fMCCuts->IsMCEventSelected(event);//AODs
205   } else {
206     IsEventMCSelected = fMCCuts->IsMCEventSelected(fMCEvent);//ESDs
207   }
208
209   // pass to the manager the event info to the cuts that need it 
210   if(!isAOD) fCFManager->SetMCEventInfo(fMCEvent);
211   else fCFManager->SetMCEventInfo(event);
212   fCFManager->SetRecEventInfo(event);
213
214   if(!IsEventMCSelected) {
215     AliDebug(3,"MC event not passing the quality criteria \n");
216     PostData(1,fHistEventsProcessed);
217     PostData(2,fCFManager->GetParticleContainer());
218     PostData(3,fQAHistList);
219     return;
220   }
221   fHistEventsProcessed->Fill(1.5); // # of Event after passing MC cuts
222        
223
224   //
225   // Step 1-3: Check the MC generated particles
226   // 
227   if(!isAOD) CheckESDMCParticles();
228   else CheckAODMCParticles();
229
230   //
231   // Step 4-7: Reconstructed event and track selection
232   //
233   Bool_t isRecoEventOk = fMCCuts->IsRecoEventSelected(event);
234
235   if(isRecoEventOk) {
236     fHistEventsProcessed->Fill(2.5); // # of Event after passing all cuts
237     CheckReconstructedParticles();
238   }
239   else AliDebug(3,"Event not passing quality criteria\n");
240
241
242   PostData(1,fHistEventsProcessed);
243   PostData(2,fCFManager->GetParticleContainer());
244   PostData(3,fQAHistList);
245
246   return;
247 }
248
249
250 //___________________________________________________________________________
251 void AliCFSingleTrackEfficiencyTask::Terminate(Option_t*)
252 {
253   //
254   // Terminate
255   //
256
257   Info("Terminate","");
258   AliAnalysisTaskSE::Terminate();
259
260   /*
261   //draw some example plots....
262   AliCFContainer *cont= dynamic_cast<AliCFContainer*> (GetOutputData(2));
263        
264   TH1D* h00 =   cont->ShowProjection(5,0) ;
265   TH1D* h01 =   cont->ShowProjection(5,1) ;
266   TH1D* h02 =   cont->ShowProjection(5,2) ;
267   TH1D* h03 =   cont->ShowProjection(5,3) ;
268   TH1D* h04 =   cont->ShowProjection(5,4) ;
269   TH1D* h05 =   cont->ShowProjection(5,5) ;
270   TH1D* h06 =   cont->ShowProjection(5,6) ;
271        
272   h00->SetMarkerStyle(23) ;
273   h01->SetMarkerStyle(24) ;
274   h02->SetMarkerStyle(25) ;
275   h03->SetMarkerStyle(26) ;
276   h04->SetMarkerStyle(27) ;
277   h05->SetMarkerStyle(28) ;
278   h06->SetMarkerStyle(28) ;
279        
280        
281        
282   TCanvas * c =new TCanvas("c","",1400,800);
283   c->Divide(4,2);
284        
285   c->cd(1);
286   h00->Draw("p");
287   c->cd(2);
288   h01->Draw("p");
289   c->cd(3);
290   h02->Draw("p");
291   c->cd(4);
292   h03->Draw("p");
293   c->cd(5);
294   h04->Draw("p");
295   c->cd(6);
296   h05->Draw("p");
297   c->cd(7);
298   h06->Draw("p");
299        
300   c->SaveAs("plots.eps");
301   */
302        
303 }
304
305
306 //___________________________________________________________________________
307 void AliCFSingleTrackEfficiencyTask::UserCreateOutputObjects() 
308 {
309   //
310   // UserCreateOutputObjects
311   //
312   
313   Info("CreateOutputObjects","CreateOutputObjects of task %s", GetName());
314
315   //slot #1
316   OpenFile(1);
317        
318   const char* nameoutput=GetOutputSlot(1)->GetContainer()->GetName();
319   //       fHistEventsProcessed = new TH1I(nameoutput"fHistEventsProcessed","fHistEventsProcessed",3,0,3) ;
320   fHistEventsProcessed = new TH1I(nameoutput,"fHistEventsProcessed",3,0,3) ;
321   fHistEventsProcessed->GetXaxis()->SetBinLabel(1,"All events");
322   fHistEventsProcessed->GetXaxis()->SetBinLabel(2,"Good MC events");
323   fHistEventsProcessed->GetXaxis()->SetBinLabel(3,"Good Reconstructed events");
324        
325   PostData(1,fHistEventsProcessed) ;
326   PostData(2,fCFManager->GetParticleContainer()) ;
327   PostData(3,fQAHistList) ;
328        
329   return;
330 }
331
332
333 //_________________________________________________________________________
334 void AliCFSingleTrackEfficiencyTask::CheckESDMCParticles()
335 {
336   //
337   // Check ESD generated particles
338   //
339   if (!fMCEvent) {
340     AliFatal("NO MC INFO FOUND");
341     return;
342   }
343
344   Double_t containerInput[6] = { 0., 0., 0., 0., 0., 0. };
345
346   TArrayF vtxPos(3);
347   AliGenEventHeader *genHeader = NULL;
348   genHeader = fMCEvent->GenEventHeader();
349   genHeader->PrimaryVertex(vtxPos);
350   containerInput[4] = vtxPos[2]; // z-vtx position
351
352   Double_t multiplicity = (Double_t)GetNumberOfTrackletsInEtaRange(-1.0,1.0);
353   containerInput[5] = multiplicity; //reconstructed number of tracklets
354
355   // loop on the MC Generated Particles
356   for (Int_t ipart=0; ipart<fMCEvent->GetNumberOfTracks(); ipart++) {
357
358     AliMCParticle *mcPart = (AliMCParticle*)fMCEvent->GetTrack(ipart);
359     containerInput[0] = (Float_t)mcPart->Pt();
360     containerInput[1] = mcPart->Eta();
361     containerInput[2] = mcPart->Phi();
362     containerInput[3] = mcPart->Theta();
363
364     // Step 1. Particle passing through Generation criteria and filling
365     if( !fMCCuts->IsMCParticleGenerated(mcPart) ) {
366       AliDebug(3,"MC Particle not passing through genetations criteria\n");
367       continue;
368     }
369     fCFManager->GetParticleContainer()->Fill(containerInput,kStepMCGenCut);
370
371     // Step 2. Particle passing through Kinematic criteria and filling
372     if( !fMCCuts->IsMCParticleInKineAcceptance(mcPart) ) {
373       AliDebug(3,"MC Particle not in the kine acceptance\n");
374       continue;
375     }
376     fCFManager->GetParticleContainer()->Fill(containerInput,kStepMCKineCut);
377
378     // Step 3. Particle passing through Track ref criteria and filling
379     // did leave signal (enough clusters) on the detector
380     if( !fMCCuts->IsMCParticleInReconstructable(mcPart) ) {
381       AliDebug(3,"MC Particle not in the reconstructible\n");
382       continue;
383     }
384     fCFManager->GetParticleContainer()->Fill(containerInput,kStepMCAccpCut);
385
386   }// end of particle loop
387
388   return;
389 }
390
391
392 //________________________________________________________________________
393 void AliCFSingleTrackEfficiencyTask::CheckAODMCParticles()
394 {
395   //
396   // Check AOD generated particles
397   //
398   if (!fInputEvent) {
399     AliFatal("NO EVENT FOUND!");
400     return;
401   }
402
403   AliAODEvent* event = dynamic_cast<AliAODEvent*>(fInputEvent);
404   if(!event && AODEvent() && IsStandardAOD()) {
405     // In case there is an AOD handler writing a standard AOD, use the AOD 
406     // event in memory rather than the input (ESD) event.    
407     event = dynamic_cast<AliAODEvent*> (AODEvent());
408   }
409   TClonesArray* mcArray = dynamic_cast<TClonesArray*>(event->FindListObject(AliAODMCParticle::StdBranchName()));
410   if (!mcArray) {
411     AliError("Could not find Monte-Carlo in AOD");
412     return;
413   }
414   AliAODMCHeader *mcHeader=NULL;
415   mcHeader = dynamic_cast<AliAODMCHeader*>(event->GetList()->FindObject(AliAODMCHeader::StdBranchName()));
416   if (!mcHeader) {
417     AliError("Could not find MC Header in AOD");
418     return;
419   }
420
421   Double_t containerInput[6] = { 0., 0., 0., 0., 0., 0. };
422   // Set the z-vertex position
423   containerInput[4]  = mcHeader->GetVtxZ();
424   // Multiplicity of the event defined as Ntracklets |eta|<1.0
425   Double_t multiplicity = (Double_t)GetNumberOfTrackletsInEtaRange(-1.0,1.0);
426   containerInput[5] = multiplicity;
427
428   for (Int_t ipart=0; ipart<mcArray->GetEntriesFast(); ipart++) {
429
430     AliAODMCParticle* mcPart = dynamic_cast<AliAODMCParticle*>(mcArray->At(ipart));
431     if (!mcPart){
432       AliError("Failed casting particle from MC array!, Skipping particle");
433       continue;
434     }
435     containerInput[0] = (Float_t)mcPart->Pt();
436     containerInput[1] = mcPart->Eta();
437     containerInput[2] = mcPart->Phi();
438     containerInput[3] = mcPart->Theta();
439
440     // Step 1. Particle passing through Generation criteria and filling
441     if( !fMCCuts->IsMCParticleGenerated(mcPart) ) {
442       AliDebug(3,"MC Particle not passing quality criteria\n");
443       continue;
444     }
445     fCFManager->GetParticleContainer()->Fill(containerInput,kStepMCGenCut);
446
447     // Step 2. Particle passing through Kinematic criteria and filling
448     if( !fMCCuts->IsMCParticleInKineAcceptance(mcPart) ) {
449       AliDebug(3,"MC Particle not in the acceptance\n");
450       continue;
451     }
452     fCFManager->GetParticleContainer()->Fill(containerInput,kStepMCKineCut);
453
454     // Step 3. Particle passing through Track Ref criteria and filling
455     // but no info available for Track ref in AOD fillng same as above
456     fCFManager->GetParticleContainer()->Fill(containerInput,kStepMCAccpCut);
457   }
458
459   return;
460 }
461
462
463 //_______________________________________________________________________________
464 AliESDtrack * AliCFSingleTrackEfficiencyTask::ConvertTrack(AliAODTrack *track)
465 {
466   //
467   // Convert an AOD track to an  ESD track to apply ESDtrackCuts
468   //
469
470   AliAODEvent *aodEvent = (AliAODEvent*)fInputEvent;
471   const AliAODVertex *primary = aodEvent->GetPrimaryVertex();
472   Double_t pos[3],cov[6];
473   primary->GetXYZ(pos);
474   primary->GetCovarianceMatrix(cov);
475   const AliESDVertex vESD(pos,cov,100.,100);
476
477   AliESDtrack *esdTrack =  new AliESDtrack(track);
478   // set the TPC cluster info
479   esdTrack->SetTPCClusterMap(track->GetTPCClusterMap());
480   esdTrack->SetTPCSharedMap(track->GetTPCSharedMap());
481   esdTrack->SetTPCPointsF(track->GetTPCNclsF());
482   // needed to calculate the impact parameters
483   esdTrack->RelateToVertex(&vESD,0.,3.);
484   //  std::cout << " primary vtx "<< primary << std::endl;
485   //  std::cout << " esdvtx "<< vESD.GetName() << std::endl;
486   //  std::cout<< " esdtrack pt "<< esdTrack.Pt() << " and status " << esdTrack.GetStatus() <<endl;
487   //  std::cout << " aod track "<< track<< " and status " << track->GetStatus() << std::endl;
488   //  std::cout << " esd track "<< esdTrack<< " and status " << esdTrack->GetStatus() << std::endl;
489   return esdTrack;
490 }
491
492
493 //___________________________________________________________________________
494 void AliCFSingleTrackEfficiencyTask::CheckReconstructedParticles()
495 {
496   //
497   // Check reconstructed particles
498   //
499
500   AliVEvent* event = fInputEvent;
501   Bool_t isAOD = fInputEvent->IsA()->InheritsFrom("AliAODEvent");
502   if(!event && AODEvent() && IsStandardAOD()) {
503    event  = dynamic_cast<AliAODEvent*> (AODEvent());
504   }
505   if (!event) return;
506
507   Double_t containerInput[6] = { 0, 0, 0, 0, 0, 0};   // contains reconstructed quantities
508   Double_t containerInputMC[6] = { 0, 0, 0, 0, 0, 0}; // contains generated quantities
509   
510   const AliVVertex *vertex = event->GetPrimaryVertex();
511   // set the z-vertex position
512   containerInput[4] = vertex->GetZ();
513   containerInputMC[4] = containerInput[4];
514   // set the event multiplicity as Ntracklets in |eta|<1.0
515   Double_t multiplicity = (Double_t)GetNumberOfTrackletsInEtaRange(-1.0,1.0);
516   containerInput[5] = multiplicity;
517   containerInputMC[5] = multiplicity;
518
519   // Reco tracks track loop
520   AliVParticle* track = NULL;
521   for (Int_t iTrack = 0; iTrack<event->GetNumberOfTracks(); iTrack++) {
522
523     track = event->GetTrack(iTrack);
524     if(!track) {
525       AliDebug(3,Form("Track %d not found",iTrack));
526       continue;
527     }
528
529     Double_t mom[3];
530     track->PxPyPz(mom);
531     Double_t pt = TMath::Sqrt(mom[0]*mom[0]+mom[1]*mom[1]);
532     containerInput[0] = pt;
533     containerInput[1] = track->Eta();
534     containerInput[2] = track->Phi();
535     containerInput[3] = track->Theta();
536
537     //
538     // Step 4. Track that are recostructed, filling
539     //
540     fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructed);
541
542     //
543     // Step 5. Track that are recostructed and pass acceptance cuts, filling
544     //
545     if (!fMCCuts->IsRecoParticleKineAcceptance(track)) continue;
546     fCFManager->GetParticleContainer()->Fill(containerInput,kStepRecoKineCuts);
547
548     // is track associated to particle ? if yes + implimenting the physical primary..
549     Int_t label = TMath::Abs(track->GetLabel());
550     if (label<=0) {
551       AliDebug(3,"Particle not matching MC label \n");
552       continue;
553     }
554
555     //
556     // Step 6-7. Track that are recostructed + Kine + Quality criteria, filling
557     //
558
559     // check particle selections at MC level
560     AliVParticle *mcPart  = (AliVParticle*)fMCEvent->GetTrack(label);
561     if(!mcPart) continue;
562     containerInputMC[0] = (Float_t)mcPart->Pt();
563     containerInputMC[1] = mcPart->Eta();
564     containerInputMC[2] = mcPart->Phi();
565     containerInputMC[3] = mcPart->Theta();
566
567     if (!fMCCuts->IsMCParticleGenerated(mcPart)) continue;
568     //    cout<< "MC matching did work"<<endl;
569
570
571     // for filter bit selection
572     AliAODTrack *aodTrack = dynamic_cast<AliAODTrack*>(track);
573     if(isAOD && !aodTrack) continue;
574     if(isAOD && fSetFilterBit) if (!aodTrack->TestFilterMask(BIT(fbit))) continue;
575     //    cout<<" Filter bit check passed"<<endl;
576
577     Bool_t isESDtrack = track->IsA()->InheritsFrom("AliESDtrack");
578     AliESDtrack *tmptrack = NULL;
579     if(isESDtrack) {
580       tmptrack = dynamic_cast<AliESDtrack*>(track); // ESDs
581     } else {
582       tmptrack = ConvertTrack(aodTrack); // AODs
583     }
584     if (!tmptrack) continue;
585
586     // exclude global constrained and TPC only tracks (filter bits 128 and 512)
587     Int_t id = tmptrack->GetID();
588     if(isAOD && id<0) {
589       AliDebug(3,"Track removed bc corresponts to either filter bit 128 or 512 (TPC only tracks)\n");
590       continue;
591     }
592
593     // Apply ESD track cuts
594     if( !fTrackCuts->IsSelected(tmptrack) ){
595       AliDebug(3,"Reconstructed track not passing quality criteria\n");
596       if(isAOD) delete tmptrack;
597       continue;
598     }
599     fCFManager->GetParticleContainer()->Fill(containerInputMC,kStepReconstructedMC);
600     fCFManager->GetParticleContainer()->Fill(containerInput,kStepRecoQualityCuts);
601     if(isAOD) delete tmptrack;
602
603     //
604     // Step8, PID check
605     //
606     if( !fMCCuts->IsRecoParticlePID(track) ){
607       AliDebug(3,"Reconstructed track not passing PID criteria\n");
608       continue;
609     }
610     fCFManager->GetParticleContainer()->Fill(containerInput,kStepRecoPID);
611
612   }
613   return;
614 }
615
616 //______________________________________________________________________
617 Int_t AliCFSingleTrackEfficiencyTask::GetNumberOfTrackletsInEtaRange(Double_t mineta, Double_t maxeta)
618 {
619   //
620   // counts tracklets in given eta range
621   //
622
623   AliAODEvent* event = dynamic_cast<AliAODEvent*>(fInputEvent);
624   Bool_t isAOD = fInputEvent->IsA()->InheritsFrom("AliAODEvent");
625   if(!event && AODEvent() && IsStandardAOD()) {
626    event  = dynamic_cast<AliAODEvent*> (AODEvent());
627   }
628   if (!event) return -1;
629   Int_t count=0;
630
631   if(isAOD) {
632     AliAODTracklets* tracklets = event->GetTracklets();
633     Int_t nTr=tracklets->GetNumberOfTracklets();
634     for(Int_t iTr=0; iTr<nTr; iTr++){
635       Double_t theta=tracklets->GetTheta(iTr);
636       Double_t eta=-TMath::Log(TMath::Tan(theta/2.));
637       if(eta>mineta && eta<maxeta) count++;
638     }
639   } else {  
640     AliESDEvent *esdEvent = dynamic_cast<AliESDEvent*>(fInputEvent);
641     if (!esdEvent) return -1;
642     const AliMultiplicity *mult = esdEvent->GetMultiplicity();
643     if (mult) {
644       if (mult->GetNumberOfTracklets()>0) {     
645         for (Int_t n=0; n<mult->GetNumberOfTracklets(); n++) {
646           Double_t eta = -TMath::Log( TMath::Tan( mult->GetTheta(n) / 2.0 ) );
647           if(eta>mineta && eta<maxeta) count++;
648         }
649       }
650     }
651   }
652
653   return count;
654 }