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