]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenReaderCwn.cxx
added settings for magnetic field
[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 /* $Id$ */
18
19 // Read the old ALICE event format based on CW-ntuples
20 // http://consult.cern.ch/alice/Internal_Notes/1995/32/abstract
21 // .cwn file have to be converted to .root using h2root
22 // Use SetFileName(file) to read from "file" 
23 // Author: andreas.morsch@cern.ch
24
25 #include <TFile.h>
26 #include <TParticle.h>
27 #include <TDatabasePDG.h>
28 #include <TTree.h>
29 #include <TVirtualMC.h>
30
31 #include "AliGenReaderCwn.h"
32
33 ClassImp(AliGenReaderCwn)
34
35 AliGenReaderCwn::AliGenReaderCwn():
36     fNcurrent(0),
37     fNparticle(0),
38     fNparticleMax(0),
39     fTreeNtuple(0),
40     fNihead(0),
41     fNrhead(0),
42     fIdpart(0),
43     fTheta(0.),
44     fPhi(0.),
45     fP(0.),
46     fE(0.)
47 {
48 // Default constructor
49 }
50
51 AliGenReaderCwn::AliGenReaderCwn(const AliGenReaderCwn &reader):
52     AliGenReader(reader),
53     fNcurrent(0),
54     fNparticle(0),
55     fNparticleMax(0),
56     fTreeNtuple(0),
57     fNihead(0),
58     fNrhead(0),
59     fIdpart(0),
60     fTheta(0.),
61     fPhi(0.),
62     fP(0.),
63     fE(0.)
64 {
65     // Copy constructor
66     reader.Copy(*this);
67 }
68
69
70 AliGenReaderCwn::~AliGenReaderCwn()
71 {
72     delete fTreeNtuple;
73 }
74
75 void AliGenReaderCwn::Init() 
76 {
77 //
78 // reset the existing file environment and open a new root file if
79 // the pointer to the Fluka tree is null
80     
81     TFile *pFile=0;
82     if (!pFile) {
83         pFile = new TFile(fFileName);
84         pFile->cd();
85         printf("\n I have opened %s file \n", fFileName);
86     }
87 // get the tree address in the Fluka boundary source file
88     fTreeNtuple = (TTree*)gDirectory->Get("h888");
89
90     TTree *h2=fTreeNtuple;
91 //Set branch addresses
92     h2->SetBranchAddress("Nihead",&fNihead);
93     h2->SetBranchAddress("Ihead",fIhead);
94     h2->SetBranchAddress("Nrhead",&fNrhead);
95     h2->SetBranchAddress("Rhead",fRhead);
96     h2->SetBranchAddress("Idpart",&fIdpart);
97     h2->SetBranchAddress("Theta",&fTheta);
98     h2->SetBranchAddress("Phi",&fPhi);
99     h2->SetBranchAddress("P",&fP);
100     h2->SetBranchAddress("E",&fE);
101 }
102
103 Int_t AliGenReaderCwn::NextEvent() 
104 {
105 // Read the next event  
106     Int_t nTracks;
107     fNparticle = 0;
108     TFile* pFile = fTreeNtuple->GetCurrentFile();
109     pFile->cd();
110
111     Int_t nentries = (Int_t) fTreeNtuple->GetEntries();
112     if (fNcurrent < nentries) {
113         fNcurrent++;
114         
115         Int_t i5=fIhead[4];
116         Int_t i6=fIhead[5];
117         if (i5==0) {
118             printf("\n This should never happen !\n");
119             nTracks = 0;
120         } else {
121             printf("\n Next event contains %d tracks! \n", i6);
122             nTracks = i6;
123         }    
124         fNparticleMax = nTracks;
125         return nTracks;
126     }
127
128     return 0;
129 }
130
131 TParticle* AliGenReaderCwn::NextParticle() 
132 {
133 //
134 //  
135     Float_t prwn;
136     Float_t p[4];
137 // Read the next particle
138     if (fCode == kGEANT3) fIdpart=gMC->PDGFromId(fIdpart);
139     Double_t amass = TDatabasePDG::Instance()->GetParticle(fIdpart)->Mass();
140     if(fE<=amass) {
141         Warning("Generate","Particle %d  E = %f mass = %f %f %f \n",
142                 fIdpart,fE,amass, fPhi, fTheta);
143         prwn=0;
144     } else {
145         prwn=sqrt((fE+amass)*(fE-amass));
146     }
147
148     fTheta *= TMath::Pi()/180.;
149     fPhi    = (fPhi-180)*TMath::Pi()/180.;      
150     p[0] = prwn*TMath::Sin(fTheta)*TMath::Cos(fPhi);
151     p[1] = prwn*TMath::Sin(fTheta)*TMath::Sin(fPhi);      
152     p[2] = prwn*TMath::Cos(fTheta);
153     p[3] = fE;
154     TParticle* particle = new TParticle(fIdpart, 0, -1, -1, -1, -1, p[0], p[1], p[2], p[3], 0., 0., 0., 0.);
155     fNcurrent++;
156     fNparticle++;
157     return particle;
158 }
159
160
161
162 AliGenReaderCwn& AliGenReaderCwn::operator=(const  AliGenReaderCwn& rhs)
163 {
164 // Assignment operator
165     rhs.Copy(*this);
166     return *this;
167 }
168
169 void AliGenReaderCwn::Copy(TObject&) const
170 {
171     //
172     // Copy 
173     //
174     Fatal("Copy","Not implemented!\n");
175 }
176
177
178
179