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