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 // 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
26 #include "AliGenFLUKAsource.h"
34 ClassImp(AliGenFLUKAsource)
36 AliGenFLUKAsource::AliGenFLUKAsource()
41 fTitle="FLUKA Boundary Source";
42 // Read in all particle types by default
44 // Set maximum admitted age of particles to 1.e-05 by default
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
54 // Set the default file
58 fTreeChain = new TChain("h1");
64 AliGenFLUKAsource::AliGenFLUKAsource(Int_t npart)
69 fTitle = "FLUKA Boundary Source";
70 // Read in all particle types by default
72 // Set maximum admitted age of particles to 1.e-05 by default
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
82 // Set the default file
86 fTreeChain = new TChain("h1");
90 AliGenFLUKAsource::AliGenFLUKAsource(const AliGenFLUKAsource & FLUKAsource):
91 AliGenerator(FLUKAsource)
94 FLUKAsource.Copy(*this);
98 //____________________________________________________________
99 AliGenFLUKAsource::~AliGenFLUKAsource()
102 if (fTreeFluka) delete fTreeFluka;
103 if (fTreeChain) delete fTreeChain;
106 void AliGenFLUKAsource::AddFile(const Text_t *filename)
108 // Add a file to the chain
109 fTreeChain->Add(filename);
114 //____________________________________________________________
115 void AliGenFLUKAsource::FlukaInit()
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);
139 //____________________________________________________________
140 void AliGenFLUKAsource::Generate()
142 // Generate one event
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,
150 Float_t polar[3]= {0,0,0};
158 Float_t amass, charge, tlife;
161 Int_t i, j, part, nt;
167 TChain *h2=fTreeChain;
168 Int_t nentries = (Int_t) h2->GetEntries();
169 if (fNpart == -1) fNpart=Int_t(nentries*fFrac);
172 // loop over number of particles
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");
181 // Get AliRun object or create it
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");
187 TTree *fAli=gAlice->TreeK();
188 if (fAli) pFile =fAli->GetCurrentFile();
190 printf("Generate - I'm out \n");
194 Int_t ifip = Int_t(fIp);
197 if (fSourceId != -1 && fIgas !=fSourceId) {
202 if (ifip > 28 || ifip < 0) {
207 if ((ifip != fIkine && fIkine != kAll && fIkine != kCharged
208 && fIkine != 10) || fAge > fAgeMax){
211 } else if (fIkine == kCharged) {
212 if (ifip == 7 || ifip == 8 || fAge > fAgeMax) {
216 } else if (fIkine == kNoNeutron) {
217 if (ifip == 8 || fAge > fAgeMax) {
226 // PDG code from FLUKA particle type (ifip)
227 part=kIfluge[int(ifip)-1];
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);
243 //handle particle weight correctly
244 wgt = (part == 13) ? fWgt*fAddWeight : fWgt;
246 fwgt=wgt-Float_t(iwgt);
248 if (random[0] < fwgt) iwgt++;
249 if (part==1 && iwgt>100) iwgt=100;
251 for (j=0; j<iwgt; j++) {
252 PushTrack(fTrackIt,-1,part,p,origin,polar,fAge,kPPrimary,nt);
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);
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);
265 if (nstack == 0) continue;
269 // Get AliRun object or create it
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");
275 TTree *fAli=gAlice->TreeK();
276 if (fAli) pFile =fAli->GetCurrentFile();
281 AliGenFLUKAsource& AliGenFLUKAsource::operator=(const AliGenFLUKAsource& rhs)
283 // Assignment operator
289 void AliGenFLUKAsource::Copy(AliGenFLUKAsource &) const
291 Fatal("Copy","Not implemented!\n");