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