]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenSTRANGElib.cxx
Updates (B. Vulpescu)
[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 /* $Id$ */
38
39 #include "TMath.h"
40 #include "TRandom.h"
41
42 #include "AliGenSTRANGElib.h"
43
44 ClassImp(AliGenSTRANGElib)
45
46 //============================================================= 
47 //
48  Double_t AliGenSTRANGElib::PtScal(Double_t pt, Int_t np)
49 {
50 // Mt-scaling
51 // Function for the calculation of the Pt distribution for a 
52 // given particle np, from the pion Pt distribution using the 
53 // mt scaling. This function was taken from AliGenMUONlib 
54 // aliroot version 3.01, and was extended for hyperons.
55 // np = 1=>Pions 2=>Kaons 3=>Etas 4=>Omegas 5=>ETA' 6=>PHI
56 //      7=>BARYONS-BARYONBARS
57 //      8=>Lambda-antiLambda
58 //      9=>Xi-antiXi
59 //     10=>Omega-antiOmega
60
61   //    MASS SCALING RESPECT TO PIONS
62   //    MASS                1=>PI,  2=>K, 3=>ETA,4=>OMEGA,5=>ETA',6=>PHI 
63   const Double_t khm[11] = {0.1396, 0.494,0.547, 0.782,   0.957,  1.02, 
64   //    MASS               7=>BARYON-BARYONBAR  
65                                  0.938, 
66   //    MASS               8=>Lambda-antiLambda  
67                                   1.1157,
68   //    MASS               9=>Xi-antiXi  
69                                   1.3213, 
70   //    MASS              10=>Omega-antiOmega  
71                                   1.6725,
72   //    MASS              11=>Lambda(1520)
73                                   1.5195};
74   //     VALUE MESON/PI AT 5 GEV
75   const Double_t kfmax[11]={1., 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 *ran)
229 {
230 //                 particle composition
231 //                 generation of fixed type of particle
232 //
233     Float_t random = ran->Rndm();
234     if (random < 0.5) {       
235       return  3122;   //   Lambda 
236     } else {  
237       return -3122;   //   Anti-Lambda
238     }
239 }
240 // End Lambda
241 //============================================================================
242 //    XIminus
243  Double_t AliGenSTRANGElib::PtXiMinus( Double_t *px, Double_t *)
244 {
245 // Xi
246 //                pt-distribution
247 //____________________________________________________________
248
249   return PtScal(*px,9);  //  9==> Xi-antiXi in the PtScal function
250 }
251
252  Double_t AliGenSTRANGElib::YXiMinus( Double_t *py, Double_t *)
253 {
254 // y-distribution
255 //____________________________________________________________
256
257   const Double_t ka    = 1000.;
258   const Double_t kdy   = 4.;
259
260
261   Double_t y=TMath::Abs(*py);
262   //
263   Double_t ex = y*y/(2*kdy*kdy);
264   return ka*TMath::Exp(-ex);
265 }
266
267  Int_t AliGenSTRANGElib::IpXiMinus(TRandom *ran)
268 {
269 //                 particle composition
270 //                 generation of fixed type of particle
271 //
272     Float_t random = ran->Rndm();
273     if (random < 0.5) {       
274       return  3312;   //   Xi- 
275     } else {  
276       return -3312;   //   Xi+
277     }
278 }
279 // End Ximinus
280 //============================================================================
281 //    Omegaminus
282  Double_t AliGenSTRANGElib::PtOmegaMinus( Double_t *px, Double_t *)
283 {
284 // Omega
285 //                pt-distribution
286 //____________________________________________________________
287
288   return PtScal(*px,10);  //  10==> Omega-antiOmega in the PtScal function
289 }
290
291  Double_t AliGenSTRANGElib::YOmegaMinus( Double_t *py, Double_t *)
292 {
293 // y-distribution
294 //____________________________________________________________
295
296   const Double_t ka    = 1000.;
297   const Double_t kdy   = 4.;
298
299
300   Double_t y=TMath::Abs(*py);
301   //
302   Double_t ex = y*y/(2*kdy*kdy);
303   return ka*TMath::Exp(-ex);
304 }
305
306  Int_t AliGenSTRANGElib::IpOmegaMinus(TRandom * ran)
307 {
308 //                 particle composition
309 //                 generation of fixed type of particle
310 //
311
312     Float_t random = ran->Rndm();
313     if (random < 0.5) {       
314       return  3334;   //   Omega- 
315     } else {  
316       return -3334;   //   Omega+
317     }
318 }
319 // End Omegaminus
320 //============================================================================
321 //     Lambda(1520)
322 Double_t AliGenSTRANGElib::PtLambda1520( Double_t *px, Double_t *)
323 {
324 // Lambda(1520)
325 //                  pt-distribution
326 //____________________________________________________________
327
328   return PtScal(*px,11);   //   11=> Lambda(1520) in the PtScal function
329 }
330
331 Double_t AliGenSTRANGElib::YLambda1520( Double_t *py, Double_t *)
332 {
333 // y-distribution
334 //____________________________________________________________
335
336   const Double_t ka   = 1000.;
337   const Double_t kdy  = 4.;
338
339   
340   Double_t y=TMath::Abs(*py);
341   //
342   Double_t ex = y*y/(2*kdy*kdy);
343   return ka*TMath::Exp(-ex);
344 }
345
346 Int_t AliGenSTRANGElib::IpLambda1520(TRandom * ran)
347 {
348 //                 particle composition
349 //                 generation of fixed type of particle
350 //
351
352    Float_t random = ran->Rndm();
353    if (random < 0.5) {       
354      return  3124;   //   Lambda(1520) 
355    } else {  
356      return -3124;   //   antiLambda(1520)
357    }
358 }
359 // End Lambda(1520)
360 //============================================================================
361
362 typedef Double_t (*GenFunc) (Double_t*,  Double_t*);
363  GenFunc AliGenSTRANGElib::GetPt(Int_t param, const char* /*tname*/) const
364 {
365 // Return pinter to pT parameterisation
366     GenFunc func;
367     
368     switch (param)
369     {
370     case kKaon:
371         func=PtKaon;
372         break;
373     case kPhi:
374         func=PtPhi;
375         break;
376     case kLambda:
377         func=PtLambda;
378         break;
379     case kXiMinus:
380         func=PtXiMinus;
381         break;
382     case kOmegaMinus:
383         func=PtOmegaMinus;
384         break;
385     case kLambda1520:
386         func=PtLambda1520;
387         break;
388     default:
389         func=0;
390         printf("<AliGenSTRANGElib::GetPt> unknown parametrisationn");
391     }
392     return func;
393 }
394
395  GenFunc AliGenSTRANGElib::GetY(Int_t param, const char* /*tname*/) const
396 {
397 // Return pointer to Y parameterisation
398     GenFunc func;
399     switch (param)
400     {
401     case kKaon:
402         func=YKaon;
403         break;
404     case kPhi:
405         func=YPhi;
406         break;
407     case kLambda:
408         func=YLambda;
409         break;
410     case kXiMinus:
411         func=YXiMinus;
412         break;
413     case kOmegaMinus:
414         func=YOmegaMinus;
415         break;
416     case kLambda1520:
417         func=YLambda1520;
418         break;
419     default:
420         func=0;
421         printf("<AliGenSTRANGElib::GetY> unknown parametrisationn");
422     }
423     return func;
424 }
425 typedef Int_t (*GenFuncIp) (TRandom *);
426  GenFuncIp AliGenSTRANGElib::GetIp(Int_t param,  const char* /*tname*/) const
427 {
428 // Return pointer to particle composition
429     GenFuncIp func;
430     switch (param)
431     {
432     case kKaon:
433         func=IpKaon;
434         break;
435     case kPhi:
436         func=IpPhi;
437         break;
438     case kLambda:
439         func=IpLambda;
440         break;
441     case kXiMinus:
442         func=IpXiMinus;
443         break;
444     case kOmegaMinus:
445         func=IpOmegaMinus;
446         break;
447     case kLambda1520:
448         func=IpLambda1520;
449         break;
450     default:
451         func=0;
452         printf("<AliGenSTRANGElib::GetIp> unknown parametrisationn");
453     }
454     return func;
455 }
456