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