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>
33 #include <TParticle.h>
36 #include "AliCollisionGeometry.h"
40 #include "AliGenSlowNucleons.h"
41 #include "AliSlowNucleonModel.h"
43 ClassImp(AliGenSlowNucleons)
46 AliGenSlowNucleons::AliGenSlowNucleons()
68 fBeamCrossingAngle(0.),
74 // Default constructor
75 fCollisionGeometry = 0;
78 AliGenSlowNucleons::AliGenSlowNucleons(Int_t npart)
100 fBeamCrossingAngle(0.),
104 fSlowNucleonModel(new AliSlowNucleonModel())
108 fName = "SlowNucleons";
109 fTitle = "Generator for gray particles in pA collisions";
110 fCollisionGeometry = 0;
113 //____________________________________________________________
114 AliGenSlowNucleons::~AliGenSlowNucleons()
117 delete fSlowNucleonModel;
120 void AliGenSlowNucleons::SetProtonDirection(Float_t dir) {
121 // Set direction of the proton to change between pA (1) and Ap (-1)
122 fProtonDirection = dir / TMath::Abs(dir);
125 void AliGenSlowNucleons::Init()
130 Double_t kMass = TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
131 fMomentum = fCMS/2. * Float_t(fZTarget) / Float_t(fATarget);
132 fBeta = fMomentum / TMath::Sqrt(kMass * kMass + fMomentum * fMomentum);
133 //printf(" fMomentum %f fBeta %1.10f\n",fMomentum, fBeta);
135 fDebugHist1 = new TH2F("DebugHist1", "nu vs N_slow", 100, 0., 100., 20, 0., 20.);
136 fDebugHist2 = new TH2F("DebugHist2", "b vs N_slow", 100, 0., 100., 15, 0., 15.);
137 fCosThetaGrayHist = new TH1F("fCosThetaGrayHist", "Gray particles angle", 100, -1., 1.);
140 // non-uniform cos(theta) distribution
142 if(fThetaDistribution != 0) {
143 fCosTheta = new TF1("fCosTheta",
144 "(2./3.14159265358979312)/(exp(2./3.14159265358979312)-exp(-2./3.14159265358979312))*exp(2.*x/3.14159265358979312)",
148 printf("\n AliGenSlowNucleons: applying crossing angle %f mrad to slow nucleons\n",fBeamCrossingAngle*1000.);
151 void AliGenSlowNucleons::FinishRun()
154 // Show histogram for debugging if requested.
156 TCanvas *c = new TCanvas("c","Canvas 1",400,10,600,700);
159 fDebugHist1->Draw("colz");
163 fCosThetaGrayHist->Draw();
168 void AliGenSlowNucleons::Generate()
171 // Generate one event
174 // Communication with Gray Particle Model
176 if (fCollisionGeometry) {
177 Float_t b = fCollisionGeometry->ImpactParameter();
178 // Int_t nn = fCollisionGeometry->NN();
179 // Int_t nwn = fCollisionGeometry->NwN();
180 // Int_t nnw = fCollisionGeometry->NNw();
181 // Int_t nwnw = fCollisionGeometry->NwNw();
184 if(fSmearMode==0) fSlowNucleonModel->GetNumberOfSlowNucleons(fCollisionGeometry, fNgp, fNgn, fNbp, fNbn);
185 // (2) Model inspired on exp. data at lower energy (Gallio-Oppedisano)
186 // --- smearing the Ncoll fron generator used as input
187 else if(fSmearMode==1) fSlowNucleonModel->GetNumberOfSlowNucleons2(fCollisionGeometry, fNgp, fNgn, fNbp, fNbn);
188 // --- smearing directly Nslow
189 else if(fSmearMode==2) fSlowNucleonModel->GetNumberOfSlowNucleons2s(fCollisionGeometry, fNgp, fNgn, fNbp, fNbn);
191 //printf("Collision Geometry %f %d %d %d %d\n", b, nn, nwn, nnw, nwnw);
192 printf("Slow nucleons: %d grayp %d grayn %d blackp %d blackn \n", fNgp, fNgn, fNbp, fNbn);
193 fDebugHist1->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), fCollisionGeometry->NNw(), 1.);
194 fDebugHist2->Fill(Float_t(fNgp + fNgn + fNbp + fNbn), b, 1.);
200 Float_t p[3] = {0., 0., 0.}, theta=0;
201 Float_t origin[3] = {0., 0., 0.};
203 Float_t polar [3] = {0., 0., 0.};
207 // Extracting 1 value per event for the divergence angle
208 Double_t rvec = gRandom->Gaus(0.0, 1.0);
209 fBeamDivEvent = fBeamDivergence * TMath::Abs(rvec);
210 printf("\n AliGenSlowNucleons: applying beam divergence %f mrad to slow nucleons\n",fBeamDivEvent*1000.);
212 if(fVertexSmear == kPerEvent) {
214 for (j=0; j < 3; j++) origin[j] = fVertex[j];
222 for(i = 0; i < fNgp; i++) {
223 GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p, theta);
224 if (fDebug) fCosThetaGrayHist->Fill(TMath::Cos(theta));
225 PushTrack(fTrackIt, -1, kf, p, origin, polar,
226 time, kPNoProcess, nt, 1.,-2);
228 SetProcessID(nt,kGrayProcess);
235 for(i = 0; i < fNgn; i++) {
236 GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p, theta);
237 if (fDebug) fCosThetaGrayHist->Fill(TMath::Cos(theta));
238 PushTrack(fTrackIt, -1, kf, p, origin, polar,
239 time, kPNoProcess, nt, 1.,-2);
241 SetProcessID(nt,kGrayProcess);
248 for(i = 0; i < fNbp; i++) {
249 GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p, theta);
250 PushTrack(fTrackIt, -1, kf, p, origin, polar,
251 time, kPNoProcess, nt, 1.,-1);
253 SetProcessID(nt,kBlackProcess);
260 for(i = 0; i < fNbn; i++) {
261 GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p, theta);
262 PushTrack(fTrackIt, -1, kf, p, origin, polar,
263 time, kPNoProcess, nt, 1.,-1);
265 SetProcessID(nt,kBlackProcess);
272 void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T,
273 Double_t beta, Float_t* q, Float_t &theta)
277 Emit a slow nucleon with "temperature" T [GeV],
278 from a source moving with velocity beta
279 Three-momentum [GeV/c] is given back in q[3]
282 //printf("Generating slow nuc. with: charge %d. temp. %1.4f, beta %f \n",charge,T,beta);
284 Double_t m, pmax, p, f, phi;
285 TDatabasePDG * pdg = TDatabasePDG::Instance();
286 const Double_t kMassProton = pdg->GetParticle(kProton) ->Mass();
287 const Double_t kMassNeutron = pdg->GetParticle(kNeutron)->Mass();
289 /* Select nucleon type */
290 if(charge == 0) m = kMassNeutron;
291 else m = kMassProton;
293 /* Momentum at maximum of Maxwell-distribution */
295 pmax = TMath::Sqrt(2*T*(T+TMath::Sqrt(T*T+m*m)));
297 /* Try until proper momentum */
298 /* for lack of primitive function of the Maxwell-distribution */
299 /* a brute force trial-accept loop, normalized at pmax */
304 f = Maxwell(m, p, T) / Maxwell(m , pmax, T);
308 /* Spherical symmetric emission for black particles (beta=0)*/
309 if(beta==0 || fThetaDistribution==0) theta = TMath::ACos(2. * Rndm() - 1.);
310 /* cos theta distributed according to experimental results for gray particles (beta=0.05)*/
311 else if(fThetaDistribution!=0){
312 Double_t costheta = fCosTheta->GetRandom();
313 theta = TMath::ACos(costheta);
316 phi = 2. * TMath::Pi() * Rndm();
319 /* Determine momentum components in system of the moving source */
320 q[0] = p * TMath::Sin(theta) * TMath::Cos(phi);
321 q[1] = p * TMath::Sin(theta) * TMath::Sin(phi);
322 q[2] = p * TMath::Cos(theta);
323 //if(fDebug==1) printf("\n Momentum in RS of the moving source: p = (%f, %f, %f)\n",q[0],q[1],q[2]);
326 /* Transform to system of the target nucleus */
327 /* beta is passed as negative, because the gray nucleons are slowed down */
328 Lorentz(m, -beta, q);
329 //if(fDebug==1) printf(" Momentum in RS of the target nucleus: p = (%f, %f, %f)\n",q[0],q[1],q[2]);
331 /* Transform to laboratory system */
332 Lorentz(m, fBeta, q);
333 q[2] *= fProtonDirection;
334 if(fDebug==1)printf("\n Momentum after LHC boost: p = (%f, %f, %f)\n",q[0],q[1],q[2]);
336 if(fBeamCrossingAngle>0.) BeamCrossDivergence(1, q); // applying crossing angle
337 if(fBeamDivergence>0.) BeamCrossDivergence(2, q); // applying divergence
341 Double_t AliGenSlowNucleons::Maxwell(Double_t m, Double_t p, Double_t T)
343 /* Relativistic Maxwell-distribution */
345 ekin = TMath::Sqrt(p*p+m*m)-m;
346 return (p*p * exp(-ekin/T));
350 //_____________________________________________________________________________
351 void AliGenSlowNucleons::Lorentz(Double_t m, Double_t beta, Float_t* q)
353 /* Lorentz transform in the direction of q[2] */
355 Double_t gamma = 1./TMath::Sqrt((1.-beta)*(1.+beta));
356 Double_t energy = TMath::Sqrt(m*m + q[0]*q[0] + q[1]*q[1] + q[2]*q[2]);
357 q[2] = gamma * (q[2] + beta*energy);
358 //printf(" \t beta %1.10f gamma %f energy %f -> p_z = %f\n",beta, gamma, energy,q[2]);
361 //_____________________________________________________________________________
362 void AliGenSlowNucleons::BeamCrossDivergence(Int_t iwhat, Float_t *pLab)
364 // Applying beam divergence and crossing angle
366 Double_t pmod = TMath::Sqrt(pLab[0]*pLab[0]+pLab[1]*pLab[1]+pLab[2]*pLab[2]);
368 Double_t tetdiv = 0.;
371 tetdiv = fBeamCrossingAngle;
375 tetdiv = fBeamDivEvent;
376 fidiv = (gRandom->Rndm())*k2PI;
379 Double_t tetpart = TMath::ATan2(TMath::Sqrt(pLab[0]*pLab[0]+pLab[1]*pLab[1]), pLab[2]);
381 if(TMath::Abs(pLab[1])>0. || TMath::Abs(pLab[0])>0.) fipart = TMath::ATan2(pLab[1], pLab[0]);
382 if(fipart<0.) {fipart = fipart+k2PI;}
383 tetdiv = tetdiv*kRaddeg;
384 fidiv = fidiv*kRaddeg;
385 tetpart = tetpart*kRaddeg;
386 fipart = fipart*kRaddeg;
388 Double_t angleSum[2]={0., 0.};
389 AddAngle(tetpart,fipart,tetdiv,fidiv,angleSum);
391 Double_t tetsum = angleSum[0];
392 Double_t fisum = angleSum[1];
393 //printf("tetpart %f fipart %f tetdiv %f fidiv %f angleSum %f %f\n",tetpart,fipart,tetdiv,fidiv,angleSum[0],angleSum[1]);
394 tetsum = tetsum*kDegrad;
395 fisum = fisum*kDegrad;
397 pLab[0] = pmod*TMath::Sin(tetsum)*TMath::Cos(fisum);
398 pLab[1] = pmod*TMath::Sin(tetsum)*TMath::Sin(fisum);
399 pLab[2] = pmod*TMath::Cos(tetsum);
401 if(iwhat==1) printf(" Beam crossing angle %f mrad ", fBeamCrossingAngle*1000.);
402 else if(iwhat==2) printf(" Beam divergence %f mrad ", fBeamDivEvent*1000.);
403 printf(" p = (%f, %f, %f)\n",pLab[0],pLab[1],pLab[2]);
407 //_____________________________________________________________________________
408 void AliGenSlowNucleons::AddAngle(Double_t theta1, Double_t phi1,
409 Double_t theta2, Double_t phi2, Double_t *angleSum)
411 // Calculating the sum of 2 angles
413 Double_t conv = 180./TMath::ACos(temp);
415 Double_t ct1 = TMath::Cos(theta1/conv);
416 Double_t st1 = TMath::Sin(theta1/conv);
417 Double_t cp1 = TMath::Cos(phi1/conv);
418 Double_t sp1 = TMath::Sin(phi1/conv);
419 Double_t ct2 = TMath::Cos(theta2/conv);
420 Double_t st2 = TMath::Sin(theta2/conv);
421 Double_t cp2 = TMath::Cos(phi2/conv);
422 Double_t sp2 = TMath::Sin(phi2/conv);
423 Double_t cx = ct1*cp1*st2*cp2+st1*cp1*ct2-sp1*st2*sp2;
424 Double_t cy = ct1*sp1*st2*cp2+st1*sp1*ct2+cp1*st2*sp2;
425 Double_t cz = ct1*ct2-st1*st2*cp2;
427 Double_t rtetsum = TMath::ACos(cz);
428 Double_t tetsum = conv*rtetsum;
429 if(TMath::Abs(tetsum)<1e-4 || tetsum==180.) return;
431 temp = cx/TMath::Sin(rtetsum);
433 if(temp<-1.) temp=-1.;
434 Double_t fisum = conv*TMath::ACos(temp);
435 if(cy<0) {fisum = 360.-fisum;}
437 angleSum[0] = tetsum;
441 //_____________________________________________________________________________
442 void AliGenSlowNucleons::SetProcessID(Int_t nt, UInt_t process)
444 // Tag the particle as
447 fStack->Particle(nt)->SetUniqueID(process);
449 gAlice->GetMCApp()->Particle(nt)->SetUniqueID(process);