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