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