]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenReaderCwn.cxx
Warnings cleaned up (Alpha)
[u/mrichter/AliRoot.git] / EVGEN / AliGenReaderCwn.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
17 /*
18 $Log$
19 Revision 1.3  2002/10/21 09:20:51  alibrary
20 Introduce Riostream.h and remove unused variables
21
22 Revision 1.2  2001/11/12 14:31:00  morsch
23 Memory leaks fixed. (M. Bondila)
24
25 Revision 1.1  2001/11/09 09:10:46  morsch
26 Realisation of AliGenReader that reads the old cwn event format.
27
28 */
29
30 // Read the old ALICE event format based on CW-ntuples
31 // http://consult.cern.ch/alice/Internal_Notes/1995/32/abstract
32 // .cwn file have to be converted to .root using h2root
33 // Use SetFileName(file) to read from "file" 
34 // Author: andreas.morsch@cern.ch
35
36 #include <TFile.h>
37 #include <TTree.h>
38 #include <TParticle.h>
39
40 #include "AliGenReaderCwn.h"
41 #include "AliMC.h"
42 ClassImp(AliGenReaderCwn);
43
44
45 AliGenReaderCwn::AliGenReaderCwn() 
46 {
47 // Default constructor
48     fNcurrent   = 0;
49     fTreeNtuple = 0;
50 }
51
52 AliGenReaderCwn::~AliGenReaderCwn()
53 {
54     delete fTreeNtuple;
55 }
56
57 void AliGenReaderCwn::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("h888");
71
72     TTree *h2=fTreeNtuple;
73 //Set branch addresses
74     h2->SetBranchAddress("Nihead",&fNihead);
75     h2->SetBranchAddress("Ihead",fIhead);
76     h2->SetBranchAddress("Nrhead",&fNrhead);
77     h2->SetBranchAddress("Rhead",fRhead);
78     h2->SetBranchAddress("Idpart",&fIdpart);
79     h2->SetBranchAddress("Theta",&fTheta);
80     h2->SetBranchAddress("Phi",&fPhi);
81     h2->SetBranchAddress("P",&fP);
82     h2->SetBranchAddress("E",&fE);
83 }
84
85 Int_t AliGenReaderCwn::NextEvent() 
86 {
87 // Read the next event  
88     Int_t nTracks;
89     fNparticle = 0;
90     TFile* pFile = fTreeNtuple->GetCurrentFile();
91     pFile->cd();
92
93     Int_t nentries = (Int_t) fTreeNtuple->GetEntries();
94     if (fNcurrent < nentries) {
95         fNcurrent++;
96         
97         Int_t i5=fIhead[4];
98         Int_t i6=fIhead[5];
99         if (i5==0) {
100             printf("\n This should never happen !\n");
101             nTracks = 0;
102         } else {
103             printf("\n Next event contains %d tracks! \n", i6);
104             nTracks = i6;
105         }    
106         fNparticleMax = nTracks;
107         return nTracks;
108     }
109
110     return 0;
111 }
112
113 TParticle* AliGenReaderCwn::NextParticle() 
114 {
115 //
116 //  
117     Float_t prwn;
118     Float_t p[4];
119 // Read the next particle
120     if (fCode == kGEANT3) fIdpart=gMC->PDGFromId(fIdpart);
121     Double_t amass = TDatabasePDG::Instance()->GetParticle(fIdpart)->Mass();
122     if(fE<=amass) {
123         Warning("Generate","Particle %d  E = %f mass = %f %f %f \n",
124                 fIdpart,fE,amass, fPhi, fTheta);
125         prwn=0;
126     } else {
127         prwn=sqrt((fE+amass)*(fE-amass));
128     }
129
130     fTheta *= TMath::Pi()/180.;
131     fPhi    = (fPhi-180)*TMath::Pi()/180.;      
132     p[0] = prwn*TMath::Sin(fTheta)*TMath::Cos(fPhi);
133     p[1] = prwn*TMath::Sin(fTheta)*TMath::Sin(fPhi);      
134     p[2] = prwn*TMath::Cos(fTheta);
135     p[3] = fE;
136     TParticle* particle = new TParticle(fIdpart, 0, -1, -1, -1, -1, p[0], p[1], p[2], p[3], 0., 0., 0., 0.);
137     fNcurrent++;
138     fNparticle++;
139     return particle;
140 }
141
142
143
144 AliGenReaderCwn& AliGenReaderCwn::operator=(const  AliGenReaderCwn& rhs)
145 {
146 // Assignment operator
147     return *this;
148 }
149
150
151