]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - STEER/AliAODEvent.h
Fixing ompilation warnings
[u/mrichter/AliRoot.git] / STEER / AliAODEvent.h
... / ...
CommitLineData
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 "AliAODcascade.h"
26#include "AliAODTracklets.h"
27#include "AliAODJet.h"
28#include "AliAODCaloCells.h"
29#include "AliAODCaloCluster.h"
30#include "AliAODPmdCluster.h"
31#include "AliAODFmdCluster.h"
32
33class TTree;
34class TFolder;
35
36class AliAODEvent : public AliVEvent {
37
38 public :
39 enum AODListIndex_t {kAODHeader,
40 kAODTracks,
41 kAODVertices,
42 kAODv0,
43 kAODcascade,
44 kAODTracklets,
45 kAODJets,
46 kAODEmcalCells,
47 kAODPhosCells,
48 kAODCaloClusters,
49 kAODFmdClusters,
50 kAODPmdClusters,
51 kAODListN
52 };
53
54 AliAODEvent();
55 virtual ~AliAODEvent();
56
57 AliAODEvent(const AliAODEvent& aodevent);
58 AliAODEvent& operator=(const AliAODEvent& aodevent);
59
60 void AddObject(TObject *obj);
61 void RemoveObject(TObject *obj);
62 TObject *FindListObject(const char *objName) ;
63 TList *GetList() const { return fAODObjects; }
64
65 // -- Header
66 AliAODHeader *GetHeader() const { return fHeader; }
67 void AddHeader(const AliAODHeader* hdx)
68 {
69 delete fHeader; fHeader = new AliAODHeader(*hdx);
70 (fAODObjects->FirstLink())->SetObject(fHeader);
71 }
72
73 // setters and getters for header information
74 void SetRunNumber(Int_t n) {if (fHeader) fHeader->SetRunNumber(n);}
75 void SetPeriodNumber(UInt_t n){if (fHeader) fHeader->SetPeriodNumber(n);}
76 void SetOrbitNumber(UInt_t n) {if (fHeader) fHeader->SetOrbitNumber(n);}
77 void SetBunchCrossNumber(UShort_t n) {if (fHeader) fHeader->SetBunchCrossNumber(n);}
78 void SetMagneticField(Double_t mf){if (fHeader) fHeader->SetMagneticField(mf);}
79 void SetDiamond(Float_t xy[2],Float_t cov[3]){if (fHeader) fHeader->SetDiamond(xy,cov);}
80
81 Int_t GetRunNumber() const {return fHeader ? fHeader->GetRunNumber() : -999;}
82 UInt_t GetPeriodNumber() const {return fHeader ? fHeader->GetPeriodNumber() : 0;}
83 UInt_t GetOrbitNumber() const {return fHeader ? fHeader->GetOrbitNumber() : 0;}
84 UShort_t GetBunchCrossNumber() const {return fHeader ? fHeader->GetBunchCrossNumber() : 0;}
85 Double_t GetMagneticField() const {return fHeader ? fHeader->GetMagneticField() : -999.;}
86 Double_t GetDiamondX() const {return fHeader ? fHeader->GetDiamondX() : -999.;}
87 Double_t GetDiamondY() const {return fHeader ? fHeader->GetDiamondY() : -999.;}
88 void GetDiamondCovXY(Float_t cov[3]) const {cov[0]=-999.; if(fHeader) fHeader->GetDiamondCovXY(cov);}
89
90 void SetEventType(UInt_t eventType){fHeader->SetEventType(eventType);}
91 void SetTriggerMask(ULong64_t n) {fHeader->SetTriggerMask(n);}
92 void SetTriggerCluster(UChar_t n) {fHeader->SetTriggerCluster(n);}
93
94 UInt_t GetEventType() const { return fHeader ? fHeader->GetEventType() : 0;}
95 ULong64_t GetTriggerMask() const {return fHeader ? fHeader->GetTriggerMask() : 0;}
96 UChar_t GetTriggerCluster() const {return fHeader ? fHeader->GetTriggerCluster() : 0;}
97
98 Double_t GetZDCN1Energy() const { return fHeader ? fHeader->GetZDCN1Energy() : -999.; }
99 Double_t GetZDCP1Energy() const { return fHeader ? fHeader->GetZDCP1Energy() : -999.; }
100 Double_t GetZDCN2Energy() const { return fHeader ? fHeader->GetZDCN2Energy() : -999.; }
101 Double_t GetZDCP2Energy() const { return fHeader ? fHeader->GetZDCP2Energy() : -999.; }
102 Double_t GetZDCEMEnergy(Int_t i) const { return fHeader ? fHeader->GetZDCEMEnergy(i) : -999.; }
103
104 // -- Tracks
105 TClonesArray *GetTracks() const { return fTracks; }
106 Int_t GetNTracks() const { return fTracks->GetEntriesFast(); }
107 Int_t GetNumberOfTracks() const { return GetNTracks(); }
108 AliAODTrack *GetTrack(Int_t nTrack) const { return (AliAODTrack*)fTracks->UncheckedAt(nTrack); }
109 Int_t AddTrack(const AliAODTrack* trk)
110 {new((*fTracks)[fTracks->GetEntriesFast()]) AliAODTrack(*trk); return fTracks->GetEntriesFast()-1;}
111 Int_t GetMuonTracks(TRefArray *muonTracks) const;
112
113 // -- Vertex
114 TClonesArray *GetVertices() const { return fVertices; }
115 Int_t GetNumberOfVertices() const { return fVertices->GetEntriesFast(); }
116 AliAODVertex *GetVertex(Int_t nVertex) const { return (AliAODVertex*)fVertices->UncheckedAt(nVertex); }
117 Int_t AddVertex(const AliAODVertex* vtx)
118 {new((*fVertices)[fVertices->GetEntriesFast()]) AliAODVertex(*vtx); return fVertices->GetEntriesFast()-1;}
119
120 // primary vertex
121 virtual AliAODVertex *GetPrimaryVertex() const { return GetVertex(0); }
122
123 // V0
124 TClonesArray *GetV0s() const { return fV0s; }
125 Int_t GetNumberOfV0s() const { return fV0s->GetEntriesFast(); }
126 AliAODv0 *GetV0(Int_t nV0) const { return (AliAODv0*)fV0s->UncheckedAt(nV0); }
127 Int_t AddV0(const AliAODv0* v0)
128 {new((*fV0s)[fV0s->GetEntriesFast()]) AliAODv0(*v0); return fV0s->GetEntriesFast()-1;}
129
130 // Cascades
131 TClonesArray *GetCascades() const { return fCascades; }
132 Int_t GetNumberOfCascades() const { return fCascades->GetEntriesFast(); }
133 AliAODcascade *GetCascade(Int_t nCasc) const { return (AliAODcascade*)fCascades->UncheckedAt(nCasc); }
134 Int_t AddCascade(const AliAODcascade* cascade)
135 {new((*fCascades)[fCascades->GetEntriesFast()]) AliAODcascade(*cascade); return fCascades->GetEntriesFast()-1;}
136
137 // -- EMCAL and PHOS Cluster
138 TClonesArray *GetCaloClusters() const { return fCaloClusters; }
139 Int_t GetNCaloClusters() const { return fCaloClusters->GetEntriesFast(); }
140 AliAODCaloCluster *GetCaloCluster(Int_t nCluster) const { return (AliAODCaloCluster*)fCaloClusters->UncheckedAt(nCluster); }
141 Int_t AddCaloCluster(const AliAODCaloCluster* clus)
142 {new((*fCaloClusters)[fCaloClusters->GetEntriesFast()]) AliAODCaloCluster(*clus); return fCaloClusters->GetEntriesFast()-1;}
143
144 Int_t GetEMCALClusters(TRefArray *clusters) const;
145 Int_t GetPHOSClusters(TRefArray *clusters) const;
146
147
148 // -- FMD Cluster
149 TClonesArray *GetFmdClusters() const { return fFmdClusters; }
150 Int_t GetNFmdClusters() const { return fFmdClusters->GetEntriesFast(); }
151 AliAODFmdCluster *GetFmdCluster(Int_t nCluster) const { return (AliAODFmdCluster*)fFmdClusters->UncheckedAt(nCluster); }
152 Int_t AddFmdCluster(const AliAODFmdCluster* clus)
153 {new((*fFmdClusters)[fFmdClusters->GetEntriesFast()]) AliAODFmdCluster(*clus); return fFmdClusters->GetEntriesFast()-1;}
154
155 // -- PMD Cluster
156 TClonesArray *GetPmdClusters() const { return fPmdClusters; }
157 Int_t GetNPmdClusters() const { return fPmdClusters->GetEntriesFast(); }
158 AliAODPmdCluster *GetPmdCluster(Int_t nCluster) const { return (AliAODPmdCluster*)fPmdClusters->UncheckedAt(nCluster); }
159 Int_t AddPmdCluster(const AliAODPmdCluster* clus)
160 {new((*fPmdClusters)[fPmdClusters->GetEntriesFast()]) AliAODPmdCluster(*clus); return fPmdClusters->GetEntriesFast()-1;}
161
162 // -- Jet
163 TClonesArray *GetJets() const { return fJets; }
164 Int_t GetNJets() const { return fJets->GetEntriesFast(); }
165 AliAODJet *GetJet(Int_t nJet) const { return (AliAODJet*)fJets->UncheckedAt(nJet); }
166 Int_t AddJet(const AliAODJet* vtx)
167 {new((*fJets)[fJets->GetEntriesFast()]) AliAODJet(*vtx); return fJets->GetEntriesFast()-1;}
168
169 // -- Tracklets
170 AliAODTracklets *GetTracklets() const { return fTracklets; }
171
172 // -- Calorimeter Cells
173 AliAODCaloCells *GetEMCALCells() const { return fEmcalCells; }
174 AliAODCaloCells *GetPHOSCells() const { return fPhosCells; }
175
176 // -- Services
177 void CreateStdContent();
178 void SetStdNames();
179 void GetStdContent();
180 void CreateStdFolders();
181 void ResetStd(Int_t trkArrSize = 0,
182 Int_t vtxArrSize = 0,
183 Int_t v0ArrSize = 0,
184 Int_t cascadeArrSize = 0,
185 Int_t jetSize = 0,
186 Int_t caloClusSize = 0,
187 Int_t fmdClusSize = 0,
188 Int_t pmdClusSize = 0
189 );
190 void ClearStd();
191 void ReadFromTree(TTree *tree, Option_t* opt = "");
192 void WriteToTree(TTree* tree) const {tree->Branch(fAODObjects);}
193
194 void Print(Option_t *option="") const;
195 void MakeEntriesReferencable();
196 static void AssignIDtoCollection(TCollection* col);
197 private :
198
199 TList *fAODObjects; // list of AODObjects
200 TFolder *fAODFolder; // folder structure of branches
201 Bool_t fConnected; //! flag if leaves are alreday connected
202 // standard content
203 AliAODHeader *fHeader; //! event information
204 TClonesArray *fTracks; //! charged tracks
205 TClonesArray *fVertices; //! vertices
206 TClonesArray *fV0s; //! V0s
207 TClonesArray *fCascades; //! Cascades
208 AliAODTracklets *fTracklets; //! SPD tracklets
209 TClonesArray *fJets; //! jets
210 AliAODCaloCells *fEmcalCells; //! EMCAL calorimenter cells
211 AliAODCaloCells *fPhosCells; //! PHOS calorimenter cells
212 TClonesArray *fCaloClusters; //! calorimeter clusters
213 TClonesArray *fFmdClusters; //! FMDclusters
214 TClonesArray *fPmdClusters; //! PMDclusters
215
216 static const char* fAODListName[kAODListN]; //!
217
218 ClassDef(AliAODEvent, 5);
219};
220
221#endif