]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/AliRsnDaughter.h
update on v2 calculation
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / AliRsnDaughter.h
1 #ifndef ALIRSNDAUGHTER_H
2 #define ALIRSNDAUGHTER_H
3
4 //
5 // Interface to single daughter candidate.
6 // Used for old-version analysis and for
7 // selecting tracks in mini package
8 //
9
10 #include <TMath.h>
11 #include <TLorentzVector.h>
12
13 #include "AliPID.h"
14 #include "AliVTrack.h"
15 #include "AliESDtrack.h"
16 #include "AliESDv0.h"
17 #include "AliESDcascade.h"
18 #include "AliMCParticle.h"
19 #include "AliAODTrack.h"
20 #include "AliAODv0.h"
21 #include "AliAODcascade.h"
22 #include "AliAODMCParticle.h"
23
24 class AliRsnEvent;
25
26 typedef AliPID::EParticleType EPARTYPE;
27
28 class AliRsnDaughter : public TObject {
29 public:
30
31    enum ERefType {
32       kTrack,
33       kV0,
34       kCascade,
35       kNoType
36    };
37
38    enum ESpecies {
39       kElectron,
40       kMuon,
41       kPion,
42       kKaon,
43       kProton,
44       kKaon0,
45       kLambda,
46       kXi,
47       kOmega,
48       kUnknown
49    };
50
51    AliRsnDaughter() : fOK(kFALSE), fLabel(-1), fMotherPDG(0), fRsnID(-1),
52       fPrec(), fPsim(), fRef(0x0), fRefMC(0x0), fOwnerEvent(0x0) { }
53    AliRsnDaughter(const AliRsnDaughter &copy);
54    AliRsnDaughter &operator= (const AliRsnDaughter &copy);
55    virtual ~AliRsnDaughter() { /*empty, since pointers must not be deleted*/ }
56
57    // basic getters
58    Bool_t          IsOK()          const {return fOK;}
59    Int_t           GetLabel()      const {return fLabel;}
60    Int_t           GetMotherPDG()  const {return fMotherPDG;}
61    Int_t           GetRsnID()      const {return fRsnID;}
62    AliVParticle   *GetRef()              {return fRef;}
63    AliVParticle   *GetRefMC()            {return fRefMC;}
64    AliRsnEvent    *GetOwnerEvent()       {return fOwnerEvent;}
65    TLorentzVector &P(Bool_t mc)          {return (mc ? fPsim : fPrec);}
66    TLorentzVector &Prec()                {return fPrec;}
67    TLorentzVector &Psim()                {return fPsim;}
68    Int_t           GetPDG();
69    Int_t           GetPDGAbs()           {return TMath::Abs(GetPDG());}
70    Int_t           GetID();
71    Int_t           GetMother();
72
73    // basic setters (for data members)
74    void  Reset();
75    void  SetOK(Bool_t ok)              {fOK = ok;}
76    void  SetBad()                      {fOK = kFALSE;}
77    void  SetGood()                     {fOK = kTRUE;}
78    void  SetLabel(Int_t label)         {fLabel = label;}
79    void  SetMotherPDG(Int_t value)     {fMotherPDG = value;}
80    void  SetRsnID(Int_t id)            {fRsnID = id;}
81    void  SetRef(AliVParticle *p)       {fRef = p;}
82    void  SetRefMC(AliVParticle *p)     {fRefMC = p;}
83    void  SetOwnerEvent(AliRsnEvent *e) {fOwnerEvent = e;}
84
85    // additional functions
86    void  FillP(Double_t mass);
87    void  Print(Option_t *o = "") const;
88
89    // getters related to charge
90    Bool_t   IsPos()      const {if (fRef) return (fRef->Charge() > 0); return kFALSE;}
91    Bool_t   IsNeg()      const {if (fRef) return (fRef->Charge() < 0); return kFALSE;}
92    Bool_t   IsCharged()  const {if (fRef) return (IsPos() || IsNeg()); return kFALSE;}
93    Bool_t   IsNeutral()  const {if (fRef) return (!IsCharged());       return kFALSE;}
94    Short_t  ChargeS()    const {if (IsPos()) return  1 ; else if (IsNeg()) return -1 ; else return  0 ;}
95    Char_t   ChargeC()    const {if (IsPos()) return '+'; else if (IsNeg()) return '-'; else return '0';}
96
97    // getters which automatically convert refs into allowed types
98    static Bool_t     Match(AliVParticle *p, TClass *ref) {if (p) return (p->InheritsFrom(ref)); return kFALSE;}
99    Bool_t            MatchRef(TClass *ref)               {return Match(fRef, ref);}
100    Bool_t            MatchRefMC(TClass *ref)             {return Match(fRefMC, ref);}
101    ERefType          RefType();
102    Bool_t            IsESD();
103    Bool_t            IsAOD();
104
105    AliVTrack        *Ref2Vtrack()        {if (Match(fRef, AliVTrack    ::Class())) return static_cast<AliVTrack *>       (fRef); return 0x0;}
106    AliESDtrack      *Ref2ESDtrack()      {if (Match(fRef, AliESDtrack  ::Class())) return static_cast<AliESDtrack *>     (fRef); return 0x0;}
107    AliAODTrack      *Ref2AODtrack()      {if (Match(fRef, AliAODTrack  ::Class())) return static_cast<AliAODTrack *>     (fRef); return 0x0;}
108    AliMCParticle    *Ref2MCparticle()    {if (Match(fRef, AliMCParticle::Class())) return static_cast<AliMCParticle *>   (fRef); return 0x0;}
109    AliAODMCParticle *Ref2AODMCparticle() {if (Match(fRef, AliMCParticle::Class())) return static_cast<AliAODMCParticle *>(fRef); return 0x0;}
110
111    AliESDv0         *Ref2ESDv0()         {if (Match(fRef, AliESDv0     ::Class())) return static_cast<AliESDv0 *>(fRef); return 0x0;}
112    AliAODv0         *Ref2AODv0()         {if (Match(fRef, AliAODv0     ::Class())) return static_cast<AliAODv0 *>(fRef); return 0x0;}
113
114    AliESDcascade    *Ref2ESDcascade()    {if (Match(fRef, AliESDcascade::Class())) return static_cast<AliESDcascade *>(fRef); return 0x0;}
115    AliAODcascade    *Ref2AODcascade()    {if (Match(fRef, AliAODcascade::Class())) return static_cast<AliAODcascade *>(fRef); return 0x0;}
116
117    AliMCParticle    *RefMC2ESD()         {if (Match(fRefMC, AliMCParticle   ::Class())) return static_cast<AliMCParticle *>   (fRef)  ; return 0x0;}
118    AliAODMCParticle *RefMC2AOD()         {if (Match(fRefMC, AliAODMCParticle::Class())) return static_cast<AliAODMCParticle *>(fRefMC); return 0x0;}
119
120    // static functions related to internal ESpecies enum
121    static ERefType    RefType(ESpecies species);
122    static Bool_t      IsCharged(ESpecies species) {return (species <= kProton);}
123    static const char *SpeciesName(ESpecies species);
124    static Int_t       SpeciesPDG(ESpecies species);
125    static Double_t    SpeciesMass(ESpecies species);
126    static EPARTYPE    ToAliPID(ESpecies species);
127    static ESpecies    FromAliPID(EPARTYPE species);
128
129 private:
130
131    Bool_t         fOK;          // internal utility flag which is kFALSE when this object should not be used
132    Int_t          fLabel;       // index of MC particle
133    Int_t          fMotherPDG;   // PDG code of mother (makes sense only if fRefMC is defined)
134    Int_t          fRsnID;       // internal ID for monitoring purposes
135
136    TLorentzVector fPrec;        // 4-momentum for rec
137    TLorentzVector fPsim;        // 4-momentum for MC
138
139    AliVParticle  *fRef;         // reference to reconstructed object
140    AliVParticle  *fRefMC;       // reference to corresponding MC particle
141    AliRsnEvent   *fOwnerEvent;  // pointer to owner event
142
143    ClassDef(AliRsnDaughter, 12)
144 };
145
146 //__________________________________________________________________________________________________
147 inline AliRsnDaughter::ERefType AliRsnDaughter::RefType()
148 {
149 //
150 // Returns the enum value corresponding to the real nature
151 // of the object pointed by the fRef data member
152 //
153
154    if (Match(fRef, AliESDtrack  ::Class())) return kTrack;
155    if (Match(fRef, AliESDv0     ::Class())) return kV0;
156    if (Match(fRef, AliESDcascade::Class())) return kCascade;
157    if (Match(fRef, AliAODTrack  ::Class())) return kTrack;
158    if (Match(fRef, AliAODv0     ::Class())) return kV0;
159    if (Match(fRef, AliAODcascade::Class())) return kCascade;
160    if (Match(fRef, AliMCParticle::Class())) return kTrack;
161
162    return kNoType;
163 }
164
165 //__________________________________________________________________________________________________
166 inline Bool_t AliRsnDaughter::IsESD()
167 {
168 //
169 // Tells if the object pointed by fRef data member is ESD
170 // NOTE: it is true even when fRef is the MC corresponding
171 //       object (= AliMCParticle)
172 //
173
174    if (Match(fRef, AliESDtrack  ::Class())) return kTRUE;
175    if (Match(fRef, AliESDv0     ::Class())) return kTRUE;
176    if (Match(fRef, AliESDcascade::Class())) return kTRUE;
177    if (Match(fRef, AliMCParticle::Class())) return kTRUE;
178
179    return kFALSE;
180 }
181
182 //__________________________________________________________________________________________________
183 inline Bool_t AliRsnDaughter::IsAOD()
184 {
185 //
186 // Tells if the object pointed by fRef data member is AOD
187 // NOTE: it is true even when fRef is the MC corresponding
188 //       object (= AliAODMCParticle)
189 //
190
191    if (Match(fRef, AliAODTrack     ::Class())) return kTRUE;
192    if (Match(fRef, AliAODv0        ::Class())) return kTRUE;
193    if (Match(fRef, AliAODcascade   ::Class())) return kTRUE;
194    if (Match(fRef, AliAODMCParticle::Class())) return kTRUE;
195
196    return kFALSE;
197 }
198
199 //__________________________________________________________________________________________________
200 inline AliRsnDaughter::ERefType AliRsnDaughter::RefType(ESpecies species)
201 {
202 //
203 // Returns the expected object type for a candidate daughter
204 // of the given species.
205 //
206
207    switch (species) {
208       case kElectron:
209       case kMuon:
210       case kPion:
211       case kKaon:
212       case kProton:
213          return kTrack;
214       case kKaon0:
215       case kLambda:
216          return kV0;
217       case kXi:
218       case kOmega:
219          return kCascade;
220       default:
221          return kNoType;
222    }
223 }
224
225 //__________________________________________________________________________________________________
226 inline void AliRsnDaughter::FillP(Double_t mass)
227 {
228 //
229 // Fills the 4-momentum data member using the values in the
230 // AliVParticle pointer data members, choosing according to arguments
231 //
232
233    if (fRefMC) fPsim.SetXYZM(fRefMC->Px(), fRefMC->Py(), fRefMC->Pz(), mass);
234    if (fRef)   fPrec.SetXYZM(fRef  ->Px(), fRef  ->Py(), fRef  ->Pz(), mass);
235 }
236
237 #endif