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