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