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