]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliGenGSIlib.cxx
Corrected path to mapping libraries.
[u/mrichter/AliRoot.git] / EVGEN / AliGenGSIlib.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.7  2002/10/14 14:55:35  hristov
19 Merging the VirtualMC branch to the main development branch (HEAD)
20
21 Revision 1.5.10.1  2002/06/10 14:57:41  hristov
22 Merged with v3-08-02
23
24 Revision 1.6  2002/04/16 15:53:19  morsch
25 Coding rule violations corrected.
26
27 Revision 1.5  2001/05/15 15:43:23  morsch
28 Update of library for low mass resonances (Yiota Foka)
29
30 Revision 1.4  2001/03/09 13:01:41  morsch
31 - enum constants for paramterisation type (particle family) moved to AliGen*lib.h
32 - use AliGenGSIlib::kUpsilon, AliGenPHOSlib::kEtaPrime to access the constants
33
34 Revision 1.3  2000/12/21 16:24:06  morsch
35 Coding convention clean-up
36
37 Revision 1.2  2000/11/30 07:12:50  alibrary
38 Introducing new Rndm and QA classes
39
40 Revision 1.1  2000/06/15 08:48:43  morsch
41 AliGenGSIlib with parametersations for GSI physics simulation added (YF, MI)
42
43 Revision 1.7  2000/05/02 08:12:13  morsch
44 Coding rule violations corrected.
45
46 Revision 1.6  1999/09/29 09:24:14  fca
47 Introduction of the Copyright and cvs Log
48
49 */
50
51 /////////////////////////////////////////////////////////////////////////////
52 //                                                                         //
53 // Implementation of AliGenlib to collect parametrisations used for        //
54 // GSI simulations.                                                        //
55 // It is an extension of AliMUONLib providing in addition the option       //
56 // for different parametrisations of pt, y and ip for every particle type  //
57 //                                                                         //
58 // Responsible: Andres.Sandoval@cern.ch                                    //
59 //                                                                         //
60 /////////////////////////////////////////////////////////////////////////////
61
62 #include "TMath.h"
63 #include "TRandom.h"
64 #include "TString.h"
65 #include "AliGenGSIlib.h"
66
67
68 ClassImp(AliGenGSIlib)
69
70 //==========================================================================
71 //
72 //              Definition of Particle Distributions
73 //                    
74 //==========================================================================
75 //
76 //                         Upsilon 
77 //
78 //--------------------------------------------------------------------------
79 //
80 //                  upsilon particle composition
81 //
82 //--------------------------------------------------------------------------
83 Int_t AliGenGSIlib::IpUpsilon(TRandom *)
84 {
85 // Return upsilon pdg code
86
87   return 553;     
88
89 }
90 Double_t AliGenGSIlib::PtUpsilonFlat( Double_t *px, Double_t *dummy )
91 {
92 //--------------------------------------------------------------------------
93 //
94 //                   upsilon pt-distribution FLAT
95 //
96 //____________________________________________________________--------------
97     
98   const Double_t kptmin = 0.0;
99   const Double_t kptmax  = 15.0;
100   Double_t x=*px;
101   Double_t weight = 0.;
102
103   if (kptmin < x < kptmax) weight = 1.;
104
105   return weight;
106    
107 }
108 Double_t AliGenGSIlib::YUpsilonFlat(Double_t *py, Double_t *dummy)
109 {
110 //--------------------------------------------------------------------------
111 //
112 //                    upsilon y-distribution FLAT
113 //
114 //--------------------------------------------------------------------------
115
116   const Double_t ky0 = 1.5;
117   const Double_t kb=1.;
118   Double_t yu;
119   Double_t y=TMath::Abs(*py);
120
121   if (y < ky0)
122     yu=kb;
123   else
124     yu = 0.;
125
126   return yu;
127
128 }
129 Double_t AliGenGSIlib::PtUpsilonRitman( Double_t *px, Double_t *dummy )
130 {
131 //--------------------------------------------------------------------------
132 //
133 //                 upsilon pt-distribution RITMAN
134 //
135 //--------------------------------------------------------------------------
136
137   const Double_t kpt0 = 4.7;
138   const Double_t kxn  = 3.5;
139   Double_t x=*px;
140
141   Double_t pass1 = 1.+((x*x)/(kpt0*kpt0));
142
143   return x/TMath::Power(pass1,kxn);
144    
145 }
146 Double_t AliGenGSIlib::YUpsilonRitman(Double_t *py, Double_t *dummy)
147 {
148 //--------------------------------------------------------------------------
149 //
150 //                  upsilon y-distribution RITMAN
151 //
152 //--------------------------------------------------------------------------
153
154   const Double_t ky0 = 3.;
155   const Double_t kb=1.;
156   Double_t yu;
157   Double_t y=TMath::Abs(*py);
158
159   if (y < ky0)
160     yu=kb;
161   else
162     yu=kb*TMath::Exp(-(y-ky0)*(y-ky0)/2);
163
164   return yu;
165    
166 }
167 Double_t AliGenGSIlib::PtUpsilonKarel( Double_t *px, Double_t *dummy )
168 {
169 //--------------------------------------------------------------------------
170 //
171 //                 upsilon pt-distribution kAREL
172 //
173 //--------------------------------------------------------------------------
174 // to implement
175
176   return 0.1;   
177
178 }
179 Double_t AliGenGSIlib::YUpsilonKarel(Double_t *py, Double_t *dummy)
180 {
181 //--------------------------------------------------------------------------
182 //
183 //                  upsilon y-distribution KAREL
184 //
185 //--------------------------------------------------------------------------
186   
187   //to implement
188
189   return 0.2;  
190
191 }
192 Double_t AliGenGSIlib::PtUpsilonMUON( Double_t *px, Double_t *dummy )
193 {
194 //--------------------------------------------------------------------------
195 //
196 //                 upsilon pt-distribution MUONlib
197 //
198 //--------------------------------------------------------------------------
199
200   const Double_t kpt0 = 5.3;
201   const Double_t kxn  = 2.5;
202   Double_t x=*px;
203
204   Double_t pass1 = 1.+(x/kpt0)*(x/kpt0);
205
206   return x/TMath::Power(pass1,kxn);
207
208 }
209 Double_t AliGenGSIlib::YUpsilonMUON(Double_t *py, Double_t *dummy)
210 {
211 //--------------------------------------------------------------------------
212 //
213 //                   upsilon y-distribution MUONlib
214 //
215 //--------------------------------------------------------------------------
216
217   const Double_t ky0 = 3.;
218   const Double_t kb=1.;
219   Double_t yu;
220   Double_t y=TMath::Abs(*py);
221
222   if (y < ky0)
223     yu=kb;
224   else
225     yu=kb*TMath::Exp(-(y-ky0)*(y-ky0)/2);
226
227   return yu;
228
229 }
230 //--------------------------------------------------------------------------
231 //
232 //                             J/Psi
233 //
234 Int_t AliGenGSIlib::IpJpsi(TRandom *)
235 {
236 //--------------------------------------------------------------------------
237 //
238 //                    J/Psi particle composition
239 //
240 //--------------------------------------------------------------------------
241
242   return 443;     
243
244 }
245 Double_t AliGenGSIlib::PtJpsiFlat( Double_t *px, Double_t *dummy )
246 {
247 //--------------------------------------------------------------------------
248 //
249 //                   J/Psi pt-distribution FLAT
250 //
251 //--------------------------------------------------------------------------
252
253   const Double_t kptmin = 0.0;
254   const Double_t kptmax  = 15.0;
255   Double_t x=*px;
256   Double_t weight = 0.;
257
258   if (kptmin < x < kptmax) weight = 1.;
259
260   return weight;
261    
262 }
263 Double_t AliGenGSIlib::YJpsiFlat(Double_t *py, Double_t *dummy)
264 {
265 //--------------------------------------------------------------------------
266 //
267 //                    J/Psi y-distribution FLAT
268 //
269 //--------------------------------------------------------------------------
270
271   const Double_t ky0 = 1.5;
272   const Double_t kb=1.;
273   Double_t yu;
274   Double_t y=TMath::Abs(*py);
275
276   if (y < ky0)
277     yu=kb;
278   else
279     yu = 0.;
280
281   return yu;
282
283 }
284 Double_t AliGenGSIlib::PtJpsiMUON( Double_t *px, Double_t *dummy )
285 {
286 //--------------------------------------------------------------------------
287 //
288 //                   J/Psi pt-distribution MUONlib
289 //
290 //--------------------------------------------------------------------------
291
292   const Double_t kpt0 = 4.;
293   const Double_t kxn  = 3.6;
294   Double_t x=*px;
295
296   Double_t pass1 = 1.+(x/kpt0)*(x/kpt0);
297   return x/TMath::Power(pass1,kxn);
298    
299 }
300 Double_t AliGenGSIlib::PtJpsiRitman( Double_t *px, Double_t *dummy )
301 {
302 //--------------------------------------------------------------------------
303 //
304 //                   J/Psi pt-distribution Ritman
305 //
306 //--------------------------------------------------------------------------
307
308   const Double_t kpt0 = 2.3;
309   const Double_t kxn  = 3.5;
310   Double_t x=*px;
311
312   Double_t pass1 = 1.+(x/kpt0)*(x/kpt0);
313
314   return x/TMath::Power(pass1,kxn);
315    
316 }
317 Double_t AliGenGSIlib::YJpsiMUON(Double_t *py, Double_t *dummy)
318 {
319 //--------------------------------------------------------------------------
320 //
321 //                    J/Psi y-distribution MUONlib
322 //
323 //--------------------------------------------------------------------------
324
325   const Double_t ky0 = 4.;
326   const Double_t kb=1.;
327   Double_t yj;
328   Double_t y=TMath::Abs(*py);
329
330   if (y < ky0)
331     yj=kb;
332   else
333     yj=kb*TMath::Exp(-(y-ky0)*(y-ky0)/2);
334   return yj;
335
336 }
337 //--------------------------------------------------------------------------
338 //
339 //                  J/Psi pt-distribution by Sergei
340 //
341 //--------------------------------------------------------------------------
342 //Double_t AliGenGSIlib::PtJpsi( Double_t *px, Double_t *dummy )
343 //{
344
345 //  return x = gRandom->Rndm()*10.;
346
347 //}
348 //--------------------------------------------------------------------------
349 //
350 //                  J/Psi y-distribution by Sergei
351 //
352 //--------------------------------------------------------------------------
353 /*Double_t AliGenGSIlib::YJpsi(Double_t *py, Double_t *dummy)
354 {
355
356   const Double_t ky0 = 4.;
357   const Double_t kb=1.;
358   Double_t yj;
359   Double_t y=TMath::Abs(*py);
360   //
361   if (y < ky0)
362     yj=kb;
363   else
364     yj=kb*TMath::Exp(-(y-ky0)*(y-ky0)/2);
365   return yj;
366
367 }
368 */
369 //--------------------------------------------------------------------------
370 //
371 //                        Charm
372 //
373 //--------------------------------------------------------------------------
374 Int_t AliGenGSIlib::IpCharm(TRandom *ran)
375 {
376 //
377 //                    charm particle composition
378 //
379 //--------------------------------------------------------------------------
380   
381     Float_t random;
382     Int_t ip;
383     // 411,421,431,4122
384     random = ran->Rndm();
385     if (random < 0.5) {
386         ip=411;
387     } else if (random < 0.75) {
388         ip=421;
389     } else if (random < 0.90) {
390         ip=431;
391     } else {
392         ip=4122;
393     }
394     if (ran->Rndm() < 0.5) {ip=-ip;}
395     
396     return ip;
397
398 }
399 Double_t AliGenGSIlib::PtCharmFlat( Double_t *px, Double_t *dummy)
400 {
401 //--------------------------------------------------------------------------
402 //
403 //                    charm pt-distribution, FLAT
404 //
405 //--------------------------------------------------------------------------
406
407   Double_t x=*px;
408
409   if (x>10.)  x = 0.;
410   else x=1.;
411   return x ;
412
413 }
414 Double_t AliGenGSIlib::PtCharmGSI( Double_t *px, Double_t *dummy)
415 {
416 //--------------------------------------------------------------------------
417 //
418 //                    charm pt-distribution, from Dariuzs Miskowiec
419 //
420 //--------------------------------------------------------------------------
421
422   //Taken from PYTHIA with MRS D-' (3031 from PDFLIB), K=3.0
423   const Double_t kp1 = 1.3;
424   const Double_t kp2  = 0.39;
425   const Double_t kp3  = 0.018;
426   const Double_t kp4  = 0.91;
427   Double_t x=*px;
428
429   Double_t pass1 =TMath::Exp(-x/kp2)  ;
430   Double_t pass2 =TMath::Exp(-x/kp4)  ;
431   return TMath::Power(x,kp1) * (pass1 + kp3 * pass2);
432
433 }
434 Double_t AliGenGSIlib::PtCharmMUON( Double_t *px, Double_t *dummy)
435 {
436 //--------------------------------------------------------------------------
437 //
438 //                    charm pt-distribution, from MUONlib
439 //
440 //--------------------------------------------------------------------------
441
442   const Double_t kpt0 = 4.08;
443   const Double_t kxn  = 9.40;
444   Double_t x=*px;
445
446   Double_t pass1 = 1.+(x/kpt0)*(x/kpt0);
447
448   return x/TMath::Power(pass1,kxn);
449
450 }
451 Double_t AliGenGSIlib::YCharm( Double_t *px, Double_t *dummy)
452 {
453 //--------------------------------------------------------------------------
454 //
455 //                    charm y-distribution
456 //
457 //--------------------------------------------------------------------------
458
459     Double_t *dum=0;
460
461     return YJpsiMUON(px,dum);
462
463 }
464 //--------------------------------------------------------------------------
465 //
466 //                        Beauty
467 //
468 //--------------------------------------------------------------------------
469 Int_t AliGenGSIlib::IpBeauty(TRandom *ran)
470 {  
471 //
472 //                    beauty particle composition
473 //
474 //--------------------------------------------------------------------------
475
476     Float_t random;
477     Int_t ip;
478     random = ran->Rndm();
479     if (random < 0.5) {
480         ip=511;
481     } else if (random < 0.75) {
482         ip=521;
483     } else if (random < 0.90) {
484         ip=531;
485     } else {
486         ip=5122;
487     }
488     if (ran->Rndm() < 0.5) {ip=-ip;}
489     
490     return ip;
491 }
492 Double_t AliGenGSIlib::PtBeautyFlat( Double_t *px, Double_t *dummy)
493 {
494 //--------------------------------------------------------------------------
495 //
496 //                    beauty pt-distribution, FLAT
497 //
498 //--------------------------------------------------------------------------
499
500   Double_t x=*px;
501
502   if (x>10.) x=0.;
503   else x = 1.;
504   return x ;
505
506 }
507 Double_t AliGenGSIlib::PtBeautyGSI( Double_t *px, Double_t *dummy)
508 {
509 //--------------------------------------------------------------------------
510 //
511 //
512 //                    beauty pt-distribution, from D. Miskowiec
513 //
514 //--------------------------------------------------------------------------
515
516   //Taken from PYTHIA with MRS D-' (3031 from PDFLIB), K=3.0
517   const Double_t kp1 = 1.3;
518   const Double_t kp2  = 1.78;
519   const Double_t kp3  = 0.0096;
520   const Double_t kp4  = 4.16;
521   Double_t x=*px;
522
523   Double_t pass1 =TMath::Exp(-x/kp2)  ;
524   Double_t pass2 =TMath::Exp(-x/kp4)  ;
525
526   return TMath::Power(x,kp1) * (pass1 + kp3 * pass2);
527
528 }
529 Double_t AliGenGSIlib::PtBeautyMUON( Double_t *px, Double_t *dummy)
530 {
531 //--------------------------------------------------------------------------
532 //
533 //                    beauty pt-distribution, from MUONlib
534 //
535 //--------------------------------------------------------------------------
536
537   const Double_t kpt0 = 4.;
538   const Double_t kxn  = 3.6;
539   Double_t x=*px;
540
541   Double_t pass1 = 1.+(x/kpt0)*(x/kpt0);
542
543   return x/TMath::Power(pass1,kxn);
544
545 }
546 Double_t AliGenGSIlib::YBeauty( Double_t *px, Double_t *dummy)
547 {
548 //--------------------------------------------------------------------------
549 //
550 //                    beauty y-distribution
551 //
552 //--------------------------------------------------------------------------
553
554     Double_t *dum=0;
555
556     return YJpsiMUON(px,dum);
557
558 }
559 //--------------------------------------------------------------------------
560 //
561 //                          Eta
562 //
563 //--------------------------------------------------------------------------
564 Int_t AliGenGSIlib::IpEta(TRandom *)
565 {
566 //
567 //                 eta particle composition 
568 //
569 //--------------------------------------------------------------------------
570
571   return 221;     
572
573 }
574 Double_t AliGenGSIlib::PtEtaPHOS( Double_t *px, Double_t *dummy )
575 {
576 //--------------------------------------------------------------------------
577 //
578 //                  eta pt-distribution
579 //
580 //____________________________________________________________--------------
581
582   return PtScal(*px,3);  //  3==> Eta in the PtScal function
583    
584 }
585 Double_t AliGenGSIlib::YEtaPHOS(Double_t *py, Double_t *dummy)
586 {
587 //--------------------------------------------------------------------------
588 //
589 //                   eta y-distribution 
590 //
591 //--------------------------------------------------------------------------
592
593   const Double_t ka    = 1000.;
594   const Double_t kdy   = 4.;
595
596   Double_t y=TMath::Abs(*py);
597
598   Double_t ex = y*y/(2*kdy*kdy);
599
600   return ka*TMath::Exp(-ex);
601
602 }
603 //--------------------------------------------------------------------------
604 //
605 //                       Etaprime
606 //
607 //--------------------------------------------------------------------------
608 Int_t AliGenGSIlib::IpEtaprime(TRandom *)
609 {
610 //
611 //                 etaprime particle composition 
612 //
613 //--------------------------------------------------------------------------
614
615   return 331;     
616
617 }
618 Double_t AliGenGSIlib::PtEtaprimePHOS( Double_t *px, Double_t *dummy )
619 {
620 //--------------------------------------------------------------------------
621 //
622 //                 etaprime pt-distribution
623 //
624 //____________________________________________________________--------------
625
626   return PtScal(*px,5);  //  5==> Etaprime in the PtScal function
627    
628 }
629 Double_t AliGenGSIlib::YEtaprimePHOS(Double_t *py, Double_t *dummy)
630 {
631 //--------------------------------------------------------------------------
632 //
633 //                  etaprime y-distribution 
634 //
635 //--------------------------------------------------------------------------
636
637   const Double_t ka    = 1000.;
638   const Double_t kdy   = 4.;
639
640   Double_t y=TMath::Abs(*py);
641
642   Double_t ex = y*y/(2*kdy*kdy);
643
644   return ka*TMath::Exp(-ex);
645
646 }
647 //--------------------------------------------------------------------------
648 //
649 //                       omega 
650 //
651 //--------------------------------------------------------------------------
652 Int_t AliGenGSIlib::IpOmega(TRandom *)
653 {
654 //
655 //                 omega particle composition 
656 //
657 //--------------------------------------------------------------------------
658
659   return 223;     
660
661 }
662 Double_t AliGenGSIlib::PtOmega( Double_t *px, Double_t *dummy )
663 {
664 //--------------------------------------------------------------------------
665 //
666 //                  omega pt-distribution
667 //
668 //____________________________________________________________--------------
669
670   return PtScal(*px,4);  //  4==> Omega in the PtScal function
671    
672 }
673 Double_t AliGenGSIlib::YOmega(Double_t *py, Double_t *dummy)
674 {
675 //--------------------------------------------------------------------------
676 //
677 //                   omega y-distribution 
678 //
679 //--------------------------------------------------------------------------
680
681   const Double_t ka    = 1000.;
682   const Double_t kdy   = 4.;
683
684
685   Double_t y=TMath::Abs(*py);
686
687   Double_t ex = y*y/(2*kdy*kdy);
688
689   return ka*TMath::Exp(-ex);
690
691 }
692 //--------------------------------------------------------------------------
693 //
694 //                       Rho 
695 //
696 //--------------------------------------------------------------------------
697
698 Int_t AliGenGSIlib::IpRho(TRandom *)
699 {
700 //
701 //                 rho particle composition 
702 //
703 //--------------------------------------------------------------------------
704
705   return 113;     
706
707 }
708 Double_t AliGenGSIlib::PtRho( Double_t *px, Double_t *dummy )
709 {
710 //--------------------------------------------------------------------------
711 //
712 //                  rho pt-distribution
713 //
714 //____________________________________________________________--------------
715
716   return PtScal(*px,11);  //  11==> Rho in the PtScal function
717    
718 }
719 Double_t AliGenGSIlib::YRho(Double_t *py, Double_t *dummy)
720 {
721 //--------------------------------------------------------------------------
722 //
723 //                   rho y-distribution 
724 //
725 //--------------------------------------------------------------------------
726
727   const Double_t ka    = 1000.;
728   const Double_t kdy   = 4.;
729
730
731   Double_t y=TMath::Abs(*py);
732
733   Double_t ex = y*y/(2*kdy*kdy);
734
735   return ka*TMath::Exp(-ex);
736
737 }
738 //--------------------------------------------------------------------------
739 //
740 //                              Pion
741 //
742 //--------------------------------------------------------------------------
743 Int_t AliGenGSIlib::IpPionPHOS(TRandom *ran)
744 {
745 //
746 //                 particle composition  pi+, pi0, pi-
747 //
748 //--------------------------------------------------------------------------
749
750     Float_t random = ran->Rndm();
751
752     if ( (3.*random)  < 1. ) 
753     {
754           return 211 ;
755     } 
756     else
757     {  
758       if ( (3.*random) >= 2.)
759       {
760          return -211 ;
761       }
762       else 
763       {
764         return 111  ;
765       }
766     }
767 }
768 Double_t AliGenGSIlib::PtPion( Double_t *px, Double_t *dummy )
769 {
770 //--------------------------------------------------------------------------
771 //
772 //                  pion pt-distribution
773 //
774 //       Pion transverse momentum distribtuion as in AliGenMUONlib class, 
775 //       version 3.01 of aliroot:
776 //       PT-PARAMETERIZATION CDF, PRL 61(88) 1819
777 //       POWER LAW FOR PT > 500 MEV
778 //       MT SCALING BELOW (T=160 MEV)
779 //
780 //____________________________________________________________--------------
781
782   const Double_t kp0 = 1.3;
783   const Double_t kxn = 8.28;
784   const Double_t kxlim=0.5;
785   const Double_t kt=0.160;
786   const Double_t kxmpi=0.139;
787   const Double_t kb=1.;
788   Double_t y, y1, kxmpi2, ynorm, a;
789   Double_t x=*px;
790
791   y1=TMath::Power(kp0/(kp0+kxlim),kxn);
792   kxmpi2=kxmpi*kxmpi;
793   ynorm=kb*(TMath::Exp(-sqrt(kxlim*kxlim+kxmpi2)/kt));
794   a=ynorm/y1;
795   if (x > kxlim)
796     y=a*TMath::Power(kp0/(kp0+x),kxn);
797   else
798     y=kb*TMath::Exp(-sqrt(x*x+kxmpi2)/kt);
799   return y*x;
800    
801 }
802 Double_t AliGenGSIlib::YPion(Double_t *py, Double_t *dummy)
803 {
804 //--------------------------------------------------------------------------
805 //
806 //                    pion y-distribution 
807 //
808 //--------------------------------------------------------------------------
809
810   const Double_t ka    = 7000.;   
811   const Double_t kdy   = 4.;
812
813   Double_t y=TMath::Abs(*py);
814
815   Double_t ex = y*y/(2*kdy*kdy);
816
817   return ka*TMath::Exp(-ex);
818
819 }
820 Int_t AliGenGSIlib::IpKaonPHOS(TRandom *ran)
821 {
822 //--------------------------------------------------------------------------
823 //
824 //
825 //                        Kaon
826 //--------------------------------------------------------------------------
827 //
828 //                kaon particle composition K+, K-, Ko_short, Ko_long
829 //
830 //--------------------------------------------------------------------------
831
832     Float_t random = ran->Rndm();
833     Float_t random2 = ran->Rndm();
834     if (random2 < 0.5) 
835     {
836       if (random < 0.5) {       
837         return  321;   //   K+
838       } else {
839         return -321;   // K-
840       }
841     }
842     else
843     {  
844       if (random < 0.5) {       
845         return  130;   // K^0 short
846       } else {  
847         return  310;   // K^0 long
848       }
849     }
850 }
851 Double_t AliGenGSIlib::PtKaonPHOS( Double_t *px, Double_t *dummy )
852 {
853 //--------------------------------------------------------------------------
854 //
855 //                   kaon pt-distribution
856 //
857 //____________________________________________________________--------------
858
859   return PtScal(*px,2);  //  2==> Kaon in the PtScal function
860    
861 }
862 Double_t AliGenGSIlib::YKaonPHOS(Double_t *py, Double_t *dummy)
863 {
864 //--------------------------------------------------------------------------
865 //
866 //                    kaon y-distribution 
867 //
868 //--------------------------------------------------------------------------
869
870   const Double_t ka    = 1000.;
871   const Double_t kdy   = 4.;
872
873   Double_t y=TMath::Abs(*py);
874
875   Double_t ex = y*y/(2*kdy*kdy);
876
877   return ka*TMath::Exp(-ex);
878
879 }
880 //--------------------------------------------------------------------------
881 //
882 //                        Phi  
883 //
884 Int_t AliGenGSIlib::IpPhi(TRandom *)
885 {
886 //--------------------------------------------------------------------------
887 //
888 //                 particle composition 
889 //
890 //--------------------------------------------------------------------------
891
892   return 333;     
893
894 }
895 Double_t AliGenGSIlib::PtPhiPHOS( Double_t *px, Double_t *dummy )
896 {
897 //--------------------------------------------------------------------------
898 //
899 //                   phi pt-distribution
900 //
901 //____________________________________________________________--------------
902
903   return PtScal(*px,6);  //  6==> Phi in the PtScal function
904    
905 }
906 Double_t AliGenGSIlib::YPhiPHOS(Double_t *py, Double_t *dummy)
907 {
908 //--------------------------------------------------------------------------
909 //
910 //                    phi y-distribution 
911 //
912 //--------------------------------------------------------------------------
913
914   const Double_t ka    = 1000.;
915   const Double_t kdy   = 4.;
916
917
918   Double_t y=TMath::Abs(*py);
919  
920   Double_t ex = y*y/(2*kdy*kdy);
921
922   return ka*TMath::Exp(-ex);
923
924 }
925 Int_t AliGenGSIlib::IpBaryons(TRandom *ran)
926 {
927 //--------------------------------------------------------------------------
928 //
929 //                          Baryons
930 //
931 //--------------------------------------------------------------------------
932 //
933 //                 baryons particle composition p, pbar, n, nbar
934 //
935 //--------------------------------------------------------------------------
936
937     Float_t random = ran->Rndm();
938     Float_t random2 = ran->Rndm();
939     if (random2 < 0.5) 
940     {
941       if (random < 0.5) {       
942         return  2212;   //   p
943       } else {
944         return -2212;   // pbar
945       }
946     }
947     else
948     {  
949       if (random < 0.5) {       
950         return  2112;   // n
951       } else {  
952         return -2112;   // n bar
953       }
954     }
955 }
956 Double_t AliGenGSIlib::PtBaryons( Double_t *px, Double_t *dummy )
957 {
958 //--------------------------------------------------------------------------
959 //
960 //                  baryons pt-distribution
961 //
962 //____________________________________________________________--------------
963
964   return PtScal(*px,7);  //  7==> Baryon in the PtScal function
965    
966 }
967 Double_t AliGenGSIlib::YBaryons(Double_t *py, Double_t *dummy)
968 {
969 //--------------------------------------------------------------------------
970 //
971 //                   baryons y-distribution 
972 //
973 //--------------------------------------------------------------------------
974
975   const Double_t ka    = 1000.;
976   const Double_t kdy   = 4.;
977
978   Double_t y=TMath::Abs(*py);
979
980   Double_t ex = y*y/(2*kdy*kdy);
981
982   return ka*TMath::Exp(-ex);
983
984 }
985 //============================================================= 
986 //
987 //                    Mt-scaling as in AliGenPHOSlib  
988 //
989 //============================================================= 
990 //
991  Double_t AliGenGSIlib::PtScal(Double_t pt, Int_t np)
992 {
993 // Function for the calculation of the Pt distribution for a 
994 // given particle np, from the pion Pt distribution using the 
995 // mt scaling. 
996
997 // It was taken from AliGenPHOSlib aliroot version 3.04, which 
998 // is an update of the one in AliGenMUONlib aliroot version 3.01
999 // with an extension for Baryons but supressing Rhos
1000 // np = 1=>Pions 2=>Kaons 3=>Etas 4=>Omegas 5=>ETA' 6=>PHI
1001 //      7=>BARYONS-BARYONBARS
1002
1003 // The present adds the Rhos
1004
1005 // MASS   1=>PI, 2=>K, 3=>ETA, 4=>OMEGA, 5=>ETA', 6=>PHI 
1006 //        7=>BARYON-BARYONBAR, 11==>RHO
1007
1008   const Double_t khm[11] = {0.1396, 0.494,  0.547,    0.782,   0.957,   1.02, 
1009                                          0.938, 0. , 0., 0., 0.769};
1010
1011   //     VALUE MESON/PI AT 5 GEV
1012
1013   const Double_t kfmax[11]={1., 1., 1., 1., 1., 1., 1., 1., 1., 1., 1.};
1014
1015   np--;
1016   Double_t f5=TMath::Power(((sqrt(100.018215)+2.)/(sqrt(100.+khm[np]*khm[np])+2.0)),12.3);
1017   Double_t kfmax2=f5/kfmax[np];
1018   // PIONS
1019   Double_t ptpion=100.*PtPion(&pt, (Double_t*) 0);
1020   Double_t fmtscal=TMath::Power(((sqrt(pt*pt+0.018215)+2.)/
1021                                  (sqrt(pt*pt+khm[np]*khm[np])+2.0)),12.3)/ kfmax2;
1022   return fmtscal*ptpion;
1023 }
1024
1025 //==========================================================================
1026 //
1027 //                     Set Getters 
1028 //    
1029 //==========================================================================
1030
1031 typedef Double_t (*GenFunc) (Double_t*,  Double_t*);
1032
1033 typedef Int_t (*GenFuncIp) (TRandom *);
1034
1035 GenFunc AliGenGSIlib::GetPt(Int_t param, const char * tname) const
1036 {
1037 // Return pointer to pT parameterisation
1038    GenFunc func=0;
1039    TString sname(tname);
1040
1041    switch (param) 
1042     {
1043     case kUpsilon:
1044       if (sname=="FLAT"){
1045         func= PtUpsilonFlat;
1046         break;
1047       }
1048       if (sname=="MUON"){
1049         func= PtUpsilonMUON;
1050         break;
1051       }
1052       if (sname=="RITMAN"){
1053         func=PtUpsilonRitman;
1054         break;
1055       }
1056       if (sname=="KAREL"){
1057         func=PtUpsilonKarel;
1058         break;
1059       }
1060       func=0;
1061       printf("<AliGenGSIlib::GetPt> unknown parametrisation\n");
1062       break;
1063
1064     case kJPsi:
1065       if (sname=="FLAT"){
1066         func= PtJpsiFlat;
1067         break;
1068       }
1069       if (sname=="MUON"){
1070         func= PtJpsiMUON;
1071         break;
1072       }
1073       //      if (sname=="SERGEI"){
1074       //        func= PtJpsi;
1075       //        break;
1076       //     }
1077       func=0;
1078       printf("<AliGenGSIlib::GetPt> unknown parametrisation\n");
1079       break;
1080
1081     case kCharm: 
1082       if (sname=="FLAT"){
1083         func= PtCharmFlat;
1084         break;
1085       }
1086
1087       if (sname=="MUON"){
1088         func= PtCharmMUON;
1089         break;
1090       }
1091
1092       if (sname=="GSI"){
1093         func= PtCharmGSI;
1094         break;
1095       }
1096       func=0;
1097       printf("<AliGenGSIlib::GetPt> unknown parametrisation\n");
1098       break;
1099
1100     case kBeauty: 
1101       if (sname=="FLAT"){
1102         func= PtBeautyFlat;
1103         break;
1104       }
1105       if (sname=="MUON"){
1106         func= PtBeautyMUON;
1107         break;
1108       }
1109       if (sname=="GSI"){
1110         func= PtBeautyGSI;
1111         break;
1112       }
1113       func=0;
1114       printf("<AliGenGSIlib::GetPt> unknown parametrisation\n");
1115       break;
1116
1117
1118     case kEta:
1119       func=PtEtaPHOS;
1120       break;
1121
1122     case kEtaprime:
1123       func=PtEtaprimePHOS;
1124       break;
1125
1126     case kOmega:
1127       func=PtOmega;
1128       break;
1129
1130     case kRho:
1131       func=PtRho;
1132       break;
1133
1134     case kKaon:
1135       func=PtKaonPHOS;
1136       break;
1137
1138     case kPion:
1139       func=PtPion;
1140       break;
1141
1142     case kPhi:
1143       func=PtPhiPHOS;
1144       break;
1145
1146          //    case kLambda:
1147          //         func=PtLambda;
1148          //         break;
1149
1150     case kBaryons:
1151       func=PtBaryons;
1152       break;
1153
1154     default:
1155       func=0;
1156       printf("<AliGenGSIlib::GetPt> unknown parametrisation\n");
1157     }
1158    return func;
1159 }
1160
1161 GenFunc AliGenGSIlib::GetY(Int_t param, const char * tname) const
1162 {
1163 // Return pointer to y- parameterisation
1164    GenFunc func=0;
1165    TString sname(tname);
1166
1167    switch (param) 
1168     {
1169     case kUpsilon:
1170       if (sname=="FLAT"){
1171         func= YUpsilonFlat;
1172         break;
1173       }
1174       if (sname=="MUON"){
1175         func= YUpsilonMUON;
1176         break;
1177       }
1178       if (sname=="RITMAN"){
1179         func=YUpsilonRitman;
1180         break;
1181       }
1182       if (sname=="KAREL"){
1183         func=YUpsilonKarel;
1184         break;
1185       }
1186       func=0;
1187       printf("<AliGenGSIlib::GetY> unknown parametrisation\n");
1188       break;
1189
1190     case kJPsi:
1191       if (sname=="FLAT"){
1192         func= YJpsiFlat;
1193         break;
1194       }
1195       if (sname=="MUON"){
1196         func= YJpsiMUON;
1197         break;
1198       }
1199       //      if (sname=="SERGEI"){
1200       //        func= YJpsi;
1201       //        break;
1202       //      }
1203    
1204       func=0;
1205       printf("<AliGenGSIlib::GetY> unknown parametrisation\n");
1206       break;
1207      
1208     case kCharm: 
1209         func= YCharm;
1210         break;
1211
1212     case kBeauty: 
1213         func= YBeauty;
1214         break;
1215
1216     case kEta:
1217          func=YEtaPHOS;
1218          break;
1219
1220     case kEtaprime:
1221          func=YEtaprimePHOS;
1222          break;
1223
1224     case kOmega:
1225          func=YOmega;
1226          break;
1227
1228     case kRho:
1229          func=YRho;
1230          break;
1231
1232     case kKaon:
1233          func=YKaonPHOS;
1234          break;
1235
1236     case kPion:
1237          func=YPion;
1238          break;
1239
1240     case kPhi:
1241          func=YPhiPHOS;
1242          break;
1243
1244          //    case kLambda:
1245          //         func=YLambda;
1246          //         break;
1247
1248     case kBaryons:
1249          func=YBaryons;
1250          break;
1251
1252     default:
1253         func=0;
1254         printf("<AliGenGSIlib::GetY> unknown parametrisation\n");
1255     }
1256     return func;
1257 }
1258
1259 GenFuncIp AliGenGSIlib::GetIp(Int_t param, const char * tname) const
1260 {
1261 // Return pointer to particle type parameterisation
1262    GenFuncIp func=0;
1263    TString sname(tname);
1264
1265    switch (param) 
1266     {
1267     case kUpsilon:
1268         func= IpUpsilon;
1269         break;
1270
1271     case kJPsi:
1272         func= IpJpsi;
1273         break;
1274
1275     case kCharm: 
1276         func= IpCharm;
1277         break;
1278
1279     case kBeauty: 
1280         func= IpBeauty;
1281         break;
1282
1283     case kEta:
1284          func=IpEta;
1285          break;
1286
1287     case kEtaprime:
1288          func=IpEtaprime;
1289          break;
1290
1291     case kOmega:
1292          func=IpOmega;
1293          break;
1294
1295     case kRho:
1296          func=IpRho;
1297          break;
1298
1299     case kKaon:
1300          func=IpKaonPHOS;
1301          break;
1302
1303     case kPion:
1304          func=IpPionPHOS;
1305          break;
1306
1307     case kPhi:
1308          func=IpPhi;
1309          break;
1310
1311          //    case kLambda:
1312          //         func=IpLambda;
1313          //         break;
1314
1315     case kBaryons:
1316          func=IpBaryons;
1317          break;
1318
1319     default:
1320         func=0;
1321         printf("<AliGenGSIlib::GetIp> unknown parametrisation\n");
1322     }
1323     return func;
1324 }
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336