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