]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnAnalysisEffSE.cxx
Added cut for checking multiplicity and done some adaptments for AOD analysis with...
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnAnalysisEffSE.cxx
1 //
2 // Class AliRsnAnalysisEffSE
3 //
4 // TODO
5 // TODO
6 //
7 // authors: Alberto Pulvirenti (alberto.pulvirenti@ct.infn.it)
8 //          Martin Vala (martin.vala@cern.ch)
9 //
10 #include <Riostream.h>
11 #include "AliStack.h"
12 #include "AliESDEvent.h"
13 #include "AliMCEvent.h"
14
15 #include "AliCFContainer.h"
16 #include "AliRsnCutManager.h"
17 #include "AliRsnValue.h"
18 #include "AliRsnAnalysisEffSE.h"
19 #include "AliRsnPairDef.h"
20 #include "AliRsnCutSet.h"
21
22 ClassImp(AliRsnAnalysisEffSE)
23
24 //_____________________________________________________________________________
25 AliRsnAnalysisEffSE::AliRsnAnalysisEffSE(const char *name) :
26   AliRsnVAnalysisTaskSE(name),
27   fUseITSSA(kTRUE),
28   fUseGlobal(kTRUE),
29   fStepListMC(0),
30   fStepListESD(0),
31   fAxisList("AliRsnValue", 0),
32   fPairDefList(0),
33   fContainerList(0x0),
34   fOutList(0x0),
35   fVar(0),
36   fPair(),
37   fEventCuts("eventCuts", AliRsnCut::kEvent)
38 {
39 //
40 // Default constructor.
41 //
42
43   AliDebug(AliLog::kDebug+2,"<-");
44
45   DefineOutput(2, TList::Class());
46
47   AliDebug(AliLog::kDebug+2,"->");
48 }
49
50 //_____________________________________________________________________________
51 AliRsnAnalysisEffSE::AliRsnAnalysisEffSE(const AliRsnAnalysisEffSE& copy) :
52   AliRsnVAnalysisTaskSE(copy),
53   fUseITSSA(copy.fUseITSSA),
54   fUseGlobal(copy.fUseGlobal),
55   fStepListMC(copy.fStepListMC),
56   fStepListESD(copy.fStepListESD),
57   fAxisList(copy.fAxisList),
58   fPairDefList(copy.fPairDefList),
59   fContainerList(copy.fContainerList),
60   fOutList(0x0),
61   fVar(0),
62   fPair(),
63   fEventCuts(copy.fEventCuts)
64 {
65 //
66 // Copy constrtuctor
67 //
68 }
69
70 //_____________________________________________________________________________
71 void AliRsnAnalysisEffSE::AddStepMC(AliRsnCutManager *mgr)
72 {
73 //
74 // Add a step on montecarlo
75 //
76
77   fStepListMC.AddLast(mgr);
78 }
79
80 //_____________________________________________________________________________
81 void AliRsnAnalysisEffSE::AddStepESD(AliRsnCutManager *mgr) 
82 {
83 //
84 // Add a step on ESD
85 //
86
87   fStepListESD.AddLast(mgr);
88 }
89     
90 //_____________________________________________________________________________
91 void AliRsnAnalysisEffSE::AddAxis(AliRsnValue *axis) 
92 {
93 //
94 // Add a new axis
95 //
96
97   Int_t n = fAxisList.GetEntries();
98   new (fAxisList[n]) AliRsnValue(*axis);
99 }
100
101 //_____________________________________________________________________________
102 void AliRsnAnalysisEffSE::RsnUserCreateOutputObjects()
103 {
104 //
105 // Creation of output objects.
106 // These are created through the utility methods in the analysis manager,
107 // which produces a list of histograms for each specified set of pairs.
108 // Each of these lists is added to the main list of this task.
109 //
110
111   AliDebug(AliLog::kDebug+2,"<-");
112
113   // get number of steps and axes
114   Int_t iaxis  = 0;
115   Int_t nAxes  = fAxisList.GetEntries();
116   Int_t nSteps = (Int_t)fStepListMC.GetEntries() + (Int_t)fStepListESD.GetEntries();
117
118   if (!nSteps) {
119     AliError("No steps defined");
120     return;
121   }
122   if (!nAxes) {
123     AliError("No axes defined");
124     return;
125   }
126
127   // initialize variable list
128   fVar.Set(nAxes);
129
130   // retrieve number of bins for each axis
131   Int_t *nBins = new Int_t[nAxes];
132   for (iaxis = 0; iaxis < nAxes; iaxis++) 
133   {
134     AliRsnValue *fcnAxis = (AliRsnValue*)fAxisList.At(iaxis);
135     nBins[iaxis] = fcnAxis->GetArray().GetSize() - 1;
136   }
137
138   // create ouput list of containers
139   fContainerList = new TList;
140   fContainerList->SetOwner();
141   fContainerList->SetName(Form("%s_containers", GetName()));
142
143   // initialize output list
144   OpenFile(2);
145   fOutList = new TList();
146   fOutList->SetOwner();
147
148   // create the containers
149   Int_t i = 0, nDef = (Int_t)fPairDefList.GetEntries();
150   for (i = 0; i < nDef; i++) 
151   {
152     AliRsnPairDef *def = (AliRsnPairDef*)fPairDefList[i];
153     AliCFContainer *cont = new AliCFContainer(Form("%s", def->GetPairName()), "", nSteps, nAxes, nBins);
154     // set the bin limits for each axis
155     for (iaxis = 0; iaxis < nAxes; iaxis++) 
156     {
157       AliRsnValue *fcnAxis = (AliRsnValue*)fAxisList.At(iaxis);
158       cont->SetBinLimits(iaxis, fcnAxis->GetArray().GetArray());
159     }
160     // add the container to output list
161     fContainerList->Add(cont);
162   }
163
164   fOutList->Add(fContainerList);
165   fOutList->Print();
166
167   PostData(2, fOutList);
168   
169   // clear heap
170   delete [] nBins;
171
172   AliDebug(AliLog::kDebug+2,"->");
173 }
174
175 //_____________________________________________________________________________
176 void AliRsnAnalysisEffSE::RsnUserExec(Option_t*)
177 {
178 //
179 // Execution of the analysis task.
180 // Recovers the input event and processes it with all included pair objects.
181 // In this case, we NEED to have ESD and MC, otherwise we cannod do anything.
182 //ULTIMO UNDO
183 /*
184   AliRsnDaughter trk;
185   for (Int_t i = 0; i <= AliPID::kSPECIES; i++)
186   {
187     cout << AliPID::ParticleName((AliPID::EParticleType)i) << ": " << endl;
188     for (Int_t m = 0; m < AliRsnDaughter::kMethods; m++)
189     {
190       cout << "-- method: " << AliRsnDaughter::MethodName((AliRsnDaughter::EPIDMethod)m) << endl;
191       Char_t   sign[2] = {'+', '-'};
192       for (Int_t s = 0; s < 2; s++)
193       {
194         TArrayI *a = fRsnPIDIndex.GetTracksArray((AliRsnDaughter::EPIDMethod)m, sign[s], (AliPID::EParticleType)i);
195         Int_t n = a->GetSize();
196         for (Int_t j = 0; j < n; j++)
197         {
198           Int_t k = a->At(j);
199           cout << "-- -- track " << Form("%4d ", k) << ": ";
200           fRsnEvent.SetDaughter(trk, k);
201           cout << "charge = " << (trk.IsPos() ? "+ " : (trk.IsNeg() ? "- " : "0 "));
202           cout << "truePID = " << Form("%10s ", AliPID::ParticleName(trk.PerfectPID()));
203           cout << "realPID = " << Form("%10s ", AliPID::ParticleName(trk.RealisticPID()));
204           cout << endl;
205           cout << "-- -- -- weights (computed): ";
206           for (Int_t q = 0; q < AliPID::kSPECIES; q++)
207             cout << Form("%15.12f", trk.ComputedWeights()[q]) << ' ';
208           cout << endl;
209           cout << "-- -- -- weights (original): ";
210           for (Int_t q = 0; q < AliPID::kSPECIES; q++)
211             cout << Form("%15.12f", trk.GetRef()->PID()[q]) << ' ';
212           cout << endl;
213         }
214       }
215     }
216   } return;
217   */
218
219   AliDebug(AliLog::kDebug+2,"<-");
220   if (!fESDEvent || !fMCEvent) {
221     AliError("This task can process only ESD + MC events");
222     return;
223   }
224   fRsnEvent.SetRef(fESDEvent, fMCEvent);
225
226   // if general event cuts are added to the task (recommended)
227   // they are checked here on the RSN event interface and,
228   // if the event does not pass them, it is skipped and ProcessInfo
229   // is updated accordingly
230   if (!fEventCuts.IsSelected(&fRsnEvent)) {
231     fTaskInfo.SetEventUsed(kFALSE);
232     return;
233   }
234   
235   // if cuts are passed or not cuts were defined,
236   // update the task info before processing the event
237   fTaskInfo.SetEventUsed(kTRUE);
238
239   // process first MC steps and then ESD steps
240   AliRsnPairDef *pairDef = 0;
241   TObjArrayIter iter(&fPairDefList);
242   while ( (pairDef = (AliRsnPairDef*)iter.Next()) )
243   {
244     //ProcessEventMC(pairDef);
245     //ProcessEventESD(pairDef);
246     ProcessEvent(pairDef);
247   }
248
249   // Post the data
250   PostData(2, fOutList);
251
252   AliDebug(AliLog::kDebug+2,"->");
253 }
254
255 //_____________________________________________________________________________
256 void AliRsnAnalysisEffSE::ProcessEvent(AliRsnPairDef *pairDef)
257 {
258 //
259 // Process current event with the definitions of the specified step in MC list
260 // and store results in the container slot defined in second argument.
261 // It is associated with the AliCFContainer with the name of the pair.
262 //
263
264   AliStack      *stack = fRsnEvent.GetRefMC()->Stack();
265   AliESDEvent   *esd   = fRsnEvent.GetRefESD();
266   AliMCEvent    *mc    = fRsnEvent.GetRefMC();
267   AliMCParticle *mother;
268
269   if (!pairDef) return;
270   AliCFContainer *cont = (AliCFContainer*)fContainerList->FindObject(pairDef->GetPairName());
271   if (!cont) return;
272   
273   // get informations from pairDef
274   Int_t pdgM = 0, pdgD[2] = {0, 0};
275   Short_t chargeD[2] = {0, 0};
276   pdgM    = pairDef->GetMotherPDG();
277   pdgD[0] = AliPID::ParticleCode(pairDef->GetPID(0));
278   pdgD[1] = AliPID::ParticleCode(pairDef->GetPID(1));
279   chargeD[0] = pairDef->GetChargeShort(0);
280   chargeD[1] = pairDef->GetChargeShort(1);
281
282   // other utility variables
283   Int_t      first, j, ipart;
284   Int_t      label[2] = {-1, -1};
285   Short_t    charge[2] = {0, 0};
286   Short_t    pairDefMatch[2] = {-1, -1};
287   Int_t      esdIndex[2] = {-1, -1};
288   TParticle *part[2] = {0, 0};
289
290   // in this case, we first take the resonance from MC
291   // and then we find its daughters and compute cuts on them
292   for (ipart = 0; ipart < stack->GetNprimary(); ipart++) 
293   {
294     // take a track from the MC event
295     mother = (AliMCParticle*) fMCEvent->GetTrack(ipart);
296     
297     // check that it is a binary decay and the PDG code matches
298     if (mother->Particle()->GetNDaughters() != 2) continue;
299     if (mother->Particle()->GetPdgCode() != pdgM) continue;
300
301     // store the labels of the two daughters
302     label[0] = mother->Particle()->GetFirstDaughter();
303     label[1] = mother->Particle()->GetLastDaughter();
304     
305     // retrieve the particles and other stuff
306     // check if they match the order in the pairDef
307     for (j = 0; j < 2; j++)
308     {
309       if (label[j] < 0) continue;
310       part[j]   = stack->Particle(label[j]);
311       pdgD[j]   = TMath::Abs(part[j]->GetPdgCode());
312       charge[j] = (Short_t)(part[j]->GetPDG()->Charge() / 3);
313       if (pdgD[j] == AliPID::ParticleCode(pairDef->GetPID(0)) && charge[j] == pairDef->GetChargeShort(0))
314         pairDefMatch[j] = 0;
315       else if (pdgD[j] == AliPID::ParticleCode(pairDef->GetPID(1)) && charge[j] == pairDef->GetChargeShort(1))
316         pairDefMatch[j] = 1;
317       else
318         pairDefMatch[j] = -1;
319         
320       // find corresponding ESD particle: first try rejecting fakes,
321       // and in case of failure, try accepting fakes
322       esdIndex[j] = FindESDtrack(label[j], esd, kTRUE);
323       //TArrayI idx = FindESDtracks(label[j], esd);
324       //for (Int_t kk = 0; kk < idx.GetSize(); kk++) cout << "DAUGHTER " << j << " --> FOUND INDEX: " << idx[kk] << endl;
325       if (esdIndex[j] < 0) esdIndex[j] = FindESDtrack(label[j], esd, kFALSE);
326       //cout << "DAUGHTER " << j << " SINGLE FOUND INDEX = " << esdIndex[j] << endl;
327     }
328     
329     // since each candidate good resonance is taken once, we must check
330     // that it matches the pair definition in any order, and reject in case
331     // in none of them the pair is OK
332     // anyway, we must associate the correct daughter to the correct data member
333     if (pairDefMatch[0] == 0 && pairDefMatch[1] == 1)
334     {
335       // 1st track --> 1st member of PairDef
336       fDaughter[0].SetRef(mc->GetTrack(label[0]));
337       fDaughter[0].SetRefMC((AliMCParticle*)mc->GetTrack(label[0]));
338       fDaughter[0].SetGood();
339       // 2nd track --> 2nd member of PairDef
340       fDaughter[1].SetRef(mc->GetTrack(label[1]));
341       fDaughter[1].SetRefMC((AliMCParticle*)mc->GetTrack(label[1]));
342       fDaughter[1].SetGood();
343     }
344     else if ((pairDefMatch[0] == 1 && pairDefMatch[1] == 0))
345     {
346       // 1st track --> 2nd member of PairDef
347       fDaughter[0].SetRef(mc->GetTrack(label[1]));
348       fDaughter[0].SetRefMC((AliMCParticle*)mc->GetTrack(label[1]));
349       fDaughter[0].SetGood();
350       // 2nd track --> 1st member of PairDef
351       fDaughter[1].SetRef(mc->GetTrack(label[0]));
352       fDaughter[1].SetRefMC((AliMCParticle*)mc->GetTrack(label[0]));
353       fDaughter[1].SetGood();
354     }
355     else
356     {
357       fDaughter[0].SetBad();
358       fDaughter[1].SetBad();
359     }
360     
361     // reject the pair if the matching was unsuccessful
362     if (!fDaughter[0].IsOK() || !fDaughter[1].IsOK()) continue;
363     
364     // first, we set the internal AliRsnMother object to
365     // the MC particles and then operate the selections on MC
366     fPair.SetDaughters(&fDaughter[0], pairDef->GetMass(0), &fDaughter[1], pairDef->GetMass(1));
367     FillContainer(cont, &fStepListMC, pairDef, 0);
368     
369     // then, if both particles found a good match in the ESD
370     // reassociate the ESD tracks to the pair and fill ESD containers
371     if (esdIndex[0] < 0 || esdIndex[1] < 0) continue;
372     if (pairDefMatch[0] == 0 && pairDefMatch[1] == 1)
373     {
374       // 1st track --> 1st member of PairDef
375       fDaughter[0].SetRef(esd->GetTrack(esdIndex[0]));
376       // 2nd track --> 2nd member of PairDef
377       fDaughter[1].SetRef(esd->GetTrack(esdIndex[1]));
378       //cout << "****** MATCHING SCHEME 1" << endl;
379     }
380     else if ((pairDefMatch[0] == 1 && pairDefMatch[1] == 0))
381     {
382       // 1st track --> 2nd member of PairDef
383       fDaughter[0].SetRef(esd->GetTrack(esdIndex[1]));
384       // 2nd track --> 1st member of PairDef
385       fDaughter[1].SetRef(esd->GetTrack(esdIndex[0]));
386       //cout << "****** MATCHING SCHEME 2" << endl;
387     }
388     //cout << "****** IDs = " << fDaughter[0].GetID() << ' ' << fDaughter[1].GetID() << endl;
389     // here we must remember how many steps were already filled
390     first = (Int_t)fStepListMC.GetEntries();
391     FillContainer(cont, &fStepListESD, pairDef, first);
392   }
393 }
394
395 //_____________________________________________________________________________
396 void AliRsnAnalysisEffSE::ProcessEventMC(AliRsnPairDef */*pairDef*/)
397 {
398 //
399 // Process current event with the definitions of the specified step in MC list
400 // and store results in the container slot defined in second argument
401 //
402
403   /*
404   AliStack      *stack = fMCEvent->Stack();
405   AliMCParticle *mother, *daughter;
406
407   if (!pairDef) return;
408   AliCFContainer *cont = (AliCFContainer*)fContainerList->FindObject(pairDef->GetPairName().Data());
409   if (!cont) return;
410
411   // other utility variables
412   Int_t i[2], j, ipart;
413
414   // in this case, we first take the resonance from MC
415   // and then we find its daughters and compute cuts on them
416   for (ipart = 0; ipart < stack->GetNprimary(); ipart++) 
417   {
418     mother = (AliMCParticle*) fMCEvent->GetTrack(ipart);
419     if (mother->Particle()->GetNDaughters() != 2) continue;
420
421     i[0] = mother->Particle()->GetFirstDaughter();
422     i[1] = mother->Particle()->GetLastDaughter();
423
424     for (j = 0; j < 2; j++) 
425     {
426       daughter = (AliMCParticle*) fMCEvent->GetTrack(i[j]);
427       fDaughter[j].SetRef(daughter);
428       fDaughter[j].SetRefMC(daughter);
429       fDaughter[j].FindMotherPDG(stack);
430     }
431
432     if (fDaughter[0].ChargeC() != pairDef->GetCharge(0)) continue;
433     if (fDaughter[1].ChargeC() != pairDef->GetCharge(1)) continue;
434     if (fDaughter[0].PerfectPID() != pairDef->GetType(0)) continue;
435     if (fDaughter[1].PerfectPID() != pairDef->GetType(1)) continue;
436
437     fPair.SetDaughters(&fDaughter[0], &fDaughter[1]);
438
439     // create pair
440     FillContainer(cont, &fStepListMC, pairDef, 0);
441   }
442   */
443 }
444
445 //_____________________________________________________________________________
446 void AliRsnAnalysisEffSE::ProcessEventESD(AliRsnPairDef */*pairDef*/)
447 {
448 //
449 // Process current event with the definitions of the specified step in ESD list
450 // and store results in the container slot defined in second argument
451 //
452
453   /*
454   Int_t i0, i1, first = (Int_t)fStepListMC.GetEntries();
455
456   if (!pairDef) return;
457   AliCFContainer *cont = (AliCFContainer*)fContainerList->FindObject(pairDef->GetPairName().Data());
458   if (!cont) return;
459
460   // external loop on tracks
461   for (i0 = 0; i0 < a0->GetSize(); i0++) {
462     // connect interface
463     fRsnEvent.SetDaughter(fDaughter[0], a0->At(i0));
464     if (!fDaughter[0].IsOK()) continue;
465     fDaughter[0].SetRequiredPID(pairDef->GetType(0));
466
467     // internal loop on tracks
468     for (i1 = 0; i1 < a1->GetSize(); i1++) {
469       // connect interface
470       fRsnEvent.SetDaughter(fDaughter[1], a1->At(i1));
471       if (!fDaughter[1].IsOK()) continue;
472       fDaughter[1].SetRequiredPID(pairDef->GetType(1));
473       // build pair
474       fPair.SetPair(&fDaughter[0], &fDaughter[1]);
475       if (TMath::Abs(fPair.CommonMother()) != pairDef->GetMotherPDG()) continue;
476       // fill containers
477       FillContainer(cont, &fStepListESD, pairDef, first);
478     }
479   }
480   */
481 }
482
483 //_____________________________________________________________________________
484 void AliRsnAnalysisEffSE::FillContainer(AliCFContainer *cont, const TObjArray *stepList, AliRsnPairDef *pd, Int_t firstOutStep)
485 {
486 //
487 // Fill the containers
488 //
489
490   Int_t iaxis, nAxes  = fAxisList.GetEntries();
491   Int_t istep, nSteps = stepList->GetEntries();
492   
493   // set daughters to pair
494   fPair.SetDaughters(&fDaughter[0], pd->GetMass(0), &fDaughter[1], pd->GetMass(1));
495
496   // compute values for all axes
497   for (iaxis = 0; iaxis < nAxes; iaxis++) 
498   {
499     AliRsnValue *fcnAxis = (AliRsnValue*)fAxisList.At(iaxis);
500     fVar[iaxis] = -1E10;
501     if (fcnAxis->Eval(&fPair, pd, &fRsnEvent)) fVar[iaxis] = (Float_t)fcnAxis->GetValue();
502   }
503
504   // fill all steps
505   for (istep = 0; istep < nSteps; istep++) 
506   {
507     AliRsnCutManager *cutMgr = (AliRsnCutManager*)stepList->At(istep);
508     cutMgr->SetEvent(&fRsnEvent);
509     if (!cutMgr->PassCommonDaughterCuts(&fDaughter[0])) break;
510     if (!cutMgr->PassCommonDaughterCuts(&fDaughter[1])) break;
511     if (!cutMgr->PassDaughter1Cuts(&fDaughter[0])) break;
512     if (!cutMgr->PassDaughter2Cuts(&fDaughter[1])) break;
513     if (!cutMgr->PassMotherCuts(&fPair)) break;
514     //cout << "**************************************** FILLING STEP " << istep << endl;
515     cont->Fill(fVar.GetArray(), istep + firstOutStep);
516   }
517 }
518
519 //_____________________________________________________________________________
520 void AliRsnAnalysisEffSE::RsnTerminate(Option_t*)
521 {
522 //
523 // Termination.
524 // Could be added some monitor histograms here.
525 //
526
527   AliDebug(AliLog::kDebug+2,"<-");
528   AliDebug(AliLog::kDebug+2,"->");
529 }
530
531 //_____________________________________________________________________________
532 void AliRsnAnalysisEffSE::AddPairDef(AliRsnPairDef* pairDef)
533 {
534 //
535 //  Adds pair definition
536 //
537   fPairDefList.AddLast(pairDef);
538 }
539
540 //_____________________________________________________________________________
541 Int_t AliRsnAnalysisEffSE::FindESDtrack(Int_t label, AliESDEvent *esd, Bool_t rejectFakes)
542 {
543 //
544 // Finds in the ESD a track whose label corresponds to that in argument.
545 // When global tracks are enabled, tries first to find a global track 
546 // satisfying that requirement.
547 // If no global tracks are found, if ITS-SA are enable, tries to search among them
548 // otherwise return a negative number.
549 // If global tracks are disabled, search only among ITS SA
550 //
551
552   Int_t   i = 0;
553   Int_t   ntracks = esd->GetNumberOfTracks();
554   ULong_t status;
555   Bool_t  isTPC;
556   Bool_t  isITSSA;
557   
558   // loop for global tracks
559   if (fUseGlobal)
560   {
561     for (i = 0; i < ntracks; i++)
562     {
563       AliESDtrack *track = esd->GetTrack(i);
564       status  = (ULong_t)track->GetStatus();
565       isTPC   = ((status & AliESDtrack::kTPCin)  != 0);
566       if (!isTPC) continue;
567       
568       // check that label match
569       if (TMath::Abs(track->GetLabel()) != label) continue;
570       
571       // if required, reject fakes
572       if (rejectFakes && track->GetLabel() < 0) continue;
573       
574       // if all checks are passed and we are searching among global
575       // this means that thie track is a global one with the right label
576       // then, the return value is set to this, and returned
577       return i;
578     }
579   }
580   
581   // loop for ITS-SA tracks (this happens only if no global tracks were found
582   // or searching among globals is disabled)
583   if (fUseITSSA)
584   {
585     for (i = 0; i < ntracks; i++)
586     {
587       AliESDtrack *track = esd->GetTrack(i);
588       status  = (ULong_t)track->GetStatus();
589       isITSSA = ((status & AliESDtrack::kTPCin)  == 0 && (status & AliESDtrack::kITSrefit) != 0 && (status & AliESDtrack::kITSpureSA) == 0 && (status & AliESDtrack::kITSpid) != 0);
590       if (!isITSSA) continue;
591       
592       // check that label match
593       if (TMath::Abs(track->GetLabel()) != label) continue;
594             
595       // if required, reject fakes
596       if (rejectFakes && track->GetLabel() < 0) continue;
597       
598       // if all checks are passed and we are searching among global
599       // this means that thie track is a global one with the right label
600       // then, the return value is set to this, and returned
601       return i;
602     }
603   }
604   
605   // if we reach this point, no match were found
606   return -1;
607 }
608
609 //_____________________________________________________________________________
610 TArrayI AliRsnAnalysisEffSE::FindESDtracks(Int_t label, AliESDEvent *esd)
611 {
612 //
613 // Finds in the ESD a track whose label corresponds to that in argument.
614 // When global tracks are enabled, tries first to find a global track 
615 // satisfying that requirement.
616 // If no global tracks are found, if ITS-SA are enable, tries to search among them
617 // otherwise return a negative number.
618 // If global tracks are disabled, search only among ITS SA
619 //
620
621   Int_t   i = 0;
622   Int_t   ntracks = esd->GetNumberOfTracks();
623   ULong_t status;
624   Bool_t  isTPC;
625   Bool_t  isITSSA;
626   TArrayI array(100);
627   Int_t   nfound = 0;
628   
629   // loop for global tracks
630   if (fUseGlobal)
631   {
632     for (i = 0; i < ntracks; i++)
633     {
634       AliESDtrack *track = esd->GetTrack(i);
635       status  = (ULong_t)track->GetStatus();
636       isTPC   = ((status & AliESDtrack::kTPCin)  != 0);
637       if (!isTPC) continue;
638       
639       // check that label match
640       if (TMath::Abs(track->GetLabel()) != label) continue;
641       
642       array[nfound++] = i;
643     }
644   }
645   
646   // loop for ITS-SA tracks (this happens only if no global tracks were found
647   // or searching among globals is disabled)
648   if (fUseITSSA)
649   {
650     for (i = 0; i < ntracks; i++)
651     {
652       AliESDtrack *track = esd->GetTrack(i);
653       status  = (ULong_t)track->GetStatus();
654       isITSSA = ((status & AliESDtrack::kTPCin)  == 0 && (status & AliESDtrack::kITSrefit) != 0 && (status & AliESDtrack::kITSpureSA) == 0 && (status & AliESDtrack::kITSpid) != 0);
655       if (!isITSSA) continue;
656       
657       // check that label match
658       if (TMath::Abs(track->GetLabel()) != label) continue;
659             
660       array[nfound++] = i;
661     }
662   }
663   
664   array.Set(nfound);
665   return array;
666 }