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