]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenSTRANGElib.cxx
Option kPyCharmPbMNR added. Produce charm pairs in agreement with MNR
[u/mrichter/AliRoot.git] / EVGEN / AliGenSTRANGElib.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
17 //======================================================================
18 //  AliGenSTRANGElib class contains parameterizations of the
19 //  kaon, phi and hyperon (Lambda, Anti-Lambda, Xi, anti-Xi, Omega,
20 //  anti-Omega)  for the PPR study of the strange particle production. 
21 //  These parameterizations are used by the 
22 //  AliGenParam  class:
23 //  AliGenParam(npar, param,  AliGenSTRANGElib::GetPt(param),
24 //                            AliGenSTRANGElib::GetY(param),
25 //                            AliGenSTRANGElib::GetIp(param) )
26 //  param represents the particle to be simulated. 
27 //  ?????????
28 //  Pt distributions are calculated from the transverse mass scaling 
29 //  with Pions, using the PtScal function taken from AliGenMUONlib 
30 //  version aliroot 3.01
31 //
32 //     Rocco CALIANDRO. Rosa Anna FINI, Tiziano VIRGILI
33 //     Rocco.Caliandro@cern.ch Rosanna.Fini@ba.infn.it, 
34 //     Tiziano.Virgili@roma1.infn.it
35 //======================================================================
36
37 /*
38 $Log$
39 */
40
41 #include "TMath.h"
42 #include "TRandom.h"
43
44 #include "AliGenSTRANGElib.h"
45
46 ClassImp(AliGenSTRANGElib)
47
48 //============================================================= 
49 //
50  Double_t AliGenSTRANGElib::PtScal(Double_t pt, Int_t np)
51 {
52 // Mt-scaling
53 // Function for the calculation of the Pt distribution for a 
54 // given particle np, from the pion Pt distribution using the 
55 // mt scaling. This function was taken from AliGenMUONlib 
56 // aliroot version 3.01, and was extended for hyperons.
57 // np = 1=>Pions 2=>Kaons 3=>Etas 4=>Omegas 5=>ETA' 6=>PHI
58 //      7=>BARYONS-BARYONBARS
59 //      8=>Lambda-antiLambda
60 //      9=>Xi-antiXi
61 //     10=>Omega-antiOmega
62
63   //    MASS SCALING RESPECT TO PIONS
64   //    MASS                1=>PI,  2=>K, 3=>ETA,4=>OMEGA,5=>ETA',6=>PHI 
65   const Double_t khm[10] = {0.1396, 0.494,0.547, 0.782,   0.957,  1.02, 
66   //    MASS               7=>BARYON-BARYONBAR  
67                                  0.938, 
68   //    MASS               8=>Lambda-antiLambda  
69                                   1.1157,
70   //    MASS               9=>Xi-antiXi  
71                                   1.3213, 
72   //    MASS              10=>Omega-antiOmega  
73                                   1.6725};
74   //     VALUE MESON/PI AT 5 GEV
75   const Double_t kfmax[10]={1., 1., 1., 1., 1., 1., 1., 1., 1., 1.};
76   np--;
77   Double_t f5=TMath::Power(((sqrt(100.018215)+2.)/(sqrt(100.+khm[np]*khm[np])+2.0)),12.3);
78   Double_t kfmax2=f5/kfmax[np];
79   // PIONS
80   Double_t ptpion=100.*PtPion(&pt, (Double_t*) 0);
81   Double_t fmtscal=TMath::Power(((sqrt(pt*pt+0.018215)+2.)/
82                                  (sqrt(pt*pt+khm[np]*khm[np])+2.0)),12.3)/ kfmax2;
83   return fmtscal*ptpion;
84 }
85 //============================================================= 
86 //
87  Double_t AliGenSTRANGElib::PtPion(Double_t *px, Double_t *)
88 {
89 //     Pion transverse momentum distribtuion taken 
90 //     from AliGenMUONlib class, version 3.01 of aliroot
91 //     PT-PARAMETERIZATION CDF, PRL 61(88) 1819
92 //     POWER LAW FOR PT > 500 MEV
93 //     MT SCALING BELOW (T=160 MEV)
94 //
95   const Double_t kp0 = 1.3;
96   const Double_t kxn = 8.28;
97   const Double_t kxlim=0.5;
98   const Double_t kt=0.160;
99   const Double_t kxmpi=0.139;
100   const Double_t kb=1.;
101   Double_t y, y1, kxmpi2, ynorm, a;
102   Double_t x=*px;
103   //
104   y1=TMath::Power(kp0/(kp0+kxlim),kxn);
105   kxmpi2=kxmpi*kxmpi;
106   ynorm=kb*(TMath::Exp(-sqrt(kxlim*kxlim+kxmpi2)/kt));
107   a=ynorm/y1;
108   if (x > kxlim)
109     y=a*TMath::Power(kp0/(kp0+x),kxn);
110   else
111     y=kb*TMath::Exp(-sqrt(x*x+kxmpi2)/kt);
112   return y*x;
113 }
114 // End Scaling
115 //============================================================================
116 //    K  A  O  N  
117  Double_t AliGenSTRANGElib::PtKaon( Double_t *px, Double_t *)
118 {
119 //                kaon
120 //                pt-distribution
121 //____________________________________________________________
122
123   return PtScal(*px,2);  //  2==> Kaon in the PtScal function
124 }
125
126  Double_t AliGenSTRANGElib::YKaon( Double_t *py, Double_t *)
127 {
128 // y-distribution
129 //____________________________________________________________
130
131   const Double_t ka    = 1000.;
132   const Double_t kdy   = 4.;
133
134
135   Double_t y=TMath::Abs(*py);
136   //
137   Double_t ex = y*y/(2*kdy*kdy);
138   return ka*TMath::Exp(-ex);
139 }
140
141  Int_t AliGenSTRANGElib::IpKaon(TRandom *ran)
142 {
143 //                 particle composition
144 //
145
146     Float_t random = ran->Rndm();
147     Float_t random2 = ran->Rndm();
148     if (random2 < 0.5) 
149     {
150       if (random < 0.5) {       
151         return  321;   //   K+
152       } else {
153         return -321;   // K-
154       }
155     }
156     else
157     {  
158       if (random < 0.5) {       
159         return  130;   // K^0 short
160       } else {  
161         return  310;   // K^0 long
162       }
163     }
164 }
165 // End Kaons
166 //============================================================================
167 //============================================================================
168 //    P  H  I   
169  Double_t AliGenSTRANGElib::PtPhi( Double_t *px, Double_t *)
170 {
171 // phi
172 //                pt-distribution
173 //____________________________________________________________
174
175   return PtScal(*px,6);  //  6==> Phi in the PtScal function
176 }
177
178  Double_t AliGenSTRANGElib::YPhi( Double_t *py, Double_t *)
179 {
180 // y-distribution
181 //____________________________________________________________
182
183   const Double_t ka    = 1000.;
184   const Double_t kdy   = 4.;
185
186
187   Double_t y=TMath::Abs(*py);
188   //
189   Double_t ex = y*y/(2*kdy*kdy);
190   return ka*TMath::Exp(-ex);
191 }
192
193  Int_t AliGenSTRANGElib::IpPhi(TRandom *)
194 {
195 //                 particle composition
196 //
197     
198         return  333;   //   Phi      
199 }
200 // End Phis
201 //===================================================================
202 //============================================================================
203 //    Lambda
204  Double_t AliGenSTRANGElib::PtLambda( Double_t *px, Double_t *)
205 {
206 // Lambda
207 //                pt-distribution
208 //____________________________________________________________
209
210   return PtScal(*px,8);  //  8==> Lambda-antiLambda in the PtScal function
211 }
212
213  Double_t AliGenSTRANGElib::YLambda( Double_t *py, Double_t *)
214 {
215 // y-distribution
216 //____________________________________________________________
217
218   const Double_t ka    = 1000.;
219   const Double_t kdy   = 4.;
220
221
222   Double_t y=TMath::Abs(*py);
223   //
224   Double_t ex = y*y/(2*kdy*kdy);
225   return ka*TMath::Exp(-ex);
226 }
227
228  Int_t AliGenSTRANGElib::IpLambda(TRandom *)
229 {
230 //                 particle composition
231 //                 generation of fixed type of particle
232 //
233
234         return  3122; //   Lambda 
235 }
236 // End Lambda
237 //============================================================================
238 //    XIminus
239  Double_t AliGenSTRANGElib::PtXiMinus( Double_t *px, Double_t *)
240 {
241 // Xi
242 //                pt-distribution
243 //____________________________________________________________
244
245   return PtScal(*px,9);  //  9==> Xi-antiXi in the PtScal function
246 }
247
248  Double_t AliGenSTRANGElib::YXiMinus( Double_t *py, Double_t *)
249 {
250 // y-distribution
251 //____________________________________________________________
252
253   const Double_t ka    = 1000.;
254   const Double_t kdy   = 4.;
255
256
257   Double_t y=TMath::Abs(*py);
258   //
259   Double_t ex = y*y/(2*kdy*kdy);
260   return ka*TMath::Exp(-ex);
261 }
262
263  Int_t AliGenSTRANGElib::IpXiMinus(TRandom *)
264 {
265 //                 particle composition
266 //                 generation of fixed type of particle
267 //
268
269         return  3312; //   Xi- (only)
270 //        return  -3312; //   Xi+
271 }
272 // End Ximinus
273 //============================================================================
274 //    Omegaminus
275  Double_t AliGenSTRANGElib::PtOmegaMinus( Double_t *px, Double_t *)
276 {
277 // Omega
278 //                pt-distribution
279 //____________________________________________________________
280
281   return PtScal(*px,10);  //  10==> Omega-antiOmega in the PtScal function
282 }
283
284  Double_t AliGenSTRANGElib::YOmegaMinus( Double_t *py, Double_t *)
285 {
286 // y-distribution
287 //____________________________________________________________
288
289   const Double_t ka    = 1000.;
290   const Double_t kdy   = 4.;
291
292
293   Double_t y=TMath::Abs(*py);
294   //
295   Double_t ex = y*y/(2*kdy*kdy);
296   return ka*TMath::Exp(-ex);
297 }
298
299  Int_t AliGenSTRANGElib::IpOmegaMinus(TRandom *)
300 {
301 //                 particle composition
302 //                 generation of fixed type of particle
303 //
304
305         return  3334; //   Omega- 
306 }
307 // End Omegaminus
308 //============================================================================
309
310
311 typedef Double_t (*GenFunc) (Double_t*,  Double_t*);
312  GenFunc AliGenSTRANGElib::GetPt(Int_t param, const char* tname)
313 {
314 // Return pinter to pT parameterisation
315     GenFunc func;
316     
317     switch (param)
318     {
319     case kKaon:
320         func=PtKaon;
321         break;
322     case kPhi:
323         func=PtPhi;
324         break;
325     case kLambda:
326         func=PtLambda;
327         break;
328     case kXiMinus:
329         func=PtXiMinus;
330         break;
331     case kOmegaMinus:
332         func=PtOmegaMinus;
333         break;
334     default:
335         func=0;
336         printf("<AliGenSTRANGElib::GetPt> unknown parametrisationn");
337     }
338     return func;
339 }
340
341  GenFunc AliGenSTRANGElib::GetY(Int_t param, const char* tname)
342 {
343 // Return pointer to Y parameterisation
344     GenFunc func;
345     switch (param)
346     {
347     case kKaon:
348         func=YKaon;
349         break;
350     case kPhi:
351         func=YPhi;
352         break;
353     case kLambda:
354         func=YLambda;
355         break;
356     case kXiMinus:
357         func=YXiMinus;
358         break;
359     case kOmegaMinus:
360         func=YOmegaMinus;
361         break;
362     default:
363         func=0;
364         printf("<AliGenSTRANGElib::GetY> unknown parametrisationn");
365     }
366     return func;
367 }
368 typedef Int_t (*GenFuncIp) (TRandom *);
369  GenFuncIp AliGenSTRANGElib::GetIp(Int_t param,  const char* tname)
370 {
371 // Return pointer to particle composition
372     GenFuncIp func;
373     switch (param)
374     {
375     case kKaon:
376         func=IpKaon;
377         break;
378     case kPhi:
379         func=IpPhi;
380         break;
381     case kLambda:
382         func=IpLambda;
383         break;
384     case kXiMinus:
385         func=IpXiMinus;
386         break;
387     case kOmegaMinus:
388         func=IpOmegaMinus;
389         break;
390     default:
391         func=0;
392         printf("<AliGenSTRANGElib::GetIp> unknown parametrisationn");
393     }
394     return func;
395 }
396