]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TUHKMgen/UHKM/Particle.h
updated macros for making PPR plots
[u/mrichter/AliRoot.git] / TUHKMgen / UHKM / Particle.h
1 /*                                                                            
2                                                                             
3         Nikolai Amelin, Ludmila Malinina, Timur Pocheptsov (C) JINR/Dubna
4       amelin@sunhe.jinr.ru, malinina@sunhe.jinr.ru, pocheptsov@sunhe.jinr.ru 
5                            November. 2, 2005                                
6
7 */
8
9 #ifndef PARTICLE_INCLUDED
10 #define PARTICLE_INCLUDED
11
12 #include <list>
13
14 #include <TLorentzRotation.h>
15 #include <TLorentzVector.h>
16 #include <TVector3.h>
17 #ifndef PARTICLE_PDG
18 #include "ParticlePDG.h"
19 #endif
20 #include <iostream>
21 using namespace std;
22 //class ParticlePDG;
23
24 class Particle {
25  public:
26   Particle(const TLorentzVector &, const TLorentzVector &);
27   Particle(const Particle& copy);
28   Particle& operator=(const Particle&);
29   virtual ~Particle() {};
30   Particle(ParticlePDG *pdg = 0);
31   Particle(ParticlePDG *pdg, const TLorentzVector &pos, const TLorentzVector &mom,
32            Double_t lastInterTime = 0., Int_t lastInterNum = 0, Int_t type=0);
33   Particle(ParticlePDG *pdg, const TLorentzVector &pos, const TLorentzVector &mom,
34            Double_t lastInterTime, Int_t lastInterNum, Int_t type, Int_t motherPdg, 
35            const TLorentzVector &motherPos, const TLorentzVector &motherMom);
36
37   Double_t X()const{return fPosition.X();}
38   Double_t X(Double_t val){fPosition.SetX(val); return val;}
39   Double_t Y()const{return fPosition.Y();}
40   Double_t Y(Double_t val){fPosition.SetY(val); return val;}
41   Double_t Z()const{return fPosition.Z();}
42   Double_t Z(Double_t val){fPosition.SetZ(val); return val;}
43   Double_t T()const{return fPosition.T();}
44   Double_t T(Double_t val){fPosition.SetT(val); return val;}
45   Double_t Px()const{return fMomentum.Px();}
46   Double_t Px(Double_t val){fMomentum.SetPx(val); return val;}
47   Double_t Py()const{return fMomentum.Py();}
48   Double_t Py(Double_t val){fMomentum.SetPy(val); return val;}
49   Double_t Pz()const{return fMomentum.Pz();}
50   Double_t Pz(Double_t val){fMomentum.SetPz(val); return val;}
51   Double_t E()const{return fMomentum.E();}
52   Double_t E(Double_t val){fMomentum.SetE(val); return val;}
53
54   TLorentzVector &Pos(){return fPosition;}
55   const TLorentzVector &Pos()const{return fPosition;}
56   TLorentzVector &Pos(const TLorentzVector &val){return fPosition = val;}
57   TLorentzVector &Mom(){return fMomentum;}
58   const TLorentzVector &Mom()const{return fMomentum;}
59   TLorentzVector &Mom(const TLorentzVector &val){return fMomentum = val;}
60
61   void SetDecayed() {fDecayed = kTRUE;}
62   Bool_t GetDecayed() const {return fDecayed;}
63                 
64   void Boost(const TVector3 &val){fMomentum.Boost(val);}
65   void Boost(const TLorentzVector &val){fMomentum.Boost(val.BoostVector());}
66   void TransformMomentum(const TRotation &rotator){fMomentum *= rotator;}
67   void TransformPosition(const TRotation &rotator){fPosition *= rotator;}
68   void Shift(const TVector3 &val){fPosition += TLorentzVector(val, 0.);}
69
70   //Pseudorapidity
71   Double_t Eta ()const;
72   //Rapidity
73   Double_t Rapidity()const;
74   Double_t Phi()const;
75   Double_t Theta()const;
76   Double_t Pt()const;
77
78   Int_t Encoding() const;
79   Double_t TableMass() const;
80   ParticlePDG *Def() const {return fParticleProperties;}
81   ParticlePDG *Def(ParticlePDG *newProp) {return fParticleProperties = newProp;}
82   //mother   
83   void SetLastMotherPdg(Int_t value){fLastMotherPdg = value;}
84   Int_t GetLastMotherPdg() const {return fLastMotherPdg;}
85
86   // aic(2008/08/08): functions added in order to enable tracking of mother/daughter particles by a unique index
87   // The index coincides with the position of the particle in the secondaries list.
88   Int_t SetIndex() {
89     fIndex = ++fLastIndex; 
90     return fIndex;
91   }
92   Int_t GetIndex() const {return fIndex;}
93   static Int_t GetLastIndex() {return fLastIndex;}
94   static void InitIndexing() {
95     fLastIndex = -1;
96   }
97   void SetMother(Int_t value) {fMotherIndex = value;}
98   Int_t GetMother() const {return fMotherIndex;}
99   void SetFirstDaughterIndex(Int_t index) {fFirstDaughterIndex = index;}
100   void SetLastDaughterIndex(Int_t index) {fLastDaughterIndex = index;}
101   void SetPythiaStatusCode(Int_t code) {fPythiaStatusCode = code;}
102   Int_t GetPythiaStatusCode() const {return fPythiaStatusCode;}
103
104   //  Int_t GetNDaughters() const {return fNDaughters;}
105   Int_t GetNDaughters() const {
106     if(fFirstDaughterIndex==-1 || fLastDaughterIndex==-1) 
107       return 0;
108     else
109       return fLastDaughterIndex-fFirstDaughterIndex+1;
110   }
111   Int_t GetFirstDaughterIndex() const {return fFirstDaughterIndex;}
112   Int_t GetLastDaughterIndex() const {return fLastDaughterIndex;}
113
114   TLorentzVector &SetLastMotherDecayCoor(const TLorentzVector &val){return fLastMotherDecayCoor = val;}
115   const TLorentzVector &GetLastMotherDecayCoor()const{return fLastMotherDecayCoor;}
116   TLorentzVector &SetLastMotherDecayMom(const TLorentzVector &val){return fLastMotherDecayMom = val;}
117   const TLorentzVector &GetLastMotherDecayMom()const{return fLastMotherDecayMom;}
118
119   void SetLastInterTime(Double_t value){fLastInteractionTime = value;}
120   Double_t GetLastInterTime()const{return fLastInteractionTime;}
121   void SetLastInterNumber(Int_t value){fInteractionNumber = value;}
122   Int_t GetLastInterNumber()const{return fInteractionNumber;}
123   void IncInter(){++fInteractionNumber;}
124
125   void SetType(Int_t value){fType = value;}
126   Int_t GetType()const{return fType;}
127
128  protected:
129   TLorentzVector   fPosition;
130   TLorentzVector   fMomentum;
131   TLorentzVector   fLastMotherDecayCoor;
132   TLorentzVector   fLastMotherDecayMom;
133   ParticlePDG     *fParticleProperties;
134   Double_t         fLastInteractionTime;
135   Int_t            fInteractionNumber;
136   Int_t            fPythiaStatusCode;
137   Int_t            fLastMotherPdg;
138   Int_t            fType; //0-hydro, 1-jets
139   Int_t            fIndex;                    // index (0 based) of particle in the final particle list which will contain both primaries and secondaries
140   Int_t            fMotherIndex;              // index of the mother (-1 if its a primary particle)
141   Int_t            fNDaughters;               // number of daughter particles (0 if the particle had not decayed)
142   Int_t            fFirstDaughterIndex;       // index for the first daughter particle (-1 if non-existing)
143   Int_t            fLastDaughterIndex;        // index for the last daughter particle (-1 if non-existing)
144   Bool_t           fDecayed;                  // true if the decay procedure already applied
145   static Int_t     fLastIndex;                // the last index assigned
146 };
147
148 Double_t S(const TLorentzVector &, const TLorentzVector &);
149 Double_t T(const TLorentzVector &, const TLorentzVector &);
150
151 typedef std::list<Particle> List_t;
152 typedef std::list<Particle>::iterator LPIT_t;
153
154 class ParticleAllocator {
155  public:
156   ParticleAllocator() {};
157   void AddParticle(const Particle & particle, List_t & list);
158   void FreeListNode(List_t & list, LPIT_t it);
159   void FreeList(List_t & list);
160
161  private:
162   List_t fFreeNodes;
163 };
164
165 #endif