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