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