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