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