1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.14 2001/03/21 11:28:20 morsch
19 Use enum constants for particle selection.
21 Revision 1.13 2000/12/21 16:24:06 morsch
22 Coding convention clean-up
24 Revision 1.12 2000/11/30 07:12:50 alibrary
25 Introducing new Rndm and QA classes
27 Revision 1.11 2000/06/14 15:20:40 morsch
30 Revision 1.10 2000/06/09 20:31:34 morsch
31 All coding rule violations except RS3 corrected
33 Revision 1.9 2000/03/07 13:52:54 morsch
36 Revision 1.8 2000/02/14 14:49:38 morsch
37 Correct particle type for gamma and neutrons
38 More consistent calculation of momentum from kin. energy and mass
40 Revision 1.7 1999/11/03 17:43:20 fca
41 New version from G.Martinez & A.Morsch
43 Revision 1.6 1999/09/29 09:24:12 fca
44 Introduction of the Copyright and cvs Log
50 // Read background particles from a FLUKA boundary source file
51 // This is a very special generator that works for background studies for the muon-spectrometer.
52 // The input files come from FLUKA simulations.
53 // Files can be chained.
54 // Author: andreas.morsch@cern.ch
56 #include "AliGenFLUKAsource.h"
66 ClassImp(AliGenFLUKAsource)
67 AliGenFLUKAsource::AliGenFLUKAsource()
72 fTitle="FLUKA Boundary Source";
73 // Read in all particle types by default
75 // Set maximum admitted age of particles to 1.e-05 by default
79 // Shift the z-coordinate of the impact point by 4.5 cm only if it reads
80 // from specific boundary source to the chamber (fZshift=4.5;),else there
81 // is no need to shift as it reads boundary source for the whole volume of
82 // the Muon Arm; the default value corresponds to boundary source for the
83 // whole volume of the MUON Arm
85 // Set the default file
89 fTreeChain = new TChain("h1");
95 AliGenFLUKAsource::AliGenFLUKAsource(Int_t npart)
100 fTitle = "FLUKA Boundary Source";
101 // Read in all particle types by default
103 // Set maximum admitted age of particles to 1.e-05 by default
107 // Shift the z-coordinate of the impact point by 4.5 cm only if it reads
108 // from specific boundary source to the chamber (fZshift=4.5;),else there
109 // is no need to shift as it reads boundary source for the whole volume of
110 // the Muon Arm; the default value corresponds to boundary source for the
111 // whole volume of the MUON Arm
113 // Set the default file
117 fTreeChain = new TChain("h1");
121 AliGenFLUKAsource::AliGenFLUKAsource(const AliGenFLUKAsource & FLUKAsource)
127 //____________________________________________________________
128 AliGenFLUKAsource::~AliGenFLUKAsource()
131 if (fTreeFluka) delete fTreeFluka;
132 if (fTreeChain) delete fTreeChain;
135 void AliGenFLUKAsource::AddFile(const Text_t *filename)
137 // Add a file to the chain
138 fTreeChain->Add(filename);
143 //____________________________________________________________
144 void AliGenFLUKAsource::FlukaInit()
146 // Set branch addresses of data entries
147 TChain *h2=fTreeChain;
148 //Set branch addresses
149 h2->SetBranchAddress("Ip",&fIp);
150 h2->SetBranchAddress("Ipp",&fIpp);
151 h2->SetBranchAddress("Xi",&fXi);
152 h2->SetBranchAddress("Yi",&fYi);
153 h2->SetBranchAddress("Zi",&fZi);
154 h2->SetBranchAddress("Px",&fPx);
155 h2->SetBranchAddress("Py",&fPy);
156 h2->SetBranchAddress("Pz",&fPz);
157 h2->SetBranchAddress("Ekin",&fEkin);
158 h2->SetBranchAddress("Zv",&fZv);
159 h2->SetBranchAddress("Rv",&fRv);
160 h2->SetBranchAddress("Itra",&fItra);
161 h2->SetBranchAddress("Igas",&fIgas);
162 h2->SetBranchAddress("Wgt",&fWgt);
163 h2->SetBranchAddress("Etag",&fEtag);
164 h2->SetBranchAddress("Ptg",&fPtg);
165 h2->SetBranchAddress("Age",&fAge);
168 //____________________________________________________________
169 void AliGenFLUKAsource::Generate()
171 // Generate one event
172 AliMC* gMC = AliMC::GetMC();
174 const Int_t kIfluge[28]={kProton, kProtonBar, kElectron, kPositron,
175 kNuE, kNuEBar, kGamma, kNeutron, kNeutronBar,
176 kMuonPlus, kMuonMinus, kK0Long , kPiPlus, kPiMinus,
177 kKPlus, kKMinus, kLambda0, kLambda0Bar, kK0Short,
178 kSigmaMinus, kSigmaPlus, kSigma0, kPi0, kK0, kK0Bar,
180 Float_t polar[3]= {0,0,0};
188 Float_t amass, charge, tlife;
191 Int_t i, j, part, nt;
197 TChain *h2=fTreeChain;
198 Int_t nentries = (Int_t) h2->GetEntries();
199 if (fNpart == -1) fNpart=Int_t(nentries*fFrac);
202 // loop over number of particles
204 Int_t ev=gMC->CurrentEvent();
205 for (i=0; i<fNpart;i++) {
206 Int_t entry=fNpart*(ev-1)+i;
207 nb = (Int_t)h2->GetEvent(entry);
208 if (irwn > nentries) {
209 printf("No more entries in the FLUKA boundary source file\n");
211 // Get AliRun object or create it
213 gAlice = (AliRun*)pFile->Get("gAlice");
214 if (gAlice) printf("AliRun object found on file\n");
215 if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
217 TTree *fAli=gAlice->TreeK();
218 if (fAli) pFile =fAli->GetCurrentFile();
220 printf("Generate - I'm out \n");
224 Int_t ifip = Int_t(fIp);
227 if (fSourceId != -1 && fIgas !=fSourceId) {
232 if (ifip > 28 || ifip < 0) {
237 if ((ifip != fIkine && fIkine != kAll && fIkine != kCharged
238 && fIkine != 10) || fAge > fAgeMax){
241 } else if (fIkine == kCharged) {
242 if (ifip == 7 || ifip == 8 || fAge > fAgeMax) {
246 } else if (fIkine == kNoNeutron) {
247 if (ifip == 8 || fAge > fAgeMax) {
256 // PDG code from FLUKA particle type (ifip)
257 part=kIfluge[int(ifip)-1];
259 // Calculate momentum from kinetic energy and mass of the particle
260 gMC->Gfpart(part, name, itrtyp,
261 amass, charge, tlife);
262 prwn=fEkin*sqrt(1. + 2.*amass/fEkin);
273 //handle particle weight correctly
274 wgt = (part == 13) ? fWgt*fAddWeight : fWgt;
276 fwgt=wgt-Float_t(iwgt);
278 if (random[0] < fwgt) iwgt++;
279 if (part==1 && iwgt>100) iwgt=100;
281 for (j=0; j<iwgt; j++) {
282 SetTrack(fTrackIt,-1,part,p,origin,polar,fAge,kPPrimary,nt);
284 phi=2*random[1]*TMath::Pi();
285 Float_t pn1=p[0]*TMath::Sin(phi) - p[1]*TMath::Cos(phi);
286 Float_t pn2=p[0]*TMath::Cos(phi) + p[1]*TMath::Sin(phi);
289 Float_t on1=origin[0]*TMath::Sin(phi)-origin[1]*TMath::Cos(phi);
290 Float_t on2=origin[0]*TMath::Cos(phi)+origin[1]*TMath::Sin(phi);
295 if (nstack == 0) continue;
299 // Get AliRun object or create it
301 gAlice = (AliRun*)pFile->Get("gAlice");
302 if (gAlice) printf("AliRun object found on file\n");
303 if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
305 TTree *fAli=gAlice->TreeK();
306 if (fAli) pFile =fAli->GetCurrentFile();
311 AliGenFLUKAsource& AliGenFLUKAsource::operator=(const AliGenFLUKAsource& rhs)
313 // Assignment operator