]>
Commit | Line | Data |
---|---|---|
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 | ||
88cb7938 | 16 | /* $Id$ */ |
93de1f83 | 17 | // |
18 | // Realisation of AliGenReader to be used with AliGenExtFile | |
19 | // It reads Hijing events from a ntuple like event structure. | |
20 | // The event format is defined in Init() | |
21 | // NextEvent() is used to loop over events and NextParticle() to loop over particles. | |
22 | // Author: andreas.morsch@cern.ch | |
23 | // | |
b152a071 | 24 | #include <TFile.h> |
b152a071 | 25 | #include <TParticle.h> |
88cb7938 | 26 | #include <TTree.h> |
b152a071 | 27 | |
28 | #include "AliGenReaderEcalHijing.h" | |
88cb7938 | 29 | |
706938e6 | 30 | ClassImp(AliGenReaderEcalHijing) |
b152a071 | 31 | |
1c56e311 | 32 | AliGenReaderEcalHijing::AliGenReaderEcalHijing(): |
33 | fNcurrent(0), | |
34 | fNparticle(0), | |
35 | fTreeNtuple(0), | |
36 | fNjatt(0), | |
37 | fNahij(0), | |
38 | fNphij(0) | |
b152a071 | 39 | { |
40 | // Default constructor | |
1c56e311 | 41 | } |
42 | ||
43 | AliGenReaderEcalHijing::AliGenReaderEcalHijing(const AliGenReaderEcalHijing &reader): | |
44 | AliGenReader(reader), | |
45 | fNcurrent(0), | |
46 | fNparticle(0), | |
47 | fTreeNtuple(0), | |
48 | fNjatt(0), | |
49 | fNahij(0), | |
50 | fNphij(0) | |
51 | { | |
52 | // Copy constructor | |
53 | reader.Copy(*this); | |
b152a071 | 54 | } |
55 | ||
56 | void AliGenReaderEcalHijing::Init() | |
57 | { | |
58 | // | |
59 | // reset the existing file environment and open a new root file if | |
60 | // the pointer to the Fluka tree is null | |
61 | ||
62 | TFile *pFile=0; | |
63 | if (!pFile) { | |
64 | pFile = new TFile(fFileName); | |
65 | pFile->cd(); | |
66 | printf("\n I have opened %s file \n", fFileName); | |
67 | } | |
68 | // get the tree address in the Fluka boundary source file | |
69 | fTreeNtuple = (TTree*)gDirectory->Get("h2"); | |
70 | TTree *h2=fTreeNtuple; | |
71 | h2->SetMakeClass(1); | |
72 | //Set branch addresses | |
73 | h2->SetBranchAddress("njatt", &fNjatt); | |
74 | h2->SetBranchAddress("nahij", &fNahij); | |
75 | h2->SetBranchAddress("nphij", &fNphij); | |
76 | h2->SetBranchAddress("khij", fKhij) ; | |
77 | h2->SetBranchAddress("pxhij", fPxhij); | |
78 | h2->SetBranchAddress("pyhij", fPyhij); | |
79 | h2->SetBranchAddress("pzhij", fPzhij); | |
80 | h2->SetBranchAddress("ehij", fEhij) ; | |
81 | } | |
82 | ||
83 | Int_t AliGenReaderEcalHijing::NextEvent() | |
84 | { | |
85 | // Read the next event | |
58323221 | 86 | Int_t nTracks=0, nread=0; |
b152a071 | 87 | |
88 | TFile* pFile = fTreeNtuple->GetCurrentFile(); | |
89 | pFile->cd(); | |
90 | ||
91 | Int_t nentries = (Int_t) fTreeNtuple->GetEntries(); | |
92 | if (fNcurrent < nentries) { | |
93 | Int_t nb = (Int_t)fTreeNtuple->GetEvent(fNcurrent); | |
94 | nread += nb; | |
95 | fNcurrent++; | |
96 | printf("\n Next event contains %d tracks! \n", fNphij); | |
97 | nTracks = fNphij; | |
98 | fNparticle = 0; | |
99 | return nTracks; | |
100 | } | |
101 | return 0; | |
102 | } | |
103 | ||
104 | TParticle* AliGenReaderEcalHijing::NextParticle() | |
105 | { | |
b152a071 | 106 | // Read the next particle |
93de1f83 | 107 | |
108 | Float_t p[4]; | |
b152a071 | 109 | Int_t ipart = fKhij[fNparticle]; |
110 | p[0] = fPxhij[fNparticle]; | |
111 | p[1] = fPyhij[fNparticle]; | |
112 | p[2] = fPzhij[fNparticle]; | |
113 | p[3] = fEhij[fNparticle]; | |
114 | ||
115 | Double_t amass = TDatabasePDG::Instance()->GetParticle(ipart)->Mass(); | |
116 | ||
117 | if(p[3] <= amass) { | |
118 | Warning("Generate","Particle %d E = %f mass = %f \n", | |
119 | ipart, p[3], amass); | |
120 | } | |
121 | TParticle* particle = | |
122 | new TParticle(ipart, 0, -1, -1, -1, -1, p[0], p[1], p[2], p[3], | |
123 | 0., 0., 0., 0.); | |
124 | fNparticle++; | |
125 | return particle; | |
126 | } | |
127 | ||
128 | ||
129 | ||
130 | AliGenReaderEcalHijing& AliGenReaderEcalHijing::operator=(const AliGenReaderEcalHijing& rhs) | |
131 | { | |
132 | // Assignment operator | |
198bb1c7 | 133 | rhs.Copy(*this); |
134 | return (*this); | |
135 | } | |
136 | ||
dc1d768c | 137 | void AliGenReaderEcalHijing::Copy(TObject&) const |
198bb1c7 | 138 | { |
139 | // | |
140 | // Copy | |
141 | // | |
142 | Fatal("Copy","Not implemented!\n"); | |
b152a071 | 143 | } |
144 | ||
145 | ||
146 | ||
147 | ||
148 | ||
149 |