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