]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnDaughter.cxx
aa80cd4e6c2d9000e769dfda454995984e298086
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnDaughter.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 //  This class works as generic interface to each candidate resonance daughter.
18 //  Its main purpose is to provide a unique reference which includes all the
19 //  facilities available in the AliVParticle generic base class, plus all info
20 //  which could be needed during analysis, which are not in AliVParticle but
21 //  need to be accessed from ESD or AOD objects, usually in different ways.
22 //  When MC is available, AliRsnDaughter matches each reconstructed object with
23 //  its corresponding MC particle.
24 //  
25 //  Currently, this interface can point to all kinds of single-particle object
26 //  which one can have in the reconstructed event: charged tracks, V0s and 
27 //  cascades. It is care of the user to treat each of them in the correct way,
28 //  regarding cuts, functions to be computed, etc.
29 //
30 //  authors: A. Pulvirenti (alberto.pulvirenti@ct.infn.it)
31 //           M. Vala (martin.vala@cern.ch)
32 //
33
34 #include <TParticle.h>
35 #include <TDatabasePDG.h>
36 #include "AliAODVertex.h"
37
38 #include "AliRsnDaughter.h"
39
40 ClassImp(AliRsnDaughter)
41
42 //_____________________________________________________________________________
43 AliRsnDaughter::AliRsnDaughter() :
44    fOK(kFALSE),
45    fLabel(-1),
46    fMotherPDG(0),
47    fRsnID(-1),
48    fPrec(0.0, 0.0, 0.0, 0.0),
49    fPsim(0.0, 0.0, 0.0, 0.0),
50    fRef(0x0),
51    fRefMC(0x0),
52    fOwnerEvent(0x0)
53 {
54 //
55 // Default constructor.
56 // Initializes all data members to the same values
57 // Which will be given by Reset() method.
58 //
59 }
60
61 //_____________________________________________________________________________
62 AliRsnDaughter::AliRsnDaughter(const AliRsnDaughter &copy) :
63    TObject(copy),
64    fOK(copy.fOK),
65    fLabel(copy.fLabel),
66    fMotherPDG(copy.fMotherPDG),
67    fRsnID(copy.fRsnID),
68    fPrec(copy.fPrec),
69    fPsim(copy.fPsim),
70    fRef(copy.fRef),
71    fRefMC(copy.fRefMC),
72    fOwnerEvent(copy.fOwnerEvent)
73 {
74 //
75 // Copy constructor.
76 // Pointers are NOT duplicated, since they don't come from a 'new'
77 // statement, but from just referencing something in the data source.
78 //
79 }
80
81 //_____________________________________________________________________________
82 AliRsnDaughter& AliRsnDaughter::operator=(const AliRsnDaughter &copy)
83 {
84 //
85 // Assignment operator.
86 // Pointers are NOT duplicated, since they don't come from a 'new'
87 // statement, but from just referencing something in the data source.
88 //
89
90    fOK         = copy.fOK;
91    fLabel      = copy.fLabel;
92    fMotherPDG  = copy.fMotherPDG;
93    fRsnID      = copy.fRsnID;
94    fPrec       = copy.fPrec;
95    fPsim       = copy.fPsim;
96    fRef        = copy.fRef;
97    fRefMC      = copy.fRefMC;
98    fOwnerEvent = copy.fOwnerEvent;
99
100    return (*this);
101 }
102
103 //_____________________________________________________________________________
104 void AliRsnDaughter::SetRef(AliVParticle *p)
105 {
106 //
107 // Set the pointer to reference reconstructed VParticle
108 // and copies its momentum in the 4-vector.
109 // Mass is set to 0 (must be assigned by SetMass()).
110 //
111
112    fPrec.SetXYZT(0.0, 0.0, 0.0, 0.0);
113    fRef = p;
114    
115    if (fRef) {
116       fPrec.SetX(fRef->Px());
117       fPrec.SetY(fRef->Py());
118       fPrec.SetZ(fRef->Pz());
119    }
120 }
121
122 //_____________________________________________________________________________
123 void AliRsnDaughter::SetRefMC(AliVParticle *p)
124 {
125 //
126 // Set the pointer to reference MonteCarlo VParticle
127 // and copies its momentum in the 4-vector.
128 // Mass is set to 0 (must be assigned by SetMass()).
129 //
130
131    fPsim.SetXYZT(0.0, 0.0, 0.0, 0.0);
132    fRefMC = p;
133    
134    if (fRefMC) {
135       fPsim.SetX(fRefMC->Px());
136       fPsim.SetY(fRefMC->Py());
137       fPsim.SetZ(fRefMC->Pz());
138    }
139 }
140
141 //_____________________________________________________________________________
142 void AliRsnDaughter::Reset()
143 {
144 //
145 // Reset this track to meaningless values and to a 'bad' status.
146 // After this has been done, this object should not be used
147 // for analysis unless initialized properly.
148 //
149
150    fOK        = kFALSE;
151    fLabel     = -1;
152    fMotherPDG =  0;
153    fRsnID     = -1;
154
155    SetRef(NULL);
156    SetRefMC(NULL);
157 }
158
159 //_____________________________________________________________________________
160 Int_t AliRsnDaughter::GetPDG(Bool_t abs)
161 {
162 //
163 // Return the PDG code of the particle from MC ref (if any).
164 // If argument is kTRUE, returns its absolute value.
165 //
166
167    Int_t pdg = 0;
168    if (!fRefMC) return pdg;
169
170    // ESD
171    AliMCParticle *esd = GetRefMCESD();
172    if (esd) pdg = esd->Particle()->GetPdgCode();
173
174    // AOD
175    AliAODMCParticle *aod = GetRefMCAOD();
176    if (aod) pdg = aod->GetPdgCode();
177
178    // abs value if required
179    if (abs) pdg = TMath::Abs(pdg);
180    return pdg;
181 }
182
183 //_____________________________________________________________________________
184 Int_t AliRsnDaughter::GetID()
185 {
186 //
187 // Return reference index, using the "GetID" method
188 // of the possible source object.
189 // When this method is unsuccessful (e.g.: V0s), return the label.
190 //
191
192    // ESD tracks
193    AliESDtrack *esd = GetRefESDtrack();
194    if (esd) return esd->GetID();
195
196    // AOD tracks
197    AliAODTrack *aod = GetRefAODtrack();
198    if (aod) return aod->GetID();
199
200    // whatever else
201    return GetLabel();
202 }
203
204 //_____________________________________________________________________________
205 Bool_t AliRsnDaughter::IsKinkDaughter()
206 {
207 //
208 // Checks if this track is a kink daughter.
209 // this information is important for some cuts, in some cases
210 // and it is retrieved differently from ESDs and AODs, so
211 // this is done here in order to have a unique outcome.
212 //
213
214    AliESDtrack *etrack = GetRefESDtrack();
215    AliAODTrack *atrack = GetRefAODtrack();
216
217    if (etrack) {
218       return (etrack->GetKinkIndex(0) > 0);
219    } else if (atrack) {
220       AliAODVertex *vertex = atrack->GetProdVertex();
221       if (vertex) if (vertex->GetType() == AliAODVertex::kKink) return kTRUE;
222    }
223
224    return kFALSE;
225 }
226
227 //______________________________________________________________________________
228 AliRsnDaughter::ERefType AliRsnDaughter::RefType(ESpecies species)
229 {
230 //
231 // Returns the expected object type for a candidate daughter
232 // of the given species.
233 //
234
235    switch (species) {
236       case kElectron:
237       case kMuon:
238       case kPion:
239       case kKaon:
240       case kProton:
241          return kTrack;
242       case kKaon0:
243       case kLambda:
244          return kV0;
245       case kXi:
246       case kOmega:
247          return kCascade;
248       default:
249          return kNoType;
250    }
251 }
252
253 //______________________________________________________________________________
254 void AliRsnDaughter::Print(Option_t *) const
255 {
256 //
257 // Override of TObject::Print()
258 //
259
260    AliInfo("=== DAUGHTER INFO ======================================================================");
261    //if (fRef) {
262    //   AliInfo(Form(" Ref  : %x (%15s) with px,py,pz = %6.2f %6.2f %6.2f", (UInt_t)fRef  , fRef  ->ClassName(), fPrec.X(), fPrec.Y(), fPrec.Z()));
263    //} else {
264    //   AliInfo(" Ref  : NULL");
265    //}
266    //if (fRefMC) {
267    //   AliInfo(Form(" RefMC: %x (%15s) with px,py,pz = %6.2f %6.2f %6.2f", (UInt_t)fRefMC, fRefMC->ClassName(), fPsim.X(), fPsim.Y(), fPsim.Z()));
268    //} else {
269    //   AliInfo(" RefMC: NULL");
270    //}
271    AliInfo(Form(" OK, RsnID, Label, MotherPDG = %s, %5d, %5d, %4d", (fOK ? "true " : "false"), fRsnID, fLabel, fMotherPDG));
272    AliInfo("========================================================================================");
273 }
274
275 //______________________________________________________________________________
276 const char* AliRsnDaughter::SpeciesName(ESpecies species)
277 {
278 //
279 // Return a string with the short name of the particle
280 //
281
282    switch (species) {
283       case kElectron: return "E";
284       case kMuon:     return "Mu";
285       case kPion:     return "Pi";
286       case kKaon:     return "K";
287       case kProton:   return "P";
288       case kKaon0:    return "K0s";
289       case kLambda:   return "Lambda";
290       case kXi:       return "Xi";
291       case kOmega:    return "Omega";
292       default:        return "Undef";
293    }
294 }
295
296 //______________________________________________________________________________
297 Int_t AliRsnDaughter::SpeciesPDG(ESpecies species)
298 {
299 //
300 // Return the PDG code of a particle species (abs value)
301 //
302
303    switch (species) {
304       case kElectron: return 11;
305       case kMuon:     return 13;
306       case kPion:     return 211;
307       case kKaon:     return 321;
308       case kProton:   return 2212;
309       case kKaon0:    return 310;
310       case kLambda:   return 3122;
311       case kXi:       return 3312;
312       case kOmega:    return 3334;
313       default:        return 0;
314    }
315 }
316
317 //______________________________________________________________________________
318 Double_t AliRsnDaughter::SpeciesMass(ESpecies species)
319 {
320 //
321 // Return the mass of a particle species
322 //
323
324    TDatabasePDG *db = TDatabasePDG::Instance();
325    TParticlePDG *part = 0x0;
326    
327    Int_t pdg = SpeciesPDG(species);
328    if (pdg) {
329       part = db->GetParticle(pdg);
330       return part->Mass();
331    }
332    else
333       return 0.0;
334 }
335
336 //______________________________________________________________________________
337 EPARTYPE AliRsnDaughter::ToAliPID(ESpecies species)
338 {
339 //
340 // Convert an enum element from this object
341 // into the enumeration of AliPID.
342 // If no match are cound 'kUnknown' is returned.
343 //
344
345    switch (species) {
346       case kElectron: return AliPID::kElectron;
347       case kMuon:     return AliPID::kMuon;
348       case kPion:     return AliPID::kPion;
349       case kKaon:     return AliPID::kKaon;
350       case kProton:   return AliPID::kProton;
351       case kKaon0:    return AliPID::kKaon0;
352       default:        return AliPID::kUnknown;
353    }
354 }
355
356 //______________________________________________________________________________
357 AliRsnDaughter::ESpecies AliRsnDaughter::FromAliPID(EPARTYPE pid)
358 {
359 //
360 // Convert an enum element from AliPID
361 // into the enumeration of this object.
362 // If no match are cound 'kUnknown' is returned.
363 //
364
365    switch (pid) {
366       case AliPID::kElectron: return kElectron;
367       case AliPID::kMuon:     return kMuon;
368       case AliPID::kPion:     return kPion;
369       case AliPID::kKaon:     return kKaon;
370       case AliPID::kProton:   return kProton;
371       case AliPID::kKaon0:    return kKaon0;
372       default:                return kUnknown;
373    }
374 }