]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenSlowNucleons.cxx
Purity check is corrected
[u/mrichter/AliRoot.git] / EVGEN / AliGenSlowNucleons.cxx
CommitLineData
c9a8628a 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
16/*
17$Log$
18*/
19
20/*
21 Generator for slow nucluons in pA interactions.
22 Source is modelled by a relativistic Maxwell distributions.
23 Original code by Ferenc Sikler <sikler@rmki.kfki.hu>
24 */
25
26#include <TDatabasePDG.h>
27#include <TPDGCode.h>
28#include <TH2F.h>
29
30#include "AliCollisionGeometry.h"
31#include "AliGenSlowNucleons.h"
32#include "AliSlowNucleonModel.h"
33
34 ClassImp(AliGenSlowNucleons)
35
36 AliGenSlowNucleons::AliGenSlowNucleons():AliGenerator(-1)
37{
38// Default constructor
39 fSlowNucleonModel = 0;
40 fCollisionGeometry = 0;
41}
42
43AliGenSlowNucleons::AliGenSlowNucleons(Int_t npart)
44 :AliGenerator(npart)
45{
46// Constructor
47 fName = "SlowNucleons";
48 fTitle = "Generator for gray particles in pA collisions";
49 SetPmax();
50 SetTarget();
51 SetNominalCmsEnergy();
52 SetCharge();
53 SetTemperature();
54 SetBetaSource();
55 fSlowNucleonModel = new AliSlowNucleonModel();
56 fCollisionGeometry = 0;
57 fDebug = 0;
58}
59
60//____________________________________________________________
61AliGenSlowNucleons::~AliGenSlowNucleons()
62{
63// Destructor
64 delete fSlowNucleonModel;
65}
66
67
68void AliGenSlowNucleons::Init()
69{
70 //
71 // Initialization
72 //
73 Float_t kMass = TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
74 fMomentum = fCMS/2. * fZTarget / fATarget;
75 fBeta = fMomentum / TMath::Sqrt(kMass * kMass + fMomentum * fMomentum);
76 if (fDebug) {
77 fDebugHist = new TH2F("DebugHist", "N_heavy vs nu", 50, 0., 50., 50, 0., 50.);
78 }
79}
80
81void AliGenSlowNucleons::FinishRun()
82{
83 if (fDebug) {
84 fDebugHist->Draw();
85 }
86}
87
88
89void AliGenSlowNucleons::Generate()
90{
91 //
92 // Generate one event
93 //
94 //
95 // Communication with Gray Particle Model
96 //
97 Int_t ngp, ngn, nbp, nbn;
98
99 Float_t b = fCollisionGeometry->ImpactParameter();
100 Int_t nn = fCollisionGeometry->NN();
101 Int_t nwn = fCollisionGeometry->NwN();
102 Int_t nnw = fCollisionGeometry->NNw();
103 Int_t nwnw = fCollisionGeometry->NwNw();
104
105 printf("AliGenSlowNucleons: Impact parameter from Collision Geometry %f %d %d %d %d\n",
106 b, nn, nwn, nnw, nwnw);
107
108 fSlowNucleonModel->GetNumberOfSlowNucleons(fCollisionGeometry, ngp, ngn, nbp, nbn);
109
110 if (fDebug) {
111 printf("nucleons %d %d %d %d \n", ngp, ngn, nbp, nbn);
112
113 fDebugHist->Fill(Float_t(ngp + ngn + nbp + nbn), fCollisionGeometry->NwN(), 1.);
114 }
115
116
117 //
118 Float_t p[3];
119 Float_t origin[3] = {0., 0., 0.};
120 Float_t polar [3] = {0., 0., 0.};
121 Int_t nt, i;
122 Int_t kf;
123//
124// Gray protons
125//
126 fCharge = 1;
127 kf = kProton;
128 for(i = 0; i < ngp; i++) {
129 GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p);
130 SetTrack(fTrackIt, -1, kf, p, origin, polar,
131 0., kPNoProcess, nt, 1.);
132 KeepTrack(nt);
133 }
134//
135// Gray neutrons
136//
137 fCharge = 0;
138 kf = kNeutron;
139 for(i = 0; i < ngn; i++) {
140 GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p);
141 SetTrack(fTrackIt, -1, kf, p, origin, polar,
142 0., kPNoProcess, nt, 1.);
143 KeepTrack(nt);
144 }
145//
146// Black protons
147//
148 fCharge = 1;
149 kf = kProton;
150 for(i = 0; i < nbp; i++) {
151 GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p);
152 SetTrack(fTrackIt, -1, kf, p, origin, polar,
153 0., kPNoProcess, nt, 1.);
154 KeepTrack(nt);
155 }
156//
157// Black neutrons
158//
159 fCharge = 0;
160 kf = kNeutron;
161 for(i = 0; i < nbn; i++) {
162 GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p);
163 SetTrack(fTrackIt, -1, kf, p, origin, polar,
164 0., kPNoProcess, nt, 1.);
165 KeepTrack(nt);
166 }
167}
168
169
170
171
172void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T, Double_t beta, Float_t* q)
173/*
174 Emit a slow nucleon with "temperature" T [GeV],
175 from a source moving with velocity beta
176 Three-momentum [GeV/c] is given back in q[3]
177*/
178
179{
180 Double_t m, pmax, p, f, theta, phi;
181
182 TDatabasePDG * pdg = TDatabasePDG::Instance();
183 const Double_t kMassProton = pdg->GetParticle(kProton) ->Mass();
184 const Double_t kMassNeutron = pdg->GetParticle(kNeutron)->Mass();
185
186 /* Select nucleon type */
187 if(charge == 0) m = kMassNeutron;
188 else m = kMassProton;
189
190 /* Momentum at maximum of Maxwell-distribution */
191
192 pmax = TMath::Sqrt(2*T*(T+sqrt(T*T+m*m)));
193
194 /* Try until proper momentum */
195 /* for lack of primitive function of the Maxwell-distribution */
196 /* a brute force trial-accept loop, normalized at pmax */
197
198 do
199 {
200 p = Rndm() * fPmax;
201 f = Maxwell(m, p, T) / Maxwell(m , pmax, T);
202 }
203 while(f < Rndm());
204
205 /* Spherical symmetric emission */
206 theta = TMath::ACos(2. * Rndm() - 1.);
207 phi = 2. * TMath::Pi() * Rndm();
208
209 /* Determine momentum components in system of the moving source */
210 q[0] = p * TMath::Sin(theta) * TMath::Cos(phi);
211 q[1] = p * TMath::Sin(theta) * TMath::Sin(phi);
212 q[2] = p * TMath::Cos(theta);
213
214 /* Transform to system of the target nucleus */
215 /* beta is passed as negative, because the gray nucleons are slowed down */
216 Lorentz(m, -beta, q);
217
218 /* Transform to laboratory system */
219 Lorentz(m, fBeta, q);
220}
221
222Double_t AliGenSlowNucleons::Maxwell(Double_t m, Double_t p, Double_t T)
223{
224/* Relativistic Maxwell-distribution */
225 Double_t ekin;
226 ekin = TMath::Sqrt(p*p+m*m)-m;
227 return (p*p * exp(-ekin/T));
228}
229
230
231void AliGenSlowNucleons::Lorentz(Double_t m, Double_t beta, Float_t* q)
232{
233/* Lorentz transform in the direction of q[2] */
234
235 Double_t gamma = 1/sqrt(1-beta*beta);
236 Double_t energy = sqrt(m*m + q[0]*q[0] + q[1]*q[1] + q[2]*q[2]);
237 q[2] = gamma * (q[2] + beta*energy);
238}
239
240
241
242
243
244
245