]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAODEvent.h
M AliTPCcalibDButil.cxx - Protection agains missing Offset calibration...
[u/mrichter/AliRoot.git] / STEER / 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
34 class TTree;
35 class TFolder;
36
37 class AliAODEvent : public AliVEvent {
38
39  public :
40   enum AODListIndex_t {kAODHeader,
41                        kAODTracks,
42                        kAODVertices,
43                        kAODv0,
44                        kAODcascade,
45                        kAODTracklets,
46                        kAODJets,
47                        kAODEmcalCells,
48                        kAODPhosCells,
49                        kAODCaloClusters,
50                        kAODFmdClusters,
51                        kAODPmdClusters,
52                        kAODDimuons,
53                        kAODListN
54   };
55
56   AliAODEvent();
57   virtual ~AliAODEvent();
58
59   AliAODEvent(const AliAODEvent& aodevent);
60   AliAODEvent& operator=(const AliAODEvent& aodevent);
61
62   void          AddObject(TObject *obj);
63   void          RemoveObject(TObject *obj);
64   TObject      *FindListObject(const char *objName) ;
65   TList        *GetList()                const { return fAODObjects; }
66
67   // -- Header
68   AliAODHeader *GetHeader()              const { return fHeader; }
69   void          AddHeader(const AliAODHeader* hdx)
70     {
71         delete fHeader; fHeader = new AliAODHeader(*hdx);
72         (fAODObjects->FirstLink())->SetObject(fHeader);
73     }
74
75   // setters and getters for header information
76   void     SetRunNumber(Int_t n) {if (fHeader) fHeader->SetRunNumber(n);}
77   void     SetPeriodNumber(UInt_t n){if (fHeader) fHeader->SetPeriodNumber(n);}
78   void     SetOrbitNumber(UInt_t n) {if (fHeader) fHeader->SetOrbitNumber(n);}
79   void     SetBunchCrossNumber(UShort_t n) {if (fHeader) fHeader->SetBunchCrossNumber(n);}
80   void     SetMagneticField(Double_t mf){if (fHeader) fHeader->SetMagneticField(mf);}
81   void     SetMuonMagFieldScale(Double_t mf){if (fHeader) fHeader->SetMuonMagFieldScale(mf);}
82   void     SetDiamond(Float_t xy[2],Float_t cov[3]){if (fHeader) fHeader->SetDiamond(xy,cov);}
83
84   Int_t    GetRunNumber() const {return fHeader ? fHeader->GetRunNumber() : -999;}
85   UInt_t   GetPeriodNumber() const {return fHeader ? fHeader->GetPeriodNumber() : 0;}
86   UInt_t   GetOrbitNumber() const {return fHeader ? fHeader->GetOrbitNumber() : 0;}
87   UShort_t GetBunchCrossNumber() const {return fHeader ? fHeader->GetBunchCrossNumber() : 0;}
88   Double_t GetMagneticField() const {return fHeader ? fHeader->GetMagneticField() : -999.;}
89   Double_t GetMuonMagFieldScale() const {return fHeader ? fHeader->GetMuonMagFieldScale() : -999.;}
90   Double_t GetDiamondX() const {return fHeader ? fHeader->GetDiamondX() : -999.;}
91   Double_t GetDiamondY() const {return fHeader ? fHeader->GetDiamondY() : -999.;}
92   void     GetDiamondCovXY(Float_t cov[3]) const {cov[0]=-999.; if(fHeader) fHeader->GetDiamondCovXY(cov);}
93   
94   void      SetEventType(UInt_t eventType){fHeader->SetEventType(eventType);}
95   void      SetTriggerMask(ULong64_t n) {fHeader->SetTriggerMask(n);}
96   void      SetTriggerCluster(UChar_t n) {fHeader->SetTriggerCluster(n);}
97
98   UInt_t    GetEventType()          const { return fHeader ? fHeader->GetEventType() : 0;}
99   ULong64_t GetTriggerMask()        const { return fHeader ? fHeader->GetTriggerMask() : 0;}
100   UChar_t   GetTriggerCluster()     const { return fHeader ? fHeader->GetTriggerCluster() : 0;}
101   TString   GetFiredTriggerClasses()const { return fHeader->GetFiredTriggerClasses();};
102   Double_t  GetZDCN1Energy()        const { return fHeader ? fHeader->GetZDCN1Energy() : -999.; }
103   Double_t  GetZDCP1Energy()        const { return fHeader ? fHeader->GetZDCP1Energy() : -999.; }
104   Double_t  GetZDCN2Energy()        const { return fHeader ? fHeader->GetZDCN2Energy() : -999.; }
105   Double_t  GetZDCP2Energy()        const { return fHeader ? fHeader->GetZDCP2Energy() : -999.; }
106   Double_t  GetZDCEMEnergy(Int_t i) const { return fHeader ? fHeader->GetZDCEMEnergy(i) : -999.; }
107
108   // -- Tracks
109   TClonesArray *GetTracks()              const { return fTracks; }
110   Int_t         GetNTracks()             const { return fTracks->GetEntriesFast(); }
111   Int_t         GetNumberOfTracks()      const { return GetNTracks(); }
112   AliAODTrack  *GetTrack(Int_t nTrack)   const { return (AliAODTrack*)fTracks->UncheckedAt(nTrack); }
113   Int_t         AddTrack(const AliAODTrack* trk)
114     {new((*fTracks)[fTracks->GetEntriesFast()]) AliAODTrack(*trk); return fTracks->GetEntriesFast()-1;}
115   Int_t         GetMuonTracks(TRefArray *muonTracks) const;
116
117   // -- Vertex
118   TClonesArray *GetVertices()            const { return fVertices; }
119   Int_t         GetNumberOfVertices()    const { return fVertices->GetEntriesFast(); }
120   AliAODVertex *GetVertex(Int_t nVertex) const { return (AliAODVertex*)fVertices->UncheckedAt(nVertex); }
121   Int_t         AddVertex(const AliAODVertex* vtx)
122   {new((*fVertices)[fVertices->GetEntriesFast()]) AliAODVertex(*vtx); return fVertices->GetEntriesFast()-1;}
123   
124   // primary vertex
125   virtual AliAODVertex *GetPrimaryVertex() const { return GetVertex(0); }
126
127   // V0
128   TClonesArray *GetV0s()                 const { return fV0s; }
129   Int_t         GetNumberOfV0s()         const { return fV0s->GetEntriesFast(); }
130   AliAODv0     *GetV0(Int_t nV0)         const { return (AliAODv0*)fV0s->UncheckedAt(nV0); }
131   Int_t         AddV0(const AliAODv0* v0)
132   {new((*fV0s)[fV0s->GetEntriesFast()]) AliAODv0(*v0); return fV0s->GetEntriesFast()-1;}
133
134   // Cascades
135   TClonesArray  *GetCascades()            const { return fCascades; }
136   Int_t          GetNumberOfCascades()    const { return fCascades->GetEntriesFast(); }
137   AliAODcascade *GetCascade(Int_t nCasc)  const { return (AliAODcascade*)fCascades->UncheckedAt(nCasc); }
138   Int_t          AddCascade(const AliAODcascade* cascade)
139   {new((*fCascades)[fCascades->GetEntriesFast()]) AliAODcascade(*cascade); return fCascades->GetEntriesFast()-1;}
140
141   // -- EMCAL and PHOS Cluster
142   TClonesArray *GetCaloClusters()        const { return fCaloClusters; }
143   Int_t         GetNCaloClusters()       const { return fCaloClusters->GetEntriesFast(); }
144   AliAODCaloCluster *GetCaloCluster(Int_t nCluster) const { return (AliAODCaloCluster*)fCaloClusters->UncheckedAt(nCluster); }
145   Int_t         AddCaloCluster(const AliAODCaloCluster* clus)
146   {new((*fCaloClusters)[fCaloClusters->GetEntriesFast()]) AliAODCaloCluster(*clus); return fCaloClusters->GetEntriesFast()-1;}
147
148   Int_t GetEMCALClusters(TRefArray *clusters) const;
149   Int_t GetPHOSClusters(TRefArray *clusters) const;
150
151
152   // -- FMD Cluster
153   TClonesArray *GetFmdClusters()        const { return fFmdClusters; }
154   Int_t         GetNFmdClusters()       const { return fFmdClusters->GetEntriesFast(); }
155   AliAODFmdCluster *GetFmdCluster(Int_t nCluster) const { return (AliAODFmdCluster*)fFmdClusters->UncheckedAt(nCluster); }
156   Int_t         AddFmdCluster(const AliAODFmdCluster* clus)
157   {new((*fFmdClusters)[fFmdClusters->GetEntriesFast()]) AliAODFmdCluster(*clus); return fFmdClusters->GetEntriesFast()-1;}
158
159   // -- PMD Cluster
160   TClonesArray *GetPmdClusters()        const { return fPmdClusters; }
161   Int_t         GetNPmdClusters()       const { return fPmdClusters->GetEntriesFast(); }
162   AliAODPmdCluster *GetPmdCluster(Int_t nCluster) const { return (AliAODPmdCluster*)fPmdClusters->UncheckedAt(nCluster); }
163   Int_t         AddPmdCluster(const AliAODPmdCluster* clus)
164   {new((*fPmdClusters)[fPmdClusters->GetEntriesFast()]) AliAODPmdCluster(*clus); return fPmdClusters->GetEntriesFast()-1;}
165
166   // -- Jet
167   TClonesArray *GetJets()            const { return fJets; }
168   Int_t         GetNJets()           const { return fJets->GetEntriesFast(); }
169   AliAODJet    *GetJet(Int_t nJet) const { return (AliAODJet*)fJets->UncheckedAt(nJet); }
170   Int_t         AddJet(const AliAODJet* vtx)
171     {new((*fJets)[fJets->GetEntriesFast()]) AliAODJet(*vtx); return fJets->GetEntriesFast()-1;}
172
173   // -- Tracklets
174   AliAODTracklets *GetTracklets() const { return fTracklets; }
175
176   // -- Calorimeter Cells
177   AliAODCaloCells *GetEMCALCells() const { return fEmcalCells; }
178   AliAODCaloCells *GetPHOSCells() const { return fPhosCells; }
179
180   // -- Dimuons
181   TClonesArray *GetDimuons()              const { return fDimuons; }
182   Int_t         GetNDimuons()             const { return fDimuons->GetEntriesFast(); }
183   Int_t         GetNumberOfDimuons()      const { return GetNDimuons(); }
184   AliAODDimuon *GetDimuon(Int_t nDimu)    const { return (AliAODDimuon*)fDimuons->UncheckedAt(nDimu); }
185   Int_t         AddDimuon(const AliAODDimuon* dimu)
186     {new((*fDimuons)[fDimuons->GetEntriesFast()]) AliAODDimuon(*dimu); return fDimuons->GetEntriesFast()-1;}
187   
188   // -- Services
189   void    CreateStdContent();
190   void    SetStdNames();
191   void    GetStdContent();
192   void    CreateStdFolders();
193   void    ResetStd(Int_t trkArrSize = 0, 
194                    Int_t vtxArrSize = 0, 
195                    Int_t v0ArrSize = 0, 
196                    Int_t cascadeArrSize = 0,
197                    Int_t jetSize = 0, 
198                    Int_t caloClusSize = 0, 
199                    Int_t fmdClusSize = 0, 
200                    Int_t pmdClusSize = 0,
201                    Int_t dimuonArrsize =0
202                    );
203   void    ClearStd();
204   void    Reset() {ClearStd();} 
205   void    ReadFromTree(TTree *tree, Option_t* opt = "");
206   void    WriteToTree(TTree* tree) const {tree->Branch(fAODObjects);}
207
208   void  Print(Option_t *option="") const;
209   void  MakeEntriesReferencable();
210   static void AssignIDtoCollection(TCollection* col);
211  private :
212
213   TList   *fAODObjects; //  list of AODObjects
214   TFolder *fAODFolder;  //  folder structure of branches
215   Bool_t   fConnected;  //! flag if leaves are alreday connected 
216   // standard content
217   AliAODHeader    *fHeader;       //! event information
218   TClonesArray    *fTracks;       //! charged tracks
219   TClonesArray    *fVertices;     //! vertices
220   TClonesArray    *fV0s;          //! V0s
221   TClonesArray    *fCascades;     //! Cascades
222   AliAODTracklets *fTracklets;    //! SPD tracklets
223   TClonesArray    *fJets;         //! jets
224   AliAODCaloCells *fEmcalCells;   //! EMCAL calorimenter cells
225   AliAODCaloCells *fPhosCells;    //! PHOS calorimenter cells
226   TClonesArray    *fCaloClusters; //! calorimeter clusters
227   TClonesArray    *fFmdClusters;  //! FMDclusters
228   TClonesArray    *fPmdClusters;  //! PMDclusters
229   TClonesArray    *fDimuons;      //! dimuons
230
231   static const char* fAODListName[kAODListN]; //!
232
233   ClassDef(AliAODEvent, 6);
234 };
235
236 #endif