]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/hfe/AliAnalysisTaskDisplacedElectrons.cxx
Update of the HFE package
[u/mrichter/AliRoot.git] / PWG3 / hfe / AliAnalysisTaskDisplacedElectrons.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15 //
16 // The analysis task:
17 // study displaced electrons from beauty and charm 
18 // with cut on impact parameters in various pT bins
19 // 
20 // 
21 // Authors:
22 //  Hongyan Yang <hongyan@physi.uni-heidelberg.de>
23 //
24
25 #include <TChain.h>
26 #include <TFile.h>
27 #include <TH1F.h>
28 #include <TH1I.h>
29 #include "AliLog.h"
30
31 #include <TList.h>
32 #include <TMath.h>
33 #include <TObjArray.h>
34 #include <TParticle.h>
35 #include <TString.h>
36
37 #include <TCanvas.h>
38
39 #include "AliCFManager.h"
40 #include "AliMCEvent.h"
41 #include "AliMCEventHandler.h"
42 #include "AliESDEvent.h"
43 #include "AliESDInputHandler.h"
44 #include "AliESDtrack.h"
45 #include "AliAnalysisManager.h"
46 #include "AliMCParticle.h"
47
48
49 #include "AliHFEpid.h"
50 #include "AliHFEcuts.h"
51 #include "AliHFEdisplacedElectrons.h"
52
53 #include "AliAnalysisTaskDisplacedElectrons.h"
54
55
56 //____________________________________________________________
57 AliAnalysisTaskDisplacedElectrons::AliAnalysisTaskDisplacedElectrons():
58   AliAnalysisTaskSE("Task for displaced electron study")
59   , fDebugLevel(0)
60   , fPIDdetectors("")
61   , fPIDstrategy(0)
62   , fPlugins(0)
63   , fESD(0x0)
64   , fMC(0x0)
65   , fCuts(0x0)
66   , fPID(0x0)
67   , fCFM(0x0)
68   , fNEvents(0x0)
69   , fElectronsPt(0x0)
70   , fOutput(0x0)
71   , fCorrection(0x0)
72   , fDisplacedElectrons(0x0)
73   , fHistDisplacedElectrons(0x0)
74 {
75   //
76   // Dummy constructor
77   //
78   DefineInput(0, TChain::Class());
79   DefineOutput(2, TList::Class());  // output
80   DefineOutput(1, TList::Class());  // output
81   
82   SetHasMCData();
83
84   // Initialize pid
85   fPID = new AliHFEpid;
86   
87
88 }
89
90 //____________________________________________________________
91 AliAnalysisTaskDisplacedElectrons::AliAnalysisTaskDisplacedElectrons(const char * name):
92   AliAnalysisTaskSE(name)
93   , fDebugLevel(0)
94   , fPIDdetectors("")
95   , fPIDstrategy(0)
96   , fPlugins(0)
97   , fESD(0x0)
98   , fMC(0x0)
99   , fCuts(0x0)
100   , fPID(0x0)
101   , fCFM(0x0)
102   , fNEvents(0x0)
103   , fElectronsPt(0x0)
104   , fOutput(0x0)
105   , fCorrection(0x0)
106   , fDisplacedElectrons(0x0)
107   , fHistDisplacedElectrons(0x0)
108 {
109   //
110   // Default constructor
111   //
112   DefineInput(0, TChain::Class());
113   DefineOutput(2, TList::Class());
114   DefineOutput(1, TList::Class());
115
116   SetHasMCData();
117 }
118
119 //____________________________________________________________
120 AliAnalysisTaskDisplacedElectrons::AliAnalysisTaskDisplacedElectrons(const AliAnalysisTaskDisplacedElectrons &ref):
121   AliAnalysisTaskSE(ref)
122   , fDebugLevel(ref.fDebugLevel)
123   , fPIDdetectors(ref.fPIDdetectors)
124   , fPIDstrategy(ref.fPIDstrategy)
125   , fPlugins(ref.fPlugins)
126   , fESD(ref.fESD)
127   , fMC(ref.fMC)
128   , fCuts(ref.fCuts)
129   , fPID(ref.fPID)
130   , fCFM(ref.fCFM)
131   , fNEvents(ref.fNEvents)
132   , fElectronsPt(ref.fElectronsPt)
133   , fOutput(ref.fOutput)
134   , fCorrection(ref.fCorrection)
135   , fDisplacedElectrons(ref.fDisplacedElectrons)
136   , fHistDisplacedElectrons(ref.fHistDisplacedElectrons)
137 {
138   //
139   // Copy Constructor
140   //
141 }
142
143 //____________________________________________________________
144 AliAnalysisTaskDisplacedElectrons &AliAnalysisTaskDisplacedElectrons::operator=(const AliAnalysisTaskDisplacedElectrons &ref){
145   //
146   // Assignment operator
147   //
148   if(this == &ref) return *this;
149   AliAnalysisTask::operator=(ref);
150   fDebugLevel = ref.fDebugLevel;
151   fPIDdetectors = ref.fPIDdetectors;
152   fPIDstrategy = ref.fPIDstrategy;
153   fPlugins = ref.fPlugins;
154   fESD = ref.fESD;
155   fMC = ref.fMC;
156   fPID = ref.fPID;
157   fCuts = ref.fCuts;
158   fCFM = ref.fCFM;
159   fNEvents = ref.fNEvents;
160   fOutput = ref.fOutput;
161   fCorrection = ref.fCorrection;
162   fDisplacedElectrons = ref.fDisplacedElectrons;
163   fHistDisplacedElectrons = ref.fHistDisplacedElectrons;
164
165   return *this;
166 }
167
168 //____________________________________________________________
169 AliAnalysisTaskDisplacedElectrons::~AliAnalysisTaskDisplacedElectrons(){
170   //
171   // Destructor
172   //
173
174   if(fPID) delete fPID;
175   
176   if(fESD) delete fESD;
177   if(fMC) delete fMC;
178
179   if(fCFM) delete fCFM;
180
181   if(fNEvents) delete fNEvents;
182   
183   if(fElectronsPt) delete fElectronsPt;
184
185   if(fOutput){ 
186     fOutput->Clear();
187     delete fOutput;
188   }
189   
190
191   if(fCorrection){
192     fCorrection->Clear();
193     delete fCorrection;
194   }
195  
196   if(fDisplacedElectrons) delete fDisplacedElectrons;  
197
198   if(fHistDisplacedElectrons){ 
199     fHistDisplacedElectrons->Clear();  
200     delete fHistDisplacedElectrons;  
201   }
202    
203   
204 }
205
206 //____________________________________________________________
207 void AliAnalysisTaskDisplacedElectrons::UserCreateOutputObjects(){
208   // create output objects
209   // fNEvents
210   // MC and Data containers
211
212   fNEvents = new TH1I("nEvents", "Number of Events in the Analysis", 2, 0, 2); 
213   const Int_t nBins = 14;
214   const Float_t ptBins[nBins] = {0.0,0.5,1.0,1.5,2.0,2.5,3.0,4.0,5.0,7.0,9.0,12.0,16.0,20.0};
215   fElectronsPt = new TH1F("esdPt", "p_{T} distribution of identified electrons (HFEpid);p_{T} (GeV/c);dN/dp_{T};", nBins-1, ptBins); 
216   if(!fOutput) fOutput = new TList;
217   fOutput->SetName("results");
218   fOutput->AddAt(fElectronsPt,0);
219   fOutput->AddAt(fNEvents,1);
220   // Initialize correction Framework and Cuts
221   fCFM = new AliCFManager;
222   MakeEventContainer();
223   MakeParticleContainer();
224  
225   if(!fCorrection) fCorrection = new TList();
226   fCorrection->SetName("corrections");
227   fCorrection->AddAt(fCFM->GetEventContainer(), 0);
228   fCorrection->AddAt(fCFM->GetParticleContainer(), 1);
229   fCorrection->Print();
230
231   // Temporary fix: Initialize particle cuts with 0x0
232   for(Int_t istep = 0; istep < fCFM->GetEventContainer()->GetNStep(); istep++)
233     fCFM->SetEventCutsList(istep, 0x0);
234   for(Int_t istep = 0; istep < fCFM->GetParticleContainer()->GetNStep(); istep++)
235     fCFM->SetParticleCutsList(istep, 0x0);
236   if(!fCuts){
237     AliWarning("Cuts not available. Default cuts will be used");
238     fCuts = new AliHFEcuts;
239     fCuts->CreateStandardCuts();
240
241   }
242   
243   fCuts->SetCutITSpixel(AliHFEextraCuts::kBoth);
244   
245   fCuts->Initialize(fCFM);
246     
247   fPID->SetHasMCData(HasMCData());
248   if(!fPIDdetectors.Length() && ! fPIDstrategy) AddPIDdetector("TPC");
249   if(fPIDstrategy)
250     fPID->InitializePID(Form("Strategy%d", fPIDstrategy));
251   else
252     fPID->InitializePID(fPIDdetectors.Data());     // Only restrictions to TPC allowed 
253
254   // displaced electron study----------------------------------
255   if(GetPlugin(kDisplacedElectrons)){
256     
257     fDisplacedElectrons = new AliHFEdisplacedElectrons;
258     fDisplacedElectrons->SetDebugLevel(fDebugLevel);
259     fDisplacedElectrons->SetHasMCData(HasMCData());
260     
261     if(!fHistDisplacedElectrons) fHistDisplacedElectrons = new TList();
262     
263     fDisplacedElectrons->CreateOutputs(fHistDisplacedElectrons);
264     
265     fOutput->AddAt(fHistDisplacedElectrons, 2);
266   }
267
268
269 }
270
271
272
273 //____________________________________________________________
274 void AliAnalysisTaskDisplacedElectrons::UserExec(Option_t *){
275   //
276   // Run the analysis
277   // 
278
279   AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler *>(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
280   if(!esdH){      
281     AliError("No ESD input handler");
282     return;
283   } else {
284     fESD = esdH->GetEvent();
285   }
286  
287   AliMCEventHandler *mcH = dynamic_cast<AliMCEventHandler *>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
288   if(!mcH){       
289     AliError("No MC truth handler");
290     return;
291   } else {
292     fMC = mcH->MCEvent();
293   }
294   
295   if(fDebugLevel>=5)AliInfo("Processing ESD events");
296   if(!fESD){
297     AliError("No ESD Event");
298     return;
299   }
300   
301   if(HasMCData()){
302     if(fDebugLevel>=5)AliInfo(Form("MC Event: %p", fMC));
303     if(!fMC){
304       AliError("No MC Event, but MC Data required");
305       return;
306     }
307   }
308   if(!fCuts){
309     AliError("HFE cuts not available");
310     return;
311   }
312   
313   // process data: ESD tracks with MC information
314   Double_t container[8];   // container for the output in THnSparse
315   memset(container, 0, sizeof(Double_t) * 8);
316   
317   Bool_t signal = kTRUE;
318   Bool_t alreadyseen = kFALSE;
319   LabelContainer cont(fESD->GetNumberOfTracks());
320
321   Int_t nElectrons=0;
322   
323   AliESDtrack *track = 0x0;    
324   
325   Double_t nContrib = fESD->GetPrimaryVertex()->GetNContributors();
326   if(!fCFM->CheckEventCuts(AliHFEcuts::kEventStepReconstructed, fESD)) return;
327   if(GetPlugin(kCorrection)){
328     fCFM->GetEventContainer()->Fill(&nContrib, AliHFEcuts::kEventStepReconstructed);
329   }
330   fCFM->SetRecEventInfo(fESD);
331
332   for(Int_t itrack = 0; itrack < fESD->GetNumberOfTracks(); itrack++){
333     track = fESD->GetTrack(itrack);
334     
335     //
336     //  track quality cut: 1st step: ITS and TPC cut; 2nd step: Rec prim; 3rd step: hfe cut on ITS pixel layer
337     //
338     // require within rapidity range +/-0.9
339
340     //   if((TMath::Abs(track->Eta()))>0.9) continue;
341
342         
343     if(GetPlugin(kDisplacedElectrons)) {
344
345       // 1st cut
346       // RecKine: ITSTPC cuts : ITS & TPC refit, covmatrix: (2, 2, 0.5, 0.5, 2); min_tpccls: 50, chi2_tpccls: 3.5
347       if(!fCFM->CheckParticleCuts(AliHFEcuts::kStepRecKineITSTPC, track)){
348         signal = kFALSE;
349         continue;
350       }
351       
352       container[0] = track->Pt();
353       container[1] = track->Eta();
354       container[2] = track->Phi();
355       container[3] = track->Charge();
356       
357       AliStack *stack = fMC->Stack();
358       if(!stack) continue;      
359       
360       AliMCParticle *mctrack = NULL;
361       if(HasMCData() ){ 
362         mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(track->GetLabel())));
363         if(!mctrack)  continue;                 
364         
365         container[4] = mctrack->Pt();
366         container[5] = mctrack->Eta();
367         container[6] = mctrack->Phi();
368         container[7] = mctrack->Charge();
369
370         //      if(!fCFM->CheckParticleCuts(AliHFEcuts::kStepMCGenerated, mctrack)) signal = kFALSE;
371         if(TMath::Abs(mctrack->Eta())>0.9) {
372           signal = kFALSE;
373           continue;
374         }  // cut on kinematics
375
376       }  // has MC data
377
378       
379       if(signal && GetPlugin(kCorrection)){
380         alreadyseen = cont.Find(TMath::Abs(track->GetLabel()));  // double counted track
381         cont.Append(TMath::Abs(track->GetLabel()));
382         fCFM->GetParticleContainer()->Fill(&container[4], AliHFEcuts::kStepRecKineITSTPC);
383         fCFM->GetParticleContainer()->Fill(&container[0], AliHFEcuts::kStepRecKineITSTPC + 2*AliHFEcuts::kNcutStepsESDtrack);
384         if(alreadyseen) 
385           fCFM->GetParticleContainer()->Fill(&container[4], AliHFEcuts::kStepRecKineITSTPC + AliHFEcuts::kNcutStepsESDtrack);
386         
387       }  // fill correction --- 1st
388     
389       // second cut
390       // RecPrim: cut on track quality : DCA to vertex max: 3cm and 10cm; reject kink daughters
391       if(!fCFM->CheckParticleCuts(AliHFEcuts::kStepRecPrim, track)) { 
392         signal = kFALSE;
393         continue;
394       }
395
396       if(signal) {
397         alreadyseen = cont.Find(TMath::Abs(track->GetLabel()));
398         cont.Append(TMath::Abs(track->GetLabel()));
399         
400         fCFM->GetParticleContainer()->Fill(&container[4], AliHFEcuts::kStepRecPrim);
401         fCFM->GetParticleContainer()->Fill(&container[0], AliHFEcuts::kStepRecPrim + 2*AliHFEcuts::kNcutStepsESDtrack);
402         if(alreadyseen) {
403           fCFM->GetParticleContainer()->Fill(&container[4], AliHFEcuts::kStepRecPrim + AliHFEcuts::kNcutStepsESDtrack);
404         }
405       }
406       
407       //  third cut
408       // HFEcuts: ITS layers cuts: ITS pixel layer: kFirst, kSecond, kBoth, kNone or kAny
409       if(!fCFM->CheckParticleCuts(AliHFEcuts::kStepHFEcutsITS, track)){
410         signal = kFALSE;
411         continue;
412       }
413
414       if(signal) {
415         alreadyseen = cont.Find(TMath::Abs(track->GetLabel()));
416         cont.Append(TMath::Abs(track->GetLabel()));
417         
418         fCFM->GetParticleContainer()->Fill(&container[4], AliHFEcuts::kStepHFEcutsITS);
419         fCFM->GetParticleContainer()->Fill(&container[0], AliHFEcuts::kStepHFEcutsITS + 2*AliHFEcuts::kNcutStepsESDtrack);
420         if(alreadyseen) {
421           fCFM->GetParticleContainer()->Fill(&container[4], AliHFEcuts::kStepHFEcutsITS + AliHFEcuts::kNcutStepsESDtrack);
422         }
423       }
424       
425       fDisplacedElectrons->FillMCOutput(fESD, track, stack);
426       
427       // track accepted, do PID --> only electron candidate will be processed
428       AliHFEpidObject hfetrack;
429       hfetrack.fAnalysisType = AliHFEpidObject::kESDanalysis;
430       hfetrack.fRecTrack = track;
431       if((fPID->IsSelected(&hfetrack)==kTRUE))
432         if(fDebugLevel>=10)
433           AliInfo(Form("ESD info: this particle is %s identified as electron by HFEpid method \n", (fPID->IsSelected(&hfetrack)==kTRUE)?" ":" NOT "));
434       if(!fPID->IsSelected(&hfetrack))   continue;
435       
436          // Fill Containers
437
438       nElectrons++;
439
440       fElectronsPt->Fill(track->Pt());
441
442       if(signal) {
443         fCFM->GetParticleContainer()->Fill(container, AliHFEcuts::kStepPID + 2*AliHFEcuts::kNcutStepsESDtrack);
444         fCFM->GetParticleContainer()->Fill(&container[4], AliHFEcuts::kStepPID);
445         if(alreadyseen) {
446           fCFM->GetParticleContainer()->Fill(&container[4], (AliHFEcuts::kStepPID + (AliHFEcuts::kNcutStepsESDtrack)));
447         }
448       }
449       
450       fDisplacedElectrons->FillESDOutput(fESD, track);
451       
452     } // analyze displaced electrons plugin switched on
453   } // track loop  
454
455   if(fDebugLevel>=5)
456     AliInfo(Form("ESD info: number of electrons found in this event: %d\n", nElectrons));
457
458
459   fNEvents->Fill(1);
460   
461   PostData(1, fOutput);
462   PostData(2, fCorrection);
463
464 }
465
466 //____________________________________________________________
467 void AliAnalysisTaskDisplacedElectrons::Terminate(Option_t *){
468   //
469   // Terminate not implemented at the moment
470   //  
471   
472   if(GetPlugin(kDisplacedElectrons))
473     {
474       
475       fOutput = dynamic_cast<TList *>(GetOutputData(1));
476       fCorrection= dynamic_cast<TList *>(GetOutputData(2));
477
478       if(!fOutput || !fCorrection){
479         if(!fCorrection) AliError("correction list not available\n");
480         if(!fOutput) AliError("output list not available\n");
481         
482         return;
483       }
484
485       fOutput->Print();
486       fCorrection->Print();
487       
488       AliInfo("analysis done!\n");
489       
490     }
491 }
492
493
494
495
496
497
498 //____________________________________________________________
499 void AliAnalysisTaskDisplacedElectrons::PrintStatus() const {
500   
501   //
502   // Print Analysis status
503   //
504   printf("\n\tAnalysis Settings\n\t========================================\n");
505   printf("\tdisplaced electrons' analysis is %s\n", GetPlugin(kDisplacedElectrons)?"ON":"OFF");
506   printf("\tcorrection container  is %s\n", GetPlugin(kCorrection)?"ON":"OFF");
507   printf("\tpost processing  is %s\n", GetPlugin(kPostProcess)?"ON":"OFF");
508   printf("\tcuts: %s\n", (fCuts != NULL) ? "YES" : "NO");
509   printf("\t ");
510   printf("\n");
511 }
512
513 //__________________________________________                                                  
514 void AliAnalysisTaskDisplacedElectrons::SwitchOnPlugin(Int_t plug){
515   //                                            
516   // Switch on Plugin          
517   // Available:                                  
518   //  - analyze impact parameter
519   //  - Post Processing                                                                      
520   
521   switch(plug)
522     {
523     case kDisplacedElectrons: 
524       SETBIT(fPlugins, plug); 
525       break;
526     case kPostProcess: 
527       SETBIT(fPlugins, plug); 
528       break;
529     case kCorrection:
530       SETBIT(fPlugins, plug); 
531       break;
532     default: 
533       AliError("Unknown Plugin");
534     };
535 }
536
537
538 //____________________________________________________________
539 void AliAnalysisTaskDisplacedElectrons::MakeParticleContainer(){
540   //
541   // Create the particle container (borrowed from AliAnalysisTaskHFE)
542   //
543   const Int_t kNvar   = 4; //number of variables on the grid:pt,eta, phi
544   const Double_t kPtmin = 0.1, kPtmax = 20.;   // used for fixed pt bins
545   //  const Float_t ptBins[14] = {0.0,0.5,1.0,1.5,2.0,2.5,3.0,4.0,5.0,7.0,9.0,12.0,16.0,20.0};
546   const Double_t kEtamin = -0.9, kEtamax = 0.9;
547   const Double_t kPhimin = 0., kPhimax = 2. * TMath::Pi();
548
549   //arrays for the number of bins in each dimension
550   Int_t iBin[kNvar];
551   iBin[0] = 40; //bins in pt   // used for fixed pt bins
552   //  iBin[0] = 13; //bins in pt
553   iBin[1] =  8; //bins in eta 
554   iBin[2] = 18; // bins in phi
555   iBin[3] =  2; // bins in charge
556   
557   //arrays for lower bounds :
558   Double_t* binEdges[kNvar];
559   for(Int_t ivar = 0; ivar < kNvar; ivar++)
560     binEdges[ivar] = new Double_t[iBin[ivar] + 1];
561
562   //values for bin lower bounds
563   //  for(Int_t i=0; i<=iBin[0]; i++) binEdges[0][i]=ptBins[i];  // using variable bins
564   for(Int_t i=0; i<=iBin[0]; i++) 
565      binEdges[0][i]=(Double_t)TMath::Power(10,TMath::Log10(kPtmin) + (TMath::Log10(kPtmax)-TMath::Log10(kPtmin))/iBin[0]*(Double_t)i);   // fixed pt bin
566   for(Int_t i=0; i<=iBin[1]; i++) binEdges[1][i]=(Double_t)kEtamin  + (kEtamax-kEtamin)/iBin[1]*(Double_t)i;
567   for(Int_t i=0; i<=iBin[2]; i++) binEdges[2][i]=(Double_t)kPhimin  + (kPhimax-kPhimin)/iBin[2]*(Double_t)i;
568   for(Int_t i=0; i<=iBin[3]; i++) binEdges[3][i]=1.1*i-1.1; // Numeric precision
569   
570   //one "container" for MC
571   AliCFContainer* container = new AliCFContainer("container","container for tracks", 
572                                                  (AliHFEcuts::kNcutStepsTrack + 1 + 2*(AliHFEcuts::kNcutStepsESDtrack + 1)), 
573                                                  kNvar, iBin);
574
575   //setting the bin limits
576   for(Int_t ivar = 0; ivar < kNvar; ivar++)
577     container -> SetBinLimits(ivar, binEdges[ivar]);
578   fCFM->SetParticleContainer(container);
579   
580   //create correlation matrix for unfolding
581   Int_t thnDim[2*kNvar];
582   for (int k=0; k<kNvar; k++) {
583     //first half  : reconstructed 
584     //second half : MC
585     thnDim[k]      = iBin[k];
586     thnDim[k+kNvar] = iBin[k];
587   }
588
589   // add more containers
590
591
592 }
593
594 //____________________________________________________________
595 void AliAnalysisTaskDisplacedElectrons::MakeEventContainer(){
596   //
597   // Create the event container for the correction framework and link it
598   //
599   const Int_t kNvar = 1;  // number of variables on the grid: number of tracks per event
600   const Double_t kNTrackBound[2] = {-0.5, 200.5};
601   const Int_t kNBins = 201;
602
603   AliCFContainer *evCont = new AliCFContainer("eventContainer", "Container for events", AliHFEcuts::kNcutStepsEvent, kNvar, &kNBins);
604
605   Double_t trackBins[kNBins];
606   for(Int_t ibin = 0; ibin < kNBins; ibin++) trackBins[ibin] = kNTrackBound[0] + static_cast<Double_t>(ibin);
607   evCont->SetBinLimits(0,trackBins);
608
609   fCFM->SetEventContainer(evCont);
610
611 }
612
613
614
615 //__________________________________________________________
616 void AliAnalysisTaskDisplacedElectrons::AddPIDdetector(TString detector){
617   //
618   // Adding PID detector to the task
619   //
620   if(!fPIDdetectors.Length()) 
621     fPIDdetectors = detector;
622   else
623     fPIDdetectors += ":" + detector;
624 }
625
626
627
628 //____________________________________________________________
629 AliAnalysisTaskDisplacedElectrons::LabelContainer::LabelContainer(Int_t capacity):
630   fContainer(NULL),
631   fBegin(NULL),
632   fEnd(NULL),
633   fLast(NULL),
634   fCurrent(NULL)
635 {
636   //
637   // Default constructor
638   //
639   fContainer = new Int_t[capacity];
640   fBegin = &fContainer[0];
641   fEnd = &fContainer[capacity - 1];
642   fLast = fCurrent = fBegin;
643 }
644
645 //____________________________________________________________
646 Bool_t AliAnalysisTaskDisplacedElectrons::LabelContainer::Append(Int_t label){
647   //
648   // Add Label to the container
649   //
650   if(fLast > fEnd) return kFALSE;
651   *fLast++ = label;
652   return kTRUE;
653 }
654
655 //____________________________________________________________
656 Bool_t AliAnalysisTaskDisplacedElectrons::LabelContainer::Find(Int_t label) const {
657   //
658   // Find track in the list of labels
659   //
660   for(Int_t *entry = fBegin; entry <= fLast; entry++) 
661     if(*entry == label) return kTRUE;
662   return kFALSE;
663 }
664
665 //____________________________________________________________
666 Int_t AliAnalysisTaskDisplacedElectrons::LabelContainer::Next()  { 
667   //
668   // Mimic iterator
669   //
670   if(fCurrent > fLast) return -1; 
671   return *fCurrent++;
672 }