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 Revision 1.7 2000/01/19 17:17:40 fca
20 Revision 1.6 1999/09/29 09:24:35 fca
21 Introduction of the Copyright and cvs Log
25 #include <TLorentzVector.h>
28 #include "AliGenZDC.h"
35 //_____________________________________________________________________________
36 AliGenZDC::AliGenZDC()
40 // Default constructor
45 //_____________________________________________________________________________
46 AliGenZDC::AliGenZDC(Int_t npart)
50 // Standard constructor
53 fTitle = "Generation of Test Particles for ZDCs";
60 // LHC values for beam divergence and crossing angle
62 fBeamCrossAngle = 0.0001;
66 //_____________________________________________________________________________
67 void AliGenZDC::Init()
69 printf(" Initializing AliGenZDC\n");
70 printf(" Fermi flag = %d, Beam Divergence = %f, Crossing Angle "
71 "= %f, Crossing Plane = %d\n\n", fFermiflag, fBeamDiv, fBeamCrossAngle,
73 //Initialize Fermi momentum distributions for Pb-Pb
74 FermiTwoGaussian(207.,82.,fPp,fProbintp,fProbintn);
77 //_____________________________________________________________________________
78 void AliGenZDC::Generate()
81 // Generate one trigger (n or p)
83 Double_t mass, pLab[3], balp0, balp[3], ddp[3], dddp0, dddp[3];
87 if(fPseudoRapidity==0.){
93 Float_t scang = 2*TMath::ATan(TMath::Exp(-(fPseudoRapidity)));
94 pLab[0] = -ptot*TMath::Sin(scang);
96 pLab[2] = ptot*TMath::Cos(scang);
98 for(Int_t i=0; i<=2; i++){
102 // Beam divergence and crossing angle
103 if(fBeamDiv!=0.) {BeamDivCross(0,fBeamDiv,fBeamCrossAngle,fBeamCrossPlane,pLab);}
104 if(fBeamCrossAngle!=0.) {BeamDivCross(1,fBeamDiv,fBeamCrossAngle,fBeamCrossPlane,pLab);}
106 // If required apply the Fermi momentum
108 if((fIpart==kProton) || (fIpart==kNeutron)){
109 ExtractFermi(fIpart,fPp,fProbintp,fProbintn,ddp);
111 if(fIpart==kProton) {mass = 0.93956563;}
112 if(fIpart==kNeutron) {mass = 0.93827231;}
113 // printf(" pLABx = %f pLABy = %f pLABz = %f \n",pLab[0],pLab[1],pLab[2]);
114 for(Int_t i=0; i<=2; i++){
117 balp0 = TMath::Sqrt(pLab[0]*pLab[0]+pLab[1]*pLab[1]+pLab[2]*pLab[2]+mass*mass);
118 for(Int_t i=0; i<=2; i++){
121 dddp0 = TMath::Sqrt(dddp[0]*dddp[0]+dddp[1]*dddp[1]+dddp[2]*dddp[2]+mass*mass);
123 TVector3 b(balp[0]/balp0, balp[1]/balp0, balp[2]/balp0);
124 TLorentzVector pFermi(dddp[0], dddp[1], dddp[2], dddp0);
126 // printf(" pmu -> pLABx = %f pLABy = %f pLABz = %f E = %f\n",
127 // balp[0],balp[1],balp[2],balp0);
128 // printf(" Beta -> bx = %f by = %f bz = %f\n", b[0], b[1], b[2]);
129 // printf(" pFermi -> px = %f, py = %f, pz = %f\n", pFermi[0], pFermi[1], pFermi[2]);
133 // printf(" Boosted momentum -> px = %f, py = %f, pz = %f\n",
134 // pFermi[0], pFermi[1], pFermi[2]);
135 for(Int_t i=0; i<=2; i++){
136 fBoostP[i] = pFermi[i];
141 Float_t polar[3] = {0,0,0};
142 gAlice->SetTrack(fTrackIt,-1,fIpart,fBoostP,fOrigin.GetArray(),polar,0,
146 //_____________________________________________________________________________
147 void AliGenZDC::FermiTwoGaussian(Double_t A, Float_t Z, Double_t* fPp, Double_t*
148 fProbintp, Double_t* fProbintn)
151 // Momenta distributions according to the "double-gaussian"
152 // distribution (Ilinov) - equal for protons and neutrons
154 // printf(" Initialization of Fermi momenta distribution\n");
157 Double_t sig1 = 0.113;
158 Double_t sig2 = 0.250;
159 Double_t alfa = 0.18*(TMath::Power((A/12.),(Float_t)1/3));
160 Double_t xk = (2*k2PI)/((1.+alfa)*(TMath::Power(k2PI,1.5)));
162 for(Int_t i=1; i<=200; i++){
163 Double_t p = i*0.005;
165 // printf(" fPp[%d] = %f\n",i,fPp[i]);
166 Double_t e1 = (p*p)/(2.*sig1*sig1);
167 Double_t e2 = (p*p)/(2.*sig2*sig2);
168 Double_t f1 = TMath::Exp(-(e1));
169 Double_t f2 = TMath::Exp(-(e2));
170 Double_t probp = xk*p*p*(f1/(TMath::Power(sig1,3.))+
171 alfa*f2/(TMath::Power(sig2,3.)))*0.005;
172 // printf(" probp = %f\n",probp);
173 fProbintp[i] = fProbintp[i-1] + probp;
174 fProbintn[i] = fProbintp[i];
175 // printf(" fProbintp[%d] = %f, fProbintp[%d] = %f\n",i,fProbintp[i],i,fProbintn[i]);
178 //_____________________________________________________________________________
179 void AliGenZDC::ExtractFermi(Int_t id, Double_t* fPp, Double_t* fProbintp,
180 Double_t* fProbintn, Double_t* ddp)
183 // Compute Fermi momentum for spectator nucleons
186 Float_t xx = gRandom->Rndm();
188 for(i=0; i<=200; i++){
189 if((xx>=fProbintp[i-1]) && (xx<fProbintp[i])) break;
192 else if(id==kNeutron){
193 for(i=0; i<=200; i++){
194 if((xx>=fProbintn[i-1]) && (xx<fProbintn[i])) break;
197 Float_t pext = fPp[i]+0.001;
198 Float_t phi = k2PI*(gRandom->Rndm());
199 Float_t cost = (1.-2.*(gRandom->Rndm()));
200 Float_t tet = TMath::ACos(cost);
201 ddp[0] = pext*TMath::Sin(tet)*TMath::Cos(phi);
202 ddp[1] = pext*TMath::Sin(tet)*TMath::Sin(phi);
206 //_____________________________________________________________________________
207 void AliGenZDC::BeamDivCross(Int_t icross, Float_t fBeamDiv, Float_t fBeamCrossAngle,
208 Int_t fBeamCrossPlane, Double_t* pLab)
210 Double_t tetpart, fipart, tetdiv, fidiv, angleSum[2], tetsum, fisum, dplab[3];
214 for(int i=0; i<=2; i++){
216 pmq = pmq+pLab[i]*pLab[i];
218 Double_t pmod = TMath::Sqrt(pmq);
219 // printf(" pmod = %f\n",pmod);
221 // printf(" icross = %d, fBeamDiv = %f\n",icross,fBeamDiv);
223 rvec = gRandom->Gaus(0.0,1.0);
224 tetdiv = fBeamDiv * TMath::Abs(rvec);
225 fidiv = (gRandom->Rndm())*k2PI;
228 if(fBeamCrossPlane==0.){
232 else if(fBeamCrossPlane==1.){
233 tetdiv = fBeamCrossAngle;
236 else if(fBeamCrossPlane==2.){
237 tetdiv = fBeamCrossAngle;
241 // printf(" tetdiv = %f, fidiv = %f\n",tetdiv,fidiv);
242 tetpart = TMath::ATan(TMath::Sqrt(dplab[0]*dplab[0]+dplab[1]*dplab[1])/dplab[2]);
243 if(dplab[1]!=0. || dplab[0]!=0.){
244 fipart = TMath::ATan2(dplab[1],dplab[0]);
249 if(fipart<0.) {fipart = fipart+k2PI;}
250 // printf(" tetpart = %f, fipart = %f\n",tetpart,fipart);
251 tetdiv = tetdiv*kRaddeg;
252 fidiv = fidiv*kRaddeg;
253 tetpart = tetpart*kRaddeg;
254 fipart = fipart*kRaddeg;
255 AddAngle(tetpart,fipart,tetdiv,fidiv,angleSum);
256 tetsum = angleSum[0];
258 // printf(" tetsum = %f, fisum = %f\n",tetsum,fisum);
259 tetsum = tetsum*kDegrad;
260 fisum = fisum*kDegrad;
261 pLab[0] = pmod*TMath::Sin(tetsum)*TMath::Cos(fisum);
262 pLab[1] = pmod*TMath::Sin(tetsum)*TMath::Sin(fisum);
263 pLab[2] = pmod*TMath::Cos(tetsum);
264 // printf(" pLab[0] = %f pLab[1] = %f pLab[2] = %f \n\n",
265 // pLab[0],pLab[1],pLab[2]);
266 for(Int_t i=0; i<=2; i++){
271 //_____________________________________________________________________________
272 void AliGenZDC::AddAngle(Double_t theta1, Double_t phi1, Double_t theta2,
273 Double_t phi2, Double_t* angleSum)
275 Double_t temp, conv, cx, cy, cz, ct1, st1, ct2, st2, cp1, sp1, cp2, sp2;
276 Double_t rtetsum, tetsum, fisum;
279 conv = 180./TMath::ACos(temp);
281 ct1 = TMath::Cos(theta1/conv);
282 st1 = TMath::Sin(theta1/conv);
283 cp1 = TMath::Cos(phi1/conv);
284 sp1 = TMath::Sin(phi1/conv);
285 ct2 = TMath::Cos(theta2/conv);
286 st2 = TMath::Sin(theta2/conv);
287 cp2 = TMath::Cos(phi2/conv);
288 sp2 = TMath::Sin(phi2/conv);
289 cx = ct1*cp1*st2*cp2+st1*cp1*ct2-sp1*st2*sp2;
290 cy = ct1*sp1*st2*cp2+st1*sp1*ct2+cp1*st2*sp2;
291 cz = ct1*ct2-st1*st2*cp2;
293 rtetsum = TMath::ACos(cz);
294 tetsum = conv*rtetsum;
295 if(tetsum==0. || tetsum==180.){
299 temp = cx/TMath::Sin(rtetsum);
301 if(temp<-1.) temp=-1.;
302 fisum = conv*TMath::ACos(temp);
303 if(cy<0) {fisum = 360.-fisum;}
304 // printf(" AddAngle -> tetsum = %f, fisum = %f\n",tetsum, fisum);
305 angleSum[0] = tetsum;