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