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 **************************************************************************/
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>
27 #include <TDatabasePDG.h>
34 #include "AliCollisionGeometry.h"
35 #include "AliGenSlowNucleons.h"
36 #include "AliSlowNucleonModel.h"
38 ClassImp(AliGenSlowNucleons)
40 AliGenSlowNucleons::AliGenSlowNucleons():AliGenerator(-1)
42 // Default constructor
43 fSlowNucleonModel = 0;
44 fCollisionGeometry = 0;
49 AliGenSlowNucleons::AliGenSlowNucleons(Int_t npart)
53 fName = "SlowNucleons";
54 fTitle = "Generator for gray particles in pA collisions";
57 SetNominalCmsEnergy();
61 fSlowNucleonModel = new AliSlowNucleonModel();
62 fCollisionGeometry = 0;
66 //____________________________________________________________
67 AliGenSlowNucleons::AliGenSlowNucleons(const AliGenSlowNucleons & sn):
74 //____________________________________________________________
75 AliGenSlowNucleons::~AliGenSlowNucleons()
78 delete fSlowNucleonModel;
81 void 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);
86 void AliGenSlowNucleons::Init()
91 Float_t kMass = TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
92 fMomentum = fCMS/2. * fZTarget / fATarget;
93 fBeta = fMomentum / TMath::Sqrt(kMass * kMass + fMomentum * fMomentum);
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.);
97 fCosThetaGrayHist = new TH1F("fCosThetaGrayHist", "Gray particles angle", 100, -1., 1.);
100 // non-uniform cos(theta) distribution
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)",
109 void AliGenSlowNucleons::FinishRun()
112 // Show histogram for debugging if requested.
114 TCanvas *c = new TCanvas("c","Canvas 1",400,10,600,700);
121 fCosThetaGrayHist->Draw();
126 void AliGenSlowNucleons::Generate()
129 // Generate one event
132 // Communication with Gray Particle Model
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();
141 fSlowNucleonModel->GetNumberOfSlowNucleons(fCollisionGeometry, fNgp, fNgn, fNbp, fNbn);
143 printf("Nucleons %d %d %d %d \n", fNgp, fNgn, fNbp, fNbn);
144 fDebugHist1->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), fCollisionGeometry->NwN(), 1.);
145 fDebugHist2->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), b, 1.);
146 printf("AliGenSlowNucleons: Impact parameter from Collision Geometry %f %d %d %d %d\n",
147 b, nn, nwn, nnw, nwnw);
152 Float_t p[3], theta=0;
153 Float_t origin[3] = {0., 0., 0.};
154 Float_t polar [3] = {0., 0., 0.};
158 if(fVertexSmear == kPerEvent) {
160 for (j=0; j < 3; j++) origin[j] = fVertex[j];
167 for(i = 0; i < fNgp; i++) {
168 GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p, theta);
169 if (fDebug) fCosThetaGrayHist->Fill(TMath::Cos(theta));
170 PushTrack(fTrackIt, -1, kf, p, origin, polar,
171 0., kPNoProcess, nt, 1.);
179 for(i = 0; i < fNgn; i++) {
180 GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p, theta);
181 if (fDebug) fCosThetaGrayHist->Fill(TMath::Cos(theta));
182 PushTrack(fTrackIt, -1, kf, p, origin, polar,
183 0., kPNoProcess, nt, 1.);
191 for(i = 0; i < fNbp; i++) {
192 GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p, theta);
193 PushTrack(fTrackIt, -1, kf, p, origin, polar,
194 0., kPNoProcess, nt, 1.);
202 for(i = 0; i < fNbn; i++) {
203 GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p, theta);
204 PushTrack(fTrackIt, -1, kf, p, origin, polar,
205 0., kPNoProcess, nt, 1.);
213 void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T,
214 Double_t beta, Float_t* q, Float_t &theta)
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]
223 Double_t m, pmax, p, f, phi;
224 TDatabasePDG * pdg = TDatabasePDG::Instance();
225 const Double_t kMassProton = pdg->GetParticle(kProton) ->Mass();
226 const Double_t kMassNeutron = pdg->GetParticle(kNeutron)->Mass();
228 /* Select nucleon type */
229 if(charge == 0) m = kMassNeutron;
230 else m = kMassProton;
232 /* Momentum at maximum of Maxwell-distribution */
234 pmax = TMath::Sqrt(2*T*(T+sqrt(T*T+m*m)));
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 */
243 f = Maxwell(m, p, T) / Maxwell(m , pmax, T);
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);
255 phi = 2. * TMath::Pi() * Rndm();
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);
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);
266 /* Transform to laboratory system */
267 Lorentz(m, fBeta, q);
268 q[2] *= fProtonDirection;
271 Double_t AliGenSlowNucleons::Maxwell(Double_t m, Double_t p, Double_t T)
273 /* Relativistic Maxwell-distribution */
275 ekin = TMath::Sqrt(p*p+m*m)-m;
276 return (p*p * exp(-ekin/T));
280 void AliGenSlowNucleons::Lorentz(Double_t m, Double_t beta, Float_t* q)
282 /* Lorentz transform in the direction of q[2] */
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);
290 AliGenSlowNucleons& AliGenSlowNucleons::operator=(const AliGenSlowNucleons& rhs)
292 // Assignment operator
297 void AliGenSlowNucleons::Copy(TObject&) const
302 Fatal("Copy","Not implemented!\n");