]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenReaderEcalHijing.cxx
Wagon macros, update cut on PHOS clusters close to border.
[u/mrichter/AliRoot.git] / EVGEN / AliGenReaderEcalHijing.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 // 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 //
24 #include <TFile.h>
25 #include <TParticle.h>
26 #include <TDatabasePDG.h>
27 #include <TTree.h>
28
29 #include "AliGenReaderEcalHijing.h"
30
31 ClassImp(AliGenReaderEcalHijing)
32
33 AliGenReaderEcalHijing::AliGenReaderEcalHijing():
34     fNcurrent(0),
35     fNparticle(0),
36     fTreeNtuple(0),
37     fNjatt(0),
38     fNahij(0),
39     fNphij(0)
40 {
41 // Default constructor
42 }
43
44 AliGenReaderEcalHijing::AliGenReaderEcalHijing(const AliGenReaderEcalHijing &reader):
45     AliGenReader(reader),
46     fNcurrent(0),
47     fNparticle(0),
48     fTreeNtuple(0),
49     fNjatt(0),
50     fNahij(0),
51     fNphij(0)
52 {
53     // Copy constructor
54     reader.Copy(*this);
55 }
56
57 void AliGenReaderEcalHijing::Init() 
58 {
59 //
60 // reset the existing file environment and open a new root file if
61 // the pointer to the Fluka tree is null
62     
63     TFile *pFile=0;
64     if (!pFile) {
65         pFile = new TFile(fFileName);
66         pFile->cd();
67         printf("\n I have opened %s file \n", fFileName);
68     }
69 // get the tree address in the Fluka boundary source file
70     fTreeNtuple = (TTree*)gDirectory->Get("h2");
71     TTree *h2=fTreeNtuple;
72     h2->SetMakeClass(1);
73 //Set branch addresses
74     h2->SetBranchAddress("njatt", &fNjatt);
75     h2->SetBranchAddress("nahij", &fNahij);
76     h2->SetBranchAddress("nphij", &fNphij);
77     h2->SetBranchAddress("khij",   fKhij) ;
78     h2->SetBranchAddress("pxhij",  fPxhij);
79     h2->SetBranchAddress("pyhij",  fPyhij);
80     h2->SetBranchAddress("pzhij",  fPzhij);
81     h2->SetBranchAddress("ehij",   fEhij) ;
82 }
83
84 Int_t AliGenReaderEcalHijing::NextEvent() 
85 {
86 // Read the next event  
87     Int_t nTracks=0, nread=0;
88     
89     TFile* pFile = fTreeNtuple->GetCurrentFile();
90     pFile->cd();
91
92     Int_t nentries = (Int_t) fTreeNtuple->GetEntries();
93     if (fNcurrent < nentries) {
94         Int_t nb = (Int_t)fTreeNtuple->GetEvent(fNcurrent);
95         nread += nb;
96         fNcurrent++;
97         printf("\n Next event contains %d tracks! \n", fNphij);
98         nTracks    = fNphij;
99         fNparticle = 0;
100         return nTracks;
101     }
102     return 0;
103 }
104
105 TParticle* AliGenReaderEcalHijing::NextParticle() 
106 {
107 // Read the next particle
108
109     Float_t p[4];
110     Int_t ipart = fKhij[fNparticle];
111     p[0] = fPxhij[fNparticle];
112     p[1] = fPyhij[fNparticle];      
113     p[2] = fPzhij[fNparticle];
114     p[3] = fEhij[fNparticle];
115     
116     Double_t amass = TDatabasePDG::Instance()->GetParticle(ipart)->Mass();
117
118     if(p[3] <= amass) {
119         Warning("Generate","Particle %d  E = %f mass = %f \n",
120                 ipart, p[3], amass);
121     } 
122     TParticle* particle = 
123         new TParticle(ipart, 0, -1, -1, -1, -1, p[0], p[1], p[2], p[3], 
124                       0., 0., 0., 0.);
125     fNparticle++;
126     return particle;
127 }
128
129
130
131 AliGenReaderEcalHijing& AliGenReaderEcalHijing::operator=(const  AliGenReaderEcalHijing& rhs)
132 {
133 // Assignment operator
134     rhs.Copy(*this);
135     return (*this);
136 }
137
138 void AliGenReaderEcalHijing::Copy(TObject&) const
139 {
140     //
141     // Copy 
142     //
143     Fatal("Copy","Not implemented!\n");
144 }
145
146
147
148
149
150