]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenPHOSlib.cxx
flistTreeFrame attribute added; fCanvasWindow removed
[u/mrichter/AliRoot.git] / EVGEN / AliGenPHOSlib.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 $Log$
18 Revision 1.6  2000/11/30 07:12:50  alibrary
19 Introducing new Rndm and QA classes
20
21 Revision 1.5  2000/06/29 21:08:27  morsch
22 All paramatrisation libraries derive from the pure virtual base class AliGenLib.
23 This allows to pass a pointer to a library directly to AliGenParam and avoids the
24 use of function pointers in Config.C.
25
26 Revision 1.4  2000/06/14 15:21:05  morsch
27 Include clean-up (IH)
28
29 Revision 1.3  2000/06/09 20:32:54  morsch
30 All coding rule violations except RS3 corrected
31
32 Revision 1.2  1999/11/04 11:30:48  fca
33 Improve comments
34
35 Revision 1.1  1999/11/03 17:43:20  fca
36 New version from G.Martinez & A.Morsch
37
38 */
39
40 //======================================================================
41 //  AliGenPHOSlib class contains parameterizations of the
42 //  pion, kaon, eta, omega, etaprime, phi and baryon (proton, 
43 //  antiproton, neutron and anti-neutron) particles for the 
44 //  study of the neutral background in PHOS detector. 
45 //  These parameterizations are used by the 
46 //  AliGenParam  class:
47 //  AliGenParam(npar, param,  AliGenPHOSlib::GetPt(param),
48 //                            AliGenPHOSlib::GetY(param),
49 //                            AliGenPHOSlib::GetIp(param) )
50 //  param represents the particle to be simulated : 
51 //  Pion, Kaon, Eta, Omega, Etaprime, Phi or Baryon    
52 //  Pt distributions are calculated from the transverse mass scaling 
53 //  with Pions, using the PtScal function taken from AliGenMUONlib 
54 //  version aliroot 3.01
55 //
56 //     Gines MARTINEZ. Laurent APHECETCHE and Yves SCHUTZ
57 //      GPS @ SUBATECH,  Nantes , France  (October 1999)
58 //     http://www-subatech.in2p3.fr/~photons/subatech
59 //     martinez@subatech.in2p3.fr
60 //======================================================================
61
62 #include "TMath.h"
63 #include "TRandom.h"
64
65 #include "AliGenPHOSlib.h"
66
67 ClassImp(AliGenPHOSlib)
68
69 //======================================================================
70 //    P  I  O  N  S
71 //    (From GetPt, GetY and GetIp as param = Pion)
72 //    Transverse momentum distribution" PtPion 
73 //    Rapidity distribution YPion
74 //    Particle distribution IdPion  111, 211 and -211 (pi0, pi+ and pi-)
75 //
76  Double_t AliGenPHOSlib::PtPion(Double_t *px, Double_t *)
77 {
78 //     Pion transverse momentum distribtuion taken 
79 //     from AliGenMUONlib class, version 3.01 of aliroot
80 //     PT-PARAMETERIZATION CDF, PRL 61(88) 1819
81 //     POWER LAW FOR PT > 500 MEV
82 //     MT SCALING BELOW (T=160 MEV)
83 //
84   const Double_t kp0 = 1.3;
85   const Double_t kxn = 8.28;
86   const Double_t kxlim=0.5;
87   const Double_t kt=0.160;
88   const Double_t kxmpi=0.139;
89   const Double_t kb=1.;
90   Double_t y, y1, kxmpi2, ynorm, a;
91   Double_t x=*px;
92   //
93   y1=TMath::Power(kp0/(kp0+kxlim),kxn);
94   kxmpi2=kxmpi*kxmpi;
95   ynorm=kb*(TMath::Exp(-sqrt(kxlim*kxlim+kxmpi2)/kt));
96   a=ynorm/y1;
97   if (x > kxlim)
98     y=a*TMath::Power(kp0/(kp0+x),kxn);
99   else
100     y=kb*TMath::Exp(-sqrt(x*x+kxmpi2)/kt);
101   return y*x;
102 }
103  Double_t AliGenPHOSlib::YPion( Double_t *py, Double_t *)
104 {
105 //
106 // pion y-distribution
107 //
108
109   const Double_t ka    = 7000.;   
110   const Double_t kdy   = 4.;
111
112   Double_t y=TMath::Abs(*py);
113   //
114   Double_t ex = y*y/(2*kdy*kdy);
115   return ka*TMath::Exp(-ex);
116 }
117
118  Int_t AliGenPHOSlib::IpPion(TRandom *ran)
119 {
120 //                 particle composition pi+, pi0, pi-
121 //
122
123     Float_t random = ran->Rndm();
124
125     if ( (3.*random)  < 1. ) 
126     {
127           return 211 ;
128     } 
129     else
130     {  
131       if ( (3.*random) >= 2.)
132       {
133          return -211 ;
134       }
135       else 
136       {
137         return 111  ;
138       }
139     }
140 }
141 // End Pions
142 //============================================================= 
143 //
144  Double_t AliGenPHOSlib::PtScal(Double_t pt, Int_t np)
145 {
146 // Mt-scaling
147 // Fonction for the calculation of the Pt distribution for a 
148 // given particle np, from the pion Pt distribution using the 
149 // mt scaling. This function was taken from AliGenMUONlib 
150 // aliroot version 3.01, and was extended for baryons
151 // np = 1=>Pions 2=>Kaons 3=>Etas 4=>Omegas 5=>ETA' 6=>PHI
152 //      7=>BARYONS-BARYONBARS
153
154   //    SCALING EN MASSE PAR RAPPORT A PTPI
155   //    MASS                1=>PI,  2=>K, 3=>ETA, 4=>OMEGA,  5=>ETA',6=>PHI
156   const Double_t khm[10] = {0.1396, 0.494,  0.547,    0.782,   0.957,   1.02, 
157   //    MASS               7=>BARYON-BARYONBAR  
158                                          0.938, 0. , 0., 0.};
159   //     VALUE MESON/PI AT 5 GEV
160   const Double_t kfmax[10]={1., 1., 1., 1., 1., 1., 1., 1., 1., 1.};
161   np--;
162   Double_t f5=TMath::Power(((sqrt(100.018215)+2.)/(sqrt(100.+khm[np]*khm[np])+2.0)),12.3);
163   Double_t kfmax2=f5/kfmax[np];
164   // PIONS
165   Double_t ptpion=100.*PtPion(&pt, (Double_t*) 0);
166   Double_t fmtscal=TMath::Power(((sqrt(pt*pt+0.018215)+2.)/
167                                  (sqrt(pt*pt+khm[np]*khm[np])+2.0)),12.3)/ kfmax2;
168   return fmtscal*ptpion;
169 }
170 // End Scaling
171 //============================================================================
172 //    K  A  O  N  S
173  Double_t AliGenPHOSlib::PtKaon( Double_t *px, Double_t *)
174 {
175 //                kaon
176 //                pt-distribution
177 //____________________________________________________________
178
179   return PtScal(*px,2);  //  2==> Kaon in the PtScal function
180 }
181
182  Double_t AliGenPHOSlib::YKaon( Double_t *py, Double_t *)
183 {
184 // y-distribution
185 //____________________________________________________________
186
187   const Double_t ka    = 1000.;
188   const Double_t kdy   = 4.;
189
190
191   Double_t y=TMath::Abs(*py);
192   //
193   Double_t ex = y*y/(2*kdy*kdy);
194   return ka*TMath::Exp(-ex);
195 }
196
197  Int_t AliGenPHOSlib::IpKaon(TRandom *ran)
198 {
199 //                 particle composition
200 //
201
202     Float_t random = ran->Rndm();
203     Float_t random2 = ran->Rndm();
204     if (random2 < 0.5) 
205     {
206       if (random < 0.5) {       
207         return  321;   //   K+
208       } else {
209         return -321;   // K-
210       }
211     }
212     else
213     {  
214       if (random < 0.5) {       
215         return  130;   // K^0 short
216       } else {  
217         return  310;   // K^0 long
218       }
219     }
220 }
221 // End Kaons
222 //============================================================================
223 //============================================================================
224 //   E  T  A  S
225  Double_t AliGenPHOSlib::PtEta( Double_t *px, Double_t *)
226 {
227 //                etas
228 //                pt-distribution
229 //____________________________________________________________
230
231   return PtScal(*px,3);  //  3==> Eta in the PtScal function
232 }
233
234  Double_t AliGenPHOSlib::YEta( Double_t *py, Double_t *)
235 {
236 // y-distribution
237 //____________________________________________________________
238
239   const Double_t ka    = 1000.;
240   const Double_t kdy   = 4.;
241
242
243   Double_t y=TMath::Abs(*py);
244   //
245   Double_t ex = y*y/(2*kdy*kdy);
246   return ka*TMath::Exp(-ex);
247 }
248
249  Int_t AliGenPHOSlib::IpEta(TRandom *)
250 {
251 //                 particle composition
252 //
253
254         return  221;   //   eta
255 }
256 // End Etas
257 //============================================================================
258 //============================================================================
259 //    O  M  E  G  A  S
260  Double_t AliGenPHOSlib::PtOmega( Double_t *px, Double_t *)
261 {
262 // omegas
263 //                pt-distribution
264 //____________________________________________________________
265
266   return PtScal(*px,4);  //  4==> Omega in the PtScal function
267 }
268
269  Double_t AliGenPHOSlib::YOmega( Double_t *py, Double_t *)
270 {
271 // y-distribution
272 //____________________________________________________________
273
274   const Double_t ka    = 1000.;
275   const Double_t kdy   = 4.;
276
277
278   Double_t y=TMath::Abs(*py);
279   //
280   Double_t ex = y*y/(2*kdy*kdy);
281   return ka*TMath::Exp(-ex);
282 }
283
284  Int_t AliGenPHOSlib::IpOmega(TRandom *)
285 {
286 //                 particle composition
287 //
288
289         return  223;   // Omega
290 }
291 // End Omega
292 //============================================================================
293 //============================================================================
294 //    E  T  A  P  R  I  M  E
295  Double_t AliGenPHOSlib::PtEtaprime( Double_t *px, Double_t *)
296 {
297 // etaprime
298 //                pt-distribution
299 //____________________________________________________________
300
301   return PtScal(*px,5);  //  5==> Etaprime in the PtScal function
302 }
303
304  Double_t AliGenPHOSlib::YEtaprime( Double_t *py, Double_t *)
305 {
306 // y-distribution
307 //____________________________________________________________
308
309   const Double_t ka    = 1000.;
310   const Double_t kdy   = 4.;
311
312
313   Double_t y=TMath::Abs(*py);
314   //
315   Double_t ex = y*y/(2*kdy*kdy);
316   return ka*TMath::Exp(-ex);
317 }
318
319  Int_t AliGenPHOSlib::IpEtaprime(TRandom *)
320 {
321 //                 particle composition
322 //
323
324         return  331;   //   Etaprime
325 }
326 // End EtaPrime
327 //===================================================================
328 //============================================================================
329 //    P  H  I   S
330  Double_t AliGenPHOSlib::PtPhi( Double_t *px, Double_t *)
331 {
332 // phi
333 //                pt-distribution
334 //____________________________________________________________
335
336   return PtScal(*px,6);  //  6==> Phi in the PtScal function
337 }
338
339  Double_t AliGenPHOSlib::YPhi( Double_t *py, Double_t *)
340 {
341 // y-distribution
342 //____________________________________________________________
343
344   const Double_t ka    = 1000.;
345   const Double_t kdy   = 4.;
346
347
348   Double_t y=TMath::Abs(*py);
349   //
350   Double_t ex = y*y/(2*kdy*kdy);
351   return ka*TMath::Exp(-ex);
352 }
353
354  Int_t AliGenPHOSlib::IpPhi(TRandom *)
355 {
356 //                 particle composition
357 //
358     
359         return  333;   //   Phi      
360 }
361 // End Phis
362 //===================================================================
363 //============================================================================
364 //    B  A  R  Y  O  N  S  == protons, protonsbar, neutrons, and neutronsbars
365  Double_t AliGenPHOSlib::PtBaryon( Double_t *px, Double_t *)
366 {
367 // baryons
368 //                pt-distribution
369 //____________________________________________________________
370
371   return PtScal(*px,7);  //  7==> Baryon in the PtScal function
372 }
373
374  Double_t AliGenPHOSlib::YBaryon( Double_t *py, Double_t *)
375 {
376 // y-distribution
377 //____________________________________________________________
378
379   const Double_t ka    = 1000.;
380   const Double_t kdy   = 4.;
381
382
383   Double_t y=TMath::Abs(*py);
384   //
385   Double_t ex = y*y/(2*kdy*kdy);
386   return ka*TMath::Exp(-ex);
387 }
388
389  Int_t AliGenPHOSlib::IpBaryon(TRandom *ran)
390 {
391 //                 particle composition
392 //
393
394     Float_t random = ran->Rndm();
395     Float_t random2 = ran->Rndm();
396     if (random2 < 0.5) 
397     {
398       if (random < 0.5) {       
399         return  2212;   //   p
400       } else {
401         return -2212;   // pbar
402       }
403     }
404     else
405     {  
406       if (random < 0.5) {       
407         return  2112;   // n
408       } else {  
409         return -2112;   // n bar
410       }
411     }
412 }
413 // End Baryons
414 //===================================================================
415
416
417 typedef Double_t (*GenFunc) (Double_t*,  Double_t*);
418  GenFunc AliGenPHOSlib::GetPt(Int_t param, const char* tname)
419 {
420 // Return pinter to pT parameterisation
421     GenFunc func;
422     
423     switch (param)
424     {
425     case kPion:     
426         func=PtPion;
427         break;
428     case kKaon:
429         func=PtKaon;
430         break;
431     case kEta:
432         func=PtEta;
433         break;
434     case kOmega:
435         func=PtOmega;
436         break;
437     case kEtaPrime:
438         func=PtEtaprime;
439         break;
440     case kBaryon:
441         func=PtBaryon;
442         break;
443     default:
444         func=0;
445         printf("<AliGenPHOSlib::GetPt> unknown parametrisationn");
446     }
447     return func;
448 }
449
450  GenFunc AliGenPHOSlib::GetY(Int_t param, const char* tname)
451 {
452 // Return pointer to Y parameterisation
453     GenFunc func;
454     switch (param)
455     {
456     case kPion:
457         func=YPion;
458         break;
459     case kKaon:
460         func=YKaon;
461         break;
462     case kEta:
463         func=YEta;
464         break;
465     case kOmega:
466         func=YOmega;
467         break;
468     case kEtaPrime:
469         func=YEtaprime;
470         break;
471     case kPhi:
472         func=YPhi;
473         break;
474     case kBaryon:
475         func=YBaryon;
476         break;
477     default:
478         func=0;
479         printf("<AliGenPHOSlib::GetY> unknown parametrisationn");
480     }
481     return func;
482 }
483 typedef Int_t (*GenFuncIp) (TRandom *);
484  GenFuncIp AliGenPHOSlib::GetIp(Int_t param,  const char* tname)
485 {
486 // Return pointer to particle composition
487     GenFuncIp func;
488     switch (param)
489     {
490     case kPion:
491         
492         func=IpPion;
493         break;
494     case kKaon:
495         func=IpKaon;
496         break;
497     case kEta:
498         func=IpEta;
499         break;
500     case kOmega:
501         func=IpOmega;
502         break;
503     case kEtaPrime:
504         func=IpEtaprime;
505         break;
506     case kPhi:
507         func=IpPhi;
508         break;
509     case kBaryon:
510         func=IpBaryon;
511         break;
512     default:
513         func=0;
514         printf("<AliGenPHOSlib::GetIp> unknown parametrisationn");
515     }
516     return func;
517 }
518