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