]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/hfe/AliHFEreducedEvent.h
Fix of sigmaZ for crossing tracklets from Alex
[u/mrichter/AliRoot.git] / PWGHF / hfe / AliHFEreducedEvent.h
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 // Debug event to look at the distribution of the variable we are cutting on
17 //
18 //
19 #ifndef ALIHFEREDUCEDEVENT_H
20 #define ALIHFEREDUCEDEVENT_H
21
22 #include <TObject.h>
23
24 class TObjArray;
25 class AliHFEreducedTrack;
26 class AliHFEreducedMCParticle;
27
28 class AliHFEreducedEvent : public TObject{
29  public:
30   AliHFEreducedEvent();
31   AliHFEreducedEvent(const AliHFEreducedEvent &ref);
32   AliHFEreducedEvent &operator=(const AliHFEreducedEvent &ref);
33   ~AliHFEreducedEvent();
34   
35   void AddTrack(const AliHFEreducedTrack *track);
36   const AliHFEreducedTrack *GetTrack(int itrk) const;
37   Int_t GetNumberOfTracks() const { return fNtracks; }
38   void AddMCParticle(const AliHFEreducedMCParticle *mctrack);
39   const AliHFEreducedMCParticle *GetMCParticle(int itrk) const;
40   Int_t GetNumberOfMCParticles() const { return fNmcparticles; }
41   
42   Float_t GetVX() const { return fVX[0]; }
43   Float_t GetVY() const { return fVY[0]; }
44   Float_t GetVZ() const { return fVZ[0]; }
45   Float_t GetVXSPD() const { return fVX[1]; }
46   Float_t GetVYSPD() const { return fVY[1]; }
47   Float_t GetVZSPD() const { return fVZ[1]; }
48   Double_t GetVXMC() const { return fVMC[0]; }
49   Double_t GetVYMC() const { return fVMC[1]; }
50   Double_t GetVZMC() const { return fVMC[2]; }
51   Int_t GetNContribVertex() const { return fNContrib[0]; }
52   Int_t GetNContribVertexSPD() const { return fNContrib[1]; }
53   Bool_t HasPrimaryVertex() const { return fNContrib[0] > 0; }
54   Bool_t HasPrimaryVertexSPD() const { return fNContrib[1] > 0; }
55   Float_t GetVertexZResolution() const { return fVertexResolution[0]; };
56   Float_t GetVertexZResolutionSPD() const { return fVertexResolution[1]; };
57   Float_t GetVertexDispersion()  const { return fVertexDispersion[0]; };
58   Float_t GetVertexDispersionSPD()  const { return fVertexDispersion[1]; };
59   Int_t GetRunNumber() const { return fRunNumber; }
60   Double_t GetCentrality() const { return fCentrality[0]; }
61   Double_t GetCentralityV0M() const { return fCentrality[0]; }
62   Double_t GetCentralityV0A() const { return fCentrality[1]; }
63   Double_t GetCentralityV0C() const { return fCentrality[2]; }
64   Double_t GetCentralityTracklets() const { return fCentrality[3]; }
65   Double_t GetCentralityTracks() const { return fCentrality[4]; }
66   Double_t GetCentralityZNA() const { return fCentrality[5]; }
67   Double_t GetCentralityZNC() const { return fCentrality[6]; }
68   Double_t GetCentralityCL0() const { return fCentrality[7]; }
69   Double_t GetCentralityCL1() const { return fCentrality[8]; }
70   Double_t GetCentralityCND() const { return fCentrality[9]; }
71   Float_t GetV0AMultiplicity() const { return fV0Multiplicity[0]; }
72   Float_t GetV0CMultiplicity() const { return fV0Multiplicity[1]; }
73   Float_t GetV0MMultiplicity() const { return fV0Multiplicity[0] + fV0Multiplicity[1]; }
74   Float_t GetZNAEnergy() const { return fZDCEnergy[0]; }
75   Float_t GetZNCEnergy() const { return fZDCEnergy[1]; }
76   Float_t GetZPAEnergy() const { return fZDCEnergy[2]; }
77   Float_t GetZPCEnergy() const { return fZDCEnergy[3]; }
78   Float_t GetZDCNEnergySum() const { return fZDCEnergy[0] + fZDCEnergy[1]; }
79   Float_t GetZDCNEnergyDifference() const { return fZDCEnergy[0] - fZDCEnergy[1]; }
80   Float_t GetZDCNEnergyAsymmetry() const { return fZDCEnergy[0] + fZDCEnergy[1] != 0 ? (fZDCEnergy[0] - fZDCEnergy[1])/(fZDCEnergy[0] + fZDCEnergy[1]) : 1.; }
81   Float_t GetZDCPEnergySum() const { return fZDCEnergy[2] + fZDCEnergy[3]; }
82   Float_t GetZDCPEnergyDifference() const { return fZDCEnergy[2] - fZDCEnergy[3]; }
83   Float_t GetZDCPEnergyAsymmetry() const { return fZDCEnergy[2] + fZDCEnergy[3] != 0 ? (fZDCEnergy[2] - fZDCEnergy[3])/(fZDCEnergy[2] + fZDCEnergy[3]) : 1.; }
84   Int_t   GetSPDMultiplicity() const { return fSPDMultiplicity; }
85   
86   void SetVX(Float_t vx) { fVX[0] = vx; }
87   void SetVY(Float_t vy) { fVY[0] = vy; }
88   void SetVZ(Float_t vz) { fVZ[0] = vz; }
89   void SetVXSPD(Float_t vx) { fVX[1] = vx; }
90   void SetVYSPD(Float_t vy) { fVY[1] = vy; }
91   void SetVZSPD(Float_t vz) { fVZ[1] = vz; }
92   void SetVMC(Double_t vx, Double_t vy, Double_t vz){
93       fVMC[0] = vx;
94       fVMC[1] = vy;
95       fVMC[2] = vz;
96   }
97   void SetRunNumber(Int_t runnumber) { fRunNumber = runnumber; }
98   void SetPileupFlag() { fPileupFlag = kTRUE; }
99   inline void SetCentrality(
100         Float_t centV0M, 
101         Float_t centV0A, 
102         Float_t centV0C, 
103         Float_t centTLS, 
104         Float_t centTrks, 
105         Float_t centZNA, 
106         Float_t centZNC,
107         Float_t centCL0,
108         Float_t centCL1,
109         Float_t centCND
110   );
111   void SetNContribVertex(Int_t ncontrib) { fNContrib[0] = ncontrib; }
112   void SetNContribVertexSPD(Int_t ncontrib) { fNContrib[1] = ncontrib; }
113   void SetVertexResolution(Float_t res) { fVertexResolution[0] = res; }
114   void SetVertexResolutionSPD(Float_t res) { fVertexResolution[1] = res; }
115   void SetVertexDispersion(Float_t dis) { fVertexDispersion[0] = dis; }
116   void SetVertexDispersionSPD(Float_t dis) { fVertexDispersion[1] = dis; }
117   
118   Bool_t IsMBTrigger() const { return TESTBIT(fTrigger, kMB); }
119   Bool_t IsSemiCentralTrigger() const { return TESTBIT(fTrigger, kSemiCentral); }
120   Bool_t IsCentralTrigger() const { return TESTBIT(fTrigger, kCentral); }
121   Bool_t IsEMCalTrigger() const { return TESTBIT(fTrigger, kEMCAL); }
122   Bool_t IsTRDSETrigger() const { return TESTBIT(fTrigger, kTRDSE); }
123   Bool_t IsTRDDQTrigger() const { return TESTBIT(fTrigger, kTRDDQ); }
124   Bool_t IsINTTrigger() const { return TESTBIT(fTrigger, kINTTRG); }
125   Bool_t HasPileupFlag() const { return fPileupFlag; }
126   
127   void SetMBTrigger() { SETBIT(fTrigger, kMB); }
128   void SetSemiCentralTrigger() { SETBIT(fTrigger, kSemiCentral); }
129   void SetCentralTrigger() { SETBIT(fTrigger, kCentral); }
130   void SetEMCALTrigger() { SETBIT(fTrigger, kEMCAL); }
131   void SetTRDSETrigger() { SETBIT(fTrigger, kTRDSE); }
132   void SetTRDDQTrigger() { SETBIT(fTrigger, kTRDDQ); }
133   void SetINTTrigger() { SETBIT(fTrigger, kINTTRG); }
134
135   void SetV0Multiplicity(Float_t v0A, Float_t v0C) {
136     fV0Multiplicity[0] = v0A;
137     fV0Multiplicity[1] = v0C;
138   }
139   inline void SetZDCEnergy(Float_t zna, Float_t znc, Float_t zpa, Float_t zpc);
140   void SetSPDMultiplicity(Int_t mult) { fSPDMultiplicity = mult; }
141   
142  private:
143   typedef enum{
144     kMB = 0,
145     kSemiCentral = 1,
146     kCentral = 2,
147     kEMCAL = 3,
148     kTRDSE = 4,
149     kTRDDQ = 5,
150     kINTTRG = 6
151   } Trigger_t;
152   enum{
153     kCentBuff = 15
154   };
155   TObjArray *fTracks;           // Array with reconstructed tracks
156   TObjArray *fMCparticles;      // Array with MC particles
157   Int_t fNtracks;               // Number of tracks
158   Int_t fNmcparticles;          // Number of MC Particles
159   Int_t fRunNumber;             // Run Number
160   Float_t  fCentrality[kCentBuff];     // Centrality (V0M, V0A, V0C, TLS, TRK, ZNA, ZNC, CL0, CL1, CND)
161   Int_t fTrigger;               // Trigger bits
162   Float_t fVX[2];               // Vertex X
163   Float_t fVY[2];               // Vertex Y
164   Float_t fVZ[2];               // Vertex Z
165   Double_t fVMC[3];              // Position of the MC Vertex
166   Int_t    fNContrib[2];        // Number of vertex contributors
167   Float_t  fVertexResolution[2];// z-Vertex resolution 
168   Float_t  fVertexDispersion[2];// z-Vertex dispersion
169   Float_t  fV0Multiplicity[2];  // V0 multiplicity
170   Float_t  fZDCEnergy[4];       // ZDC Energy (n,p)
171   Int_t    fSPDMultiplicity;    // SPD tracklet multiplicity
172   Bool_t   fPileupFlag;         // Flag for Pileup event
173   
174   ClassDef(AliHFEreducedEvent, 5)
175 };
176
177 //____________________________________________________________
178 void AliHFEreducedEvent::SetCentrality(
179         Float_t centV0M, 
180         Float_t centV0A, 
181         Float_t centV0C, 
182         Float_t centTLS, 
183         Float_t centTrks, 
184         Float_t centZNA, 
185         Float_t centZNC,
186         Float_t centCL0,
187         Float_t centCL1,
188         Float_t centCND
189   ) 
190
191     fCentrality[0] = centV0M; 
192     fCentrality[1] = centV0A; 
193     fCentrality[2] = centV0C; 
194     fCentrality[3] = centTLS; 
195     fCentrality[4] = centTrks; 
196     fCentrality[5] = centZNA; 
197     fCentrality[6] = centZNC; 
198     fCentrality[7] = centCL0; 
199     fCentrality[8] = centCL1; 
200     fCentrality[9] = centCND; 
201 }
202
203 //_________________________________________________________
204 void AliHFEreducedEvent::SetZDCEnergy(Float_t zna, Float_t znc, Float_t zpa, Float_t zpc){
205   fZDCEnergy[0] = zna;
206   fZDCEnergy[1] = znc;
207   fZDCEnergy[2] = zpa;
208   fZDCEnergy[3] = zpc;
209 }
210 #endif