]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenExtFile.cxx
Removal of useless dependecies via forward declarations
[u/mrichter/AliRoot.git] / EVGEN / AliGenExtFile.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 $Log$
18 Revision 1.10  2000/07/11 18:24:55  fca
19 Coding convention corrections + few minor bug fixes
20
21 Revision 1.9  2000/06/14 15:20:09  morsch
22 Include clean-up (IH)
23
24 Revision 1.8  2000/06/09 20:36:44  morsch
25 All coding rule violations except RS3 corrected
26
27 Revision 1.7  2000/02/16 14:56:27  morsch
28 Convert geant particle code into pdg code before putting particle on the stack.
29
30 Revision 1.6  1999/11/09 07:38:48  fca
31 Changes for compatibility with version 2.23 of ROOT
32
33 Revision 1.5  1999/09/29 09:24:12  fca
34 Introduction of the Copyright and cvs Log
35
36 */
37
38 #include <iostream.h>
39
40 #include "AliGenExtFile.h"
41 #include "AliMC.h"
42 #include "AliRun.h"
43
44 #include <TDirectory.h>
45 #include <TDatabasePDG.h>
46 #include <TFile.h>
47 #include <stdlib.h>
48
49  ClassImp(AliGenExtFile)
50      AliGenExtFile::AliGenExtFile()
51          :AliGenerator(-1)
52 {
53 //  Constructor
54     fName="ExtFile";
55     fTitle="Primaries from ext. File";
56     fFileName="dtujet93.root";
57     fTreeNtuple=0;
58     fNcurrent=0;
59 //
60 //  Read all particles
61     fNpart=-1;
62 }
63
64 AliGenExtFile::AliGenExtFile(Int_t npart)
65     :AliGenerator(npart)
66 {
67 //  Constructor
68     fName="ExtFile";
69     fTitle="Primaries from ext. File";
70     fFileName="dtujet93.root";
71     fTreeNtuple=0;
72     fNcurrent=0;
73 }
74
75 AliGenExtFile::AliGenExtFile(const AliGenExtFile & ExtFile)
76 {
77 // copy constructor
78 }
79 //____________________________________________________________
80 AliGenExtFile::~AliGenExtFile()
81 {
82 // Destructor
83     delete fTreeNtuple;
84 }
85
86 //____________________________________________________________
87 void AliGenExtFile::NtupleInit() 
88 {
89 //
90 // reset the existing file environment and open a new root file if
91 // the pointer to the Fluka tree is null
92     
93     TFile *pFile=0;
94     if (fTreeNtuple==0) {
95         if (!pFile) {
96             pFile = new TFile(fFileName);
97             pFile->cd();
98             cout<<"I have opened "<<fFileName<<" file "<<endl;
99         }
100 // get the tree address in the Fluka boundary source file
101         fTreeNtuple = (TTree*)gDirectory->Get("h888");
102     } else {
103         pFile = fTreeNtuple->GetCurrentFile();
104         pFile->cd();
105     }
106
107     TTree *h2=fTreeNtuple;
108 //Set branch addresses
109 //Set branch addresses
110     h2->SetBranchAddress("Nihead",&fNihead);
111     h2->SetBranchAddress("Ihead",fIhead);
112     h2->SetBranchAddress("Nrhead",&fNrhead);
113     h2->SetBranchAddress("Rhead",fRhead);
114     h2->SetBranchAddress("Idpart",&fIdpart);
115     h2->SetBranchAddress("Theta",&fTheta);
116     h2->SetBranchAddress("Phi",&fPhi);
117     h2->SetBranchAddress("P",&fP);
118     h2->SetBranchAddress("E",&fE);
119 }
120
121
122 //____________________________________________________________
123 void AliGenExtFile::Generate()
124 {
125 // Generate particles
126
127   Float_t polar[3]= {0,0,0};
128   //
129   Float_t origin[3]={0,0,0};
130   Float_t p[3];
131   Float_t random[6];
132   Float_t prwn;
133   Int_t i, j, nt, nTracks=0;
134   //
135   NtupleInit();
136   TTree *h2=fTreeNtuple;
137   Int_t nentries = (Int_t) h2->GetEntries();
138   // loop over number of particles
139   Int_t nb = (Int_t)h2->GetEvent(fNcurrent);
140   Int_t i5=fIhead[4];
141   Int_t i6=fIhead[5];
142
143   for (j=0;j<3;j++) origin[j]=fOrigin[j];
144   if(fVertexSmear==kPerTrack) {
145     gMC->Rndm(random,6);
146     for (j=0;j<3;j++) {
147         origin[j]+=fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
148             TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
149     }
150   }
151
152   if (fNcurrent >= nentries) {
153       printf("\n No more entries !!! !\n");
154       return;
155   }
156   
157           
158   if (i5==0) {
159       printf("\n This should never happen !\n");
160   } else {
161       printf("\n Next event contains %d tracks! \n", i6);
162       nTracks=i6;
163   }
164   for (i=0; i<nTracks; i++) {
165       fIdpart=gMC->PDGFromId(fIdpart);
166       Double_t amass = TDatabasePDG::Instance()->GetParticle(fIdpart)->Mass();
167       if(fE<=amass) {
168         Warning("Generate","Particle %d no %d E = %f mass = %f\n",
169                 fIdpart,i,fE,amass);
170         prwn=0;
171       } else {
172         prwn=sqrt((fE+amass)*(fE-amass));
173       }
174
175       fTheta *= TMath::Pi()/180.;
176       fPhi    = (fPhi-180)*TMath::Pi()/180.;      
177       if(fTheta<fThetaMin || fTheta>fThetaMax ||
178          fPhi<fPhiMin || fPhi>fPhiMax         ||
179          prwn<fPMin || prwn>fPMax)          
180       {
181           ;
182       } else {
183           p[0]=prwn*TMath::Sin(fTheta)*TMath::Cos(fPhi);
184           p[1]=prwn*TMath::Sin(fTheta)*TMath::Sin(fPhi);      
185           p[2]=prwn*TMath::Cos(fTheta);
186           
187           if(fVertexSmear==kPerTrack) {
188               gMC->Rndm(random,6);
189               for (j=0;j<3;j++) {
190                   origin[j]=fOrigin[j]
191                       +fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
192                       TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
193               }
194           }
195           gAlice->SetTrack(fTrackIt,-1,fIdpart,p,origin,polar,0,"Primary",nt);
196       }
197       fNcurrent++;
198       nb = (Int_t)h2->GetEvent(fNcurrent); 
199   }
200  
201     TFile *pFile=0;
202 // Get AliRun object or create it 
203     if (!gAlice) {
204         gAlice = (AliRun*)pFile->Get("gAlice");
205         if (gAlice) printf("AliRun object found on file\n");
206         if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
207     }
208     TTree *fAli=gAlice->TreeK();
209     if (fAli) pFile =fAli->GetCurrentFile();
210     pFile->cd();
211 }
212
213
214 AliGenExtFile& AliGenExtFile::operator=(const  AliGenExtFile& rhs)
215 {
216 // Assignment operator
217     return *this;
218 }
219
220
221
222
223
224
225
226