]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenReaderEcalHijing.cxx
Information about trigger jets added.
[u/mrichter/AliRoot.git] / EVGEN / AliGenReaderEcalHijing.cxx
CommitLineData
b152a071 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/*
17$Log$
58323221 18Revision 1.1 2002/01/08 09:59:34 morsch
19Readers for EMCAL primary particle input.
20
b152a071 21*/
22
23#include <TFile.h>
24#include <TTree.h>
25#include <TParticle.h>
26
27#include "AliGenReaderEcalHijing.h"
28#include "AliMC.h"
29ClassImp(AliGenReaderEcalHijing)
30
31
32AliGenReaderEcalHijing::AliGenReaderEcalHijing()
33{
34// Default constructor
35 fNcurrent = 0;
36 fTreeNtuple = 0;
37}
38
39void AliGenReaderEcalHijing::Init()
40{
41//
42// reset the existing file environment and open a new root file if
43// the pointer to the Fluka tree is null
44
45 TFile *pFile=0;
46 if (!pFile) {
47 pFile = new TFile(fFileName);
48 pFile->cd();
49 printf("\n I have opened %s file \n", fFileName);
50 }
51// get the tree address in the Fluka boundary source file
52 fTreeNtuple = (TTree*)gDirectory->Get("h2");
53 TTree *h2=fTreeNtuple;
54 h2->SetMakeClass(1);
55//Set branch addresses
56 h2->SetBranchAddress("njatt", &fNjatt);
57 h2->SetBranchAddress("nahij", &fNahij);
58 h2->SetBranchAddress("nphij", &fNphij);
59 h2->SetBranchAddress("khij", fKhij) ;
60 h2->SetBranchAddress("pxhij", fPxhij);
61 h2->SetBranchAddress("pyhij", fPyhij);
62 h2->SetBranchAddress("pzhij", fPzhij);
63 h2->SetBranchAddress("ehij", fEhij) ;
64}
65
66Int_t AliGenReaderEcalHijing::NextEvent()
67{
68// Read the next event
58323221 69 Int_t nTracks=0, nread=0;
b152a071 70
71 TFile* pFile = fTreeNtuple->GetCurrentFile();
72 pFile->cd();
73
74 Int_t nentries = (Int_t) fTreeNtuple->GetEntries();
75 if (fNcurrent < nentries) {
76 Int_t nb = (Int_t)fTreeNtuple->GetEvent(fNcurrent);
77 nread += nb;
78 fNcurrent++;
79 printf("\n Next event contains %d tracks! \n", fNphij);
80 nTracks = fNphij;
81 fNparticle = 0;
82 return nTracks;
83 }
84 return 0;
85}
86
87TParticle* AliGenReaderEcalHijing::NextParticle()
88{
89 Float_t p[4];
90// Read the next particle
91 Int_t ipart = fKhij[fNparticle];
92 p[0] = fPxhij[fNparticle];
93 p[1] = fPyhij[fNparticle];
94 p[2] = fPzhij[fNparticle];
95 p[3] = fEhij[fNparticle];
96
97 Double_t amass = TDatabasePDG::Instance()->GetParticle(ipart)->Mass();
98
99 if(p[3] <= amass) {
100 Warning("Generate","Particle %d E = %f mass = %f \n",
101 ipart, p[3], amass);
102 }
103 TParticle* particle =
104 new TParticle(ipart, 0, -1, -1, -1, -1, p[0], p[1], p[2], p[3],
105 0., 0., 0., 0.);
106 fNparticle++;
107 return particle;
108}
109
110
111
112AliGenReaderEcalHijing& AliGenReaderEcalHijing::operator=(const AliGenReaderEcalHijing& rhs)
113{
114// Assignment operator
115 return *this;
116}
117
118
119
120
121
122