]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVGEN/AliGenPHOSlib.cxx
Remove fit of gaussian on inv mass histogram
[u/mrichter/AliRoot.git] / EVGEN / AliGenPHOSlib.cxx
CommitLineData
886b6f73 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
88cb7938 16/* $Id$ */
886b6f73 17
18//======================================================================
19// AliGenPHOSlib class contains parameterizations of the
20// pion, kaon, eta, omega, etaprime, phi and baryon (proton,
21// antiproton, neutron and anti-neutron) particles for the
22// study of the neutral background in PHOS detector.
23// These parameterizations are used by the
24// AliGenParam class:
25// AliGenParam(npar, param, AliGenPHOSlib::GetPt(param),
26// AliGenPHOSlib::GetY(param),
27// AliGenPHOSlib::GetIp(param) )
28// param represents the particle to be simulated :
29// Pion, Kaon, Eta, Omega, Etaprime, Phi or Baryon
30// Pt distributions are calculated from the transverse mass scaling
9ae59e17 31// with Pions, using the PtScal function taken from AliGenMUONlib
886b6f73 32// version aliroot 3.01
33//
9ae59e17 34// Gines MARTINEZ. Laurent APHECETCHE and Yves SCHUTZ
35// GPS @ SUBATECH, Nantes , France (October 1999)
886b6f73 36// http://www-subatech.in2p3.fr/~photons/subatech
37// martinez@subatech.in2p3.fr
35ead6c1 38// Additional particle species simulation options has been added:
39// Charged Pion, Charged Kaons, KLong Proton, Anti-Proton, Neutron,
40// Anti-Neutron --> Changes made by Gustavo Conesa in November 2004
8fc6910b 41// Add flat Omega(782) distribution in Nov. 2010 by Renzhuo WAN
886b6f73 42//======================================================================
43
65fb704d 44#include "TMath.h"
45#include "TRandom.h"
46
886b6f73 47#include "AliGenPHOSlib.h"
886b6f73 48
49ClassImp(AliGenPHOSlib)
50
51//======================================================================
52// P I O N S
53// (From GetPt, GetY and GetIp as param = Pion)
54// Transverse momentum distribution" PtPion
55// Rapidity distribution YPion
56// Particle distribution IdPion 111, 211 and -211 (pi0, pi+ and pi-)
57//
75e0cc59 58 Double_t AliGenPHOSlib::PtPion(const Double_t *px, const Double_t *)
886b6f73 59{
60// Pion transverse momentum distribtuion taken
9ae59e17 61// from AliGenMUONlib class, version 3.01 of aliroot
886b6f73 62// PT-PARAMETERIZATION CDF, PRL 61(88) 1819
63// POWER LAW FOR PT > 500 MEV
64// MT SCALING BELOW (T=160 MEV)
65//
f87cfe57 66 const Double_t kp0 = 1.3;
67 const Double_t kxn = 8.28;
68 const Double_t kxlim=0.5;
69 const Double_t kt=0.160;
70 const Double_t kxmpi=0.139;
71 const Double_t kb=1.;
72 Double_t y, y1, kxmpi2, ynorm, a;
886b6f73 73 Double_t x=*px;
74 //
f87cfe57 75 y1=TMath::Power(kp0/(kp0+kxlim),kxn);
76 kxmpi2=kxmpi*kxmpi;
77 ynorm=kb*(TMath::Exp(-sqrt(kxlim*kxlim+kxmpi2)/kt));
886b6f73 78 a=ynorm/y1;
f87cfe57 79 if (x > kxlim)
80 y=a*TMath::Power(kp0/(kp0+x),kxn);
886b6f73 81 else
f87cfe57 82 y=kb*TMath::Exp(-sqrt(x*x+kxmpi2)/kt);
886b6f73 83 return y*x;
84}
75e0cc59 85 Double_t AliGenPHOSlib::YPion( const Double_t *py, const Double_t *)
886b6f73 86{
f87cfe57 87//
88// pion y-distribution
89//
90
91 const Double_t ka = 7000.;
92 const Double_t kdy = 4.;
886b6f73 93
94 Double_t y=TMath::Abs(*py);
95 //
f87cfe57 96 Double_t ex = y*y/(2*kdy*kdy);
97 return ka*TMath::Exp(-ex);
886b6f73 98}
f87cfe57 99
65fb704d 100 Int_t AliGenPHOSlib::IpPion(TRandom *ran)
886b6f73 101{
f87cfe57 102// particle composition pi+, pi0, pi-
103//
104
886b6f73 105 return 111 ;
886b6f73 106}
35ead6c1 107 Int_t AliGenPHOSlib::IpChargedPion(TRandom *ran)
108{
109// particle composition pi+, pi0, pi-
110//
111
112 Float_t random = ran->Rndm();
113
114 if ( (2.*random) < 1. )
115 {
116 return 211 ;
117 }
118 else
119 {
120 return -211 ;
121 }
122}
74c62c73 123
124//End Pions
125//======================================================================
126// Pi 0 Flat Distribution
127// Transverse momentum distribution PtPi0Flat
128// Rapidity distribution YPi0Flat
129// Particle distribution IdPi0Flat 111 (pi0)
130//
131
75e0cc59 132Double_t AliGenPHOSlib::PtPi0Flat(const Double_t */*px*/, const Double_t *)
74c62c73 133{
134// Pion transverse momentum flat distribution
135
136return 1;
137
138}
139
75e0cc59 140Double_t AliGenPHOSlib::YPi0Flat( const Double_t */*py*/, const Double_t *)
74c62c73 141{
142
143// pion y-distribution
144//
145 return 1.;
146}
147
148 Int_t AliGenPHOSlib::IpPi0Flat(TRandom *)
149{
150
151// particle composition pi0
152//
153 return 111 ;
154}
155// End Pi0Flat
886b6f73 156//=============================================================
157//
f87cfe57 158 Double_t AliGenPHOSlib::PtScal(Double_t pt, Int_t np)
159{
886b6f73 160// Mt-scaling
161// Fonction for the calculation of the Pt distribution for a
162// given particle np, from the pion Pt distribution using the
163// mt scaling. This function was taken from AliGenMUONlib
164// aliroot version 3.01, and was extended for baryons
165// np = 1=>Pions 2=>Kaons 3=>Etas 4=>Omegas 5=>ETA' 6=>PHI
166// 7=>BARYONS-BARYONBARS
f87cfe57 167
886b6f73 168 // SCALING EN MASSE PAR RAPPORT A PTPI
169 // MASS 1=>PI, 2=>K, 3=>ETA, 4=>OMEGA, 5=>ETA',6=>PHI
f87cfe57 170 const Double_t khm[10] = {0.1396, 0.494, 0.547, 0.782, 0.957, 1.02,
886b6f73 171 // MASS 7=>BARYON-BARYONBAR
172 0.938, 0. , 0., 0.};
173 // VALUE MESON/PI AT 5 GEV
f87cfe57 174 const Double_t kfmax[10]={1., 1., 1., 1., 1., 1., 1., 1., 1., 1.};
886b6f73 175 np--;
f87cfe57 176 Double_t f5=TMath::Power(((sqrt(100.018215)+2.)/(sqrt(100.+khm[np]*khm[np])+2.0)),12.3);
177 Double_t kfmax2=f5/kfmax[np];
886b6f73 178 // PIONS
179 Double_t ptpion=100.*PtPion(&pt, (Double_t*) 0);
180 Double_t fmtscal=TMath::Power(((sqrt(pt*pt+0.018215)+2.)/
f87cfe57 181 (sqrt(pt*pt+khm[np]*khm[np])+2.0)),12.3)/ kfmax2;
886b6f73 182 return fmtscal*ptpion;
74c62c73 183
886b6f73 184}
185// End Scaling
186//============================================================================
187// K A O N S
75e0cc59 188 Double_t AliGenPHOSlib::PtKaon( const Double_t *px, const Double_t *)
886b6f73 189{
f87cfe57 190// kaon
191// pt-distribution
192//____________________________________________________________
193
886b6f73 194 return PtScal(*px,2); // 2==> Kaon in the PtScal function
195}
196
75e0cc59 197 Double_t AliGenPHOSlib::YKaon( const Double_t *py, const Double_t *)
886b6f73 198{
f87cfe57 199// y-distribution
200//____________________________________________________________
201
202 const Double_t ka = 1000.;
203 const Double_t kdy = 4.;
886b6f73 204
205
206 Double_t y=TMath::Abs(*py);
207 //
f87cfe57 208 Double_t ex = y*y/(2*kdy*kdy);
209 return ka*TMath::Exp(-ex);
886b6f73 210}
211
65fb704d 212 Int_t AliGenPHOSlib::IpKaon(TRandom *ran)
886b6f73 213{
f87cfe57 214// particle composition
215//
216
65fb704d 217 Float_t random = ran->Rndm();
218 Float_t random2 = ran->Rndm();
219 if (random2 < 0.5)
886b6f73 220 {
65fb704d 221 if (random < 0.5) {
886b6f73 222 return 321; // K+
223 } else {
224 return -321; // K-
225 }
226 }
227 else
228 {
65fb704d 229 if (random < 0.5) {
9ae59e17 230 return 130; // K^0 short
886b6f73 231 } else {
9ae59e17 232 return 310; // K^0 long
886b6f73 233 }
234 }
235}
35ead6c1 236
237 Int_t AliGenPHOSlib::IpChargedKaon(TRandom *ran)
238{
239// particle composition
240//
241
242 Float_t random = ran->Rndm();
243
244 if (random < 0.5) {
245 return 321; // K+
246 } else {
247 return -321; // K-
248 }
249
250
251}
252Int_t AliGenPHOSlib::IpKaon0L(TRandom *)
253{
254 // particle composition
255 //
256
257 return 130; // K^0 long
258}
886b6f73 259// End Kaons
260//============================================================================
261//============================================================================
262// E T A S
75e0cc59 263 Double_t AliGenPHOSlib::PtEta( const Double_t *px, const Double_t *)
886b6f73 264{
f87cfe57 265// etas
266// pt-distribution
267//____________________________________________________________
268
886b6f73 269 return PtScal(*px,3); // 3==> Eta in the PtScal function
270}
271
75e0cc59 272 Double_t AliGenPHOSlib::YEta( const Double_t *py, const Double_t *)
886b6f73 273{
f87cfe57 274// y-distribution
275//____________________________________________________________
276
277 const Double_t ka = 1000.;
278 const Double_t kdy = 4.;
886b6f73 279
280
281 Double_t y=TMath::Abs(*py);
282 //
f87cfe57 283 Double_t ex = y*y/(2*kdy*kdy);
284 return ka*TMath::Exp(-ex);
886b6f73 285}
286
65fb704d 287 Int_t AliGenPHOSlib::IpEta(TRandom *)
886b6f73 288{
f87cfe57 289// particle composition
290//
291
886b6f73 292 return 221; // eta
293}
294// End Etas
74c62c73 295
296//======================================================================
297// Eta Flat Distribution
298// Transverse momentum distribution PtEtaFlat
299// Rapidity distribution YEtaFlat
300// Particle distribution IdEtaFlat 111 (pi0)
301//
302
75e0cc59 303Double_t AliGenPHOSlib::PtEtaFlat(const Double_t */*px*/, const Double_t *)
74c62c73 304{
305// Eta transverse momentum flat distribution
306
307 return 1;
308
309}
310
75e0cc59 311Double_t AliGenPHOSlib::YEtaFlat( const Double_t */*py*/, const Double_t *)
74c62c73 312{
313//
314// pion y-distribution
315//
316 return 1.;
317}
318
319 Int_t AliGenPHOSlib::IpEtaFlat(TRandom *)
320{
321//
322// particle composition eta
323//
324 return 221 ;
325}
35ead6c1 326// End EtaFlat
886b6f73 327//============================================================================
328//============================================================================
329// O M E G A S
75e0cc59 330 Double_t AliGenPHOSlib::PtOmega( const Double_t *px, const Double_t *)
f87cfe57 331{
886b6f73 332// omegas
333// pt-distribution
334//____________________________________________________________
f87cfe57 335
886b6f73 336 return PtScal(*px,4); // 4==> Omega in the PtScal function
337}
338
75e0cc59 339 Double_t AliGenPHOSlib::YOmega( const Double_t *py, const Double_t *)
886b6f73 340{
f87cfe57 341// y-distribution
342//____________________________________________________________
343
344 const Double_t ka = 1000.;
345 const Double_t kdy = 4.;
886b6f73 346
347
348 Double_t y=TMath::Abs(*py);
349 //
f87cfe57 350 Double_t ex = y*y/(2*kdy*kdy);
351 return ka*TMath::Exp(-ex);
886b6f73 352}
353
65fb704d 354 Int_t AliGenPHOSlib::IpOmega(TRandom *)
886b6f73 355{
f87cfe57 356// particle composition
357//
358
886b6f73 359 return 223; // Omega
360}
361// End Omega
362//============================================================================
8fc6910b 363//======================================================================
364// Omega(782) Flat Distribution
365// Transverse momentum distribution PtOmegaFlat
366// Rapidity distribution YOmegaFlat
367// Particle distribution IdOmegaFlat 223(0mega)
368//
369
370Double_t AliGenPHOSlib::PtOmegaFlat(const Double_t */*px*/, const Double_t *)
371{
372// omega transverse momentum flat distribution
373
374return 1;
375
376}
377
378Double_t AliGenPHOSlib::YOmegaFlat( const Double_t */*py*/, const Double_t *)
379{
380
381// omega y-distribution
382//
383 return 1.;
384}
385
386 Int_t AliGenPHOSlib::IpOmegaFlat(TRandom *)
387{
388
389// particle composition omega
390//
391 return 223 ;
392}
393// End OmegaFlat
394
395
886b6f73 396//============================================================================
397// E T A P R I M E
75e0cc59 398 Double_t AliGenPHOSlib::PtEtaprime( const Double_t *px, const Double_t *)
f87cfe57 399{
886b6f73 400// etaprime
401// pt-distribution
402//____________________________________________________________
f87cfe57 403
886b6f73 404 return PtScal(*px,5); // 5==> Etaprime in the PtScal function
405}
406
75e0cc59 407 Double_t AliGenPHOSlib::YEtaprime( const Double_t *py, const Double_t *)
886b6f73 408{
f87cfe57 409// y-distribution
410//____________________________________________________________
411
412 const Double_t ka = 1000.;
413 const Double_t kdy = 4.;
886b6f73 414
415
416 Double_t y=TMath::Abs(*py);
417 //
f87cfe57 418 Double_t ex = y*y/(2*kdy*kdy);
419 return ka*TMath::Exp(-ex);
886b6f73 420}
421
65fb704d 422 Int_t AliGenPHOSlib::IpEtaprime(TRandom *)
886b6f73 423{
f87cfe57 424// particle composition
425//
426
886b6f73 427 return 331; // Etaprime
428}
429// End EtaPrime
430//===================================================================
431//============================================================================
432// P H I S
75e0cc59 433 Double_t AliGenPHOSlib::PtPhi( const Double_t *px, const Double_t *)
f87cfe57 434{
886b6f73 435// phi
436// pt-distribution
437//____________________________________________________________
f87cfe57 438
886b6f73 439 return PtScal(*px,6); // 6==> Phi in the PtScal function
440}
441
75e0cc59 442 Double_t AliGenPHOSlib::YPhi( const Double_t *py, const Double_t *)
886b6f73 443{
f87cfe57 444// y-distribution
445//____________________________________________________________
446
447 const Double_t ka = 1000.;
448 const Double_t kdy = 4.;
886b6f73 449
450
451 Double_t y=TMath::Abs(*py);
452 //
f87cfe57 453 Double_t ex = y*y/(2*kdy*kdy);
454 return ka*TMath::Exp(-ex);
886b6f73 455}
456
65fb704d 457 Int_t AliGenPHOSlib::IpPhi(TRandom *)
f87cfe57 458{
886b6f73 459// particle composition
460//
f87cfe57 461
886b6f73 462 return 333; // Phi
463}
464// End Phis
465//===================================================================
466//============================================================================
467// B A R Y O N S == protons, protonsbar, neutrons, and neutronsbars
75e0cc59 468 Double_t AliGenPHOSlib::PtBaryon( const Double_t *px, const Double_t *)
f87cfe57 469{
886b6f73 470// baryons
471// pt-distribution
472//____________________________________________________________
f87cfe57 473
886b6f73 474 return PtScal(*px,7); // 7==> Baryon in the PtScal function
475}
476
75e0cc59 477 Double_t AliGenPHOSlib::YBaryon( const Double_t *py, const Double_t *)
886b6f73 478{
f87cfe57 479// y-distribution
480//____________________________________________________________
481
482 const Double_t ka = 1000.;
483 const Double_t kdy = 4.;
886b6f73 484
485
486 Double_t y=TMath::Abs(*py);
487 //
f87cfe57 488 Double_t ex = y*y/(2*kdy*kdy);
489 return ka*TMath::Exp(-ex);
886b6f73 490}
491
65fb704d 492 Int_t AliGenPHOSlib::IpBaryon(TRandom *ran)
886b6f73 493{
f87cfe57 494// particle composition
495//
496
65fb704d 497 Float_t random = ran->Rndm();
498 Float_t random2 = ran->Rndm();
499 if (random2 < 0.5)
886b6f73 500 {
65fb704d 501 if (random < 0.5) {
886b6f73 502 return 2212; // p
503 } else {
504 return -2212; // pbar
505 }
506 }
507 else
508 {
65fb704d 509 if (random < 0.5) {
886b6f73 510 return 2112; // n
511 } else {
512 return -2112; // n bar
513 }
514 }
35ead6c1 515}
516
517 Int_t AliGenPHOSlib::IpProton(TRandom *)
518{
519// particle composition
520//
521 return 2212; // p
522
523}
524 Int_t AliGenPHOSlib::IpAProton(TRandom *)
525{
526// particle composition
527//
528 return -2212; // p bar
529
530}
531
532 Int_t AliGenPHOSlib::IpNeutron(TRandom *)
533{
534// particle composition
535//
536 return 2112; // n
537
538}
539 Int_t AliGenPHOSlib::IpANeutron(TRandom *)
540{
541// particle composition
542//
543 return -2112; // n
544
886b6f73 545}
546// End Baryons
547//===================================================================
548
75e0cc59 549typedef Double_t (*GenFunc) (const Double_t*, const Double_t*);
198bb1c7 550GenFunc AliGenPHOSlib::GetPt(Int_t param, const char* /*tname*/) const
886b6f73 551{
f87cfe57 552// Return pinter to pT parameterisation
886b6f73 553 GenFunc func;
554
555 switch (param)
35ead6c1 556 {
557 case kPion:
886b6f73 558 func=PtPion;
559 break;
35ead6c1 560 case kPi0Flat:
74c62c73 561 func=PtPi0Flat;
562 break;
35ead6c1 563 case kKaon:
886b6f73 564 func=PtKaon;
565 break;
35ead6c1 566 case kEta:
886b6f73 567 func=PtEta;
568 break;
35ead6c1 569 case kEtaFlat:
74c62c73 570 func=PtEtaFlat;
571 break;
35ead6c1 572 case kOmega:
886b6f73 573 func=PtOmega;
574 break;
8fc6910b 575 case kOmegaFlat:
576 func=PtOmegaFlat;
577 break;
35ead6c1 578 case kEtaPrime:
886b6f73 579 func=PtEtaprime;
580 break;
35ead6c1 581 case kBaryon:
886b6f73 582 func=PtBaryon;
583 break;
35ead6c1 584 default:
886b6f73 585 func=0;
586 printf("<AliGenPHOSlib::GetPt> unknown parametrisationn");
35ead6c1 587 }
886b6f73 588 return func;
589}
590
198bb1c7 591GenFunc AliGenPHOSlib::GetY(Int_t param, const char* /*tname*/) const
886b6f73 592{
35ead6c1 593 // Return pointer to Y parameterisation
594 GenFunc func;
595 switch (param)
886b6f73 596 {
34f60c01 597 case kPion:
35ead6c1 598 func=YPion;
599 break;
74c62c73 600 case kPi0Flat:
35ead6c1 601 func=YPi0Flat;
602 break;
34f60c01 603 case kKaon:
35ead6c1 604 func=YKaon;
605 break;
34f60c01 606 case kEta:
35ead6c1 607 func=YEta;
608 break;
609 case kEtaFlat:
610 func=YEtaFlat;
611 break;
34f60c01 612 case kOmega:
35ead6c1 613 func=YOmega;
614 break;
8fc6910b 615 case kOmegaFlat:
616 func=YOmegaFlat;
617 break;
34f60c01 618 case kEtaPrime:
35ead6c1 619 func=YEtaprime;
620 break;
34f60c01 621 case kPhi:
35ead6c1 622 func=YPhi;
623 break;
34f60c01 624 case kBaryon:
35ead6c1 625 func=YBaryon;
626 break;
886b6f73 627 default:
35ead6c1 628 func=0;
629 printf("<AliGenPHOSlib::GetY> unknown parametrisationn");
886b6f73 630 }
35ead6c1 631 return func;
886b6f73 632}
65fb704d 633typedef Int_t (*GenFuncIp) (TRandom *);
198bb1c7 634GenFuncIp AliGenPHOSlib::GetIp(Int_t param, const char* /*tname*/) const
886b6f73 635{
35ead6c1 636 // Return pointer to particle composition
637 GenFuncIp func;
638 switch (param)
886b6f73 639 {
74c62c73 640 case kPion:
35ead6c1 641 func=IpPion;
642 break;
643 case kChargedPion:
644 func=IpChargedPion;
645 break;
74c62c73 646 case kPi0Flat:
35ead6c1 647 func=IpPi0Flat;
648 break;
34f60c01 649 case kKaon:
35ead6c1 650 func=IpKaon;
651 break;
652 case kChargedKaon:
653 func=IpChargedKaon;
654 break;
655 case kKaon0L:
656 func=IpKaon0L;
657 break;
34f60c01 658 case kEta:
35ead6c1 659 func=IpEta;
660 break;
74c62c73 661 case kEtaFlat:
35ead6c1 662 func=IpEtaFlat;
663 break;
664
34f60c01 665 case kOmega:
35ead6c1 666 func=IpOmega;
667 break;
8fc6910b 668 case kOmegaFlat:
669 func=IpOmegaFlat;
670 break;
34f60c01 671 case kEtaPrime:
35ead6c1 672 func=IpEtaprime;
673 break;
34f60c01 674 case kPhi:
35ead6c1 675 func=IpPhi;
676 break;
34f60c01 677 case kBaryon:
35ead6c1 678 func=IpBaryon;
679 break;
680 case kProton:
681 func=IpProton;
682 break;
683 case kAProton:
684 func=IpAProton;
685 break;
686 case kNeutron:
687 func=IpNeutron;
688 break;
689 case kANeutron:
690 func=IpANeutron;
691 break;
692
886b6f73 693 default:
35ead6c1 694 func=0;
695 printf("<AliGenPHOSlib::GetIp> unknown parametrisationn");
886b6f73 696 }
35ead6c1 697 return func;
886b6f73 698}
699