]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnDaughter.h
fix for Coverity (B.Hippolyte)
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnDaughter.h
1 #ifndef ALIRSNDAUGHTER_H
2 #define ALIRSNDAUGHTER_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6  
7 ////////////////////////////////////////////////////////////////////////////////
8 //
9 //  Interface to single daughter candidate.
10 //
11 ////////////////////////////////////////////////////////////////////////////////
12
13 #include <TLorentzVector.h>
14
15 #include "AliPID.h"
16 #include "AliVTrack.h"
17 #include "AliESDtrack.h"
18 #include "AliESDv0.h"
19 #include "AliESDcascade.h"
20 #include "AliMCParticle.h"
21 #include "AliAODTrack.h"
22 #include "AliAODv0.h"
23 #include "AliAODcascade.h"
24 #include "AliAODMCParticle.h"
25
26 typedef AliPID::EParticleType EPARTYPE;
27
28 class AliRsnDaughter : public TObject {
29 public:
30
31    enum ERefType {
32       kTrack = 0,
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();
52    AliRsnDaughter(const AliRsnDaughter &copy);
53    AliRsnDaughter& operator= (const AliRsnDaughter& copy);
54    virtual ~AliRsnDaughter() { /*empty, since pointers must not be deleted*/ }
55
56    // basic getters (for data members)
57    Bool_t          IsOK()         const {return fOK;}
58    Int_t           GetLabel()     const {return fLabel;}
59    Int_t           GetMotherPDG() const {return fMotherPDG;}
60    Int_t           GetRsnID()     const {return fRsnID;}
61    TLorentzVector& Prec()               {return fPrec;}
62    TLorentzVector& Psim()               {return fPsim;}
63    TLorentzVector& P(Bool_t mc)         {return (mc ? fPsim : fPrec);}
64    AliVParticle*   GetRef()             {return fRef;}
65    AliVParticle*   GetRefMC()           {return fRefMC;}
66
67    // basic setters (for data members)
68    void    SetBad()                  {fOK = kFALSE;}
69    void    SetGood()                 {fOK = kTRUE;}
70    void    SetLabel(Int_t label)     {fLabel = label;}
71    void    SetRsnID(Int_t id)        {fRsnID = id;}
72    void    SetMotherPDG(Int_t value) {fMotherPDG = value;}
73    void    SetRef(AliVParticle *p);
74    void    SetRefMC(AliVParticle *p);
75
76    // charge checkers
77    Bool_t  IsPos()             const {return (fRef->Charge() > 0);}
78    Bool_t  IsNeg()             const {return (fRef->Charge() < 0);}
79    Bool_t  IsNeutral()         const {return (!IsPos() && !IsNeg());}
80    Bool_t  IsSign(Char_t sign) const {if (sign == '+') return IsPos(); else if (sign == '-') return IsNeg(); else return IsNeutral();}
81    Short_t ChargeS()           const {if (IsPos()) return  1 ; else if (IsNeg()) return -1 ; else return  0 ;}
82    Char_t  ChargeC()           const {if (IsPos()) return '+'; else if (IsNeg()) return '-'; else return '0';}
83
84    // getters which automatically convert refs into allowed types
85    AliVTrack*        GetRefVtrack()     {if (classMatchRef  (AliVTrack       ::Class())) return static_cast<AliVTrack*>       (fRef)  ; return 0x0;}
86    AliESDtrack*      GetRefESDtrack()   {if (classMatchRef  (AliESDtrack     ::Class())) return static_cast<AliESDtrack*>     (fRef)  ; return 0x0;}
87    AliESDv0*         GetRefESDv0()      {if (classMatchRef  (AliESDv0        ::Class())) return static_cast<AliESDv0*>        (fRef)  ; return 0x0;}
88    AliESDcascade*    GetRefESDcascade() {if (classMatchRef  (AliESDcascade   ::Class())) return static_cast<AliESDcascade*>   (fRef)  ; return 0x0;}
89    AliAODTrack*      GetRefAODtrack()   {if (classMatchRef  (AliAODTrack     ::Class())) return static_cast<AliAODTrack*>     (fRef)  ; return 0x0;}
90    AliAODv0*         GetRefAODv0()      {if (classMatchRef  (AliAODv0        ::Class())) return static_cast<AliAODv0*>        (fRef)  ; return 0x0;}
91    AliAODcascade*    GetRefAODcascade() {if (classMatchRef  (AliAODcascade   ::Class())) return static_cast<AliAODcascade*>   (fRef)  ; return 0x0;}
92    AliMCParticle*    GetRefMCtrack()    {if (classMatchRef  (AliMCParticle   ::Class())) return static_cast<AliMCParticle*>   (fRef)  ; return 0x0;}
93    AliMCParticle*    GetRefMCESD()      {if (classMatchRefMC(AliMCParticle   ::Class())) return static_cast<AliMCParticle*>   (fRefMC); return 0x0;}
94    AliAODMCParticle* GetRefMCAOD()      {if (classMatchRefMC(AliAODMCParticle::Class())) return static_cast<AliAODMCParticle*>(fRefMC); return 0x0;}
95
96    // check the input type
97    Bool_t    IsMC()      {if (GetRefMCtrack()) return kTRUE; return kFALSE;}
98    Bool_t    IsAOD()     {if (GetRefAODtrack() || GetRefAODv0() || GetRefESDcascade()) return kTRUE; return kFALSE;}
99    Bool_t    IsESD()     {if (GetRefESDtrack() || GetRefESDv0() || GetRefAODcascade()) return kTRUE; return kFALSE;}
100    Bool_t    IsTrack()   {if (GetRefESDtrack() || GetRefAODtrack() || GetRefMCtrack()) return kTRUE; return kFALSE;}
101    Bool_t    IsV0()      {if (GetRefESDv0() || GetRefAODv0()) return kTRUE; return kFALSE;}
102    Bool_t    IsCascade() {if (GetRefESDcascade() || GetRefAODcascade()) return kTRUE; return kFALSE;}
103    ERefType  RefType()   {if (IsTrack()) return kTrack; if (IsV0()) return kV0; if (IsCascade()) return kCascade; return kNoType;}
104    
105    // utilities
106    void     Reset();
107    Int_t    GetPDG(Bool_t abs = kTRUE);
108    Int_t    GetID();
109    Bool_t   SetMass(Double_t mass);
110    Bool_t   IsKinkDaughter();
111    
112    // static utilities
113    static ERefType    RefType(ESpecies species);
114    static Bool_t      IsCharged(ESpecies species) {return (species <= kProton);}
115    static const char* SpeciesName(ESpecies species);
116    static Int_t       SpeciesPDG(ESpecies species);
117    static Double_t    SpeciesMass(ESpecies species);
118    static EPARTYPE    ToAliPID(ESpecies species);
119    static ESpecies    FromAliPID(EPARTYPE species);
120
121 private:
122
123    Bool_t classMatchRef  (TClass *ref) {if (fRef  ) return (fRef  ->InheritsFrom(ref)); return kFALSE;}
124    Bool_t classMatchRefMC(TClass *ref) {if (fRefMC) return (fRefMC->InheritsFrom(ref)); return kFALSE;}
125
126    Bool_t         fOK;          // internal utility flag which is kFALSE when this object should not be used
127    Int_t          fLabel;       // GEANT label of corresponding MC particle
128    Int_t          fMotherPDG;   // PDG code of mother (makes sense only if fRefMC is defined)
129    Int_t          fRsnID;       // internal ID for monitoring purposes
130
131    TLorentzVector fPrec;        // 4-vector filled with track info from default ref (if present)
132    TLorentzVector fPsim;        // 4-vector filled with track info from MC ref (if present)
133
134    AliVParticle  *fRef;         // reference to reconstructed object
135    AliVParticle  *fRefMC;       // reference to corresponding MC particle
136
137    ClassDef(AliRsnDaughter, 9)
138 };
139
140
141
142 #endif