]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliVEvent.h
Completely reengineered version of CMake build system (Johny)
[u/mrichter/AliRoot.git] / STEER / AliVEvent.h
CommitLineData
6bc03c45 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>
c8fe2783 18#include <TGeoMatrix.h>
6bc03c45 19#include "AliVHeader.h"
20#include "AliVParticle.h"
5303c167 21#include "AliVVertex.h"
c8fe2783 22#include "AliVCluster.h"
23#include "AliVCaloCells.h"
24#include "TRefArray.h"
6bc03c45 25
26class AliVEvent : public TObject {
27
28public:
0c6c629b 29 enum EOfflineTriggerTypes {
1209509c 30 kMB = BIT(0), // Minimum bias trigger, i.e. interaction trigger, offline SPD or V0 selection
31 // empty slot
0c6c629b 32 kMUON = BIT(2), // Muon trigger, offline SPD or V0 selection
33 kHighMult = BIT(3), // High-multiplicity trigger (threshold defined online), offline SPD or V0 selection
34 kUserDefined = BIT(31), // Set when custom trigger classes are set in AliPhysicsSelection, offline SPD or V0 selection
35 kAny = 0xffffffff // to accept any trigger
36 };
6bc03c45 37
38 AliVEvent() { }
39 virtual ~AliVEvent() { }
40 AliVEvent(const AliVEvent& vEvnt);
41 AliVEvent& operator=(const AliVEvent& vEvnt);
42
43 // Services
44 virtual void AddObject(TObject* obj) = 0;
2811495d 45 virtual TObject* FindListObject(const char *name) const = 0;
6bc03c45 46 virtual TList* GetList() const = 0;
47
48 virtual void CreateStdContent() = 0;
49 virtual void GetStdContent() = 0;
50
1d0dd492 51 virtual void ReadFromTree(TTree *tree, Option_t* opt) = 0;
f12d42ce 52 virtual void WriteToTree(TTree* tree) const = 0;
6bc03c45 53
ec8aded8 54 virtual void Reset() = 0;
6bc03c45 55 //virtual void ResetStdContent() = 0;
56 virtual void SetStdNames() = 0;
57
58 virtual void Print(Option_t *option="") const = 0;
59
60 // Header
61 virtual AliVHeader* GetHeader() const = 0;
62
63 // Delegated methods for fESDRun or AODHeader
64
65 virtual void SetRunNumber(Int_t n) = 0;
66 virtual void SetPeriodNumber(UInt_t n) = 0;
67 virtual void SetMagneticField(Double_t mf) = 0;
68
69 virtual Int_t GetRunNumber() const = 0;
70 virtual UInt_t GetPeriodNumber() const = 0;
71 virtual Double_t GetMagneticField() const = 0;
72
5303c167 73 virtual Double_t GetDiamondX() const {return -999.;}
74 virtual Double_t GetDiamondY() const {return -999.;}
75 virtual void GetDiamondCovXY(Float_t cov[3]) const
76 {cov[0]=-999.; return;}
77
6bc03c45 78 // Delegated methods for fHeader
79 virtual void SetOrbitNumber(UInt_t n) = 0;
80 virtual void SetBunchCrossNumber(UShort_t n) = 0;
81 virtual void SetEventType(UInt_t eventType)= 0;
82 virtual void SetTriggerMask(ULong64_t n) = 0;
83 virtual void SetTriggerCluster(UChar_t n) = 0;
84
85 virtual UInt_t GetOrbitNumber() const = 0;
86 virtual UShort_t GetBunchCrossNumber() const = 0;
87 virtual UInt_t GetEventType() const = 0;
88 virtual ULong64_t GetTriggerMask() const = 0;
89 virtual UChar_t GetTriggerCluster() const = 0;
90
91 virtual Double_t GetZDCN1Energy() const = 0;
92 virtual Double_t GetZDCP1Energy() const = 0;
93 virtual Double_t GetZDCN2Energy() const = 0;
94 virtual Double_t GetZDCP2Energy() const = 0;
a85132e7 95 virtual Double_t GetZDCEMEnergy(Int_t i) const = 0;
6bc03c45 96
97 // Tracks
98 virtual AliVParticle *GetTrack(Int_t i) const = 0;
99 //virtual Int_t AddTrack(const AliVParticle *t) = 0;
100 virtual Int_t GetNumberOfTracks() const = 0;
c1b93f02 101 virtual Int_t GetNumberOfV0s() const = 0;
bc0379e5 102 virtual Int_t GetNumberOfCascades() const = 0;
6bc03c45 103
c8fe2783 104 // Calorimeter Clusters/Cells
105 virtual AliVCluster *GetCaloCluster(Int_t) const {return 0;}
106 virtual Int_t GetNumberOfCaloClusters() const {return 0;}
107 virtual Int_t GetEMCALClusters(TRefArray *) const {return 0;}
108 virtual Int_t GetPHOSClusters (TRefArray *) const {return 0;}
109 virtual AliVCaloCells *GetEMCALCells() const {return 0;}
110 virtual AliVCaloCells *GetPHOSCells() const {return 0;}
111 const TGeoHMatrix* GetPHOSMatrix(Int_t /*i*/) const {return NULL;}
112 const TGeoHMatrix* GetEMCALMatrix(Int_t /*i*/) const {return NULL;}
113
114
5303c167 115 // Primary vertex
116 virtual const AliVVertex *GetPrimaryVertex() const {return 0x0;}
a98c78e5 117 virtual Bool_t IsPileupFromSPD(Int_t /*minContributors*/,
118 Double_t /*minZdist*/,
119 Double_t /*nSigmaZdist*/,
120 Double_t /*nSigmaDiamXY*/,
121 Double_t /*nSigmaDiamZ*/)
122 const{
123 return kFALSE;
124 }
b46ff4b0 125
126 virtual Bool_t IsPileupFromSPDInMultBins() const {
127 return kFALSE;
128 }
129
c8fe2783 130 virtual Int_t EventIndex(Int_t itrack) const = 0;
131 virtual Int_t EventIndexForCaloCluster(Int_t iclu) const= 0;
132 virtual Int_t EventIndexForPHOSCell(Int_t icell) const= 0;
133 virtual Int_t EventIndexForEMCALCell(Int_t icell) const= 0;
134
6bc03c45 135 //---------- end of new stuff
136
137
138 /* to be considered to go in here be implemented
139
140 void SetPrimaryVertex(const AliESDVertex *vertex) {
141 *fPrimaryVertex = *vertex;
142 fPrimaryVertex->SetName("PrimaryVertex");// error prone use class wide names?
143 }
6bc03c45 144
145 void SetMultiplicity(const AliMultiplicity *mul) {
146 *fSPDMult = *mul;
147 // CKB
148 // new (&fSPDMult) AliMultiplicity(*mul);
149 }
150 const AliMultiplicity *GetMultiplicity() const {return fSPDMult;}
151
152
153 AliESDMuonTrack *GetMuonTrack(Int_t i) const {
154 return (AliESDMuonTrack *)fMuonTracks->UncheckedAt(i);
155 }
156 void AddMuonTrack(const AliESDMuonTrack *t) {
157 TClonesArray &fmu = *fMuonTracks;
158 new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
159 }
160
161 AliESDv0 *GetV0(Int_t i) const {
162 return (AliESDv0*)fV0s->UncheckedAt(i);
163 }
164 Int_t AddV0(const AliESDv0 *v);
165
166 AliESDcascade *GetCascade(Int_t i) const {
167 return (AliESDcascade *)fCascades->UncheckedAt(i);
168 }
169 void AddCascade(const AliESDcascade *c) {
170 TClonesArray &fc = *fCascades;
171 new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
172 }
173
174 AliESDkink *GetKink(Int_t i) const {
175 return (AliESDkink *)fKinks->UncheckedAt(i);
176 }
177 Int_t AddKink(const AliESDkink *c);
178
179 AliESDCaloCluster *GetCaloCluster(Int_t i) const {
180 return (AliESDCaloCluster *)fCaloClusters->UncheckedAt(i);
181 }
182 Int_t AddCaloCluster(const AliESDCaloCluster *c);
183
184 Int_t GetNumberOfMuonTracks() const {return fMuonTracks->GetEntriesFast();}
185 Int_t GetNumberOfV0s() const {return fV0s->GetEntriesFast();}
186 Int_t GetNumberOfCascades() const {return fCascades->GetEntriesFast();}
187 Int_t GetNumberOfKinks() const {return fKinks->GetEntriesFast();}
188 Int_t GetNumberOfCaloClusters() const {return fCaloClusters->GetEntriesFast();}
189
190 */
191
a8e24688 192 ClassDef(AliVEvent,1) // base class for AliEvent data
6bc03c45 193};
194#endif
195