]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenSlowNucleons.cxx
Moving the functions used to initialize TF1 and TF2 to the pivate part of the class
[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
88cb7938 16/* $Id$ */
c9a8628a 17
f687abfe 18//
19// Generator for slow nucleons in pA interactions.
20// Source is modelled by a relativistic Maxwell distributions.
21// This class cooparates with AliCollisionGeometry if used inside AliGenCocktail.
22// In this case the number of slow nucleons is determined from the number of wounded nuclei
23// using a realisation of AliSlowNucleonModel.
24// Original code by Ferenc Sikler <sikler@rmki.kfki.hu>
25//
c9a8628a 26
27#include <TDatabasePDG.h>
28#include <TPDGCode.h>
29#include <TH2F.h>
c8f7f6f9 30#include <TCanvas.h>
c9a8628a 31
32#include "AliCollisionGeometry.h"
33#include "AliGenSlowNucleons.h"
34#include "AliSlowNucleonModel.h"
35
36 ClassImp(AliGenSlowNucleons)
37
38 AliGenSlowNucleons::AliGenSlowNucleons():AliGenerator(-1)
39{
40// Default constructor
41 fSlowNucleonModel = 0;
42 fCollisionGeometry = 0;
43}
44
45AliGenSlowNucleons::AliGenSlowNucleons(Int_t npart)
46 :AliGenerator(npart)
47{
48// Constructor
49 fName = "SlowNucleons";
50 fTitle = "Generator for gray particles in pA collisions";
51 SetPmax();
52 SetTarget();
53 SetNominalCmsEnergy();
54 SetCharge();
55 SetTemperature();
56 SetBetaSource();
57 fSlowNucleonModel = new AliSlowNucleonModel();
58 fCollisionGeometry = 0;
59 fDebug = 0;
60}
61
f687abfe 62AliGenSlowNucleons::AliGenSlowNucleons(const AliGenSlowNucleons & sn):
63 AliGenerator(sn)
64{
65// Copy constructor
66 sn.Copy(*this);
67}
68
c9a8628a 69//____________________________________________________________
70AliGenSlowNucleons::~AliGenSlowNucleons()
71{
72// Destructor
73 delete fSlowNucleonModel;
74}
75
83b3cd6c 76void AliGenSlowNucleons::SetProtonDirection(Float_t dir) {
77// Set direction of the proton to change between pA (1) and Ap (-1)
78 fProtonDirection = dir / TMath::Abs(dir);
79}
c9a8628a 80
81void AliGenSlowNucleons::Init()
82{
83 //
84 // Initialization
85 //
86 Float_t kMass = TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
87 fMomentum = fCMS/2. * fZTarget / fATarget;
88 fBeta = fMomentum / TMath::Sqrt(kMass * kMass + fMomentum * fMomentum);
89 if (fDebug) {
c8f7f6f9 90 fDebugHist1 = new TH2F("DebugHist1", "nu vs N_slow", 100, 0., 100., 20, 0., 20.);
91 fDebugHist2 = new TH2F("DebugHist2", "b vs N_slow", 100, 0., 100., 15, 0., 15.);
c9a8628a 92 }
93}
94
95void AliGenSlowNucleons::FinishRun()
96{
f687abfe 97// End of run action
98// Show histogram for debugging if requested.
c9a8628a 99 if (fDebug) {
c8f7f6f9 100 TCanvas *c = new TCanvas("c","Canvas 1",400,10,600,700);
101 c->Divide(2,1);
102 c->cd(1);
103 fDebugHist1->Draw();
104 c->cd(2);
105 fDebugHist2->Draw();
c9a8628a 106 }
107}
108
109
110void AliGenSlowNucleons::Generate()
111{
112 //
113 // Generate one event
114 //
115 //
116 // Communication with Gray Particle Model
117 //
1c22a81e 118 if (fCollisionGeometry) {
119 Float_t b = fCollisionGeometry->ImpactParameter();
120 Int_t nn = fCollisionGeometry->NN();
121 Int_t nwn = fCollisionGeometry->NwN();
122 Int_t nnw = fCollisionGeometry->NNw();
123 Int_t nwnw = fCollisionGeometry->NwNw();
124
99f5114f 125 fSlowNucleonModel->GetNumberOfSlowNucleons(fCollisionGeometry, fNgp, fNgn, fNbp, fNbn);
1c22a81e 126 if (fDebug) {
127 printf("Nucleons %d %d %d %d \n", fNgp, fNgn, fNbp, fNbn);
c8f7f6f9 128 fDebugHist1->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), fCollisionGeometry->NwN(), 1.);
129 fDebugHist2->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), b, 1.);
1c22a81e 130 printf("AliGenSlowNucleons: Impact parameter from Collision Geometry %f %d %d %d %d\n",
131 b, nn, nwn, nnw, nwnw);
132 }
133 }
c9a8628a 134
135 //
136 Float_t p[3];
137 Float_t origin[3] = {0., 0., 0.};
138 Float_t polar [3] = {0., 0., 0.};
99f5114f 139 Int_t nt, i, j;
c9a8628a 140 Int_t kf;
99f5114f 141
142 if(fVertexSmear == kPerEvent) {
143 Vertex();
144 for (j=0; j < 3; j++) origin[j] = fVertex[j];
145 } // if kPerEvent
c9a8628a 146//
147// Gray protons
148//
149 fCharge = 1;
150 kf = kProton;
1c22a81e 151 for(i = 0; i < fNgp; i++) {
c9a8628a 152 GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p);
642f15cf 153 PushTrack(fTrackIt, -1, kf, p, origin, polar,
c9a8628a 154 0., kPNoProcess, nt, 1.);
155 KeepTrack(nt);
156 }
157//
158// Gray neutrons
159//
160 fCharge = 0;
161 kf = kNeutron;
1c22a81e 162 for(i = 0; i < fNgn; i++) {
c9a8628a 163 GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p);
642f15cf 164 PushTrack(fTrackIt, -1, kf, p, origin, polar,
c9a8628a 165 0., kPNoProcess, nt, 1.);
166 KeepTrack(nt);
167 }
168//
169// Black protons
170//
171 fCharge = 1;
172 kf = kProton;
1c22a81e 173 for(i = 0; i < fNbp; i++) {
c9a8628a 174 GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p);
642f15cf 175 PushTrack(fTrackIt, -1, kf, p, origin, polar,
c9a8628a 176 0., kPNoProcess, nt, 1.);
177 KeepTrack(nt);
178 }
179//
180// Black neutrons
181//
182 fCharge = 0;
183 kf = kNeutron;
1c22a81e 184 for(i = 0; i < fNbn; i++) {
c9a8628a 185 GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p);
642f15cf 186 PushTrack(fTrackIt, -1, kf, p, origin, polar,
c9a8628a 187 0., kPNoProcess, nt, 1.);
188 KeepTrack(nt);
189 }
190}
191
192
193
194
195void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T, Double_t beta, Float_t* q)
f687abfe 196
197{
c9a8628a 198/*
199 Emit a slow nucleon with "temperature" T [GeV],
200 from a source moving with velocity beta
201 Three-momentum [GeV/c] is given back in q[3]
202*/
203
c9a8628a 204 Double_t m, pmax, p, f, theta, phi;
c9a8628a 205 TDatabasePDG * pdg = TDatabasePDG::Instance();
206 const Double_t kMassProton = pdg->GetParticle(kProton) ->Mass();
207 const Double_t kMassNeutron = pdg->GetParticle(kNeutron)->Mass();
208
209 /* Select nucleon type */
210 if(charge == 0) m = kMassNeutron;
211 else m = kMassProton;
212
213 /* Momentum at maximum of Maxwell-distribution */
214
215 pmax = TMath::Sqrt(2*T*(T+sqrt(T*T+m*m)));
216
217 /* Try until proper momentum */
218 /* for lack of primitive function of the Maxwell-distribution */
219 /* a brute force trial-accept loop, normalized at pmax */
220
221 do
222 {
223 p = Rndm() * fPmax;
224 f = Maxwell(m, p, T) / Maxwell(m , pmax, T);
225 }
226 while(f < Rndm());
227
228 /* Spherical symmetric emission */
229 theta = TMath::ACos(2. * Rndm() - 1.);
230 phi = 2. * TMath::Pi() * Rndm();
231
232 /* Determine momentum components in system of the moving source */
233 q[0] = p * TMath::Sin(theta) * TMath::Cos(phi);
234 q[1] = p * TMath::Sin(theta) * TMath::Sin(phi);
235 q[2] = p * TMath::Cos(theta);
236
237 /* Transform to system of the target nucleus */
238 /* beta is passed as negative, because the gray nucleons are slowed down */
239 Lorentz(m, -beta, q);
240
241 /* Transform to laboratory system */
242 Lorentz(m, fBeta, q);
83b3cd6c 243 q[2] *= fProtonDirection;
c9a8628a 244}
245
246Double_t AliGenSlowNucleons::Maxwell(Double_t m, Double_t p, Double_t T)
247{
248/* Relativistic Maxwell-distribution */
249 Double_t ekin;
250 ekin = TMath::Sqrt(p*p+m*m)-m;
251 return (p*p * exp(-ekin/T));
252}
253
254
255void AliGenSlowNucleons::Lorentz(Double_t m, Double_t beta, Float_t* q)
256{
257/* Lorentz transform in the direction of q[2] */
258
259 Double_t gamma = 1/sqrt(1-beta*beta);
260 Double_t energy = sqrt(m*m + q[0]*q[0] + q[1]*q[1] + q[2]*q[2]);
261 q[2] = gamma * (q[2] + beta*energy);
262}
263
f687abfe 264
265AliGenSlowNucleons& AliGenSlowNucleons::operator=(const AliGenSlowNucleons& rhs)
266{
267// Assignment operator
268 rhs.Copy(*this);
269 return *this;
270}
271
dc1d768c 272void AliGenSlowNucleons::Copy(TObject&) const
f687abfe 273{
274 //
275 // Copy
276 //
277 Fatal("Copy","Not implemented!\n");
278}
279
c9a8628a 280
281
282
283
284
285