]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnAnalysisEffSE.cxx
Corrections after survey of Coverity tests
[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   label[2] = {-1, -1}, first, j, ipart;
284   Short_t charge[2] = {0, 0};
285   Short_t pairDefMatch[2] = {-1, -1};
286   Int_t   esdIndex[2];
287   TParticle *part[2] = {0, 0};
288
289   // in this case, we first take the resonance from MC
290   // and then we find its daughters and compute cuts on them
291   for (ipart = 0; ipart < stack->GetNprimary(); ipart++) 
292   {
293     // take a track from the MC event
294     mother = (AliMCParticle*) fMCEvent->GetTrack(ipart);
295     
296     // check that it is a binary decay and the PDG code matches
297     if (mother->Particle()->GetNDaughters() != 2) continue;
298     if (mother->Particle()->GetPdgCode() != pdgM) continue;
299
300     // store the labels of the two daughters
301     label[0] = mother->Particle()->GetFirstDaughter();
302     label[1] = mother->Particle()->GetLastDaughter();
303     
304     // retrieve the particles and other stuff
305     // check if they match the order in the pairDef
306     for (j = 0; j < 2; j++)
307     {
308       if (label[j] < 0) continue;
309       part[j]   = stack->Particle(label[j]);
310       pdgD[j]   = TMath::Abs(part[j]->GetPdgCode());
311       charge[j] = (Short_t)(part[j]->GetPDG()->Charge() / 3);
312       if (pdgD[j] == AliPID::ParticleCode(pairDef->GetPID(0)) && charge[j] == pairDef->GetChargeShort(0))
313         pairDefMatch[j] = 0;
314       else if (pdgD[j] == AliPID::ParticleCode(pairDef->GetPID(1)) && charge[j] == pairDef->GetChargeShort(1))
315         pairDefMatch[j] = 1;
316       else
317         pairDefMatch[j] = -1;
318         
319       // find corresponding ESD particle: first try rejecting fakes,
320       // and in case of failure, try accepting fakes
321       esdIndex[j] = FindESDtrack(label[j], esd, kTRUE);
322       //TArrayI idx = FindESDtracks(label[j], esd);
323       //for (Int_t kk = 0; kk < idx.GetSize(); kk++) cout << "DAUGHTER " << j << " --> FOUND INDEX: " << idx[kk] << endl;
324       if (esdIndex[j] < 0) esdIndex[j] = FindESDtrack(label[j], esd, kFALSE);
325       //cout << "DAUGHTER " << j << " SINGLE FOUND INDEX = " << esdIndex[j] << endl;
326     }
327     
328     // since each candidate good resonance is taken once, we must check
329     // that it matches the pair definition in any order, and reject in case
330     // in none of them the pair is OK
331     // anyway, we must associate the correct daughter to the correct data member
332     if (pairDefMatch[0] == 0 && pairDefMatch[1] == 1)
333     {
334       // 1st track --> 1st member of PairDef
335       fDaughter[0].SetRef(mc->GetTrack(label[0]));
336       fDaughter[0].SetRefMC((AliMCParticle*)mc->GetTrack(label[0]));
337       fDaughter[0].SetGood();
338       // 2nd track --> 2nd member of PairDef
339       fDaughter[1].SetRef(mc->GetTrack(label[1]));
340       fDaughter[1].SetRefMC((AliMCParticle*)mc->GetTrack(label[1]));
341       fDaughter[1].SetGood();
342     }
343     else if ((pairDefMatch[0] == 1 && pairDefMatch[1] == 0))
344     {
345       // 1st track --> 2nd member of PairDef
346       fDaughter[0].SetRef(mc->GetTrack(label[1]));
347       fDaughter[0].SetRefMC((AliMCParticle*)mc->GetTrack(label[1]));
348       fDaughter[0].SetGood();
349       // 2nd track --> 1st member of PairDef
350       fDaughter[1].SetRef(mc->GetTrack(label[0]));
351       fDaughter[1].SetRefMC((AliMCParticle*)mc->GetTrack(label[0]));
352       fDaughter[1].SetGood();
353     }
354     else
355     {
356       fDaughter[0].SetBad();
357       fDaughter[1].SetBad();
358     }
359     
360     // reject the pair if the matching was unsuccessful
361     if (!fDaughter[0].IsOK() || !fDaughter[1].IsOK()) continue;
362     
363     // first, we set the internal AliRsnMother object to
364     // the MC particles and then operate the selections on MC
365     fPair.SetDaughters(&fDaughter[0], pairDef->GetMass(0), &fDaughter[1], pairDef->GetMass(1));
366     FillContainer(cont, &fStepListMC, pairDef, 0);
367     
368     // then, if both particles found a good match in the ESD
369     // reassociate the ESD tracks to the pair and fill ESD containers
370     if (esdIndex[0] < 0 || esdIndex[1] < 0) continue;
371     if (pairDefMatch[0] == 0 && pairDefMatch[1] == 1)
372     {
373       // 1st track --> 1st member of PairDef
374       fDaughter[0].SetRef(esd->GetTrack(esdIndex[0]));
375       // 2nd track --> 2nd member of PairDef
376       fDaughter[1].SetRef(esd->GetTrack(esdIndex[1]));
377       //cout << "****** MATCHING SCHEME 1" << endl;
378     }
379     else if ((pairDefMatch[0] == 1 && pairDefMatch[1] == 0))
380     {
381       // 1st track --> 2nd member of PairDef
382       fDaughter[0].SetRef(esd->GetTrack(esdIndex[1]));
383       // 2nd track --> 1st member of PairDef
384       fDaughter[1].SetRef(esd->GetTrack(esdIndex[0]));
385       //cout << "****** MATCHING SCHEME 2" << endl;
386     }
387     //cout << "****** IDs = " << fDaughter[0].GetID() << ' ' << fDaughter[1].GetID() << endl;
388     // here we must remember how many steps were already filled
389     first = (Int_t)fStepListMC.GetEntries();
390     FillContainer(cont, &fStepListESD, pairDef, first);
391   }
392 }
393
394 //_____________________________________________________________________________
395 void AliRsnAnalysisEffSE::ProcessEventMC(AliRsnPairDef */*pairDef*/)
396 {
397 //
398 // Process current event with the definitions of the specified step in MC list
399 // and store results in the container slot defined in second argument
400 //
401
402   /*
403   AliStack      *stack = fMCEvent->Stack();
404   AliMCParticle *mother, *daughter;
405
406   if (!pairDef) return;
407   AliCFContainer *cont = (AliCFContainer*)fContainerList->FindObject(pairDef->GetPairName().Data());
408   if (!cont) return;
409
410   // other utility variables
411   Int_t i[2], j, ipart;
412
413   // in this case, we first take the resonance from MC
414   // and then we find its daughters and compute cuts on them
415   for (ipart = 0; ipart < stack->GetNprimary(); ipart++) 
416   {
417     mother = (AliMCParticle*) fMCEvent->GetTrack(ipart);
418     if (mother->Particle()->GetNDaughters() != 2) continue;
419
420     i[0] = mother->Particle()->GetFirstDaughter();
421     i[1] = mother->Particle()->GetLastDaughter();
422
423     for (j = 0; j < 2; j++) 
424     {
425       daughter = (AliMCParticle*) fMCEvent->GetTrack(i[j]);
426       fDaughter[j].SetRef(daughter);
427       fDaughter[j].SetRefMC(daughter);
428       fDaughter[j].FindMotherPDG(stack);
429     }
430
431     if (fDaughter[0].ChargeC() != pairDef->GetCharge(0)) continue;
432     if (fDaughter[1].ChargeC() != pairDef->GetCharge(1)) continue;
433     if (fDaughter[0].PerfectPID() != pairDef->GetType(0)) continue;
434     if (fDaughter[1].PerfectPID() != pairDef->GetType(1)) continue;
435
436     fPair.SetDaughters(&fDaughter[0], &fDaughter[1]);
437
438     // create pair
439     FillContainer(cont, &fStepListMC, pairDef, 0);
440   }
441   */
442 }
443
444 //_____________________________________________________________________________
445 void AliRsnAnalysisEffSE::ProcessEventESD(AliRsnPairDef */*pairDef*/)
446 {
447 //
448 // Process current event with the definitions of the specified step in ESD list
449 // and store results in the container slot defined in second argument
450 //
451
452   /*
453   Int_t i0, i1, first = (Int_t)fStepListMC.GetEntries();
454
455   if (!pairDef) return;
456   AliCFContainer *cont = (AliCFContainer*)fContainerList->FindObject(pairDef->GetPairName().Data());
457   if (!cont) return;
458
459   // external loop on tracks
460   for (i0 = 0; i0 < a0->GetSize(); i0++) {
461     // connect interface
462     fRsnEvent.SetDaughter(fDaughter[0], a0->At(i0));
463     if (!fDaughter[0].IsOK()) continue;
464     fDaughter[0].SetRequiredPID(pairDef->GetType(0));
465
466     // internal loop on tracks
467     for (i1 = 0; i1 < a1->GetSize(); i1++) {
468       // connect interface
469       fRsnEvent.SetDaughter(fDaughter[1], a1->At(i1));
470       if (!fDaughter[1].IsOK()) continue;
471       fDaughter[1].SetRequiredPID(pairDef->GetType(1));
472       // build pair
473       fPair.SetPair(&fDaughter[0], &fDaughter[1]);
474       if (TMath::Abs(fPair.CommonMother()) != pairDef->GetMotherPDG()) continue;
475       // fill containers
476       FillContainer(cont, &fStepListESD, pairDef, first);
477     }
478   }
479   */
480 }
481
482 //_____________________________________________________________________________
483 void AliRsnAnalysisEffSE::FillContainer(AliCFContainer *cont, const TObjArray *stepList, AliRsnPairDef *pd, Int_t firstOutStep)
484 {
485 //
486 // Fill the containers
487 //
488
489   Int_t iaxis, nAxes  = fAxisList.GetEntries();
490   Int_t istep, nSteps = stepList->GetEntries();
491   
492   // set daughters to pair
493   fPair.SetDaughters(&fDaughter[0], pd->GetMass(0), &fDaughter[1], pd->GetMass(1));
494
495   // compute values for all axes
496   for (iaxis = 0; iaxis < nAxes; iaxis++) 
497   {
498     AliRsnValue *fcnAxis = (AliRsnValue*)fAxisList.At(iaxis);
499     fVar[iaxis] = -1E10;
500     if (fcnAxis->Eval(&fPair, pd, &fRsnEvent)) fVar[iaxis] = (Float_t)fcnAxis->GetValue();
501   }
502
503   // fill all steps
504   for (istep = 0; istep < nSteps; istep++) 
505   {
506     AliRsnCutManager *cutMgr = (AliRsnCutManager*)stepList->At(istep);
507     cutMgr->SetEvent(&fRsnEvent);
508     if (!cutMgr->PassCommonDaughterCuts(&fDaughter[0])) break;
509     if (!cutMgr->PassCommonDaughterCuts(&fDaughter[1])) break;
510     if (!cutMgr->PassDaughter1Cuts(&fDaughter[0])) break;
511     if (!cutMgr->PassDaughter2Cuts(&fDaughter[1])) break;
512     if (!cutMgr->PassMotherCuts(&fPair)) break;
513     //cout << "**************************************** FILLING STEP " << istep << endl;
514     cont->Fill(fVar.GetArray(), istep + firstOutStep);
515   }
516 }
517
518 //_____________________________________________________________________________
519 void AliRsnAnalysisEffSE::RsnTerminate(Option_t*)
520 {
521 //
522 // Termination.
523 // Could be added some monitor histograms here.
524 //
525
526   AliDebug(AliLog::kDebug+2,"<-");
527   AliDebug(AliLog::kDebug+2,"->");
528 }
529
530 //_____________________________________________________________________________
531 void AliRsnAnalysisEffSE::AddPairDef(AliRsnPairDef* pairDef)
532 {
533 //
534 //  Adds pair definition
535 //
536   fPairDefList.AddLast(pairDef);
537 }
538
539 //_____________________________________________________________________________
540 Int_t AliRsnAnalysisEffSE::FindESDtrack(Int_t label, AliESDEvent *esd, Bool_t rejectFakes)
541 {
542 //
543 // Finds in the ESD a track whose label corresponds to that in argument.
544 // When global tracks are enabled, tries first to find a global track 
545 // satisfying that requirement.
546 // If no global tracks are found, if ITS-SA are enable, tries to search among them
547 // otherwise return a negative number.
548 // If global tracks are disabled, search only among ITS SA
549 //
550
551   Int_t   i = 0;
552   Int_t   ntracks = esd->GetNumberOfTracks();
553   ULong_t status;
554   Bool_t  isTPC;
555   Bool_t  isITSSA;
556   
557   // loop for global tracks
558   if (fUseGlobal)
559   {
560     for (i = 0; i < ntracks; i++)
561     {
562       AliESDtrack *track = esd->GetTrack(i);
563       status  = (ULong_t)track->GetStatus();
564       isTPC   = ((status & AliESDtrack::kTPCin)  != 0);
565       if (!isTPC) continue;
566       
567       // check that label match
568       if (TMath::Abs(track->GetLabel()) != label) continue;
569       
570       // if required, reject fakes
571       if (rejectFakes && track->GetLabel() < 0) continue;
572       
573       // if all checks are passed and we are searching among global
574       // this means that thie track is a global one with the right label
575       // then, the return value is set to this, and returned
576       return i;
577     }
578   }
579   
580   // loop for ITS-SA tracks (this happens only if no global tracks were found
581   // or searching among globals is disabled)
582   if (fUseITSSA)
583   {
584     for (i = 0; i < ntracks; i++)
585     {
586       AliESDtrack *track = esd->GetTrack(i);
587       status  = (ULong_t)track->GetStatus();
588       isITSSA = ((status & AliESDtrack::kTPCin)  == 0 && (status & AliESDtrack::kITSrefit) != 0 && (status & AliESDtrack::kITSpureSA) == 0 && (status & AliESDtrack::kITSpid) != 0);
589       if (!isITSSA) continue;
590       
591       // check that label match
592       if (TMath::Abs(track->GetLabel()) != label) continue;
593             
594       // if required, reject fakes
595       if (rejectFakes && track->GetLabel() < 0) continue;
596       
597       // if all checks are passed and we are searching among global
598       // this means that thie track is a global one with the right label
599       // then, the return value is set to this, and returned
600       return i;
601     }
602   }
603   
604   // if we reach this point, no match were found
605   return -1;
606 }
607
608 //_____________________________________________________________________________
609 TArrayI AliRsnAnalysisEffSE::FindESDtracks(Int_t label, AliESDEvent *esd)
610 {
611 //
612 // Finds in the ESD a track whose label corresponds to that in argument.
613 // When global tracks are enabled, tries first to find a global track 
614 // satisfying that requirement.
615 // If no global tracks are found, if ITS-SA are enable, tries to search among them
616 // otherwise return a negative number.
617 // If global tracks are disabled, search only among ITS SA
618 //
619
620   Int_t   i = 0;
621   Int_t   ntracks = esd->GetNumberOfTracks();
622   ULong_t status;
623   Bool_t  isTPC;
624   Bool_t  isITSSA;
625   TArrayI array(100);
626   Int_t   nfound = 0;
627   
628   // loop for global tracks
629   if (fUseGlobal)
630   {
631     for (i = 0; i < ntracks; i++)
632     {
633       AliESDtrack *track = esd->GetTrack(i);
634       status  = (ULong_t)track->GetStatus();
635       isTPC   = ((status & AliESDtrack::kTPCin)  != 0);
636       if (!isTPC) continue;
637       
638       // check that label match
639       if (TMath::Abs(track->GetLabel()) != label) continue;
640       
641       array[nfound++] = i;
642     }
643   }
644   
645   // loop for ITS-SA tracks (this happens only if no global tracks were found
646   // or searching among globals is disabled)
647   if (fUseITSSA)
648   {
649     for (i = 0; i < ntracks; i++)
650     {
651       AliESDtrack *track = esd->GetTrack(i);
652       status  = (ULong_t)track->GetStatus();
653       isITSSA = ((status & AliESDtrack::kTPCin)  == 0 && (status & AliESDtrack::kITSrefit) != 0 && (status & AliESDtrack::kITSpureSA) == 0 && (status & AliESDtrack::kITSpid) != 0);
654       if (!isITSSA) continue;
655       
656       // check that label match
657       if (TMath::Abs(track->GetLabel()) != label) continue;
658             
659       array[nfound++] = i;
660     }
661   }
662   
663   array.Set(nfound);
664   return array;
665 }