]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenReaderEcalJets.cxx
Additional initialization. Do not remove fHeader, it is done in AliRunLoader
[u/mrichter/AliRoot.git] / EVGEN / AliGenReaderEcalJets.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 Pythia Jet 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 <TTree.h>
27
28 #include "AliGenReaderEcalJets.h"
29
30 ClassImp(AliGenReaderEcalJets)
31
32
33 AliGenReaderEcalJets::AliGenReaderEcalJets() 
34 {
35 // Default constructor
36     fNcurrent   = 0;
37     fTreeNtuple = 0;
38 }
39
40 void AliGenReaderEcalJets::Init() 
41 {
42 //
43 // reset the existing file environment and open a new root file if
44 // the pointer to the Fluka tree is null
45     
46     TFile *pFile=0;
47     if (!pFile) {
48         pFile = new TFile(fFileName);
49         pFile->cd();
50         printf("\n I have opened %s file \n", fFileName);
51     }
52 // get the tree address in the Fluka boundary source file
53     fTreeNtuple = (TTree*)gDirectory->Get("h1");
54     TTree *h1=fTreeNtuple;
55     h1->SetMakeClass(1);
56 //Set branch addresses
57     h1->SetBranchAddress("nev",   &fNev);
58     h1->SetBranchAddress("x",      fX);
59     h1->SetBranchAddress("xtyp",   fXtyp);
60     h1->SetBranchAddress("npart", &fNpart);
61     h1->SetBranchAddress("xpt",    fXpt);
62     h1->SetBranchAddress("xeta",   fXeta);
63     h1->SetBranchAddress("xphi",   fXphi);
64     h1->SetBranchAddress("xid",    fXid);
65     h1->SetBranchAddress("njet",  &fNjet);
66     h1->SetBranchAddress("jet",    fJet);
67     h1->SetBranchAddress("jeta",   fJeta);
68     h1->SetBranchAddress("jphi",   fJphi);
69     h1->SetBranchAddress("nsjet", &fNsjet);
70     h1->SetBranchAddress("jset",   fJset);
71     h1->SetBranchAddress("jseta",  fJseta);
72     h1->SetBranchAddress("jsphi",  fJsphi);
73     h1->SetBranchAddress("npjet", &fNpjet);
74     h1->SetBranchAddress("jpet",   fJpet);
75     h1->SetBranchAddress("jpeta",  fJpeta);
76     h1->SetBranchAddress("jpphi",  fJpphi);
77 }
78
79 Int_t AliGenReaderEcalJets::NextEvent() 
80 {
81 // Read the next event  
82     Int_t nTracks=0, nread=0;
83     
84     TFile* pFile = fTreeNtuple->GetCurrentFile();
85     pFile->cd();
86
87     Int_t nentries = (Int_t) fTreeNtuple->GetEntries();
88     if (fNcurrent < nentries) {
89         Int_t nb = (Int_t)fTreeNtuple->GetEvent(fNcurrent);
90         nread += nb;
91         fNcurrent++;
92         printf("\n Next event contains %d tracks! \n", fNpjet);
93         nTracks    = fNpjet;
94         fNparticle = 0;
95         return nTracks;
96     }
97     return 0;
98 }
99
100 TParticle* AliGenReaderEcalJets::NextParticle() 
101 {
102 // Read the next particle
103
104     Float_t p[4];
105     Int_t    ipart  = fXid[fNparticle];
106     Float_t  pt     = fXpt[fNparticle];
107     Float_t  eta    = fXeta[fNparticle];
108     Float_t  phi    = fXphi[fNparticle];
109     Float_t  theta  = 2.*TMath::ATan(TMath::Exp(-eta));
110     Double_t amass  = TDatabasePDG::Instance()->GetParticle(ipart)->Mass();
111
112     p[0] = pt*TMath::Sin(phi);
113     p[1] = pt*TMath::Cos(phi);      
114     p[2] = pt/TMath::Cos(theta);
115     p[3] = TMath::Sqrt(pt*pt+p[2]*p[2]+amass*amass);
116     
117
118     TParticle* particle = 
119         new TParticle(ipart, 0, -1, -1, -1, -1, p[0], p[1], p[2], p[3], 
120                       0., 0., 0., 0.);
121     fNparticle++;
122     return particle;
123 }
124
125
126
127 AliGenReaderEcalJets& AliGenReaderEcalJets::operator=(const  AliGenReaderEcalJets& rhs)
128 {
129 // Assignment operator
130     rhs.Copy(*this);
131     return (*this);
132 }
133
134 void AliGenReaderEcalJets::Copy(TObject&) const
135 {
136     //
137     // Copy 
138     //
139     Fatal("Copy","Not implemented!\n");
140 }
141
142
143
144
145
146
147