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 <TDatabasePDG.h>
32 #include <TVirtualMC.h>
34 #include "AliGenFLUKAsource.h"
37 ClassImp(AliGenFLUKAsource)
39 AliGenFLUKAsource::AliGenFLUKAsource()
44 fTitle="FLUKA Boundary Source";
45 // Read in all particle types by default
47 // Set maximum admitted age of particles to 1.e-05 by default
51 // Shift the z-coordinate of the impact point by 4.5 cm only if it reads
52 // from specific boundary source to the chamber (fZshift=4.5;),else there
53 // is no need to shift as it reads boundary source for the whole volume of
54 // the Muon Arm; the default value corresponds to boundary source for the
55 // whole volume of the MUON Arm
57 // Set the default file
61 fTreeChain = new TChain("h1");
67 AliGenFLUKAsource::AliGenFLUKAsource(Int_t npart)
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");
93 AliGenFLUKAsource::AliGenFLUKAsource(const AliGenFLUKAsource & FLUKAsource):
94 AliGenerator(FLUKAsource)
97 FLUKAsource.Copy(*this);
101 //____________________________________________________________
102 AliGenFLUKAsource::~AliGenFLUKAsource()
105 if (fTreeFluka) delete fTreeFluka;
106 if (fTreeChain) delete fTreeChain;
109 void AliGenFLUKAsource::AddFile(const Text_t *filename)
111 // Add a file to the chain
112 fTreeChain->Add(filename);
117 //____________________________________________________________
118 void AliGenFLUKAsource::FlukaInit()
120 // Set branch addresses of data entries
121 TChain *h2=fTreeChain;
122 //Set branch addresses
123 h2->SetBranchAddress("Ip",&fIp);
124 h2->SetBranchAddress("Ipp",&fIpp);
125 h2->SetBranchAddress("Xi",&fXi);
126 h2->SetBranchAddress("Yi",&fYi);
127 h2->SetBranchAddress("Zi",&fZi);
128 h2->SetBranchAddress("Px",&fPx);
129 h2->SetBranchAddress("Py",&fPy);
130 h2->SetBranchAddress("Pz",&fPz);
131 h2->SetBranchAddress("Ekin",&fEkin);
132 h2->SetBranchAddress("Zv",&fZv);
133 h2->SetBranchAddress("Rv",&fRv);
134 h2->SetBranchAddress("Itra",&fItra);
135 h2->SetBranchAddress("Igas",&fIgas);
136 h2->SetBranchAddress("Wgt",&fWgt);
137 h2->SetBranchAddress("Etag",&fEtag);
138 h2->SetBranchAddress("Ptg",&fPtg);
139 h2->SetBranchAddress("Age",&fAge);
142 //____________________________________________________________
143 void AliGenFLUKAsource::Generate()
145 // Generate one event
147 const Int_t kIfluge[28]={kProton, kProtonBar, kElectron, kPositron,
148 kNuE, kNuEBar, kGamma, kNeutron, kNeutronBar,
149 kMuonPlus, kMuonMinus, kK0Long , kPiPlus, kPiMinus,
150 kKPlus, kKMinus, kLambda0, kLambda0Bar, kK0Short,
151 kSigmaMinus, kSigmaPlus, kSigma0, kPi0, kK0, kK0Bar,
153 Float_t polar[3]= {0,0,0};
162 Int_t i, j, part, nt;
168 TChain *h2=fTreeChain;
169 Int_t nentries = (Int_t) h2->GetEntries();
170 if (fNpart == -1) fNpart=Int_t(nentries*fFrac);
173 // loop over number of particles
175 Int_t ev=gMC->CurrentEvent();
176 for (i=0; i<fNpart;i++) {
177 Int_t entry=fNpart*(ev-1)+i;
178 nb = (Int_t)h2->GetEvent(entry);
179 if (irwn > nentries) {
180 printf("No more entries in the FLUKA boundary source file\n");
182 // Get AliRun object or create it
184 gAlice = (AliRun*)pFile->Get("gAlice");
185 if (gAlice) printf("AliRun object found on file\n");
186 if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
188 TTree *fAli=gAlice->TreeK();
189 if (fAli) pFile =fAli->GetCurrentFile();
191 printf("Generate - I'm out \n");
195 Int_t ifip = Int_t(fIp);
198 if (fSourceId != -1 && fIgas !=fSourceId) {
203 if (ifip > 28 || ifip < 0) {
208 if ((ifip != fIkine && fIkine != kAll && fIkine != kCharged
209 && fIkine != 10) || fAge > fAgeMax){
212 } else if (fIkine == kCharged) {
213 if (ifip == 7 || ifip == 8 || fAge > fAgeMax) {
217 } else if (fIkine == kNoNeutron) {
218 if (ifip == 8 || fAge > fAgeMax) {
227 // PDG code from FLUKA particle type (ifip)
228 part=kIfluge[int(ifip)-1];
230 // Calculate momentum from kinetic energy and mass of the particle
231 #if ROOT_VERSION_CODE > 197895
232 amass = gMC->ParticleMass(part);
234 amass = (TDatabasePDG::Instance())->GetParticle(part)->Mass();
236 prwn=fEkin*sqrt(1. + 2.*amass/fEkin);
247 //handle particle weight correctly
248 wgt = (part == 13) ? fWgt*fAddWeight : fWgt;
250 fwgt=wgt-Float_t(iwgt);
252 if (random[0] < fwgt) iwgt++;
253 if (part==1 && iwgt>100) iwgt=100;
255 for (j=0; j<iwgt; j++) {
256 PushTrack(fTrackIt,-1,part,p,origin,polar,fAge,kPPrimary,nt);
258 phi=2*random[1]*TMath::Pi();
259 Float_t pn1=p[0]*TMath::Sin(phi) - p[1]*TMath::Cos(phi);
260 Float_t pn2=p[0]*TMath::Cos(phi) + p[1]*TMath::Sin(phi);
263 Float_t on1=origin[0]*TMath::Sin(phi)-origin[1]*TMath::Cos(phi);
264 Float_t on2=origin[0]*TMath::Cos(phi)+origin[1]*TMath::Sin(phi);
269 if (nstack == 0) continue;
273 // Get AliRun object or create it
275 gAlice = (AliRun*)pFile->Get("gAlice");
276 if (gAlice) printf("AliRun object found on file\n");
277 if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
279 TTree *fAli=gAlice->TreeK();
280 if (fAli) pFile =fAli->GetCurrentFile();
285 AliGenFLUKAsource& AliGenFLUKAsource::operator=(const AliGenFLUKAsource& rhs)
287 // Assignment operator
293 void AliGenFLUKAsource::Copy(TObject &) const
295 Fatal("Copy","Not implemented!\n");