8 #include "Pythia8/Pythia.h"
25 const double Pythia::VERSIONNUMBERCODE = 8.186;
33 const int Pythia::NTRY = 10;
36 const int Pythia::SUBRUNDEFAULT = -999;
42 Pythia::Pythia(
string xmlDir,
bool printBanner) {
47 useNewPdfHard =
false;
48 useNewPdfPomA =
false;
49 useNewPdfPomB =
false;
63 couplingsPtr = &couplings;
73 hasMergingHooks =
false;
74 hasOwnMergingHooks =
false;
78 useNewBeamShape =
false;
91 const char* PYTHIA8DATA =
"PYTHIA8DATA";
92 char* envPath = getenv(PYTHIA8DATA);
93 if (envPath != 0 && *envPath !=
'\0') {
95 while (*(envPath+i) !=
'\0') xmlPath += *(envPath+(i++));
97 else xmlPath = xmlDir;
98 if (xmlPath[ xmlPath.length() - 1 ] !=
'/') xmlPath +=
"/";
101 settings.initPtr( &info);
102 string initFile = xmlPath +
"Index.xml";
103 isConstructed = settings.init( initFile);
104 if (!isConstructed) {
105 info.errorMsg(
"Abort from Pythia::Pythia: settings unavailable");
110 double versionNumberXML = parm(
"Pythia:versionNumber");
111 isConstructed = (abs(versionNumberXML - VERSIONNUMBERCODE) < 0.0005);
112 if (!isConstructed) {
113 ostringstream errCode;
114 errCode << fixed << setprecision(3) <<
": in code " << VERSIONNUMBERCODE
115 <<
" but in XML " << versionNumberXML;
116 info.errorMsg(
"Abort from Pythia::Pythia: unmatched version numbers",
122 particleData.initPtr( &info, &settings, &rndm, couplingsPtr);
123 string dataFile = xmlPath +
"ParticleData.xml";
124 isConstructed = particleData.init( dataFile);
125 if (!isConstructed) {
126 info.errorMsg(
"Abort from Pythia::Pythia: particle data unavailable");
131 if (printBanner) banner();
146 if (useNewPdfHard && pdfHardAPtr != pdfAPtr)
delete pdfHardAPtr;
147 if (useNewPdfHard && pdfHardBPtr != pdfBPtr)
delete pdfHardBPtr;
148 if (useNewPdfA)
delete pdfAPtr;
149 if (useNewPdfB)
delete pdfBPtr;
150 if (useNewPdfPomA)
delete pdfPomAPtr;
151 if (useNewPdfPomB)
delete pdfPomBPtr;
154 if (useNewLHA)
delete lhaUpPtr;
157 if (hasOwnMergingHooks)
delete mergingHooksPtr;
160 if (useNewBeamShape)
delete beamShapePtr;
163 if (useNewTimes)
delete timesPtr;
164 if (useNewSpace)
delete spacePtr;
172 bool Pythia::readString(
string line,
bool warn) {
175 if (!isConstructed)
return false;
178 if (line.find_first_not_of(
" \n\t\v\b\r\f\a") == string::npos)
return true;
181 int firstChar = line.find_first_not_of(
" \n\t\v\b\r\f\a");
182 if (!isalnum(line[firstChar]))
return true;
185 if (isdigit(line[firstChar]))
186 return particleData.readString(line, warn);
189 return settings.readString(line, warn);
197 bool Pythia::readFile(
string fileName,
bool warn,
int subrun) {
200 if (!isConstructed)
return false;
203 const char* cstring = fileName.c_str();
204 ifstream is(cstring);
206 info.errorMsg(
"Error in Pythia::readFile: did not find file", fileName);
211 return readFile( is, warn, subrun);
220 bool Pythia::readFile(istream& is,
bool warn,
int subrun) {
223 if (!isConstructed)
return false;
227 bool isCommented =
false;
228 bool accepted =
true;
229 int subrunNow = SUBRUNDEFAULT;
230 while ( getline(is, line) ) {
233 int commentLine = readCommented( line);
234 if (commentLine == +1) isCommented =
true;
235 else if (commentLine == -1) isCommented =
false;
236 else if (isCommented) ;
240 int subrunLine = readSubrun( line, warn);
241 if (subrunLine >= 0) subrunNow = subrunLine;
244 if ( (subrunNow == subrun || subrunNow == SUBRUNDEFAULT)
245 && !readString( line, warn) ) accepted =
false;
258 bool Pythia::setPDFPtr( PDF* pdfAPtrIn, PDF* pdfBPtrIn, PDF* pdfHardAPtrIn,
259 PDF* pdfHardBPtrIn, PDF* pdfPomAPtrIn, PDF* pdfPomBPtrIn) {
262 if (useNewPdfHard && pdfHardAPtr != pdfAPtr)
delete pdfHardAPtr;
263 if (useNewPdfHard && pdfHardBPtr != pdfBPtr)
delete pdfHardBPtr;
264 if (useNewPdfA)
delete pdfAPtr;
265 if (useNewPdfB)
delete pdfBPtr;
266 if (useNewPdfPomA)
delete pdfPomAPtr;
267 if (useNewPdfPomB)
delete pdfPomBPtr;
272 useNewPdfHard =
false;
273 useNewPdfPomA =
false;
274 useNewPdfPomB =
false;
283 if (pdfAPtrIn == 0 && pdfBPtrIn == 0)
return true;
286 if (pdfAPtrIn == pdfBPtrIn)
return false;
293 pdfHardAPtr = pdfAPtrIn;
294 pdfHardBPtr = pdfBPtrIn;
297 if (pdfHardAPtrIn != 0 && pdfHardBPtrIn != 0) {
298 if (pdfHardAPtrIn == pdfHardBPtrIn)
return false;
299 pdfHardAPtr = pdfHardAPtrIn;
300 pdfHardBPtr = pdfHardBPtrIn;
304 if (pdfPomAPtrIn != 0 && pdfPomBPtrIn != 0) {
305 if (pdfPomAPtrIn == pdfPomBPtrIn)
return false;
306 pdfPomAPtr = pdfPomAPtrIn;
307 pdfPomBPtr = pdfPomBPtrIn;
318 bool Pythia::init() {
322 if (!isConstructed) {
323 info.errorMsg(
"Abort from Pythia::init: constructor "
324 "initialization failed");
329 doProcessLevel = settings.flag(
"ProcessLevel:all");
332 if (settings.readingFailed()) {
333 info.errorMsg(
"Abort from Pythia::init: some user settings "
334 "did not make sense");
337 if (particleData.readingFailed()) {
338 info.errorMsg(
"Abort from Pythia::init: some user particle data "
339 "did not make sense");
345 frameType = mode(
"Beams:frameType");
348 if (frameType < 4 ) {
350 boostType = frameType;
351 idA = mode(
"Beams:idA");
352 idB = mode(
"Beams:idB");
353 eCM = parm(
"Beams:eCM");
354 eA = parm(
"Beams:eA");
355 eB = parm(
"Beams:eB");
356 pxA = parm(
"Beams:pxA");
357 pyA = parm(
"Beams:pyA");
358 pzA = parm(
"Beams:pzA");
359 pxB = parm(
"Beams:pxB");
360 pyB = parm(
"Beams:pyB");
361 pzB = parm(
"Beams:pzB");
367 string lhef = word(
"Beams:LHEF");
368 string lhefHeader = word(
"Beams:LHEFheader");
369 bool readHeaders = flag(
"Beams:readLHEFheaders");
370 bool skipInit = flag(
"Beams:newLHEFsameInit");
371 int nSkipAtInit = mode(
"Beams:nSkipLHEFatInit");
374 if (frameType == 4) {
375 const char* cstring1 = lhef.c_str();
376 if (useNewLHA && skipInit) lhaUpPtr->newEventFile(cstring1);
378 if (useNewLHA)
delete lhaUpPtr;
380 const char* cstring2 = (lhefHeader ==
"void") ?
381 NULL : lhefHeader.c_str();
382 lhaUpPtr =
new LHAupLHEF(cstring1, cstring2, readHeaders);
387 if (!lhaUpPtr->fileFound()) {
388 info.errorMsg(
"Abort from Pythia::init: "
389 "Les Houches Event File not found");
396 info.errorMsg(
"Abort from Pythia::init: "
397 "LHAup object not found");
402 if (!lhaUpPtr->fileFound()) {
403 info.errorMsg(
"Abort from Pythia::init: "
404 "LHAup initialisation error");
410 lhaUpPtr->setPtr( &info);
411 processLevel.setLHAPtr( lhaUpPtr);
418 if (nSkipAtInit > 0) lhaUpPtr->skipEvent(nSkipAtInit);
423 if (frameType == 4 || frameType == 5) {
426 bool doUserMerging = settings.flag(
"Merging:doUserMerging");
427 bool doMGMerging = settings.flag(
"Merging:doMGMerging");
428 bool doKTMerging = settings.flag(
"Merging:doKTMerging");
429 bool doPTLundMerging = settings.flag(
"Merging:doPTLundMerging");
430 bool doCutBasedMerging = settings.flag(
"Merging:doCutBasedMerging");
432 bool doUMEPSTree = settings.flag(
"Merging:doUMEPSTree");
433 bool doUMEPSSubt = settings.flag(
"Merging:doUMEPSSubt");
435 bool doNL3Tree = settings.flag(
"Merging:doNL3Tree");
436 bool doNL3Loop = settings.flag(
"Merging:doNL3Loop");
437 bool doNL3Subt = settings.flag(
"Merging:doNL3Subt");
439 bool doUNLOPSTree = settings.flag(
"Merging:doUNLOPSTree");
440 bool doUNLOPSLoop = settings.flag(
"Merging:doUNLOPSLoop");
441 bool doUNLOPSSubt = settings.flag(
"Merging:doUNLOPSSubt");
442 bool doUNLOPSSubtNLO = settings.flag(
"Merging:doUNLOPSSubtNLO");
443 bool doXSectionEst = settings.flag(
"Merging:doXSectionEstimate");
444 doMerging = doUserMerging || doMGMerging || doKTMerging
445 || doPTLundMerging || doCutBasedMerging || doUMEPSTree || doUMEPSSubt
446 || doNL3Tree || doNL3Loop || doNL3Subt || doUNLOPSTree
447 || doUNLOPSLoop || doUNLOPSSubt || doUNLOPSSubtNLO || doXSectionEst;
450 bool inputMergingHooks = (mergingHooksPtr != 0);
451 if (doMerging && !inputMergingHooks){
452 if (hasOwnMergingHooks && mergingHooksPtr)
delete mergingHooksPtr;
453 mergingHooksPtr =
new MergingHooks();
454 hasOwnMergingHooks =
true;
457 hasMergingHooks = (mergingHooksPtr != 0);
460 if (doMerging && !hasMergingHooks) {
461 info.errorMsg(
"Abort from Pythia::init: "
462 "no merging hooks object has been provided");
464 }
else if (doMerging) {
465 string lhefIn = (frameType == 4) ? lhef :
"";
466 mergingHooksPtr->setLHEInputFile( lhefIn);
470 info.setCounter(41,0);
474 if ( !lhaUpPtr->setInit()) {
475 info.errorMsg(
"Abort from Pythia::init: "
476 "Les Houches initialization failed");
481 idA = lhaUpPtr->idBeamA();
482 idB = lhaUpPtr->idBeamB();
483 int idRenameBeams = settings.mode(
"LesHouches:idRenameBeams");
484 if (abs(idA) == idRenameBeams) idA = 16;
485 if (abs(idB) == idRenameBeams) idB = -16;
486 if (idA == 0 || idB == 0) doProcessLevel =
false;
487 eA = lhaUpPtr->eBeamA();
488 eB = lhaUpPtr->eBeamB();
491 if (nSkipAtInit > 0) lhaUpPtr->skipEvent(nSkipAtInit);
495 hasUserHooks = (userHooksPtr != 0);
496 doVetoProcess =
false;
497 doVetoPartons =
false;
498 retryPartonLevel =
false;
500 userHooksPtr->initPtr( &info, &settings, &particleData, &rndm, &beamA,
501 &beamB, &beamPomA, &beamPomB, couplingsPtr, &partonSystems, &sigmaTot);
502 if (!userHooksPtr->initAfterBeams()) {
503 info.errorMsg(
"Abort from Pythia::init: could not initialise UserHooks");
506 doVetoProcess = userHooksPtr->canVetoProcessLevel();
507 doVetoPartons = userHooksPtr->canVetoPartonLevel();
508 retryPartonLevel = userHooksPtr->retryPartonLevel();
513 info.setTooLowPTmin(
false);
517 doPartonLevel = settings.flag(
"PartonLevel:all");
518 doHadronLevel = settings.flag(
"HadronLevel:all");
519 doDiffraction = settings.flag(
"SoftQCD:all")
520 || settings.flag(
"SoftQCD:singleDiffractive")
521 || settings.flag(
"SoftQCD:doubleDiffractive")
522 || settings.flag(
"SoftQCD:centralDiffractive")
523 || settings.flag(
"SoftQCD:inelastic");
524 doResDec = settings.flag(
"ProcessLevel:resonanceDecays");
525 doFSRinRes = doPartonLevel && settings.flag(
"PartonLevel:FSR")
526 && settings.flag(
"PartonLevel:FSRinResonances");
527 decayRHadrons = settings.flag(
"RHadrons:allowDecay");
528 doMomentumSpread = settings.flag(
"Beams:allowMomentumSpread");
529 doVertexSpread = settings.flag(
"Beams:allowVertexSpread");
530 abortIfVeto = settings.flag(
"Check:abortIfVeto");
531 checkEvent = settings.flag(
"Check:event");
532 checkHistory = settings.flag(
"Check:history");
533 nErrList = settings.mode(
"Check:nErrList");
534 epTolErr = settings.parm(
"Check:epTolErr");
535 epTolWarn = settings.parm(
"Check:epTolWarn");
538 if ( doMerging && (hasMergingHooks || hasOwnMergingHooks) )
539 mergingHooksPtr->init( settings, &info, &particleData, &partonSystems );
542 if ( settings.flag(
"Random:setSeed") )
543 rndm.init( settings.mode(
"Random:seed") );
549 couplingsPtr->init( settings, &rndm );
552 bool useSLHAcouplings =
false;
553 slhaInterface.setPtr( &info );
554 slhaInterface.init( settings, &rndm, couplingsPtr, &particleData,
556 if (useSLHAcouplings) couplingsPtr = slhaInterface.couplingsPtr;
559 particleData.initPtr( &info, &settings, &rndm, couplingsPtr);
560 if (hasUserHooks) userHooksPtr->initPtr( &info, &settings, &particleData,
561 &rndm, &beamA, &beamB, &beamPomA, &beamPomB, couplingsPtr,
562 &partonSystems, &sigmaTot);
565 int startColTag = settings.mode(
"Event:startColTag");
566 process.init(
"(hard process)", &particleData, startColTag);
567 event.init(
"(complete event)", &particleData, startColTag);
570 particleData.initWidths( resonancePtrs);
573 rHadrons.init( &info, settings, &particleData, &rndm);
576 if (timesDecPtr == 0 || timesPtr == 0) {
577 TimeShower* timesNow =
new TimeShower();
578 if (timesDecPtr == 0) timesDecPtr = timesNow;
579 if (timesPtr == 0) timesPtr = timesNow;
583 spacePtr =
new SpaceShower();
588 timesPtr->initPtr( &info, &settings, &particleData, &rndm, couplingsPtr,
589 &partonSystems, userHooksPtr, mergingHooksPtr);
590 timesDecPtr->initPtr( &info, &settings, &particleData, &rndm, couplingsPtr,
591 &partonSystems, userHooksPtr, mergingHooksPtr);
592 spacePtr->initPtr( &info, &settings, &particleData, &rndm, couplingsPtr,
593 &partonSystems, userHooksPtr, mergingHooksPtr);
594 timesDecPtr->init( 0, 0);
597 if (beamShapePtr == 0) {
598 beamShapePtr =
new BeamShape();
599 useNewBeamShape =
true;
601 beamShapePtr->init( settings, &rndm);
605 info.errorMsg(
"Abort from Pythia::init: "
606 "checkBeams initialization failed");
611 if (!doProcessLevel) boostType = 1;
615 if (!initKinematics()) {
616 info.errorMsg(
"Abort from Pythia::init: "
617 "kinematics initialization failed");
623 info.errorMsg(
"Abort from Pythia::init: PDF initialization failed");
628 StringFlav* flavSelPtr = hadronLevel.getStringFlavPtr();
629 beamA.init( idA, pzAcm, eA, mA, &info, settings, &particleData, &rndm,
630 pdfAPtr, pdfHardAPtr, isUnresolvedA, flavSelPtr);
631 beamB.init( idB, pzBcm, eB, mB, &info, settings, &particleData, &rndm,
632 pdfBPtr, pdfHardBPtr, isUnresolvedB, flavSelPtr);
635 if ( doDiffraction) {
636 beamPomA.init( 990, 0.5 * eCM, 0.5 * eCM, 0., &info, settings,
637 &particleData, &rndm, pdfPomAPtr, pdfPomAPtr,
false, flavSelPtr);
638 beamPomB.init( 990, -0.5 * eCM, 0.5 * eCM, 0., &info, settings,
639 &particleData, &rndm, pdfPomBPtr, pdfPomBPtr,
false, flavSelPtr);
644 if ( doProcessLevel && !processLevel.init( &info, settings, &particleData,
645 &rndm, &beamA, &beamB, couplingsPtr, &sigmaTot, doLHA, &slhaInterface,
646 userHooksPtr, sigmaPtrs) ) {
647 info.errorMsg(
"Abort from Pythia::init: "
648 "processLevel initialization failed");
653 if ( !doProcessLevel) processLevel.initDecays( &info, &particleData,
657 if ( doPartonLevel && doProcessLevel && !partonLevel.init( &info, settings,
658 &particleData, &rndm, &beamA, &beamB, &beamPomA, &beamPomB, couplingsPtr,
659 &partonSystems, &sigmaTot, timesDecPtr, timesPtr, spacePtr, &rHadrons,
660 userHooksPtr, mergingHooksPtr,
false) ) {
661 info.errorMsg(
"Abort from Pythia::init: "
662 "partonLevel initialization failed" );
667 if ( !doProcessLevel || !doPartonLevel) partonLevel.init( &info, settings,
668 &particleData, &rndm, 0, 0, 0, 0, couplingsPtr, &partonSystems, 0,
669 timesDecPtr, 0, 0, &rHadrons, 0, 0,
false);
672 if ( doMerging && !trialPartonLevel.init( &info, settings, &particleData,
673 &rndm, &beamA, &beamB, &beamPomA, &beamPomB, couplingsPtr,
674 &partonSystems, &sigmaTot, timesDecPtr, timesPtr, spacePtr, &rHadrons,
675 NULL, mergingHooksPtr,
true) ) {
676 info.errorMsg(
"Abort from Pythia::init: "
677 "trialPartonLevel initialization failed");
682 if (doMerging ) merging.init( &settings, &info, &particleData, &rndm,
683 &beamA, &beamB, mergingHooksPtr, &trialPartonLevel );
687 if ( !hadronLevel.init( &info, settings, &particleData, &rndm,
688 couplingsPtr, timesDecPtr, &rHadrons, decayHandlePtr,
689 handledParticles) ) {
690 info.errorMsg(
"Abort from Pythia::init: "
691 "hadronLevel initialization failed");
696 if ( settings.flag(
"Check:particleData") )
697 particleData.checkTable( settings.mode(
"Check:levelParticleData") );
700 bool showCS = settings.flag(
"Init:showChangedSettings");
701 bool showAS = settings.flag(
"Init:showAllSettings");
702 bool showCPD = settings.flag(
"Init:showChangedParticleData");
703 bool showCRD = settings.flag(
"Init:showChangedResonanceData");
704 bool showAPD = settings.flag(
"Init:showAllParticleData");
705 int show1PD = settings.mode(
"Init:showOneParticleData");
706 bool showPro = settings.flag(
"Init:showProcesses");
707 if (showCS) settings.listChanged();
708 if (showAS) settings.listAll();
709 if (show1PD > 0) particleData.list(show1PD);
710 if (showCPD) particleData.listChanged(showCRD);
711 if (showAPD) particleData.listAll();
714 nCount = settings.mode(
"Next:numberCount");
715 nShowLHA = settings.mode(
"Next:numberShowLHA");
716 nShowInfo = settings.mode(
"Next:numberShowInfo");
717 nShowProc = settings.mode(
"Next:numberShowProcess");
718 nShowEvt = settings.mode(
"Next:numberShowEvent");
719 showSaV = settings.flag(
"Next:showScaleAndVertex");
720 showMaD = settings.flag(
"Next:showMothersAndDaughters");
725 if (useNewLHA && showPro) lhaUpPtr->listInit();
734 bool Pythia::init(
int idAin,
int idBin,
double eCMin) {
737 settings.mode(
"Beams:idA", idAin);
738 settings.mode(
"Beams:idB", idBin);
739 settings.mode(
"Beams:frameType", 1);
740 settings.parm(
"Beams:eCM", eCMin);
751 bool Pythia::init(
int idAin,
int idBin,
double eAin,
double eBin) {
754 settings.mode(
"Beams:idA", idAin);
755 settings.mode(
"Beams:idB", idBin);
756 settings.mode(
"Beams:frameType", 2);
757 settings.parm(
"Beams:eA", eAin);
758 settings.parm(
"Beams:eB", eBin);
769 bool Pythia::init(
int idAin,
int idBin,
double pxAin,
double pyAin,
770 double pzAin,
double pxBin,
double pyBin,
double pzBin) {
773 settings.mode(
"Beams:idA", idAin);
774 settings.mode(
"Beams:idB", idBin);
775 settings.mode(
"Beams:frameType", 3);
776 settings.parm(
"Beams:pxA", pxAin);
777 settings.parm(
"Beams:pyA", pyAin);
778 settings.parm(
"Beams:pzA", pzAin);
779 settings.parm(
"Beams:pxB", pxBin);
780 settings.parm(
"Beams:pyB", pyBin);
781 settings.parm(
"Beams:pzB", pzBin);
792 bool Pythia::init(
string LesHouchesEventFile,
bool skipInit) {
795 settings.mode(
"Beams:frameType", 4);
796 settings.word(
"Beams:LHEF", LesHouchesEventFile);
797 settings.flag(
"Beams:newLHEFsameInit", skipInit);
808 bool Pythia::init( LHAup* lhaUpPtrIn) {
811 setLHAupPtr( lhaUpPtrIn);
812 settings.mode(
"Beams:frameType", 5);
813 string lhaFileName = lhaUpPtr->getFileName();
814 if (lhaFileName !=
"void") settings.word(
"Beams:LHEF", lhaFileName);
825 void Pythia::checkSettings() {
828 if ((settings.flag(
"PartonLevel:ISR") || settings.flag(
"PartonLevel:FSR"))
829 && settings.flag(
"MultipartonInteractions:allowDoubleRescatter")) {
830 info.errorMsg(
"Warning in Pythia::checkSettings: "
831 "double rescattering switched off since showering is on");
832 settings.flag(
"MultipartonInteractions:allowDoubleRescatter",
false);
841 bool Pythia::checkBeams() {
844 int idAabs = abs(idA);
845 int idBabs = abs(idB);
846 if (!doProcessLevel)
return true;
849 bool isLeptonA = (idAabs > 10 && idAabs < 17);
850 bool isLeptonB = (idBabs > 10 && idBabs < 17);
851 bool isUnresLep = !settings.flag(
"PDF:lepton");
852 isUnresolvedA = isLeptonA && (idAabs%2 == 0 || isUnresLep);
853 isUnresolvedB = isLeptonB && (idBabs%2 == 0 || isUnresLep);
856 if (isLeptonA && isLeptonB && isUnresolvedA == isUnresolvedB)
return true;
859 int PomFlux = settings.mode(
"Diffraction:PomFlux");
861 bool ispp = (idAabs == 2212 && idBabs == 2212);
862 bool ispbarpbar = (idA == -2212 && idB == -2212);
863 if (ispp && !ispbarpbar)
return true;
864 info.errorMsg(
"Error in Pythia::init: cannot handle this beam combination"
865 " with PomFlux == 5");
870 bool isHadronA = (idAabs == 2212) || (idAabs == 2112) || (idA == 111)
871 || (idAabs == 211) || (idA == 990);
872 bool isHadronB = (idBabs == 2212) || (idBabs == 2112) || (idB == 111)
873 || (idBabs == 211) || (idB == 990);
874 if (isHadronA && isHadronB)
return true;
877 info.errorMsg(
"Error in Pythia::init: cannot handle this beam combination");
886 bool Pythia::initKinematics() {
889 mA = particleData.m0(idA);
890 mB = particleData.m0(idB);
895 if (boostType == 2) {
898 pzA = sqrt(eA*eA - mA*mA);
899 pzB = -sqrt(eB*eB - mB*mB);
900 pAinit = Vec4( 0., 0., pzA, eA);
901 pBinit = Vec4( 0., 0., pzB, eB);
902 eCM = sqrt( pow2(eA + eB) - pow2(pzA + pzB) );
905 betaZ = (pzA + pzB) / (eA + eB);
906 gammaZ = (eA + eB) / eCM;
907 if (abs(betaZ) < 1e-10) boostType = 1;
911 else if (boostType == 3) {
912 eA = sqrt( pxA*pxA + pyA*pyA + pzA*pzA + mA*mA);
913 eB = sqrt( pxB*pxB + pyB*pyB + pzB*pzB + mB*mB);
914 pAinit = Vec4( pxA, pyA, pzA, eA);
915 pBinit = Vec4( pxB, pyB, pzB, eB);
916 eCM = (pAinit + pBinit).mCalc();
920 MfromCM.fromCMframe( pAinit, pBinit);
927 info.errorMsg(
"Error in Pythia::initKinematics: too low energy");
932 pzAcm = 0.5 * sqrtpos( (eCM + mA + mB) * (eCM - mA - mB)
933 * (eCM - mA + mB) * (eCM + mA - mB) ) / eCM;
935 eA = sqrt(mA*mA + pzAcm*pzAcm);
936 eB = sqrt(mB*mB + pzBcm*pzBcm);
939 if (boostType != 2 && boostType != 3) {
940 pAinit = Vec4( 0., 0., pzAcm, eA);
941 pBinit = Vec4( 0., 0., pzBcm, eB);
945 info.setBeamA( idA, pzAcm, eA, mA);
946 info.setBeamB( idB, pzBcm, eB, mB);
950 if (doMomentumSpread) boostType = 3;
961 bool Pythia::initPDFs() {
965 if (pdfHardAPtr != pdfAPtr) {
969 if (pdfHardBPtr != pdfBPtr) {
973 useNewPdfHard =
false;
987 useNewPdfPomA =
false;
992 useNewPdfPomB =
false;
998 pdfAPtr = getPDFPtr(idA);
999 if (pdfAPtr == 0 || !pdfAPtr->isSetup()) {
1000 info.errorMsg(
"Error in Pythia::init: "
1001 "could not set up PDF for beam A");
1004 pdfHardAPtr = pdfAPtr;
1008 pdfBPtr = getPDFPtr(idB);
1009 if (pdfBPtr == 0 || !pdfBPtr->isSetup()) {
1010 info.errorMsg(
"Error in Pythia::init: "
1011 "could not set up PDF for beam B");
1014 pdfHardBPtr = pdfBPtr;
1019 if (settings.flag(
"PDF:useHard") && useNewPdfA && useNewPdfB) {
1020 pdfHardAPtr = getPDFPtr(idA, 2);
1021 if (!pdfHardAPtr->isSetup())
return false;
1022 pdfHardBPtr = getPDFPtr(idB, 2);
1023 if (!pdfHardBPtr->isSetup())
return false;
1024 useNewPdfHard =
true;
1028 if ( doDiffraction) {
1029 if (pdfPomAPtr == 0) {
1030 pdfPomAPtr = getPDFPtr(990);
1031 useNewPdfPomA =
true;
1033 if (pdfPomBPtr == 0) {
1034 pdfPomBPtr = getPDFPtr(990);
1035 useNewPdfPomB =
true;
1048 bool Pythia::next() {
1051 if (!isConstructed)
return false;
1054 int nPrevious = info.getCounter(3);
1055 if (nCount > 0 && nPrevious > 0 && nPrevious%nCount == 0)
1056 cout <<
"\n Pythia::next(): " << nPrevious
1057 <<
" events have been generated " << endl;
1061 for (
int i = 10; i < 13; ++i) info.setCounter(i);
1064 if (!doProcessLevel) {
1067 if (doLHA && !processLevel.nextLHAdec( event)) {
1068 if (info.atEndOfFile()) info.errorMsg(
"Abort from Pythia::next:"
1069 " reached end of Les Houches Events File");
1078 for (
int i = 1; i <
event.size(); ++i)
1079 if (event[i].isFinal()) pSum += event[i].p();
1081 event[0].m( pSum.mCalc() );
1084 bool status = forceHadronLevel();
1085 if (status) info.addCounter(4);
1086 if (status && nPrevious < nShowEvt)
event.list(showSaV, showMaD);
1094 partonSystems.clear();
1102 info.errorMsg(
"Abort from Pythia::next: "
1103 "not properly initialized so cannot generate events");
1108 if (doMomentumSpread || doVertexSpread) beamShapePtr->pick();
1111 if (doMomentumSpread) nextKinematics();
1116 info.addCounter(10);
1117 bool hasVetoed =
false;
1123 if ( !processLevel.next( process) ) {
1124 if (doLHA && info.atEndOfFile()) info.errorMsg(
"Abort from "
1125 "Pythia::next: reached end of Les Houches Events File");
1126 else info.errorMsg(
"Abort from Pythia::next: "
1127 "processLevel failed; giving up");
1131 info.addCounter(11);
1134 if (doVetoProcess) {
1135 hasVetoed = userHooksPtr->doVetoProcessLevel( process);
1137 if (abortIfVeto)
return false;
1144 int veto = merging.mergeProcess( process );
1148 if (abortIfVeto)
return false;
1151 }
else if ( veto == 0 )
break;
1155 if (veto == 2 && doResDec) processLevel.nextDecays( process);
1159 if (!doPartonLevel) {
1160 boostAndVertex(
true,
true);
1161 processLevel.accumulate();
1163 if (doLHA && nPrevious < nShowLHA) lhaUpPtr->listEvent();
1164 if (nPrevious < nShowInfo) info.list();
1165 if (nPrevious < nShowProc) process.list(showSaV, showMaD);
1170 Event processSave = process;
1171 int sizeMPI = info.sizeMPIarrays();
1172 info.addCounter(12);
1173 for (
int i = 14; i < 19; ++i) info.setCounter(i);
1176 bool physical =
true;
1177 for (
int iTry = 0; iTry < NTRY; ++iTry) {
1179 info.addCounter(14);
1184 if (iTry > 0) process = processSave;
1185 if (iTry > 0) info.resizeMPIarrays( sizeMPI);
1193 partonSystems.clear();
1196 if ( !partonLevel.next( process, event) ) {
1199 hasVetoed = partonLevel.hasVetoed();
1201 if (retryPartonLevel) {
1205 if (abortIfVeto)
return false;
1209 info.errorMsg(
"Error in Pythia::next: "
1210 "partonLevel failed; try again");
1214 info.addCounter(15);
1217 if (doVetoPartons) {
1218 hasVetoed = userHooksPtr->doVetoPartonLevel( event);
1220 if (abortIfVeto)
return false;
1226 boostAndVertex(
true,
true);
1229 if (!doHadronLevel) {
1230 processLevel.accumulate();
1231 partonLevel.accumulate();
1233 if (checkEvent && !check()) {
1234 info.errorMsg(
"Abort from Pythia::next: "
1235 "check of event revealed problems");
1239 if (doLHA && nPrevious < nShowLHA) lhaUpPtr->listEvent();
1240 if (nPrevious < nShowInfo) info.list();
1241 if (nPrevious < nShowProc) process.list(showSaV, showMaD);
1242 if (nPrevious < nShowEvt)
event.list(showSaV, showMaD);
1247 info.addCounter(16);
1248 if ( !hadronLevel.next( event) ) {
1249 info.errorMsg(
"Error in Pythia::next: "
1250 "hadronLevel failed; try again");
1256 if (decayRHadrons && rHadrons.exist() && !doRHadronDecays()) {
1257 info.errorMsg(
"Error in Pythia::next: "
1258 "decayRHadrons failed; try again");
1262 info.addCounter(17);
1265 if (checkEvent && !check()) {
1266 info.errorMsg(
"Error in Pythia::next: "
1267 "check of event revealed problems");
1273 info.addCounter(18);
1279 if (abortIfVeto)
return false;
1285 info.errorMsg(
"Abort from Pythia::next: "
1286 "parton+hadronLevel failed; giving up");
1291 processLevel.accumulate();
1292 partonLevel.accumulate();
1293 event.scale( process.scale() );
1296 info.addCounter(13);
1301 if (doLHA && nPrevious < nShowLHA) lhaUpPtr->listEvent();
1302 if (nPrevious < nShowInfo) info.list();
1303 if (nPrevious < nShowProc) process.list(showSaV,showMaD);
1304 if (nPrevious < nShowEvt)
event.list(showSaV, showMaD);
1317 bool Pythia::forceHadronLevel(
bool findJunctions) {
1321 info.errorMsg(
"Abort from Pythia::forceHadronLevel: "
1322 "not properly initialized so cannot generate events");
1328 if (findJunctions) {
1329 event.clearJunctions();
1330 for (
int i = 0; i <
event.size(); ++i)
1331 if (event[i].isFinal()
1332 && (
event[i].col() != 0 ||
event[i].acol() != 0)) {
1333 processLevel.findJunctions( event);
1339 Event spareEvent = event;
1342 bool physical =
true;
1343 for (
int iTry = 0; iTry < NTRY; ++ iTry) {
1349 processLevel.nextDecays( process);
1352 if (process.size() >
event.size()) {
1354 partonLevel.setupShowerSys( process, event);
1355 partonLevel.resonanceShowers( process, event,
false);
1356 }
else event = process;
1361 if (hadronLevel.next( event))
break;
1364 info.errorMsg(
"Error in Pythia::forceHadronLevel: "
1365 "hadronLevel failed; try again");
1372 info.errorMsg(
"Abort from Pythia::forceHadronLevel: "
1373 "hadronLevel failed; giving up");
1378 if (checkEvent && !check()) {
1379 info.errorMsg(
"Abort from Pythia::forceHadronLevel: "
1380 "check of event revealed problems");
1393 void Pythia::nextKinematics() {
1396 pAnow = pAinit + beamShapePtr->deltaPA();
1397 pAnow.e( sqrt(pAnow.pAbs2() + mA * mA) );
1398 pBnow = pBinit + beamShapePtr->deltaPB();
1399 pBnow.e( sqrt(pBnow.pAbs2() + mB * mB) );
1402 eCM = (pAnow + pBnow).mCalc();
1403 pzAcm = 0.5 * sqrtpos( (eCM + mA + mB) * (eCM - mA - mB)
1404 * (eCM - mA + mB) * (eCM + mA - mB) ) / eCM;
1406 eA = sqrt(mA*mA + pzAcm*pzAcm);
1407 eB = sqrt(mB*mB + pzBcm*pzBcm);
1410 info.setBeamA( idA, pzAcm, eA, mA);
1411 info.setBeamB( idB, pzBcm, eB, mB);
1413 beamA.newPzE( pzAcm, eA);
1414 beamB.newPzE( pzBcm, eB);
1418 MfromCM.fromCMframe( pAnow, pBnow);
1428 void Pythia::boostAndVertex(
bool toLab,
bool setVertex) {
1432 if (boostType == 2) process.bst(0., 0., betaZ, gammaZ);
1433 else if (boostType == 3) process.rotbst(MfromCM);
1436 if (event.size() > 0) {
1437 if (boostType == 2)
event.bst(0., 0., betaZ, gammaZ);
1438 else if (boostType == 3)
event.rotbst(MfromCM);
1443 if (boostType == 2) process.bst(0., 0., -betaZ, gammaZ);
1444 else if (boostType == 3) process.rotbst(MtoCM);
1447 if (event.size() > 0) {
1448 if (boostType == 2)
event.bst(0., 0., -betaZ, gammaZ);
1449 else if (boostType == 3)
event.rotbst(MtoCM);
1454 if (setVertex && doVertexSpread) {
1455 Vec4
vertex = beamShapePtr->vertex();
1456 for (
int i = 0; i < process.size(); ++i) process[i].vProd( vertex);
1457 for (
int i = 0; i <
event.size(); ++i) event[i].vProd( vertex);
1466 bool Pythia::doRHadronDecays( ) {
1469 if ( !rHadrons.exist() )
return true;
1472 if ( !rHadrons.decay( event) )
return false;
1475 if ( !partonLevel.resonanceShowers( process, event,
false) )
return false;
1478 if ( !hadronLevel.next( event) )
return false;
1489 void Pythia::stat() {
1492 bool showPrL = settings.flag(
"Stat:showProcessLevel");
1493 bool showPaL = settings.flag(
"Stat:showPartonLevel");
1494 bool showErr = settings.flag(
"Stat:showErrors");
1495 bool reset = settings.flag(
"Stat:reset");
1498 if (doProcessLevel) {
1499 if (showPrL) processLevel.statistics(
false);
1500 if (reset) processLevel.resetStatistics();
1504 if (showPaL) partonLevel.statistics(
false);
1505 if (reset) partonLevel.resetStatistics();
1508 if (doMerging) merging.statistics();
1511 if (showErr) info.errorStatistics();
1512 if (reset) info.errorReset();
1520 void Pythia::statistics(
bool all,
bool reset) {
1523 if (doProcessLevel) processLevel.statistics(reset);
1526 if (all) partonLevel.statistics(reset);
1529 if (doMerging) merging.statistics();
1532 info.errorStatistics();
1533 if (reset) info.errorReset();
1541 void Pythia::banner(ostream& os) {
1544 double versionNumber = settings.parm(
"Pythia:versionNumber");
1545 int versionDate = settings.mode(
"Pythia:versionDate");
1546 string month[12] = {
"Jan",
"Feb",
"Mar",
"Apr",
"May",
"Jun",
1547 "Jul",
"Aug",
"Sep",
"Oct",
"Nov",
"Dec"};
1552 strftime(dateNow,12,
"%d %b %Y",localtime(&t));
1554 strftime(timeNow,9,
"%H:%M:%S",localtime(&t));
1557 <<
" *-------------------------------------------"
1558 <<
"-----------------------------------------* \n"
1561 <<
" | *----------------------------------------"
1562 <<
"--------------------------------------* | \n"
1567 <<
" | | PPP Y Y TTTTT H H III A "
1568 <<
" Welcome to the Lund Monte Carlo! | | \n"
1569 <<
" | | P P Y Y T H H I A A "
1570 <<
" This is PYTHIA version " << fixed << setprecision(3)
1571 << setw(5) << versionNumber <<
" | | \n"
1572 <<
" | | PPP Y T HHHHH I AAAAA"
1573 <<
" Last date of change: " << setw(2) << versionDate%100
1574 <<
" " << month[ (versionDate/100)%100 - 1 ]
1575 <<
" " << setw(4) << versionDate/10000 <<
" | | \n"
1576 <<
" | | P Y T H H I A A"
1578 <<
" | | P Y T H H III A A"
1579 <<
" Now is " << dateNow <<
" at " << timeNow <<
" | | \n"
1582 <<
" | | Torbjorn Sjostrand; Department of As"
1583 <<
"tronomy and Theoretical Physics, | | \n"
1584 <<
" | | Lund University, Solvegatan 14A, S"
1585 <<
"E-223 62 Lund, Sweden; | | \n"
1586 <<
" | | phone: + 46 - 46 - 222 48 16; e-ma"
1587 <<
"il: torbjorn@thep.lu.se | | \n"
1588 <<
" | | Jesper Roy Christiansen; Department "
1589 <<
"of Astronomy and Theoretical Physics, | | \n"
1590 <<
" | | Lund University, Solvegatan 14A, S"
1591 <<
"E-223 62 Lund, Sweden; | | \n"
1592 <<
" | | e-mail: Jesper.Roy.Christiansen@th"
1593 <<
"ep.lu.se | | \n"
1594 <<
" | | Nishita Desai; Institut fuer Theoret"
1595 <<
"ische Physik, | | \n"
1596 <<
" | | Universitaet Heidelberg, Philosophe"
1597 <<
"nweg 16, D-69120 Heidelberg, Germany; | | \n"
1598 <<
" | | phone: +49 - 6221 54 9424; e-mail:"
1599 <<
" Nishita.Desai@cern.ch | | \n"
1600 <<
" | | Philip Ilten; Massachusetts Institut"
1601 <<
"e of Technology, | | \n"
1602 <<
" | | stationed at CERN, CH-1211 Geneva "
1603 <<
"23, Switzerland; | | \n"
1604 <<
" | | e-mail: philten@cern.ch "
1606 <<
" | | Stephen Mrenna; Computing Division, "
1607 <<
"Simulations Group, | | \n"
1608 <<
" | | Fermi National Accelerator Laborat"
1609 <<
"ory, MS 234, Batavia, IL 60510, USA; | | \n"
1610 <<
" | | phone: + 1 - 630 - 840 - 2556; e-m"
1611 <<
"ail: mrenna@fnal.gov | | \n"
1612 <<
" | | Stefan Prestel; Theory Group, DESY, "
1614 <<
" | | Notkestrasse 85, D-22607 Hamburg, "
1615 <<
"Germany; | | \n"
1616 <<
" | | phone: + 49 - 40 - 8998-4250; e-ma"
1617 <<
"il: stefan.prestel@thep.lu.se | | \n"
1618 <<
" | | Peter Skands; Theoretical Physics, C"
1619 <<
"ERN, CH-1211 Geneva 23, Switzerland; | | \n"
1620 <<
" | | phone: + 41 - 22 - 767 2447; e-mai"
1621 <<
"l: peter.skands@cern.ch | | \n"
1622 <<
" | | Stefan Ask; former author; e-mail: as"
1623 <<
"k.stefan@gmail.com | | \n"
1624 <<
" | | Richard Corke; former author; e-mail:"
1625 <<
" r.corke@errno.net | | \n"
1628 <<
" | | The main program reference is the 'Br"
1629 <<
"ief Introduction to PYTHIA 8.1', | | \n"
1630 <<
" | | T. Sjostrand, S. Mrenna and P. Skands"
1631 <<
", Comput. Phys. Comm. 178 (2008) 85 | | \n"
1632 <<
" | | [arXiv:0710.3820] "
1636 <<
" | | The main physics reference is the 'PY"
1637 <<
"THIA 6.4 Physics and Manual', | | \n"
1638 <<
" | | T. Sjostrand, S. Mrenna and P. Skands"
1639 <<
", JHEP05 (2006) 026 [hep-ph/0603175]. | | \n"
1642 <<
" | | An archive of program versions and do"
1643 <<
"cumentation is found on the web: | | \n"
1644 <<
" | | http://www.thep.lu.se/~torbjorn/Pythi"
1648 <<
" | | This program is released under the GN"
1649 <<
"U General Public Licence version 2. | | \n"
1650 <<
" | | Please respect the MCnet Guidelines f"
1651 <<
"or Event Generator Authors and Users. | | \n"
1654 <<
" | | Disclaimer: this program comes withou"
1655 <<
"t any guarantees. | | \n"
1656 <<
" | | Beware of errors and use common sense"
1657 <<
" when interpreting results. | | \n"
1660 <<
" | | Copyright (C) 2014 Torbjorn Sjostrand"
1666 <<
" | *----------------------------------------"
1667 <<
"--------------------------------------* | \n"
1670 <<
" *-------------------------------------------"
1671 <<
"-----------------------------------------* \n" << endl;
1679 int Pythia::readSubrun(
string line,
bool warn, ostream& os) {
1682 int subrunLine = SUBRUNDEFAULT;
1683 if (line.find_first_not_of(
" \n\t\v\b\r\f\a") == string::npos)
1687 string lineNow = line;
1688 int firstChar = lineNow.find_first_not_of(
" \n\t\v\b\r\f\a");
1689 if (!isalpha(lineNow[firstChar]))
return subrunLine;
1692 while (lineNow.find(
"=") != string::npos) {
1693 int firstEqual = lineNow.find_first_of(
"=");
1694 lineNow.replace(firstEqual, 1,
" ");
1698 istringstream splitLine(lineNow);
1703 while (name.find(
"::") != string::npos) {
1704 int firstColonColon = name.find_first_of(
"::");
1705 name.replace(firstColonColon, 2,
":");
1710 for (
int i = 0; i < int(name.length()); ++i) name[i] = tolower(name[i]);
1713 if (name !=
"main:subrun")
return subrunLine;
1716 splitLine >> subrunLine;
1718 if (warn) os <<
"\n PYTHIA Warning: Main:subrun number not"
1719 <<
" recognized; skip:\n " << line << endl;
1720 subrunLine = SUBRUNDEFAULT;
1731 int Pythia::readCommented(
string line) {
1734 if (line.find_first_not_of(
" \n\t\v\b\r\f\a") == string::npos)
return 0;
1735 int firstChar = line.find_first_not_of(
" \n\t\v\b\r\f\a");
1736 if (
int(line.size()) < firstChar + 2)
return 0;
1739 if (line.substr(firstChar, 2) ==
"/*")
return +1;
1740 if (line.substr(firstChar, 2) ==
"*/")
return -1;
1752 bool Pythia::check(ostream& os) {
1755 bool physical =
true;
1756 bool listVertices =
false;
1757 bool listHistory =
false;
1758 bool listSystems =
false;
1759 bool listBeams =
false;
1764 iErrNanVtx.resize(0);
1766 double chargeSum = 0.;
1769 if (doProcessLevel) {
1770 pSum = - (
event[1].p() +
event[2].p());
1771 chargeSum = - (
event[1].charge() +
event[2].charge());
1774 }
else if (process.size() > 0) {
1775 pSum = - process[0].p();
1776 for (
int i = 0; i < process.size(); ++i)
1777 if (process[i].isFinal()) chargeSum -= process[i].charge();
1781 pSum = -
event[0].p();
1782 for (
int i = 1; i <
event.size(); ++i)
1783 if (event[i].statusAbs() < 10 ||
event[i].statusAbs() == 23)
1784 chargeSum -= event[i].charge();
1786 double eLab = abs(pSum.e());
1789 for (
int i = 0; i <
event.size(); ++i) {
1792 int id =
event[i].id();
1793 if (
id == 0 || !particleData.isParticle(
id)) {
1794 ostringstream errCode;
1795 errCode <<
", i = " << i <<
", id = " << id;
1796 info.errorMsg(
"Error in Pythia::check: "
1797 "unknown particle code", errCode.str());
1799 iErrId.push_back(i);
1803 int colType =
event[i].colType();
1804 int col =
event[i].col();
1805 int acol =
event[i].acol();
1806 if ( (colType == 0 && (col > 0 || acol > 0))
1807 || (colType == 1 && (col <= 0 || acol > 0))
1808 || (colType == -1 && (col > 0 || acol <= 0))
1809 || (colType == 2 && (col <= 0 || acol <= 0)) ) {
1810 ostringstream errCode;
1811 errCode <<
", i = " << i <<
", id = " <<
id <<
" cols = " << col
1813 info.errorMsg(
"Error in Pythia::check: "
1814 "incorrect colours", errCode.str());
1816 iErrCol.push_back(i);
1821 if (abs(event[i].px()) >= 0. && abs(event[i].py()) >= 0.
1822 && abs(event[i].pz()) >= 0. && abs(event[i].e()) >= 0.
1823 && abs(event[i].m()) >= 0.) {
1824 if (abs(event[i].mCalc() - event[i].m()) > epTolErr * event[i].e()) {
1825 info.errorMsg(
"Error in Pythia::check: "
1826 "unmatched particle energy/momentum/mass");
1828 iErrEpm.push_back(i);
1831 info.errorMsg(
"Error in Pythia::check: "
1832 "not-a-number energy/momentum/mass");
1834 iErrNan.push_back(i);
1838 if (abs(event[i].xProd()) >= 0. && abs(event[i].yProd()) >= 0.
1839 && abs(event[i].zProd()) >= 0. && abs(event[i].tProd()) >= 0.
1840 && abs(event[i].tau()) >= 0.) ;
1842 info.errorMsg(
"Error in Pythia::check: "
1843 "not-a-number vertex/lifetime");
1845 listVertices =
true;
1846 iErrNanVtx.push_back(i);
1850 if (event[i].isFinal()) {
1851 pSum +=
event[i].p();
1852 chargeSum +=
event[i].charge();
1859 double epDev = abs(pSum.e()) + abs(pSum.px()) + abs(pSum.py())
1861 if (epDev > epTolErr * eLab) {
1862 info.errorMsg(
"Error in Pythia::check: energy-momentum not conserved");
1864 }
else if (epDev > epTolWarn * eLab) {
1865 info.errorMsg(
"Warning in Pythia::check: "
1866 "energy-momentum not quite conserved");
1868 if (abs(chargeSum) > 0.1) {
1869 info.errorMsg(
"Error in Pythia::check: charge not conserved");
1875 if (info.isResolved())
1876 for (
int iSys = 0; iSys < beamA.sizeInit(); ++iSys) {
1877 int eventANw = partonSystems.getInA(iSys);
1878 int eventBNw = partonSystems.getInB(iSys);
1879 int beamANw = beamA[iSys].iPos();
1880 int beamBNw = beamB[iSys].iPos();
1881 if (eventANw != beamANw || eventBNw != beamBNw) {
1882 info.errorMsg(
"Error in Pythia::check: "
1883 "event and beams records disagree");
1893 vector< pair<int,int> > noMotDau;
1897 bool hasBeams =
false;
1898 for (
int i = 0; i <
event.size(); ++i) {
1899 int status =
event[i].status();
1900 if (abs(status) == 12) hasBeams =
true;
1903 vector<int> mList =
event[i].motherList();
1904 vector<int> dList =
event[i].daughterList();
1905 if (mList.size() == 0 && abs(status) != 11 && abs(status) != 12)
1907 if (dList.size() == 0 && status < 0 && status != -11)
1911 for (
int j = 0; j < int(mList.size()); ++j) {
1912 if ( event[mList[j]].daughter1() <= i
1913 &&
event[mList[j]].daughter2() >= i )
continue;
1914 vector<int> dmList =
event[mList[j]].daughterList();
1915 bool foundMatch =
false;
1916 for (
int k = 0; k < int(dmList.size()); ++k)
1917 if (dmList[k] == i) {
1921 if (!hasBeams && mList.size() == 1 && mList[0] == 0) foundMatch =
true;
1923 bool oldPair =
false;
1924 for (
int k = 0; k < int(noMotDau.size()); ++k)
1925 if (noMotDau[k].first == mList[j] && noMotDau[k].second == i) {
1929 if (!oldPair) noMotDau.push_back( make_pair( mList[j], i) );
1934 for (
int j = 0; j < int(dList.size()); ++j) {
1935 if ( event[dList[j]].statusAbs() > 80
1936 &&
event[dList[j]].statusAbs() < 90
1937 &&
event[dList[j]].mother1() <= i
1938 &&
event[dList[j]].mother2() >= i)
continue;
1939 vector<int> mdList =
event[dList[j]].motherList();
1940 bool foundMatch =
false;
1941 for (
int k = 0; k < int(mdList.size()); ++k)
1942 if (mdList[k] == i) {
1947 bool oldPair =
false;
1948 for (
int k = 0; k < int(noMotDau.size()); ++k)
1949 if (noMotDau[k].first == i && noMotDau[k].second == dList[j]) {
1953 if (!oldPair) noMotDau.push_back( make_pair( i, dList[j]) );
1959 if (noMot.size() > 0 || noDau.size() > 0 || noMotDau.size() > 0) {
1960 info.errorMsg(
"Error in Pythia::check: "
1961 "mismatch in daughter and mother lists");
1968 if (physical)
return true;
1971 if (nErrEvent < nErrList) {
1972 os <<
"\n PYTHIA erroneous event info: \n";
1973 if (iErrId.size() > 0) {
1974 os <<
" unknown particle codes in lines ";
1975 for (
int i = 0; i < int(iErrId.size()); ++i)
1976 os << iErrId[i] <<
" ";
1979 if (iErrCol.size() > 0) {
1980 os <<
" incorrect colour assignments in lines ";
1981 for (
int i = 0; i < int(iErrCol.size()); ++i)
1982 os << iErrCol[i] <<
" ";
1985 if (iErrEpm.size() > 0) {
1986 os <<
" mismatch between energy/momentum/mass in lines ";
1987 for (
int i = 0; i < int(iErrEpm.size()); ++i)
1988 os << iErrEpm[i] <<
" ";
1991 if (iErrNan.size() > 0) {
1992 os <<
" not-a-number energy/momentum/mass in lines ";
1993 for (
int i = 0; i < int(iErrNan.size()); ++i)
1994 os << iErrNan[i] <<
" ";
1997 if (iErrNanVtx.size() > 0) {
1998 os <<
" not-a-number vertex/lifetime in lines ";
1999 for (
int i = 0; i < int(iErrNanVtx.size()); ++i)
2000 os << iErrNanVtx[i] <<
" ";
2003 if (epDev > epTolErr * eLab) os << scientific << setprecision(3)
2004 <<
" total energy-momentum non-conservation = " << epDev <<
"\n";
2005 if (abs(chargeSum) > 0.1) os << fixed << setprecision(2)
2006 <<
" total charge non-conservation = " << chargeSum <<
"\n";
2007 if (noMot.size() > 0) {
2008 os <<
" missing mothers for particles ";
2009 for (
int i = 0; i < int(noMot.size()); ++i) os << noMot[i] <<
" ";
2012 if (noDau.size() > 0) {
2013 os <<
" missing daughters for particles ";
2014 for (
int i = 0; i < int(noDau.size()); ++i) os << noDau[i] <<
" ";
2017 if (noMotDau.size() > 0) {
2018 os <<
" inconsistent history for (mother,daughter) pairs ";
2019 for (
int i = 0; i < int(noMotDau.size()); ++i)
2020 os <<
"(" << noMotDau[i].first <<
"," << noMotDau[i].second <<
") ";
2026 event.list(listVertices, listHistory);
2027 if (listSystems) partonSystems.list();
2028 if (listBeams) beamA.list();
2029 if (listBeams) beamB.list();
2042 PDF* Pythia::getPDFPtr(
int idIn,
int sequence) {
2045 PDF* tempPDFPtr = 0;
2048 if (idIn == 990 && settings.mode(
"PDF:PomSet") == 2) idIn = 111;
2051 if ((abs(idIn) == 2212 || abs(idIn) == 2112) && sequence == 1) {
2052 int pSet = settings.mode(
"PDF:pSet");
2053 bool useLHAPDF = settings.flag(
"PDF:useLHAPDF");
2057 if (pSet == 1) tempPDFPtr =
new GRV94L(idIn);
2058 else if (pSet == 2) tempPDFPtr =
new CTEQ5L(idIn);
2060 tempPDFPtr =
new MSTWpdf(idIn, pSet - 2, xmlPath, &info);
2061 else if (pSet <= 12)
2062 tempPDFPtr =
new CTEQ6pdf(idIn, pSet - 6, xmlPath, &info);
2063 else tempPDFPtr =
new NNPDF(idIn, pSet - 12, xmlPath, &info);
2068 string LHAPDFset = settings.word(
"PDF:LHAPDFset");
2069 int LHAPDFmember = settings.mode(
"PDF:LHAPDFmember");
2070 int nSet = LHAPDF::findNSet(LHAPDFset, LHAPDFmember);
2071 if (nSet == -1) nSet = LHAPDF::freeNSet();
2072 tempPDFPtr =
new LHAPDF(idIn, LHAPDFset, LHAPDFmember, nSet, &info);
2075 tempPDFPtr->setExtrapolate( settings.flag(
"PDF:extrapolateLHAPDF") );
2080 else if (abs(idIn) == 2212 || abs(idIn) == 2112) {
2081 int pSet = settings.mode(
"PDF:pHardSet");
2082 bool useLHAPDF = settings.flag(
"PDF:useHardLHAPDF");
2086 if (pSet == 1) tempPDFPtr =
new GRV94L(idIn);
2087 else if (pSet == 2) tempPDFPtr =
new CTEQ5L(idIn);
2089 tempPDFPtr =
new MSTWpdf(idIn, pSet - 2, xmlPath, &info);
2090 else if (pSet <= 12)
2091 tempPDFPtr =
new CTEQ6pdf(idIn, pSet - 6, xmlPath, &info);
2092 else tempPDFPtr =
new NNPDF(idIn, pSet - 12, xmlPath, &info);
2097 string LHAPDFset = settings.word(
"PDF:hardLHAPDFset");
2098 int LHAPDFmember = settings.mode(
"PDF:hardLHAPDFmember");
2099 int nSet = LHAPDF::findNSet(LHAPDFset, LHAPDFmember);
2100 if (nSet == -1) nSet = LHAPDF::freeNSet();
2101 tempPDFPtr =
new LHAPDF(idIn, LHAPDFset, LHAPDFmember, nSet, &info);
2104 tempPDFPtr->setExtrapolate( settings.flag(
"PDF:extrapolateLHAPDF") );
2109 else if (abs(idIn) == 211 || idIn == 111) {
2110 bool useLHAPDF = settings.flag(
"PDF:piUseLHAPDF");
2114 tempPDFPtr =
new GRVpiL(idIn);
2119 string LHAPDFset = settings.word(
"PDF:piLHAPDFset");
2120 int LHAPDFmember = settings.mode(
"PDF:piLHAPDFmember");
2121 tempPDFPtr =
new LHAPDF(idIn, LHAPDFset, LHAPDFmember, 1, &info);
2124 tempPDFPtr->setExtrapolate( settings.flag(
"PDF:extrapolateLHAPDF") );
2129 else if (idIn == 990) {
2130 int pomSet = settings.mode(
"PDF:PomSet");
2131 double rescale = settings.parm(
"PDF:PomRescale");
2135 double gluonA = settings.parm(
"PDF:PomGluonA");
2136 double gluonB = settings.parm(
"PDF:PomGluonB");
2137 double quarkA = settings.parm(
"PDF:PomQuarkA");
2138 double quarkB = settings.parm(
"PDF:PomQuarkB");
2139 double quarkFrac = settings.parm(
"PDF:PomQuarkFrac");
2140 double strangeSupp = settings.parm(
"PDF:PomStrangeSupp");
2141 tempPDFPtr =
new PomFix( 990, gluonA, gluonB, quarkA, quarkB,
2142 quarkFrac, strangeSupp);
2146 else if (pomSet == 3 || pomSet == 4)
2147 tempPDFPtr =
new PomH1FitAB( 990, pomSet - 2, rescale, xmlPath, &info);
2148 else if (pomSet == 5)
2149 tempPDFPtr =
new PomH1Jets( 990, rescale, xmlPath, &info);
2150 else if (pomSet == 6)
2151 tempPDFPtr =
new PomH1FitAB( 990, 3, rescale, xmlPath, &info);
2155 else if (abs(idIn) > 10 && abs(idIn) < 17) {
2156 if (abs(idIn)%2 == 0) tempPDFPtr =
new NeutrinoPoint(idIn);
2157 else if (settings.flag(
"PDF:lepton")) tempPDFPtr =
new Lepton(idIn);
2158 else tempPDFPtr =
new LeptonPoint(idIn);
2162 else tempPDFPtr = 0;