]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAODEvent.h
Possibility to create trigger configuration from a custom file
[u/mrichter/AliRoot.git] / STEER / AliAODEvent.h
1 #ifndef AliAODEvent_H
2 #define AliAODEvent_H
3 /* Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id$ */
7
8 //-------------------------------------------------------------------------
9 //     AOD base class
10 //     Author: Markus Oldenburg, CERN
11 //-------------------------------------------------------------------------
12
13 #include <TBuffer.h>
14 #include <TClonesArray.h>
15 #include <TList.h>
16 #include <TTree.h>
17 #include <TNamed.h>
18
19 #include "AliVEvent.h"
20 #include "AliVParticle.h"
21 #include "AliAODHeader.h"
22 #include "AliAODTrack.h"
23 #include "AliAODVertex.h"
24 #include "AliAODv0.h"
25 #include "AliAODTracklets.h"
26 #include "AliAODJet.h"
27 #include "AliAODCaloCells.h"
28 #include "AliAODCaloCluster.h"
29 #include "AliAODPmdCluster.h"
30 #include "AliAODFmdCluster.h"
31
32 class TTree;
33
34 class AliAODEvent : public AliVEvent {
35
36  public :
37   enum AODListIndex_t {kAODHeader,
38                        kAODTracks,
39                        kAODVertices,
40                        kAODv0,
41                        kAODTracklets,
42                        kAODJets,
43                        kAODEmcalCells,
44                        kAODPhosCells,
45                        kAODCaloClusters,
46                        kAODFmdClusters,
47                        kAODPmdClusters,
48                        kAODListN
49   };
50
51   AliAODEvent();
52   virtual ~AliAODEvent();
53
54   //AliAODEvent(const AliAODEvent& aodevent);  // not implemented
55   //AliAODEvent& operator=(const AliAODEvent& aodevent);  // not implemented
56
57   void          AddObject(TObject *obj);
58   void          RemoveObject(TObject *obj);
59   TObject      *FindListObject(const char *objName) ;
60   TList        *GetList()                const { return fAODObjects; }
61
62   // -- Header
63   AliAODHeader *GetHeader()              const { return fHeader; }
64   void          AddHeader(const AliAODHeader* hdx)
65     {
66         delete fHeader; fHeader = new AliAODHeader(*hdx);
67         (fAODObjects->FirstLink())->SetObject(fHeader);
68     }
69
70   // setters and getters for header information
71   void     SetRunNumber(Int_t n) {if (fHeader) fHeader->SetRunNumber(n);}
72   void     SetPeriodNumber(UInt_t n){if (fHeader) fHeader->SetPeriodNumber(n);}
73   void     SetOrbitNumber(UInt_t n) {if (fHeader) fHeader->SetOrbitNumber(n);}
74   void     SetBunchCrossNumber(UShort_t n) {if (fHeader) fHeader->SetBunchCrossNumber(n);}
75   void     SetMagneticField(Double_t mf){if (fHeader) fHeader->SetMagneticField(mf);}
76
77   Int_t    GetRunNumber() const {return fHeader ? fHeader->GetRunNumber() : -999;}
78   UInt_t   GetPeriodNumber() const {return fHeader ? fHeader->GetPeriodNumber() : 0;}
79   UInt_t   GetOrbitNumber() const {return fHeader ? fHeader->GetOrbitNumber() : 0;}
80   UShort_t GetBunchCrossNumber() const {return fHeader ? fHeader->GetBunchCrossNumber() : 0;}
81   Double_t GetMagneticField() const {return fHeader ? fHeader->GetMagneticField() : -999.;}
82   
83   void      SetEventType(UInt_t eventType){fHeader->SetEventType(eventType);}
84   void      SetTriggerMask(ULong64_t n) {fHeader->SetTriggerMask(n);}
85   void      SetTriggerCluster(UChar_t n) {fHeader->SetTriggerCluster(n);}
86
87   UInt_t    GetEventType()  const { return fHeader ? fHeader->GetEventType() : 0;}
88   ULong64_t GetTriggerMask() const {return fHeader ? fHeader->GetTriggerMask() : 0;}
89   UChar_t   GetTriggerCluster() const {return fHeader ? fHeader->GetTriggerCluster() : 0;}
90
91   Double_t  GetZDCN1Energy()        const { return fHeader ? fHeader->GetZDCN1Energy() : -999.; }
92   Double_t  GetZDCP1Energy()        const { return fHeader ? fHeader->GetZDCP1Energy() : -999.; }
93   Double_t  GetZDCN2Energy()        const { return fHeader ? fHeader->GetZDCN2Energy() : -999.; }
94   Double_t  GetZDCP2Energy()        const { return fHeader ? fHeader->GetZDCP2Energy() : -999.; }
95   Double_t  GetZDCEMEnergy(Int_t i) const { return fHeader ? fHeader->GetZDCEMEnergy(i) : -999.; }
96
97   // -- Tracks
98   TClonesArray *GetTracks()              const { return fTracks; }
99   Int_t         GetNTracks()             const { return fTracks->GetEntriesFast(); }
100   Int_t         GetNumberOfTracks()      const { return GetNTracks(); }
101   AliAODTrack  *GetTrack(Int_t nTrack)   const { return (AliAODTrack*)fTracks->UncheckedAt(nTrack); }
102   Int_t         AddTrack(const AliAODTrack* trk)
103     {new((*fTracks)[fTracks->GetEntriesFast()]) AliAODTrack(*trk); return fTracks->GetEntriesFast()-1;}
104   Int_t         GetMuonTracks(TRefArray *muonTracks) const;
105
106   // -- Vertex
107   TClonesArray *GetVertices()            const { return fVertices; }
108   Int_t         GetNVertices()           const { return fVertices->GetEntriesFast(); }
109   AliAODVertex *GetVertex(Int_t nVertex) const { return (AliAODVertex*)fVertices->UncheckedAt(nVertex); }
110   Int_t         AddVertex(const AliAODVertex* vtx)
111   {new((*fVertices)[fVertices->GetEntriesFast()]) AliAODVertex(*vtx); return fVertices->GetEntriesFast()-1;}
112   
113   // primary vertex
114   virtual AliAODVertex *GetPrimaryVertex() const { return GetVertex(0); }
115
116   // V0
117   TClonesArray *GetV0s()                 const { return fV0s; }
118   Int_t         GetNV0s()                const { return fV0s->GetEntriesFast(); }
119   AliAODv0     *GetV0(Int_t nV0)         const { return (AliAODv0*)fV0s->UncheckedAt(nV0); }
120   Int_t         AddV0(const AliAODv0* v0)
121   {new((*fV0s)[fV0s->GetEntriesFast()]) AliAODv0(*v0); return fV0s->GetEntriesFast()-1;}
122
123   // -- EMCAL and PHOS Cluster
124   TClonesArray *GetCaloClusters()        const { return fCaloClusters; }
125   Int_t         GetNCaloClusters()       const { return fCaloClusters->GetEntriesFast(); }
126   AliAODCaloCluster *GetCaloCluster(Int_t nCluster) const { return (AliAODCaloCluster*)fCaloClusters->UncheckedAt(nCluster); }
127   Int_t         AddCaloCluster(const AliAODCaloCluster* clus)
128   {new((*fCaloClusters)[fCaloClusters->GetEntriesFast()]) AliAODCaloCluster(*clus); return fCaloClusters->GetEntriesFast()-1;}
129
130   // -- FMD Cluster
131   TClonesArray *GetFmdClusters()        const { return fFmdClusters; }
132   Int_t         GetNFmdClusters()       const { return fFmdClusters->GetEntriesFast(); }
133   AliAODFmdCluster *GetFmdCluster(Int_t nCluster) const { return (AliAODFmdCluster*)fFmdClusters->UncheckedAt(nCluster); }
134   Int_t         AddFmdCluster(const AliAODFmdCluster* clus)
135   {new((*fFmdClusters)[fFmdClusters->GetEntriesFast()]) AliAODFmdCluster(*clus); return fFmdClusters->GetEntriesFast()-1;}
136
137   // -- PMD Cluster
138   TClonesArray *GetPmdClusters()        const { return fPmdClusters; }
139   Int_t         GetNPmdClusters()       const { return fPmdClusters->GetEntriesFast(); }
140   AliAODPmdCluster *GetPmdCluster(Int_t nCluster) const { return (AliAODPmdCluster*)fPmdClusters->UncheckedAt(nCluster); }
141   Int_t         AddPmdCluster(const AliAODPmdCluster* clus)
142   {new((*fPmdClusters)[fPmdClusters->GetEntriesFast()]) AliAODPmdCluster(*clus); return fPmdClusters->GetEntriesFast()-1;}
143
144   // -- Jet
145   TClonesArray *GetJets()            const { return fJets; }
146   Int_t         GetNJets()           const { return fJets->GetEntriesFast(); }
147   AliAODJet    *GetJet(Int_t nJet) const { return (AliAODJet*)fJets->UncheckedAt(nJet); }
148   Int_t         AddJet(const AliAODJet* vtx)
149     {new((*fJets)[fJets->GetEntriesFast()]) AliAODJet(*vtx); return fJets->GetEntriesFast()-1;}
150
151   // -- Tracklets
152   AliAODTracklets *GetTracklets() const { return fTracklets; }
153
154   // -- Calorimeter Cells
155   AliAODCaloCells *GetEMCALCells() const { return fEmcalCells; }
156   AliAODCaloCells *GetPHOSCells() const { return fPhosCells; }
157
158   // -- Services
159   void    CreateStdContent();
160   void    SetStdNames();
161   void    GetStdContent();
162   void    ResetStd(Int_t trkArrSize = 0, 
163                    Int_t vtxArrSize = 0, 
164                    Int_t v0ArrSize = 0, 
165                    Int_t jetSize = 0, 
166                    Int_t caloClusSize = 0, 
167                    Int_t fmdClusSize = 0, 
168                    Int_t pmdClusSize = 0);
169   void    ClearStd();
170   void    ReadFromTree(TTree *tree);
171   const void WriteToTree(TTree* tree) const {tree->Branch(fAODObjects);}
172
173   void  Print(Option_t *option="") const;
174
175  private :
176
177   TList *fAODObjects; //  list of AODObjects
178   
179   // standard content
180   AliAODHeader    *fHeader;       //! event information
181   TClonesArray    *fTracks;       //! charged tracks
182   TClonesArray    *fVertices;     //! vertices
183   TClonesArray    *fV0s;          //! V0s
184   AliAODTracklets *fTracklets;    //! SPD tracklets
185   TClonesArray    *fJets;         //! jets
186   AliAODCaloCells *fEmcalCells;   //! EMCAL calorimenter cells
187   AliAODCaloCells *fPhosCells;    //! PHOS calorimenter cells
188   TClonesArray    *fCaloClusters; //! calorimeter clusters
189   TClonesArray    *fFmdClusters;  //! FMDclusters
190   TClonesArray    *fPmdClusters;  //! PMDclusters
191
192   static const char* fAODListName[kAODListN]; //!
193
194   ClassDef(AliAODEvent,3);
195 };
196
197 #endif