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 **************************************************************************/
18 //////////////////////////////////////////////////////////////////////
20 // Generator of spectator nucleons (either protons or neutrons)//
21 // computes beam crossing and divergence and Fermi momentum //
23 /////////////////////////////////////////////////////////////////////
27 #include <TDatabasePDG.h>
28 #include <TLorentzVector.h>
29 #include <TMCProcess.h>
35 #include "AliGenZDC.h"
41 //_____________________________________________________________________________
42 AliGenZDC::AliGenZDC()
56 // Default constructor
58 for(Int_t i=0; i<201; i++){
65 //_____________________________________________________________________________
66 AliGenZDC::AliGenZDC(Int_t npart)
75 fBeamCrossAngle(0.0001),
80 // Standard constructor
83 fTitle = "Generation of Test Particles for ZDCs";
85 for(Int_t i=0; i<201; i++){
92 //_____________________________________________________________________________
93 void AliGenZDC::Init()
95 //Initialize Fermi momentum distributions for Pb-Pb
97 printf("\n\n AliGenZDC initialization:\n");
98 printf(" Particle: %d, Track cosines: x = %f, y = %f, z = %f \n",
99 fIpart,fCosx,fCosy,fCosz);
100 printf(" Fermi flag = %d, Beam divergence = %f, Crossing angle "
101 "= %f, Crossing plane = %d\n\n", fFermiflag, fBeamDiv, fBeamCrossAngle,
104 FermiTwoGaussian(208.);
107 //_____________________________________________________________________________
108 void AliGenZDC::Generate()
111 // Generate one trigger (n or p)
115 Double_t mass, pLab[3], fP0, fP[3], fBoostP[3], ddp[3]={0.,0.,0.}, dddp0, dddp[3];
116 Float_t fPTrack[3], ptot = fPMin;
119 if(fPseudoRapidity==0.){
120 pLab[0] = ptot*fCosx;
121 pLab[1] = ptot*fCosy;
122 pLab[2] = ptot*fCosz;
125 Float_t scang = 2*TMath::ATan(TMath::Exp(-(fPseudoRapidity)));
126 pLab[0] = -ptot*TMath::Sin(scang);
128 pLab[2] = ptot*TMath::Cos(scang);
130 for(i=0; i<=2; i++) fP[i] = pLab[i];
132 printf("\n\n Particle momentum before divergence and crossing\n");
133 for(i=0; i<=2; i++)printf(" pLab[%d] = %f\n",i,pLab[i]);
136 // Beam divergence and crossing angle
137 if(fBeamCrossAngle!=0.) {
138 BeamDivCross(1, pLab);
139 for(i=0; i<=2; i++) fP[i] = pLab[i];
142 BeamDivCross(0, pLab);
143 for(i=0; i<=2; i++) fP[i] = pLab[i];
146 // If required apply the Fermi momentum
148 if((fIpart==kProton) || (fIpart==kNeutron))
149 ExtractFermi(fIpart, ddp);
150 mass=TDatabasePDG::Instance()->GetParticle(fIpart)->Mass();
151 fP0 = TMath::Sqrt(fP[0]*fP[0]+fP[1]*fP[1]+fP[2]*fP[2]+mass*mass);
152 for(i=0; i<=2; i++) dddp[i] = ddp[i];
153 dddp0 = TMath::Sqrt(dddp[0]*dddp[0]+dddp[1]*dddp[1]+dddp[2]*dddp[2]+mass*mass);
155 TVector3 b(fP[0]/fP0, fP[1]/fP0, fP[2]/fP0);
156 TLorentzVector pFermi(dddp[0], dddp[1], dddp[2], dddp0);
160 fBoostP[i] = pFermi[i];
166 for(i=0; i<=2; i++) fPTrack[i] = fP[i];
168 Float_t polar[3] = {0,0,0};
169 gAlice->GetMCApp()->PushTrack(fTrackIt,-1,fIpart,fPTrack,fOrigin.GetArray(),polar,0,
171 // -----------------------------------------------------------------------
173 printf("\n\n Track momentum:\n");
174 printf("\n fPTrack = %f, %f, %f \n",fPTrack[0],fPTrack[1],fPTrack[2]);
176 else if(fDebugOpt == 2){
178 if((file = fopen("SpectMomentum.dat","a")) == NULL){
179 printf("Cannot open file SpectMomentum.dat\n");
182 fprintf(file," %f \t %f \t %f \n",fPTrack[0],fPTrack[1],fPTrack[2]);
188 //_____________________________________________________________________________
189 void AliGenZDC::FermiTwoGaussian(Float_t A)
192 // Momenta distributions according to the "double-gaussian"
193 // distribution (Ilinov) - equal for protons and neutrons
196 Double_t sig1 = 0.113;
197 Double_t sig2 = 0.250;
198 Double_t alfa = 0.18*(TMath::Power((A/12.),(Float_t)1/3));
199 Double_t xk = (2*k2PI)/((1.+alfa)*(TMath::Power(k2PI,1.5)));
201 for(Int_t i=1; i<=200; i++){
202 Double_t p = i*0.005;
204 Double_t e1 = (p*p)/(2.*sig1*sig1);
205 Double_t e2 = (p*p)/(2.*sig2*sig2);
206 Double_t f1 = TMath::Exp(-(e1));
207 Double_t f2 = TMath::Exp(-(e2));
208 Double_t probp = xk*p*p*(f1/(TMath::Power(sig1,3.))+
209 alfa*f2/(TMath::Power(sig2,3.)))*0.005;
210 fProbintp[i] = fProbintp[i-1] + probp;
211 fProbintn[i] = fProbintp[i];
214 printf("\n\n Initialization of Fermi momenta distribution \n");
215 //for(Int_t i=0; i<=200; i++)
216 // printf(" fProbintp[%d] = %f, fProbintn[%d] = %f\n",i,fProbintp[i],i,fProbintn[i]);
219 //_____________________________________________________________________________
220 void AliGenZDC::ExtractFermi(Int_t id, Double_t *ddp)
223 // Compute Fermi momentum for spectator nucleons
227 Float_t xx = gRandom->Rndm();
228 assert ( id==kProton || id==kNeutron );
230 for(Int_t i=1; i<=200; i++){
231 if((xx>=fProbintp[i-1]) && (xx<fProbintp[i])) break;
236 for(Int_t i=1; i<=200; i++){
237 if((xx>=fProbintn[i-1]) && (xx<fProbintn[i])) break;
241 Float_t pext = fPp[index]+0.001;
242 Float_t phi = k2PI*(gRandom->Rndm());
243 Float_t cost = (1.-2.*(gRandom->Rndm()));
244 Float_t tet = TMath::ACos(cost);
245 ddp[0] = pext*TMath::Sin(tet)*TMath::Cos(phi);
246 ddp[1] = pext*TMath::Sin(tet)*TMath::Sin(phi);
250 printf("\n\n Extraction of Fermi momentum\n");
251 printf("\n pxFermi = %f pyFermi = %f pzFermi = %f \n",ddp[0],ddp[1],ddp[2]);
255 //_____________________________________________________________________________
256 void AliGenZDC::BeamDivCross(Int_t icross, Double_t *pLab)
258 // Applying beam divergence and crossing angle
260 Double_t tetpart, fipart, tetdiv=0, fidiv=0, angleSum[2], tetsum, fisum;
265 for(i=0; i<=2; i++) pmq = pmq+pLab[i]*pLab[i];
266 Double_t pmod = TMath::Sqrt(pmq);
268 if(icross==0){ // ##### Beam divergence
269 rvec = gRandom->Gaus(0.0,1.0);
270 tetdiv = fBeamDiv * TMath::Abs(rvec);
271 fidiv = (gRandom->Rndm())*k2PI;
273 else if(icross==1){ // ##### Crossing angle
274 if(fBeamCrossPlane==0){
278 else if(fBeamCrossPlane==1){ // Horizontal crossing plane
279 tetdiv = fBeamCrossAngle;
282 else if(fBeamCrossPlane==2){ // Vertical crossing plane
283 tetdiv = fBeamCrossAngle;
288 tetpart = TMath::ATan2(TMath::Sqrt(pLab[0]*pLab[0]+pLab[1]*pLab[1]),pLab[2]);
289 if(pLab[1]!=0. || pLab[0]!=0.) fipart = TMath::ATan2(pLab[1],pLab[0]);
291 if(fipart<0.) {fipart = fipart+k2PI;}
292 tetdiv = tetdiv*kRaddeg;
293 fidiv = fidiv*kRaddeg;
294 tetpart = tetpart*kRaddeg;
295 fipart = fipart*kRaddeg;
296 AddAngle(tetpart,fipart,tetdiv,fidiv,angleSum);
297 tetsum = angleSum[0];
299 tetsum = tetsum*kDegrad;
300 fisum = fisum*kDegrad;
301 pLab[0] = pmod*TMath::Sin(tetsum)*TMath::Cos(fisum);
302 pLab[1] = pmod*TMath::Sin(tetsum)*TMath::Sin(fisum);
303 pLab[2] = pmod*TMath::Cos(tetsum);
305 if(icross==0) printf("\n\n Beam divergence \n");
306 else printf("\n\n Beam crossing \n");
307 for(i=0; i<=2; i++)printf(" pLab[%d] = %f\n",i,pLab[i]);
311 //_____________________________________________________________________________
312 void AliGenZDC::AddAngle(Double_t theta1, Double_t phi1, Double_t theta2,
313 Double_t phi2, Double_t *angleSum)
315 // Calculating the sum of 2 angles
316 Double_t temp, conv, cx, cy, cz, ct1, st1, ct2, st2, cp1, sp1, cp2, sp2;
317 Double_t rtetsum, tetsum, fisum;
320 conv = 180./TMath::ACos(temp);
322 ct1 = TMath::Cos(theta1/conv);
323 st1 = TMath::Sin(theta1/conv);
324 cp1 = TMath::Cos(phi1/conv);
325 sp1 = TMath::Sin(phi1/conv);
326 ct2 = TMath::Cos(theta2/conv);
327 st2 = TMath::Sin(theta2/conv);
328 cp2 = TMath::Cos(phi2/conv);
329 sp2 = TMath::Sin(phi2/conv);
330 cx = ct1*cp1*st2*cp2+st1*cp1*ct2-sp1*st2*sp2;
331 cy = ct1*sp1*st2*cp2+st1*sp1*ct2+cp1*st2*sp2;
332 cz = ct1*ct2-st1*st2*cp2;
334 rtetsum = TMath::ACos(cz);
335 tetsum = conv*rtetsum;
336 if(tetsum==0. || tetsum==180.){
340 temp = cx/TMath::Sin(rtetsum);
342 if(temp<-1.) temp=-1.;
343 fisum = conv*TMath::ACos(temp);
344 if(cy<0) {fisum = 360.-fisum;}
345 angleSum[0] = tetsum;