]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ZDC/AliGenZDC.cxx
Using getter instead of global constant
[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
88cb7938 16/* $Id$ */
116cbefd 17
94de3818 18#include <assert.h>
19
116cbefd 20#include <TDatabasePDG.h>
68ca986e 21#include <TLorentzVector.h>
116cbefd 22#include <TMCProcess.h>
23#include <TPDGCode.h>
24#include <TRandom.h>
68ca986e 25#include <TVector3.h>
26
68ca986e 27#include "AliConst.h"
116cbefd 28#include "AliGenZDC.h"
68ca986e 29#include "AliRun.h"
5d12ce38 30#include "AliMC.h"
68ca986e 31
32ClassImp(AliGenZDC)
33
34//_____________________________________________________________________________
35AliGenZDC::AliGenZDC()
36 :AliGenerator()
37{
38 //
39 // Default constructor
40 //
41 fIpart = 0;
42}
43
44//_____________________________________________________________________________
45AliGenZDC::AliGenZDC(Int_t npart)
46 :AliGenerator(npart)
47{
48 //
49 // Standard constructor
50 //
51 fName = "AliGenZDC";
52 fTitle = "Generation of Test Particles for ZDCs";
53 fIpart = kNeutron;
54 fCosx = 0.;
55 fCosy = 0.;
56 fCosz = 1.;
57 fPseudoRapidity = 0.;
5a881c97 58
68ca986e 59 fFermiflag = 1;
60 // LHC values for beam divergence and crossing angle
61 fBeamDiv = 0.000032;
62 fBeamCrossAngle = 0.0001;
63 fBeamCrossPlane = 2;
5a881c97 64
65 Int_t i, j;
66 for(i=0; i<201; i++){
67 fProbintp[i] = 0;
68 fProbintn[i] = 0;
69 }
c294be39 70 for(j=0; j<3; j++) fPp[i] = 0;
5a881c97 71 fDebugOpt = 0;
68ca986e 72}
73
74//_____________________________________________________________________________
75void AliGenZDC::Init()
76{
c294be39 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,
68ca986e 82 fBeamCrossPlane);
5a881c97 83
68ca986e 84 //Initialize Fermi momentum distributions for Pb-Pb
f5cb71ad 85 FermiTwoGaussian(207.,fPp,fProbintp,fProbintn);
68ca986e 86}
87
88//_____________________________________________________________________________
89void AliGenZDC::Generate()
90{
91 //
92 // Generate one trigger (n or p)
93 //
c0ceba4c 94 Int_t i;
95
5a881c97 96 Double_t Mass, pLab[3], fP0, fP[3], fBoostP[3], ddp[3], dddp0, dddp[3];
97 Float_t fPTrack[3], ptot = fPMin;
68ca986e 98 Int_t nt;
99
866ab5a2 100 if(fPseudoRapidity==0.){
68ca986e 101 pLab[0] = ptot*fCosx;
102 pLab[1] = ptot*fCosy;
103 pLab[2] = ptot*fCosz;
104 }
105 else{
106 Float_t scang = 2*TMath::ATan(TMath::Exp(-(fPseudoRapidity)));
107 pLab[0] = -ptot*TMath::Sin(scang);
108 pLab[1] = 0.;
109 pLab[2] = ptot*TMath::Cos(scang);
110 }
c294be39 111 for(i=0; i<=2; i++) fP[i] = pLab[i];
112 if(fDebugOpt == 1){
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]);
68ca986e 115 }
116
117 // Beam divergence and crossing angle
866ab5a2 118 if(fBeamCrossAngle!=0.) {
119 BeamDivCross(1,fBeamDiv,fBeamCrossAngle,fBeamCrossPlane,pLab);
c294be39 120 for(i=0; i<=2; i++) fP[i] = pLab[i];
866ab5a2 121 }
122 if(fBeamDiv!=0.) {
123 BeamDivCross(0,fBeamDiv,fBeamCrossAngle,fBeamCrossPlane,pLab);
c294be39 124 for(i=0; i<=2; i++) fP[i] = pLab[i];
866ab5a2 125 }
126
68ca986e 127 // If required apply the Fermi momentum
128 if(fFermiflag==1){
c294be39 129 if((fIpart==kProton) || (fIpart==kNeutron))
68ca986e 130 ExtractFermi(fIpart,fPp,fProbintp,fProbintn,ddp);
5a881c97 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);
c294be39 133 for(i=0; i<=2; i++) dddp[i] = ddp[i];
5a881c97 134 dddp0 = TMath::Sqrt(dddp[0]*dddp[0]+dddp[1]*dddp[1]+dddp[2]*dddp[2]+Mass*Mass);
68ca986e 135
866ab5a2 136 TVector3 b(fP[0]/fP0, fP[1]/fP0, fP[2]/fP0);
68ca986e 137 TLorentzVector pFermi(dddp[0], dddp[1], dddp[2], dddp0);
138
68ca986e 139 pFermi.Boost(b);
c0ceba4c 140 for(i=0; i<=2; i++){
68ca986e 141 fBoostP[i] = pFermi[i];
866ab5a2 142 fP[i] = pFermi[i];
68ca986e 143 }
144
145 }
866ab5a2 146
c294be39 147 for(i=0; i<=2; i++) fPTrack[i] = fP[i];
866ab5a2 148
68ca986e 149 Float_t polar[3] = {0,0,0};
5d12ce38 150 gAlice->GetMCApp()->PushTrack(fTrackIt,-1,fIpart,fPTrack,fOrigin.GetArray(),polar,0,
1de555dc 151 kPPrimary,nt);
5a881c97 152 if(fDebugOpt == 1){
153 printf("\n\n Track momentum:\n");
154 printf("\n fPTrack = %f, %f, %f \n",fPTrack[0],fPTrack[1],fPTrack[2]);
155 }
68ca986e 156}
157
158//_____________________________________________________________________________
f5cb71ad 159void AliGenZDC::FermiTwoGaussian(Float_t A, Double_t *fPp,
5a881c97 160 Double_t *fProbintp, Double_t *fProbintn)
68ca986e 161{
162//
163// Momenta distributions according to the "double-gaussian"
164// distribution (Ilinov) - equal for protons and neutrons
165//
5a881c97 166
68ca986e 167 fProbintp[0] = 0;
168 fProbintn[0] = 0;
169 Double_t sig1 = 0.113;
170 Double_t sig2 = 0.250;
171 Double_t alfa = 0.18*(TMath::Power((A/12.),(Float_t)1/3));
172 Double_t xk = (2*k2PI)/((1.+alfa)*(TMath::Power(k2PI,1.5)));
173
174 for(Int_t i=1; i<=200; i++){
175 Double_t p = i*0.005;
176 fPp[i] = p;
68ca986e 177 Double_t e1 = (p*p)/(2.*sig1*sig1);
178 Double_t e2 = (p*p)/(2.*sig2*sig2);
179 Double_t f1 = TMath::Exp(-(e1));
180 Double_t f2 = TMath::Exp(-(e2));
181 Double_t probp = xk*p*p*(f1/(TMath::Power(sig1,3.))+
182 alfa*f2/(TMath::Power(sig2,3.)))*0.005;
68ca986e 183 fProbintp[i] = fProbintp[i-1] + probp;
184 fProbintn[i] = fProbintp[i];
5a881c97 185 }
186 if(fDebugOpt == 1){
187 printf("\n\n Initialization of Fermi momenta distribution \n");
c294be39 188 //for(Int_t i=0; i<=200; i++)
189 // printf(" fProbintp[%d] = %f, fProbintn[%d] = %f\n",i,fProbintp[i],i,fProbintn[i]);
68ca986e 190 }
191}
192//_____________________________________________________________________________
5a881c97 193void AliGenZDC::ExtractFermi(Int_t id, Double_t *fPp, Double_t *fProbintp,
194 Double_t *fProbintn, Double_t *ddp)
68ca986e 195{
196//
197// Compute Fermi momentum for spectator nucleons
198//
5a881c97 199
68ca986e 200 Int_t i;
201 Float_t xx = gRandom->Rndm();
699b37ac 202 assert ( id==kProton || id==kNeutron );
68ca986e 203 if(id==kProton){
0ff3ad02 204 for(i=1; i<=200; i++){
68ca986e 205 if((xx>=fProbintp[i-1]) && (xx<fProbintp[i])) break;
206 }
207 }
94de3818 208 else {
68ca986e 209 for(i=0; i<=200; i++){
210 if((xx>=fProbintn[i-1]) && (xx<fProbintn[i])) break;
211 }
212 }
213 Float_t pext = fPp[i]+0.001;
214 Float_t phi = k2PI*(gRandom->Rndm());
215 Float_t cost = (1.-2.*(gRandom->Rndm()));
216 Float_t tet = TMath::ACos(cost);
217 ddp[0] = pext*TMath::Sin(tet)*TMath::Cos(phi);
218 ddp[1] = pext*TMath::Sin(tet)*TMath::Sin(phi);
219 ddp[2] = pext*cost;
5a881c97 220
221 if(fDebugOpt == 1){
222 printf("\n\n Extraction of Fermi momentum\n");
223 printf("\n pxFermi = %f pyFermi = %f pzFermi = %f \n",ddp[0],ddp[1],ddp[2]);
224 }
68ca986e 225}
226
227//_____________________________________________________________________________
228void AliGenZDC::BeamDivCross(Int_t icross, Float_t fBeamDiv, Float_t fBeamCrossAngle,
5a881c97 229 Int_t fBeamCrossPlane, Double_t *pLab)
68ca986e 230{
866ab5a2 231 Double_t tetpart, fipart, tetdiv=0, fidiv=0, angleSum[2], tetsum, fisum;
68ca986e 232 Double_t rvec;
c0ceba4c 233
c294be39 234 Int_t sign=0;
235 if(pLab[2]>=0.) sign=1;
236 else sign=-1;
68ca986e 237 Double_t pmq = 0.;
c294be39 238 Int_t i;
239 for(i=0; i<=2; i++) pmq = pmq+pLab[i]*pLab[i];
68ca986e 240 Double_t pmod = TMath::Sqrt(pmq);
68ca986e 241
68ca986e 242 if(icross==0){
243 rvec = gRandom->Gaus(0.0,1.0);
244 tetdiv = fBeamDiv * TMath::Abs(rvec);
245 fidiv = (gRandom->Rndm())*k2PI;
246 }
247 else if(icross==1){
248 if(fBeamCrossPlane==0.){
249 tetdiv = 0.;
250 fidiv = 0.;
251 }
252 else if(fBeamCrossPlane==1.){
253 tetdiv = fBeamCrossAngle;
254 fidiv = 0.;
255 }
256 else if(fBeamCrossPlane==2.){
257 tetdiv = fBeamCrossAngle;
736c9b58 258 fidiv = k2PI/4.;
68ca986e 259 }
260 }
866ab5a2 261
262 tetpart = TMath::ATan(TMath::Sqrt(pLab[0]*pLab[0]+pLab[1]*pLab[1])/pLab[2]);
c294be39 263 if(pLab[1]!=0. || pLab[0]!=0.) fipart = TMath::ATan2(pLab[1],pLab[0]);
264 else fipart = 0.;
68ca986e 265 if(fipart<0.) {fipart = fipart+k2PI;}
68ca986e 266 tetdiv = tetdiv*kRaddeg;
267 fidiv = fidiv*kRaddeg;
268 tetpart = tetpart*kRaddeg;
269 fipart = fipart*kRaddeg;
270 AddAngle(tetpart,fipart,tetdiv,fidiv,angleSum);
271 tetsum = angleSum[0];
272 fisum = angleSum[1];
68ca986e 273 tetsum = tetsum*kDegrad;
274 fisum = fisum*kDegrad;
275 pLab[0] = pmod*TMath::Sin(tetsum)*TMath::Cos(fisum);
276 pLab[1] = pmod*TMath::Sin(tetsum)*TMath::Sin(fisum);
c294be39 277 if(sign==1) pLab[2] = pmod*TMath::Cos(tetsum);
278 else pLab[2] = -pmod*TMath::Cos(tetsum);
5a881c97 279 if(fDebugOpt == 1){
280 printf("\n\n Beam divergence and crossing angle\n");
c294be39 281 for(i=0; i<=2; i++)printf(" pLab[%d] = %f\n",i,pLab[i]);
68ca986e 282 }
283}
284
285//_____________________________________________________________________________
286void AliGenZDC::AddAngle(Double_t theta1, Double_t phi1, Double_t theta2,
5a881c97 287 Double_t phi2, Double_t *angleSum)
68ca986e 288{
289 Double_t temp, conv, cx, cy, cz, ct1, st1, ct2, st2, cp1, sp1, cp2, sp2;
290 Double_t rtetsum, tetsum, fisum;
291
292 temp = -1.;
293 conv = 180./TMath::ACos(temp);
294
295 ct1 = TMath::Cos(theta1/conv);
296 st1 = TMath::Sin(theta1/conv);
297 cp1 = TMath::Cos(phi1/conv);
298 sp1 = TMath::Sin(phi1/conv);
299 ct2 = TMath::Cos(theta2/conv);
300 st2 = TMath::Sin(theta2/conv);
301 cp2 = TMath::Cos(phi2/conv);
302 sp2 = TMath::Sin(phi2/conv);
303 cx = ct1*cp1*st2*cp2+st1*cp1*ct2-sp1*st2*sp2;
304 cy = ct1*sp1*st2*cp2+st1*sp1*ct2+cp1*st2*sp2;
305 cz = ct1*ct2-st1*st2*cp2;
306
307 rtetsum = TMath::ACos(cz);
308 tetsum = conv*rtetsum;
309 if(tetsum==0. || tetsum==180.){
310 fisum = 0.;
311 return;
312 }
313 temp = cx/TMath::Sin(rtetsum);
314 if(temp>1.) temp=1.;
315 if(temp<-1.) temp=-1.;
316 fisum = conv*TMath::ACos(temp);
317 if(cy<0) {fisum = 360.-fisum;}
68ca986e 318 angleSum[0] = tetsum;
319 angleSum[1] = fisum;
320}
321