]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliMCEvent.h
c846ee63abb401d8b58943217f881eec453bb320
[u/mrichter/AliRoot.git] / STEER / AliMCEvent.h
1 // -*- mode: C++ -*- 
2 #ifndef ALIMCEVENT_H
3 #define ALIMCEVENT_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 AliMCEvent
12 //      
13 // Origin: Andreas.Morsch, CERN, andreas.morsch@cern.ch 
14 //-------------------------------------------------------------------------
15
16
17 #include <TTree.h>
18 #include <TRefArray.h>
19 #include <TClonesArray.h>
20
21 #include <AliVEvent.h>
22 #include "AliVHeader.h"
23 #include "AliVParticle.h"
24 #include "AliMCParticle.h"
25
26 class AliStack;
27 class AliHeader;
28 class AliGenEventHeader;
29
30 class TClonesArray;
31 class TList;
32
33 class AliMCEvent : public AliVEvent {
34
35 public:
36
37     AliMCEvent();
38     virtual ~AliMCEvent() {;} 
39     AliMCEvent(const AliMCEvent& mcEvnt); 
40     AliMCEvent& operator=(const AliMCEvent& mcEvnt);
41     //
42     // Methods implementing the interface
43     //
44     // Services
45     virtual void AddObject(TObject* /*obj*/)               {;}
46     virtual TObject* FindListObject(const char */*name*/)  {return 0;}
47     virtual TList* GetList() const                         {return 0;}
48     virtual void CreateStdContent()                        {;} 
49     virtual void GetStdContent()                           {;}
50     virtual void ReadFromTree(TTree * /*tree*/, Option_t* /*opt*/) {;}
51     virtual void WriteToTree(TTree* /*tree*/)  const       {;}
52
53     virtual void SetStdNames()                             {;}
54     virtual void Print(Option_t */*option=""*/)  const     {;}
55     virtual void PreReadAll();
56     virtual void Reset()                                   {;}
57
58     // Header
59     virtual AliVHeader* GetHeader()          const         {return 0;}
60
61     // Delegated methods for fESDRun or AODHeader
62   
63     virtual void     SetRunNumber(Int_t /*n*/)             {;}
64     virtual void     SetPeriodNumber(UInt_t /*n*/)         {;}
65     virtual void     SetMagneticField(Double_t /*mf*/)     {;}
66     
67   
68     virtual Int_t    GetRunNumber()          const         {return 0;}
69     virtual UInt_t   GetPeriodNumber()       const         {return 0;}
70     virtual Double_t GetMagneticField()      const         {return 0.;}
71
72     // Setters not needed
73     virtual void      SetOrbitNumber(UInt_t /*n*/)         {;}
74     virtual void      SetBunchCrossNumber(UShort_t /*n*/)  {;}
75     virtual void      SetEventType(UInt_t /*eventType*/)   {;}
76     virtual void      SetTriggerMask(ULong64_t /*n*/)      {;}
77     virtual void      SetTriggerCluster(UChar_t /*n*/)     {;} 
78
79     virtual UInt_t    GetOrbitNumber()        const {return 0;}
80     virtual UShort_t  GetBunchCrossNumber()   const {return 0;}
81     
82     virtual UInt_t    GetEventType()          const {return 0;}
83
84     virtual ULong64_t GetTriggerMask()        const {return 0;}
85     virtual UChar_t   GetTriggerCluster()     const {return 0;}
86     virtual Double_t  GetZDCN1Energy()        const {return 0.;}
87     virtual Double_t  GetZDCP1Energy()        const {return 0.;}
88     virtual Double_t  GetZDCN2Energy()        const {return 0.;}
89     virtual Double_t  GetZDCP2Energy()        const {return 0.;}
90     virtual Double_t  GetZDCEMEnergy(Int_t /*i*/) 
91                                               const {return 0.;}
92     // Tracks
93     virtual AliVParticle *GetTrack(Int_t i) const;
94     virtual Int_t     GetNumberOfTracks()    const {return fNparticles;}
95     virtual Int_t     GetNumberOfV0s()       const {return -1;}
96     virtual Int_t     GetNumberOfCascades()  const {return -1;}
97
98     //
99     // MC Specific methods
100     //
101     // Getters
102     AliStack*    Stack()   {return fStack;}
103     AliHeader*   Header()  {return fHeader;}
104     AliGenEventHeader* GenEventHeader();
105     // Services
106     virtual void      ConnectTreeE (TTree* tree);
107     virtual void      ConnectTreeK (TTree* tree);
108     virtual void      ConnectTreeTR(TTree* tree);
109     virtual void      Clean();
110     virtual void      InitEvent();
111     virtual void      FinishEvent();
112     virtual Int_t     GetParticleAndTR(Int_t i, TParticle*& particle, TClonesArray*& trefs);
113     virtual void      DrawCheck(Int_t i, Int_t search);
114     virtual void      AddSubsidiaryEvent(AliMCEvent* event);
115     virtual Int_t     GetNumberOfPrimaries() {return fNprimaries;}
116     virtual Int_t     GetPrimaryOffset()    const {return fPrimaryOffset;}
117     virtual Int_t     GetSecondaryOffset()  const {return fSecondaryOffset;}    
118     virtual void      SetPrimaryOffset(Int_t ioff)    {fPrimaryOffset = ioff;}
119     virtual void      SetSecondaryOffset(Int_t ioff)  {fSecondaryOffset = ioff;}    
120     virtual Bool_t    IsPhysicalPrimary(Int_t i);
121     virtual Int_t     BgLabelToIndex(Int_t label);
122     static  Int_t     BgLabelOffset() {return fgkBgLabelOffset;}
123     // External particle array
124     virtual void      SetParticleArray(TClonesArray* mcParticles) 
125         {fMCParticles = mcParticles; fNparticles = fMCParticles->GetEntries(); fExternal = kTRUE;}
126     
127      
128 private:
129     virtual void      ReorderAndExpandTreeTR();
130     virtual Int_t     FindIndexAndEvent(Int_t oldidx, AliMCEvent*& event) const;
131 private: 
132     // Stanndard implementation for ESD production
133     AliStack         *fStack;            // Current pointer to stack
134     TClonesArray     *fMCParticles;      // Pointer to list of particles
135     TRefArray        *fMCParticleMap;    // Map of MC Particles
136     AliHeader        *fHeader;           // Current pointer to header
137     TClonesArray     *fTRBuffer;         // Track reference buffer    
138     TClonesArray     *fTrackReferences;  // Array of track references
139     TTree            *fTreeTR;           // Pointer to Track Reference Tree
140     TTree            *fTmpTreeTR;        // Temporary tree TR to read old format
141     TFile            *fTmpFileTR;        // Temporary file with TreeTR to read old format
142     Int_t             fNprimaries;       // Number of primaries
143     Int_t             fNparticles;       // Number of particles
144     TList            *fSubsidiaryEvents; // List of possible subsidiary events (for example merged underlying event) 
145     Int_t             fPrimaryOffset;    // Offset for primaries
146     Int_t             fSecondaryOffset;  // Offset for secondaries
147     Bool_t            fExternal;         // True if external particle array
148     static Int_t      fgkBgLabelOffset;  // Standard branch name    
149     ClassDef(AliMCEvent, 1)  // AliVEvent realisation for MC data
150 };
151
152
153 #endif 
154