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