]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDEvent.h
Next step towards supporting multiple primary vertices in ESD. Every primary vertex...
[u/mrichter/AliRoot.git] / STEER / AliESDEvent.h
1 // -*- mode: C++ -*- 
2 #ifndef ALIESDEVENT_H
3 #define ALIESDEVENT_H
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7
8 /* $Id$ */
9
10 //-------------------------------------------------------------------------
11 //                          Class AliESDEvent
12 //   This is the class to deal with during the physics analysis of data.
13 //   It also ensures the backward compatibility with the old ESD format.
14 //      
15 // Origin: Christian Klein-Boesing, CERN, Christian.Klein-Boesing@cern.ch 
16 //-------------------------------------------------------------------------
17
18 #include <TClonesArray.h>
19 #include <TObject.h>
20 #include <TTree.h>
21 #include <TArrayF.h>
22
23 class TList;
24
25 #include "AliVEvent.h"
26 // some includes for delegated methods
27 #include "AliESDCaloTrigger.h"
28 #include "AliESDRun.h"
29 #include "AliESDHeader.h"
30 #include "AliESDTZERO.h"
31 #include "AliESDZDC.h"
32 #include "AliESDACORDE.h"
33
34 // AliESDtrack has to be included so that the compiler 
35 // knows its inheritance tree (= that it is a AliVParticle).
36 #include "AliESDtrack.h"
37 // same for AliESDVertex (which is a AliVVertex)
38 #include "AliESDVertex.h"
39
40 class AliESDfriend;
41 class AliESDVZERO;
42 class AliESDHLTtrack;
43 class AliESDVertex;
44 class AliESDPmdTrack;
45 class AliESDFMD;
46 class AliESDkink;
47 class AliESDCaloCluster;
48 class AliESDCaloCells;
49 class AliESDv0;
50 class AliMultiplicity;
51 class AliRawDataErrorLog;
52 class AliESDRun;
53 class AliESDTrdTrack;
54 class AliESDMuonTrack;
55 class AliESD;
56 class AliESDcascade;
57 class TRefArray;
58 class AliESDACORDE;
59
60 class AliESDEvent : public AliVEvent {
61 public:
62
63
64   enum ESDListIndex   {kESDRun,
65                        kHeader,
66                        kESDZDC,
67                        kESDFMD,
68                        kESDVZERO,
69                        kESDTZERO,
70                        kTPCVertex,
71                        kSPDVertex,
72                        kPrimaryVertex,
73                        kSPDMult,
74                        kPHOSTrigger,
75                        kEMCALTrigger,
76                        kSPDPileupVertices,
77                        kTrkPileupVertices,
78                        kTracks,
79                        kMuonTracks,
80                        kPmdTracks,
81                        kTrdTracks,
82                        kV0s,
83                        kCascades,
84                        kKinks,
85                        kCaloClusters,
86                        kEMCALCells,
87                        kPHOSCells,
88                        kErrorLogs,
89                        kESDACORDE,
90                        kESDListN
91   };
92
93   AliESDEvent();
94   virtual ~AliESDEvent();
95   AliESDEvent &operator=(const AliESDEvent& source); // or make private and use only copy? 
96   virtual void Copy(TObject& obj) const;
97
98   // RUN
99   // move this to the UserData!!!
100   const AliESDRun*    GetESDRun() const {return fESDRun;}
101
102   // Delegated methods for fESDRun
103   void     SetRunNumber(Int_t n) {fESDRun->SetRunNumber(n);}
104   Int_t    GetRunNumber() const {return fESDRun->GetRunNumber();}
105   void     SetPeriodNumber(UInt_t n){fESDRun->SetPeriodNumber(n);}
106   UInt_t   GetPeriodNumber() const {return fESDRun->GetPeriodNumber();}
107   void     SetMagneticField(Double_t mf){fESDRun->SetMagneticField(mf);}
108   Double_t GetMagneticField() const {return fESDRun->GetMagneticField();}
109   void     SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
110   Double_t  GetDiamondX() const {return fESDRun->GetDiamondX();}
111   Double_t  GetDiamondY() const {return fESDRun->GetDiamondY();}
112   Double_t  GetSigma2DiamondX() const {return  fESDRun->GetSigma2DiamondX();}
113   Double_t  GetSigma2DiamondY() const {return  fESDRun->GetSigma2DiamondY();}
114   void      GetDiamondCovXY(Float_t cov[3]) const {fESDRun->GetDiamondCovXY(cov);}   
115   void     SetTriggerClass(const char*name, Int_t index) {fESDRun->SetTriggerClass(name,index);}
116   
117
118   // HEADER
119   AliESDHeader* GetHeader() const {return fHeader;}
120
121   // Delegated methods for fHeader
122   void      SetTriggerMask(ULong64_t n) {fHeader->SetTriggerMask(n);}
123   void      SetOrbitNumber(UInt_t n) {fHeader->SetOrbitNumber(n);}
124   void      SetTimeStamp(UInt_t timeStamp){fHeader->SetTimeStamp(timeStamp);}
125   void      SetEventType(UInt_t eventType){fHeader->SetEventType(eventType);}
126   void      SetEventNumberInFile(Int_t n) {fHeader->SetEventNumberInFile(n);}
127   //  void     SetRunNumber(Int_t n) {fHeader->SetRunNumber(n);}
128   void      SetBunchCrossNumber(UShort_t n) {fHeader->SetBunchCrossNumber(n);}
129   void      SetTriggerCluster(UChar_t n) {fHeader->SetTriggerCluster(n);}
130   
131   ULong64_t GetTriggerMask() const {return fHeader->GetTriggerMask();}
132   TString   GetFiredTriggerClasses() const {return fESDRun->GetFiredTriggerClasses(fHeader->GetTriggerMask());}
133   Bool_t    IsTriggerClassFired(const char *name) const {return fESDRun->IsTriggerClassFired(fHeader->GetTriggerMask(),name);}
134   UInt_t    GetOrbitNumber() const {return fHeader->GetOrbitNumber();}
135   UInt_t    GetTimeStamp()  const { return fHeader->GetTimeStamp();}
136   UInt_t    GetEventType()  const { return fHeader->GetEventType();}
137   Int_t     GetEventNumberInFile() const {return fHeader->GetEventNumberInFile();}
138   UShort_t  GetBunchCrossNumber() const {return fHeader->GetBunchCrossNumber();}
139   UChar_t   GetTriggerCluster() const {return fHeader->GetTriggerCluster();}
140
141   // ZDC CKB: put this in the header?
142   AliESDZDC*    GetESDZDC() const {return fESDZDC;}
143
144   // Delegated methods for fESDZDC
145   Double_t GetZDCN1Energy() const {return fESDZDC->GetZDCN1Energy();}
146   Double_t GetZDCP1Energy() const {return fESDZDC->GetZDCP1Energy();}
147   Double_t GetZDCN2Energy() const {return fESDZDC->GetZDCN2Energy();}
148   Double_t GetZDCP2Energy() const {return fESDZDC->GetZDCP2Energy();}
149   Double_t GetZDCEMEnergy(Int_t i=0) const {return fESDZDC->GetZDCEMEnergy(i);}
150   Int_t   GetZDCParticipants() const {return fESDZDC->GetZDCParticipants();}
151   Int_t   GetZDCParticipants2() const {return fESDZDC->GetZDCParticipants2();}
152   void    SetZDC(Float_t n1Energy, Float_t p1Energy, Float_t em1Energy, Float_t em2Energy,
153                  Float_t n2Energy, Float_t p2Energy, 
154                  Int_t participants, Int_t participants2)
155   {fESDZDC->SetZDC(n1Energy, p1Energy, em1Energy, em2Energy, n2Energy, p2Energy, 
156    participants, participants2);}
157
158
159   // FMD
160   void SetFMDData(AliESDFMD * obj);
161   AliESDFMD *GetFMDData(){ return fESDFMD; }
162
163
164   // TZERO CKB: put this in the header?
165   const AliESDTZERO*    GetESDTZERO() const {return fESDTZERO;}
166   // delegetated methods for fESDTZERO
167
168   Double_t GetT0zVertex() const {return fESDTZERO->GetT0zVertex();}
169   void SetT0zVertex(Float_t z) {fESDTZERO->SetT0zVertex(z);}
170   Double_t GetT0() const {return fESDTZERO->GetT0();}
171   void SetT0(Float_t timeStart) {fESDTZERO->SetT0(timeStart);}
172   const Double_t * GetT0time() const {return fESDTZERO->GetT0time();}
173   void SetT0time(Float_t time[24]) {fESDTZERO->SetT0time(time);}
174   const Double_t * GetT0amplitude() const {return fESDTZERO->GetT0amplitude();}
175   void SetT0amplitude(Float_t amp[24]){fESDTZERO->SetT0amplitude(amp);}
176
177   // VZERO 
178   AliESDVZERO *GetVZEROData() const { return fESDVZERO; }
179   void SetVZEROData(AliESDVZERO * obj);
180
181  // ACORDE
182   AliESDACORDE *GetACORDEData() const { return fESDACORDE;}
183   void SetACORDEData(AliESDACORDE * obj);
184
185   void SetESDfriend(const AliESDfriend *f) const;
186   void GetESDfriend(AliESDfriend *f) const;
187
188
189
190   void SetPrimaryVertexTPC(const AliESDVertex *vertex); 
191   const AliESDVertex *GetPrimaryVertexTPC() const {return fTPCVertex;}
192
193   void SetPrimaryVertexSPD(const AliESDVertex *vertex); 
194   const AliESDVertex *GetPrimaryVertexSPD() const {return fSPDVertex;}
195   const AliESDVertex *GetVertex() const {
196     //For the backward compatibily only
197      return GetPrimaryVertexSPD();
198   }
199
200   void SetPrimaryVertexTracks(const AliESDVertex *vertex);
201   const AliESDVertex *GetPrimaryVertexTracks() const {return fPrimaryVertex;}
202
203   const AliESDVertex *GetPrimaryVertex() const;
204
205   void SetMultiplicity(const AliMultiplicity *mul);
206
207   const AliMultiplicity *GetMultiplicity() const {return fSPDMult;}
208
209
210   Bool_t Clean(Float_t *cleanPars);
211   Bool_t RemoveKink(Int_t i)   const;
212   Bool_t RemoveV0(Int_t i)     const;
213   Bool_t RemoveTrack(Int_t i)  const;
214
215   const AliESDVertex *GetPileupVertexSPD(Int_t i) const {
216     return (const AliESDVertex *)fSPDPileupVertices->UncheckedAt(i);
217   }
218   Char_t  AddPileupVertexSPD(const AliESDVertex *vtx);
219
220   const AliESDVertex *GetPileupVertexTracks(Int_t i) const {
221     return (const AliESDVertex *)fTrkPileupVertices->UncheckedAt(i);
222   }
223   Char_t  AddPileupVertexTracks(const AliESDVertex *vtx);
224
225   AliESDtrack *GetTrack(Int_t i) const {
226     return (AliESDtrack *)fTracks->UncheckedAt(i);
227   }
228   Int_t  AddTrack(const AliESDtrack *t);
229
230   
231   AliESDHLTtrack *GetHLTConfMapTrack(Int_t /*i*/) const {
232     //    return (AliESDHLTtrack *)fHLTConfMapTracks->UncheckedAt(i);
233     return 0;
234   }
235   void AddHLTConfMapTrack(const AliESDHLTtrack */*t*/) {
236     printf("ESD:: AddHLTConfMapTrack do nothing \n");
237     //    TClonesArray &fhlt = *fHLTConfMapTracks;
238     //  new(fhlt[fHLTConfMapTracks->GetEntriesFast()]) AliESDHLTtrack(*t);
239   }
240   
241
242   AliESDHLTtrack *GetHLTHoughTrack(Int_t /*i*/) const {
243     //    return (AliESDHLTtrack *)fHLTHoughTracks->UncheckedAt(i);
244     return 0;
245   }
246   void AddHLTHoughTrack(const AliESDHLTtrack */*t*/) {
247     printf("ESD:: AddHLTHoughTrack do nothing \n");
248     //    TClonesArray &fhlt = *fHLTHoughTracks;
249     //     new(fhlt[fHLTHoughTracks->GetEntriesFast()]) AliESDHLTtrack(*t);
250   }
251   
252   AliESDMuonTrack *GetMuonTrack(Int_t i) const {
253     return (AliESDMuonTrack *)fMuonTracks->UncheckedAt(i);
254   }
255
256   void AddMuonTrack(const AliESDMuonTrack *t);
257
258   AliESDPmdTrack *GetPmdTrack(Int_t i) const {
259     return (AliESDPmdTrack *)fPmdTracks->UncheckedAt(i);
260   }
261
262   void AddPmdTrack(const AliESDPmdTrack *t);
263
264
265   AliESDTrdTrack *GetTrdTrack(Int_t i) const {
266     return (AliESDTrdTrack *)fTrdTracks->UncheckedAt(i);
267   }
268
269   
270   void AddTrdTrack(const AliESDTrdTrack *t);
271
272   AliESDv0 *GetV0(Int_t i) const {
273     return (AliESDv0*)fV0s->UncheckedAt(i);
274   }
275   Int_t AddV0(const AliESDv0 *v);
276
277   AliESDcascade *GetCascade(Int_t i) const {
278     return (AliESDcascade *)fCascades->UncheckedAt(i);
279   }
280
281   void AddCascade(const AliESDcascade *c);
282
283   AliESDkink *GetKink(Int_t i) const {
284     return (AliESDkink *)fKinks->UncheckedAt(i);
285   }
286   Int_t AddKink(const AliESDkink *c);
287
288   AliESDCaloCluster *GetCaloCluster(Int_t i) const {
289     return (AliESDCaloCluster *)fCaloClusters->UncheckedAt(i);
290   }
291
292   Int_t AddCaloCluster(const AliESDCaloCluster *c);
293
294   AliESDCaloCells *GetEMCALCells() const {return fEMCALCells; }  
295   AliESDCaloCells *GetPHOSCells() const {return fPHOSCells; }  
296
297   AliRawDataErrorLog *GetErrorLog(Int_t i) const {
298     return (AliRawDataErrorLog *)fErrorLogs->UncheckedAt(i);
299   }
300   void  AddRawDataErrorLog(const AliRawDataErrorLog *log) const;
301
302   Int_t GetNumberOfErrorLogs()   const {return fErrorLogs->GetEntriesFast();}
303
304     
305   void AddPHOSTriggerPosition(TArrayF array)   { fPHOSTrigger->AddTriggerPosition(array); }
306   void AddPHOSTriggerAmplitudes(TArrayF array) { fPHOSTrigger->AddTriggerAmplitudes(array);}
307   void AddEMCALTriggerPosition(TArrayF array)  { fEMCALTrigger->AddTriggerPosition(array); }
308   void AddEMCALTriggerAmplitudes(TArrayF array){ fEMCALTrigger->AddTriggerAmplitudes(array); }
309
310   Int_t GetNumberOfPileupVerticesSPD() const {
311      return fSPDPileupVertices->GetEntriesFast();
312   }
313   Int_t GetNumberOfPileupVerticesTracks() const {
314      return fTrkPileupVertices->GetEntriesFast();
315   }
316   Int_t GetNumberOfTracks()     const {return fTracks->GetEntriesFast();}
317   Int_t GetNumberOfHLTConfMapTracks()     const {return 0;} 
318   // fHLTConfMapTracks->GetEntriesFast();}
319   Int_t GetNumberOfHLTHoughTracks()     const {return  0;  }
320   //  fHLTHoughTracks->GetEntriesFast();  }
321
322   Int_t GetNumberOfMuonTracks() const {return fMuonTracks->GetEntriesFast();}
323   Int_t GetNumberOfPmdTracks() const {return fPmdTracks->GetEntriesFast();}
324   Int_t GetNumberOfTrdTracks() const {return fTrdTracks->GetEntriesFast();}
325   Int_t GetNumberOfV0s()      const {return fV0s->GetEntriesFast();}
326   Int_t GetNumberOfCascades() const {return fCascades->GetEntriesFast();}
327   Int_t GetNumberOfKinks() const {return fKinks->GetEntriesFast();}
328   
329   Int_t GetEMCALClusters(TRefArray *clusters) const;
330   Int_t GetPHOSClusters(TRefArray *clusters) const;
331   Int_t GetNumberOfCaloClusters() const {return fCaloClusters->GetEntriesFast();}
332
333   void SetUseOwnList(Bool_t b){fUseOwnList = b;}
334   Bool_t GetUseOwnList(){return fUseOwnList;}
335   
336   // Remove this stuff CKB?
337   //---------------------------------------------------
338   Int_t GetNumberOfEMCALClusters() const {return fEMCALClusters;}
339   void  SetNumberOfEMCALClusters(Int_t clus) {fEMCALClusters = clus;}
340   Int_t GetFirstEMCALCluster() const {return fFirstEMCALCluster;}
341   void  SetFirstEMCALCluster(Int_t index) {fFirstEMCALCluster = index;}
342  
343   Int_t GetNumberOfPHOSClusters() const {return fPHOSClusters;}
344   void  SetNumberOfPHOSClusters(Int_t part) { fPHOSClusters = part ; }
345   void  SetFirstPHOSCluster(Int_t index) { fFirstPHOSCluster = index ; } 
346   Int_t GetFirstPHOSCluster() const  { return fFirstPHOSCluster ; }
347   //-------------------------------------------------------
348
349   TArrayF *GetEMCALTriggerPosition() const {return  fEMCALTrigger->GetTriggerPosition();}
350   TArrayF *GetEMCALTriggerAmplitudes() const {return  fEMCALTrigger->GetTriggerAmplitudes();}
351   TArrayF *GetPHOSTriggerPosition() const {return  fPHOSTrigger->GetTriggerPosition();}
352   TArrayF *GetPHOSTriggerAmplitudes() const {return  fPHOSTrigger->GetTriggerAmplitudes();}
353
354   void ResetV0s() { fV0s->Clear(); }
355   void ResetCascades() { fCascades->Clear(); }
356   void Reset();
357
358   void  Print(Option_t *option="") const;
359
360   void AddObject(TObject* obj);
361   void ReadFromTree(TTree *tree, Option_t* opt = "");
362   TObject* FindListObject(const char *name);
363   AliESD *GetAliESDOld(){return fESDOld;}
364   void WriteToTree(TTree* tree) const;
365   void GetStdContent();
366   void ResetStdContent();
367   void CreateStdContent();
368   void CreateStdContent(Bool_t bUseThisList);
369   void SetStdNames();
370   void CopyFromOldESD();
371   TList* GetList() const {return fESDObjects;}
372
373 protected:
374   AliESDEvent(const AliESDEvent&);
375   static Bool_t ResetWithPlacementNew(TObject *pObject);
376
377   TList *fESDObjects;             // List of esd Objects
378
379   AliESDRun       *fESDRun;           //! Run information tmp put in the Userdata
380   AliESDHeader    *fHeader;           //! ESD Event Header
381   AliESDZDC       *fESDZDC;           //! ZDC information
382   AliESDFMD       *fESDFMD;           //! FMD object containing rough multiplicity
383   AliESDVZERO     *fESDVZERO;         //! VZERO object containing rough multiplicity
384   AliESDTZERO     *fESDTZERO;         //! TZEROObject
385   AliESDVertex    *fTPCVertex;        //! Primary vertex estimated by the TPC
386   AliESDVertex    *fSPDVertex;        //! Primary vertex estimated by the SPD
387   AliESDVertex    *fPrimaryVertex;    //! Primary vertex estimated using ESD tracks
388   AliMultiplicity *fSPDMult;          //! SPD tracklet multiplicity
389   AliESDCaloTrigger* fPHOSTrigger;     //! PHOS Trigger information
390   AliESDCaloTrigger* fEMCALTrigger;    //! PHOS Trigger information
391   AliESDACORDE    *fESDACORDE;        //! ACORDE ESD object caontaining bit pattern
392
393   TClonesArray *fSPDPileupVertices;//! Pileup primary vertices reconstructed by SPD 
394   TClonesArray *fTrkPileupVertices;//! Pileup primary vertices reconstructed using the tracks 
395   TClonesArray *fTracks;           //! ESD tracks 
396   TClonesArray *fMuonTracks;       //! MUON ESD tracks
397   TClonesArray *fPmdTracks;        //! PMD ESD tracks
398   TClonesArray *fTrdTracks;        //! TRD ESD tracks (triggered)
399   TClonesArray *fV0s;              //! V0 vertices
400   TClonesArray *fCascades;         //! Cascade vertices
401   TClonesArray *fKinks;            //! Kinks
402   TClonesArray *fCaloClusters;     //! Calorimeter clusters for PHOS/EMCAL
403   AliESDCaloCells *fEMCALCells;     //! EMCAL cell info
404   AliESDCaloCells *fPHOSCells;     //! PHOS cell info
405   TClonesArray *fErrorLogs;        //! Raw-data reading error messages
406  
407
408
409   AliESD       *fESDOld;           //! Old esd Structure
410   AliESDfriend *fESDFriendOld;     //! Old friend esd Structure
411   Bool_t    fConnected;            //! flag if leaves are alreday connected
412   Bool_t    fUseOwnList;           //! Do not use the list from the esdTree but use the one created by this class 
413
414   static const char* fgkESDListName[kESDListN]; //!
415
416   // Remove this stuff CKB
417   Int_t        fEMCALClusters;   // Number of EMCAL clusters (subset of caloclusters)
418   Int_t        fFirstEMCALCluster; // First EMCAL cluster in the fCaloClusters list 
419
420   Int_t        fPHOSClusters;     // Number of PHOS clusters (subset of caloclusters)
421   Int_t        fFirstPHOSCluster; // First PHOS cluster in the fCaloClusters list 
422
423   ClassDef(AliESDEvent,10)  //ESDEvent class 
424 };
425 #endif 
426