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