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