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
60 //_____________________________________________________________________________
61 AliGenZDC::AliGenZDC(Int_t npart)
70 fBeamCrossAngle(0.0001),
75 // Standard constructor
78 fTitle = "Generation of Test Particles for ZDCs";
80 for(Int_t i=0; i<201; i++){
87 //_____________________________________________________________________________
88 void AliGenZDC::Init()
90 //Initialize Fermi momentum distributions for Pb-Pb
92 printf("\n\n AliGenZDC initialization:\n");
93 printf(" Particle: %d, Track cosines: x = %f, y = %f, z = %f \n",
94 fIpart,fCosx,fCosy,fCosz);
95 printf(" Fermi flag = %d, Beam divergence = %f, Crossing angle "
96 "= %f, Crossing plane = %d\n\n", fFermiflag, fBeamDiv, fBeamCrossAngle,
99 FermiTwoGaussian(208.);
102 //_____________________________________________________________________________
103 void AliGenZDC::Generate()
106 // Generate one trigger (n or p)
110 Double_t mass, pLab[3], fP0, fP[3], fBoostP[3], ddp[3], dddp0, dddp[3];
111 Float_t fPTrack[3], ptot = fPMin;
114 if(fPseudoRapidity==0.){
115 pLab[0] = ptot*fCosx;
116 pLab[1] = ptot*fCosy;
117 pLab[2] = ptot*fCosz;
120 Float_t scang = 2*TMath::ATan(TMath::Exp(-(fPseudoRapidity)));
121 pLab[0] = -ptot*TMath::Sin(scang);
123 pLab[2] = ptot*TMath::Cos(scang);
125 for(i=0; i<=2; i++) fP[i] = pLab[i];
127 printf("\n\n Particle momentum before divergence and crossing\n");
128 for(i=0; i<=2; i++)printf(" pLab[%d] = %f\n",i,pLab[i]);
131 // Beam divergence and crossing angle
132 if(fBeamCrossAngle!=0.) {
133 BeamDivCross(1, pLab);
134 for(i=0; i<=2; i++) fP[i] = pLab[i];
137 BeamDivCross(0, pLab);
138 for(i=0; i<=2; i++) fP[i] = pLab[i];
141 // If required apply the Fermi momentum
143 if((fIpart==kProton) || (fIpart==kNeutron))
144 ExtractFermi(fIpart, ddp);
145 mass=TDatabasePDG::Instance()->GetParticle(fIpart)->Mass();
146 fP0 = TMath::Sqrt(fP[0]*fP[0]+fP[1]*fP[1]+fP[2]*fP[2]+mass*mass);
147 for(i=0; i<=2; i++) dddp[i] = ddp[i];
148 dddp0 = TMath::Sqrt(dddp[0]*dddp[0]+dddp[1]*dddp[1]+dddp[2]*dddp[2]+mass*mass);
150 TVector3 b(fP[0]/fP0, fP[1]/fP0, fP[2]/fP0);
151 TLorentzVector pFermi(dddp[0], dddp[1], dddp[2], dddp0);
155 fBoostP[i] = pFermi[i];
161 for(i=0; i<=2; i++) fPTrack[i] = fP[i];
163 Float_t polar[3] = {0,0,0};
164 gAlice->GetMCApp()->PushTrack(fTrackIt,-1,fIpart,fPTrack,fOrigin.GetArray(),polar,0,
166 // -----------------------------------------------------------------------
168 printf("\n\n Track momentum:\n");
169 printf("\n fPTrack = %f, %f, %f \n",fPTrack[0],fPTrack[1],fPTrack[2]);
171 else if(fDebugOpt == 2){
173 if((file = fopen("SpectMomentum.dat","a")) == NULL){
174 printf("Cannot open file SpectMomentum.dat\n");
177 fprintf(file," %f \t %f \t %f \n",fPTrack[0],fPTrack[1],fPTrack[2]);
183 //_____________________________________________________________________________
184 void AliGenZDC::FermiTwoGaussian(Float_t A)
187 // Momenta distributions according to the "double-gaussian"
188 // distribution (Ilinov) - equal for protons and neutrons
191 Double_t sig1 = 0.113;
192 Double_t sig2 = 0.250;
193 Double_t alfa = 0.18*(TMath::Power((A/12.),(Float_t)1/3));
194 Double_t xk = (2*k2PI)/((1.+alfa)*(TMath::Power(k2PI,1.5)));
196 for(Int_t i=1; i<=200; i++){
197 Double_t p = i*0.005;
199 Double_t e1 = (p*p)/(2.*sig1*sig1);
200 Double_t e2 = (p*p)/(2.*sig2*sig2);
201 Double_t f1 = TMath::Exp(-(e1));
202 Double_t f2 = TMath::Exp(-(e2));
203 Double_t probp = xk*p*p*(f1/(TMath::Power(sig1,3.))+
204 alfa*f2/(TMath::Power(sig2,3.)))*0.005;
205 fProbintp[i] = fProbintp[i-1] + probp;
206 fProbintn[i] = fProbintp[i];
209 printf("\n\n Initialization of Fermi momenta distribution \n");
210 //for(Int_t i=0; i<=200; i++)
211 // printf(" fProbintp[%d] = %f, fProbintn[%d] = %f\n",i,fProbintp[i],i,fProbintn[i]);
214 //_____________________________________________________________________________
215 void AliGenZDC::ExtractFermi(Int_t id, Double_t *ddp)
218 // Compute Fermi momentum for spectator nucleons
222 Float_t xx = gRandom->Rndm();
223 assert ( id==kProton || id==kNeutron );
225 for(i=1; i<=200; i++){
226 if((xx>=fProbintp[i-1]) && (xx<fProbintp[i])) break;
230 for(i=0; i<=200; i++){
231 if((xx>=fProbintn[i-1]) && (xx<fProbintn[i])) break;
234 Float_t pext = fPp[i]+0.001;
235 Float_t phi = k2PI*(gRandom->Rndm());
236 Float_t cost = (1.-2.*(gRandom->Rndm()));
237 Float_t tet = TMath::ACos(cost);
238 ddp[0] = pext*TMath::Sin(tet)*TMath::Cos(phi);
239 ddp[1] = pext*TMath::Sin(tet)*TMath::Sin(phi);
243 printf("\n\n Extraction of Fermi momentum\n");
244 printf("\n pxFermi = %f pyFermi = %f pzFermi = %f \n",ddp[0],ddp[1],ddp[2]);
248 //_____________________________________________________________________________
249 void AliGenZDC::BeamDivCross(Int_t icross, Double_t *pLab)
251 // Applying beam divergence and crossing angle
253 Double_t tetpart, fipart, tetdiv=0, fidiv=0, angleSum[2], tetsum, fisum;
258 for(i=0; i<=2; i++) pmq = pmq+pLab[i]*pLab[i];
259 Double_t pmod = TMath::Sqrt(pmq);
261 if(icross==0){ // ##### Beam divergence
262 rvec = gRandom->Gaus(0.0,1.0);
263 tetdiv = fBeamDiv * TMath::Abs(rvec);
264 fidiv = (gRandom->Rndm())*k2PI;
266 else if(icross==1){ // ##### Crossing angle
267 if(fBeamCrossPlane==0){
271 else if(fBeamCrossPlane==1){ // Horizontal crossing plane
272 tetdiv = fBeamCrossAngle;
275 else if(fBeamCrossPlane==2){ // Vertical crossing plane
276 tetdiv = fBeamCrossAngle;
281 tetpart = TMath::ATan2(TMath::Sqrt(pLab[0]*pLab[0]+pLab[1]*pLab[1]),pLab[2]);
282 if(pLab[1]!=0. || pLab[0]!=0.) fipart = TMath::ATan2(pLab[1],pLab[0]);
284 if(fipart<0.) {fipart = fipart+k2PI;}
285 tetdiv = tetdiv*kRaddeg;
286 fidiv = fidiv*kRaddeg;
287 tetpart = tetpart*kRaddeg;
288 fipart = fipart*kRaddeg;
289 AddAngle(tetpart,fipart,tetdiv,fidiv,angleSum);
290 tetsum = angleSum[0];
292 tetsum = tetsum*kDegrad;
293 fisum = fisum*kDegrad;
294 pLab[0] = pmod*TMath::Sin(tetsum)*TMath::Cos(fisum);
295 pLab[1] = pmod*TMath::Sin(tetsum)*TMath::Sin(fisum);
296 pLab[2] = pmod*TMath::Cos(tetsum);
298 if(icross==0) printf("\n\n Beam divergence \n");
299 else printf("\n\n Beam crossing \n");
300 for(i=0; i<=2; i++)printf(" pLab[%d] = %f\n",i,pLab[i]);
304 //_____________________________________________________________________________
305 void AliGenZDC::AddAngle(Double_t theta1, Double_t phi1, Double_t theta2,
306 Double_t phi2, Double_t *angleSum)
308 // Calculating the sum of 2 angles
309 Double_t temp, conv, cx, cy, cz, ct1, st1, ct2, st2, cp1, sp1, cp2, sp2;
310 Double_t rtetsum, tetsum, fisum;
313 conv = 180./TMath::ACos(temp);
315 ct1 = TMath::Cos(theta1/conv);
316 st1 = TMath::Sin(theta1/conv);
317 cp1 = TMath::Cos(phi1/conv);
318 sp1 = TMath::Sin(phi1/conv);
319 ct2 = TMath::Cos(theta2/conv);
320 st2 = TMath::Sin(theta2/conv);
321 cp2 = TMath::Cos(phi2/conv);
322 sp2 = TMath::Sin(phi2/conv);
323 cx = ct1*cp1*st2*cp2+st1*cp1*ct2-sp1*st2*sp2;
324 cy = ct1*sp1*st2*cp2+st1*sp1*ct2+cp1*st2*sp2;
325 cz = ct1*ct2-st1*st2*cp2;
327 rtetsum = TMath::ACos(cz);
328 tetsum = conv*rtetsum;
329 if(tetsum==0. || tetsum==180.){
333 temp = cx/TMath::Sin(rtetsum);
335 if(temp<-1.) temp=-1.;
336 fisum = conv*TMath::ACos(temp);
337 if(cy<0) {fisum = 360.-fisum;}
338 angleSum[0] = tetsum;