AliRoot Core  edcc906 (edcc906)
AliVEvent.h
Go to the documentation of this file.
1 // -*- mode: C++ -*-
2 #ifndef ALIVEVENT_H
3 #define ALIVEVENT_H
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice */
6 
7 
8 /* $Id$ */
9 
10 //-------------------------------------------------------------------------
11 // Class AliVEvent
12 //
13 // Origin: Markus Oldenburg, CERN, Markus.Oldenburg@cern.ch
14 //-------------------------------------------------------------------------
15 
16 #include <TObject.h>
17 #include <TTree.h>
18 #include <TGeoMatrix.h>
19 #include "AliVHeader.h"
20 #include "AliVParticle.h"
21 #include "AliVVertex.h"
22 #include "AliVCluster.h"
23 #include "AliVCaloCells.h"
24 #include "AliVCaloTrigger.h"
25 #include "TRefArray.h"
26 #include "AliTOFHeader.h"
27 #include "AliVTrdTrack.h"
28 #include "AliVMultiplicity.h"
29 class AliVfriendEvent;
30 class AliCentrality;
31 class AliEventplane;
32 class AliVVZERO;
33 class AliVZDC;
34 class AliVMFT; // AU
35 class AliESDkink;
36 class AliESDv0;
37 class AliESDVertex;
38 class AliESDVZERO;
39 class AliMultiplicity;
40 class AliVTrack;
41 class AliVAD;
42 
43 class AliVEvent : public TObject {
44 
45 public:
48  kMB = BIT( 0), // Minimum bias trigger in PbPb 2010-11
49  kINT1 = BIT( 0), // V0A | V0C | SPD minimum bias trigger
50  kINT7 = BIT( 1), // V0AND minimum bias trigger
51  kMUON = BIT( 2), // Single muon trigger in pp2010-11, INT1 suite
52  kHighMult = BIT( 3), // High-multiplicity SPD trigger
53  kHighMultSPD = BIT( 3), // High-multiplicity SPD trigger
54  kEMC1 = BIT( 4), // EMCAL trigger in pp2011, INT1 suite
55  kCINT5 = BIT( 5), // V0OR minimum bias trigger
56  kINT5 = BIT( 5), // V0OR minimum bias trigger
57  kCMUS5 = BIT( 6), // Single muon trigger, INT5 suite
58  kMUSPB = BIT( 6), // Single muon trigger in PbPb 2011
59  kINT7inMUON = BIT( 6), // INT7 in MUON or MUFAST cluster
60  kMuonSingleHighPt7 = BIT( 7), // Single muon high-pt, INT7 suite
61  kMUSH7 = BIT( 7), // Single muon high-pt, INT7 suite
62  kMUSHPB = BIT( 7), // Single muon high-pt in PbPb 2011
63  kMuonLikeLowPt7 = BIT( 8), // Like-sign dimuon low-pt, INT7 suite
64  kMUL7 = BIT( 8), // Like-sign dimuon low-pt, INT7 suite
65  kMuonLikePB = BIT( 8), // Like-sign dimuon low-pt in PbPb 2011
66  kMuonUnlikeLowPt7 = BIT( 9), // Unlike-sign dimuon low-pt, INT7 suite
67  kMUU7 = BIT( 9), // Unlike-sign dimuon low-pt, INT7 suite
68  kMuonUnlikePB = BIT( 9), // Unlike-sign dimuon low-pt in PbPb 2011
69  kEMC7 = BIT(10), // EMCAL/DCAL L0 trigger, INT7 suite
70  kEMC8 = BIT(10), // EMCAL/DCAL L0 trigger, INT8 suite
71  kMUS7 = BIT(11), // Single muon low-pt, INT7 suite
72  kMuonSingleLowPt7 = BIT(11), // Single muon low-pt, INT7 suite
73  kPHI1 = BIT(12), // PHOS L0 trigger in pp2011, INT1 suite
74  kPHI7 = BIT(13), // PHOS trigger, INT7 suite
75  kPHI8 = BIT(13), // PHOS trigger, INT8 suite
76  kPHOSPb = BIT(13), // PHOS trigger in PbPb 2011
77  kEMCEJE = BIT(14), // EMCAL/DCAL L1 jet trigger
78  kEMCEGA = BIT(15), // EMCAL/DCAL L1 gamma trigger
79  kHighMultV0 = BIT(16), // High-multiplicity V0 trigger
80  kCentral = BIT(16), // Central trigger in PbPb 2011
81  kSemiCentral = BIT(17), // Semicentral trigger in PbPb 2011
82  kDG = BIT(18), // Double gap diffractive
83  kDG5 = BIT(18), // Double gap diffractive
84  kZED = BIT(19), // ZDC electromagnetic dissociation
85  kSPI7 = BIT(20), // Power interaction trigger
86  kSPI = BIT(20), // Power interaction trigger
87  kINT8 = BIT(21), // 0TVX trigger
88  kMuonSingleLowPt8 = BIT(22), // Single muon low-pt, INT8 suite
89  kMuonSingleHighPt8 = BIT(23), // Single muon high-pt, INT8 suite
90  kMuonLikeLowPt8 = BIT(24), // Like-sign dimuon low-pt, INT8 suite
91  kMuonUnlikeLowPt8 = BIT(25), // Unlike-sign dimuon low-pt, INT8 suite
92  kMuonUnlikeLowPt0 = BIT(26), // Unlike-sign dimuon low-pt, no additional L0 requirement
93  kUserDefined = BIT(27), // Set when custom trigger classes are set in AliPhysicsSelection
94  kTRD = BIT(28), // TRD trigger
95  kMuonCalo = BIT(29), // Muon-calo triggers
96  // Bits 30 and above are reserved for FLAGS
97  kFastOnly = BIT(30), // The fast cluster fired. This bit is set in to addition another trigger bit, e.g. kMB
98  kAny = 0xffffffff, // to accept any defined trigger
99  kAnyINT = kMB | kINT7 | kINT5 | kINT8 | kSPI7 // to accept any interaction (aka minimum bias) trigger
100  };
101 
102  AliVEvent() { }
103  virtual ~AliVEvent() { }
104  AliVEvent(const AliVEvent& vEvnt);
105  AliVEvent& operator=(const AliVEvent& vEvnt);
106 
107  // Services
108  virtual void AddObject(TObject* obj) = 0;
109  virtual TObject* FindListObject(const char *name) const = 0;
110  virtual TList* GetList() const = 0;
111 
112  virtual void CreateStdContent() = 0;
113  virtual void GetStdContent() = 0;
114 
115  virtual void ReadFromTree(TTree *tree, Option_t* opt) = 0;
116  virtual void WriteToTree(TTree* tree) const = 0;
117 
118  virtual void Reset() = 0;
119  //virtual void ResetStdContent() = 0;
120  virtual void SetStdNames() = 0;
121 
122  virtual void Print(Option_t *option="") const = 0;
123 
124  // Header
125  virtual AliVHeader* GetHeader() const = 0;
126  //
127  // field initialization
128  virtual Bool_t InitMagneticField() const {return kFALSE;}
129 
130  // Delegated methods for fESDRun or AODHeader
131 
132  virtual void SetRunNumber(Int_t n) = 0;
133  virtual void SetPeriodNumber(UInt_t n) = 0;
134  virtual void SetMagneticField(Double_t mf) = 0;
135 
136  virtual Int_t GetRunNumber() const = 0;
137  virtual UInt_t GetPeriodNumber() const = 0;
138  virtual Double_t GetMagneticField() const = 0;
139 
140  virtual Double_t GetDiamondX() const {return -999.;}
141  virtual Double_t GetDiamondY() const {return -999.;}
142  virtual void GetDiamondCovXY(Float_t cov[3]) const
143  {cov[0]=-999.; return;}
144 
145  // Delegated methods for fHeader
146  virtual void SetOrbitNumber(UInt_t n) = 0;
147  virtual void SetBunchCrossNumber(UShort_t n) = 0;
148  virtual void SetEventType(UInt_t eventType)= 0;
149  virtual void SetTriggerMask(ULong64_t n) = 0;
150  virtual void SetTriggerCluster(UChar_t n) = 0;
151 
152  virtual UInt_t GetOrbitNumber() const = 0;
153  virtual UShort_t GetBunchCrossNumber() const = 0;
154  virtual UInt_t GetEventType() const = 0;
155  virtual ULong64_t GetTriggerMask() const = 0;
156  virtual ULong64_t GetTriggerMaskNext50() const {return 0;}
157  virtual UChar_t GetTriggerCluster() const = 0;
158  virtual TString GetFiredTriggerClasses() const = 0;
159  virtual Bool_t IsTriggerClassFired(const char* /*name*/) const {return 0;}
160  virtual Double_t GetZDCN1Energy() const = 0;
161  virtual Double_t GetZDCP1Energy() const = 0;
162  virtual Double_t GetZDCN2Energy() const = 0;
163  virtual Double_t GetZDCP2Energy() const = 0;
164  virtual Double_t GetZDCEMEnergy(Int_t i) const = 0;
165 
166  // Tracks
167  virtual AliVParticle *GetTrack(Int_t i) const = 0;
168  virtual AliVTrack *GetVTrack(Int_t /*i*/) const {return NULL;}
169  //virtual AliVTrack *GetVTrack(Int_t /*i*/) {return NULL;}
170  //virtual Int_t AddTrack(const AliVParticle *t) = 0;
171  virtual Int_t GetNumberOfTracks() const = 0;
172  virtual Int_t GetNumberOfV0s() const = 0;
173  virtual Int_t GetNumberOfCascades() const = 0;
174 
175  // TOF header and T0 methods
176  virtual const AliTOFHeader *GetTOFHeader() const {return NULL;}
177  virtual Float_t GetEventTimeSpread() const {return 0.;}
178  virtual Float_t GetTOFTimeResolution() const {return 0.;}
179  virtual Double32_t GetT0TOF(Int_t icase) const {return 0.0*icase;}
180  virtual const Double32_t * GetT0TOF() const {return NULL;}
181  virtual Float_t GetT0spread(Int_t /*i*/) const {return 0.;}
182 
183  // Calorimeter Clusters/Cells
184  virtual AliVCluster *GetCaloCluster(Int_t) const {return 0;}
185  virtual Int_t GetNumberOfCaloClusters() const {return 0;}
186  virtual Int_t GetEMCALClusters(TRefArray *) const {return 0;}
187  virtual Int_t GetPHOSClusters (TRefArray *) const {return 0;}
188  virtual AliVCaloCells *GetEMCALCells() const {return 0;}
189  virtual AliVCaloCells *GetPHOSCells() const {return 0;}
190  const TGeoHMatrix* GetPHOSMatrix(Int_t /*i*/) const {return NULL;}
191  const TGeoHMatrix* GetEMCALMatrix(Int_t /*i*/) const {return NULL;}
192  virtual AliVCaloTrigger *GetCaloTrigger(TString /*calo*/) const {return NULL;}
193 
194 
195  // Primary vertex
196  virtual Bool_t IsPileupFromSPD(Int_t /*minContributors*/,
197  Double_t /*minZdist*/,
198  Double_t /*nSigmaZdist*/,
199  Double_t /*nSigmaDiamXY*/,
200  Double_t /*nSigmaDiamZ*/)
201  const{
202  return kFALSE;
203  }
204 
205  // Tracklets
206  virtual AliVMultiplicity* GetMultiplicity() const {return 0;}
207  virtual Int_t GetNumberOfITSClusters(Int_t) const {return 0;}
208 
209  virtual Bool_t IsPileupFromSPDInMultBins() const {
210  return kFALSE;
211  }
212  virtual AliCentrality* GetCentrality() = 0;
213  virtual AliEventplane* GetEventplane() = 0;
214  virtual Int_t EventIndex(Int_t itrack) const = 0;
215  virtual Int_t EventIndexForCaloCluster(Int_t iclu) const = 0;
216  virtual Int_t EventIndexForPHOSCell(Int_t icell) const = 0;
217  virtual Int_t EventIndexForEMCALCell(Int_t icell) const = 0;
218 
219  virtual AliVVZERO *GetVZEROData() const = 0;
220  virtual const Float_t* GetVZEROEqFactors() const {return NULL;}
221  virtual Float_t GetVZEROEqMultiplicity(Int_t /* i */) const {return -1;}
222  virtual void SetVZEROEqFactors(Float_t /* factors */[64]) const {return;}
223  virtual AliVZDC *GetZDCData() const = 0;
224 
225  virtual AliVAD *GetADData() const { return NULL;}
226 
227  virtual Int_t GetNumberOfTrdTracks() const { return 0; }
228  virtual AliVTrdTrack* GetTrdTrack(Int_t /* iTrack */) const { return 0x0; }
229 
230  virtual Int_t GetNumberOfESDTracks() const { return 0; }
231  virtual Int_t GetEventNumberInFile() const {return 0;}
232 
233  //used in calibration:
234  virtual Int_t GetV0(AliESDv0&, Int_t /*iv0*/) const {return 0;}
235  virtual UInt_t GetTimeStamp() const { return 0; }
236  virtual AliVfriendEvent* FindFriend() const { return 0; }
237  virtual void SetFriendEvent( AliVfriendEvent* ) {}
238  virtual UInt_t GetEventSpecie() const { return 0; }
239  virtual AliESDkink* GetKink(Int_t /*i*/) const { return NULL; }
240  virtual Int_t GetNumberOfKinks() const { return 0; }
241 
242  virtual Int_t GetVZEROData( AliESDVZERO & ) const {return -1;}
243  virtual Int_t GetMultiplicity( AliMultiplicity & ) const {return -1;}
244 
245  // Primary vertex
246  virtual const AliVVertex *GetPrimaryVertex() const {return 0x0;}
247  virtual const AliVVertex *GetPrimaryVertexSPD() const {return 0x0;}
248  virtual const AliVVertex *GetPrimaryVertexTPC() const {return 0x0;}
249  virtual const AliVVertex *GetPrimaryVertexTracks() const {return 0x0;}
250 
251  virtual Int_t GetPrimaryVertex( AliESDVertex & ) const {return 0;}
252  virtual Int_t GetPrimaryVertexTPC( AliESDVertex & ) const {return 0;}
253  virtual Int_t GetPrimaryVertexSPD( AliESDVertex & ) const {return 0;}
254  virtual Int_t GetPrimaryVertexTracks( AliESDVertex & ) const {return 0;}
255 
256  // event status
257  virtual Bool_t IsIncompleteDAQ() {return kFALSE;}
258 
259 
260  virtual void ConnectTracks() {}
261  virtual EDataLayoutType GetDataLayoutType() const = 0;
262  const char* Whoami();
263  virtual ULong64_t GetSize() const {return 0;}
264 
265  virtual void AdjustMCLabels(const AliVEvent */*mctruth*/) {return;}
266 
267  ClassDef(AliVEvent, 3) // base class for AliEvent data
268 };
269 #endif
270 
virtual Int_t GetV0(AliESDv0 &, Int_t) const
Definition: AliVEvent.h:234
virtual Int_t GetEventNumberInFile() const
Definition: AliVEvent.h:231
virtual void Reset()=0
virtual UInt_t GetEventSpecie() const
Definition: AliVEvent.h:238
virtual Int_t GetNumberOfCascades() const =0
virtual Int_t GetNumberOfTracks() const =0
virtual const AliVVertex * GetPrimaryVertex() const
Definition: AliVEvent.h:246
virtual const AliTOFHeader * GetTOFHeader() const
Definition: AliVEvent.h:176
virtual const AliVVertex * GetPrimaryVertexTPC() const
Definition: AliVEvent.h:248
virtual void CreateStdContent()=0
virtual ULong64_t GetTriggerMaskNext50() const
Definition: AliVEvent.h:156
virtual const Double32_t * GetT0TOF() const
Definition: AliVEvent.h:180
virtual AliVCaloCells * GetPHOSCells() const
Definition: AliVEvent.h:189
virtual Double_t GetZDCP2Energy() const =0
virtual Double_t GetDiamondX() const
Definition: AliVEvent.h:140
virtual void SetVZEROEqFactors(Float_t[64]) const
Definition: AliVEvent.h:222
virtual void GetDiamondCovXY(Float_t cov[3]) const
Definition: AliVEvent.h:142
virtual Double_t GetZDCN2Energy() const =0
virtual Double_t GetMagneticField() const =0
virtual UInt_t GetOrbitNumber() const =0
virtual UShort_t GetBunchCrossNumber() const =0
virtual const AliVVertex * GetPrimaryVertexTracks() const
Definition: AliVEvent.h:249
virtual ~AliVEvent()
Definition: AliVEvent.h:103
virtual AliVfriendEvent * FindFriend() const
Definition: AliVEvent.h:236
Virtual class for calorimeter cell data handling.
Definition: AliVCaloCells.h:19
virtual Int_t GetEMCALClusters(TRefArray *) const
Definition: AliVEvent.h:186
virtual Float_t GetVZEROEqMultiplicity(Int_t) const
Definition: AliVEvent.h:221
virtual Int_t GetPrimaryVertexSPD(AliESDVertex &) const
Definition: AliVEvent.h:253
virtual TList * GetList() const =0
const char * Whoami()
Definition: AliVEvent.cxx:37
const TGeoHMatrix * GetPHOSMatrix(Int_t) const
Definition: AliVEvent.h:190
virtual AliVTrack * GetVTrack(Int_t) const
Definition: AliVEvent.h:168
virtual void ReadFromTree(TTree *tree, Option_t *opt)=0
virtual void SetTriggerMask(ULong64_t n)=0
virtual Bool_t IsPileupFromSPD(Int_t, Double_t, Double_t, Double_t, Double_t) const
Definition: AliVEvent.h:196
virtual const Float_t * GetVZEROEqFactors() const
Definition: AliVEvent.h:220
virtual void ConnectTracks()
Definition: AliVEvent.h:260
virtual void SetEventType(UInt_t eventType)=0
virtual Int_t GetMultiplicity(AliMultiplicity &) const
Definition: AliVEvent.h:243
Virtual class for calorimeter cluster data handling.
Definition: AliVCluster.h:20
virtual AliVHeader * GetHeader() const =0
virtual void WriteToTree(TTree *tree) const =0
virtual Int_t GetPHOSClusters(TRefArray *) const
Definition: AliVEvent.h:187
TTree * tree
virtual Int_t GetPrimaryVertex(AliESDVertex &) const
Definition: AliVEvent.h:251
virtual AliVCluster * GetCaloCluster(Int_t) const
Definition: AliVEvent.h:184
virtual Int_t GetNumberOfKinks() const
Definition: AliVEvent.h:240
virtual Double_t GetZDCN1Energy() const =0
Virtual class to access calorimeter (EMCAL, PHOS, PMD, FMD) trigger data.
virtual void GetStdContent()=0
virtual void SetOrbitNumber(UInt_t n)=0
virtual Int_t EventIndexForEMCALCell(Int_t icell) const =0
virtual void SetFriendEvent(AliVfriendEvent *)
Definition: AliVEvent.h:237
virtual UInt_t GetTimeStamp() const
Definition: AliVEvent.h:235
virtual Int_t GetPrimaryVertexTPC(AliESDVertex &) const
Definition: AliVEvent.h:252
virtual Int_t GetNumberOfESDTracks() const
Definition: AliVEvent.h:230
EOfflineTriggerTypes
Definition: AliVEvent.h:47
virtual AliVParticle * GetTrack(Int_t i) const =0
virtual Double_t GetZDCEMEnergy(Int_t i) const =0
virtual TObject * FindListObject(const char *name) const =0
virtual void SetPeriodNumber(UInt_t n)=0
virtual Int_t EventIndexForCaloCluster(Int_t iclu) const =0
virtual UInt_t GetEventType() const =0
virtual AliEventplane * GetEventplane()=0
virtual EDataLayoutType GetDataLayoutType() const =0
virtual AliVAD * GetADData() const
Definition: AliVEvent.h:225
virtual Int_t EventIndex(Int_t itrack) const =0
const TGeoHMatrix * GetEMCALMatrix(Int_t) const
Definition: AliVEvent.h:191
virtual void SetStdNames()=0
virtual Bool_t IsPileupFromSPDInMultBins() const
Definition: AliVEvent.h:209
virtual AliVZDC * GetZDCData() const =0
virtual Float_t GetT0spread(Int_t) const
Definition: AliVEvent.h:181
virtual UChar_t GetTriggerCluster() const =0
virtual void SetTriggerCluster(UChar_t n)=0
virtual Bool_t InitMagneticField() const
Definition: AliVEvent.h:128
virtual AliCentrality * GetCentrality()=0
AliVEvent & operator=(const AliVEvent &vEvnt)
Definition: AliVEvent.cxx:29
virtual Float_t GetTOFTimeResolution() const
Definition: AliVEvent.h:178
virtual Bool_t IsIncompleteDAQ()
Definition: AliVEvent.h:257
virtual ULong64_t GetSize() const
Definition: AliVEvent.h:263
virtual Float_t GetEventTimeSpread() const
Definition: AliVEvent.h:177
virtual Int_t GetNumberOfCaloClusters() const
Definition: AliVEvent.h:185
virtual Double_t GetZDCP1Energy() const =0
virtual Double_t GetDiamondY() const
Definition: AliVEvent.h:141
virtual Int_t EventIndexForPHOSCell(Int_t icell) const =0
virtual Int_t GetRunNumber() const =0
EDataLayoutType
Definition: AliVEvent.h:46
Definition: AliVAD.h:13
virtual TString GetFiredTriggerClasses() const =0
virtual void AddObject(TObject *obj)=0
virtual Int_t GetNumberOfTrdTracks() const
Definition: AliVEvent.h:227
virtual AliVVZERO * GetVZEROData() const =0
virtual AliVTrdTrack * GetTrdTrack(Int_t) const
Definition: AliVEvent.h:228
virtual void SetRunNumber(Int_t n)=0
virtual Double32_t GetT0TOF(Int_t icase) const
Definition: AliVEvent.h:179
virtual Bool_t IsTriggerClassFired(const char *) const
Definition: AliVEvent.h:159
virtual ULong64_t GetTriggerMask() const =0
virtual AliESDkink * GetKink(Int_t) const
Definition: AliVEvent.h:239
virtual AliVMultiplicity * GetMultiplicity() const
Definition: AliVEvent.h:206
virtual void SetBunchCrossNumber(UShort_t n)=0
virtual void AdjustMCLabels(const AliVEvent *)
Definition: AliVEvent.h:265
virtual Int_t GetPrimaryVertexTracks(AliESDVertex &) const
Definition: AliVEvent.h:254
virtual Int_t GetNumberOfV0s() const =0
virtual UInt_t GetPeriodNumber() const =0
virtual Int_t GetNumberOfITSClusters(Int_t) const
Definition: AliVEvent.h:207
virtual Int_t GetVZEROData(AliESDVZERO &) const
Definition: AliVEvent.h:242
virtual const AliVVertex * GetPrimaryVertexSPD() const
Definition: AliVEvent.h:247
virtual AliVCaloTrigger * GetCaloTrigger(TString) const
Definition: AliVEvent.h:192
virtual AliVCaloCells * GetEMCALCells() const
Definition: AliVEvent.h:188
virtual void Print(Option_t *option="") const =0
virtual void SetMagneticField(Double_t mf)=0