]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliMCEvent.cxx
AliMCEvent: Access to Transport MC Truth during analysis
[u/mrichter/AliRoot.git] / STEER / AliMCEvent.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17 //---------------------------------------------------------------------------------
18 //                          Class AliMCEvent
19 // This class gives access to MC truth during the analysis.
20 // Monte Carlo truth is containe in the kinematics tree (produced particles) and 
21 // the tree of reference hits.
22 //      
23 // Origin: Andreas Morsch, CERN, andreas.morsch@cern.ch 
24 //---------------------------------------------------------------------------------
25
26
27
28 #include "AliMCEvent.h"
29 #include "AliTrackReference.h"
30 #include "AliHeader.h"
31 #include "AliStack.h"
32
33 #include <TTree.h>
34 #include <TFile.h>
35 #include <TParticle.h>
36 #include <TClonesArray.h>
37 #include <TDirectoryFile.h>
38 #include <TArrow.h>
39 #include <TMarker.h>
40 #include <TH2F.h>
41
42
43 ClassImp(AliMCEvent)
44
45 AliMCEvent::AliMCEvent() :
46     AliVirtualEventHandler(),
47     fFileE(0),
48     fFileK(0),
49     fFileTR(0),
50     fTreeE(0),
51     fTreeK(0),
52     fTreeTR(0),
53     fStack(0),
54     fHeader(0),
55     fTrackReferences(0),
56     fNEvent(-1),
57     fEvent(-1),
58     fNprimaries(-1),
59     fNparticles(-1)
60 {
61     // Default constructor
62 }
63
64 AliMCEvent::AliMCEvent(const char* name, const char* title) :
65     AliVirtualEventHandler(name, title),
66     fFileE(0),
67     fFileK(0),
68     fFileTR(0),
69     fTreeE(0),
70     fTreeK(0),
71     fTreeTR(0),
72     fStack(0),
73     fHeader(new AliHeader()),
74     fTrackReferences(new TClonesArray("AliTrackReference", 200)),
75     fNEvent(-1),
76     fEvent(-1),
77     fNprimaries(-1),
78     fNparticles(-1)
79 {
80     // Constructor
81 }
82 AliMCEvent::~AliMCEvent()
83
84     // Destructor
85     delete fFileE;
86     delete fFileK;
87     delete fFileTR;
88
89     delete fTreeE;
90     delete fTreeK;
91     delete fTreeTR;
92     
93 }
94
95 Bool_t AliMCEvent::InitIO(Option_t* /*opt*/) 
96
97     // Initialize input
98     fFileE = new TFile("galice.root");
99     fFileE->GetObject("TE", fTreeE);
100     fTreeE->SetBranchAddress("Header", &fHeader);
101     fNEvent = fTreeE->GetEntries();
102     //
103     // Tree K
104     fFileK = new TFile("Kinematics.root");
105     //
106     // Tree TR
107     fFileTR = new TFile("TrackRefs.root");
108     //
109     // Reset the event number
110     fEvent = -1;
111     printf("Number of events in this directory %5d \n", fNEvent);
112     return kTRUE;
113     
114 }
115
116 Bool_t AliMCEvent::BeginEvent()
117
118     // Read the next event
119     fEvent++;
120     char folder[20];
121     sprintf(folder, "Event%d", fEvent);
122     // TreeE
123     fTreeE->GetEntry(fEvent);
124     fStack = fHeader->Stack();
125     // Tree K
126     TDirectoryFile* dirK  = 0;
127     fFileK->GetObject(folder, dirK);
128     dirK->GetObject("TreeK", fTreeK);
129     fStack->ConnectTree(fTreeK);
130     fStack->GetEvent();
131     //Tree TR 
132     TDirectoryFile* dirTR = 0;
133     fFileTR->GetObject(folder, dirTR);
134     dirTR->GetObject("TreeTR", fTreeTR);
135     fTreeTR->SetBranchAddress("TrackReferences", &fTrackReferences);
136     //
137     fNparticles = fStack->GetNtrack();
138     fNprimaries = fStack->GetNprimary();
139     printf("Event#: %5d Number of particles %5d \n", fEvent, fNparticles);
140     
141     return kTRUE;
142 }
143
144 Int_t AliMCEvent::GetParticleAndTR(Int_t i, TParticle*& particle, TClonesArray*& trefs)
145 {
146     // Retrieve entry i
147     if (i > -1 && i < fNparticles) {
148         fTreeTR->GetEntry(fStack->TreeKEntry(i));
149     } else {
150         printf("AliMCEvent::GetEntry: Index out of range");
151     }
152     particle = fStack->Particle(i);
153     trefs    = fTrackReferences;
154     printf("Returning %5d \n", particle->GetPdgCode());
155     
156     return trefs->GetEntries();
157     
158 }
159
160 void AliMCEvent::DrawCheck(Int_t i)
161 {
162     // Retrieve entry i and draw momentum vector and hits
163     if (i > -1 && i < fNparticles) {
164         fTreeTR->GetEntry(fStack->TreeKEntry(i));
165     } else {
166         printf("AliMCEvent::GetEntry: Index out of range");
167     }
168
169     TParticle* particle = fStack->Particle(i);
170     Int_t nh =  fTrackReferences->GetEntries();
171     
172     TH2F*    h = new TH2F("", "", 100, -500, 500, 100, -500, 500);
173     Float_t x0 = particle->Vx();
174     Float_t y0 = particle->Vy();
175
176     Float_t x1 = particle->Vx() + particle->Px() * 50.;
177     Float_t y1 = particle->Vy() + particle->Py() * 50.;
178     
179     TArrow*  a = new TArrow(x0, y0, x1, y1, 0.01);
180     h->Draw();
181     a->SetLineColor(2);
182     
183     a->Draw();
184     
185     for (Int_t ih = 0; ih < nh; ih++) {
186         AliTrackReference* ref = (AliTrackReference*) fTrackReferences->At(ih);
187         TMarker* m = new TMarker(ref->X(), ref->Y(), 20);
188         m->Draw();
189         m->SetMarkerSize(0.4);
190         
191     }
192 }
193
194     
195
196                             
197 Bool_t AliMCEvent::FinishEvent()
198 {
199     // Dummy 
200     return kTRUE;
201 }
202
203 Bool_t AliMCEvent::Terminate()
204
205     // Dummy 
206     return kTRUE;
207 }
208
209 Bool_t AliMCEvent::TerminateIO()
210
211     // Dummy
212     return kTRUE;
213 }
214