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
24 #include "AliGenFLUKAsource.h"
27 #include <TDatabasePDG.h>
33 #include <TVirtualMC.h>
37 ClassImp(AliGenFLUKAsource)
39 AliGenFLUKAsource::AliGenFLUKAsource()
70 fTitle="FLUKA Boundary Source";
71 // Read in all particle types by default
72 // Set maximum admitted age of particles to 1.e-05 by default
74 // Shift the z-coordinate of the impact point by 4.5 cm only if it reads
75 // from specific boundary source to the chamber (fZshift=4.5;),else there
76 // is no need to shift as it reads boundary source for the whole volume of
77 // the Muon Arm; the default value corresponds to boundary source for the
78 // whole volume of the MUON Arm
79 // Set the default file
80 fTreeChain = new TChain("h1");
86 AliGenFLUKAsource::AliGenFLUKAsource(Int_t npart)
95 fTreeChain(new TChain("h1")),
117 fTitle = "FLUKA Boundary Source";
120 //____________________________________________________________
121 AliGenFLUKAsource::~AliGenFLUKAsource()
124 if (fTreeFluka) delete fTreeFluka;
125 if (fTreeChain) delete fTreeChain;
128 void AliGenFLUKAsource::AddFile(const Text_t *filename)
130 // Add a file to the chain
131 fTreeChain->Add(filename);
136 //____________________________________________________________
137 void AliGenFLUKAsource::FlukaInit()
139 // Set branch addresses of data entries
140 TChain *h2=fTreeChain;
141 //Set branch addresses
142 h2->SetBranchAddress("Ip",&fIp);
143 h2->SetBranchAddress("Ipp",&fIpp);
144 h2->SetBranchAddress("Xi",&fXi);
145 h2->SetBranchAddress("Yi",&fYi);
146 h2->SetBranchAddress("Zi",&fZi);
147 h2->SetBranchAddress("Px",&fPx);
148 h2->SetBranchAddress("Py",&fPy);
149 h2->SetBranchAddress("Pz",&fPz);
150 h2->SetBranchAddress("Ekin",&fEkin);
151 h2->SetBranchAddress("Zv",&fZv);
152 h2->SetBranchAddress("Rv",&fRv);
153 h2->SetBranchAddress("Itra",&fItra);
154 h2->SetBranchAddress("Igas",&fIgas);
155 h2->SetBranchAddress("Wgt",&fWgt);
156 h2->SetBranchAddress("Etag",&fEtag);
157 h2->SetBranchAddress("Ptg",&fPtg);
158 h2->SetBranchAddress("Age",&fAge);
161 //____________________________________________________________
162 void AliGenFLUKAsource::Generate()
164 // Generate one event
166 const Int_t kIfluge[28]={kProton, kProtonBar, kElectron, kPositron,
167 kNuE, kNuEBar, kGamma, kNeutron, kNeutronBar,
168 kMuonPlus, kMuonMinus, kK0Long , kPiPlus, kPiMinus,
169 kKPlus, kKMinus, kLambda0, kLambda0Bar, kK0Short,
170 kSigmaMinus, kSigmaPlus, kSigma0, kPi0, kK0, kK0Bar,
172 Float_t polar[3]= {0,0,0};
181 Int_t i, j, part, nt;
187 TChain *h2=fTreeChain;
188 Int_t nentries = (Int_t) h2->GetEntries();
189 if (fNpart == -1) fNpart=Int_t(nentries*fFrac);
192 // loop over number of particles
194 Int_t ev=gMC->CurrentEvent();
195 for (i=0; i<fNpart;i++) {
196 Int_t entry=fNpart*(ev-1)+i;
197 nb = (Int_t)h2->GetEvent(entry);
198 if (irwn > nentries) {
199 printf("No more entries in the FLUKA boundary source file\n");
201 // Get AliRun object or create it
203 gAlice = (AliRun*)pFile->Get("gAlice");
204 if (gAlice) printf("AliRun object found on file\n");
205 if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
207 TTree *fAli=gAlice->TreeK();
208 if (fAli) pFile =fAli->GetCurrentFile();
210 printf("Generate - I'm out \n");
214 Int_t ifip = Int_t(fIp);
217 if (fSourceId != -1 && fIgas !=fSourceId) {
222 if (ifip > 28 || ifip < 0) {
227 if ((ifip != fIkine && fIkine != kAll && fIkine != kCharged
228 && fIkine != 10) || fAge > fAgeMax){
231 } else if (fIkine == kCharged) {
232 if (ifip == 7 || ifip == 8 || fAge > fAgeMax) {
236 } else if (fIkine == kNoNeutron) {
237 if (ifip == 8 || fAge > fAgeMax) {
246 // PDG code from FLUKA particle type (ifip)
247 part=kIfluge[int(ifip)-1];
249 // Calculate momentum from kinetic energy and mass of the particle
250 #if ROOT_VERSION_CODE > 197895
251 amass = gMC->ParticleMass(part);
253 amass = (TDatabasePDG::Instance())->GetParticle(part)->Mass();
255 prwn=fEkin*sqrt(1. + 2.*amass/fEkin);
266 //handle particle weight correctly
267 wgt = (part == 13) ? fWgt*fAddWeight : fWgt;
269 fwgt=wgt-Float_t(iwgt);
271 if (random[0] < fwgt) iwgt++;
272 if (part==1 && iwgt>100) iwgt=100;
274 for (j=0; j<iwgt; j++) {
275 PushTrack(fTrackIt,-1,part,p,origin,polar,fAge,kPPrimary,nt);
277 phi=2*random[1]*TMath::Pi();
278 Float_t pn1=p[0]*TMath::Sin(phi) - p[1]*TMath::Cos(phi);
279 Float_t pn2=p[0]*TMath::Cos(phi) + p[1]*TMath::Sin(phi);
282 Float_t on1=origin[0]*TMath::Sin(phi)-origin[1]*TMath::Cos(phi);
283 Float_t on2=origin[0]*TMath::Cos(phi)+origin[1]*TMath::Sin(phi);
288 if (nstack == 0) continue;
292 // Get AliRun object or create it
294 gAlice = (AliRun*)pFile->Get("gAlice");
295 if (gAlice) printf("AliRun object found on file\n");
296 if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
298 TTree *fAli=gAlice->TreeK();
299 if (fAli) pFile =fAli->GetCurrentFile();