1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 // Pythia 6 interface used by AliGenPythia
17 // Some settings are done by AliGenPythia, others here :)
21 #include "AliPythia.h"
22 #include "AliPythiaRndm.h"
23 #include "AliFastGlauber.h"
24 #include "AliQuenchingWeights.h"
25 #include "AliOmegaDalitz.h"
26 #include "AliDecayerExodus.h"
29 #include "TLorentzVector.h"
30 #include "PyquenCommon.h"
35 # define pyclus pyclus_
36 # define pycell pycell_
37 # define pyshow pyshow_
38 # define pyrobo pyrobo_
39 # define pyquen pyquen_
40 # define pyevnw pyevnw_
41 # define pyshowq pyshowq_
42 # define qpygin0 qpygin0_
43 # define pytune pytune_
44 # define py2ent py2ent_
45 # define setpowwght setpowwght_
48 # define pyclus PYCLUS
49 # define pycell PYCELL
50 # define pyrobo PYROBO
51 # define pyquen PYQUEN
52 # define pyevnw PYEVNW
53 # define pyshowq PYSHOWQ
54 # define qpygin0 QPYGIN0
55 # define pytune PYTUNE
56 # define py2ent PY2ENT
57 # define setpowwght SETPOWWGHT
58 # define type_of_call _stdcall
61 extern "C" void type_of_call pyclus(Int_t & );
62 extern "C" void type_of_call pycell(Int_t & );
63 extern "C" void type_of_call pyshow(Int_t &, Int_t &, Double_t &);
64 extern "C" void type_of_call pyrobo(Int_t &, Int_t &, Double_t &, Double_t &, Double_t &, Double_t &, Double_t &);
65 extern "C" void type_of_call pyquen(Double_t &, Int_t &, Double_t &);
66 extern "C" void type_of_call pyevnw();
67 extern "C" void type_of_call pyshowq(Int_t &, Int_t &, Double_t &);
68 extern "C" void type_of_call pytune(Int_t &);
69 extern "C" void type_of_call py2ent(Int_t &, Int_t&, Int_t&, Double_t&);
70 extern "C" void type_of_call qpygin0();
71 extern "C" void type_of_call setpowwght(Double_t &);
72 //_____________________________________________________________________________
74 AliPythia* AliPythia::fgAliPythia=NULL;
76 AliPythia::AliPythia():
92 // Default Constructor
95 if (!AliPythiaRndm::GetPythiaRandom())
96 AliPythiaRndm::SetPythiaRandom(GetRandom());
98 fQuenchingWeights = 0;
100 for (i = 0; i < 501; i++) fDefMDCY[i] = 0;
101 for (i = 0; i < 2001; i++) fDefMDME[i] = 0;
102 for (i = 0; i < 4; i++) fZQuench[i] = 0;
105 AliPythia::AliPythia(const AliPythia& pythia):
118 fQuenchingWeights(0),
125 for (i = 0; i < 501; i++) fDefMDCY[i] = 0;
126 for (i = 0; i < 2001; i++) fDefMDME[i] = 0;
127 for (i = 0; i < 4; i++) fZQuench[i] = 0;
131 void AliPythia::ProcInit(Process_t process, Float_t energy, StrucFunc_t strucfunc, Int_t itune)
133 // Initialise the process to generate
134 if (!AliPythiaRndm::GetPythiaRandom())
135 AliPythiaRndm::SetPythiaRandom(GetRandom());
141 fStrucFunc = strucfunc;
142 //...Switch off decay of pi0, K0S, Lambda, Sigma+-, Xi0-, Omega-.
143 SetMDCY(Pycomp(111) ,1,0); // pi0
144 SetMDCY(Pycomp(310) ,1,0); // K0S
145 SetMDCY(Pycomp(3122),1,0); // kLambda
146 SetMDCY(Pycomp(3112),1,0); // sigma -
147 SetMDCY(Pycomp(3222),1,0); // sigma +
148 SetMDCY(Pycomp(3312),1,0); // xi -
149 SetMDCY(Pycomp(3322),1,0); // xi 0
150 SetMDCY(Pycomp(3334),1,0); // omega-
151 // Select structure function
153 SetMSTP(51, AliStructFuncType::PDFsetIndex(strucfunc));
154 // Particles produced in string fragmentation point directly to either of the two endpoints
155 // of the string (depending in the side they were generated from).
159 // Pythia initialisation for selected processes//
163 for (Int_t i=1; i<= 200; i++) {
166 // select charm production
169 case kPyOldUEQ2ordered: //Old underlying events with Q2 ordered QCD processes
170 // Multiple interactions on.
172 // Double Gaussian matter distribution.
178 // Reference energy for pT0 and energy rescaling pace.
181 // String drawing almost completely minimizes string length.
184 // ISR and FSR activity.
190 case kPyOldUEQ2ordered2:
191 // Old underlying events with Q2 ordered QCD processes
192 // Multiple interactions on.
194 // Double Gaussian matter distribution.
200 // Reference energy for pT0 and energy rescaling pace.
202 SetPARP(90,0.16); // here is the difference with kPyOldUEQ2ordered
203 // String drawing almost completely minimizes string length.
206 // ISR and FSR activity.
213 // Old production mechanism: Old Popcorn
216 // (D=2) Like MSTJ(12)=2 but added prod ofthe 1er rank baryon
218 // (D=1)see can be used to form baryons (BARYON JUNCTION)
224 // heavy quark masses
254 case kPyCharmUnforced:
263 case kPyBeautyUnforced:
273 // Minimum Bias pp-Collisions
276 // select Pythia min. bias model
278 SetMSUB(92,1); // single diffraction AB-->XB
279 SetMSUB(93,1); // single diffraction AB-->AX
280 SetMSUB(94,1); // double diffraction
281 SetMSUB(95,1); // low pt production
286 case kPyMbAtlasTuneMC09:
287 // Minimum Bias pp-Collisions
290 // select Pythia min. bias model
292 SetMSUB(92,1); // single diffraction AB-->XB
293 SetMSUB(93,1); // single diffraction AB-->AX
294 SetMSUB(94,1); // double diffraction
295 SetMSUB(95,1); // low pt production
300 case kPyMbWithDirectPhoton:
301 // Minimum Bias pp-Collisions with direct photon processes added
304 // select Pythia min. bias model
306 SetMSUB(92,1); // single diffraction AB-->XB
307 SetMSUB(93,1); // single diffraction AB-->AX
308 SetMSUB(94,1); // double diffraction
309 SetMSUB(95,1); // low pt production
322 // Minimum Bias pp-Collisions
325 // select Pythia min. bias model
327 SetMSUB(92,1); // single diffraction AB-->XB
328 SetMSUB(93,1); // single diffraction AB-->AX
329 SetMSUB(94,1); // double diffraction
330 SetMSUB(95,1); // low pt production
333 // Les Houches Working Group 05 Minimum Bias pp-Collisions: hep-ph/0604120
334 // -> Pythia 6.3 or above is needed
337 SetMSUB(92,1); // single diffraction AB-->XB
338 SetMSUB(93,1); // single diffraction AB-->AX
339 SetMSUB(94,1); // double diffraction
340 SetMSUB(95,1); // low pt production
342 SetMSTP(51,AliStructFuncType::PDFsetIndex(kCTEQ6ll)); // CTEQ6ll pdf
346 SetMSTP(81,1); // Multiple Interactions ON
347 SetMSTP(82,4); // Double Gaussian Model
350 SetPARP(82,2.3); // [GeV] PT_min at Ref. energy
351 SetPARP(83,0.5); // Core density in proton matter distribution (def.value)
352 SetPARP(84,0.5); // Core radius
353 SetPARP(85,0.9); // Regulates gluon prod. mechanism
354 SetPARP(90,0.2); // 2*epsilon (exponent in power law)
358 // Minimum Bias pp-Collisions
361 // select Pythia min. bias model
363 SetMSUB(95,1); // low pt production
370 SetMSTP(91,1);// Width (1=gaussian) primordial kT dist. inside hadrons
371 SetPARP(91,1.); // <kT^2> = PARP(91,1.)^2
372 SetPARP(93,5.); // Upper cut-off
374 SetPMAS(4,1,1.2); // Charm quark mass
375 SetPMAS(5,1,4.78); // Beauty quark mass
376 SetPARP(71,4.); // Defaut value
386 // Pythia Tune A (CDF)
389 SetPARP(67,2.5); // Regulates Initial State Radiation (value from best fit to D0 dijet analysis)
390 SetMSTP(82,4); // Double Gaussian Model
391 SetPARP(82,2.0); // [GeV] PT_min at Ref. energy
392 SetPARP(84,0.4); // Core radius
393 SetPARP(85,0.90) ; // Regulates gluon prod. mechanism
394 SetPARP(86,0.95); // Regulates gluon prod. mechanism
395 SetPARP(89,1800.); // [GeV] Ref. energy
396 SetPARP(90,0.25); // 2*epsilon (exponent in power law)
402 case kPyCharmPbPbMNR:
404 case kPyDPlusPbPbMNR:
405 case kPyDPlusStrangePbPbMNR:
406 // Tuning of Pythia parameters aimed to get a resonable agreement
407 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
408 // c-cbar single inclusive and double differential distributions.
409 // This parameter settings are meant to work with Pb-Pb collisions
410 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
411 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
412 // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
424 case kPyDPlusStrangepPbMNR:
425 // Tuning of Pythia parameters aimed to get a resonable agreement
426 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
427 // c-cbar single inclusive and double differential distributions.
428 // This parameter settings are meant to work with p-Pb collisions
429 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
430 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
431 // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
444 case kPyDPlusStrangeppMNR:
445 case kPyLambdacppMNR:
446 // Tuning of Pythia parameters aimed to get a resonable agreement
447 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
448 // c-cbar single inclusive and double differential distributions.
449 // This parameter settings are meant to work with pp collisions
450 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
451 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
452 // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
462 case kPyCharmppMNRwmi:
463 // Tuning of Pythia parameters aimed to get a resonable agreement
464 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
465 // c-cbar single inclusive and double differential distributions.
466 // This parameter settings are meant to work with pp collisions
467 // and with kCTEQ5L PDFs.
468 // Added multiple interactions according to ATLAS tune settings.
469 // To get a "reasonable" agreement with MNR results, events have to be
470 // generated with the minimum ptHard (AliGenPythia::SetPtHard)
472 // To get a "perfect" agreement with MNR results, events have to be
473 // generated in four ptHard bins with the following relative
489 case kPyBeautyPbPbMNR:
490 // Tuning of Pythia parameters aimed to get a resonable agreement
491 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
492 // b-bbar single inclusive and double differential distributions.
493 // This parameter settings are meant to work with Pb-Pb collisions
494 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
495 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
496 // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
508 case kPyBeautypPbMNR:
509 // Tuning of Pythia parameters aimed to get a resonable agreement
510 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
511 // b-bbar single inclusive and double differential distributions.
512 // This parameter settings are meant to work with p-Pb collisions
513 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
514 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
515 // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
528 // Tuning of Pythia parameters aimed to get a resonable agreement
529 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
530 // b-bbar single inclusive and double differential distributions.
531 // This parameter settings are meant to work with pp collisions
532 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
533 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
534 // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
549 case kPyBeautyppMNRwmi:
550 // Tuning of Pythia parameters aimed to get a resonable agreement
551 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
552 // b-bbar single inclusive and double differential distributions.
553 // This parameter settings are meant to work with pp collisions
554 // and with kCTEQ5L PDFs.
555 // Added multiple interactions according to ATLAS tune settings.
556 // To get a "reasonable" agreement with MNR results, events have to be
557 // generated with the minimum ptHard (AliGenPythia::SetPtHard)
559 // To get a "perfect" agreement with MNR results, events have to be
560 // generated in four ptHard bins with the following relative
583 //Inclusive production of W+/-
589 // //f fbar -> gamma W+
596 // Initial/final parton shower on (Pythia default)
597 // With parton showers on we are generating "W inclusive process"
598 SetMSTP(61,1); //Initial QCD & QED showers on
599 SetMSTP(71,1); //Final QCD & QED showers on
605 //Inclusive production of Z
610 // // f fbar -> g Z/gamma
612 // // f fbar -> gamma Z/gamma
614 // // f g -> f Z/gamma
616 // // f gamma -> f Z/gamma
619 //only Z included, not gamma
622 // Initial/final parton shower on (Pythia default)
623 // With parton showers on we are generating "Z inclusive process"
624 SetMSTP(61,1); //Initial QCD & QED showers on
625 SetMSTP(71,1); //Final QCD & QED showers on
629 //Inclusive production of Z
633 // Initial/final parton shower on (Pythia default)
634 // With parton showers on we are generating "Z inclusive process"
635 SetMSTP(61,1); //Initial QCD & QED showers on
636 SetMSTP(71,1); //Final QCD & QED showers on
638 case kPyMBRSingleDiffraction:
639 case kPyMBRDoubleDiffraction:
640 case kPyMBRCentralDiffraction:
645 // For the case of jet production the following parameter setting
646 // limits the transverse momentum of secondary scatterings, due
647 // to multiple parton interactions, to be less than that of the
648 // primary interaction (see POWHEG Dijet paper arXiv:1012.3380
649 // [hep-ph] sec. 4.1 and also the PYTHIA Manual).
652 // maximum number of errors before pythia aborts (def=10)
654 // number of warnings printed on the shell
661 // number of warnings printed on the shell
672 if (GetMSTP(192) > 1 || GetMSTP(193) > 1) {
673 AliWarning(Form("Structure function for tune %5d set to %5s\n",
674 itune, AliStructFuncType::PDFsetName(strucfunc).Data()));
676 SetMSTP(51, AliStructFuncType::PDFsetIndex(strucfunc));
680 SetMSTP(41,1); // all resonance decays switched on
681 if (process == kPyJetsPWHG || process == kPyCharmPWHG || process == kPyBeautyPWHG || process == kPyWPWHG) {
682 Initialize("USER","","",0.);
684 Initialize("CMS",fProjectile,fTarget,fEcms);
690 Int_t AliPythia::CheckedLuComp(Int_t kf)
692 // Check Lund particle code (for debugging)
694 printf("\n Lucomp kf,kc %d %d",kf,kc);
698 void AliPythia::SetNuclei(Int_t a1, Int_t a2, Int_t pdf)
700 // Treat protons as inside nuclei with mass numbers a1 and a2
701 // The MSTP array in the PYPARS common block is used to enable and
702 // select the nuclear structure functions.
703 // MSTP(52) : (D=1) choice of proton and nuclear structure-function library
704 // =1: internal PYTHIA acording to MSTP(51)
705 // =2: PDFLIB proton s.f., with MSTP(51) = 1000xNGROUP+NSET
706 // If the following mass number both not equal zero, nuclear corrections of the stf are used.
707 // MSTP(192) : Mass number of nucleus side 1
708 // MSTP(193) : Mass number of nucleus side 2
709 // MSTP(194) : Nuclear structure function: 0: EKS98 8:EPS08 9:EPS09LO 19:EPS09NLO
717 AliPythia* AliPythia::Instance()
719 // Set random number generator
723 fgAliPythia = new AliPythia();
728 void AliPythia::PrintParticles()
730 // Print list of particl properties
732 char* name = new char[16];
733 for (Int_t kf=0; kf<1000000; kf++) {
734 for (Int_t c = 1; c > -2; c-=2) {
735 Int_t kc = Pycomp(c*kf);
737 Float_t mass = GetPMAS(kc,1);
738 Float_t width = GetPMAS(kc,2);
739 Float_t tau = GetPMAS(kc,4);
745 printf("\n mass, width, tau: %6d %s %10.3f %10.3e %10.3e",
746 c*kf, name, mass, width, tau);
750 printf("\n Number of particles %d \n \n", np);
753 void AliPythia::ResetDecayTable()
755 // Set default values for pythia decay switches
757 for (i = 1; i < 501; i++) SetMDCY(i,1,fDefMDCY[i]);
758 for (i = 1; i < 2001; i++) SetMDME(i,1,fDefMDME[i]);
761 void AliPythia::SetDecayTable()
763 // Set default values for pythia decay switches
766 for (i = 1; i < 501; i++) fDefMDCY[i] = GetMDCY(i,1);
767 for (i = 1; i < 2001; i++) fDefMDME[i] = GetMDME(i,1);
770 void AliPythia::Pyclus(Int_t& njet)
772 // Call Pythia clustering algorithm
777 void AliPythia::Pycell(Int_t& njet)
779 // Call Pythia jet reconstruction algorithm
784 void AliPythia::Pyshow(Int_t ip1, Int_t ip2, Double_t qmax)
786 // Call Pythia jet reconstruction algorithm
788 pyshow(ip1, ip2, qmax);
791 void AliPythia::Pyrobo(Int_t imi, Int_t ima, Double_t the, Double_t phi, Double_t bex, Double_t bey, Double_t bez)
793 pyrobo(imi, ima, the, phi, bex, bey, bez);
796 void AliPythia::Pytune(Int_t itune)
800 C ITUNE NAME (detailed descriptions below)
801 C 0 Default : No settings changed => linked Pythia version's defaults.
802 C ====== Old UE, Q2-ordered showers ==========================================
803 C 100 A : Rick Field's CDF Tune A
804 C 101 AW : Rick Field's CDF Tune AW
805 C 102 BW : Rick Field's CDF Tune BW
806 C 103 DW : Rick Field's CDF Tune DW
807 C 104 DWT : Rick Field's CDF Tune DW with slower UE energy scaling
808 C 105 QW : Rick Field's CDF Tune QW (NB: needs CTEQ6.1M pdfs externally)
809 C 106 ATLAS-DC2: Arthur Moraes' (old) ATLAS tune (ATLAS DC2 / Rome)
810 C 107 ACR : Tune A modified with annealing CR
811 C 108 D6 : Rick Field's CDF Tune D6 (NB: needs CTEQ6L pdfs externally)
812 C 109 D6T : Rick Field's CDF Tune D6T (NB: needs CTEQ6L pdfs externally)
813 C ====== Intermediate Models =================================================
814 C 200 IM 1 : Intermediate model: new UE, Q2-ordered showers, annealing CR
815 C 201 APT : Tune A modified to use pT-ordered final-state showers
816 C ====== New UE, interleaved pT-ordered showers, annealing CR ================
817 C 300 S0 : Sandhoff-Skands Tune 0
818 C 301 S1 : Sandhoff-Skands Tune 1
819 C 302 S2 : Sandhoff-Skands Tune 2
820 C 303 S0A : S0 with "Tune A" UE energy scaling
821 C 304 NOCR : New UE "best try" without colour reconnections
822 C 305 Old : New UE, original (primitive) colour reconnections
823 C 306 ATLAS-CSC: Arthur Moraes' (new) ATLAS tune (needs CTEQ6L externally)
824 C ======= The Uppsala models =================================================
825 C ( NB! must be run with special modified Pythia 6.215 version )
826 C ( available from http://www.isv.uu.se/thep/MC/scigal/ )
827 C 400 GAL 0 : Generalized area-law model. Old parameters
828 C 401 SCI 0 : Soft-Colour-Interaction model. Old parameters
829 C 402 GAL 1 : Generalized area-law model. Tevatron MB retuned (Skands)
834 void AliPythia::Py2ent(Int_t idx, Int_t pdg1, Int_t pdg2, Double_t p){
835 // Inset 2-parton system at line idx
836 py2ent(idx, pdg1, pdg2, p);
839 void AliPythia::SetWeightPower(Double_t pow)
842 SetMSTP(142, 1); // Tell Pythia to use pyevwt to calculate event wghts
844 AliWarning("Need to set minimum p_T,hard to nonzero value for weighted event generation");
847 void AliPythia::InitQuenching(Float_t cMin, Float_t cMax, Float_t k, Int_t iECMethod, Float_t zmax, Int_t ngmax)
850 // (1) The quenching model using quenching weights according to C. Salgado and U. Wiedemann
851 // (2) The nuclear geometry using the Glauber Model
854 fGlauber = AliFastGlauber::Instance();
856 fGlauber->SetCentralityClass(cMin, cMax);
858 fQuenchingWeights = new AliQuenchingWeights();
859 fQuenchingWeights->InitMult();
860 fQuenchingWeights->SetK(k);
861 fQuenchingWeights->SetECMethod(AliQuenchingWeights::kECMethod(iECMethod));
868 void AliPythia::Quench()
872 // Simple Jet Quenching routine:
873 // =============================
874 // The jet formed by all final state partons radiated by the parton created
875 // in the hard collisions is quenched by a factor (1-z) using light cone variables in
876 // the initial parton reference frame:
877 // (E + p_z)new = (1-z) (E + p_z)old
882 // The lost momentum is first balanced by one gluon with virtuality > 0.
883 // Subsequently the gluon splits to yield two gluons with E = p.
887 static Float_t eMean = 0.;
888 static Int_t icall = 0;
893 Int_t klast[4] = {-1, -1, -1, -1};
895 Int_t numpart = fPyjets->N;
896 Double_t px = 0., py = 0., pz = 0., e = 0., m = 0., p = 0., pt = 0., theta = 0., phi = 0.;
897 Double_t pxq[4], pyq[4], pzq[4], eq[4], yq[4], mq[4], pq[4], phiq[4], thetaq[4], ptq[4];
899 Double_t wjtKick[4] = {0., 0., 0., 0.};
905 // Sore information about Primary partons
908 // 0, 1 partons from hard scattering
909 // 2, 3 partons from initial state radiation
911 for (Int_t i = 2; i <= 7; i++) {
913 // Skip gluons that participate in hard scattering
914 if (i == 4 || i == 5) continue;
915 // Gluons from hard Scattering
916 if (i == 6 || i == 7) {
918 pxq[j] = fPyjets->P[0][i];
919 pyq[j] = fPyjets->P[1][i];
920 pzq[j] = fPyjets->P[2][i];
921 eq[j] = fPyjets->P[3][i];
922 mq[j] = fPyjets->P[4][i];
924 // Gluons from initial state radiation
926 // Obtain 4-momentum vector from difference between original parton and parton after gluon
927 // radiation. Energy is calculated independently because initial state radition does not
928 // conserve strictly momentum and energy for each partonic system independently.
930 // Not very clean. Should be improved !
934 pxq[j] = fPyjets->P[0][i] - fPyjets->P[0][i+2];
935 pyq[j] = fPyjets->P[1][i] - fPyjets->P[1][i+2];
936 pzq[j] = fPyjets->P[2][i] - fPyjets->P[2][i+2];
937 mq[j] = fPyjets->P[4][i];
938 eq[j] = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j] + pzq[j] * pzq[j] + mq[j] * mq[j]);
941 // Calculate some kinematic variables
943 yq[j] = 0.5 * TMath::Log((eq[j] + pzq[j] + 1.e-14) / (eq[j] - pzq[j] + 1.e-14));
944 pq[j] = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j] + pzq[j] * pzq[j]);
945 phiq[j] = TMath::Pi()+TMath::ATan2(-pyq[j], -pxq[j]);
946 ptq[j] = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j]);
947 thetaq[j] = TMath::ATan2(ptq[j], pzq[j]);
948 qPdg[j] = fPyjets->K[1][i];
954 fGlauber->GetI0I1ForPythiaAndXY(4, phiq, int0, int1, fXJet, fYJet, 15.);
956 for (Int_t j = 0; j < 4; j++) {
958 // Quench only central jets and with E > 10.
962 Int_t itype = (qPdg[j] == 21) ? 2 : 1;
963 Double_t eloss = fQuenchingWeights->GetELossRandomKFast(itype, int0[j], int1[j], eq[j]);
965 if (TMath::Abs(yq[j]) > 2.5 || eq[j] < 10.) {
968 if (eq[j] > 40. && TMath::Abs(yq[j]) < 0.5) {
974 Double_t l = fQuenchingWeights->CalcLk(int0[j], int1[j]);
975 wjtKick[j] = TMath::Sqrt(l * fQuenchingWeights->CalcQk(int0[j], int1[j]));
977 // Fractional energy loss
978 fZQuench[j] = eloss / eq[j];
980 // Avoid complete loss
982 if (fZQuench[j] > fZmax) fZQuench[j] = fZmax;
984 // Some debug printing
987 // printf("Initial parton # %3d, Type %3d Energy %10.3f Phi %10.3f Length %10.3f Loss %10.3f Kick %10.3f Mean: %10.3f %10.3f\n",
988 // j, itype, eq[j], phiq[j], l, eloss, wjtKick[j], eMean / Float_t(icall+1), yq[j]);
990 // fZQuench[j] = 0.8;
991 // while (fZQuench[j] >= 0.95) fZQuench[j] = gRandom->Exp(0.2);
994 quenched[j] = (fZQuench[j] > 0.01);
999 Double_t pNew[1000][4];
1002 Double_t zquench[4];
1006 for (Int_t isys = 0; isys < 4; isys++) {
1007 // Skip to next system if not quenched.
1008 if (!quenched[isys]) continue;
1010 nGluon[isys] = 1 + Int_t(fZQuench[isys] / (1. - fZQuench[isys]));
1011 if (nGluon[isys] > fNGmax) nGluon[isys] = fNGmax;
1012 zquench[isys] = 1. - TMath::Power(1. - fZQuench[isys], 1./Double_t(nGluon[isys]));
1013 wjtKick[isys] = wjtKick[isys] / TMath::Sqrt(Double_t(nGluon[isys]));
1019 Double_t pg[4] = {0., 0., 0., 0.};
1022 // Loop on radiation events
1024 for (Int_t iglu = 0; iglu < nGluon[isys]; iglu++) {
1027 for (Int_t k = 0; k < 4; k++)
1033 // Loop over partons
1034 for (Int_t i = 0; i < numpart; i++)
1036 imo = fPyjets->K[2][i];
1037 kst = fPyjets->K[0][i];
1038 pdg = fPyjets->K[1][i];
1042 // Quarks and gluons only
1043 if (pdg != 21 && TMath::Abs(pdg) > 6) continue;
1044 // Particles from hard scattering only
1046 if (imo > 8 && imo < 1000) imo = fPyjets->K[2][imo - 1];
1047 Int_t imom = imo % 1000;
1048 if ((isys == 0 || isys == 1) && ((imom != (isys + 7)))) continue;
1049 if ((isys == 2 || isys == 3) && ((imom != (isys + 1)))) continue;
1052 // Skip comment lines
1053 if (kst != 1 && kst != 2) continue;
1056 px = fPyjets->P[0][i];
1057 py = fPyjets->P[1][i];
1058 pz = fPyjets->P[2][i];
1059 e = fPyjets->P[3][i];
1060 m = fPyjets->P[4][i];
1061 pt = TMath::Sqrt(px * px + py * py);
1062 p = TMath::Sqrt(px * px + py * py + pz * pz);
1063 phi = TMath::Pi() + TMath::ATan2(-py, -px);
1064 theta = TMath::ATan2(pt, pz);
1067 // Save 4-momentum sum for balancing
1078 // Fractional energy loss
1079 Double_t z = zquench[index];
1082 // Don't fully quench radiated gluons
1085 // This small factor makes sure that the gluons are not too close in phase space to avoid recombination
1090 // printf("z: %d %f\n", imo, z);
1097 // Transform into frame in which initial parton is along z-axis
1099 TVector3 v(px, py, pz);
1100 v.RotateZ(-phiq[index]); v.RotateY(-thetaq[index]);
1101 Double_t pxs = v.X(); Double_t pys = v.Y(); Double_t pl = v.Z();
1103 Double_t jt = TMath::Sqrt(pxs * pxs + pys * pys);
1104 Double_t mt2 = jt * jt + m * m;
1107 // Kinematic limit on z
1109 if (m > 0.) zmax = 1. - m / TMath::Sqrt(m * m + jt * jt);
1111 // Change light-cone kinematics rel. to initial parton
1113 Double_t eppzOld = e + pl;
1114 Double_t empzOld = e - pl;
1116 Double_t eppzNew = (1. - z) * eppzOld;
1117 Double_t empzNew = empzOld - mt2 * z / eppzOld;
1118 Double_t eNew = 0.5 * (eppzNew + empzNew);
1119 Double_t plNew = 0.5 * (eppzNew - empzNew);
1123 // if mt very small (or sometimes even < 0 for numerical reasons) set it to 0
1124 Double_t mt2New = eppzNew * empzNew;
1125 if (mt2New < 1.e-8) mt2New = 0.;
1127 if (m * m > mt2New) {
1129 // This should not happen
1131 Fatal("Quench()", "This should never happen %e %e %e!", m, eppzNew, empzNew);
1134 jtNew = TMath::Sqrt(mt2New - m * m);
1137 // If pT is to small (probably a leading massive particle) we scale only the energy
1138 // This can cause negative masses of the radiated gluon
1139 // Let's hope for the best ...
1141 eNew = TMath::Sqrt(plNew * plNew + mt2);
1145 // Calculate new px, py
1151 pxNew = jtNew / jt * pxs;
1152 pyNew = jtNew / jt * pys;
1154 // Double_t dpx = pxs - pxNew;
1155 // Double_t dpy = pys - pyNew;
1156 // Double_t dpz = pl - plNew;
1157 // Double_t de = e - eNew;
1158 // Double_t dmass2 = de * de - dpx * dpx - dpy * dpy - dpz * dpz;
1159 // printf("New mass (1) %e %e %e %e %e %e %e \n", dmass2, jt, jtNew, pl, plNew, e, eNew);
1160 // printf("New mass (2) %e %e \n", pxNew, pyNew);
1164 TVector3 w(pxNew, pyNew, plNew);
1165 w.RotateY(thetaq[index]); w.RotateZ(phiq[index]);
1166 pxNew = w.X(); pyNew = w.Y(); plNew = w.Z();
1168 p1[index][0] += pxNew;
1169 p1[index][1] += pyNew;
1170 p1[index][2] += plNew;
1171 p1[index][3] += eNew;
1173 // Updated 4-momentum vectors
1175 pNew[icount][0] = pxNew;
1176 pNew[icount][1] = pyNew;
1177 pNew[icount][2] = plNew;
1178 pNew[icount][3] = eNew;
1183 // Check if there was phase-space for quenching
1186 if (icount == 0) quenched[isys] = kFALSE;
1187 if (!quenched[isys]) break;
1189 for (Int_t j = 0; j < 4; j++)
1191 p2[isys][j] = p0[isys][j] - p1[isys][j];
1193 p2[isys][4] = p2[isys][3] * p2[isys][3] - p2[isys][0] * p2[isys][0] - p2[isys][1] * p2[isys][1] - p2[isys][2] * p2[isys][2];
1194 if (p2[isys][4] > 0.) {
1195 p2[isys][4] = TMath::Sqrt(p2[isys][4]);
1198 printf("Warning negative mass squared in system %d %f ! \n", isys, zquench[isys]);
1199 printf("4-Momentum: %10.3e %10.3e %10.3e %10.3e %10.3e \n", p2[isys][0], p2[isys][1], p2[isys][2], p2[isys][3], p2[isys][4]);
1200 if (p2[isys][4] < -0.01) {
1201 printf("Negative mass squared !\n");
1202 // Here we have to put the gluon back to mass shell
1203 // This will lead to a small energy imbalance
1205 p2[isys][3] = TMath::Sqrt(p2[isys][0] * p2[isys][0] + p2[isys][1] * p2[isys][1] + p2[isys][2] * p2[isys][2]);
1214 printf("zHeavy lowered to %f\n", zHeavy);
1215 if (zHeavy < 0.01) {
1216 printf("No success ! \n");
1218 quenched[isys] = kFALSE;
1222 } // iteration on z (while)
1224 // Update event record
1225 for (Int_t k = 0; k < icount; k++) {
1226 // printf("%6d %6d %10.3e %10.3e %10.3e %10.3e\n", k, kNew[k], pNew[k][0],pNew[k][1], pNew[k][2], pNew[k][3] );
1227 fPyjets->P[0][kNew[k]] = pNew[k][0];
1228 fPyjets->P[1][kNew[k]] = pNew[k][1];
1229 fPyjets->P[2][kNew[k]] = pNew[k][2];
1230 fPyjets->P[3][kNew[k]] = pNew[k][3];
1237 if (!quenched[isys]) continue;
1239 // Last parton from shower i
1240 Int_t in = klast[isys];
1242 // Continue if no parton in shower i selected
1243 if (in == -1) continue;
1245 // If this is the second initial parton and it is behind the first move pointer by previous ish
1246 if (isys == 1 && klast[1] > klast[0]) in += ish;
1251 // How many additional gluons will be generated
1253 if (p2[isys][4] > 0.05) ish = 2;
1255 // Position of gluons
1257 if (iglu == 0) igMin = iGlu;
1260 (fPyjets->N) += ish;
1263 fPyjets->P[0][iGlu] = p2[isys][0];
1264 fPyjets->P[1][iGlu] = p2[isys][1];
1265 fPyjets->P[2][iGlu] = p2[isys][2];
1266 fPyjets->P[3][iGlu] = p2[isys][3];
1267 fPyjets->P[4][iGlu] = p2[isys][4];
1269 fPyjets->K[0][iGlu] = 1;
1270 if (iglu == nGluon[isys] - 1) fPyjets->K[0][iGlu] = 1;
1271 fPyjets->K[1][iGlu] = 21;
1272 fPyjets->K[2][iGlu] = fPyjets->K[2][in] + 1000;
1273 fPyjets->K[3][iGlu] = -1;
1274 fPyjets->K[4][iGlu] = -1;
1276 pg[0] += p2[isys][0];
1277 pg[1] += p2[isys][1];
1278 pg[2] += p2[isys][2];
1279 pg[3] += p2[isys][3];
1282 // Split gluon in rest frame.
1284 Double_t bx = p2[isys][0] / p2[isys][3];
1285 Double_t by = p2[isys][1] / p2[isys][3];
1286 Double_t bz = p2[isys][2] / p2[isys][3];
1287 Double_t pst = p2[isys][4] / 2.;
1289 // Isotropic decay ????
1290 Double_t cost = 2. * gRandom->Rndm() - 1.;
1291 Double_t sint = TMath::Sqrt((1.-cost)*(1.+cost));
1292 Double_t phis = 2. * TMath::Pi() * gRandom->Rndm();
1294 Double_t pz1 = pst * cost;
1295 Double_t pz2 = -pst * cost;
1296 Double_t pt1 = pst * sint;
1297 Double_t pt2 = -pst * sint;
1298 Double_t px1 = pt1 * TMath::Cos(phis);
1299 Double_t py1 = pt1 * TMath::Sin(phis);
1300 Double_t px2 = pt2 * TMath::Cos(phis);
1301 Double_t py2 = pt2 * TMath::Sin(phis);
1303 fPyjets->P[0][iGlu] = px1;
1304 fPyjets->P[1][iGlu] = py1;
1305 fPyjets->P[2][iGlu] = pz1;
1306 fPyjets->P[3][iGlu] = pst;
1307 fPyjets->P[4][iGlu] = 0.;
1309 fPyjets->K[0][iGlu] = 1 ;
1310 fPyjets->K[1][iGlu] = 21;
1311 fPyjets->K[2][iGlu] = fPyjets->K[2][in] + 1000;
1312 fPyjets->K[3][iGlu] = -1;
1313 fPyjets->K[4][iGlu] = -1;
1315 fPyjets->P[0][iGlu+1] = px2;
1316 fPyjets->P[1][iGlu+1] = py2;
1317 fPyjets->P[2][iGlu+1] = pz2;
1318 fPyjets->P[3][iGlu+1] = pst;
1319 fPyjets->P[4][iGlu+1] = 0.;
1321 fPyjets->K[0][iGlu+1] = 1;
1322 if (iglu == nGluon[isys] - 1) fPyjets->K[0][iGlu+1] = 1;
1323 fPyjets->K[1][iGlu+1] = 21;
1324 fPyjets->K[2][iGlu+1] = fPyjets->K[2][in] + 1000;
1325 fPyjets->K[3][iGlu+1] = -1;
1326 fPyjets->K[4][iGlu+1] = -1;
1332 Pyrobo(iGlu + 1, iGlu + 2, 0., 0., bx, by, bz);
1335 for (Int_t ig = iGlu; ig < iGlu+ish; ig++) {
1336 Double_t px, py, pz;
1337 px = fPyjets->P[0][ig];
1338 py = fPyjets->P[1][ig];
1339 pz = fPyjets->P[2][ig];
1340 TVector3 v(px, py, pz);
1341 v.RotateZ(-phiq[isys]);
1342 v.RotateY(-thetaq[isys]);
1343 Double_t pxs = v.X(); Double_t pys = v.Y(); Double_t pzs = v.Z();
1344 Double_t r = AliPythiaRndm::GetPythiaRandom()->Rndm();
1345 Double_t jtKick = 0.3 * TMath::Sqrt(-TMath::Log(r));
1346 if (ish == 2) jtKick = wjtKick[i] * TMath::Sqrt(-TMath::Log(r)) / TMath::Sqrt(2.);
1347 Double_t phiKick = 2. * TMath::Pi() * AliPythiaRndm::GetPythiaRandom()->Rndm();
1348 pxs += jtKick * TMath::Cos(phiKick);
1349 pys += jtKick * TMath::Sin(phiKick);
1350 TVector3 w(pxs, pys, pzs);
1351 w.RotateY(thetaq[isys]);
1352 w.RotateZ(phiq[isys]);
1353 fPyjets->P[0][ig] = w.X();
1354 fPyjets->P[1][ig] = w.Y();
1355 fPyjets->P[2][ig] = w.Z();
1356 fPyjets->P[2][ig] = w.Mag();
1362 // Check energy conservation
1366 Double_t es = 14000.;
1368 for (Int_t i = 0; i < numpart; i++)
1370 kst = fPyjets->K[0][i];
1371 if (kst != 1 && kst != 2) continue;
1372 pxs += fPyjets->P[0][i];
1373 pys += fPyjets->P[1][i];
1374 pzs += fPyjets->P[2][i];
1375 es -= fPyjets->P[3][i];
1377 if (TMath::Abs(pxs) > 1.e-2 ||
1378 TMath::Abs(pys) > 1.e-2 ||
1379 TMath::Abs(pzs) > 1.e-1) {
1380 printf("%e %e %e %e\n", pxs, pys, pzs, es);
1381 // Fatal("Quench()", "4-Momentum non-conservation");
1384 } // end quenching loop (systems)
1386 for (Int_t i = 0; i < numpart; i++)
1388 imo = fPyjets->K[2][i];
1390 fPyjets->K[2][i] = fPyjets->K[2][i] % 1000;
1397 void AliPythia::Pyquen(Double_t a, Int_t ibf, Double_t b)
1399 // Igor Lokthine's quenching routine
1400 // http://lokhtin.web.cern.ch/lokhtin/pyquen/pyquen.txt
1405 void AliPythia::SetPyquenParameters(Double_t t0, Double_t tau0, Int_t nf, Int_t iengl, Int_t iangl)
1407 // Set the parameters for the PYQUEN package.
1408 // See comments in PyquenCommon.h
1414 PYQPAR.iengl = iengl;
1415 PYQPAR.iangl = iangl;
1419 void AliPythia::Pyevnw()
1421 // New multiple interaction scenario
1425 void AliPythia::Pyshowq(Int_t ip1, Int_t ip2, Double_t qmax)
1427 // Call medium-modified Pythia jet reconstruction algorithm
1429 pyshowq(ip1, ip2, qmax);
1431 void AliPythia::Qpygin0()
1433 // New multiple interaction scenario
1437 void AliPythia::GetQuenchingParameters(Double_t& xp, Double_t& yp, Double_t z[4])
1439 // Return event specific quenching parameters
1442 for (Int_t i = 0; i < 4; i++) z[i] = fZQuench[i];
1446 void AliPythia::ConfigHeavyFlavor()
1449 // Default configuration for Heavy Flavor production
1451 // All QCD processes
1457 // No multiple interactions
1462 // Initial/final parton shower on (Pythia default)
1466 // 2nd order alpha_s
1474 void AliPythia::AtlasTuning()
1477 // Configuration for the ATLAS tuning
1478 if (fItune > -1) return;
1479 printf("ATLAS TUNE \n");
1481 SetMSTP(51, AliStructFuncType::PDFsetIndex(kCTEQ5L)); // CTEQ5L pdf
1482 SetMSTP(81,1); // Multiple Interactions ON
1483 SetMSTP(82,4); // Double Gaussian Model
1484 SetPARP(81,1.9); // Min. pt for multiple interactions (default in 6.2-14)
1485 SetPARP(82,1.8); // [GeV] PT_min at Ref. energy
1486 SetPARP(89,1000.); // [GeV] Ref. energy
1487 SetPARP(90,0.16); // 2*epsilon (exponent in power law)
1488 SetPARP(83,0.5); // Core density in proton matter distribution (def.value)
1489 SetPARP(84,0.5); // Core radius
1490 SetPARP(85,0.33); // Regulates gluon prod. mechanism
1491 SetPARP(86,0.66); // Regulates gluon prod. mechanism
1492 SetPARP(67,1); // Regulates Initial State Radiation
1495 void AliPythia::AtlasTuningMC09()
1498 // Configuration for the ATLAS tuning
1499 if (fItune > -1) return;
1500 printf("ATLAS New TUNE MC09\n");
1501 SetMSTP(81,21); // treatment for MI, ISR, FSR and beam remnants: MI on, new model
1502 SetMSTP(82, 4); // Double Gaussian Model
1503 SetMSTP(52, 2); // External PDF
1504 SetMSTP(51, 20650); // MRST LO*
1507 SetMSTP(70, 0); // (was 2: def manual 1, def code 0) virtuality scale for ISR
1508 SetMSTP(72, 1); // (was 0: def 1) maximum scale for FSR
1509 SetMSTP(88, 1); // (was 0: def 1) strategy for qq junction to di-quark or baryon in beam remnant
1510 SetMSTP(90, 0); // (was 1: def 0) strategy of compensate the primordial kT
1512 SetPARP(78, 0.3); // the amount of color reconnection in the final state
1513 SetPARP(80, 0.1); // probability of color partons kicked out from beam remnant
1514 SetPARP(82, 2.3); // [GeV] PT_min at Ref. energy
1515 SetPARP(83, 0.8); // Core density in proton matter distribution (def.value)
1516 SetPARP(84, 0.7); // Core radius
1517 SetPARP(90, 0.25); // 2*epsilon (exponent in power law)
1518 SetPARJ(81, 0.29); // (was 0.14: def 0.29) Labmda value in running alpha_s for parton showers
1521 SetPARJ(41, 0.3); // a and b parameters of the symmm. Lund FF
1523 SetPARJ(46, 0.75); // mod. of the Lund FF for heavy end-point quarks
1524 SetPARP(89,1800.); // [GeV] Ref. energy
1527 AliPythia& AliPythia::operator=(const AliPythia& rhs)
1529 // Assignment operator
1534 void AliPythia::Copy(TObject&) const
1539 Fatal("Copy","Not implemented!\n");
1542 void AliPythia::DalitzDecays()
1546 // Replace all omega dalitz decays with the correct matrix element decays
1548 Int_t nt = fPyjets->N;
1549 for (Int_t i = 0; i < nt; i++) {
1550 if (fPyjets->K[1][i] != 223) continue;
1551 Int_t fd = fPyjets->K[3][i] - 1;
1552 Int_t ld = fPyjets->K[4][i] - 1;
1553 if (fd < 0) continue;
1554 if ((ld - fd) != 2) continue;
1555 if ((fPyjets->K[1][fd] != 111) ||
1556 ((TMath::Abs(fPyjets->K[1][fd+1]) != 11) && (TMath::Abs(fPyjets->K[1][fd+1]) != 13)))
1558 TLorentzVector omega(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1559 Int_t pdg = TMath::Abs(fPyjets->K[1][fd+1]);
1560 fOmegaDalitz.Decay(pdg, &omega);
1561 for (Int_t j = 0; j < 3; j++) {
1562 for (Int_t k = 0; k < 4; k++) {
1563 TLorentzVector vec = (fOmegaDalitz.Products())[2-j];
1564 fPyjets->P[k][fd+j] = vec[k];
1572 // Replace all dalitz(pi0,eta,omega,eta',phi) and resonance(rho,omega,phi,jpsi) decays with the correct matrix element decays
1573 // for di-electron cocktail calculations
1577 void AliPythia::PizeroDalitz()
1580 Int_t nt = fPyjets->N;
1581 for (Int_t i = 0; i < nt; i++) {
1582 if (fPyjets->K[1][i] != 111) continue;
1583 Int_t fd = fPyjets->K[3][i] - 1;
1584 Int_t ld = fPyjets->K[4][i] - 1;
1585 if (fd < 0) continue;
1586 if ((ld - fd) != 2) continue;
1587 if ((fPyjets->K[1][fd] != 22) || (TMath::Abs(fPyjets->K[1][fd+1]) != 11) )
1589 TLorentzVector pizero(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1590 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1591 fExodus.Decay(pdg, &pizero);
1592 for (Int_t j = 0; j < 3; j++) {
1593 for (Int_t k = 0; k < 4; k++) {
1594 TLorentzVector vec = (fExodus.Products_pion())[2-j];
1595 fPyjets->P[k][fd+j] = vec[k];
1602 void AliPythia::EtaDalitz()
1605 Int_t nt = fPyjets->N;
1606 for (Int_t i = 0; i < nt; i++) {
1607 if (fPyjets->K[1][i] != 221) continue;
1608 Int_t fd = fPyjets->K[3][i] - 1;
1609 Int_t ld = fPyjets->K[4][i] - 1;
1610 if (fd < 0) continue;
1611 if ((ld - fd) != 2) continue;
1612 if ((fPyjets->K[1][fd] != 22) || (TMath::Abs(fPyjets->K[1][fd+1]) != 11))
1614 TLorentzVector eta(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1615 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1616 fExodus.Decay(pdg, &eta);
1617 for (Int_t j = 0; j < 3; j++) {
1618 for (Int_t k = 0; k < 4; k++) {
1619 TLorentzVector vec = (fExodus.Products_eta())[2-j];
1620 fPyjets->P[k][fd+j] = vec[k];
1627 void AliPythia::RhoDirect()
1630 Int_t nt = fPyjets->N;
1631 for (Int_t i = 0; i < nt; i++) {
1632 if (fPyjets->K[1][i] != 113) continue;
1633 Int_t fd = fPyjets->K[3][i] - 1;
1634 Int_t ld = fPyjets->K[4][i] - 1;
1635 if (fd < 0) continue;
1636 if ((ld - fd) != 1) continue;
1637 if ((TMath::Abs(fPyjets->K[1][fd]) != 11))
1639 TLorentzVector rho(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1640 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1641 fExodus.Decay(pdg, &rho);
1642 for (Int_t j = 0; j < 2; j++) {
1643 for (Int_t k = 0; k < 4; k++) {
1644 TLorentzVector vec = (fExodus.Products_rho())[1-j];
1645 fPyjets->P[k][fd+j] = vec[k];
1652 void AliPythia::OmegaDalitz()
1655 Int_t nt = fPyjets->N;
1656 for (Int_t i = 0; i < nt; i++) {
1657 if (fPyjets->K[1][i] != 223) continue;
1658 Int_t fd = fPyjets->K[3][i] - 1;
1659 Int_t ld = fPyjets->K[4][i] - 1;
1660 if (fd < 0) continue;
1661 if ((ld - fd) != 2) continue;
1662 if ((fPyjets->K[1][fd] != 111) || (TMath::Abs(fPyjets->K[1][fd+1]) != 11))
1664 TLorentzVector omegadalitz(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1665 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1666 fExodus.Decay(pdg, &omegadalitz);
1667 for (Int_t j = 0; j < 3; j++) {
1668 for (Int_t k = 0; k < 4; k++) {
1669 TLorentzVector vec = (fExodus.Products_omega_dalitz())[2-j];
1670 fPyjets->P[k][fd+j] = vec[k];
1677 void AliPythia::OmegaDirect()
1680 Int_t nt = fPyjets->N;
1681 for (Int_t i = 0; i < nt; i++) {
1682 if (fPyjets->K[1][i] != 223) continue;
1683 Int_t fd = fPyjets->K[3][i] - 1;
1684 Int_t ld = fPyjets->K[4][i] - 1;
1685 if (fd < 0) continue;
1686 if ((ld - fd) != 1) continue;
1687 if ((TMath::Abs(fPyjets->K[1][fd]) != 11))
1689 TLorentzVector omegadirect(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1690 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1691 fExodus.Decay(pdg, &omegadirect);
1692 for (Int_t j = 0; j < 2; j++) {
1693 for (Int_t k = 0; k < 4; k++) {
1694 TLorentzVector vec = (fExodus.Products_omega())[1-j];
1695 fPyjets->P[k][fd+j] = vec[k];
1702 void AliPythia::EtaprimeDalitz()
1705 Int_t nt = fPyjets->N;
1706 for (Int_t i = 0; i < nt; i++) {
1707 if (fPyjets->K[1][i] != 331) continue;
1708 Int_t fd = fPyjets->K[3][i] - 1;
1709 Int_t ld = fPyjets->K[4][i] - 1;
1710 if (fd < 0) continue;
1711 if ((ld - fd) != 2) continue;
1712 if ((fPyjets->K[1][fd] != 22) || (TMath::Abs(fPyjets->K[1][fd+1]) != 11))
1714 TLorentzVector etaprime(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1715 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1716 fExodus.Decay(pdg, &etaprime);
1717 for (Int_t j = 0; j < 3; j++) {
1718 for (Int_t k = 0; k < 4; k++) {
1719 TLorentzVector vec = (fExodus.Products_etaprime())[2-j];
1720 fPyjets->P[k][fd+j] = vec[k];
1727 void AliPythia::PhiDalitz()
1730 Int_t nt = fPyjets->N;
1731 for (Int_t i = 0; i < nt; i++) {
1732 if (fPyjets->K[1][i] != 333) continue;
1733 Int_t fd = fPyjets->K[3][i] - 1;
1734 Int_t ld = fPyjets->K[4][i] - 1;
1735 if (fd < 0) continue;
1736 if ((ld - fd) != 2) continue;
1737 if ((fPyjets->K[1][fd] != 221) || (TMath::Abs(fPyjets->K[1][fd+1]) != 11))
1739 TLorentzVector phidalitz(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1740 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1741 fExodus.Decay(pdg, &phidalitz);
1742 for (Int_t j = 0; j < 3; j++) {
1743 for (Int_t k = 0; k < 4; k++) {
1744 TLorentzVector vec = (fExodus.Products_phi_dalitz())[2-j];
1745 fPyjets->P[k][fd+j] = vec[k];
1752 void AliPythia::PhiDirect()
1755 Int_t nt = fPyjets->N;
1756 for (Int_t i = 0; i < nt; i++) {
1757 if (fPyjets->K[1][i] != 333) continue;
1758 Int_t fd = fPyjets->K[3][i] - 1;
1759 Int_t ld = fPyjets->K[4][i] - 1;
1760 if (fd < 0) continue;
1761 if ((ld - fd) != 1) continue;
1762 if ((TMath::Abs(fPyjets->K[1][fd]) != 11))
1764 TLorentzVector phi(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1765 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1766 fExodus.Decay(pdg, &phi);
1767 for (Int_t j = 0; j < 2; j++) {
1768 for (Int_t k = 0; k < 4; k++) {
1769 TLorentzVector vec = (fExodus.Products_phi())[1-j];
1770 fPyjets->P[k][fd+j] = vec[k];
1776 void AliPythia::JPsiDirect()
1779 Int_t nt = fPyjets->N;
1780 for (Int_t i = 0; i < nt; i++) {
1781 if (fPyjets->K[1][i] != 443) continue;
1782 Int_t fd = fPyjets->K[3][i] - 1;
1783 Int_t ld = fPyjets->K[4][i] - 1;
1784 if (fd < 0) continue;
1785 if ((ld - fd) != 1) continue;
1786 if ((TMath::Abs(fPyjets->K[1][fd]) != 11))
1788 TLorentzVector jpsi(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1789 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1790 fExodus.Decay(pdg, &jpsi);
1791 for (Int_t j = 0; j < 2; j++) {
1792 for (Int_t k = 0; k < 4; k++) {
1793 TLorentzVector vec = (fExodus.Products_jpsi())[1-j];
1794 fPyjets->P[k][fd+j] = vec[k];