]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ZDC/AliGenZDC.cxx
Removed hit generation from step manager.
[u/mrichter/AliRoot.git] / ZDC / AliGenZDC.cxx
CommitLineData
68ca986e 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/*
17$Log$
94de3818 18Revision 1.2 2000/07/11 11:12:34 fca
19Some syntax corrections for non standard HP aCC
20
c0ceba4c 21Revision 1.1 2000/07/10 13:58:01 fca
22New version of ZDC from E.Scomparin & C.Oppedisano
23
68ca986e 24Revision 1.7 2000/01/19 17:17:40 fca
25
26Revision 1.6 1999/09/29 09:24:35 fca
27Introduction of the Copyright and cvs Log
28
29*/
94de3818 30#include <assert.h>
31
68ca986e 32#include <TRandom.h>
33#include <TLorentzVector.h>
34#include <TVector3.h>
94de3818 35#include "TDatabasePDG.h"
68ca986e 36
37#include "AliGenZDC.h"
38#include "AliConst.h"
39#include "AliPDG.h"
40#include "AliRun.h"
41
42ClassImp(AliGenZDC)
43
44//_____________________________________________________________________________
45AliGenZDC::AliGenZDC()
46 :AliGenerator()
47{
48 //
49 // Default constructor
50 //
51 fIpart = 0;
52}
53
54//_____________________________________________________________________________
55AliGenZDC::AliGenZDC(Int_t npart)
56 :AliGenerator(npart)
57{
58 //
59 // Standard constructor
60 //
61 fName = "AliGenZDC";
62 fTitle = "Generation of Test Particles for ZDCs";
63 fIpart = kNeutron;
64 fCosx = 0.;
65 fCosy = 0.;
66 fCosz = 1.;
67 fPseudoRapidity = 0.;
68 fFermiflag = 1;
69 // LHC values for beam divergence and crossing angle
70 fBeamDiv = 0.000032;
71 fBeamCrossAngle = 0.0001;
72 fBeamCrossPlane = 2;
73}
74
75//_____________________________________________________________________________
76void AliGenZDC::Init()
77{
78 printf(" Initializing AliGenZDC\n");
79 printf(" Fermi flag = %d, Beam Divergence = %f, Crossing Angle "
80 "= %f, Crossing Plane = %d\n\n", fFermiflag, fBeamDiv, fBeamCrossAngle,
81 fBeamCrossPlane);
82 //Initialize Fermi momentum distributions for Pb-Pb
83 FermiTwoGaussian(207.,82.,fPp,fProbintp,fProbintn);
84}
85
86//_____________________________________________________________________________
87void AliGenZDC::Generate()
88{
89 //
90 // Generate one trigger (n or p)
91 //
c0ceba4c 92 Int_t i;
93
68ca986e 94 Double_t mass, pLab[3], balp0, balp[3], ddp[3], dddp0, dddp[3];
95 Float_t ptot = fPMin;
96 Int_t nt;
97
98 if(fPseudoRapidity==0.){
99 pLab[0] = ptot*fCosx;
100 pLab[1] = ptot*fCosy;
101 pLab[2] = ptot*fCosz;
102 }
103 else{
104 Float_t scang = 2*TMath::ATan(TMath::Exp(-(fPseudoRapidity)));
105 pLab[0] = -ptot*TMath::Sin(scang);
106 pLab[1] = 0.;
107 pLab[2] = ptot*TMath::Cos(scang);
108 }
c0ceba4c 109 for(i=0; i<=2; i++){
68ca986e 110 fP[i] = pLab[i];
111 }
112
113 // Beam divergence and crossing angle
114 if(fBeamDiv!=0.) {BeamDivCross(0,fBeamDiv,fBeamCrossAngle,fBeamCrossPlane,pLab);}
115 if(fBeamCrossAngle!=0.) {BeamDivCross(1,fBeamDiv,fBeamCrossAngle,fBeamCrossPlane,pLab);}
116
117 // If required apply the Fermi momentum
118 if(fFermiflag==1){
119 if((fIpart==kProton) || (fIpart==kNeutron)){
120 ExtractFermi(fIpart,fPp,fProbintp,fProbintn,ddp);
121 }
94de3818 122 mass=gAlice->PDGDB()->GetParticle(fIpart)->Mass();
68ca986e 123// printf(" pLABx = %f pLABy = %f pLABz = %f \n",pLab[0],pLab[1],pLab[2]);
c0ceba4c 124 for(i=0; i<=2; i++){
68ca986e 125 balp[i] = -pLab[i];
126 }
127 balp0 = TMath::Sqrt(pLab[0]*pLab[0]+pLab[1]*pLab[1]+pLab[2]*pLab[2]+mass*mass);
c0ceba4c 128 for(i=0; i<=2; i++){
68ca986e 129 dddp[i] = ddp[i];
130 }
131 dddp0 = TMath::Sqrt(dddp[0]*dddp[0]+dddp[1]*dddp[1]+dddp[2]*dddp[2]+mass*mass);
132
133 TVector3 b(balp[0]/balp0, balp[1]/balp0, balp[2]/balp0);
134 TLorentzVector pFermi(dddp[0], dddp[1], dddp[2], dddp0);
135
136// printf(" pmu -> pLABx = %f pLABy = %f pLABz = %f E = %f\n",
137// balp[0],balp[1],balp[2],balp0);
138// printf(" Beta -> bx = %f by = %f bz = %f\n", b[0], b[1], b[2]);
139// printf(" pFermi -> px = %f, py = %f, pz = %f\n", pFermi[0], pFermi[1], pFermi[2]);
140
141 pFermi.Boost(b);
142
143// printf(" Boosted momentum -> px = %f, py = %f, pz = %f\n",
144// pFermi[0], pFermi[1], pFermi[2]);
c0ceba4c 145 for(i=0; i<=2; i++){
68ca986e 146 fBoostP[i] = pFermi[i];
147 }
148
149 }
150
151 Float_t polar[3] = {0,0,0};
152 gAlice->SetTrack(fTrackIt,-1,fIpart,fBoostP,fOrigin.GetArray(),polar,0,
153 "Primary",nt);
154}
155
156//_____________________________________________________________________________
157void AliGenZDC::FermiTwoGaussian(Double_t A, Float_t Z, Double_t* fPp, Double_t*
158 fProbintp, Double_t* fProbintn)
159{
160//
161// Momenta distributions according to the "double-gaussian"
162// distribution (Ilinov) - equal for protons and neutrons
163//
164// printf(" Initialization of Fermi momenta distribution\n");
165 fProbintp[0] = 0;
166 fProbintn[0] = 0;
167 Double_t sig1 = 0.113;
168 Double_t sig2 = 0.250;
169 Double_t alfa = 0.18*(TMath::Power((A/12.),(Float_t)1/3));
170 Double_t xk = (2*k2PI)/((1.+alfa)*(TMath::Power(k2PI,1.5)));
171
172 for(Int_t i=1; i<=200; i++){
173 Double_t p = i*0.005;
174 fPp[i] = p;
175// printf(" fPp[%d] = %f\n",i,fPp[i]);
176 Double_t e1 = (p*p)/(2.*sig1*sig1);
177 Double_t e2 = (p*p)/(2.*sig2*sig2);
178 Double_t f1 = TMath::Exp(-(e1));
179 Double_t f2 = TMath::Exp(-(e2));
180 Double_t probp = xk*p*p*(f1/(TMath::Power(sig1,3.))+
181 alfa*f2/(TMath::Power(sig2,3.)))*0.005;
182// printf(" probp = %f\n",probp);
183 fProbintp[i] = fProbintp[i-1] + probp;
184 fProbintn[i] = fProbintp[i];
185// printf(" fProbintp[%d] = %f, fProbintp[%d] = %f\n",i,fProbintp[i],i,fProbintn[i]);
186 }
187}
188//_____________________________________________________________________________
189void AliGenZDC::ExtractFermi(Int_t id, Double_t* fPp, Double_t* fProbintp,
190 Double_t* fProbintn, Double_t* ddp)
191{
192//
193// Compute Fermi momentum for spectator nucleons
194//
195 Int_t i;
196 Float_t xx = gRandom->Rndm();
94de3818 197 assert ( id==kProton && id==kNeutron );
68ca986e 198 if(id==kProton){
199 for(i=0; i<=200; i++){
200 if((xx>=fProbintp[i-1]) && (xx<fProbintp[i])) break;
201 }
202 }
94de3818 203 else {
68ca986e 204 for(i=0; i<=200; i++){
205 if((xx>=fProbintn[i-1]) && (xx<fProbintn[i])) break;
206 }
207 }
208 Float_t pext = fPp[i]+0.001;
209 Float_t phi = k2PI*(gRandom->Rndm());
210 Float_t cost = (1.-2.*(gRandom->Rndm()));
211 Float_t tet = TMath::ACos(cost);
212 ddp[0] = pext*TMath::Sin(tet)*TMath::Cos(phi);
213 ddp[1] = pext*TMath::Sin(tet)*TMath::Sin(phi);
214 ddp[2] = pext*cost;
215}
216
217//_____________________________________________________________________________
218void AliGenZDC::BeamDivCross(Int_t icross, Float_t fBeamDiv, Float_t fBeamCrossAngle,
219 Int_t fBeamCrossPlane, Double_t* pLab)
220{
94de3818 221 Double_t tetpart, fipart, tetdiv=0, fidiv=0, angleSum[2], tetsum, fisum, dplab[3];
68ca986e 222 Double_t rvec;
c0ceba4c 223
224 Int_t i;
68ca986e 225
226 Double_t pmq = 0.;
c0ceba4c 227 for(i=0; i<=2; i++){
68ca986e 228 dplab[i] = pLab[i];
229 pmq = pmq+pLab[i]*pLab[i];
230 }
231 Double_t pmod = TMath::Sqrt(pmq);
232// printf(" pmod = %f\n",pmod);
233
234// printf(" icross = %d, fBeamDiv = %f\n",icross,fBeamDiv);
235 if(icross==0){
236 rvec = gRandom->Gaus(0.0,1.0);
237 tetdiv = fBeamDiv * TMath::Abs(rvec);
238 fidiv = (gRandom->Rndm())*k2PI;
239 }
240 else if(icross==1){
241 if(fBeamCrossPlane==0.){
242 tetdiv = 0.;
243 fidiv = 0.;
244 }
245 else if(fBeamCrossPlane==1.){
246 tetdiv = fBeamCrossAngle;
247 fidiv = 0.;
248 }
249 else if(fBeamCrossPlane==2.){
250 tetdiv = fBeamCrossAngle;
251 fidiv = k2PI/4.;
252 }
253 }
254// printf(" tetdiv = %f, fidiv = %f\n",tetdiv,fidiv);
255 tetpart = TMath::ATan(TMath::Sqrt(dplab[0]*dplab[0]+dplab[1]*dplab[1])/dplab[2]);
256 if(dplab[1]!=0. || dplab[0]!=0.){
257 fipart = TMath::ATan2(dplab[1],dplab[0]);
258 }
259 else{
260 fipart = 0.;
261 }
262 if(fipart<0.) {fipart = fipart+k2PI;}
263// printf(" tetpart = %f, fipart = %f\n",tetpart,fipart);
264 tetdiv = tetdiv*kRaddeg;
265 fidiv = fidiv*kRaddeg;
266 tetpart = tetpart*kRaddeg;
267 fipart = fipart*kRaddeg;
268 AddAngle(tetpart,fipart,tetdiv,fidiv,angleSum);
269 tetsum = angleSum[0];
270 fisum = angleSum[1];
271// printf(" tetsum = %f, fisum = %f\n",tetsum,fisum);
272 tetsum = tetsum*kDegrad;
273 fisum = fisum*kDegrad;
274 pLab[0] = pmod*TMath::Sin(tetsum)*TMath::Cos(fisum);
275 pLab[1] = pmod*TMath::Sin(tetsum)*TMath::Sin(fisum);
276 pLab[2] = pmod*TMath::Cos(tetsum);
277// printf(" pLab[0] = %f pLab[1] = %f pLab[2] = %f \n\n",
278// pLab[0],pLab[1],pLab[2]);
c0ceba4c 279 for(i=0; i<=2; i++){
68ca986e 280 fDivP[i] = pLab[i];
281 }
282}
283
284//_____________________________________________________________________________
285void AliGenZDC::AddAngle(Double_t theta1, Double_t phi1, Double_t theta2,
286 Double_t phi2, Double_t* angleSum)
287{
288 Double_t temp, conv, cx, cy, cz, ct1, st1, ct2, st2, cp1, sp1, cp2, sp2;
289 Double_t rtetsum, tetsum, fisum;
290
291 temp = -1.;
292 conv = 180./TMath::ACos(temp);
293
294 ct1 = TMath::Cos(theta1/conv);
295 st1 = TMath::Sin(theta1/conv);
296 cp1 = TMath::Cos(phi1/conv);
297 sp1 = TMath::Sin(phi1/conv);
298 ct2 = TMath::Cos(theta2/conv);
299 st2 = TMath::Sin(theta2/conv);
300 cp2 = TMath::Cos(phi2/conv);
301 sp2 = TMath::Sin(phi2/conv);
302 cx = ct1*cp1*st2*cp2+st1*cp1*ct2-sp1*st2*sp2;
303 cy = ct1*sp1*st2*cp2+st1*sp1*ct2+cp1*st2*sp2;
304 cz = ct1*ct2-st1*st2*cp2;
305
306 rtetsum = TMath::ACos(cz);
307 tetsum = conv*rtetsum;
308 if(tetsum==0. || tetsum==180.){
309 fisum = 0.;
310 return;
311 }
312 temp = cx/TMath::Sin(rtetsum);
313 if(temp>1.) temp=1.;
314 if(temp<-1.) temp=-1.;
315 fisum = conv*TMath::ACos(temp);
316 if(cy<0) {fisum = 360.-fisum;}
317// printf(" AddAngle -> tetsum = %f, fisum = %f\n",tetsum, fisum);
318 angleSum[0] = tetsum;
319 angleSum[1] = fisum;
320}
321