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