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>
32 #include "AliCollisionGeometry.h"
33 #include "AliGenSlowNucleons.h"
34 #include "AliSlowNucleonModel.h"
36 ClassImp(AliGenSlowNucleons)
38 AliGenSlowNucleons::AliGenSlowNucleons():AliGenerator(-1)
40 // Default constructor
41 fSlowNucleonModel = 0;
42 fCollisionGeometry = 0;
45 AliGenSlowNucleons::AliGenSlowNucleons(Int_t npart)
49 fName = "SlowNucleons";
50 fTitle = "Generator for gray particles in pA collisions";
53 SetNominalCmsEnergy();
57 fSlowNucleonModel = new AliSlowNucleonModel();
58 fCollisionGeometry = 0;
62 AliGenSlowNucleons::AliGenSlowNucleons(const AliGenSlowNucleons & sn):
69 //____________________________________________________________
70 AliGenSlowNucleons::~AliGenSlowNucleons()
73 delete fSlowNucleonModel;
76 void 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);
81 void AliGenSlowNucleons::Init()
86 Float_t kMass = TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
87 fMomentum = fCMS/2. * fZTarget / fATarget;
88 fBeta = fMomentum / TMath::Sqrt(kMass * kMass + fMomentum * fMomentum);
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.);
95 void AliGenSlowNucleons::FinishRun()
98 // Show histogram for debugging if requested.
100 TCanvas *c = new TCanvas("c","Canvas 1",400,10,600,700);
110 void AliGenSlowNucleons::Generate()
113 // Generate one event
116 // Communication with Gray Particle Model
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();
125 fSlowNucleonModel->GetNumberOfSlowNucleons(fCollisionGeometry, fNgp, fNgn, fNbp, fNbn);
127 printf("Nucleons %d %d %d %d \n", fNgp, fNgn, fNbp, fNbn);
128 fDebugHist1->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), fCollisionGeometry->NwN(), 1.);
129 fDebugHist2->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), b, 1.);
130 printf("AliGenSlowNucleons: Impact parameter from Collision Geometry %f %d %d %d %d\n",
131 b, nn, nwn, nnw, nwnw);
137 Float_t origin[3] = {0., 0., 0.};
138 Float_t polar [3] = {0., 0., 0.};
142 if(fVertexSmear == kPerEvent) {
144 for (j=0; j < 3; j++) origin[j] = fVertex[j];
151 for(i = 0; i < fNgp; i++) {
152 GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p);
153 PushTrack(fTrackIt, -1, kf, p, origin, polar,
154 0., kPNoProcess, nt, 1.);
162 for(i = 0; i < fNgn; i++) {
163 GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p);
164 PushTrack(fTrackIt, -1, kf, p, origin, polar,
165 0., kPNoProcess, nt, 1.);
173 for(i = 0; i < fNbp; i++) {
174 GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p);
175 PushTrack(fTrackIt, -1, kf, p, origin, polar,
176 0., kPNoProcess, nt, 1.);
184 for(i = 0; i < fNbn; i++) {
185 GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p);
186 PushTrack(fTrackIt, -1, kf, p, origin, polar,
187 0., kPNoProcess, nt, 1.);
195 void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T, Double_t beta, Float_t* q)
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]
204 Double_t m, pmax, p, f, theta, phi;
205 TDatabasePDG * pdg = TDatabasePDG::Instance();
206 const Double_t kMassProton = pdg->GetParticle(kProton) ->Mass();
207 const Double_t kMassNeutron = pdg->GetParticle(kNeutron)->Mass();
209 /* Select nucleon type */
210 if(charge == 0) m = kMassNeutron;
211 else m = kMassProton;
213 /* Momentum at maximum of Maxwell-distribution */
215 pmax = TMath::Sqrt(2*T*(T+sqrt(T*T+m*m)));
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 */
224 f = Maxwell(m, p, T) / Maxwell(m , pmax, T);
228 /* Spherical symmetric emission */
229 theta = TMath::ACos(2. * Rndm() - 1.);
230 phi = 2. * TMath::Pi() * Rndm();
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);
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);
241 /* Transform to laboratory system */
242 Lorentz(m, fBeta, q);
243 q[2] *= fProtonDirection;
246 Double_t AliGenSlowNucleons::Maxwell(Double_t m, Double_t p, Double_t T)
248 /* Relativistic Maxwell-distribution */
250 ekin = TMath::Sqrt(p*p+m*m)-m;
251 return (p*p * exp(-ekin/T));
255 void AliGenSlowNucleons::Lorentz(Double_t m, Double_t beta, Float_t* q)
257 /* Lorentz transform in the direction of q[2] */
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);
265 AliGenSlowNucleons& AliGenSlowNucleons::operator=(const AliGenSlowNucleons& rhs)
267 // Assignment operator
272 void AliGenSlowNucleons::Copy(TObject&) const
277 Fatal("Copy","Not implemented!\n");