]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnDaughter.cxx
e17b09296d80b6460d832c573626be8eac196f34
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnDaughter.cxx
1 //
2 //  This class works as generic interface to each candidate resonance daughter.
3 //  Its main purpose is to provide a unique reference which includes all the
4 //  facilities available in the AliVParticle generic base class, plus all info
5 //  which could be needed during analysis, which are not in AliVParticle but
6 //  need to be accessed from ESD or AOD objects, usually in different ways.
7 //  When MC is available, AliRsnDaughter matches each reconstructed object with
8 //  its corresponding MC particle.
9 //  
10 //  Currently, this interface can point to all kinds of single-particle object
11 //  which one can have in the reconstructed event: charged tracks, V0s and 
12 //  cascades. It is care of the user to treat each of them in the correct way,
13 //  regarding cuts, functions to be computed, etc.
14 //
15 //  authors: A. Pulvirenti (alberto.pulvirenti@ct.infn.it)
16 //           M. Vala (martin.vala@cern.ch)
17 //
18
19 #include <TParticle.h>
20 #include <TDatabasePDG.h>
21
22 #include "AliRsnDaughter.h"
23
24 ClassImp(AliRsnDaughter)
25
26 //_____________________________________________________________________________
27 AliRsnDaughter::AliRsnDaughter(const AliRsnDaughter &copy) :
28    TObject(copy),
29    fOK(copy.fOK),
30    fLabel(copy.fLabel),
31    fMotherPDG(copy.fMotherPDG),
32    fRsnID(copy.fRsnID),
33    fPrec(copy.fPrec),
34    fPsim(copy.fPsim),
35    fRef(copy.fRef),
36    fRefMC(copy.fRefMC),
37    fOwnerEvent(copy.fOwnerEvent)
38 {
39 //
40 // Copy constructor.
41 // Pointers are NOT duplicated, since they don't come from a 'new'
42 // statement, but from just referencing something in the data source.
43 //
44 }
45
46 //_____________________________________________________________________________
47 AliRsnDaughter& AliRsnDaughter::operator=(const AliRsnDaughter &copy)
48 {
49 //
50 // Assignment operator.
51 // Pointers are NOT duplicated, since they don't come from a 'new'
52 // statement, but from just referencing something in the data source.
53 //
54
55    fOK         = copy.fOK;
56    fLabel      = copy.fLabel;
57    fMotherPDG  = copy.fMotherPDG;
58    fRsnID      = copy.fRsnID;
59    fPsim       = copy.fPsim;
60    fPrec       = copy.fPrec;
61    fRef        = copy.fRef;
62    fRefMC      = copy.fRefMC;
63    fOwnerEvent = copy.fOwnerEvent;
64
65    return (*this);
66 }
67
68 //_____________________________________________________________________________
69 void AliRsnDaughter::Reset()
70 {
71 //
72 // Reset this track to meaningless values and to a 'bad' status.
73 // After this has been done, this object should not be used
74 // for analysis unless initialized properly.
75 //
76
77    fOK        = kFALSE;
78    fLabel     = -1;
79    fMotherPDG =  0;
80    fRsnID     = -1;
81    
82    fPsim.SetXYZT(0.0, 0.0, 0.0, 0.0);
83    fPrec.SetXYZT(0.0, 0.0, 0.0, 0.0);
84
85    fRef = fRefMC = 0x0;
86    fOwnerEvent = 0x0;
87 }
88
89 //_____________________________________________________________________________
90 Int_t AliRsnDaughter::GetPDG()
91 {
92 //
93 // Return the PDG code of the particle from MC ref (if any).
94 // If argument is kTRUE, returns its absolute value.
95 //
96
97    if (Match(fRefMC, AliMCParticle::Class()))
98       return ((AliMCParticle*)fRefMC)->Particle()->GetPdgCode();
99    else if (Match(fRefMC, AliAODMCParticle::Class()))
100       return ((AliAODMCParticle*)fRefMC)->GetPdgCode();
101    else {
102       AliWarning("Cannot retrieve PDG");
103       return 0;
104    }
105 }
106
107 //_____________________________________________________________________________
108 Int_t AliRsnDaughter::GetID()
109 {
110 //
111 // Return reference index, using the "GetID" method
112 // of the possible source object.
113 // When this method is unsuccessful (e.g.: V0s), return the label.
114 //
115
116    // ESD tracks
117    AliESDtrack *esd = Ref2ESDtrack();
118    if (esd) return esd->GetID();
119
120    // AOD tracks
121    AliAODTrack *aod = Ref2AODtrack();
122    if (aod) return aod->GetID();
123
124    // whatever else
125    return GetLabel();
126 }
127
128 //_____________________________________________________________________________
129 Int_t AliRsnDaughter::GetMother()
130 {
131 //
132 // Return index of the first mother of the MC reference, if any.
133 // Otherwise, returns -1 (the same as for primary tracks)
134 //
135
136    if (!fRefMC) return -1;
137
138    if (fRefMC->InheritsFrom(AliMCParticle::Class())) {
139       AliMCParticle *mc = (AliMCParticle*)fRefMC;
140       return mc->Particle()->GetFirstMother();
141    } else if (fRefMC->InheritsFrom(AliAODMCParticle::Class())) {
142       AliAODMCParticle *mc = (AliAODMCParticle*)fRefMC;
143       return mc->GetMother();
144    }
145    else
146       return -1;
147 }
148    
149    
150
151 //______________________________________________________________________________
152 void AliRsnDaughter::Print(Option_t *opt) const
153 {
154 //
155 // Override of TObject::Print()
156 //
157
158    AliInfo("=== DAUGHTER INFO ======================================================================");
159    AliInfo(Form(" (sim) px,py,pz = %6.2f %6.2f %6.2f", fPsim.X(), fPsim.Y(), fPsim.Z()));
160    AliInfo(Form(" (rec) px,py,pz = %6.2f %6.2f %6.2f", fPrec.X(), fPrec.Y(), fPrec.Z()));
161    AliInfo(Form(" OK, RsnID, Label, MotherPDG = %s, %5d, %5d, %4d", (fOK ? "true " : "false"), fRsnID, fLabel, fMotherPDG));
162    AliInfo("========================================================================================");
163 }
164
165 //______________________________________________________________________________
166 const char* AliRsnDaughter::SpeciesName(ESpecies species)
167 {
168 //
169 // Return a string with the short name of the particle
170 //
171
172    switch (species) {
173       case kElectron: return "E";
174       case kMuon:     return "Mu";
175       case kPion:     return "Pi";
176       case kKaon:     return "K";
177       case kProton:   return "P";
178       case kKaon0:    return "K0s";
179       case kLambda:   return "Lambda";
180       case kXi:       return "Xi";
181       case kOmega:    return "Omega";
182       default:        return "Undef";
183    }
184 }
185
186 //______________________________________________________________________________
187 Int_t AliRsnDaughter::SpeciesPDG(ESpecies species)
188 {
189 //
190 // Return the PDG code of a particle species (abs value)
191 //
192
193    switch (species) {
194       case kElectron: return 11;
195       case kMuon:     return 13;
196       case kPion:     return 211;
197       case kKaon:     return 321;
198       case kProton:   return 2212;
199       case kKaon0:    return 310;
200       case kLambda:   return 3122;
201       case kXi:       return 3312;
202       case kOmega:    return 3334;
203       default:        return 0;
204    }
205 }
206
207 //______________________________________________________________________________
208 Double_t AliRsnDaughter::SpeciesMass(ESpecies species)
209 {
210 //
211 // Return the mass of a particle species
212 //
213
214    TDatabasePDG *db = TDatabasePDG::Instance();
215    TParticlePDG *part = 0x0;
216    
217    Int_t pdg = SpeciesPDG(species);
218    if (pdg) {
219       part = db->GetParticle(pdg);
220       return part->Mass();
221    }
222    else
223       return 0.0;
224 }
225
226 //______________________________________________________________________________
227 EPARTYPE AliRsnDaughter::ToAliPID(ESpecies species)
228 {
229 //
230 // Convert an enum element from this object
231 // into the enumeration of AliPID.
232 // If no match are cound 'kUnknown' is returned.
233 //
234
235    switch (species) {
236       case kElectron: return AliPID::kElectron;
237       case kMuon:     return AliPID::kMuon;
238       case kPion:     return AliPID::kPion;
239       case kKaon:     return AliPID::kKaon;
240       case kProton:   return AliPID::kProton;
241       case kKaon0:    return AliPID::kKaon0;
242       default:        return AliPID::kUnknown;
243    }
244 }
245
246 //______________________________________________________________________________
247 AliRsnDaughter::ESpecies AliRsnDaughter::FromAliPID(EPARTYPE pid)
248 {
249 //
250 // Convert an enum element from AliPID
251 // into the enumeration of this object.
252 // If no match are cound 'kUnknown' is returned.
253 //
254
255    switch (pid) {
256       case AliPID::kElectron: return kElectron;
257       case AliPID::kMuon:     return kMuon;
258       case AliPID::kPion:     return kPion;
259       case AliPID::kKaon:     return kKaon;
260       case AliPID::kProton:   return kProton;
261       case AliPID::kKaon0:    return kKaon0;
262       default:                return kUnknown;
263    }
264 }