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