/************************************************************************** * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ /* $Id$ */ #include "AliPythia.h" #include "AliPythiaRndm.h" ClassImp(AliPythia) #ifndef WIN32 # define pyclus pyclus_ # define pycell pycell_ # define pyshow pyshow_ # define pyrobo pyrobo_ # define type_of_call #else # define pyclus PYCLUS # define pycell PYCELL # define pyrobo PYROBO # define type_of_call _stdcall #endif extern "C" void type_of_call pyclus(Int_t & ); extern "C" void type_of_call pycell(Int_t & ); extern "C" void type_of_call pyshow(Int_t &, Int_t &, Double_t &); extern "C" void type_of_call pyrobo(Int_t &, Int_t &, Double_t &, Double_t &, Double_t &, Double_t &, Double_t &); //_____________________________________________________________________________ AliPythia* AliPythia::fgAliPythia=NULL; AliPythia::AliPythia() { // Default Constructor // // Set random number if (!AliPythiaRndm::GetPythiaRandom()) AliPythiaRndm::SetPythiaRandom(GetRandom()); } void AliPythia::ProcInit(Process_t process, Float_t energy, StrucFunc_t strucfunc) { // Initialise the process to generate if (!AliPythiaRndm::GetPythiaRandom()) AliPythiaRndm::SetPythiaRandom(GetRandom()); fProcess = process; fEcms = energy; fStrucFunc = strucfunc; // don't decay p0 SetMDCY(Pycomp(111),1,0); // select structure function SetMSTP(52,2); SetMSTP(51,strucfunc); // // Pythia initialisation for selected processes// // // Make MSEL clean // for (Int_t i=1; i<= 200; i++) { SetMSUB(i,0); } // select charm production switch (process) { case kPyCharm: SetMSEL(4); // // heavy quark masses SetPMAS(4,1,1.2); SetMSTU(16,2); // // primordial pT SetMSTP(91,1); SetPARP(91,1.); SetPARP(93,5.); // break; case kPyBeauty: SetMSEL(5); SetPMAS(5,1,4.75); SetMSTU(16,2); break; case kPyJpsi: SetMSEL(0); // gg->J/Psi g SetMSUB(86,1); break; case kPyJpsiChi: SetMSEL(0); // gg->J/Psi g SetMSUB(86,1); // gg-> chi_0c g SetMSUB(87,1); // gg-> chi_1c g SetMSUB(88,1); // gg-> chi_2c g SetMSUB(89,1); break; case kPyCharmUnforced: SetMSEL(0); // gq->qg SetMSUB(28,1); // gg->qq SetMSUB(53,1); // gg->gg SetMSUB(68,1); break; case kPyBeautyUnforced: SetMSEL(0); // gq->qg SetMSUB(28,1); // gg->qq SetMSUB(53,1); // gg->gg SetMSUB(68,1); break; case kPyMb: // Minimum Bias pp-Collisions // // // select Pythia min. bias model SetMSEL(0); SetMSUB(92,1); // single diffraction AB-->XB SetMSUB(93,1); // single diffraction AB-->AX SetMSUB(94,1); // double diffraction SetMSUB(95,1); // low pt production SetMSTP(81,1); // multiple interactions switched on SetMSTP(82,3); // model with varying impact param. & a single Gaussian SetPARP(82,3.47); // set value pT_0 for turn-off of the cross section of // multiple interaction at a reference energy = 14000 GeV SetPARP(89,14000.); // reference energy for the above parameter SetPARP(90,0.174); // set exponent for energy dependence of pT_0 case kPyMbNonDiffr: // Minimum Bias pp-Collisions // // // select Pythia min. bias model SetMSEL(0); SetMSUB(95,1); // low pt production SetMSTP(81,1); // multiple interactions switched on SetMSTP(82,3); // model with varying impact param. & a single Gaussian SetPARP(82,3.47); // set value pT_0 for turn-off of the cross section of // multiple interaction at a reference energy = 14000 GeV SetPARP(89,14000.); // reference energy for the above parameter SetPARP(90,0.174); // set exponent for energy dependence of pT_0 break; case kPyJets: // // QCD Jets // SetMSEL(1); break; case kPyDirectGamma: SetMSEL(10); break; case kPyCharmPbPbMNR: case kPyD0PbPbMNR: // Tuning of Pythia parameters aimed to get a resonable agreement // between with the NLO calculation by Mangano, Nason, Ridolfi for the // c-cbar single inclusive and double differential distributions. // This parameter settings are meant to work with Pb-Pb collisions // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs. // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard) // has to be set to 2.1GeV. Example in ConfigCharmPPR.C. // All QCD processes SetMSEL(1); // No multiple interactions SetMSTP(81,0); SetPARP(81,0.0); SetPARP(82,0.0); // Initial/final parton shower on (Pythia default) SetMSTP(61,1); SetMSTP(71,1); // 2nd order alpha_s SetMSTP(2,2); // QCD scales SetMSTP(32,2); SetPARP(34,1.0); // Intrinsic SetMSTP(91,1); SetPARP(91,1.304); SetPARP(93,6.52); // Set c-quark mass SetPMAS(4,1,1.2); break; case kPyCharmpPbMNR: case kPyD0pPbMNR: // Tuning of Pythia parameters aimed to get a resonable agreement // between with the NLO calculation by Mangano, Nason, Ridolfi for the // c-cbar single inclusive and double differential distributions. // This parameter settings are meant to work with p-Pb collisions // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs. // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard) // has to be set to 2.1GeV. Example in ConfigCharmPPR.C. // All QCD processes SetMSEL(1); // No multiple interactions SetMSTP(81,0); SetPARP(81,0.0); SetPARP(82,0.0); // Initial/final parton shower on (Pythia default) SetMSTP(61,1); SetMSTP(71,1); // 2nd order alpha_s SetMSTP(2,2); // QCD scales SetMSTP(32,2); SetPARP(34,1.0); // Intrinsic SetMSTP(91,1); SetPARP(91,1.16); SetPARP(93,5.8); // Set c-quark mass SetPMAS(4,1,1.2); break; case kPyCharmppMNR: case kPyD0ppMNR: // Tuning of Pythia parameters aimed to get a resonable agreement // between with the NLO calculation by Mangano, Nason, Ridolfi for the // c-cbar single inclusive and double differential distributions. // This parameter settings are meant to work with pp collisions // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs. // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard) // has to be set to 2.1GeV. Example in ConfigCharmPPR.C. // All QCD processes SetMSEL(1); // No multiple interactions SetMSTP(81,0); SetPARP(81,0.0); SetPARP(82,0.0); // Initial/final parton shower on (Pythia default) SetMSTP(61,1); SetMSTP(71,1); // 2nd order alpha_s SetMSTP(2,2); // QCD scales SetMSTP(32,2); SetPARP(34,1.0); // Intrinsic SetMSTP(91,1); SetPARP(91,1.); SetPARP(93,5.); // Set c-quark mass SetPMAS(4,1,1.2); break; case kPyBeautyPbPbMNR: // Tuning of Pythia parameters aimed to get a resonable agreement // between with the NLO calculation by Mangano, Nason, Ridolfi for the // b-bbar single inclusive and double differential distributions. // This parameter settings are meant to work with Pb-Pb collisions // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs. // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard) // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C. // All QCD processes SetMSEL(1); // No multiple interactions SetMSTP(81,0); SetPARP(81,0.0); SetPARP(82,0.0); // Initial/final parton shower on (Pythia default) SetMSTP(61,1); SetMSTP(71,1); // 2nd order alpha_s SetMSTP(2,2); // QCD scales SetMSTP(32,2); SetPARP(34,1.0); SetPARP(67,1.0); SetPARP(71,1.0); // Intrinsic SetMSTP(91,1); SetPARP(91,2.035); SetPARP(93,10.17); // Set b-quark mass SetPMAS(5,1,4.75); break; case kPyBeautypPbMNR: // Tuning of Pythia parameters aimed to get a resonable agreement // between with the NLO calculation by Mangano, Nason, Ridolfi for the // b-bbar single inclusive and double differential distributions. // This parameter settings are meant to work with p-Pb collisions // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs. // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard) // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C. // All QCD processes SetMSEL(1); // No multiple interactions SetMSTP(81,0); SetPARP(81,0.0); SetPARP(82,0.0); // Initial/final parton shower on (Pythia default) SetMSTP(61,1); SetMSTP(71,1); // 2nd order alpha_s SetMSTP(2,2); // QCD scales SetMSTP(32,2); SetPARP(34,1.0); SetPARP(67,1.0); SetPARP(71,1.0); // Intrinsic SetMSTP(91,1); SetPARP(91,1.60); SetPARP(93,8.00); // Set b-quark mass SetPMAS(5,1,4.75); break; case kPyBeautyppMNR: // Tuning of Pythia parameters aimed to get a resonable agreement // between with the NLO calculation by Mangano, Nason, Ridolfi for the // b-bbar single inclusive and double differential distributions. // This parameter settings are meant to work with pp collisions // (AliGenPythia::SetNuclei) and with kCTEQ4L PDFs. // To get a good agreement the minimum ptHard (AliGenPythia::SetPtHard) // has to be set to 2.75GeV. Example in ConfigBeautyPPR.C. // All QCD processes SetMSEL(1); // No multiple interactions SetMSTP(81,0); SetPARP(81,0.0); SetPARP(82,0.0); // Initial/final parton shower on (Pythia default) SetMSTP(61,1); SetMSTP(71,1); // 2nd order alpha_s SetMSTP(2,2); // QCD scales SetMSTP(32,2); SetPARP(34,1.0); SetPARP(67,1.0); SetPARP(71,1.0); // Intrinsic SetMSTP(91,1); SetPARP(91,1.); SetPARP(93,5.); // Set b-quark mass SetPMAS(5,1,4.75); break; } // // Initialize PYTHIA SetMSTP(41,1); // all resonance decays switched on Initialize("CMS","p","p",fEcms); } Int_t AliPythia::CheckedLuComp(Int_t kf) { // Check Lund particle code (for debugging) Int_t kc=Pycomp(kf); printf("\n Lucomp kf,kc %d %d",kf,kc); return kc; } void AliPythia::SetNuclei(Int_t a1, Int_t a2) { // Treat protons as inside nuclei with mass numbers a1 and a2 // The MSTP array in the PYPARS common block is used to enable and // select the nuclear structure functions. // MSTP(52) : (D=1) choice of proton and nuclear structure-function library // =1: internal PYTHIA acording to MSTP(51) // =2: PDFLIB proton s.f., with MSTP(51) = 1000xNGROUP+NSET // If the following mass number both not equal zero, nuclear corrections of the stf are used. // MSTP(192) : Mass number of nucleus side 1 // MSTP(193) : Mass number of nucleus side 2 SetMSTP(52,2); SetMSTP(192, a1); SetMSTP(193, a2); } AliPythia* AliPythia::Instance() { // Set random number generator if (fgAliPythia) { return fgAliPythia; } else { fgAliPythia = new AliPythia(); return fgAliPythia; } } void AliPythia::PrintParticles() { // Print list of particl properties Int_t np = 0; char* name = new char[16]; for (Int_t kf=0; kf<1000000; kf++) { for (Int_t c = 1; c > -2; c-=2) { Int_t kc = Pycomp(c*kf); if (kc) { Float_t mass = GetPMAS(kc,1); Float_t width = GetPMAS(kc,2); Float_t tau = GetPMAS(kc,4); Pyname(kf,name); np++; printf("\n mass, width, tau: %6d %s %10.3f %10.3e %10.3e", c*kf, name, mass, width, tau); } } } printf("\n Number of particles %d \n \n", np); } void AliPythia::ResetDecayTable() { // Set default values for pythia decay switches Int_t i; for (i = 1; i < 501; i++) SetMDCY(i,1,fDefMDCY[i]); for (i = 1; i < 2001; i++) SetMDME(i,1,fDefMDME[i]); } void AliPythia::SetDecayTable() { // Set default values for pythia decay switches // Int_t i; for (i = 1; i < 501; i++) fDefMDCY[i] = GetMDCY(i,1); for (i = 1; i < 2001; i++) fDefMDME[i] = GetMDME(i,1); } void AliPythia::Pyclus(Int_t& njet) { // Call Pythia clustering algorithm // pyclus(njet); } void AliPythia::Pycell(Int_t& njet) { // Call Pythia jet reconstruction algorithm // pycell(njet); } void AliPythia::Pyshow(Int_t ip1, Int_t ip2, Double_t qmax) { // Call Pythia jet reconstruction algorithm // Int_t numpart = fPyjets->N; for (Int_t i = 0; i < numpart; i++) { if (fPyjets->K[2][i] == 7) ip1 = i+1; if (fPyjets->K[2][i] == 8) ip2 = i+1; } qmax = 2. * GetVINT(51); printf("Pyshow %d %d %f", ip1, ip2, qmax); pyshow(ip1, ip2, qmax); } void AliPythia::Pyrobo(Int_t imi, Int_t ima, Double_t the, Double_t phi, Double_t bex, Double_t bey, Double_t bez) { pyrobo(imi, ima, the, phi, bex, bey, bez); } void AliPythia::Quench() { // // // Simple Jet Quenching routine: // ============================= // The jet formed by all final state partons radiated by the parton created // in the hard collisions is quenched by a factor z using: // (E + p_z)new = (1-z) (E + p_z)old // // The lost momentum is first balanced by one gluon with virtuality > 0. // Subsequently the gluon splits to yield two gluons with E = p. // Float_t p0[2][5]; Float_t p1[2][5]; Float_t p2[2][5]; Int_t klast[2] = {-1, -1}; Int_t kglu[2]; for (Int_t i = 0; i < 4; i++) { p0[0][i] = 0.; p0[1][i] = 0.; p1[0][i] = 0.; p1[1][i] = 0.; p2[0][i] = 0.; p2[1][i] = 0.; } Int_t numpart = fPyjets->N; for (Int_t i = 0; i < numpart; i++) { Int_t imo = fPyjets->K[2][i]; Int_t kst = fPyjets->K[0][i]; Int_t pdg = fPyjets->K[1][i]; // Quarks and gluons only if (pdg != 21 && TMath::Abs(pdg) > 6) continue; // Particles from hard scattering only Float_t px = fPyjets->P[0][i]; Float_t py = fPyjets->P[1][i]; Float_t pz = fPyjets->P[2][i]; Float_t e = fPyjets->P[3][i]; Float_t m = fPyjets->P[4][i]; Float_t pt = TMath::Sqrt(px * px + py * py); // Skip comment lines if (kst != 1 && kst != 2) continue; Float_t mt = TMath::Sqrt(px * px + py * py + m * m); // // Some cuts to be in a save kinematic region // if (imo != 7 && imo != 8) continue; Int_t index = imo - 7; klast[index] = i; p0[index][0] += px; p0[index][1] += py; p0[index][2] += pz; p0[index][3] += e; // // Fix z // Float_t z = 0.2; Float_t eppzOld = e + pz; Float_t empzOld = e - pz; // // Kinematics of the original parton // Float_t eppzNew = (1. - z) * eppzOld; Float_t empzNew = empzOld - mt * mt * z / eppzOld; Float_t eNew0 = 0.5 * (eppzNew + empzNew); Float_t pzNew0 = 0.5 * (eppzNew - empzNew); // // Skip if pt too small // if (m * m > eppzNew * empzNew) continue; Float_t ptNew = TMath::Sqrt(eppzNew * empzNew - m * m); Float_t pxNew0 = ptNew / pt * px; Float_t pyNew0 = ptNew / pt * py; p1[index][0] += pxNew0; p1[index][1] += pyNew0; p1[index][2] += pzNew0; p1[index][3] += eNew0; // // Update event record // fPyjets->P[0][i] = pxNew0; fPyjets->P[1][i] = pyNew0; fPyjets->P[2][i] = pzNew0; fPyjets->P[3][i] = eNew0; } // // Gluons // for (Int_t k = 0; k < 2; k++) { for (Int_t j = 0; j < 4; j++) { p2[k][j] = p0[k][j] - p1[k][j]; } p2[k][4] = p2[k][3] * p2[k][3] - p2[k][0] * p2[k][0] - p2[k][1] * p2[k][1] - p2[k][2] * p2[k][2]; if (p2[k][4] > 0.) { // // Bring gluon back to mass shell via momentum scaling // (momentum will not be conserved, but energy) // // not used anymore /* Float_t psq = p2[k][0] * p2[k][0] + p2[k][1] * p2[k][1] + p2[k][2] * p2[k][2]; Float_t fact = TMath::Sqrt(1. + p2[k][4] / psq); p2[k][0] *= fact; p2[k][1] *= fact; p2[k][2] *= fact; p2[k][3] = TMath::Sqrt(psq) * fact; p2[k][4] = 0.; */ } } if (p2[0][4] > 0.) { p2[0][4] = TMath::Sqrt(p2[0][4]); } else { printf("Warning negative mass squared ! \n"); } if (p2[1][4] > 0.) { p2[1][4] = TMath::Sqrt(p2[1][4]); } else { printf("Warning negative mass squared ! \n"); } // // Add the gluons // for (Int_t i = 0; i < 2; i++) { Int_t ish, jmin, jmax, iGlu, iNew; Int_t in = klast[i]; ish = 0; if (in == -1) continue; if (i == 1 && klast[1] > klast[0]) in += ish; jmin = in - 1; ish = 1; if (p2[i][4] > 0) ish = 2; iGlu = in; iNew = in + ish; jmax = numpart + ish - 1; if (fPyjets->K[0][in-1] == 1 || fPyjets->K[0][in-1] == 21 || fPyjets->K[0][in-1] == 11) { jmin = in; iGlu = in + 1; iNew = in; } kglu[i] = iGlu; for (Int_t j = jmax; j > jmin; j--) { for (Int_t k = 0; k < 5; k++) { fPyjets->K[k][j] = fPyjets->K[k][j-ish]; fPyjets->P[k][j] = fPyjets->P[k][j-ish]; fPyjets->V[k][j] = fPyjets->V[k][j-ish]; } } // end shifting numpart += ish; (fPyjets->N) += ish; if (ish == 1) { fPyjets->P[0][iGlu] = p2[i][0]; fPyjets->P[1][iGlu] = p2[i][1]; fPyjets->P[2][iGlu] = p2[i][2]; fPyjets->P[3][iGlu] = p2[i][3]; fPyjets->P[4][iGlu] = p2[i][4]; fPyjets->K[0][iGlu] = 2; fPyjets->K[1][iGlu] = 21; fPyjets->K[2][iGlu] = fPyjets->K[2][iNew]; fPyjets->K[3][iGlu] = -1; fPyjets->K[4][iGlu] = -1; } else { // // Split gluon in rest frame. // Double_t bx = p2[i][0] / p2[i][3]; Double_t by = p2[i][1] / p2[i][3]; Double_t bz = p2[i][2] / p2[i][3]; Float_t pst = p2[i][4] / 2.; Float_t cost = 2. * gRandom->Rndm() - 1.; Float_t sint = TMath::Sqrt(1. - cost * cost); Float_t phi = 2. * TMath::Pi() * gRandom->Rndm(); Float_t pz1 = pst * cost; Float_t pz2 = -pst * cost; Float_t pt1 = pst * sint; Float_t pt2 = -pst * sint; Float_t px1 = pt1 * TMath::Cos(phi); Float_t py1 = pt1 * TMath::Sin(phi); Float_t px2 = pt2 * TMath::Cos(phi); Float_t py2 = pt2 * TMath::Sin(phi); fPyjets->P[0][iGlu] = px1; fPyjets->P[1][iGlu] = py1; fPyjets->P[2][iGlu] = pz1; fPyjets->P[3][iGlu] = pst; fPyjets->P[4][iGlu] = 0.; fPyjets->K[0][iGlu] = 2; fPyjets->K[1][iGlu] = 21; fPyjets->K[2][iGlu] = fPyjets->K[2][iNew]; fPyjets->K[3][iGlu] = -1; fPyjets->K[4][iGlu] = -1; fPyjets->P[0][iGlu+1] = px2; fPyjets->P[1][iGlu+1] = py2; fPyjets->P[2][iGlu+1] = pz2; fPyjets->P[3][iGlu+1] = pst; fPyjets->P[4][iGlu+1] = 0.; fPyjets->K[0][iGlu+1] = 2; fPyjets->K[1][iGlu+1] = 21; fPyjets->K[2][iGlu+1] = fPyjets->K[2][iNew]; fPyjets->K[3][iGlu+1] = -1; fPyjets->K[4][iGlu+1] = -1; SetMSTU(1,0); SetMSTU(2,0); // // Boost back // Pyrobo(iGlu + 1, iGlu + 2, 0., 0., bx, by, bz); } } // end adding gluons } // end quench