1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 //---------------------------------------------------------------------------------
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.
23 // Origin: Andreas Morsch, CERN, andreas.morsch@cern.ch
24 //---------------------------------------------------------------------------------
28 #include "AliMCEvent.h"
29 #include "AliTrackReference.h"
30 #include "AliHeader.h"
35 #include <TParticle.h>
36 #include <TClonesArray.h>
37 #include <TDirectoryFile.h>
45 AliMCEvent::AliMCEvent() :
46 AliVirtualEventHandler(),
61 // Default constructor
64 AliMCEvent::AliMCEvent(const char* name, const char* title) :
65 AliVirtualEventHandler(name, title),
73 fHeader(new AliHeader()),
74 fTrackReferences(new TClonesArray("AliTrackReference", 200)),
82 AliMCEvent::~AliMCEvent()
95 Bool_t AliMCEvent::InitIO(Option_t* /*opt*/)
98 fFileE = new TFile("galice.root");
99 fFileE->GetObject("TE", fTreeE);
100 fTreeE->SetBranchAddress("Header", &fHeader);
101 fNEvent = fTreeE->GetEntries();
104 fFileK = new TFile("Kinematics.root");
107 fFileTR = new TFile("TrackRefs.root");
109 // Reset the event number
111 printf("Number of events in this directory %5d \n", fNEvent);
116 Bool_t AliMCEvent::BeginEvent()
118 // Read the next event
121 sprintf(folder, "Event%d", fEvent);
123 fTreeE->GetEntry(fEvent);
124 fStack = fHeader->Stack();
126 TDirectoryFile* dirK = 0;
127 fFileK->GetObject(folder, dirK);
128 dirK->GetObject("TreeK", fTreeK);
129 fStack->ConnectTree(fTreeK);
132 TDirectoryFile* dirTR = 0;
133 fFileTR->GetObject(folder, dirTR);
134 dirTR->GetObject("TreeTR", fTreeTR);
135 fTreeTR->SetBranchAddress("TrackReferences", &fTrackReferences);
137 fNparticles = fStack->GetNtrack();
138 fNprimaries = fStack->GetNprimary();
139 printf("Event#: %5d Number of particles %5d \n", fEvent, fNparticles);
144 Int_t AliMCEvent::GetParticleAndTR(Int_t i, TParticle*& particle, TClonesArray*& trefs)
147 if (i > -1 && i < fNparticles) {
148 fTreeTR->GetEntry(fStack->TreeKEntry(i));
150 printf("AliMCEvent::GetEntry: Index out of range");
152 particle = fStack->Particle(i);
153 trefs = fTrackReferences;
154 printf("Returning %5d \n", particle->GetPdgCode());
156 return trefs->GetEntries();
160 void AliMCEvent::DrawCheck(Int_t i)
162 // Retrieve entry i and draw momentum vector and hits
163 if (i > -1 && i < fNparticles) {
164 fTreeTR->GetEntry(fStack->TreeKEntry(i));
166 printf("AliMCEvent::GetEntry: Index out of range");
169 TParticle* particle = fStack->Particle(i);
170 Int_t nh = fTrackReferences->GetEntries();
172 TH2F* h = new TH2F("", "", 100, -500, 500, 100, -500, 500);
173 Float_t x0 = particle->Vx();
174 Float_t y0 = particle->Vy();
176 Float_t x1 = particle->Vx() + particle->Px() * 50.;
177 Float_t y1 = particle->Vy() + particle->Py() * 50.;
179 TArrow* a = new TArrow(x0, y0, x1, y1, 0.01);
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);
189 m->SetMarkerSize(0.4);
197 Bool_t AliMCEvent::FinishEvent()
203 Bool_t AliMCEvent::Terminate()
209 Bool_t AliMCEvent::TerminateIO()