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 **************************************************************************/
18 #include "AliPythia.h"
19 #include "AliPythiaRndm.h"
20 #include "AliFastGlauber.h"
21 #include "AliQuenchingWeights.h"
22 #include "AliOmegaDalitz.h"
23 #include "AliDecayerExodus.h"
26 #include "TLorentzVector.h"
27 #include "PyquenCommon.h"
32 # define pyclus pyclus_
33 # define pycell pycell_
34 # define pyshow pyshow_
35 # define pyrobo pyrobo_
36 # define pyquen pyquen_
37 # define pyevnw pyevnw_
38 # define pyshowq pyshowq_
39 # define qpygin0 qpygin0_
40 # define pytune pytune_
41 # define py2ent py2ent_
44 # define pyclus PYCLUS
45 # define pycell PYCELL
46 # define pyrobo PYROBO
47 # define pyquen PYQUEN
48 # define pyevnw PYEVNW
49 # define pyshowq PYSHOWQ
50 # define qpygin0 QPYGIN0
51 # define pytune PYTUNE
52 # define py2ent PY2ENT
53 # define type_of_call _stdcall
56 extern "C" void type_of_call pyclus(Int_t & );
57 extern "C" void type_of_call pycell(Int_t & );
58 extern "C" void type_of_call pyshow(Int_t &, Int_t &, Double_t &);
59 extern "C" void type_of_call pyrobo(Int_t &, Int_t &, Double_t &, Double_t &, Double_t &, Double_t &, Double_t &);
60 extern "C" void type_of_call pyquen(Double_t &, Int_t &, Double_t &);
61 extern "C" void type_of_call pyevnw();
62 extern "C" void type_of_call pyshowq(Int_t &, Int_t &, Double_t &);
63 extern "C" void type_of_call pytune(Int_t &);
64 extern "C" void type_of_call py2ent(Int_t &, Int_t&, Int_t&, Double_t&);
65 extern "C" void type_of_call qpygin0();
66 //_____________________________________________________________________________
68 AliPythia* AliPythia::fgAliPythia=NULL;
70 AliPythia::AliPythia():
86 // Default Constructor
89 if (!AliPythiaRndm::GetPythiaRandom())
90 AliPythiaRndm::SetPythiaRandom(GetRandom());
92 fQuenchingWeights = 0;
94 for (i = 0; i < 501; i++) fDefMDCY[i] = 0;
95 for (i = 0; i < 2001; i++) fDefMDME[i] = 0;
96 for (i = 0; i < 4; i++) fZQuench[i] = 0;
99 AliPythia::AliPythia(const AliPythia& pythia):
112 fQuenchingWeights(0),
119 for (i = 0; i < 501; i++) fDefMDCY[i] = 0;
120 for (i = 0; i < 2001; i++) fDefMDME[i] = 0;
121 for (i = 0; i < 4; i++) fZQuench[i] = 0;
125 void AliPythia::ProcInit(Process_t process, Float_t energy, StrucFunc_t strucfunc, Int_t itune)
127 // Initialise the process to generate
128 if (!AliPythiaRndm::GetPythiaRandom())
129 AliPythiaRndm::SetPythiaRandom(GetRandom());
135 fStrucFunc = strucfunc;
136 //...Switch off decay of pi0, K0S, Lambda, Sigma+-, Xi0-, Omega-.
137 SetMDCY(Pycomp(111) ,1,0); // pi0
138 SetMDCY(Pycomp(310) ,1,0); // K0S
139 SetMDCY(Pycomp(3122),1,0); // kLambda
140 SetMDCY(Pycomp(3112),1,0); // sigma -
141 SetMDCY(Pycomp(3222),1,0); // sigma +
142 SetMDCY(Pycomp(3312),1,0); // xi -
143 SetMDCY(Pycomp(3322),1,0); // xi 0
144 SetMDCY(Pycomp(3334),1,0); // omega-
145 // Select structure function
147 SetMSTP(51, AliStructFuncType::PDFsetIndex(strucfunc));
148 // Particles produced in string fragmentation point directly to either of the two endpoints
149 // of the string (depending in the side they were generated from).
153 // Pythia initialisation for selected processes//
157 for (Int_t i=1; i<= 200; i++) {
160 // select charm production
163 case kPyOldUEQ2ordered: //Old underlying events with Q2 ordered QCD processes
164 // Multiple interactions on.
166 // Double Gaussian matter distribution.
172 // Reference energy for pT0 and energy rescaling pace.
175 // String drawing almost completely minimizes string length.
178 // ISR and FSR activity.
184 case kPyOldUEQ2ordered2:
185 // Old underlying events with Q2 ordered QCD processes
186 // Multiple interactions on.
188 // Double Gaussian matter distribution.
194 // Reference energy for pT0 and energy rescaling pace.
196 SetPARP(90,0.16); // here is the difference with kPyOldUEQ2ordered
197 // String drawing almost completely minimizes string length.
200 // ISR and FSR activity.
207 // Old production mechanism: Old Popcorn
210 // (D=2) Like MSTJ(12)=2 but added prod ofthe 1er rank baryon
212 // (D=1)see can be used to form baryons (BARYON JUNCTION)
218 // heavy quark masses
248 case kPyCharmUnforced:
257 case kPyBeautyUnforced:
267 // Minimum Bias pp-Collisions
270 // select Pythia min. bias model
272 SetMSUB(92,1); // single diffraction AB-->XB
273 SetMSUB(93,1); // single diffraction AB-->AX
274 SetMSUB(94,1); // double diffraction
275 SetMSUB(95,1); // low pt production
280 case kPyMbAtlasTuneMC09:
281 // Minimum Bias pp-Collisions
284 // select Pythia min. bias model
286 SetMSUB(92,1); // single diffraction AB-->XB
287 SetMSUB(93,1); // single diffraction AB-->AX
288 SetMSUB(94,1); // double diffraction
289 SetMSUB(95,1); // low pt production
294 case kPyMbWithDirectPhoton:
295 // Minimum Bias pp-Collisions with direct photon processes added
298 // select Pythia min. bias model
300 SetMSUB(92,1); // single diffraction AB-->XB
301 SetMSUB(93,1); // single diffraction AB-->AX
302 SetMSUB(94,1); // double diffraction
303 SetMSUB(95,1); // low pt production
316 // Minimum Bias pp-Collisions
319 // select Pythia min. bias model
321 SetMSUB(92,1); // single diffraction AB-->XB
322 SetMSUB(93,1); // single diffraction AB-->AX
323 SetMSUB(94,1); // double diffraction
324 SetMSUB(95,1); // low pt production
327 // Les Houches Working Group 05 Minimum Bias pp-Collisions: hep-ph/0604120
328 // -> Pythia 6.3 or above is needed
331 SetMSUB(92,1); // single diffraction AB-->XB
332 SetMSUB(93,1); // single diffraction AB-->AX
333 SetMSUB(94,1); // double diffraction
334 SetMSUB(95,1); // low pt production
336 SetMSTP(51,AliStructFuncType::PDFsetIndex(kCTEQ6ll)); // CTEQ6ll pdf
340 SetMSTP(81,1); // Multiple Interactions ON
341 SetMSTP(82,4); // Double Gaussian Model
344 SetPARP(82,2.3); // [GeV] PT_min at Ref. energy
345 SetPARP(83,0.5); // Core density in proton matter distribution (def.value)
346 SetPARP(84,0.5); // Core radius
347 SetPARP(85,0.9); // Regulates gluon prod. mechanism
348 SetPARP(90,0.2); // 2*epsilon (exponent in power law)
352 // Minimum Bias pp-Collisions
355 // select Pythia min. bias model
357 SetMSUB(95,1); // low pt production
364 SetMSTP(91,1);// Width (1=gaussian) primordial kT dist. inside hadrons
365 SetPARP(91,1.); // <kT^2> = PARP(91,1.)^2
366 SetPARP(93,5.); // Upper cut-off
368 SetPMAS(4,1,1.2); // Charm quark mass
369 SetPMAS(5,1,4.78); // Beauty quark mass
370 SetPARP(71,4.); // Defaut value
379 // Pythia Tune A (CDF)
382 SetPARP(67,2.5); // Regulates Initial State Radiation (value from best fit to D0 dijet analysis)
383 SetMSTP(82,4); // Double Gaussian Model
384 SetPARP(82,2.0); // [GeV] PT_min at Ref. energy
385 SetPARP(84,0.4); // Core radius
386 SetPARP(85,0.90) ; // Regulates gluon prod. mechanism
387 SetPARP(86,0.95); // Regulates gluon prod. mechanism
388 SetPARP(89,1800.); // [GeV] Ref. energy
389 SetPARP(90,0.25); // 2*epsilon (exponent in power law)
395 case kPyCharmPbPbMNR:
397 case kPyDPlusPbPbMNR:
398 case kPyDPlusStrangePbPbMNR:
399 // Tuning of Pythia parameters aimed to get a resonable agreement
400 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
401 // c-cbar single inclusive and double differential distributions.
402 // This parameter settings are meant to work with Pb-Pb collisions
403 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
404 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
405 // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
417 case kPyDPlusStrangepPbMNR:
418 // Tuning of Pythia parameters aimed to get a resonable agreement
419 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
420 // c-cbar single inclusive and double differential distributions.
421 // This parameter settings are meant to work with p-Pb collisions
422 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
423 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
424 // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
437 case kPyDPlusStrangeppMNR:
438 case kPyLambdacppMNR:
439 // Tuning of Pythia parameters aimed to get a resonable agreement
440 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
441 // c-cbar single inclusive and double differential distributions.
442 // This parameter settings are meant to work with pp collisions
443 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
444 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
445 // has to be set to 2.1GeV. Example in ConfigCharmPPR.C.
455 case kPyCharmppMNRwmi:
456 // Tuning of Pythia parameters aimed to get a resonable agreement
457 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
458 // c-cbar single inclusive and double differential distributions.
459 // This parameter settings are meant to work with pp collisions
460 // and with kCTEQ5L PDFs.
461 // Added multiple interactions according to ATLAS tune settings.
462 // To get a "reasonable" agreement with MNR results, events have to be
463 // generated with the minimum ptHard (AliGenPythia::SetPtHard)
465 // To get a "perfect" agreement with MNR results, events have to be
466 // generated in four ptHard bins with the following relative
482 case kPyBeautyPbPbMNR:
483 // Tuning of Pythia parameters aimed to get a resonable agreement
484 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
485 // b-bbar single inclusive and double differential distributions.
486 // This parameter settings are meant to work with Pb-Pb collisions
487 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
488 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
489 // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
501 case kPyBeautypPbMNR:
502 // Tuning of Pythia parameters aimed to get a resonable agreement
503 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
504 // b-bbar single inclusive and double differential distributions.
505 // This parameter settings are meant to work with p-Pb collisions
506 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
507 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
508 // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
521 // Tuning of Pythia parameters aimed to get a resonable agreement
522 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
523 // b-bbar single inclusive and double differential distributions.
524 // This parameter settings are meant to work with pp collisions
525 // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs.
526 // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard)
527 // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C.
542 case kPyBeautyppMNRwmi:
543 // Tuning of Pythia parameters aimed to get a resonable agreement
544 // between with the NLO calculation by Mangano, Nason, Ridolfi for the
545 // b-bbar single inclusive and double differential distributions.
546 // This parameter settings are meant to work with pp collisions
547 // and with kCTEQ5L PDFs.
548 // Added multiple interactions according to ATLAS tune settings.
549 // To get a "reasonable" agreement with MNR results, events have to be
550 // generated with the minimum ptHard (AliGenPythia::SetPtHard)
552 // To get a "perfect" agreement with MNR results, events have to be
553 // generated in four ptHard bins with the following relative
576 //Inclusive production of W+/-
582 // //f fbar -> gamma W+
589 // Initial/final parton shower on (Pythia default)
590 // With parton showers on we are generating "W inclusive process"
591 SetMSTP(61,1); //Initial QCD & QED showers on
592 SetMSTP(71,1); //Final QCD & QED showers on
598 //Inclusive production of Z
603 // // f fbar -> g Z/gamma
605 // // f fbar -> gamma Z/gamma
607 // // f g -> f Z/gamma
609 // // f gamma -> f Z/gamma
612 //only Z included, not gamma
615 // Initial/final parton shower on (Pythia default)
616 // With parton showers on we are generating "Z inclusive process"
617 SetMSTP(61,1); //Initial QCD & QED showers on
618 SetMSTP(71,1); //Final QCD & QED showers on
622 //Inclusive production of Z
626 // Initial/final parton shower on (Pythia default)
627 // With parton showers on we are generating "Z inclusive process"
628 SetMSTP(61,1); //Initial QCD & QED showers on
629 SetMSTP(71,1); //Final QCD & QED showers on
631 case kPyMBRSingleDiffraction:
632 case kPyMBRDoubleDiffraction:
633 case kPyMBRCentralDiffraction:
638 // For the case of jet production the following parameter setting
639 // limits the transverse momentum of secondary scatterings, due
640 // to multiple parton interactions, to be less than that of the
641 // primary interaction (see POWHEG Dijet paper arXiv:1012.3380
642 // [hep-ph] sec. 4.1 and also the PYTHIA Manual).
645 // maximum number of errors before pythia aborts (def=10)
647 // number of warnings printed on the shell
654 // number of warnings printed on the shell
665 if (GetMSTP(192) > 1 || GetMSTP(193) > 1) {
666 AliWarning(Form("Structure function for tune %5d set to %5s\n",
667 itune, AliStructFuncType::PDFsetName(strucfunc).Data()));
669 SetMSTP(51, AliStructFuncType::PDFsetIndex(strucfunc));
673 SetMSTP(41,1); // all resonance decays switched on
674 if (process == kPyJetsPWHG || process == kPyCharmPWHG || process == kPyBeautyPWHG || process == kPyWPWHG) {
675 Initialize("USER","","",0.);
677 Initialize("CMS",fProjectile,fTarget,fEcms);
683 Int_t AliPythia::CheckedLuComp(Int_t kf)
685 // Check Lund particle code (for debugging)
687 printf("\n Lucomp kf,kc %d %d",kf,kc);
691 void AliPythia::SetNuclei(Int_t a1, Int_t a2, Int_t pdf)
693 // Treat protons as inside nuclei with mass numbers a1 and a2
694 // The MSTP array in the PYPARS common block is used to enable and
695 // select the nuclear structure functions.
696 // MSTP(52) : (D=1) choice of proton and nuclear structure-function library
697 // =1: internal PYTHIA acording to MSTP(51)
698 // =2: PDFLIB proton s.f., with MSTP(51) = 1000xNGROUP+NSET
699 // If the following mass number both not equal zero, nuclear corrections of the stf are used.
700 // MSTP(192) : Mass number of nucleus side 1
701 // MSTP(193) : Mass number of nucleus side 2
702 // MSTP(194) : Nuclear structure function: 0: EKS98 8:EPS08 9:EPS09LO 19:EPS09NLO
710 AliPythia* AliPythia::Instance()
712 // Set random number generator
716 fgAliPythia = new AliPythia();
721 void AliPythia::PrintParticles()
723 // Print list of particl properties
725 char* name = new char[16];
726 for (Int_t kf=0; kf<1000000; kf++) {
727 for (Int_t c = 1; c > -2; c-=2) {
728 Int_t kc = Pycomp(c*kf);
730 Float_t mass = GetPMAS(kc,1);
731 Float_t width = GetPMAS(kc,2);
732 Float_t tau = GetPMAS(kc,4);
738 printf("\n mass, width, tau: %6d %s %10.3f %10.3e %10.3e",
739 c*kf, name, mass, width, tau);
743 printf("\n Number of particles %d \n \n", np);
746 void AliPythia::ResetDecayTable()
748 // Set default values for pythia decay switches
750 for (i = 1; i < 501; i++) SetMDCY(i,1,fDefMDCY[i]);
751 for (i = 1; i < 2001; i++) SetMDME(i,1,fDefMDME[i]);
754 void AliPythia::SetDecayTable()
756 // Set default values for pythia decay switches
759 for (i = 1; i < 501; i++) fDefMDCY[i] = GetMDCY(i,1);
760 for (i = 1; i < 2001; i++) fDefMDME[i] = GetMDME(i,1);
763 void AliPythia::Pyclus(Int_t& njet)
765 // Call Pythia clustering algorithm
770 void AliPythia::Pycell(Int_t& njet)
772 // Call Pythia jet reconstruction algorithm
777 void AliPythia::Pyshow(Int_t ip1, Int_t ip2, Double_t qmax)
779 // Call Pythia jet reconstruction algorithm
781 pyshow(ip1, ip2, qmax);
784 void AliPythia::Pyrobo(Int_t imi, Int_t ima, Double_t the, Double_t phi, Double_t bex, Double_t bey, Double_t bez)
786 pyrobo(imi, ima, the, phi, bex, bey, bez);
789 void AliPythia::Pytune(Int_t itune)
793 C ITUNE NAME (detailed descriptions below)
794 C 0 Default : No settings changed => linked Pythia version's defaults.
795 C ====== Old UE, Q2-ordered showers ==========================================
796 C 100 A : Rick Field's CDF Tune A
797 C 101 AW : Rick Field's CDF Tune AW
798 C 102 BW : Rick Field's CDF Tune BW
799 C 103 DW : Rick Field's CDF Tune DW
800 C 104 DWT : Rick Field's CDF Tune DW with slower UE energy scaling
801 C 105 QW : Rick Field's CDF Tune QW (NB: needs CTEQ6.1M pdfs externally)
802 C 106 ATLAS-DC2: Arthur Moraes' (old) ATLAS tune (ATLAS DC2 / Rome)
803 C 107 ACR : Tune A modified with annealing CR
804 C 108 D6 : Rick Field's CDF Tune D6 (NB: needs CTEQ6L pdfs externally)
805 C 109 D6T : Rick Field's CDF Tune D6T (NB: needs CTEQ6L pdfs externally)
806 C ====== Intermediate Models =================================================
807 C 200 IM 1 : Intermediate model: new UE, Q2-ordered showers, annealing CR
808 C 201 APT : Tune A modified to use pT-ordered final-state showers
809 C ====== New UE, interleaved pT-ordered showers, annealing CR ================
810 C 300 S0 : Sandhoff-Skands Tune 0
811 C 301 S1 : Sandhoff-Skands Tune 1
812 C 302 S2 : Sandhoff-Skands Tune 2
813 C 303 S0A : S0 with "Tune A" UE energy scaling
814 C 304 NOCR : New UE "best try" without colour reconnections
815 C 305 Old : New UE, original (primitive) colour reconnections
816 C 306 ATLAS-CSC: Arthur Moraes' (new) ATLAS tune (needs CTEQ6L externally)
817 C ======= The Uppsala models =================================================
818 C ( NB! must be run with special modified Pythia 6.215 version )
819 C ( available from http://www.isv.uu.se/thep/MC/scigal/ )
820 C 400 GAL 0 : Generalized area-law model. Old parameters
821 C 401 SCI 0 : Soft-Colour-Interaction model. Old parameters
822 C 402 GAL 1 : Generalized area-law model. Tevatron MB retuned (Skands)
827 void AliPythia::Py2ent(Int_t idx, Int_t pdg1, Int_t pdg2, Double_t p){
828 // Inset 2-parton system at line idx
829 py2ent(idx, pdg1, pdg2, p);
833 void AliPythia::InitQuenching(Float_t cMin, Float_t cMax, Float_t k, Int_t iECMethod, Float_t zmax, Int_t ngmax)
836 // (1) The quenching model using quenching weights according to C. Salgado and U. Wiedemann
837 // (2) The nuclear geometry using the Glauber Model
840 fGlauber = AliFastGlauber::Instance();
842 fGlauber->SetCentralityClass(cMin, cMax);
844 fQuenchingWeights = new AliQuenchingWeights();
845 fQuenchingWeights->InitMult();
846 fQuenchingWeights->SetK(k);
847 fQuenchingWeights->SetECMethod(AliQuenchingWeights::kECMethod(iECMethod));
854 void AliPythia::Quench()
858 // Simple Jet Quenching routine:
859 // =============================
860 // The jet formed by all final state partons radiated by the parton created
861 // in the hard collisions is quenched by a factor (1-z) using light cone variables in
862 // the initial parton reference frame:
863 // (E + p_z)new = (1-z) (E + p_z)old
868 // The lost momentum is first balanced by one gluon with virtuality > 0.
869 // Subsequently the gluon splits to yield two gluons with E = p.
873 static Float_t eMean = 0.;
874 static Int_t icall = 0;
879 Int_t klast[4] = {-1, -1, -1, -1};
881 Int_t numpart = fPyjets->N;
882 Double_t px = 0., py = 0., pz = 0., e = 0., m = 0., p = 0., pt = 0., theta = 0., phi = 0.;
883 Double_t pxq[4], pyq[4], pzq[4], eq[4], yq[4], mq[4], pq[4], phiq[4], thetaq[4], ptq[4];
885 Double_t wjtKick[4] = {0., 0., 0., 0.};
891 // Sore information about Primary partons
894 // 0, 1 partons from hard scattering
895 // 2, 3 partons from initial state radiation
897 for (Int_t i = 2; i <= 7; i++) {
899 // Skip gluons that participate in hard scattering
900 if (i == 4 || i == 5) continue;
901 // Gluons from hard Scattering
902 if (i == 6 || i == 7) {
904 pxq[j] = fPyjets->P[0][i];
905 pyq[j] = fPyjets->P[1][i];
906 pzq[j] = fPyjets->P[2][i];
907 eq[j] = fPyjets->P[3][i];
908 mq[j] = fPyjets->P[4][i];
910 // Gluons from initial state radiation
912 // Obtain 4-momentum vector from difference between original parton and parton after gluon
913 // radiation. Energy is calculated independently because initial state radition does not
914 // conserve strictly momentum and energy for each partonic system independently.
916 // Not very clean. Should be improved !
920 pxq[j] = fPyjets->P[0][i] - fPyjets->P[0][i+2];
921 pyq[j] = fPyjets->P[1][i] - fPyjets->P[1][i+2];
922 pzq[j] = fPyjets->P[2][i] - fPyjets->P[2][i+2];
923 mq[j] = fPyjets->P[4][i];
924 eq[j] = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j] + pzq[j] * pzq[j] + mq[j] * mq[j]);
927 // Calculate some kinematic variables
929 yq[j] = 0.5 * TMath::Log((eq[j] + pzq[j] + 1.e-14) / (eq[j] - pzq[j] + 1.e-14));
930 pq[j] = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j] + pzq[j] * pzq[j]);
931 phiq[j] = TMath::Pi()+TMath::ATan2(-pyq[j], -pxq[j]);
932 ptq[j] = TMath::Sqrt(pxq[j] * pxq[j] + pyq[j] * pyq[j]);
933 thetaq[j] = TMath::ATan2(ptq[j], pzq[j]);
934 qPdg[j] = fPyjets->K[1][i];
940 fGlauber->GetI0I1ForPythiaAndXY(4, phiq, int0, int1, fXJet, fYJet, 15.);
942 for (Int_t j = 0; j < 4; j++) {
944 // Quench only central jets and with E > 10.
948 Int_t itype = (qPdg[j] == 21) ? 2 : 1;
949 Double_t eloss = fQuenchingWeights->GetELossRandomKFast(itype, int0[j], int1[j], eq[j]);
951 if (TMath::Abs(yq[j]) > 2.5 || eq[j] < 10.) {
954 if (eq[j] > 40. && TMath::Abs(yq[j]) < 0.5) {
960 Double_t l = fQuenchingWeights->CalcLk(int0[j], int1[j]);
961 wjtKick[j] = TMath::Sqrt(l * fQuenchingWeights->CalcQk(int0[j], int1[j]));
963 // Fractional energy loss
964 fZQuench[j] = eloss / eq[j];
966 // Avoid complete loss
968 if (fZQuench[j] > fZmax) fZQuench[j] = fZmax;
970 // Some debug printing
973 // 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",
974 // j, itype, eq[j], phiq[j], l, eloss, wjtKick[j], eMean / Float_t(icall+1), yq[j]);
976 // fZQuench[j] = 0.8;
977 // while (fZQuench[j] >= 0.95) fZQuench[j] = gRandom->Exp(0.2);
980 quenched[j] = (fZQuench[j] > 0.01);
985 Double_t pNew[1000][4];
992 for (Int_t isys = 0; isys < 4; isys++) {
993 // Skip to next system if not quenched.
994 if (!quenched[isys]) continue;
996 nGluon[isys] = 1 + Int_t(fZQuench[isys] / (1. - fZQuench[isys]));
997 if (nGluon[isys] > fNGmax) nGluon[isys] = fNGmax;
998 zquench[isys] = 1. - TMath::Power(1. - fZQuench[isys], 1./Double_t(nGluon[isys]));
999 wjtKick[isys] = wjtKick[isys] / TMath::Sqrt(Double_t(nGluon[isys]));
1005 Double_t pg[4] = {0., 0., 0., 0.};
1008 // Loop on radiation events
1010 for (Int_t iglu = 0; iglu < nGluon[isys]; iglu++) {
1013 for (Int_t k = 0; k < 4; k++)
1019 // Loop over partons
1020 for (Int_t i = 0; i < numpart; i++)
1022 imo = fPyjets->K[2][i];
1023 kst = fPyjets->K[0][i];
1024 pdg = fPyjets->K[1][i];
1028 // Quarks and gluons only
1029 if (pdg != 21 && TMath::Abs(pdg) > 6) continue;
1030 // Particles from hard scattering only
1032 if (imo > 8 && imo < 1000) imo = fPyjets->K[2][imo - 1];
1033 Int_t imom = imo % 1000;
1034 if ((isys == 0 || isys == 1) && ((imom != (isys + 7)))) continue;
1035 if ((isys == 2 || isys == 3) && ((imom != (isys + 1)))) continue;
1038 // Skip comment lines
1039 if (kst != 1 && kst != 2) continue;
1042 px = fPyjets->P[0][i];
1043 py = fPyjets->P[1][i];
1044 pz = fPyjets->P[2][i];
1045 e = fPyjets->P[3][i];
1046 m = fPyjets->P[4][i];
1047 pt = TMath::Sqrt(px * px + py * py);
1048 p = TMath::Sqrt(px * px + py * py + pz * pz);
1049 phi = TMath::Pi() + TMath::ATan2(-py, -px);
1050 theta = TMath::ATan2(pt, pz);
1053 // Save 4-momentum sum for balancing
1064 // Fractional energy loss
1065 Double_t z = zquench[index];
1068 // Don't fully quench radiated gluons
1071 // This small factor makes sure that the gluons are not too close in phase space to avoid recombination
1076 // printf("z: %d %f\n", imo, z);
1083 // Transform into frame in which initial parton is along z-axis
1085 TVector3 v(px, py, pz);
1086 v.RotateZ(-phiq[index]); v.RotateY(-thetaq[index]);
1087 Double_t pxs = v.X(); Double_t pys = v.Y(); Double_t pl = v.Z();
1089 Double_t jt = TMath::Sqrt(pxs * pxs + pys * pys);
1090 Double_t mt2 = jt * jt + m * m;
1093 // Kinematic limit on z
1095 if (m > 0.) zmax = 1. - m / TMath::Sqrt(m * m + jt * jt);
1097 // Change light-cone kinematics rel. to initial parton
1099 Double_t eppzOld = e + pl;
1100 Double_t empzOld = e - pl;
1102 Double_t eppzNew = (1. - z) * eppzOld;
1103 Double_t empzNew = empzOld - mt2 * z / eppzOld;
1104 Double_t eNew = 0.5 * (eppzNew + empzNew);
1105 Double_t plNew = 0.5 * (eppzNew - empzNew);
1109 // if mt very small (or sometimes even < 0 for numerical reasons) set it to 0
1110 Double_t mt2New = eppzNew * empzNew;
1111 if (mt2New < 1.e-8) mt2New = 0.;
1113 if (m * m > mt2New) {
1115 // This should not happen
1117 Fatal("Quench()", "This should never happen %e %e %e!", m, eppzNew, empzNew);
1120 jtNew = TMath::Sqrt(mt2New - m * m);
1123 // If pT is to small (probably a leading massive particle) we scale only the energy
1124 // This can cause negative masses of the radiated gluon
1125 // Let's hope for the best ...
1127 eNew = TMath::Sqrt(plNew * plNew + mt2);
1131 // Calculate new px, py
1137 pxNew = jtNew / jt * pxs;
1138 pyNew = jtNew / jt * pys;
1140 // Double_t dpx = pxs - pxNew;
1141 // Double_t dpy = pys - pyNew;
1142 // Double_t dpz = pl - plNew;
1143 // Double_t de = e - eNew;
1144 // Double_t dmass2 = de * de - dpx * dpx - dpy * dpy - dpz * dpz;
1145 // printf("New mass (1) %e %e %e %e %e %e %e \n", dmass2, jt, jtNew, pl, plNew, e, eNew);
1146 // printf("New mass (2) %e %e \n", pxNew, pyNew);
1150 TVector3 w(pxNew, pyNew, plNew);
1151 w.RotateY(thetaq[index]); w.RotateZ(phiq[index]);
1152 pxNew = w.X(); pyNew = w.Y(); plNew = w.Z();
1154 p1[index][0] += pxNew;
1155 p1[index][1] += pyNew;
1156 p1[index][2] += plNew;
1157 p1[index][3] += eNew;
1159 // Updated 4-momentum vectors
1161 pNew[icount][0] = pxNew;
1162 pNew[icount][1] = pyNew;
1163 pNew[icount][2] = plNew;
1164 pNew[icount][3] = eNew;
1169 // Check if there was phase-space for quenching
1172 if (icount == 0) quenched[isys] = kFALSE;
1173 if (!quenched[isys]) break;
1175 for (Int_t j = 0; j < 4; j++)
1177 p2[isys][j] = p0[isys][j] - p1[isys][j];
1179 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];
1180 if (p2[isys][4] > 0.) {
1181 p2[isys][4] = TMath::Sqrt(p2[isys][4]);
1184 printf("Warning negative mass squared in system %d %f ! \n", isys, zquench[isys]);
1185 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]);
1186 if (p2[isys][4] < -0.01) {
1187 printf("Negative mass squared !\n");
1188 // Here we have to put the gluon back to mass shell
1189 // This will lead to a small energy imbalance
1191 p2[isys][3] = TMath::Sqrt(p2[isys][0] * p2[isys][0] + p2[isys][1] * p2[isys][1] + p2[isys][2] * p2[isys][2]);
1200 printf("zHeavy lowered to %f\n", zHeavy);
1201 if (zHeavy < 0.01) {
1202 printf("No success ! \n");
1204 quenched[isys] = kFALSE;
1208 } // iteration on z (while)
1210 // Update event record
1211 for (Int_t k = 0; k < icount; k++) {
1212 // 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] );
1213 fPyjets->P[0][kNew[k]] = pNew[k][0];
1214 fPyjets->P[1][kNew[k]] = pNew[k][1];
1215 fPyjets->P[2][kNew[k]] = pNew[k][2];
1216 fPyjets->P[3][kNew[k]] = pNew[k][3];
1223 if (!quenched[isys]) continue;
1225 // Last parton from shower i
1226 Int_t in = klast[isys];
1228 // Continue if no parton in shower i selected
1229 if (in == -1) continue;
1231 // If this is the second initial parton and it is behind the first move pointer by previous ish
1232 if (isys == 1 && klast[1] > klast[0]) in += ish;
1237 // How many additional gluons will be generated
1239 if (p2[isys][4] > 0.05) ish = 2;
1241 // Position of gluons
1243 if (iglu == 0) igMin = iGlu;
1246 (fPyjets->N) += ish;
1249 fPyjets->P[0][iGlu] = p2[isys][0];
1250 fPyjets->P[1][iGlu] = p2[isys][1];
1251 fPyjets->P[2][iGlu] = p2[isys][2];
1252 fPyjets->P[3][iGlu] = p2[isys][3];
1253 fPyjets->P[4][iGlu] = p2[isys][4];
1255 fPyjets->K[0][iGlu] = 1;
1256 if (iglu == nGluon[isys] - 1) fPyjets->K[0][iGlu] = 1;
1257 fPyjets->K[1][iGlu] = 21;
1258 fPyjets->K[2][iGlu] = fPyjets->K[2][in] + 1000;
1259 fPyjets->K[3][iGlu] = -1;
1260 fPyjets->K[4][iGlu] = -1;
1262 pg[0] += p2[isys][0];
1263 pg[1] += p2[isys][1];
1264 pg[2] += p2[isys][2];
1265 pg[3] += p2[isys][3];
1268 // Split gluon in rest frame.
1270 Double_t bx = p2[isys][0] / p2[isys][3];
1271 Double_t by = p2[isys][1] / p2[isys][3];
1272 Double_t bz = p2[isys][2] / p2[isys][3];
1273 Double_t pst = p2[isys][4] / 2.;
1275 // Isotropic decay ????
1276 Double_t cost = 2. * gRandom->Rndm() - 1.;
1277 Double_t sint = TMath::Sqrt((1.-cost)*(1.+cost));
1278 Double_t phis = 2. * TMath::Pi() * gRandom->Rndm();
1280 Double_t pz1 = pst * cost;
1281 Double_t pz2 = -pst * cost;
1282 Double_t pt1 = pst * sint;
1283 Double_t pt2 = -pst * sint;
1284 Double_t px1 = pt1 * TMath::Cos(phis);
1285 Double_t py1 = pt1 * TMath::Sin(phis);
1286 Double_t px2 = pt2 * TMath::Cos(phis);
1287 Double_t py2 = pt2 * TMath::Sin(phis);
1289 fPyjets->P[0][iGlu] = px1;
1290 fPyjets->P[1][iGlu] = py1;
1291 fPyjets->P[2][iGlu] = pz1;
1292 fPyjets->P[3][iGlu] = pst;
1293 fPyjets->P[4][iGlu] = 0.;
1295 fPyjets->K[0][iGlu] = 1 ;
1296 fPyjets->K[1][iGlu] = 21;
1297 fPyjets->K[2][iGlu] = fPyjets->K[2][in] + 1000;
1298 fPyjets->K[3][iGlu] = -1;
1299 fPyjets->K[4][iGlu] = -1;
1301 fPyjets->P[0][iGlu+1] = px2;
1302 fPyjets->P[1][iGlu+1] = py2;
1303 fPyjets->P[2][iGlu+1] = pz2;
1304 fPyjets->P[3][iGlu+1] = pst;
1305 fPyjets->P[4][iGlu+1] = 0.;
1307 fPyjets->K[0][iGlu+1] = 1;
1308 if (iglu == nGluon[isys] - 1) fPyjets->K[0][iGlu+1] = 1;
1309 fPyjets->K[1][iGlu+1] = 21;
1310 fPyjets->K[2][iGlu+1] = fPyjets->K[2][in] + 1000;
1311 fPyjets->K[3][iGlu+1] = -1;
1312 fPyjets->K[4][iGlu+1] = -1;
1318 Pyrobo(iGlu + 1, iGlu + 2, 0., 0., bx, by, bz);
1321 for (Int_t ig = iGlu; ig < iGlu+ish; ig++) {
1322 Double_t px, py, pz;
1323 px = fPyjets->P[0][ig];
1324 py = fPyjets->P[1][ig];
1325 pz = fPyjets->P[2][ig];
1326 TVector3 v(px, py, pz);
1327 v.RotateZ(-phiq[isys]);
1328 v.RotateY(-thetaq[isys]);
1329 Double_t pxs = v.X(); Double_t pys = v.Y(); Double_t pzs = v.Z();
1330 Double_t r = AliPythiaRndm::GetPythiaRandom()->Rndm();
1331 Double_t jtKick = 0.3 * TMath::Sqrt(-TMath::Log(r));
1332 if (ish == 2) jtKick = wjtKick[i] * TMath::Sqrt(-TMath::Log(r)) / TMath::Sqrt(2.);
1333 Double_t phiKick = 2. * TMath::Pi() * AliPythiaRndm::GetPythiaRandom()->Rndm();
1334 pxs += jtKick * TMath::Cos(phiKick);
1335 pys += jtKick * TMath::Sin(phiKick);
1336 TVector3 w(pxs, pys, pzs);
1337 w.RotateY(thetaq[isys]);
1338 w.RotateZ(phiq[isys]);
1339 fPyjets->P[0][ig] = w.X();
1340 fPyjets->P[1][ig] = w.Y();
1341 fPyjets->P[2][ig] = w.Z();
1342 fPyjets->P[2][ig] = w.Mag();
1348 // Check energy conservation
1352 Double_t es = 14000.;
1354 for (Int_t i = 0; i < numpart; i++)
1356 kst = fPyjets->K[0][i];
1357 if (kst != 1 && kst != 2) continue;
1358 pxs += fPyjets->P[0][i];
1359 pys += fPyjets->P[1][i];
1360 pzs += fPyjets->P[2][i];
1361 es -= fPyjets->P[3][i];
1363 if (TMath::Abs(pxs) > 1.e-2 ||
1364 TMath::Abs(pys) > 1.e-2 ||
1365 TMath::Abs(pzs) > 1.e-1) {
1366 printf("%e %e %e %e\n", pxs, pys, pzs, es);
1367 // Fatal("Quench()", "4-Momentum non-conservation");
1370 } // end quenching loop (systems)
1372 for (Int_t i = 0; i < numpart; i++)
1374 imo = fPyjets->K[2][i];
1376 fPyjets->K[2][i] = fPyjets->K[2][i] % 1000;
1383 void AliPythia::Pyquen(Double_t a, Int_t ibf, Double_t b)
1385 // Igor Lokthine's quenching routine
1386 // http://lokhtin.web.cern.ch/lokhtin/pyquen/pyquen.txt
1391 void AliPythia::SetPyquenParameters(Double_t t0, Double_t tau0, Int_t nf, Int_t iengl, Int_t iangl)
1393 // Set the parameters for the PYQUEN package.
1394 // See comments in PyquenCommon.h
1400 PYQPAR.iengl = iengl;
1401 PYQPAR.iangl = iangl;
1405 void AliPythia::Pyevnw()
1407 // New multiple interaction scenario
1411 void AliPythia::Pyshowq(Int_t ip1, Int_t ip2, Double_t qmax)
1413 // Call medium-modified Pythia jet reconstruction algorithm
1415 pyshowq(ip1, ip2, qmax);
1417 void AliPythia::Qpygin0()
1419 // New multiple interaction scenario
1423 void AliPythia::GetQuenchingParameters(Double_t& xp, Double_t& yp, Double_t z[4])
1425 // Return event specific quenching parameters
1428 for (Int_t i = 0; i < 4; i++) z[i] = fZQuench[i];
1432 void AliPythia::ConfigHeavyFlavor()
1435 // Default configuration for Heavy Flavor production
1437 // All QCD processes
1443 // No multiple interactions
1448 // Initial/final parton shower on (Pythia default)
1452 // 2nd order alpha_s
1460 void AliPythia::AtlasTuning()
1463 // Configuration for the ATLAS tuning
1464 if (fItune > -1) return;
1465 printf("ATLAS TUNE \n");
1467 SetMSTP(51, AliStructFuncType::PDFsetIndex(kCTEQ5L)); // CTEQ5L pdf
1468 SetMSTP(81,1); // Multiple Interactions ON
1469 SetMSTP(82,4); // Double Gaussian Model
1470 SetPARP(81,1.9); // Min. pt for multiple interactions (default in 6.2-14)
1471 SetPARP(82,1.8); // [GeV] PT_min at Ref. energy
1472 SetPARP(89,1000.); // [GeV] Ref. energy
1473 SetPARP(90,0.16); // 2*epsilon (exponent in power law)
1474 SetPARP(83,0.5); // Core density in proton matter distribution (def.value)
1475 SetPARP(84,0.5); // Core radius
1476 SetPARP(85,0.33); // Regulates gluon prod. mechanism
1477 SetPARP(86,0.66); // Regulates gluon prod. mechanism
1478 SetPARP(67,1); // Regulates Initial State Radiation
1481 void AliPythia::AtlasTuningMC09()
1484 // Configuration for the ATLAS tuning
1485 if (fItune > -1) return;
1486 printf("ATLAS New TUNE MC09\n");
1487 SetMSTP(81,21); // treatment for MI, ISR, FSR and beam remnants: MI on, new model
1488 SetMSTP(82, 4); // Double Gaussian Model
1489 SetMSTP(52, 2); // External PDF
1490 SetMSTP(51, 20650); // MRST LO*
1493 SetMSTP(70, 0); // (was 2: def manual 1, def code 0) virtuality scale for ISR
1494 SetMSTP(72, 1); // (was 0: def 1) maximum scale for FSR
1495 SetMSTP(88, 1); // (was 0: def 1) strategy for qq junction to di-quark or baryon in beam remnant
1496 SetMSTP(90, 0); // (was 1: def 0) strategy of compensate the primordial kT
1498 SetPARP(78, 0.3); // the amount of color reconnection in the final state
1499 SetPARP(80, 0.1); // probability of color partons kicked out from beam remnant
1500 SetPARP(82, 2.3); // [GeV] PT_min at Ref. energy
1501 SetPARP(83, 0.8); // Core density in proton matter distribution (def.value)
1502 SetPARP(84, 0.7); // Core radius
1503 SetPARP(90, 0.25); // 2*epsilon (exponent in power law)
1504 SetPARJ(81, 0.29); // (was 0.14: def 0.29) Labmda value in running alpha_s for parton showers
1507 SetPARJ(41, 0.3); // a and b parameters of the symmm. Lund FF
1509 SetPARJ(46, 0.75); // mod. of the Lund FF for heavy end-point quarks
1510 SetPARP(89,1800.); // [GeV] Ref. energy
1513 AliPythia& AliPythia::operator=(const AliPythia& rhs)
1515 // Assignment operator
1520 void AliPythia::Copy(TObject&) const
1525 Fatal("Copy","Not implemented!\n");
1528 void AliPythia::DalitzDecays()
1532 // Replace all omega dalitz decays with the correct matrix element decays
1534 Int_t nt = fPyjets->N;
1535 for (Int_t i = 0; i < nt; i++) {
1536 if (fPyjets->K[1][i] != 223) continue;
1537 Int_t fd = fPyjets->K[3][i] - 1;
1538 Int_t ld = fPyjets->K[4][i] - 1;
1539 if (fd < 0) continue;
1540 if ((ld - fd) != 2) continue;
1541 if ((fPyjets->K[1][fd] != 111) ||
1542 ((TMath::Abs(fPyjets->K[1][fd+1]) != 11) && (TMath::Abs(fPyjets->K[1][fd+1]) != 13)))
1544 TLorentzVector omega(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1545 Int_t pdg = TMath::Abs(fPyjets->K[1][fd+1]);
1546 fOmegaDalitz.Decay(pdg, &omega);
1547 for (Int_t j = 0; j < 3; j++) {
1548 for (Int_t k = 0; k < 4; k++) {
1549 TLorentzVector vec = (fOmegaDalitz.Products())[2-j];
1550 fPyjets->P[k][fd+j] = vec[k];
1558 // Replace all dalitz(pi0,eta,omega,eta',phi) and resonance(rho,omega,phi,jpsi) decays with the correct matrix element decays
1559 // for di-electron cocktail calculations
1563 void AliPythia::PizeroDalitz()
1566 Int_t nt = fPyjets->N;
1567 for (Int_t i = 0; i < nt; i++) {
1568 if (fPyjets->K[1][i] != 111) continue;
1569 Int_t fd = fPyjets->K[3][i] - 1;
1570 Int_t ld = fPyjets->K[4][i] - 1;
1571 if (fd < 0) continue;
1572 if ((ld - fd) != 2) continue;
1573 if ((fPyjets->K[1][fd] != 22) || (TMath::Abs(fPyjets->K[1][fd+1]) != 11) )
1575 TLorentzVector pizero(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1576 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1577 fExodus.Decay(pdg, &pizero);
1578 for (Int_t j = 0; j < 3; j++) {
1579 for (Int_t k = 0; k < 4; k++) {
1580 TLorentzVector vec = (fExodus.Products_pion())[2-j];
1581 fPyjets->P[k][fd+j] = vec[k];
1588 void AliPythia::EtaDalitz()
1591 Int_t nt = fPyjets->N;
1592 for (Int_t i = 0; i < nt; i++) {
1593 if (fPyjets->K[1][i] != 221) continue;
1594 Int_t fd = fPyjets->K[3][i] - 1;
1595 Int_t ld = fPyjets->K[4][i] - 1;
1596 if (fd < 0) continue;
1597 if ((ld - fd) != 2) continue;
1598 if ((fPyjets->K[1][fd] != 22) || (TMath::Abs(fPyjets->K[1][fd+1]) != 11))
1600 TLorentzVector eta(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1601 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1602 fExodus.Decay(pdg, &eta);
1603 for (Int_t j = 0; j < 3; j++) {
1604 for (Int_t k = 0; k < 4; k++) {
1605 TLorentzVector vec = (fExodus.Products_eta())[2-j];
1606 fPyjets->P[k][fd+j] = vec[k];
1613 void AliPythia::RhoDirect()
1616 Int_t nt = fPyjets->N;
1617 for (Int_t i = 0; i < nt; i++) {
1618 if (fPyjets->K[1][i] != 113) continue;
1619 Int_t fd = fPyjets->K[3][i] - 1;
1620 Int_t ld = fPyjets->K[4][i] - 1;
1621 if (fd < 0) continue;
1622 if ((ld - fd) != 1) continue;
1623 if ((TMath::Abs(fPyjets->K[1][fd]) != 11))
1625 TLorentzVector rho(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1626 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1627 fExodus.Decay(pdg, &rho);
1628 for (Int_t j = 0; j < 2; j++) {
1629 for (Int_t k = 0; k < 4; k++) {
1630 TLorentzVector vec = (fExodus.Products_rho())[1-j];
1631 fPyjets->P[k][fd+j] = vec[k];
1638 void AliPythia::OmegaDalitz()
1641 Int_t nt = fPyjets->N;
1642 for (Int_t i = 0; i < nt; i++) {
1643 if (fPyjets->K[1][i] != 223) continue;
1644 Int_t fd = fPyjets->K[3][i] - 1;
1645 Int_t ld = fPyjets->K[4][i] - 1;
1646 if (fd < 0) continue;
1647 if ((ld - fd) != 2) continue;
1648 if ((fPyjets->K[1][fd] != 111) || (TMath::Abs(fPyjets->K[1][fd+1]) != 11))
1650 TLorentzVector omegadalitz(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1651 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1652 fExodus.Decay(pdg, &omegadalitz);
1653 for (Int_t j = 0; j < 3; j++) {
1654 for (Int_t k = 0; k < 4; k++) {
1655 TLorentzVector vec = (fExodus.Products_omega_dalitz())[2-j];
1656 fPyjets->P[k][fd+j] = vec[k];
1663 void AliPythia::OmegaDirect()
1666 Int_t nt = fPyjets->N;
1667 for (Int_t i = 0; i < nt; i++) {
1668 if (fPyjets->K[1][i] != 223) continue;
1669 Int_t fd = fPyjets->K[3][i] - 1;
1670 Int_t ld = fPyjets->K[4][i] - 1;
1671 if (fd < 0) continue;
1672 if ((ld - fd) != 1) continue;
1673 if ((TMath::Abs(fPyjets->K[1][fd]) != 11))
1675 TLorentzVector omegadirect(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1676 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1677 fExodus.Decay(pdg, &omegadirect);
1678 for (Int_t j = 0; j < 2; j++) {
1679 for (Int_t k = 0; k < 4; k++) {
1680 TLorentzVector vec = (fExodus.Products_omega())[1-j];
1681 fPyjets->P[k][fd+j] = vec[k];
1688 void AliPythia::EtaprimeDalitz()
1691 Int_t nt = fPyjets->N;
1692 for (Int_t i = 0; i < nt; i++) {
1693 if (fPyjets->K[1][i] != 331) continue;
1694 Int_t fd = fPyjets->K[3][i] - 1;
1695 Int_t ld = fPyjets->K[4][i] - 1;
1696 if (fd < 0) continue;
1697 if ((ld - fd) != 2) continue;
1698 if ((fPyjets->K[1][fd] != 22) || (TMath::Abs(fPyjets->K[1][fd+1]) != 11))
1700 TLorentzVector etaprime(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1701 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1702 fExodus.Decay(pdg, &etaprime);
1703 for (Int_t j = 0; j < 3; j++) {
1704 for (Int_t k = 0; k < 4; k++) {
1705 TLorentzVector vec = (fExodus.Products_etaprime())[2-j];
1706 fPyjets->P[k][fd+j] = vec[k];
1713 void AliPythia::PhiDalitz()
1716 Int_t nt = fPyjets->N;
1717 for (Int_t i = 0; i < nt; i++) {
1718 if (fPyjets->K[1][i] != 333) continue;
1719 Int_t fd = fPyjets->K[3][i] - 1;
1720 Int_t ld = fPyjets->K[4][i] - 1;
1721 if (fd < 0) continue;
1722 if ((ld - fd) != 2) continue;
1723 if ((fPyjets->K[1][fd] != 221) || (TMath::Abs(fPyjets->K[1][fd+1]) != 11))
1725 TLorentzVector phidalitz(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1726 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1727 fExodus.Decay(pdg, &phidalitz);
1728 for (Int_t j = 0; j < 3; j++) {
1729 for (Int_t k = 0; k < 4; k++) {
1730 TLorentzVector vec = (fExodus.Products_phi_dalitz())[2-j];
1731 fPyjets->P[k][fd+j] = vec[k];
1738 void AliPythia::PhiDirect()
1741 Int_t nt = fPyjets->N;
1742 for (Int_t i = 0; i < nt; i++) {
1743 if (fPyjets->K[1][i] != 333) continue;
1744 Int_t fd = fPyjets->K[3][i] - 1;
1745 Int_t ld = fPyjets->K[4][i] - 1;
1746 if (fd < 0) continue;
1747 if ((ld - fd) != 1) continue;
1748 if ((TMath::Abs(fPyjets->K[1][fd]) != 11))
1750 TLorentzVector phi(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1751 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1752 fExodus.Decay(pdg, &phi);
1753 for (Int_t j = 0; j < 2; j++) {
1754 for (Int_t k = 0; k < 4; k++) {
1755 TLorentzVector vec = (fExodus.Products_phi())[1-j];
1756 fPyjets->P[k][fd+j] = vec[k];
1762 void AliPythia::JPsiDirect()
1765 Int_t nt = fPyjets->N;
1766 for (Int_t i = 0; i < nt; i++) {
1767 if (fPyjets->K[1][i] != 443) continue;
1768 Int_t fd = fPyjets->K[3][i] - 1;
1769 Int_t ld = fPyjets->K[4][i] - 1;
1770 if (fd < 0) continue;
1771 if ((ld - fd) != 1) continue;
1772 if ((TMath::Abs(fPyjets->K[1][fd]) != 11))
1774 TLorentzVector jpsi(fPyjets->P[0][i], fPyjets->P[1][i], fPyjets->P[2][i], fPyjets->P[3][i]);
1775 Int_t pdg = TMath::Abs(fPyjets->K[1][i]);
1776 fExodus.Decay(pdg, &jpsi);
1777 for (Int_t j = 0; j < 2; j++) {
1778 for (Int_t k = 0; k < 4; k++) {
1779 TLorentzVector vec = (fExodus.Products_jpsi())[1-j];
1780 fPyjets->P[k][fd+j] = vec[k];