]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/AliRsnDaughter.cxx
fixing psi in MC from the header
[u/mrichter/AliRoot.git] / PWGLF / 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    if (this == &copy)
55       return *this;
56
57    fOK         = copy.fOK;
58    fLabel      = copy.fLabel;
59    fMotherPDG  = copy.fMotherPDG;
60    fRsnID      = copy.fRsnID;
61    fPsim       = copy.fPsim;
62    fPrec       = copy.fPrec;
63    fRef        = copy.fRef;
64    fRefMC      = copy.fRefMC;
65    fOwnerEvent = copy.fOwnerEvent;
66
67    return (*this);
68 }
69
70 //_____________________________________________________________________________
71 void AliRsnDaughter::Reset()
72 {
73 //
74 // Reset this track to meaningless values and to a 'bad' status.
75 // After this has been done, this object should not be used
76 // for analysis unless initialized properly.
77 //
78
79    fOK        = kFALSE;
80    fLabel     = -1;
81    fMotherPDG =  0;
82    fRsnID     = -1;
83
84    fPsim.SetXYZT(0.0, 0.0, 0.0, 0.0);
85    fPrec.SetXYZT(0.0, 0.0, 0.0, 0.0);
86
87    fRef = fRefMC = 0x0;
88    fOwnerEvent = 0x0;
89 }
90
91 //_____________________________________________________________________________
92 Int_t AliRsnDaughter::GetPDG()
93 {
94 //
95 // Return the PDG code of the particle from MC ref (if any).
96 // If argument is kTRUE, returns its absolute value.
97 //
98
99    if (Match(fRefMC, AliMCParticle::Class()))
100       return ((AliMCParticle *)fRefMC)->Particle()->GetPdgCode();
101    else if (Match(fRefMC, AliAODMCParticle::Class()))
102       return ((AliAODMCParticle *)fRefMC)->GetPdgCode();
103    else {
104       AliWarning("Cannot retrieve PDG");
105       return 0;
106    }
107 }
108
109 //_____________________________________________________________________________
110 Int_t AliRsnDaughter::GetID()
111 {
112 //
113 // Return reference index, using the "GetID" method
114 // of the possible source object.
115 // When this method is unsuccessful (e.g.: V0s), return the label.
116 //
117
118    // ESD tracks
119    AliESDtrack *esd = Ref2ESDtrack();
120    if (esd) return esd->GetID();
121
122    // AOD tracks
123    AliAODTrack *aod = Ref2AODtrack();
124    if (aod) return aod->GetID();
125
126    // whatever else
127    return GetLabel();
128 }
129
130 //_____________________________________________________________________________
131 Int_t AliRsnDaughter::GetMother()
132 {
133 //
134 // Return index of the first mother of the MC reference, if any.
135 // Otherwise, returns -1 (the same as for primary tracks)
136 //
137
138    if (!fRefMC) return -1;
139
140    if (fRefMC->InheritsFrom(AliMCParticle::Class())) {
141       AliMCParticle *mc = (AliMCParticle *)fRefMC;
142       return mc->Particle()->GetFirstMother();
143    } else if (fRefMC->InheritsFrom(AliAODMCParticle::Class())) {
144       AliAODMCParticle *mc = (AliAODMCParticle *)fRefMC;
145       return mc->GetMother();
146    }
147    else
148       return -1;
149 }
150
151
152
153 //______________________________________________________________________________
154 void AliRsnDaughter::Print(Option_t *) const
155 {
156 //
157 // Override of TObject::Print()
158 //
159
160    AliInfo("=== DAUGHTER INFO ======================================================================");
161    AliInfo(Form(" (sim) px,py,pz = %6.2f %6.2f %6.2f", fPsim.X(), fPsim.Y(), fPsim.Z()));
162    AliInfo(Form(" (rec) px,py,pz = %6.2f %6.2f %6.2f", fPrec.X(), fPrec.Y(), fPrec.Z()));
163    AliInfo(Form(" OK, RsnID, Label, MotherPDG = %s, %5d, %5d, %4d", (fOK ? "true " : "false"), fRsnID, fLabel, fMotherPDG));
164    AliInfo("========================================================================================");
165 }
166
167 //______________________________________________________________________________
168 const char *AliRsnDaughter::SpeciesName(ESpecies species)
169 {
170 //
171 // Return a string with the short name of the particle
172 //
173
174    switch (species) {
175       case kElectron: return "E";
176       case kMuon:     return "Mu";
177       case kPion:     return "Pi";
178       case kKaon:     return "K";
179       case kProton:   return "P";
180       case kKaon0:    return "K0s";
181       case kLambda:   return "Lambda";
182       case kXi:       return "Xi";
183       case kOmega:    return "Omega";
184       default:        return "Undef";
185    }
186 }
187
188 //______________________________________________________________________________
189 Int_t AliRsnDaughter::SpeciesPDG(ESpecies species)
190 {
191 //
192 // Return the PDG code of a particle species (abs value)
193 //
194
195    switch (species) {
196       case kElectron: return 11;
197       case kMuon:     return 13;
198       case kPion:     return 211;
199       case kKaon:     return 321;
200       case kProton:   return 2212;
201       case kKaon0:    return 310;
202       case kLambda:   return 3122;
203       case kXi:       return 3312;
204       case kOmega:    return 3334;
205       default:        return 0;
206    }
207 }
208
209 //______________________________________________________________________________
210 Double_t AliRsnDaughter::SpeciesMass(ESpecies species)
211 {
212 //
213 // Return the mass of a particle species
214 //
215
216    TDatabasePDG *db = TDatabasePDG::Instance();
217    TParticlePDG *part = 0x0;
218
219    Int_t pdg = SpeciesPDG(species);
220    if (pdg) {
221       part = db->GetParticle(pdg);
222       return part->Mass();
223    }
224    else
225       return 0.0;
226 }
227
228 //______________________________________________________________________________
229 EPARTYPE AliRsnDaughter::ToAliPID(ESpecies species)
230 {
231 //
232 // Convert an enum element from this object
233 // into the enumeration of AliPID.
234 // If no match are cound 'kUnknown' is returned.
235 //
236
237    switch (species) {
238       case kElectron: return AliPID::kElectron;
239       case kMuon:     return AliPID::kMuon;
240       case kPion:     return AliPID::kPion;
241       case kKaon:     return AliPID::kKaon;
242       case kProton:   return AliPID::kProton;
243       case kKaon0:    return AliPID::kKaon0;
244       default:        return AliPID::kUnknown;
245    }
246 }
247
248 //______________________________________________________________________________
249 AliRsnDaughter::ESpecies AliRsnDaughter::FromAliPID(EPARTYPE pid)
250 {
251 //
252 // Convert an enum element from AliPID
253 // into the enumeration of this object.
254 // If no match are cound 'kUnknown' is returned.
255 //
256
257    switch (pid) {
258       case AliPID::kElectron: return kElectron;
259       case AliPID::kMuon:     return kMuon;
260       case AliPID::kPion:     return kPion;
261       case AliPID::kKaon:     return kKaon;
262       case AliPID::kProton:   return kProton;
263       case AliPID::kKaon0:    return kKaon0;
264       default:                return kUnknown;
265    }
266 }