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