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