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