StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Pythia.cc
1 // Pythia.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2014 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the Pythia class.
7 
8 #include "Pythia8/Pythia.h"
9 
10 // Access time information.
11 #include <ctime>
12 
13 // Allow string and character manipulation.
14 #include <cctype>
15 
16 namespace Pythia8 {
17 
18 //==========================================================================
19 
20 // The Pythia class.
21 
22 //--------------------------------------------------------------------------
23 
24 // The current Pythia (sub)version number, to agree with XML version.
25 const double Pythia::VERSIONNUMBERCODE = 8.186;
26 
27 //--------------------------------------------------------------------------
28 
29 // Constants: could be changed here if desired, but normally should not.
30 // These are of technical nature, as described for each.
31 
32 // Maximum number of tries to produce parton level from given input.
33 const int Pythia::NTRY = 10;
34 
35 // Negative integer to denote that no subrun has been set.
36 const int Pythia::SUBRUNDEFAULT = -999;
37 
38 //--------------------------------------------------------------------------
39 
40 // Constructor.
41 
42 Pythia::Pythia(string xmlDir, bool printBanner) {
43 
44  // Initial values for pointers to PDF's.
45  useNewPdfA = false;
46  useNewPdfB = false;
47  useNewPdfHard = false;
48  useNewPdfPomA = false;
49  useNewPdfPomB = false;
50  pdfAPtr = 0;
51  pdfBPtr = 0;
52  pdfHardAPtr = 0;
53  pdfHardBPtr = 0;
54  pdfPomAPtr = 0;
55  pdfPomBPtr = 0;
56 
57  // Initial values for pointers to Les Houches Event objects.
58  doLHA = false;
59  useNewLHA = false;
60  lhaUpPtr = 0;
61 
62  //Initial value for couplings pointer
63  couplingsPtr = &couplings;
64 
65  // Initial value for pointer to external decay handler.
66  decayHandlePtr = 0;
67 
68  // Initial value for pointer to user hooks.
69  userHooksPtr = 0;
70 
71  // Initial value for pointer to merging hooks.
72  doMerging = false;
73  hasMergingHooks = false;
74  hasOwnMergingHooks = false;
75  mergingHooksPtr = 0;
76 
77  // Initial value for pointer to beam shape.
78  useNewBeamShape = false;
79  beamShapePtr = 0;
80 
81  // Initial values for pointers to timelike and spacelike showers.
82  useNewTimes = false;
83  useNewSpace = false;
84  timesDecPtr = 0;
85  timesPtr = 0;
86  spacePtr = 0;
87 
88  // Find path to data files, i.e. xmldoc directory location.
89  // Environment variable takes precedence, else use constructor input.
90  xmlPath = "";
91  const char* PYTHIA8DATA = "PYTHIA8DATA";
92  char* envPath = getenv(PYTHIA8DATA);
93  if (envPath != 0 && *envPath != '\0') {
94  int i = 0;
95  while (*(envPath+i) != '\0') xmlPath += *(envPath+(i++));
96  }
97  else xmlPath = xmlDir;
98  if (xmlPath[ xmlPath.length() - 1 ] != '/') xmlPath += "/";
99 
100  // Read in files with all flags, modes, parms and words.
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");
106  return;
107  }
108 
109  // Check that XML version number matches code version number.
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",
117  errCode.str());
118  return;
119  }
120 
121  // Read in files with all particle data.
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");
127  return;
128  }
129 
130  // Write the Pythia banner to output.
131  if (printBanner) banner();
132 
133  // Not initialized until at the end of the init() call.
134  isInit = false;
135  info.addCounter(0);
136 
137 }
138 
139 //--------------------------------------------------------------------------
140 
141 // Destructor.
142 
143 Pythia::~Pythia() {
144 
145  // Delete the PDF's created with new.
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;
152 
153  // Delete the Les Houches object created with new.
154  if (useNewLHA) delete lhaUpPtr;
155 
156  // Delete the MergingHooks object created with new.
157  if (hasOwnMergingHooks) delete mergingHooksPtr;
158 
159  // Delete the BeamShape object created with new.
160  if (useNewBeamShape) delete beamShapePtr;
161 
162  // Delete the timelike and spacelike showers created with new.
163  if (useNewTimes) delete timesPtr;
164  if (useNewSpace) delete spacePtr;
165 
166 }
167 
168 //--------------------------------------------------------------------------
169 
170 // Read in one update for a setting or particle data from a single line.
171 
172 bool Pythia::readString(string line, bool warn) {
173 
174  // Check that constructor worked.
175  if (!isConstructed) return false;
176 
177  // If empty line then done.
178  if (line.find_first_not_of(" \n\t\v\b\r\f\a") == string::npos) return true;
179 
180  // If first character is not a letter/digit, then taken to be a comment.
181  int firstChar = line.find_first_not_of(" \n\t\v\b\r\f\a");
182  if (!isalnum(line[firstChar])) return true;
183 
184  // Send on particle data to the ParticleData database.
185  if (isdigit(line[firstChar]))
186  return particleData.readString(line, warn);
187 
188  // Everything else sent on to Settings.
189  return settings.readString(line, warn);
190 
191 }
192 
193 //--------------------------------------------------------------------------
194 
195 // Read in updates for settings or particle data from user-defined file.
196 
197 bool Pythia::readFile(string fileName, bool warn, int subrun) {
198 
199  // Check that constructor worked.
200  if (!isConstructed) return false;
201 
202  // Open file for reading.
203  const char* cstring = fileName.c_str();
204  ifstream is(cstring);
205  if (!is.good()) {
206  info.errorMsg("Error in Pythia::readFile: did not find file", fileName);
207  return false;
208  }
209 
210  // Hand over real work to next method.
211  return readFile( is, warn, subrun);
212 
213 }
214 
215 //--------------------------------------------------------------------------
216 
217 // Read in updates for settings or particle data
218 // from user-defined stream (or file).
219 
220 bool Pythia::readFile(istream& is, bool warn, int subrun) {
221 
222  // Check that constructor worked.
223  if (!isConstructed) return false;
224 
225  // Read in one line at a time.
226  string line;
227  bool isCommented = false;
228  bool accepted = true;
229  int subrunNow = SUBRUNDEFAULT;
230  while ( getline(is, line) ) {
231 
232  // Check whether entering, leaving or inside commented-commands section.
233  int commentLine = readCommented( line);
234  if (commentLine == +1) isCommented = true;
235  else if (commentLine == -1) isCommented = false;
236  else if (isCommented) ;
237 
238  else {
239  // Check whether entered new subrun.
240  int subrunLine = readSubrun( line, warn);
241  if (subrunLine >= 0) subrunNow = subrunLine;
242 
243  // Process the line if in correct subrun.
244  if ( (subrunNow == subrun || subrunNow == SUBRUNDEFAULT)
245  && !readString( line, warn) ) accepted = false;
246  }
247 
248  // Reached end of input file.
249  };
250  return accepted;
251 
252 }
253 
254 //--------------------------------------------------------------------------
255 
256 // Routine to pass in pointers to PDF's. Usage optional.
257 
258 bool Pythia::setPDFPtr( PDF* pdfAPtrIn, PDF* pdfBPtrIn, PDF* pdfHardAPtrIn,
259  PDF* pdfHardBPtrIn, PDF* pdfPomAPtrIn, PDF* pdfPomBPtrIn) {
260 
261  // Delete any PDF's created in a previous init call.
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;
268 
269  // Reset pointers to be empty.
270  useNewPdfA = false;
271  useNewPdfB = false;
272  useNewPdfHard = false;
273  useNewPdfPomA = false;
274  useNewPdfPomB = false;
275  pdfAPtr = 0;
276  pdfBPtr = 0;
277  pdfHardAPtr = 0;
278  pdfHardBPtr = 0;
279  pdfPomAPtr = 0;
280  pdfPomBPtr = 0;
281 
282  // Switch off external PDF's by zero as input.
283  if (pdfAPtrIn == 0 && pdfBPtrIn == 0) return true;
284 
285  // The two PDF objects cannot be one and the same.
286  if (pdfAPtrIn == pdfBPtrIn) return false;
287 
288  // Save pointers.
289  pdfAPtr = pdfAPtrIn;
290  pdfBPtr = pdfBPtrIn;
291 
292  // By default same pointers for hard-process PDF's.
293  pdfHardAPtr = pdfAPtrIn;
294  pdfHardBPtr = pdfBPtrIn;
295 
296  // Optionally allow separate pointers for hard process.
297  if (pdfHardAPtrIn != 0 && pdfHardBPtrIn != 0) {
298  if (pdfHardAPtrIn == pdfHardBPtrIn) return false;
299  pdfHardAPtr = pdfHardAPtrIn;
300  pdfHardBPtr = pdfHardBPtrIn;
301  }
302 
303  // Optionally allow pointers for Pomerons in the proton.
304  if (pdfPomAPtrIn != 0 && pdfPomBPtrIn != 0) {
305  if (pdfPomAPtrIn == pdfPomBPtrIn) return false;
306  pdfPomAPtr = pdfPomAPtrIn;
307  pdfPomBPtr = pdfPomBPtrIn;
308  }
309 
310  // Done.
311  return true;
312 }
313 
314 //--------------------------------------------------------------------------
315 
316 // Routine to initialize with the variable values of the Beams kind.
317 
318 bool Pythia::init() {
319 
320  // Check that constructor worked.
321  isInit = false;
322  if (!isConstructed) {
323  info.errorMsg("Abort from Pythia::init: constructor "
324  "initialization failed");
325  return false;
326  }
327 
328  // Early readout, if return false or changed when no beams.
329  doProcessLevel = settings.flag("ProcessLevel:all");
330 
331  // Check that changes in Settings and ParticleData have worked.
332  if (settings.readingFailed()) {
333  info.errorMsg("Abort from Pythia::init: some user settings "
334  "did not make sense");
335  return false;
336  }
337  if (particleData.readingFailed()) {
338  info.errorMsg("Abort from Pythia::init: some user particle data "
339  "did not make sense");
340  return false;
341  }
342 
343  // Begin initialization. Find which frame type to use.
344  info.addCounter(1);
345  frameType = mode("Beams:frameType");
346 
347  // Initialization with internal processes: read in and set values.
348  if (frameType < 4 ) {
349  doLHA = false;
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");
362 
363  // Initialization with a Les Houches Event File or an LHAup object.
364  } else {
365  doLHA = true;
366  boostType = 2;
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");
372 
373  // For file input: renew file stream or (re)new Les Houches object.
374  if (frameType == 4) {
375  const char* cstring1 = lhef.c_str();
376  if (useNewLHA && skipInit) lhaUpPtr->newEventFile(cstring1);
377  else {
378  if (useNewLHA) delete lhaUpPtr;
379  // Header is optional, so use NULL pointer to indicate no value.
380  const char* cstring2 = (lhefHeader == "void") ?
381  NULL : lhefHeader.c_str();
382  lhaUpPtr = new LHAupLHEF(cstring1, cstring2, readHeaders);
383  useNewLHA = true;
384  }
385 
386  // Check that file was properly opened.
387  if (!lhaUpPtr->fileFound()) {
388  info.errorMsg("Abort from Pythia::init: "
389  "Les Houches Event File not found");
390  return false;
391  }
392 
393  // For object input: at least check that not null pointer.
394  } else {
395  if (lhaUpPtr == 0) {
396  info.errorMsg("Abort from Pythia::init: "
397  "LHAup object not found");
398  return false;
399  }
400 
401  // LHAup object generic abort using fileFound() routine.
402  if (!lhaUpPtr->fileFound()) {
403  info.errorMsg("Abort from Pythia::init: "
404  "LHAup initialisation error");
405  return false;
406  }
407  }
408 
409  // Send in pointer to info. Store or replace LHA pointer in other classes.
410  lhaUpPtr->setPtr( &info);
411  processLevel.setLHAPtr( lhaUpPtr);
412 
413  // If second time around, only with new file, then simplify.
414  // Optionally skip ahead a number of events at beginning of file.
415  if (skipInit) {
416  isInit = true;
417  info.addCounter(2);
418  if (nSkipAtInit > 0) lhaUpPtr->skipEvent(nSkipAtInit);
419  return true;
420  }
421 
422  // Set up values related to merging hooks.
423  if (frameType == 4 || frameType == 5) {
424 
425  // Set up values related to CKKW-L merging.
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");
431  // Set up values related to unitarised CKKW merging
432  bool doUMEPSTree = settings.flag("Merging:doUMEPSTree");
433  bool doUMEPSSubt = settings.flag("Merging:doUMEPSSubt");
434  // Set up values related to NL3 NLO merging
435  bool doNL3Tree = settings.flag("Merging:doNL3Tree");
436  bool doNL3Loop = settings.flag("Merging:doNL3Loop");
437  bool doNL3Subt = settings.flag("Merging:doNL3Subt");
438  // Set up values related to unitarised NLO merging
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;
448 
449  // Set up MergingHooks object.
450  bool inputMergingHooks = (mergingHooksPtr != 0);
451  if (doMerging && !inputMergingHooks){
452  if (hasOwnMergingHooks && mergingHooksPtr) delete mergingHooksPtr;
453  mergingHooksPtr = new MergingHooks();
454  hasOwnMergingHooks = true;
455  }
456 
457  hasMergingHooks = (mergingHooksPtr != 0);
458  // Merging hooks required for merging. If no merging hooks pointer is
459  // available, exit.
460  if (doMerging && !hasMergingHooks) {
461  info.errorMsg("Abort from Pythia::init: "
462  "no merging hooks object has been provided");
463  return false;
464  } else if (doMerging) {
465  string lhefIn = (frameType == 4) ? lhef : "";
466  mergingHooksPtr->setLHEInputFile( lhefIn);
467  }
468  // Initialise counting of Les Houches Events significantly above the
469  // merging scale.
470  info.setCounter(41,0);
471  }
472 
473  // Set LHAinit information (in some external program).
474  if ( !lhaUpPtr->setInit()) {
475  info.errorMsg("Abort from Pythia::init: "
476  "Les Houches initialization failed");
477  return false;
478  }
479 
480  // Extract beams from values set in an LHAinit object.
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();
489 
490  // Optionally skip ahead a number of events at beginning of file.
491  if (nSkipAtInit > 0) lhaUpPtr->skipEvent(nSkipAtInit);
492  }
493 
494  // Set up values related to user hooks.
495  hasUserHooks = (userHooksPtr != 0);
496  doVetoProcess = false;
497  doVetoPartons = false;
498  retryPartonLevel = false;
499  if (hasUserHooks) {
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");
504  return false;
505  }
506  doVetoProcess = userHooksPtr->canVetoProcessLevel();
507  doVetoPartons = userHooksPtr->canVetoPartonLevel();
508  retryPartonLevel = userHooksPtr->retryPartonLevel();
509  }
510 
511  // Back to common initialization. Reset error counters.
512  nErrEvent = 0;
513  info.setTooLowPTmin(false);
514  info.sigmaReset();
515 
516  // Initialize data members extracted from database.
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");
536 
537  // Initialise merging hooks.
538  if ( doMerging && (hasMergingHooks || hasOwnMergingHooks) )
539  mergingHooksPtr->init( settings, &info, &particleData, &partonSystems );
540 
541  // Initialize the random number generator.
542  if ( settings.flag("Random:setSeed") )
543  rndm.init( settings.mode("Random:seed") );
544 
545  // Check that combinations of settings are allowed; change if not.
546  checkSettings();
547 
548  // Initialize the SM couplings (needed to initialize resonances).
549  couplingsPtr->init( settings, &rndm );
550 
551  // Initialize SLHA interface (including SLHA/BSM couplings).
552  bool useSLHAcouplings = false;
553  slhaInterface.setPtr( &info );
554  slhaInterface.init( settings, &rndm, couplingsPtr, &particleData,
555  useSLHAcouplings );
556  if (useSLHAcouplings) couplingsPtr = slhaInterface.couplingsPtr;
557 
558  // Reset couplingsPtr to the correct memory address.
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);
563 
564  // Set headers to distinguish the two event listing kinds.
565  int startColTag = settings.mode("Event:startColTag");
566  process.init("(hard process)", &particleData, startColTag);
567  event.init("(complete event)", &particleData, startColTag);
568 
569  // Final setup stage of particle data, notably resonance widths.
570  particleData.initWidths( resonancePtrs);
571 
572  // Set up R-hadrons particle data, where relevant.
573  rHadrons.init( &info, settings, &particleData, &rndm);
574 
575  // Set up objects for timelike and spacelike showers.
576  if (timesDecPtr == 0 || timesPtr == 0) {
577  TimeShower* timesNow = new TimeShower();
578  if (timesDecPtr == 0) timesDecPtr = timesNow;
579  if (timesPtr == 0) timesPtr = timesNow;
580  useNewTimes = true;
581  }
582  if (spacePtr == 0) {
583  spacePtr = new SpaceShower();
584  useNewSpace = true;
585  }
586 
587  // Initialize showers, especially for simple showers in decays.
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);
595 
596  // Set up values related to beam shape.
597  if (beamShapePtr == 0) {
598  beamShapePtr = new BeamShape();
599  useNewBeamShape = true;
600  }
601  beamShapePtr->init( settings, &rndm);
602 
603  // Check that beams and beam combination can be handled.
604  if (!checkBeams()) {
605  info.errorMsg("Abort from Pythia::init: "
606  "checkBeams initialization failed");
607  return false;
608  }
609 
610  // Do not set up beam kinematics when no process level.
611  if (!doProcessLevel) boostType = 1;
612  else {
613 
614  // Set up beam kinematics.
615  if (!initKinematics()) {
616  info.errorMsg("Abort from Pythia::init: "
617  "kinematics initialization failed");
618  return false;
619  }
620 
621  // Set up pointers to PDFs.
622  if (!initPDFs()) {
623  info.errorMsg("Abort from Pythia::init: PDF initialization failed");
624  return false;
625  }
626 
627  // Set up the two beams and the common remnant system.
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);
633 
634  // Optionally set up new alternative beams for these Pomerons.
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);
640  }
641  }
642 
643  // Send info/pointers to process level for initialization.
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");
649  return false;
650  }
651 
652  // Alternatively only initialize resonance decays.
653  if ( !doProcessLevel) processLevel.initDecays( &info, &particleData,
654  &rndm, lhaUpPtr);
655 
656  // Send info/pointers to parton level for initialization.
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" );
663  return false;
664  }
665 
666  // Alternatively only initialize final-state showers in resonance decays.
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);
670 
671  // Send info/pointers to parton level for trial shower initialization.
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");
678  return false;
679  }
680 
681  // Initialise the merging wrapper class.
682  if (doMerging ) merging.init( &settings, &info, &particleData, &rndm,
683  &beamA, &beamB, mergingHooksPtr, &trialPartonLevel );
684 
685  // Send info/pointers to hadron level for initialization.
686  // Note: forceHadronLevel() can come, so we must always initialize.
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");
692  return false;
693  }
694 
695  // Optionally check particle data table for inconsistencies.
696  if ( settings.flag("Check:particleData") )
697  particleData.checkTable( settings.mode("Check:levelParticleData") );
698 
699  // Optionally show settings and particle data, changed or all.
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();
712 
713  // Listing options for the next() routine.
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");
721 
722  // Succeeded.
723  isInit = true;
724  info.addCounter(2);
725  if (useNewLHA && showPro) lhaUpPtr->listInit();
726  return true;
727 
728 }
729 
730 //--------------------------------------------------------------------------
731 
732 // Routine to initialize with CM energy.
733 
734 bool Pythia::init( int idAin, int idBin, double eCMin) {
735 
736  // Set arguments in Settings database.
737  settings.mode("Beams:idA", idAin);
738  settings.mode("Beams:idB", idBin);
739  settings.mode("Beams:frameType", 1);
740  settings.parm("Beams:eCM", eCMin);
741 
742  // Send on to common initialization.
743  return init();
744 
745 }
746 
747 //--------------------------------------------------------------------------
748 
749 // Routine to initialize with two collinear beams, energies specified.
750 
751 bool Pythia::init( int idAin, int idBin, double eAin, double eBin) {
752 
753  // Set arguments in Settings database.
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);
759 
760  // Send on to common initialization.
761  return init();
762 
763 }
764 
765 //--------------------------------------------------------------------------
766 
767 // Routine to initialize with two beams specified by three-momenta.
768 
769 bool Pythia::init( int idAin, int idBin, double pxAin, double pyAin,
770  double pzAin, double pxBin, double pyBin, double pzBin) {
771 
772  // Set arguments in Settings database.
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);
782 
783  // Send on to common initialization.
784  return init();
785 
786 }
787 
788 //--------------------------------------------------------------------------
789 
790 // Routine to initialize when all info is given in a Les Houches Event File.
791 
792 bool Pythia::init( string LesHouchesEventFile, bool skipInit) {
793 
794  // Set arguments in Settings database.
795  settings.mode("Beams:frameType", 4);
796  settings.word("Beams:LHEF", LesHouchesEventFile);
797  settings.flag("Beams:newLHEFsameInit", skipInit);
798 
799  // Send on to common initialization.
800  return init();
801 
802 }
803 
804 //--------------------------------------------------------------------------
805 
806 // Routine to initialize when beam info is given in an LHAup object.
807 
808 bool Pythia::init( LHAup* lhaUpPtrIn) {
809 
810  // Store LHAup object and set arguments in Settings database.
811  setLHAupPtr( lhaUpPtrIn);
812  settings.mode("Beams:frameType", 5);
813  string lhaFileName = lhaUpPtr->getFileName();
814  if (lhaFileName != "void") settings.word("Beams:LHEF", lhaFileName);
815 
816  // Send on to common initialization.
817  return init();
818 
819 }
820 
821 //--------------------------------------------------------------------------
822 
823 // Check that combinations of settings are allowed; change if not.
824 
825 void Pythia::checkSettings() {
826 
827  // Double rescattering not allowed if ISR or FSR.
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);
833  }
834 
835 }
836 
837 //--------------------------------------------------------------------------
838 
839 // Check that beams and beam combination can be handled. Set up unresolved.
840 
841 bool Pythia::checkBeams() {
842 
843  // Absolute flavours. If not to do process level then no check needed.
844  int idAabs = abs(idA);
845  int idBabs = abs(idB);
846  if (!doProcessLevel) return true;
847 
848  // Neutrino beams always unresolved, charged lepton ones conditionally.
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);
854 
855  // Lepton-lepton collisions OK (including neutrinos) if both (un)resolved.
856  if (isLeptonA && isLeptonB && isUnresolvedA == isUnresolvedB) return true;
857 
858  // MBR model only implemented for pp/ppbar/pbarp collisions.
859  int PomFlux = settings.mode("Diffraction:PomFlux");
860  if (PomFlux == 5) {
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");
866  return false;
867  }
868 
869  // Hadron-hadron collisions OK, with Pomeron counted as hadron.
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;
875 
876  // If no case above then failed.
877  info.errorMsg("Error in Pythia::init: cannot handle this beam combination");
878  return false;
879 
880 }
881 
882 //--------------------------------------------------------------------------
883 
884 // Calculate kinematics at initialization. Store beam four-momenta.
885 
886 bool Pythia::initKinematics() {
887 
888  // Find masses. Initial guess that we are in CM frame.
889  mA = particleData.m0(idA);
890  mB = particleData.m0(idB);
891  betaZ = 0.;
892  gammaZ = 1.;
893 
894  // Collinear beams not in CM frame: find CM energy.
895  if (boostType == 2) {
896  eA = max(eA, mA);
897  eB = max(eB, mB);
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) );
903 
904  // Find boost to rest frame.
905  betaZ = (pzA + pzB) / (eA + eB);
906  gammaZ = (eA + eB) / eCM;
907  if (abs(betaZ) < 1e-10) boostType = 1;
908  }
909 
910  // Completely general beam directions: find CM energy.
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();
917 
918  // Find boost+rotation needed to move from/to CM frame.
919  MfromCM.reset();
920  MfromCM.fromCMframe( pAinit, pBinit);
921  MtoCM = MfromCM;
922  MtoCM.invert();
923  }
924 
925  // Fail if CM energy below beam masses.
926  if (eCM < mA + mB) {
927  info.errorMsg("Error in Pythia::initKinematics: too low energy");
928  return false;
929  }
930 
931  // Set up CM-frame kinematics with beams along +-z axis.
932  pzAcm = 0.5 * sqrtpos( (eCM + mA + mB) * (eCM - mA - mB)
933  * (eCM - mA + mB) * (eCM + mA - mB) ) / eCM;
934  pzBcm = -pzAcm;
935  eA = sqrt(mA*mA + pzAcm*pzAcm);
936  eB = sqrt(mB*mB + pzBcm*pzBcm);
937 
938  // If in CM frame then store beam four-vectors (else already done above).
939  if (boostType != 2 && boostType != 3) {
940  pAinit = Vec4( 0., 0., pzAcm, eA);
941  pBinit = Vec4( 0., 0., pzBcm, eB);
942  }
943 
944  // Store main info for access in process generation.
945  info.setBeamA( idA, pzAcm, eA, mA);
946  info.setBeamB( idB, pzBcm, eB, mB);
947  info.setECM( eCM);
948 
949  // Must allow for generic boost+rotation when beam momentum spread.
950  if (doMomentumSpread) boostType = 3;
951 
952  // Done.
953  return true;
954 
955 }
956 
957 //--------------------------------------------------------------------------
958 
959 // Set up pointers to PDFs.
960 
961 bool Pythia::initPDFs() {
962 
963  // Delete any PDF's created in a previous init call.
964  if (useNewPdfHard) {
965  if (pdfHardAPtr != pdfAPtr) {
966  delete pdfHardAPtr;
967  pdfHardAPtr = 0;
968  }
969  if (pdfHardBPtr != pdfBPtr) {
970  delete pdfHardBPtr;
971  pdfHardBPtr = 0;
972  }
973  useNewPdfHard = false;
974  }
975  if (useNewPdfA) {
976  delete pdfAPtr;
977  useNewPdfA = false;
978  pdfAPtr = 0;
979  }
980  if (useNewPdfB) {
981  delete pdfBPtr;
982  useNewPdfB = false;
983  pdfBPtr = 0;
984  }
985  if (useNewPdfPomA) {
986  delete pdfPomAPtr;
987  useNewPdfPomA = false;
988  pdfPomAPtr = 0;
989  }
990  if (useNewPdfPomB) {
991  delete pdfPomBPtr;
992  useNewPdfPomB = false;
993  pdfPomBPtr = 0;
994  }
995 
996  // Set up the PDF's, if not already done.
997  if (pdfAPtr == 0) {
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");
1002  return false;
1003  }
1004  pdfHardAPtr = pdfAPtr;
1005  useNewPdfA = true;
1006  }
1007  if (pdfBPtr == 0) {
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");
1012  return false;
1013  }
1014  pdfHardBPtr = pdfBPtr;
1015  useNewPdfB = true;
1016  }
1017 
1018  // Optionally set up separate PDF's for hard process.
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;
1025  }
1026 
1027  // Optionally set up Pomeron PDF's for diffractive physics.
1028  if ( doDiffraction) {
1029  if (pdfPomAPtr == 0) {
1030  pdfPomAPtr = getPDFPtr(990);
1031  useNewPdfPomA = true;
1032  }
1033  if (pdfPomBPtr == 0) {
1034  pdfPomBPtr = getPDFPtr(990);
1035  useNewPdfPomB = true;
1036  }
1037  }
1038 
1039  // Done.
1040  return true;
1041 
1042 }
1043 
1044 //--------------------------------------------------------------------------
1045 
1046 // Main routine to generate the next event, using internal machinery.
1047 
1048 bool Pythia::next() {
1049 
1050  // Check that constructor worked.
1051  if (!isConstructed) return false;
1052 
1053  // Regularly print how many events have been generated.
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;
1058 
1059  // Set/reset info counters specific to each event.
1060  info.addCounter(3);
1061  for (int i = 10; i < 13; ++i) info.setCounter(i);
1062 
1063  // Simpler option when no hard process, i.e. mainly hadron level.
1064  if (!doProcessLevel) {
1065 
1066  // Optionally fetch in resonance decays from LHA interface.
1067  if (doLHA && !processLevel.nextLHAdec( event)) {
1068  if (info.atEndOfFile()) info.errorMsg("Abort from Pythia::next:"
1069  " reached end of Les Houches Events File");
1070  return false;
1071  }
1072 
1073  // Reset info array (while event record contains data).
1074  info.clear();
1075 
1076  // Set correct energy for system.
1077  Vec4 pSum = 0.;
1078  for (int i = 1; i < event.size(); ++i)
1079  if (event[i].isFinal()) pSum += event[i].p();
1080  event[0].p( pSum );
1081  event[0].m( pSum.mCalc() );
1082 
1083  // Generate hadronization and decays.
1084  bool status = forceHadronLevel();
1085  if (status) info.addCounter(4);
1086  if (status && nPrevious < nShowEvt) event.list(showSaV, showMaD);
1087  return status;
1088  }
1089 
1090  // Reset arrays.
1091  info.clear();
1092  process.clear();
1093  event.clear();
1094  partonSystems.clear();
1095  beamA.clear();
1096  beamB.clear();
1097  beamPomA.clear();
1098  beamPomB.clear();
1099 
1100  // Can only generate event if initialization worked.
1101  if (!isInit) {
1102  info.errorMsg("Abort from Pythia::next: "
1103  "not properly initialized so cannot generate events");
1104  return false;
1105  }
1106 
1107  // Pick beam momentum spread and beam vertex.
1108  if (doMomentumSpread || doVertexSpread) beamShapePtr->pick();
1109 
1110  // Recalculate kinematics when beam momentum spread.
1111  if (doMomentumSpread) nextKinematics();
1112 
1113  // Outer loop over hard processes; only relevant for user-set vetoes.
1114  for ( ; ; ) {
1115 
1116  info.addCounter(10);
1117  bool hasVetoed = false;
1118 
1119  // Provide the hard process that starts it off. Only one try.
1120  info.clear();
1121  process.clear();
1122 
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");
1128  return false;
1129  }
1130 
1131  info.addCounter(11);
1132 
1133  // Possibility for a user veto of the process-level event.
1134  if (doVetoProcess) {
1135  hasVetoed = userHooksPtr->doVetoProcessLevel( process);
1136  if (hasVetoed) {
1137  if (abortIfVeto) return false;
1138  continue;
1139  }
1140  }
1141 
1142  // Possibility to perform matrix element merging for this event.
1143  if (doMerging) {
1144  int veto = merging.mergeProcess( process );
1145  // Apply possible merging scale cut.
1146  if ( veto == -1 ) {
1147  hasVetoed = true;
1148  if (abortIfVeto) return false;
1149  continue;
1150  // Exit because of vanishing no-emission probability.
1151  } else if ( veto == 0 ) break;
1152 
1153  // Redo resonance decays after the merging, in case the resonance
1154  // structure has been changed because of reclusterings.
1155  if (veto == 2 && doResDec) processLevel.nextDecays( process);
1156  }
1157 
1158  // Possibility to stop the generation at this stage.
1159  if (!doPartonLevel) {
1160  boostAndVertex( true, true);
1161  processLevel.accumulate();
1162  info.addCounter(4);
1163  if (doLHA && nPrevious < nShowLHA) lhaUpPtr->listEvent();
1164  if (nPrevious < nShowInfo) info.list();
1165  if (nPrevious < nShowProc) process.list(showSaV, showMaD);
1166  return true;
1167  }
1168 
1169  // Save spare copy of process record in case of problems.
1170  Event processSave = process;
1171  int sizeMPI = info.sizeMPIarrays();
1172  info.addCounter(12);
1173  for (int i = 14; i < 19; ++i) info.setCounter(i);
1174 
1175  // Allow up to ten tries for parton- and hadron-level processing.
1176  bool physical = true;
1177  for (int iTry = 0; iTry < NTRY; ++iTry) {
1178 
1179  info.addCounter(14);
1180  physical = true;
1181  hasVetoed = false;
1182 
1183  // Restore original process record if problems.
1184  if (iTry > 0) process = processSave;
1185  if (iTry > 0) info.resizeMPIarrays( sizeMPI);
1186 
1187  // Reset event record and (extracted partons from) beam remnants.
1188  event.clear();
1189  beamA.clear();
1190  beamB.clear();
1191  beamPomA.clear();
1192  beamPomB.clear();
1193  partonSystems.clear();
1194 
1195  // Parton-level evolution: ISR, FSR, MPI.
1196  if ( !partonLevel.next( process, event) ) {
1197  // Skip to next hard process for failure owing to deliberate veto,
1198  // or alternatively retry for the same hard process.
1199  hasVetoed = partonLevel.hasVetoed();
1200  if (hasVetoed) {
1201  if (retryPartonLevel) {
1202  --iTry;
1203  continue;
1204  }
1205  if (abortIfVeto) return false;
1206  break;
1207  }
1208  // Else make a new try for other failures.
1209  info.errorMsg("Error in Pythia::next: "
1210  "partonLevel failed; try again");
1211  physical = false;
1212  continue;
1213  }
1214  info.addCounter(15);
1215 
1216  // Possibility for a user veto of the parton-level event.
1217  if (doVetoPartons) {
1218  hasVetoed = userHooksPtr->doVetoPartonLevel( event);
1219  if (hasVetoed) {
1220  if (abortIfVeto) return false;
1221  break;
1222  }
1223  }
1224 
1225  // Boost to lab frame (before decays, for vertices).
1226  boostAndVertex( true, true);
1227 
1228  // Possibility to stop the generation at this stage.
1229  if (!doHadronLevel) {
1230  processLevel.accumulate();
1231  partonLevel.accumulate();
1232  // Optionally check final event for problems.
1233  if (checkEvent && !check()) {
1234  info.errorMsg("Abort from Pythia::next: "
1235  "check of event revealed problems");
1236  return false;
1237  }
1238  info.addCounter(4);
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);
1243  return true;
1244  }
1245 
1246  // Hadron-level: hadronization, decays.
1247  info.addCounter(16);
1248  if ( !hadronLevel.next( event) ) {
1249  info.errorMsg("Error in Pythia::next: "
1250  "hadronLevel failed; try again");
1251  physical = false;
1252  continue;
1253  }
1254 
1255  // If R-hadrons have been formed, then (optionally) let them decay.
1256  if (decayRHadrons && rHadrons.exist() && !doRHadronDecays()) {
1257  info.errorMsg("Error in Pythia::next: "
1258  "decayRHadrons failed; try again");
1259  physical = false;
1260  continue;
1261  }
1262  info.addCounter(17);
1263 
1264  // Optionally check final event for problems.
1265  if (checkEvent && !check()) {
1266  info.errorMsg("Error in Pythia::next: "
1267  "check of event revealed problems");
1268  physical = false;
1269  continue;
1270  }
1271 
1272  // Stop parton- and hadron-level looping if you got this far.
1273  info.addCounter(18);
1274  break;
1275  }
1276 
1277  // If event vetoed then to make a new try.
1278  if (hasVetoed) {
1279  if (abortIfVeto) return false;
1280  continue;
1281  }
1282 
1283  // If event failed any other way (after ten tries) then give up.
1284  if (!physical) {
1285  info.errorMsg("Abort from Pythia::next: "
1286  "parton+hadronLevel failed; giving up");
1287  return false;
1288  }
1289 
1290  // Process- and parton-level statistics. Event scale.
1291  processLevel.accumulate();
1292  partonLevel.accumulate();
1293  event.scale( process.scale() );
1294 
1295  // End of outer loop over hard processes. Done with normal option.
1296  info.addCounter(13);
1297  break;
1298  }
1299 
1300  // List events.
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);
1305 
1306  // Done.
1307  info.addCounter(4);
1308  return true;
1309 
1310 }
1311 
1312 //--------------------------------------------------------------------------
1313 
1314 // Generate only the hadronization/decay stage, using internal machinery.
1315 // The "event" instance should already contain a parton-level configuration.
1316 
1317 bool Pythia::forceHadronLevel(bool findJunctions) {
1318 
1319  // Can only generate event if initialization worked.
1320  if (!isInit) {
1321  info.errorMsg("Abort from Pythia::forceHadronLevel: "
1322  "not properly initialized so cannot generate events");
1323  return false;
1324  }
1325 
1326  // Check whether any junctions in system. (Normally done in ProcessLevel.)
1327  // Avoid it if there are no final-state coloured partons.
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);
1334  break;
1335  }
1336  }
1337 
1338  // Save spare copy of event in case of failure.
1339  Event spareEvent = event;
1340 
1341  // Allow up to ten tries for hadron-level processing.
1342  bool physical = true;
1343  for (int iTry = 0; iTry < NTRY; ++ iTry) {
1344  physical = true;
1345 
1346  // Check whether any resonances need to be handled at process level.
1347  if (doResDec) {
1348  process = event;
1349  processLevel.nextDecays( process);
1350 
1351  // Allow for showers if decays happened at process level.
1352  if (process.size() > event.size()) {
1353  if (doFSRinRes) {
1354  partonLevel.setupShowerSys( process, event);
1355  partonLevel.resonanceShowers( process, event, false);
1356  } else event = process;
1357  }
1358  }
1359 
1360  // Hadron-level: hadronization, decays.
1361  if (hadronLevel.next( event)) break;
1362 
1363  // If failure then warn, restore original configuration and try again.
1364  info.errorMsg("Error in Pythia::forceHadronLevel: "
1365  "hadronLevel failed; try again");
1366  physical = false;
1367  event = spareEvent;
1368  }
1369 
1370  // Done for simpler option.
1371  if (!physical) {
1372  info.errorMsg("Abort from Pythia::forceHadronLevel: "
1373  "hadronLevel failed; giving up");
1374  return false;
1375  }
1376 
1377  // Optionally check final event for problems.
1378  if (checkEvent && !check()) {
1379  info.errorMsg("Abort from Pythia::forceHadronLevel: "
1380  "check of event revealed problems");
1381  return false;
1382  }
1383 
1384  // Done.
1385  return true;
1386 
1387 }
1388 
1389 //--------------------------------------------------------------------------
1390 
1391 // Recalculate kinematics for each event when beam momentum has a spread.
1392 
1393 void Pythia::nextKinematics() {
1394 
1395  // Read out momentum shift to give current beam momenta.
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) );
1400 
1401  // Construct CM frame kinematics.
1402  eCM = (pAnow + pBnow).mCalc();
1403  pzAcm = 0.5 * sqrtpos( (eCM + mA + mB) * (eCM - mA - mB)
1404  * (eCM - mA + mB) * (eCM + mA - mB) ) / eCM;
1405  pzBcm = -pzAcm;
1406  eA = sqrt(mA*mA + pzAcm*pzAcm);
1407  eB = sqrt(mB*mB + pzBcm*pzBcm);
1408 
1409  // Set relevant info for other classes to use.
1410  info.setBeamA( idA, pzAcm, eA, mA);
1411  info.setBeamB( idB, pzBcm, eB, mB);
1412  info.setECM( eCM);
1413  beamA.newPzE( pzAcm, eA);
1414  beamB.newPzE( pzBcm, eB);
1415 
1416  // Set boost/rotation matrices from/to CM frame.
1417  MfromCM.reset();
1418  MfromCM.fromCMframe( pAnow, pBnow);
1419  MtoCM = MfromCM;
1420  MtoCM.invert();
1421 
1422 }
1423 
1424 //--------------------------------------------------------------------------
1425 
1426 // Boost from CM frame to lab frame, or inverse. Set production vertex.
1427 
1428 void Pythia::boostAndVertex( bool toLab, bool setVertex) {
1429 
1430  // Boost process from CM frame to lab frame.
1431  if (toLab) {
1432  if (boostType == 2) process.bst(0., 0., betaZ, gammaZ);
1433  else if (boostType == 3) process.rotbst(MfromCM);
1434 
1435  // Boost nonempty event from CM frame to lab frame.
1436  if (event.size() > 0) {
1437  if (boostType == 2) event.bst(0., 0., betaZ, gammaZ);
1438  else if (boostType == 3) event.rotbst(MfromCM);
1439  }
1440 
1441  // Boost process from lab frame to CM frame.
1442  } else {
1443  if (boostType == 2) process.bst(0., 0., -betaZ, gammaZ);
1444  else if (boostType == 3) process.rotbst(MtoCM);
1445 
1446  // Boost nonempty event from lab frame to CM frame.
1447  if (event.size() > 0) {
1448  if (boostType == 2) event.bst(0., 0., -betaZ, gammaZ);
1449  else if (boostType == 3) event.rotbst(MtoCM);
1450  }
1451  }
1452 
1453  // Set production vertex; assumes particles are in lab frame and at origin.
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);
1458  }
1459 
1460 }
1461 
1462 //--------------------------------------------------------------------------
1463 
1464 // Perform R-hadron decays, either as part of normal evolution or forced.
1465 
1466 bool Pythia::doRHadronDecays( ) {
1467 
1468  // Check if R-hadrons exist to be processed.
1469  if ( !rHadrons.exist() ) return true;
1470 
1471  // Do the R-hadron decay itself.
1472  if ( !rHadrons.decay( event) ) return false;
1473 
1474  // Perform showers in resonance decay chains.
1475  if ( !partonLevel.resonanceShowers( process, event, false) ) return false;
1476 
1477  // Subsequent hadronization and decays.
1478  if ( !hadronLevel.next( event) ) return false;
1479 
1480  // Done.
1481  return true;
1482 
1483 }
1484 
1485 //--------------------------------------------------------------------------
1486 
1487 // Print statistics on event generation.
1488 
1489 void Pythia::stat() {
1490 
1491  // Read out settings for what to include.
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");
1496 
1497  // Statistics on cross section and number of events.
1498  if (doProcessLevel) {
1499  if (showPrL) processLevel.statistics(false);
1500  if (reset) processLevel.resetStatistics();
1501  }
1502 
1503  // Statistics from other classes, currently multiparton interactions.
1504  if (showPaL) partonLevel.statistics(false);
1505  if (reset) partonLevel.resetStatistics();
1506 
1507  // Merging statistics.
1508  if (doMerging) merging.statistics();
1509 
1510  // Summary of which and how many warnings/errors encountered.
1511  if (showErr) info.errorStatistics();
1512  if (reset) info.errorReset();
1513 
1514 }
1515 
1516 //--------------------------------------------------------------------------
1517 
1518 // Print statistics on event generation - deprecated version.
1519 
1520 void Pythia::statistics(bool all, bool reset) {
1521 
1522  // Statistics on cross section and number of events.
1523  if (doProcessLevel) processLevel.statistics(reset);
1524 
1525  // Statistics from other classes, e.g. multiparton interactions.
1526  if (all) partonLevel.statistics(reset);
1527 
1528  // Merging statistics.
1529  if (doMerging) merging.statistics();
1530 
1531  // Summary of which and how many warnings/errors encountered.
1532  info.errorStatistics();
1533  if (reset) info.errorReset();
1534 
1535 }
1536 
1537 //--------------------------------------------------------------------------
1538 
1539 // Write the Pythia banner, with symbol and version information.
1540 
1541 void Pythia::banner(ostream& os) {
1542 
1543  // Read in version number and last date of change.
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"};
1548 
1549  // Get date and time.
1550  time_t t = time(0);
1551  char dateNow[12];
1552  strftime(dateNow,12,"%d %b %Y",localtime(&t));
1553  char timeNow[9];
1554  strftime(timeNow,9,"%H:%M:%S",localtime(&t));
1555 
1556  os << "\n"
1557  << " *-------------------------------------------"
1558  << "-----------------------------------------* \n"
1559  << " | "
1560  << " | \n"
1561  << " | *----------------------------------------"
1562  << "--------------------------------------* | \n"
1563  << " | | "
1564  << " | | \n"
1565  << " | | "
1566  << " | | \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"
1577  << " | | \n"
1578  << " | | P Y T H H III A A"
1579  << " Now is " << dateNow << " at " << timeNow << " | | \n"
1580  << " | | "
1581  << " | | \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 "
1605  << " | | \n"
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, "
1613  << " | | \n"
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"
1626  << " | | "
1627  << " | | \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] "
1633  << " | | \n"
1634  << " | | "
1635  << " | | \n"
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"
1640  << " | | "
1641  << " | | \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"
1645  << "a.html | | \n"
1646  << " | | "
1647  << " | | \n"
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"
1652  << " | | "
1653  << " | | \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"
1658  << " | | "
1659  << " | | \n"
1660  << " | | Copyright (C) 2014 Torbjorn Sjostrand"
1661  << " | | \n"
1662  << " | | "
1663  << " | | \n"
1664  << " | | "
1665  << " | | \n"
1666  << " | *----------------------------------------"
1667  << "--------------------------------------* | \n"
1668  << " | "
1669  << " | \n"
1670  << " *-------------------------------------------"
1671  << "-----------------------------------------* \n" << endl;
1672 
1673 }
1674 
1675 //--------------------------------------------------------------------------
1676 
1677 // Check for lines in file that mark the beginning of new subrun.
1678 
1679 int Pythia::readSubrun(string line, bool warn, ostream& os) {
1680 
1681  // If empty line then done.
1682  int subrunLine = SUBRUNDEFAULT;
1683  if (line.find_first_not_of(" \n\t\v\b\r\f\a") == string::npos)
1684  return subrunLine;
1685 
1686  // If first character is not a letter, then done.
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;
1690 
1691  // Replace an equal sign by a blank to make parsing simpler.
1692  while (lineNow.find("=") != string::npos) {
1693  int firstEqual = lineNow.find_first_of("=");
1694  lineNow.replace(firstEqual, 1, " ");
1695  }
1696 
1697  // Get first word of a line.
1698  istringstream splitLine(lineNow);
1699  string name;
1700  splitLine >> name;
1701 
1702  // Replace two colons by one (:: -> :) to allow for such mistakes.
1703  while (name.find("::") != string::npos) {
1704  int firstColonColon = name.find_first_of("::");
1705  name.replace(firstColonColon, 2, ":");
1706  }
1707 
1708 
1709  // Convert to lowercase.
1710  for (int i = 0; i < int(name.length()); ++i) name[i] = tolower(name[i]);
1711 
1712  // If no match then done.
1713  if (name != "main:subrun") return subrunLine;
1714 
1715  // Else find new subrun number and return it.
1716  splitLine >> subrunLine;
1717  if (!splitLine) {
1718  if (warn) os << "\n PYTHIA Warning: Main:subrun number not"
1719  << " recognized; skip:\n " << line << endl;
1720  subrunLine = SUBRUNDEFAULT;
1721  }
1722  return subrunLine;
1723 
1724 }
1725 
1726 //--------------------------------------------------------------------------
1727 
1728 // Check for lines in file that mark the beginning or end of commented section.
1729 // Return +1 for beginning, -1 for end, 0 else.
1730 
1731 int Pythia::readCommented(string line) {
1732 
1733  // If less than two nontrivial characters on line then done.
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;
1737 
1738  // If first two nontrivial characters are /* or */ then done.
1739  if (line.substr(firstChar, 2) == "/*") return +1;
1740  if (line.substr(firstChar, 2) == "*/") return -1;
1741 
1742  // Else done.
1743  return 0;
1744 
1745 }
1746 
1747 //--------------------------------------------------------------------------
1748 
1749 // Check that the final event makes sense: no unknown id codes;
1750 // charge and energy-momentum conserved.
1751 
1752 bool Pythia::check(ostream& os) {
1753 
1754  // Reset.
1755  bool physical = true;
1756  bool listVertices = false;
1757  bool listHistory = false;
1758  bool listSystems = false;
1759  bool listBeams = false;
1760  iErrId.resize(0);
1761  iErrCol.resize(0);
1762  iErrEpm.resize(0);
1763  iErrNan.resize(0);
1764  iErrNanVtx.resize(0);
1765  Vec4 pSum;
1766  double chargeSum = 0.;
1767 
1768  // Incoming beams counted with negative momentum and charge.
1769  if (doProcessLevel) {
1770  pSum = - (event[1].p() + event[2].p());
1771  chargeSum = - (event[1].charge() + event[2].charge());
1772 
1773  // If no ProcessLevel then sum final state of process record.
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();
1778 
1779  // If process not filled, then use outgoing primary in event.
1780  } else {
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();
1785  }
1786  double eLab = abs(pSum.e());
1787 
1788  // Loop over particles in the event.
1789  for (int i = 0; i < event.size(); ++i) {
1790 
1791  // Look for any unrecognized particle codes.
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());
1798  physical = false;
1799  iErrId.push_back(i);
1800 
1801  // Check that colour assignments are the expected ones.
1802  } else {
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
1812  << " " << acol;
1813  info.errorMsg("Error in Pythia::check: "
1814  "incorrect colours", errCode.str());
1815  physical = false;
1816  iErrCol.push_back(i);
1817  }
1818  }
1819 
1820  // Look for particles with mismatched or not-a-number energy/momentum/mass.
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");
1827  physical = false;
1828  iErrEpm.push_back(i);
1829  }
1830  } else {
1831  info.errorMsg("Error in Pythia::check: "
1832  "not-a-number energy/momentum/mass");
1833  physical = false;
1834  iErrNan.push_back(i);
1835  }
1836 
1837  // Look for particles with not-a-number vertex/lifetime.
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.) ;
1841  else {
1842  info.errorMsg("Error in Pythia::check: "
1843  "not-a-number vertex/lifetime");
1844  physical = false;
1845  listVertices = true;
1846  iErrNanVtx.push_back(i);
1847  }
1848 
1849  // Add final-state four-momentum and charge.
1850  if (event[i].isFinal()) {
1851  pSum += event[i].p();
1852  chargeSum += event[i].charge();
1853  }
1854 
1855  // End of particle loop.
1856  }
1857 
1858  // Check energy-momentum/charge conservation.
1859  double epDev = abs(pSum.e()) + abs(pSum.px()) + abs(pSum.py())
1860  + abs(pSum.pz());
1861  if (epDev > epTolErr * eLab) {
1862  info.errorMsg("Error in Pythia::check: energy-momentum not conserved");
1863  physical = false;
1864  } else if (epDev > epTolWarn * eLab) {
1865  info.errorMsg("Warning in Pythia::check: "
1866  "energy-momentum not quite conserved");
1867  }
1868  if (abs(chargeSum) > 0.1) {
1869  info.errorMsg("Error in Pythia::check: charge not conserved");
1870  physical = false;
1871  }
1872 
1873  // Check that beams and event records agree on incoming partons.
1874  // Only meaningful for resolved beams.
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");
1884  physical = false;
1885  listSystems = true;
1886  listBeams = true;
1887  }
1888  }
1889 
1890  // Check that mother and daughter information match for each particle.
1891  vector<int> noMot;
1892  vector<int> noDau;
1893  vector< pair<int,int> > noMotDau;
1894  if (checkHistory) {
1895 
1896  // Loop through the event and check that there are beam particles.
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;
1901 
1902  // Check that mother and daughter lists not empty where not expected to.
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)
1906  noMot.push_back(i);
1907  if (dList.size() == 0 && status < 0 && status != -11)
1908  noDau.push_back(i);
1909 
1910  // Check that the particle appears in the daughters list of each mother.
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) {
1918  foundMatch = true;
1919  break;
1920  }
1921  if (!hasBeams && mList.size() == 1 && mList[0] == 0) foundMatch = true;
1922  if (!foundMatch) {
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) {
1926  oldPair = true;
1927  break;
1928  }
1929  if (!oldPair) noMotDau.push_back( make_pair( mList[j], i) );
1930  }
1931  }
1932 
1933  // Check that the particle appears in the mothers list of each daughter.
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) {
1943  foundMatch = true;
1944  break;
1945  }
1946  if (!foundMatch) {
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]) {
1950  oldPair = true;
1951  break;
1952  }
1953  if (!oldPair) noMotDau.push_back( make_pair( i, dList[j]) );
1954  }
1955  }
1956  }
1957 
1958  // Warn if any errors were found.
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");
1962  physical = false;
1963  listHistory = true;
1964  }
1965  }
1966 
1967  // Done for sensible events.
1968  if (physical) return true;
1969 
1970  // Print (the first few) flawed events: local info.
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] << " ";
1977  os << "\n";
1978  }
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] << " ";
1983  os << "\n";
1984  }
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] << " ";
1989  os << "\n";
1990  }
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] << " ";
1995  os << "\n";
1996  }
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] << " ";
2001  os << "\n";
2002  }
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] << " ";
2010  os << "\n";
2011  }
2012  if (noDau.size() > 0) {
2013  os << " missing daughters for particles ";
2014  for (int i = 0; i < int(noDau.size()); ++i) os << noDau[i] << " ";
2015  os << "\n";
2016  }
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 << ") ";
2021  os << "\n";
2022  }
2023 
2024  // Print (the first few) flawed events: standard listings.
2025  info.list();
2026  event.list(listVertices, listHistory);
2027  if (listSystems) partonSystems.list();
2028  if (listBeams) beamA.list();
2029  if (listBeams) beamB.list();
2030  }
2031 
2032  // Update error counter. Done also for flawed event.
2033  ++nErrEvent;
2034  return false;
2035 
2036 }
2037 
2038 //--------------------------------------------------------------------------
2039 
2040 // Routine to set up a PDF pointer.
2041 
2042 PDF* Pythia::getPDFPtr(int idIn, int sequence) {
2043 
2044  // Temporary pointer to be returned.
2045  PDF* tempPDFPtr = 0;
2046 
2047  // One option is to treat a Pomeron like a pi0.
2048  if (idIn == 990 && settings.mode("PDF:PomSet") == 2) idIn = 111;
2049 
2050  // Proton beam, normal choice. Also used for neutron.
2051  if ((abs(idIn) == 2212 || abs(idIn) == 2112) && sequence == 1) {
2052  int pSet = settings.mode("PDF:pSet");
2053  bool useLHAPDF = settings.flag("PDF:useLHAPDF");
2054 
2055  // Use internal sets.
2056  if (!useLHAPDF) {
2057  if (pSet == 1) tempPDFPtr = new GRV94L(idIn);
2058  else if (pSet == 2) tempPDFPtr = new CTEQ5L(idIn);
2059  else if (pSet <= 6)
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);
2064  }
2065 
2066  // Use sets from LHAPDF.
2067  else {
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);
2073 
2074  // Optionally allow extrapolation beyond x and Q2 limits.
2075  tempPDFPtr->setExtrapolate( settings.flag("PDF:extrapolateLHAPDF") );
2076  }
2077  }
2078 
2079  // Proton beam, special choice for the hard process. Also neutron.
2080  else if (abs(idIn) == 2212 || abs(idIn) == 2112) {
2081  int pSet = settings.mode("PDF:pHardSet");
2082  bool useLHAPDF = settings.flag("PDF:useHardLHAPDF");
2083 
2084  // Use internal sets.
2085  if (!useLHAPDF) {
2086  if (pSet == 1) tempPDFPtr = new GRV94L(idIn);
2087  else if (pSet == 2) tempPDFPtr = new CTEQ5L(idIn);
2088  else if (pSet <= 6)
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);
2093  }
2094 
2095  // Use sets from LHAPDF.
2096  else {
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);
2102 
2103  // Optionally allow extrapolation beyond x and Q2 limits.
2104  tempPDFPtr->setExtrapolate( settings.flag("PDF:extrapolateLHAPDF") );
2105  }
2106  }
2107 
2108  // Pion beam (or, in one option, Pomeron beam).
2109  else if (abs(idIn) == 211 || idIn == 111) {
2110  bool useLHAPDF = settings.flag("PDF:piUseLHAPDF");
2111 
2112  // Use internal sets.
2113  if (!useLHAPDF) {
2114  tempPDFPtr = new GRVpiL(idIn);
2115  }
2116 
2117  // Use sets from LHAPDF.
2118  else {
2119  string LHAPDFset = settings.word("PDF:piLHAPDFset");
2120  int LHAPDFmember = settings.mode("PDF:piLHAPDFmember");
2121  tempPDFPtr = new LHAPDF(idIn, LHAPDFset, LHAPDFmember, 1, &info);
2122 
2123  // Optionally allow extrapolation beyond x and Q2 limits.
2124  tempPDFPtr->setExtrapolate( settings.flag("PDF:extrapolateLHAPDF") );
2125  }
2126  }
2127 
2128  // Pomeron beam, if not treated like a pi0 beam.
2129  else if (idIn == 990) {
2130  int pomSet = settings.mode("PDF:PomSet");
2131  double rescale = settings.parm("PDF:PomRescale");
2132 
2133  // A generic Q2-independent parametrization.
2134  if (pomSet == 1) {
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);
2143  }
2144 
2145  // The H1 Q2-dependent parametrizations. Initialization requires files.
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);
2152  }
2153 
2154  // Lepton beam: neutrino, resolved charged lepton or unresolved ditto.
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);
2159  }
2160 
2161  // Failure for unrecognized particle.
2162  else tempPDFPtr = 0;
2163 
2164  // Done.
2165  return tempPDFPtr;
2166 }
2167 
2168 //==========================================================================
2169 
2170 } // end namespace Pythia8
Definition: AgUStep.h:26