]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnEvent.cxx
Restored call to CreateDigitizer (F.Prino)
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnEvent.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 ////////////////////////////////////////////////////////////////////////////////
17 //
18 //  This class works as generic interface to an event.
19 //  Its main purpose is to provide a unique reference which includes all the
20 //  facilities available in the AliVEvent generic base class, plus all info
21 //  which could be needed during analysis, which are not in AliVEvent but
22 //  need to be accessed from ESD or AOD objects, usually in different ways.
23 //  When MC is available, it is properly taken into account.
24 //  
25 //  authors: A. Pulvirenti (alberto.pulvirenti@ct.infn.it)
26 //           M. Vala (martin.vala@cern.ch)
27 //
28 ////////////////////////////////////////////////////////////////////////////////
29
30 #include <Riostream.h>
31 #include <TArrayF.h>
32
33 #include "AliLog.h"
34 #include "AliVEvent.h"
35 #include "AliMCEvent.h"
36 #include "AliStack.h"
37 #include "AliGenEventHeader.h"
38 #include "AliAODEvent.h"
39 #include "AliRsnCutPID.h"
40 #include "AliESDtrackCuts.h"
41
42 #include "AliRsnEvent.h"
43
44 ClassImp(AliRsnEvent)
45
46 AliRsnEvent* AliRsnEvent::fgRsnEvent1 = 0;
47 AliRsnEvent* AliRsnEvent::fgRsnEvent2 = 0;
48
49 //_____________________________________________________________________________
50 AliRsnEvent::AliRsnEvent(AliVEvent *ref, AliVEvent *refMC) :
51    fRef(ref),
52    fRefMC(refMC),
53    fLeading(-1),
54    fLocalID(-1)
55 {
56 //
57 // Default constructor.
58 //
59 }
60
61 //_____________________________________________________________________________
62 AliRsnEvent::AliRsnEvent(const AliRsnEvent &event) :
63    TObject(event),
64    fRef(event.fRef),
65    fRefMC(event.fRefMC),
66    fLeading(event.fLeading),
67    fLocalID(event.fLocalID)
68 {
69 //
70 // Copy constructor.
71 //
72 }
73
74 //_____________________________________________________________________________
75 AliRsnEvent& AliRsnEvent::operator= (const AliRsnEvent & event)
76 {
77 //
78 // Works in the same way as the copy constructor.
79 //
80
81    (TObject)(*this) = (TObject)event;
82    fRef             = event.fRef;
83    fRefMC           = event.fRefMC;
84    fLeading         = event.fLeading;
85    fLocalID         = event.fLocalID;
86
87    return (*this);
88 }
89
90 //_____________________________________________________________________________
91 AliRsnEvent::~AliRsnEvent()
92 {
93 //
94 // Destructor.
95 // Dereferences global pointer, if needed.
96 //
97
98    //if (gRsnCurrentEvent == this) gRsnCurrentEvent = 0;
99    //if (gRsnMixedEvent   == this) gRsnMixedEvent = 0;
100 }
101
102 //_____________________________________________________________________________
103 Bool_t AliRsnEvent::SetDaughter(AliRsnDaughter &out, Int_t i, AliRsnDaughter::ERefType type)
104 {
105 //
106 // Using the second and third arguments, retrieves the i-th object in the
107 // appropriate sample (tracks or V0s) and sets the first reference object
108 // in order to point to that.
109 // If a MonteCarlo information is provided, sets the useful informations from there,
110 // and in case of a V0, sets the 'label' data member only when the two daughters
111 // of the V0 point to the same mother.
112 // Returns kFALSE whenever the operation fails (out of range, NULL references).
113 //
114
115    Bool_t ok = kFALSE;
116
117    if (IsESD() && type == AliRsnDaughter::kTrack)   ok = SetDaughterESDtrack  (out, i);
118    if (IsAOD() && type == AliRsnDaughter::kTrack)   ok = SetDaughterAODtrack  (out, i);
119    if (IsESD() && type == AliRsnDaughter::kV0)      ok = SetDaughterESDv0     (out, i);
120    if (IsAOD() && type == AliRsnDaughter::kV0)      ok = SetDaughterAODv0     (out, i);
121    if (IsESD() && type == AliRsnDaughter::kCascade) ok = SetDaughterESDcascade(out, i);
122    if (IsAOD() && type == AliRsnDaughter::kCascade) ok = SetDaughterAODcascade(out, i);
123
124    return ok;
125 }
126
127 //_____________________________________________________________________________
128 Bool_t AliRsnEvent::SetDaughterAbs(AliRsnDaughter &out, Int_t absIndex)
129 {
130 //
131 // Sets the first argument daughter using the absolute index, which
132 // runs continuously from tracks, to V0s, to cascades.
133 // In case the conversion to real index fails, the track is flagged as bad.
134 // Additionally, sets the daughter internal 'fRsnID' member to this index.
135 //
136
137    Int_t index;
138    AliRsnDaughter::ERefType type;
139
140    if (!ConvertAbsoluteIndex(absIndex, index, type)) {
141       out.Reset();
142       out.SetBad();
143       return kFALSE;
144    }
145    
146    SetDaughter(out, index, type);
147    out.SetRsnID(absIndex);
148    return kTRUE;
149 }
150
151 //_____________________________________________________________________________
152 Bool_t AliRsnEvent::SetDaughterMC(AliRsnDaughter &out, Int_t label)
153 {
154 //
155 // Using the second argument, retrieves the i-th object in the
156 // MC sample (if present) and assigns the track using only that,
157 // so that it is considered both as main reference and MC reference.
158 // (used for MC-only analysis).
159 //
160
161    if (!fRefMC) {
162       out.SetBad();
163       return kFALSE;
164    }
165
166    // try to convert into both types
167    Int_t        imum;
168    AliMCEvent  *esd = GetRefMCESD();
169    AliAODEvent *aod = GetRefMCAOD();
170
171    // ESD
172    if (esd) {
173       // if the MC track exists, assign it
174       AliMCParticle *track = (AliMCParticle*)fRef->GetTrack(label);
175       if (!track) {
176          out.SetBad();
177          return kFALSE;
178       }
179       out.SetRef(track);
180       out.SetRefMC(track);
181       out.SetLabel(label);
182       out.SetGood();
183
184       // search for its mother in stack
185       imum = track->GetMother();
186       if (imum >= 0 && imum < esd->Stack()->GetNtrack()) {
187          TParticle *mum = esd->Stack()->Particle(imum);
188          if (mum) out.SetMotherPDG(TMath::Abs(mum->GetPdgCode()));
189       }
190    }
191
192    // AOD
193    if (aod) {
194       // checks that array of MC particles exists
195       TClonesArray *mcArray = (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
196       if (!mcArray) {
197          out.SetBad();
198          return kFALSE;
199       }
200
201       // in this case one has to loop over the sample to find the good one
202       TObjArrayIter next(mcArray);
203       AliAODMCParticle *part = 0x0;
204       while ((part = (AliAODMCParticle*)next())) {
205          if (TMath::Abs(part->GetLabel()) == label) {
206             // if the MC track exists, assign it
207             out.SetRef(part);
208             out.SetRefMC(part);
209             out.SetLabel(label);
210             out.SetGood();
211
212             // search for the mother
213             imum = part->GetMother();
214             if (imum >= 0 && imum < mcArray->GetEntriesFast()) {
215                AliAODMCParticle *mum = (AliAODMCParticle*)mcArray->At(imum);
216                if (mum) out.SetMotherPDG(TMath::Abs(mum->GetPdgCode()));
217             }
218             break;
219          }
220       }
221       return kTRUE;
222    }
223
224    return kFALSE;
225 }
226
227 //_____________________________________________________________________________
228 AliRsnDaughter AliRsnEvent::GetDaughter(Int_t i, AliRsnDaughter::ERefType type)
229 {
230 //
231 // Returns a daughter set using same criteria as SetDaughter
232 //
233
234    AliRsnDaughter d;
235    SetDaughter(d, i, type);
236    return d;
237 }
238
239 //_____________________________________________________________________________
240 AliRsnDaughter AliRsnEvent::GetDaughterAbs(Int_t absIndex)
241 {
242 //
243 // Returns a daughter set using same criteria as SetDaughter
244 //
245
246    AliRsnDaughter d;
247    SetDaughterAbs(d, absIndex);
248    return d;
249 }
250
251 //_____________________________________________________________________________
252 AliRsnDaughter AliRsnEvent::GetDaughterMC(Int_t i)
253 {
254 //
255 // Returns a daughter set using same criteria as SetDaughterMC
256 //
257
258    AliRsnDaughter d;
259    SetDaughterMC(d, i);
260    return d;
261 }
262
263 //_____________________________________________________________________________
264 Int_t AliRsnEvent::GetAbsoluteSum()
265 {
266 //
267 // Utility method that returns the sum of all daughter-like objects:
268 // tracks, V0s and cascades
269 //
270
271    Int_t total = 0;
272
273    total += fRef->GetNumberOfTracks();
274    total += fRef->GetNumberOfV0s();
275    total += fRef->GetNumberOfCascades();
276
277    return total;
278 }
279
280 //_____________________________________________________________________________
281 Bool_t AliRsnEvent::ConvertAbsoluteIndex(Int_t index, Int_t &realIndex, AliRsnDaughter::ERefType &type)
282 {
283 //
284 // Using the phylosophy of the absolute index, which loops over
285 // all tracks, V0s and cascades, returns the result of a check
286 // on it (first argument) based on this criterion:
287 // 1) if the absolute ID is smaller than number of tracks,
288 //    return itself and the type 'track'
289 // 2) if the absolute ID is larger than number of tracks, subtract it
290 //    and if the result is smaller than number of V0s,
291 //    return the corresponding V0 index and type
292 // 3) if the absolute ID is larger than number of tracks + V0s, subtract them
293 //    and if the result is smaller than number of cascades,
294 //    return the corresponding cascade index and type
295 // The results of this check are stored in the reference arguments, while the outcome of
296 // the function is kTRUE if one of these checks was successful, otherwise it is kFALSE,
297 // meaning that the absolute index reached the end.
298 //
299
300    Int_t nTracks   = fRef->GetNumberOfTracks();
301    Int_t nV0s      = fRef->GetNumberOfV0s();
302    Int_t nCascades = fRef->GetNumberOfCascades();
303
304    if (index < nTracks) {
305       realIndex = index;
306       type = AliRsnDaughter::kTrack;
307       return kTRUE;
308    } else if (index >= nTracks && index < nTracks + nV0s) {
309       realIndex = index - nTracks;
310       type = AliRsnDaughter::kV0;
311       return kTRUE;
312    } else if (index >= nTracks + nV0s && index < nTracks + nV0s + nCascades) {
313       realIndex = index - nTracks - nV0s;
314       type = AliRsnDaughter::kCascade;
315       return kTRUE;
316    }
317
318    realIndex = -1;
319    type = AliRsnDaughter::kNoType;
320    return kFALSE;
321 }
322
323 //_____________________________________________________________________________
324 Int_t AliRsnEvent::ConvertRealIndex(Int_t index, AliRsnDaughter::ERefType type)
325 {
326 //
327 // Translates a pair made by index + object type into the corresponding
328 // absolute index, which is set to -1 in case the real index overflows.
329 //
330
331    Int_t nTracks   = fRef->GetNumberOfTracks();
332    Int_t nV0s      = fRef->GetNumberOfV0s();
333    Int_t nCascades = fRef->GetNumberOfCascades();
334
335    switch (type) {
336       case AliRsnDaughter::kTrack:
337          if (index >= 0 && index < nTracks)
338             return index;
339          else
340             return -1;
341       case AliRsnDaughter::kV0:
342          if (index >= 0 && index < nV0s)
343             return nTracks + index;
344          else
345             return -1;
346       case AliRsnDaughter::kCascade:
347          if (index >= 0 && index < nCascades)
348             return nTracks + nV0s + index;
349          else
350             return -1;
351       default:
352          return -1;
353    }
354 }
355
356 //_____________________________________________________________________________
357 Int_t AliRsnEvent::GetMultiplicity(AliESDtrackCuts *cuts)
358 {
359 //
360 // Returns event multiplicity as the number of tracks.
361 // If the argument is not NULL, returns instead the
362 // number of tracks passing the cuts hereby defined.
363 //
364
365    if (!fRef) return 0;
366
367    AliESDEvent *esd = GetRefESD();
368    if (cuts && esd) return cuts->CountAcceptedTracks(esd);
369    else return fRef->GetNumberOfTracks();
370 }
371
372 //_____________________________________________________________________________
373 Int_t AliRsnEvent::GetMultiplicityMC()
374 {
375 //
376 // Returns event multiplicity as the number of primary tracks
377 // counted from the MC sample, if present.
378 //
379
380    if (!fRefMC) return 0;
381
382    if (IsESD())
383       return GetRefMCESD()->Stack()->GetNprimary();
384    else if (IsAOD())
385       return GetRefMCAOD()->GetNumberOfTracks();
386    else
387       return 0;
388 }
389
390 //_____________________________________________________________________________
391 Double_t AliRsnEvent::GetVz()
392 {
393 //
394 // Return Z coord of primary vertex
395 //
396
397    return fRef->GetPrimaryVertex()->GetZ();
398 }
399
400 //_____________________________________________________________________________
401 Int_t AliRsnEvent::SelectLeadingParticle
402 (Double_t ptMin, AliRsnCutPID *cutPID)
403 {
404 //
405 // Searches the collection of all particles with given PID type and charge,
406 // and returns the one with largest momentum, provided that it is greater than 1st argument.
407 // If one specifies AliRsnPID::kUnknown as type or AliRsnDaughter::kNoPID as method,
408 // the check is done over all particles irrespectively of their PID.
409 // If the sign argument is '+' or '-', the check is done over the particles of that charge,
410 // otherwise it is done irrespectively of the charge.
411 //
412
413    Int_t i, nTracks = fRef->GetNumberOfTracks();
414    fLeading = -1;
415    AliRsnDaughter leading;
416    leading.SetBad();
417
418    for (i = 0; i < nTracks; i++) {
419       AliRsnDaughter track = GetDaughter(i);
420       if (cutPID) if (!cutPID->IsSelected(&track)) continue;
421       const AliVParticle *ref = track.GetRef();
422       if (ref->Pt() < ptMin) continue;
423       //double pt = track.P().Perp();
424       //Printf("track %d %g", i, pt);
425       if (!leading.IsOK() || ref->Pt() > ptMin) {
426          fLeading = i;
427          //leading = track;
428          ptMin = ref->Pt();
429       }
430    }
431    return fLeading;
432 }
433
434 //_________________________________________________________________________________________________
435 Double_t AliRsnEvent::GetAverageMomentum(Int_t &count, AliRsnCutPID *cutPID)
436 {
437 //
438 // Loops on the list of tracks and computes average total momentum.
439 //
440
441    Int_t i, nTracks = fRef->GetNumberOfTracks();
442    Double_t pmean = 0.0;
443
444    for (i = 0, count = 0; i < nTracks; i++) {
445       AliRsnDaughter track = GetDaughter(i);
446       if (cutPID) if (!cutPID->IsSelected(&track)) continue;
447       pmean += track.Prec().Mag();
448       count++;
449    }
450
451    if (count > 0) pmean /= (Double_t)count;
452    else pmean = 0.0;
453
454    return pmean;
455 }
456
457 //_____________________________________________________________________________
458 Bool_t AliRsnEvent::GetAngleDistr
459 (Double_t &angleMean, Double_t &angleRMS, AliRsnDaughter *ref)
460 {
461 //
462 // Takes the leading particle and computes the mean and RMS
463 // of the distribution of directions of all other tracks
464 // with respect to the direction of leading particle.
465 //
466
467    AliRsnDaughter leading;
468    if (ref) leading = *ref;
469    else if (fLeading >= 0) SetDaughter(leading, fLeading, AliRsnDaughter::kTrack);
470    if (!leading.IsOK()) return kFALSE;
471
472    Int_t i, count, nTracks = fRef->GetNumberOfTracks();
473    Double_t angle, angle2Mean = 0.0;
474
475    angleMean = angle2Mean = 0.0;
476
477    for (i = 0, count = 0; i < nTracks; i++) {
478       AliRsnDaughter trk = GetDaughter(i);
479       if (trk.GetID() == leading.GetID()) continue;
480
481       angle = leading.Prec().Angle(trk.Prec().Vect());
482
483       angleMean += angle;
484       angle2Mean += angle * angle;
485       count++;
486    }
487
488    if (!count) return kFALSE;
489
490    angleMean /= (Double_t)count;
491    angle2Mean /= (Double_t)count;
492    angleRMS = TMath::Sqrt(angle2Mean - angleMean * angleMean);
493
494    return kTRUE;
495 }
496
497 //_____________________________________________________________________________
498 Bool_t AliRsnEvent::SetDaughterESDtrack(AliRsnDaughter &out, Int_t i)
499 {
500 //
501 // Setup the first argument to the track identified by the index.
502 // When available, adds the MC information and references.
503 // ---
504 // Version #1: ESD tracks
505 //
506
507    // check 1: index in good range
508    if (i < 0 || i >= fRef->GetNumberOfTracks()) {
509       out.SetBad();
510       return kFALSE;
511    }
512
513    // check 2: not NULL object
514    AliVTrack *track = (AliVTrack*)fRef->GetTrack(i);
515    if (!track) {
516       out.SetBad();
517       return kFALSE;
518    }
519
520    // assign references of reconstructed track
521    Int_t label = TMath::Abs(track->GetLabel());
522    out.SetRef(track);
523    out.SetLabel(label);
524    out.SetGood();
525
526    // assign MC info, if available
527    return SetMCInfoESD(out);
528 }
529
530 //_____________________________________________________________________________
531 Bool_t AliRsnEvent::SetDaughterAODtrack(AliRsnDaughter &out, Int_t i)
532 {
533 //
534 // Setup the first argument to the track identified by the index.
535 // When available, adds the MC information and references.
536 // ---
537 // Version #2: AOD tracks
538 //
539
540    // check 1: index in good range
541    if (i < 0 || i >= fRef->GetNumberOfTracks()) {
542       out.SetBad();
543       return kFALSE;
544    }
545
546    // check 2: not NULL object
547    AliVTrack *track = (AliVTrack*)fRef->GetTrack(i);
548    if (!track) {
549       out.SetBad();
550       return kFALSE;
551    }
552
553    // assign references of reconstructed track
554    Int_t label = TMath::Abs(track->GetLabel());
555    out.SetRef(track);
556    out.SetLabel(label);
557    out.SetGood();
558
559    // assign MC info, if available
560    return SetMCInfoAOD(out);
561 }
562
563 //_____________________________________________________________________________
564 Bool_t AliRsnEvent::SetDaughterESDv0(AliRsnDaughter &out, Int_t i)
565 {
566 //
567 // Setup the first argument to the track identified by the index.
568 // When available, adds the MC information and references.
569 // ---
570 // Version #3: ESD v0
571 //
572
573    // check 1: index in good range
574    if (i > fRef->GetNumberOfV0s()) {
575       out.SetBad();
576       return 1;
577    }
578
579    // check 2: not NULL object
580    AliESDEvent *ev = GetRefESD();
581    AliESDv0    *v0 = ev->GetV0(i);
582    if (!v0) {
583       out.SetBad();
584       return 2;
585    }
586
587    // assign references of reconstructed track
588    out.SetRef(v0);
589    out.SetGood();
590    out.SetLabel(-1);
591
592    // this time, assigning label is not trivial,
593    // it is done only if MC is present and both
594    // daughters come from a true particle
595    AliMCEvent  *mc = GetRefMCESD();
596    AliESDtrack *tp = ev->GetTrack(v0->GetPindex());
597    AliESDtrack *tn = ev->GetTrack(v0->GetNindex());
598    if (mc && tp && tn) {
599       Int_t        lp = TMath::Abs(tp->GetLabel());
600       Int_t        ln = TMath::Abs(tn->GetLabel());
601       TParticle   *pp = mc->Stack()->Particle(lp);
602       TParticle   *pn = mc->Stack()->Particle(ln);
603       // if their first mothers are the same, the V0 is true
604       // otherwise no label can be assigned
605       if (pp->GetFirstMother() == pn->GetFirstMother() && pp->GetFirstMother() >= 0) out.SetLabel(pp->GetFirstMother());
606    }
607
608    // assign MC info, if available
609    return SetMCInfoESD(out);
610 }
611
612 //_____________________________________________________________________________
613 Bool_t AliRsnEvent::SetDaughterAODv0(AliRsnDaughter &out, Int_t i)
614 {
615 //
616 // Setup the first argument to the track identified by the index.
617 // When available, adds the MC information and references.
618 // ---
619 // Version #4: AOD v0
620 //
621
622    // check 1: index in good range
623    if (i > fRef->GetNumberOfV0s()) {
624       out.SetBad();
625       return kFALSE;
626    }
627
628    // check 2: not NULL object
629    AliAODEvent *ev = GetRefAOD();
630    AliAODv0    *v0 = ev->GetV0(i);
631    if (!v0) {
632       out.SetBad();
633       return kFALSE;
634    }
635
636    TClonesArray *mcArray = (TClonesArray*)ev->GetList()->FindObject(AliAODMCParticle::StdBranchName());
637    if (!mcArray) {
638       out.SetBad();
639       return kFALSE;
640    }
641
642    // assign references of reconstructed track
643    out.SetRef(v0);
644    out.SetGood();
645    out.SetLabel(-1);
646
647    // this time, assigning label is not trivial,
648    // it is done only if MC is present and both
649    // daughters come from a true particle
650    AliAODTrack  *tp = ev->GetTrack(v0->GetPosID());
651    AliAODTrack  *tn = ev->GetTrack(v0->GetNegID());
652    if (mcArray && tp && tn) {
653       Int_t        lp = TMath::Abs(tp->GetLabel());
654       Int_t        ln = TMath::Abs(tn->GetLabel());
655       // loop on array to find MC daughters
656       AliAODMCParticle *pp = 0x0, *pn = 0x0;
657       Int_t ipp = -1, ipn = -1;
658       TObjArrayIter next(mcArray);
659       AliAODMCParticle *part = 0x0;
660       while ((part = (AliAODMCParticle*)next())) {
661          if (TMath::Abs(part->GetLabel()) == lp) ipp = mcArray->IndexOf(part);
662          if (TMath::Abs(part->GetLabel()) == ln) ipn = mcArray->IndexOf(part);
663       }
664       if (ipp >= 0) pp = (AliAODMCParticle*)mcArray->At(ipp);
665       if (ipn >= 0) pn = (AliAODMCParticle*)mcArray->At(ipn);
666       if (!pp || !pn) return kTRUE;
667       // assign a MC reference and a label only to true V0s
668       if (pp->GetMother() == pn->GetMother() && pp->GetMother() >= 0) out.SetLabel(pp->GetMother());
669    }
670
671    // assign MC info, if available
672    return SetMCInfoAOD(out);
673 }
674
675 //_____________________________________________________________________________
676 Bool_t AliRsnEvent::SetDaughterESDcascade(AliRsnDaughter &, Int_t)
677 {
678 //
679 // Setup the first argument to the track identified by the index.
680 // When available, adds the MC information and references.
681 // ---
682 // Version #3: ESD cascade
683 //
684
685    return kTRUE;
686 }
687
688 //_____________________________________________________________________________
689 Bool_t AliRsnEvent::SetDaughterAODcascade(AliRsnDaughter &, Int_t)
690 {
691 //
692 // Setup the first argument to the track identified by the index.
693 // When available, adds the MC information and references.
694 // ---
695 // Version #4: AOD cascade
696 //
697
698    return kTRUE;
699 }
700
701 //_____________________________________________________________________________
702 Bool_t AliRsnEvent::SetMCInfoESD(AliRsnDaughter &out)
703 {
704 //
705 // Using the label assigned to the daughter, searches for the MC informations:
706 // -- MC reference
707 // -- mother
708 //
709
710    Int_t label = out.GetLabel();
711
712    // if no MC reference is available, exit here (successfully)
713    AliMCEvent *mc = GetRefMCESD();
714    if (!mc) return kTRUE;
715    Int_t nMC = mc->GetNumberOfTracks();
716    
717    // debug message for fakes
718    if (label < 0) {
719       AliDebug(AliLog::kDebug + 1, "Fake object (fake track or false V0)");
720       return kFALSE;
721    }
722
723    // assign MC reference, being aware of eventual
724    // overflows in the array (sometimes happened)
725    if (label >= nMC) {
726       AliWarning(Form("Stack overflow: track label = %d -- stack maximum = %d", label, nMC));
727       return kFALSE;
728    }
729    AliMCParticle *mcPart = (AliMCParticle*)mc->GetTrack(label);
730    if (!mcPart) {
731       AliWarning(Form("Stack discontinuity: label %d refers to a NULL object", label));
732       return kFALSE;
733    }
734    out.SetRefMC(mcPart);
735
736    // if this is a primary track, exit here (successfully)
737    Int_t imum = mcPart->Particle()->GetFirstMother();
738    if (imum < 0) {
739       out.SetMotherPDG(0);
740       return kTRUE;
741    }
742
743    // if didn't stop there, search for the mother
744    if (imum >= nMC) {
745       AliWarning(Form("Stack overflow: track mother label = %d -- stack maximum = %d", label, nMC));
746       return kFALSE;
747    }
748    AliMCParticle *mcMother = (AliMCParticle*)mc->GetTrack(imum);
749    if (!mcMother) {
750       AliWarning(Form("Stack discontinuity: label mother %d refers to a NULL object", imum));
751       return kFALSE;
752    }
753    out.SetMotherPDG(TMath::Abs(mcMother->Particle()->GetPdgCode()));
754
755    return kTRUE;
756 }
757
758 //_____________________________________________________________________________
759 Bool_t AliRsnEvent::SetMCInfoAOD(AliRsnDaughter &out)
760 {
761 //
762 // Using the label assigned to the daughter, searches for the MC informations:
763 // -- MC reference
764 // -- mother
765 //
766
767    Int_t label = out.GetLabel();
768    
769    // debug message for fakes
770    if (label < 0) {
771       AliDebug(AliLog::kDebug + 1, "Fake object (fake track or false V0)");
772       return kFALSE;
773    }
774
775    // if no MC reference is available, exit here (successfully)
776    AliAODEvent *mc = GetRefMCAOD();
777    if (!mc) return kTRUE;
778
779    // loop on particles using the track label as reference
780    // and then assign also the mother type, if found
781    TClonesArray *mcArray = (TClonesArray*)mc->GetList()->FindObject(AliAODMCParticle::StdBranchName());
782    if (!mcArray) return kFALSE;
783    TObjArrayIter next(mcArray);
784    AliAODMCParticle *part = 0x0;
785    while ((part = (AliAODMCParticle*)next())) {
786       if (TMath::Abs(part->GetLabel()) == label) {
787          out.SetRefMC(part);
788          out.SetMotherPDG(0);
789          Int_t imum = part->GetMother();
790          if (imum >= 0 && imum <= mcArray->GetEntriesFast()) {
791             AliAODMCParticle *mum = (AliAODMCParticle*)mcArray->At(imum);
792             if (mum) out.SetMotherPDG(TMath::Abs(mum->GetPdgCode()));
793             return kTRUE;
794          } else {
795             AliWarning(Form("Array overflow: track mother label = %d -- stack maximum = %d", imum, mcArray->GetEntriesFast()));
796             return kFALSE;
797          }
798
799          // if a MC reference is found, there is no need to go on with the loop
800          break;
801       }
802    }
803
804    return kTRUE;
805 }