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 **************************************************************************/
20 #include <TDatabasePDG.h>
21 #include <TLorentzVector.h>
22 #include <TMCProcess.h>
28 #include "AliGenZDC.h"
34 //_____________________________________________________________________________
35 AliGenZDC::AliGenZDC()
39 // Default constructor
44 //_____________________________________________________________________________
45 AliGenZDC::AliGenZDC(Int_t npart)
49 // Standard constructor
52 fTitle = "Generation of Test Particles for ZDCs";
60 // LHC values for beam divergence and crossing angle
62 fBeamCrossAngle = 0.0001;
70 for(j=0; j<3; j++) fPp[i] = 0;
74 //_____________________________________________________________________________
75 void AliGenZDC::Init()
77 printf("\n\n AliGenZDC initialization:\n");
78 printf(" Particle: %d, Track cosines: x = %f, y = %f, z = %f \n",
79 fIpart,fCosx,fCosy,fCosz);
80 printf(" Fermi flag = %d, Beam divergence = %f, Crossing angle "
81 "= %f, Crossing plane = %d\n\n", fFermiflag, fBeamDiv, fBeamCrossAngle,
84 //Initialize Fermi momentum distributions for Pb-Pb
85 FermiTwoGaussian(207.,fPp,fProbintp,fProbintn);
88 //_____________________________________________________________________________
89 void AliGenZDC::Generate()
92 // Generate one trigger (n or p)
96 Double_t Mass, pLab[3], fP0, fP[3], fBoostP[3], ddp[3], dddp0, dddp[3];
97 Float_t fPTrack[3], ptot = fPMin;
100 if(fPseudoRapidity==0.){
101 pLab[0] = ptot*fCosx;
102 pLab[1] = ptot*fCosy;
103 pLab[2] = ptot*fCosz;
106 Float_t scang = 2*TMath::ATan(TMath::Exp(-(fPseudoRapidity)));
107 pLab[0] = -ptot*TMath::Sin(scang);
109 pLab[2] = ptot*TMath::Cos(scang);
111 for(i=0; i<=2; i++) fP[i] = pLab[i];
113 printf("\n\n Particle momentum before divergence and crossing\n");
114 for(i=0; i<=2; i++)printf(" pLab[%d] = %f\n",i,pLab[i]);
117 // Beam divergence and crossing angle
118 if(fBeamCrossAngle!=0.) {
119 BeamDivCross(1,fBeamDiv,fBeamCrossAngle,fBeamCrossPlane,pLab);
120 for(i=0; i<=2; i++) fP[i] = pLab[i];
123 BeamDivCross(0,fBeamDiv,fBeamCrossAngle,fBeamCrossPlane,pLab);
124 for(i=0; i<=2; i++) fP[i] = pLab[i];
127 // If required apply the Fermi momentum
129 if((fIpart==kProton) || (fIpart==kNeutron))
130 ExtractFermi(fIpart,fPp,fProbintp,fProbintn,ddp);
131 Mass=gAlice->PDGDB()->GetParticle(fIpart)->Mass();
132 fP0 = TMath::Sqrt(fP[0]*fP[0]+fP[1]*fP[1]+fP[2]*fP[2]+Mass*Mass);
133 for(i=0; i<=2; i++) dddp[i] = ddp[i];
134 dddp0 = TMath::Sqrt(dddp[0]*dddp[0]+dddp[1]*dddp[1]+dddp[2]*dddp[2]+Mass*Mass);
136 TVector3 b(fP[0]/fP0, fP[1]/fP0, fP[2]/fP0);
137 TLorentzVector pFermi(dddp[0], dddp[1], dddp[2], dddp0);
141 fBoostP[i] = pFermi[i];
147 for(i=0; i<=2; i++) fPTrack[i] = fP[i];
149 Float_t polar[3] = {0,0,0};
150 gAlice->GetMCApp()->PushTrack(fTrackIt,-1,fIpart,fPTrack,fOrigin.GetArray(),polar,0,
152 // -----------------------------------------------------------------------
154 printf("\n\n Track momentum:\n");
155 printf("\n fPTrack = %f, %f, %f \n",fPTrack[0],fPTrack[1],fPTrack[2]);
157 else if(fDebugOpt == 2){
159 if((file = fopen("SpectMomentum.dat","a")) == NULL){
160 printf("Cannot open file SpectMomentum.dat\n");
163 fprintf(file," %f \t %f \t %f \n",fPTrack[0],fPTrack[1],fPTrack[2]);
169 //_____________________________________________________________________________
170 void AliGenZDC::FermiTwoGaussian(Float_t A, Double_t *fPp,
171 Double_t *fProbintp, Double_t *fProbintn)
174 // Momenta distributions according to the "double-gaussian"
175 // distribution (Ilinov) - equal for protons and neutrons
180 Double_t sig1 = 0.113;
181 Double_t sig2 = 0.250;
182 Double_t alfa = 0.18*(TMath::Power((A/12.),(Float_t)1/3));
183 Double_t xk = (2*k2PI)/((1.+alfa)*(TMath::Power(k2PI,1.5)));
185 for(Int_t i=1; i<=200; i++){
186 Double_t p = i*0.005;
188 Double_t e1 = (p*p)/(2.*sig1*sig1);
189 Double_t e2 = (p*p)/(2.*sig2*sig2);
190 Double_t f1 = TMath::Exp(-(e1));
191 Double_t f2 = TMath::Exp(-(e2));
192 Double_t probp = xk*p*p*(f1/(TMath::Power(sig1,3.))+
193 alfa*f2/(TMath::Power(sig2,3.)))*0.005;
194 fProbintp[i] = fProbintp[i-1] + probp;
195 fProbintn[i] = fProbintp[i];
198 printf("\n\n Initialization of Fermi momenta distribution \n");
199 //for(Int_t i=0; i<=200; i++)
200 // printf(" fProbintp[%d] = %f, fProbintn[%d] = %f\n",i,fProbintp[i],i,fProbintn[i]);
203 //_____________________________________________________________________________
204 void AliGenZDC::ExtractFermi(Int_t id, Double_t *fPp, Double_t *fProbintp,
205 Double_t *fProbintn, Double_t *ddp)
208 // Compute Fermi momentum for spectator nucleons
212 Float_t xx = gRandom->Rndm();
213 assert ( id==kProton || id==kNeutron );
215 for(i=1; i<=200; i++){
216 if((xx>=fProbintp[i-1]) && (xx<fProbintp[i])) break;
220 for(i=0; i<=200; i++){
221 if((xx>=fProbintn[i-1]) && (xx<fProbintn[i])) break;
224 Float_t pext = fPp[i]+0.001;
225 Float_t phi = k2PI*(gRandom->Rndm());
226 Float_t cost = (1.-2.*(gRandom->Rndm()));
227 Float_t tet = TMath::ACos(cost);
228 ddp[0] = pext*TMath::Sin(tet)*TMath::Cos(phi);
229 ddp[1] = pext*TMath::Sin(tet)*TMath::Sin(phi);
233 printf("\n\n Extraction of Fermi momentum\n");
234 printf("\n pxFermi = %f pyFermi = %f pzFermi = %f \n",ddp[0],ddp[1],ddp[2]);
238 //_____________________________________________________________________________
239 void AliGenZDC::BeamDivCross(Int_t icross, Float_t fBeamDiv, Float_t fBeamCrossAngle,
240 Int_t fBeamCrossPlane, Double_t *pLab)
242 Double_t tetpart, fipart, tetdiv=0, fidiv=0, angleSum[2], tetsum, fisum;
247 for(i=0; i<=2; i++) pmq = pmq+pLab[i]*pLab[i];
248 Double_t pmod = TMath::Sqrt(pmq);
250 if(icross==0){ // ##### Beam divergence
251 rvec = gRandom->Gaus(0.0,1.0);
252 tetdiv = fBeamDiv * TMath::Abs(rvec);
253 fidiv = (gRandom->Rndm())*k2PI;
255 else if(icross==1){ // ##### Crossing angle
256 if(fBeamCrossPlane==0.){
260 else if(fBeamCrossPlane==1.){ // Horizontal crossing plane
261 tetdiv = fBeamCrossAngle;
264 else if(fBeamCrossPlane==2.){ // Vertical crossing plane
265 tetdiv = fBeamCrossAngle;
270 tetpart = TMath::ATan2(TMath::Sqrt(pLab[0]*pLab[0]+pLab[1]*pLab[1]),pLab[2]);
271 if(pLab[1]!=0. || pLab[0]!=0.) fipart = TMath::ATan2(pLab[1],pLab[0]);
273 if(fipart<0.) {fipart = fipart+k2PI;}
274 tetdiv = tetdiv*kRaddeg;
275 fidiv = fidiv*kRaddeg;
276 tetpart = tetpart*kRaddeg;
277 fipart = fipart*kRaddeg;
278 AddAngle(tetpart,fipart,tetdiv,fidiv,angleSum);
279 tetsum = angleSum[0];
281 tetsum = tetsum*kDegrad;
282 fisum = fisum*kDegrad;
283 pLab[0] = pmod*TMath::Sin(tetsum)*TMath::Cos(fisum);
284 pLab[1] = pmod*TMath::Sin(tetsum)*TMath::Sin(fisum);
285 pLab[2] = pmod*TMath::Cos(tetsum);
287 if(icross==0) printf("\n\n Beam divergence \n");
288 else printf("\n\n Beam crossing \n");
289 for(i=0; i<=2; i++)printf(" pLab[%d] = %f\n",i,pLab[i]);
293 //_____________________________________________________________________________
294 void AliGenZDC::AddAngle(Double_t theta1, Double_t phi1, Double_t theta2,
295 Double_t phi2, Double_t *angleSum)
297 Double_t temp, conv, cx, cy, cz, ct1, st1, ct2, st2, cp1, sp1, cp2, sp2;
298 Double_t rtetsum, tetsum, fisum;
301 conv = 180./TMath::ACos(temp);
303 ct1 = TMath::Cos(theta1/conv);
304 st1 = TMath::Sin(theta1/conv);
305 cp1 = TMath::Cos(phi1/conv);
306 sp1 = TMath::Sin(phi1/conv);
307 ct2 = TMath::Cos(theta2/conv);
308 st2 = TMath::Sin(theta2/conv);
309 cp2 = TMath::Cos(phi2/conv);
310 sp2 = TMath::Sin(phi2/conv);
311 cx = ct1*cp1*st2*cp2+st1*cp1*ct2-sp1*st2*sp2;
312 cy = ct1*sp1*st2*cp2+st1*sp1*ct2+cp1*st2*sp2;
313 cz = ct1*ct2-st1*st2*cp2;
315 rtetsum = TMath::ACos(cz);
316 tetsum = conv*rtetsum;
317 if(tetsum==0. || tetsum==180.){
321 temp = cx/TMath::Sin(rtetsum);
323 if(temp<-1.) temp=-1.;
324 fisum = conv*TMath::ACos(temp);
325 if(cy<0) {fisum = 360.-fisum;}
326 angleSum[0] = tetsum;