]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenFLUKAsource.cxx
All coding rule violations except RS3 corrected
[u/mrichter/AliRoot.git] / EVGEN / AliGenFLUKAsource.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.9  2000/03/07 13:52:54  morsch
19 static Int_t irwn=0;
20
21 Revision 1.8  2000/02/14 14:49:38  morsch
22 Correct particle type for gamma and neutrons
23 More consistent calculation of momentum from kin. energy and mass
24
25 Revision 1.7  1999/11/03 17:43:20  fca
26 New version from G.Martinez & A.Morsch
27
28 Revision 1.6  1999/09/29 09:24:12  fca
29 Introduction of the Copyright and cvs Log
30
31 */
32
33 #include "AliGenFLUKAsource.h"
34 #include "AliGenMUONlib.h"
35 #include "AliMC.h"
36 #include "AliRun.h"
37 #include "AliPDG.h"
38 #include <TDirectory.h>
39 #include <TFile.h>
40 #include <TTree.h>
41 #include <TChain.h>
42 #include <stdlib.h>
43  ClassImp(AliGenFLUKAsource)
44      AliGenFLUKAsource::AliGenFLUKAsource()
45          :AliGenerator(-1)
46 {
47     // Constructor
48     fName="FLUKA";
49     fTitle="FLUKA Boundary Source";
50     // Read in all particle types by default
51     fIkine=6;
52     // Set maximum admitted age of particles to 1.e-05 by default 
53     fAgeMax=1.e-05;
54     // Do not add weight
55     fAddWeight=1.;
56     // Shift the z-coordinate of the impact point by 4.5 cm only if it reads 
57     // from  specific boundary source to the chamber (fZshift=4.5;),else there 
58     // is no need to shift as it reads boundary source for the whole volume of 
59     // the Muon Arm; the default value corresponds to boundary source for the
60     // whole volume of the MUON Arm 
61     fZshift=0;
62     // Set the default file 
63     fFileName="flukasource.root";
64
65     fTreeFluka=0;
66     fTreeChain = new TChain("h1");
67 //
68 //  Read all particles
69     fNpart=-1;
70
71     
72 }
73
74 AliGenFLUKAsource::AliGenFLUKAsource(Int_t npart)
75     :AliGenerator(npart)
76 {
77     // Constructor
78     fName="FLUKA";
79     fTitle="FLUKA Boundary Source";
80     // Read in all particle types by default
81     fIkine=6;
82     // Set maximum admitted age of particles to 1.e-05 by default 
83     fAgeMax=1.e-05;
84     // Do not add weight
85     fAddWeight=1.;
86     // Shift the z-coordinate of the impact point by 4.5 cm only if it reads 
87     // from  specific boundary source to the chamber (fZshift=4.5;),else there 
88     // is no need to shift as it reads boundary source for the whole volume of 
89     // the Muon Arm; the default value corresponds to boundary source for the
90     // whole volume of the MUON Arm 
91     fZshift=0;
92     // Set the default file 
93     fFileName="flukasource.root";
94
95     fTreeFluka=0;
96     fTreeChain = new TChain("h1"); 
97     fSourceId=-1;
98 }
99
100 AliGenFLUKAsource::AliGenFLUKAsource(const AliGenFLUKAsource & FLUKAsource)
101 {
102 // copy constructor
103 }
104
105
106 //____________________________________________________________
107 AliGenFLUKAsource::~AliGenFLUKAsource()
108 {
109 // Destructor
110  if (fTreeFluka) delete fTreeFluka;
111  if (fTreeChain) delete fTreeChain;
112 }
113
114 void AliGenFLUKAsource::AddFile(const Text_t *filename)
115 {
116 // Add a file to the chain
117     fTreeChain->Add(filename);
118     
119 }
120
121
122 //____________________________________________________________
123 void AliGenFLUKAsource::FlukaInit() 
124 {
125 // Set branch addresses of data entries
126     TChain *h2=fTreeChain;
127 //Set branch addresses
128     h2->SetBranchAddress("Ip",&fIp);
129     h2->SetBranchAddress("Ipp",&fIpp);
130     h2->SetBranchAddress("Xi",&fXi);
131     h2->SetBranchAddress("Yi",&fYi);
132     h2->SetBranchAddress("Zi",&fZi);
133     h2->SetBranchAddress("Px",&fPx);
134     h2->SetBranchAddress("Py",&fPy);
135     h2->SetBranchAddress("Pz",&fPz);
136     h2->SetBranchAddress("Ekin",&fEkin);
137     h2->SetBranchAddress("Zv",&fZv);
138     h2->SetBranchAddress("Rv",&fRv);
139     h2->SetBranchAddress("Itra",&fItra);
140     h2->SetBranchAddress("Igas",&fIgas);
141     h2->SetBranchAddress("Wgt",&fWgt);
142     h2->SetBranchAddress("Etag",&fEtag);
143     h2->SetBranchAddress("Ptg",&fPtg);
144     h2->SetBranchAddress("Age",&fAge);
145 }
146
147 //____________________________________________________________
148 void AliGenFLUKAsource::Generate()
149 {
150 // Generate one event 
151     AliMC* gMC = AliMC::GetMC();
152
153     const Int_t kIfluge[28]={kProton, kProtonBar, kElectron, kPositron,
154                           kNuE, kNuEBar, kGamma, kNeutron, kNeutronBar,
155                           kMuonPlus, kMuonMinus, kK0Long , kPiPlus, kPiMinus,
156                           kKPlus, kKMinus, kLambda0, kLambda0Bar, kK0Short,
157                           kSigmaMinus, kSigmaPlus, kSigma0, kPi0, kK0, kK0Bar,
158                           0,kNuMu,kNuMuBar};
159     Float_t polar[3]= {0,0,0};
160   //
161     Float_t origin[3];
162     Float_t p[3];
163     Float_t prwn;
164     Float_t wgt, fwgt;
165     Float_t phi;
166     char name[100];
167     Float_t amass, charge, tlife;
168     Int_t itrtyp;
169     Int_t iwgt;
170     Int_t i, j, part, nt;
171     static Int_t irwn=0;
172     //
173     Float_t random[2];
174   //
175     FlukaInit();
176     TChain *h2=fTreeChain;
177     Int_t nentries = (Int_t) h2->GetEntries();
178     if (fNpart == -1) fNpart=Int_t(nentries*fFrac);
179   
180   // loop over number of particles
181     Int_t nb=0;
182     Int_t ev=gMC->CurrentEvent();
183     for (i=0; i<fNpart;i++) {
184         Int_t entry=fNpart*(ev-1)+i; 
185         nb = (Int_t)h2->GetEvent(entry); 
186         if (irwn > nentries) {
187             printf("No more entries in the FLUKA boundary source file\n");
188             TFile *pFile=0;
189             // Get AliRun object or create it 
190             if (!gAlice) {
191                 gAlice = (AliRun*)pFile->Get("gAlice");
192                 if (gAlice) printf("AliRun object found on file\n");
193                 if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
194             }
195             TTree *fAli=gAlice->TreeK();
196             if (fAli) pFile =fAli->GetCurrentFile();
197             pFile->cd();
198             printf("Generate - I'm out \n");
199             return;
200         }   
201
202         if (fSourceId != -1 && fIgas !=fSourceId) {
203             irwn++;
204             continue;
205         }
206         
207         if (fIp > 28 || fIp < 0) {
208             irwn++;
209             continue;
210         }
211         
212         if ((fIp != fIkine && fIkine != 6 && fIkine != 9 && fIkine != 10) || fAge > fAgeMax){
213             irwn++;
214             continue;
215         } else if (fIkine == 9) {
216             if (fIp == 7 || fIp == 8 || fAge > fAgeMax) { 
217                 irwn++;
218                 continue;
219             }
220         } else if (fIkine == 10) {
221             if (fIp == 8 || fAge > fAgeMax) { 
222                 irwn++;
223                 continue;
224             }
225         }
226     
227
228         irwn++;
229 //
230 // PDG code from FLUKA particle type (fIp)
231         part=kIfluge[int(fIp)-1];       
232 //
233 // Calculate momentum from kinetic energy and mass of the particle
234         gMC->Gfpart(part, name, itrtyp,  
235                     amass, charge, tlife); 
236         prwn=fEkin*sqrt(1. + 2.*amass/fEkin);
237
238
239         origin[0]=fYi;
240         origin[1]=fXi;
241         origin[2]=fZi;
242         
243         p[0]=fPy*prwn;
244         p[1]=fPx*prwn;
245         p[2]=fPz*prwn;
246
247         //handle particle weight correctly
248         wgt = (part == 13) ? fWgt*fAddWeight : fWgt;
249         iwgt=Int_t(wgt);
250         fwgt=wgt-Float_t(iwgt);
251         gMC->Rndm(random,2);
252         if (random[0] < fwgt) iwgt++;
253         if (part==1 && iwgt>100) iwgt=100;
254         Int_t nstack=0;
255         for (j=0; j<iwgt; j++) {
256             gAlice->SetTrack(fTrackIt,-1,part,p,origin,polar,fAge,"Primary",nt);
257             gMC->Rndm(random,2);
258             phi=2*random[1]*TMath::Pi();
259             Float_t pn1=p[0]*TMath::Sin(phi) - p[1]*TMath::Cos(phi);
260             Float_t pn2=p[0]*TMath::Cos(phi) + p[1]*TMath::Sin(phi);
261             p[0]=pn1;
262             p[1]=pn2;
263             Float_t on1=origin[0]*TMath::Sin(phi)-origin[1]*TMath::Cos(phi);
264             Float_t on2=origin[0]*TMath::Cos(phi)+origin[1]*TMath::Sin(phi);
265             origin[0]=on1;
266             origin[1]=on2;
267             nstack++;
268         }
269         if (nstack == 0) continue;
270     }
271     
272     TFile *pFile=0;
273 // Get AliRun object or create it 
274     if (!gAlice) {
275         gAlice = (AliRun*)pFile->Get("gAlice");
276         if (gAlice) printf("AliRun object found on file\n");
277         if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
278     }
279     TTree *fAli=gAlice->TreeK();
280     if (fAli) pFile =fAli->GetCurrentFile();
281     pFile->cd();
282 }
283
284
285 AliGenFLUKAsource& AliGenFLUKAsource::operator=(const  AliGenFLUKAsource& rhs)
286 {
287 // Assignment operator
288     return *this;
289 }
290
291
292
293
294
295