]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnDaughter.cxx
Where possible, replaced dynamic_cast with ROOT RTTI using TClass.
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnDaughter.cxx
1 //
2 // Class AliRsnDaughter
3 //
4 // Interface to candidate daughters of a resonance (tracks).
5 // Points to the source of information, which is generally an AliVParticle-derived object
6 // and contains few internal data-members to store "on fly" some important information
7 // for the computations required during resonance analysis.
8 // It contains a TLorentzVector data-member which, provided that a meaningful mass was assigned,
9 // eases a lot the computation of invariant masses from summing up several of these objects.
10 //
11 // authors: A. Pulvirenti (alberto.pulvirenti@ct.infn.it)
12 //          M. Vala (martin.vala@cern.ch)
13 //
14
15 #include <TParticle.h>
16 #include "AliAODVertex.h"
17
18 #include "AliRsnDaughter.h"
19
20 ClassImp(AliRsnDaughter)
21
22 //_____________________________________________________________________________
23 AliRsnDaughter::AliRsnDaughter() :
24    fOK(kFALSE),
25    fLabel(-1),
26    fMotherPDG(0),
27    fRsnID(-1),
28    fPrec(0.0, 0.0, 0.0, 0.0),
29    fPsim(0.0, 0.0, 0.0, 0.0),
30    fRef(0x0),
31    fRefMC(0x0)
32 {
33 //
34 // Default constructor.
35 // Initializes all data members to the same values
36 // Which will be given by Reset() method.
37 //
38 }
39
40 //_____________________________________________________________________________
41 AliRsnDaughter::AliRsnDaughter(const AliRsnDaughter &copy) :
42    TObject(copy),
43    fOK(copy.fOK),
44    fLabel(copy.fLabel),
45    fMotherPDG(copy.fMotherPDG),
46    fRsnID(copy.fRsnID),
47    fPrec(copy.fPrec),
48    fPsim(copy.fPsim),
49    fRef(copy.fRef),
50    fRefMC(copy.fRefMC)
51 {
52 //
53 // Copy constructor.
54 // Pointers are NOT duplicated, since they don't come from a 'new'
55 // statement, but from just referencing something in the data source.
56 //
57 }
58
59 //_____________________________________________________________________________
60 AliRsnDaughter& AliRsnDaughter::operator=(const AliRsnDaughter &copy)
61 {
62 //
63 // Assignment operator.
64 // Pointers are NOT duplicated, since they don't come from a 'new'
65 // statement, but from just referencing something in the data source.
66 //
67
68    fOK        = copy.fOK;
69    fLabel     = copy.fLabel;
70    fMotherPDG = copy.fMotherPDG;
71    fRsnID     = copy.fRsnID;
72    fPrec      = copy.fPrec;
73    fPsim      = copy.fPsim;
74    fRef       = copy.fRef;
75    fRefMC     = copy.fRefMC;
76
77    return (*this);
78 }
79
80 //_____________________________________________________________________________
81 void AliRsnDaughter::SetRef(AliVParticle *p)
82 {
83 //
84 // Set the pointer to reference VParticle
85 // and copies its momentum in the 4-vector.
86 // Mass is assigned by SetMass().
87 //
88
89    fRef = p;
90    if (p) fPrec.SetXYZT(p->Px(), p->Py(), p->Pz(), 0.0);
91    else   fPrec.SetXYZT(0.0, 0.0, 0.0, 0.0);
92 }
93
94 //_____________________________________________________________________________
95 void AliRsnDaughter::SetRefMC(AliVParticle *p)
96 {
97 //
98 // Set the pointer to reference MonteCarlo VParticle
99 // and copies its momentum in the 4-vector.
100 // Mass is assigned by SetMass().
101 //
102
103    fRefMC = p;
104    if (p) fPsim.SetXYZT(p->Px(), p->Py(), p->Pz(), 0.0);
105    else   fPrec.SetXYZT(0.0, 0.0, 0.0, 0.0);
106 }
107
108 //_____________________________________________________________________________
109 void AliRsnDaughter::Reset()
110 {
111 //
112 // Reset this track to meaningless values and to a 'bad' status.
113 // After this has been done, this object should not be used
114 // for analysis unless initialized properly.
115 //
116
117    fOK        = kFALSE;
118    fLabel     = -1;
119    fMotherPDG =  0;
120    fRsnID     = -1;
121
122    SetRef(NULL);
123    SetRefMC(NULL);
124 }
125
126
127 //_____________________________________________________________________________
128 Int_t AliRsnDaughter::GetPDG(Bool_t abs)
129 {
130 //
131 // Return the PDG code of the particle from MC ref (if any).
132 // If argument is kTRUE, returns its absolute value.
133 //
134
135    Int_t pdg = 0;
136
137    // ESD
138    AliMCParticle *esd = GetRefMCESD();
139    if (esd) pdg = esd->Particle()->GetPdgCode();
140
141    // AOD
142    AliAODMCParticle *aod = GetRefMCAOD();
143    if (aod) pdg = aod->GetPdgCode();
144
145    // abs value if required
146    if (abs) pdg = TMath::Abs(pdg);
147    return pdg;
148 }
149
150 //_____________________________________________________________________________
151 Int_t AliRsnDaughter::GetID()
152 {
153 //
154 // Return reference index, using the "GetID" method
155 // of the possible source object.
156 // When this method is unsuccessful (e.g.: V0s), return the label.
157 //
158
159    // ESD tracks
160    AliESDtrack *esd = GetRefESDtrack();
161    if (esd) return esd->GetID();
162
163    // AOD tracks
164    AliAODTrack *aod = GetRefAODtrack();
165    if (aod) return aod->GetID();
166
167    // whatever else
168    return GetLabel();
169 }
170
171 //_____________________________________________________________________________
172 Bool_t AliRsnDaughter::HasFlag(ULong_t flag)
173 {
174 //
175 // Checks that the 'status' flag of the source object has one or
176 // a combination of status flags specified in argument.
177 // Works only with objects inheriting from AliVTrack base class.
178 //
179
180    AliVTrack *track  = GetRefVtrack();
181    if (!track) return kFALSE;
182
183    ULong_t status = (ULong_t)track->GetStatus();
184
185    return ((status & flag) != 0);
186 }
187
188 //_____________________________________________________________________________
189 Bool_t AliRsnDaughter::SetMass(Double_t mass)
190 {
191 //
192 // Assign a mass hypothesis to the track.
193 // This causes the 4-momentum data members to be re-initialized
194 // using the momenta of referenced tracks/v0s and this mass.
195 // This step is fundamental for the following of the analysis.
196 //
197
198    if (mass < 0.) return kFALSE;
199
200    if (fRef)   fPrec.SetXYZM(fRef  ->Px(), fRef  ->Py(), fRef  ->Pz(), mass);
201    if (fRefMC) fPsim.SetXYZM(fRefMC->Px(), fRefMC->Py(), fRefMC->Pz(), mass);
202
203    return kTRUE;
204 }
205
206 //_____________________________________________________________________________
207 Bool_t AliRsnDaughter::IsKinkDaughter()
208 {
209 //
210 // Checks if this track is a kink daughter.
211 // this information is important for some cuts, in some cases
212 // and it is retrieved differently from ESDs and AODs, so
213 // this is done here in order to have a unique outcome.
214 //
215
216    AliESDtrack *etrack = GetRefESDtrack();
217    AliAODTrack *atrack = GetRefAODtrack();
218
219    if (etrack) {
220       return (etrack->GetKinkIndex(0) > 0);
221    } else if (atrack) {
222       AliAODVertex *vertex = atrack->GetProdVertex();
223       if (vertex) if (vertex->GetType() == AliAODVertex::kKink) return kTRUE;
224    }
225
226    return kFALSE;
227 }