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