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) 2018 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, 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 #include "Pythia8/HeavyIons.h"
10 
11 // Access time information.
12 #include <ctime>
13 
14 // Allow string and character manipulation.
15 #include <cctype>
16 
17 namespace Pythia8 {
18 
19 //==========================================================================
20 
21 // The Pythia class.
22 
23 //--------------------------------------------------------------------------
24 
25 // The current Pythia (sub)version number, to agree with XML version.
26 const double Pythia::VERSIONNUMBERHEAD = PYTHIA_VERSION;
27 const double Pythia::VERSIONNUMBERCODE = 8.235;
28 
29 //--------------------------------------------------------------------------
30 
31 // Constants: could be changed here if desired, but normally should not.
32 // These are of technical nature, as described for each.
33 
34 // Maximum number of tries to produce parton level from given input.
35 const int Pythia::NTRY = 10;
36 
37 // Negative integer to denote that no subrun has been set.
38 const int Pythia::SUBRUNDEFAULT = -999;
39 
40 //--------------------------------------------------------------------------
41 
42 // Constructor.
43 
44 Pythia::Pythia(string xmlDir, bool printBanner) {
45 
46  // Initialise / reset pointers and global variables.
47  initPtrs();
48 
49  // Find path to data files, i.e. xmldoc directory location.
50  // Environment variable takes precedence, then constructor input,
51  // and finally the pre-processor constant XMLDIR.
52  xmlPath = "";
53  const char* PYTHIA8DATA = "PYTHIA8DATA";
54  char* envPath = getenv(PYTHIA8DATA);
55  if (envPath != 0 && *envPath != '\0') {
56  int i = 0;
57  while (*(envPath+i) != '\0') xmlPath += *(envPath+(i++));
58  } else {
59  if (xmlDir[ xmlDir.length() - 1 ] != '/') xmlDir += "/";
60  xmlPath = xmlDir;
61  ifstream xmlFile((xmlPath + "Index.xml").c_str());
62  if (!xmlFile.good()) xmlPath = XMLDIR;
63  xmlFile.close();
64  }
65  if (xmlPath[ xmlPath.length() - 1 ] != '/') xmlPath += "/";
66 
67  // Read in files with all flags, modes, parms and words.
68  settings.initPtr( &info);
69  string initFile = xmlPath + "Index.xml";
70  isConstructed = settings.init( initFile);
71  if (!isConstructed) {
72  info.errorMsg("Abort from Pythia::Pythia: settings unavailable");
73  return;
74  }
75 
76  // Save XML path in settings.
77  settings.addWord( "xmlPath", xmlPath);
78 
79  // Check that XML and header version numbers match code version number.
80  if (!checkVersion()) return;
81 
82  // Read in files with all particle data.
83  particleData.initPtr( &info, &settings, &rndm, couplingsPtr);
84  string dataFile = xmlPath + "ParticleData.xml";
85  isConstructed = particleData.init( dataFile);
86  if (!isConstructed) {
87  info.errorMsg("Abort from Pythia::Pythia: particle data unavailable");
88  return;
89  }
90 
91  // Write the Pythia banner to output.
92  if (printBanner) banner();
93 
94  // Not initialized until at the end of the init() call.
95  isInit = false;
96  info.addCounter(0);
97 
99 
100 }
101 
102 //--------------------------------------------------------------------------
103 
104 // Constructor from pre-initialised ParticleData and Settings objects.
105 
106 Pythia::Pythia(Settings& settingsIn, ParticleData& particleDataIn,
107  bool printBanner) {
108 
109  // Initialise / reset pointers and global variables.
110  initPtrs();
111 
112  // Copy XML path from existing Settings database.
113  xmlPath = settingsIn.word("xmlPath");
114 
115  // Copy settings database and redirect pointers.
116  settings = settingsIn;
117  settings.initPtr( &info);
118  isConstructed = settings.getIsInit();
119  if (!isConstructed) {
120  info.errorMsg("Abort from Pythia::Pythia: settings unavailable");
121  return;
122  }
123 
124  // Check XML and header version numbers match code version number.
125  if (!checkVersion()) return;
126 
127  // Copy particleData database and redirect pointers.
128  particleData = particleDataIn;
129  particleData.initPtr( &info, &settings, &rndm, couplingsPtr);
130  isConstructed = particleData.getIsInit();
131  if (!isConstructed) {
132  info.errorMsg("Abort from Pythia::Pythia: particle data unavailable");
133  return;
134  }
135 
136  // Write the Pythia banner to output.
137  if (printBanner) banner();
138 
139  // Not initialized until at the end of the init() call.
140  isInit = false;
141  info.addCounter(0);
142 
143 }
144 
145 //--------------------------------------------------------------------------
146 
147 // Constructor from string streams.
148 
149 Pythia::Pythia( istream& settingsStrings, istream& particleDataStrings,
150  bool printBanner) {
151 
152  // Initialise / reset pointers and global variables.
153  initPtrs();
154 
155  // Copy settings database
156  settings.init( settingsStrings );
157 
158  // Reset pointers to pertain to this PYTHIA object.
159  settings.initPtr( &info);
160  isConstructed = settings.getIsInit();
161  if (!isConstructed) {
162  info.errorMsg("Abort from Pythia::Pythia: settings unavailable");
163  return;
164  }
165 
166  // Check XML and header version numbers match code version number.
167  if (!checkVersion()) return;
168 
169  // Read in files with all particle data.
170  particleData.initPtr( &info, &settings, &rndm, couplingsPtr);
171  isConstructed = particleData.init( particleDataStrings );
172  if (!isConstructed) {
173  info.errorMsg("Abort from Pythia::Pythia: particle data unavailable");
174  return;
175  }
176 
177  // Write the Pythia banner to output.
178  if (printBanner) banner();
179 
180  // Not initialized until at the end of the init() call.
181  isInit = false;
182  info.addCounter(0);
183 
184 }
185 
186 //--------------------------------------------------------------------------
187 
188 // Destructor.
189 
190 Pythia::~Pythia() {
191 
192  // Delete the PDF's created with new.
193  if (useNewPdfHard && pdfHardAPtr != pdfAPtr) delete pdfHardAPtr;
194  if (useNewPdfHard && pdfHardBPtr != pdfBPtr) delete pdfHardBPtr;
195  if (useNewPdfA) delete pdfAPtr;
196  if (useNewPdfB) delete pdfBPtr;
197  if (useNewPdfPomA) delete pdfPomAPtr;
198  if (useNewPdfPomB) delete pdfPomBPtr;
199  if (useNewPdfGamA) delete pdfGamAPtr;
200  if (useNewPdfGamB) delete pdfGamBPtr;
201  if (useNewPdfUnresA) delete pdfUnresAPtr;
202  if (useNewPdfUnresB) delete pdfUnresBPtr;
203  if (useNewPdfUnresGamA) delete pdfUnresGamAPtr;
204  if (useNewPdfUnresGamB) delete pdfUnresGamBPtr;
205  if (useNewPdfVMDA) delete pdfVMDAPtr;
206  if (useNewPdfVMDB) delete pdfVMDBPtr;
207 
208  // Delete the Les Houches object created with new.
209  if (useNewLHA) delete lhaUpPtr;
210 
211  // Delete vector of UserHooks (but not the UserHooks themselves).
212  if (hasUserHooksVector) delete userHooksPtr;
213 
214  // Delete the Merging object created with new.
215  if (hasOwnMerging) delete mergingPtr;
216 
217  // Delete the MergingHooks object created with new.
218  if (hasOwnMergingHooks) delete mergingHooksPtr;
219 
220  // Delete the HeavyIons object created with new.
221  if (hasOwnHeavyIons) delete heavyIonsPtr;
222 
223  // Delete the BeamShape object created with new.
224  if (useNewBeamShape) delete beamShapePtr;
225 
226  // Delete the timelike and spacelike showers created with new.
227  if (useNewTimesDec) delete timesDecPtr;
228  if (useNewTimes && !useNewTimesDec) delete timesPtr;
229  if (useNewSpace) delete spacePtr;
230 
231  // Delete the parton vertex object created with new.
232  if (useNewPartonVertex) delete partonVertexPtr;
233 
234 }
235 
236 //--------------------------------------------------------------------------
237 
238 // Initialise new Pythia object (common code called by constructors).
239 
240 void Pythia::initPtrs() {
241 
242  // Initial values for pointers to PDF's.
243  useNewPdfA = false;
244  useNewPdfB = false;
245  useNewPdfHard = false;
246  useNewPdfPomA = false;
247  useNewPdfPomB = false;
248  useNewPdfGamA = false;
249  useNewPdfGamB = false;
250  useNewPdfHardGamA = false;
251  useNewPdfHardGamB = false;
252  useNewPdfUnresA = false;
253  useNewPdfUnresB = false;
254  useNewPdfUnresGamA = false;
255  useNewPdfUnresGamB = false;
256  useNewPdfVMDA = false;
257  useNewPdfVMDB = false;
258  pdfAPtr = 0;
259  pdfBPtr = 0;
260  pdfHardAPtr = 0;
261  pdfHardBPtr = 0;
262  pdfPomAPtr = 0;
263  pdfPomBPtr = 0;
264  pdfGamAPtr = 0;
265  pdfGamBPtr = 0;
266  pdfHardGamAPtr = 0;
267  pdfHardGamBPtr = 0;
268  pdfUnresAPtr = 0;
269  pdfUnresBPtr = 0;
270  pdfUnresGamAPtr = 0;
271  pdfUnresGamBPtr = 0;
272  pdfVMDAPtr = 0;
273  pdfVMDBPtr = 0;
274 
275  // Initial pointers to externally provided photon fluxes.
276  pdfGamFluxAPtr = 0;
277  pdfGamFluxBPtr = 0;
278 
279  // Initial values for pointers to Les Houches Event objects.
280  doLHA = false;
281  useNewLHA = false;
282  lhaUpPtr = 0;
283 
284  //Initial value for couplings pointer
285  couplingsPtr = &couplings;
286 
287  // Initial value for pointer to external decay handler.
288  decayHandlePtr = 0;
289 
290  // Initial value for pointer to user hooks.
291  userHooksPtr = 0;
292  hasUserHooksVector = false;
293 
294  // Initial value for pointer to merging hooks.
295  doMerging = false;
296  hasMerging = false;
297  hasOwnMerging = false;
298  mergingPtr = 0;
299  hasMergingHooks = false;
300  hasOwnMergingHooks = false;
301  mergingHooksPtr = 0;
302 
303  // Initial value for pointer to HeavyIons objects.
304  doHeavyIons = false;
305  hasHeavyIons = false;
306  hasOwnHeavyIons = false;
307  heavyIonsPtr = 0;
308  hiHooksPtr = 0;
309 
310  // Initial value for pointer to beam shape.
311  useNewBeamShape = false;
312  beamShapePtr = 0;
313 
314  // Initial values for pointers to timelike and spacelike showers.
315  useNewTimesDec = false;
316  useNewTimes = false;
317  useNewSpace = false;
318  timesDecPtr = 0;
319  timesPtr = 0;
320  spacePtr = 0;
321 
322  // Initial values for pointer parton-vertex setting.
323  useNewPartonVertex = false;
324  partonVertexPtr = 0;
325 
326 }
327 
328 //--------------------------------------------------------------------------
329 
330 // Check for consistency of version numbers (called by constructors).
331 
332 bool Pythia::checkVersion() {
333 
334  // Check that XML version number matches code version number.
335  double versionNumberXML = parm("Pythia:versionNumber");
336  isConstructed = (abs(versionNumberXML - VERSIONNUMBERCODE) < 0.0005);
337  if (!isConstructed) {
338  ostringstream errCode;
339  errCode << fixed << setprecision(3) << ": in code " << VERSIONNUMBERCODE
340  << " but in XML " << versionNumberXML;
341  info.errorMsg("Abort from Pythia::Pythia: unmatched version numbers",
342  errCode.str());
343  return false;
344  }
345 
346  // Check that header version number matches code version number.
347  isConstructed = (abs(VERSIONNUMBERHEAD - VERSIONNUMBERCODE) < 0.0005);
348  if (!isConstructed) {
349  ostringstream errCode;
350  errCode << fixed << setprecision(3) << ": in code " << VERSIONNUMBERCODE
351  << " but in header " << VERSIONNUMBERHEAD;
352  info.errorMsg("Abort from Pythia::Pythia: unmatched version numbers",
353  errCode.str());
354  return false;
355  }
356 
357  // All is well that ends well.
358  return true;
359 
360 }
361 
362 //--------------------------------------------------------------------------
363 
364 // Read in one update for a setting or particle data from a single line.
365 
366 bool Pythia::readString(string line, bool warn) {
367 
368  // Check that constructor worked.
369  if (!isConstructed) return false;
370 
371  // If empty line then done.
372  if (line.find_first_not_of(" \n\t\v\b\r\f\a") == string::npos) return true;
373 
374  // If Settings input stretches over several lines then continue with it.
375  if (settings.unfinishedInput()) return settings.readString(line, warn);
376 
377  // If first character is not a letter/digit, then taken to be a comment.
378  int firstChar = line.find_first_not_of(" \n\t\v\b\r\f\a");
379  if (!isalnum(line[firstChar])) return true;
380 
381  // Send on particle data to the ParticleData database.
382  if (isdigit(line[firstChar])) {
383  bool passed = particleData.readString(line, warn);
384  if (passed) particleDataBuffer << line << endl;
385  return passed;
386  }
387 
388  // Everything else sent on to Settings.
389  return settings.readString(line, warn);
390 
391 }
392 
393 //--------------------------------------------------------------------------
394 
395 // Read in updates for settings or particle data from user-defined file.
396 
397 bool Pythia::readFile(string fileName, bool warn, int subrun) {
398 
399  // Check that constructor worked.
400  if (!isConstructed) return false;
401 
402  // Open file for reading.
403  const char* cstring = fileName.c_str();
404  ifstream is(cstring);
405  if (!is.good()) {
406  info.errorMsg("Error in Pythia::readFile: did not find file", fileName);
407  return false;
408  }
409 
410  // Hand over real work to next method.
411  return readFile( is, warn, subrun);
412 
413 }
414 
415 //--------------------------------------------------------------------------
416 
417 // Read in updates for settings or particle data
418 // from user-defined stream (or file).
419 
420 bool Pythia::readFile(istream& is, bool warn, int subrun) {
421 
422  // Check that constructor worked.
423  if (!isConstructed) return false;
424 
425  // Read in one line at a time.
426  string line;
427  bool isCommented = false;
428  bool accepted = true;
429  int subrunNow = SUBRUNDEFAULT;
430  while ( getline(is, line) ) {
431 
432  // Check whether entering, leaving or inside commented-commands section.
433  int commentLine = readCommented( line);
434  if (commentLine == +1) isCommented = true;
435  else if (commentLine == -1) isCommented = false;
436  else if (isCommented) ;
437 
438  else {
439  // Check whether entered new subrun.
440  int subrunLine = readSubrun( line, warn);
441  if (subrunLine >= 0) subrunNow = subrunLine;
442 
443  // Process the line if in correct subrun.
444  if ( (subrunNow == subrun || subrunNow == SUBRUNDEFAULT)
445  && !readString( line, warn) ) accepted = false;
446  }
447 
448  // Reached end of input file.
449  };
450  return accepted;
451 
452 }
453 
454 //--------------------------------------------------------------------------
455 
456 // Routine to pass in pointers to PDF's. Usage optional.
457 
458 bool Pythia::setPDFPtr( PDF* pdfAPtrIn, PDF* pdfBPtrIn, PDF* pdfHardAPtrIn,
459  PDF* pdfHardBPtrIn, PDF* pdfPomAPtrIn, PDF* pdfPomBPtrIn,
460  PDF* pdfGamAPtrIn, PDF* pdfGamBPtrIn, PDF* pdfHardGamAPtrIn,
461  PDF* pdfHardGamBPtrIn, PDF* pdfUnresAPtrIn, PDF* pdfUnresBPtrIn,
462  PDF* pdfUnresGamAPtrIn, PDF* pdfUnresGamBPtrIn, PDF* pdfVMDAPtrIn,
463  PDF* pdfVMDBPtrIn) {
464 
465  // Delete any PDF's created in a previous init call.
466  if (useNewPdfHard && pdfHardAPtr != pdfAPtr) delete pdfHardAPtr;
467  if (useNewPdfHard && pdfHardBPtr != pdfBPtr) delete pdfHardBPtr;
468  if (useNewPdfA) delete pdfAPtr;
469  if (useNewPdfB) delete pdfBPtr;
470  if (useNewPdfPomA) delete pdfPomAPtr;
471  if (useNewPdfPomB) delete pdfPomBPtr;
472  if (useNewPdfGamA) delete pdfGamAPtr;
473  if (useNewPdfGamB) delete pdfGamBPtr;
474  if (useNewPdfUnresA) delete pdfUnresAPtr;
475  if (useNewPdfUnresB) delete pdfUnresBPtr;
476  if (useNewPdfUnresGamA) delete pdfUnresGamAPtr;
477  if (useNewPdfUnresGamB) delete pdfUnresGamBPtr;
478  if (useNewPdfHardGamA && pdfHardGamAPtr != pdfGamAPtr) delete pdfHardGamAPtr;
479  if (useNewPdfHardGamB && pdfHardGamBPtr != pdfGamBPtr) delete pdfHardGamBPtr;
480  if (useNewPdfVMDA) delete pdfVMDAPtr;
481  if (useNewPdfVMDB) delete pdfVMDBPtr;
482 
483  // Reset pointers to be empty.
484  useNewPdfA = false;
485  useNewPdfB = false;
486  useNewPdfHard = false;
487  useNewPdfPomA = false;
488  useNewPdfPomB = false;
489  useNewPdfGamA = false;
490  useNewPdfGamB = false;
491  useNewPdfHardGamA = false;
492  useNewPdfHardGamB = false;
493  useNewPdfUnresA = false;
494  useNewPdfUnresB = false;
495  useNewPdfUnresGamA = false;
496  useNewPdfUnresGamB = false;
497  useNewPdfVMDA = false;
498  useNewPdfVMDB = false;
499  pdfAPtr = 0;
500  pdfBPtr = 0;
501  pdfHardAPtr = 0;
502  pdfHardBPtr = 0;
503  pdfPomAPtr = 0;
504  pdfPomBPtr = 0;
505  pdfGamAPtr = 0;
506  pdfGamBPtr = 0;
507  pdfHardGamAPtr = 0;
508  pdfHardGamBPtr = 0;
509  pdfUnresAPtr = 0;
510  pdfUnresBPtr = 0;
511  pdfUnresGamAPtr = 0;
512  pdfUnresGamBPtr = 0;
513  pdfVMDAPtr = 0;
514  pdfVMDBPtr = 0;
515 
516  // Switch off external PDF's by zero as input.
517  if (pdfAPtrIn == 0 && pdfBPtrIn == 0) return true;
518 
519  // The two PDF objects cannot be one and the same.
520  if (pdfAPtrIn == pdfBPtrIn) return false;
521 
522  // Save pointers.
523  pdfAPtr = pdfAPtrIn;
524  pdfBPtr = pdfBPtrIn;
525 
526  // By default same pointers for hard-process PDF's.
527  pdfHardAPtr = pdfAPtrIn;
528  pdfHardBPtr = pdfBPtrIn;
529 
530  // Optionally allow separate pointers for hard process.
531  if (pdfHardAPtrIn != 0 && pdfHardBPtrIn != 0) {
532  if (pdfHardAPtrIn == pdfHardBPtrIn) return false;
533  pdfHardAPtr = pdfHardAPtrIn;
534  pdfHardBPtr = pdfHardBPtrIn;
535  }
536 
537  // Optionally allow pointers for Pomerons in the proton.
538  if (pdfPomAPtrIn != 0 && pdfPomBPtrIn != 0) {
539  if (pdfPomAPtrIn == pdfPomBPtrIn) return false;
540  pdfPomAPtr = pdfPomAPtrIn;
541  pdfPomBPtr = pdfPomBPtrIn;
542  }
543 
544  // Optionally allow pointers for Gammas in the leptons.
545  if (pdfGamAPtrIn != 0 && pdfGamBPtrIn != 0) {
546  if (pdfGamAPtrIn == pdfGamBPtrIn) return false;
547  pdfGamAPtr = pdfGamAPtrIn;
548  pdfGamBPtr = pdfGamBPtrIn;
549  }
550 
551  // Optionally allow pointers for Hard PDFs for photons in the leptons.
552  if (pdfHardGamAPtrIn != 0 && pdfHardGamBPtrIn != 0) {
553  if (pdfHardGamAPtrIn == pdfHardGamBPtrIn) return false;
554  pdfHardGamAPtr = pdfHardGamAPtrIn;
555  pdfHardGamBPtr = pdfHardGamBPtrIn;
556  }
557 
558  // Optionally allow pointers for unresolved PDFs.
559  if (pdfUnresAPtrIn != 0 && pdfUnresBPtrIn != 0) {
560  if (pdfUnresAPtrIn == pdfUnresBPtrIn) return false;
561  pdfUnresAPtr = pdfUnresAPtrIn;
562  pdfUnresBPtr = pdfUnresBPtrIn;
563  }
564 
565  // Optionally allow pointers for unresolved PDFs for photons from leptons.
566  if (pdfUnresGamAPtrIn != 0 && pdfUnresGamBPtrIn != 0) {
567  if (pdfUnresGamAPtrIn == pdfUnresGamBPtrIn) return false;
568  pdfUnresGamAPtr = pdfUnresGamAPtrIn;
569  pdfUnresGamBPtr = pdfUnresGamBPtrIn;
570  }
571 
572  // Optionally allow pointers for VMD in the gamma.
573  if (pdfVMDAPtrIn != 0 && pdfVMDBPtrIn != 0) {
574  if (pdfVMDAPtrIn == pdfVMDBPtrIn) return false;
575  pdfVMDAPtr = pdfVMDAPtrIn;
576  pdfVMDBPtr = pdfVMDBPtrIn;
577  }
578 
579  // Done.
580  return true;
581 }
582 
583 //--------------------------------------------------------------------------
584 
585 // Routine to initialize with the variable values of the Beams kind.
586 
587 bool Pythia::init() {
588 
589  // Check that constructor worked.
590  isInit = false;
591  if (!isConstructed) {
592  info.errorMsg("Abort from Pythia::init: constructor "
593  "initialization failed");
594  return false;
595  }
596 
597  // Early catching of heavy ion mode.
598  doHeavyIons = HeavyIons::isHeavyIon(settings) ||
599  settings.mode("HeavyIon:mode") == 2;
600  if ( doHeavyIons ) {
601  if ( !heavyIonsPtr ) {
602  heavyIonsPtr = new Angantyr(*this);
603  hasOwnHeavyIons = true;
604  }
605  if ( hiHooksPtr ) heavyIonsPtr->setHIUserHooksPtr(hiHooksPtr);
606  if ( !heavyIonsPtr->init() ) doHeavyIons = false;
607  }
608 
609  // Early readout, if return false or changed when no beams.
610  doProcessLevel = settings.flag("ProcessLevel:all");
611 
612  // Check that changes in Settings and ParticleData have worked.
613  if (settings.unfinishedInput()) {
614  info.errorMsg("Abort from Pythia::init: opening { not matched by "
615  "closing }");
616  return false;
617  }
618  if (settings.readingFailed()) {
619  info.errorMsg("Abort from Pythia::init: some user settings "
620  "did not make sense");
621  return false;
622  }
623  if (particleData.readingFailed()) {
624  info.errorMsg("Abort from Pythia::init: some user particle data "
625  "did not make sense");
626  return false;
627  }
628 
629  // Initialize the random number generator.
630  if ( settings.flag("Random:setSeed") )
631  rndm.init( settings.mode("Random:seed") );
632 
633  // Find which frame type to use.
634  info.addCounter(1);
635  frameType = mode("Beams:frameType");
636 
637  // Set up values related to CKKW-L merging.
638  bool doUserMerging = settings.flag("Merging:doUserMerging");
639  bool doMGMerging = settings.flag("Merging:doMGMerging");
640  bool doKTMerging = settings.flag("Merging:doKTMerging");
641  bool doPTLundMerging = settings.flag("Merging:doPTLundMerging");
642  bool doCutBasedMerging = settings.flag("Merging:doCutBasedMerging");
643  // Set up values related to unitarised CKKW merging
644  bool doUMEPSTree = settings.flag("Merging:doUMEPSTree");
645  bool doUMEPSSubt = settings.flag("Merging:doUMEPSSubt");
646  // Set up values related to NL3 NLO merging
647  bool doNL3Tree = settings.flag("Merging:doNL3Tree");
648  bool doNL3Loop = settings.flag("Merging:doNL3Loop");
649  bool doNL3Subt = settings.flag("Merging:doNL3Subt");
650  // Set up values related to unitarised NLO merging
651  bool doUNLOPSTree = settings.flag("Merging:doUNLOPSTree");
652  bool doUNLOPSLoop = settings.flag("Merging:doUNLOPSLoop");
653  bool doUNLOPSSubt = settings.flag("Merging:doUNLOPSSubt");
654  bool doUNLOPSSubtNLO = settings.flag("Merging:doUNLOPSSubtNLO");
655  bool doXSectionEst = settings.flag("Merging:doXSectionEstimate");
656  doMerging = doUserMerging || doMGMerging || doKTMerging
657  || doPTLundMerging || doCutBasedMerging || doUMEPSTree || doUMEPSSubt
658  || doNL3Tree || doNL3Loop || doNL3Subt || doUNLOPSTree
659  || doUNLOPSLoop || doUNLOPSSubt || doUNLOPSSubtNLO || doXSectionEst;
660 
661  // Set up MergingHooks object.
662  bool inputMergingHooks = (mergingHooksPtr != 0);
663  if (doMerging && !inputMergingHooks) {
664  if (hasOwnMergingHooks && mergingHooksPtr) delete mergingHooksPtr;
665  mergingHooksPtr = new MergingHooks();
666  hasOwnMergingHooks = true;
667  }
668 
669  hasMergingHooks = (mergingHooksPtr != 0);
670  // Merging hooks required for merging. If no merging hooks pointer is
671  // available, exit.
672  if (doMerging && !hasMergingHooks) {
673  info.errorMsg("Abort from Pythia::init: "
674  "no merging hooks object has been provided");
675  return false;
676  } else if (doMerging) {
677  mergingHooksPtr->setLHEInputFile("");
678  }
679 
680  // Set up Merging object.
681  bool inputMerging = (mergingPtr != 0);
682  if (doMerging && !inputMerging) {
683  if (hasOwnMerging && mergingPtr) delete mergingPtr;
684  mergingPtr = new Merging();
685  hasOwnMerging = true;
686  }
687  hasMerging = (mergingPtr != 0);
688 
689  // Initialization with internal processes: read in and set values.
690  if (frameType < 4 ) {
691  doLHA = false;
692  boostType = frameType;
693  idA = mode("Beams:idA");
694  idB = mode("Beams:idB");
695  eCM = parm("Beams:eCM");
696  eA = parm("Beams:eA");
697  eB = parm("Beams:eB");
698  pxA = parm("Beams:pxA");
699  pyA = parm("Beams:pyA");
700  pzA = parm("Beams:pzA");
701  pxB = parm("Beams:pxB");
702  pyB = parm("Beams:pyB");
703  pzB = parm("Beams:pzB");
704 
705  // Initialization with a Les Houches Event File or an LHAup object.
706  } else {
707  doLHA = true;
708  boostType = 2;
709  string lhef = word("Beams:LHEF");
710  string lhefHeader = word("Beams:LHEFheader");
711  bool readHeaders = flag("Beams:readLHEFheaders");
712  bool setScales = flag("Beams:setProductionScalesFromLHEF");
713  bool skipInit = flag("Beams:newLHEFsameInit");
714  int nSkipAtInit = mode("Beams:nSkipLHEFatInit");
715 
716  // For file input: renew file stream or (re)new Les Houches object.
717  if (frameType == 4) {
718  const char* cstring1 = lhef.c_str();
719  bool useExternal = (lhaUpPtr && !useNewLHA && lhaUpPtr->useExternal());
720  if (!useExternal && useNewLHA && skipInit)
721  lhaUpPtr->newEventFile(cstring1);
722  else if (!useExternal) {
723  if (useNewLHA) delete lhaUpPtr;
724  // Header is optional, so use NULL pointer to indicate no value.
725  const char* cstring2 = (lhefHeader == "void")
726  ? NULL : lhefHeader.c_str();
727  lhaUpPtr = new LHAupLHEF(&info, cstring1, cstring2, readHeaders,
728  setScales);
729  useNewLHA = true;
730  }
731 
732  // Check that file was properly opened.
733  if (!lhaUpPtr->fileFound()) {
734  info.errorMsg("Abort from Pythia::init: "
735  "Les Houches Event File not found");
736  return false;
737  }
738 
739  // For object input: at least check that not null pointer.
740  } else {
741  if (lhaUpPtr == 0) {
742  info.errorMsg("Abort from Pythia::init: "
743  "LHAup object not found");
744  return false;
745  }
746 
747  // LHAup object generic abort using fileFound() routine.
748  if (!lhaUpPtr->fileFound()) {
749  info.errorMsg("Abort from Pythia::init: "
750  "LHAup initialisation error");
751  return false;
752  }
753  }
754 
755  // Send in pointer to info. Store or replace LHA pointer in other classes.
756  lhaUpPtr->setPtr( &info);
757  processLevel.setLHAPtr( lhaUpPtr);
758 
759  // If second time around, only with new file, then simplify.
760  // Optionally skip ahead a number of events at beginning of file.
761  if (skipInit) {
762  isInit = true;
763  info.addCounter(2);
764  if (nSkipAtInit > 0) lhaUpPtr->skipEvent(nSkipAtInit);
765  return true;
766  }
767 
768  // Set up values related to merging hooks.
769  if (frameType == 4 || frameType == 5) {
770 
771  // Store the name of the input LHEF for merging.
772  if (doMerging) {
773  string lhefIn = (frameType == 4) ? lhef : "";
774  mergingHooksPtr->setLHEInputFile( lhefIn);
775  }
776 
777  }
778 
779  // Set LHAinit information (in some external program).
780  if ( !lhaUpPtr->setInit()) {
781  info.errorMsg("Abort from Pythia::init: "
782  "Les Houches initialization failed");
783  return false;
784  }
785 
786  // Extract beams from values set in an LHAinit object.
787  idA = lhaUpPtr->idBeamA();
788  idB = lhaUpPtr->idBeamB();
789  int idRenameBeams = settings.mode("LesHouches:idRenameBeams");
790  if (abs(idA) == idRenameBeams) idA = 16;
791  if (abs(idB) == idRenameBeams) idB = -16;
792  if (idA == 0 || idB == 0) doProcessLevel = false;
793  eA = lhaUpPtr->eBeamA();
794  eB = lhaUpPtr->eBeamB();
795 
796  // Optionally skip ahead a number of events at beginning of file.
797  if (nSkipAtInit > 0) lhaUpPtr->skipEvent(nSkipAtInit);
798  }
799 
800  // Set up values related to user hooks.
801  hasUserHooks = (userHooksPtr != 0);
802  doVetoProcess = false;
803  doVetoPartons = false;
804  retryPartonLevel = false;
805  if (hasUserHooks) {
806  userHooksPtr->initPtr( &info, &settings, &particleData, &rndm, &beamA,
807  &beamB, &beamPomA, &beamPomB, couplingsPtr, &partonSystems, &sigmaTot);
808  if (!userHooksPtr->initAfterBeams()) {
809  info.errorMsg("Abort from Pythia::init: could not initialise UserHooks");
810  return false;
811  }
812  doVetoProcess = userHooksPtr->canVetoProcessLevel();
813  doVetoPartons = userHooksPtr->canVetoPartonLevel();
814  retryPartonLevel = userHooksPtr->retryPartonLevel();
815  }
816 
817  // Back to common initialization. Reset error counters.
818  nErrEvent = 0;
819  info.setTooLowPTmin(false);
820  info.sigmaReset();
821 
822  // Initialize data members extracted from database.
823  doPartonLevel = settings.flag("PartonLevel:all");
824  doHadronLevel = settings.flag("HadronLevel:all");
825  doCentralDiff = settings.flag("SoftQCD:centralDiffractive");
826  doSoftQCDall = settings.flag("SoftQCD:all");
827  doSoftQCDinel = settings.flag("SoftQCD:inelastic");
828  doDiffraction = settings.flag("SoftQCD:singleDiffractive")
829  || settings.flag("SoftQCD:doubleDiffractive")
830  || doSoftQCDall || doSoftQCDinel || doCentralDiff;
831  doSoftQCD = doDiffraction ||
832  settings.flag("SoftQCD:nonDiffractive");
833  doHardDiff = settings.flag("Diffraction:doHard");
834  doResDec = settings.flag("ProcessLevel:resonanceDecays");
835  doFSRinRes = doPartonLevel && settings.flag("PartonLevel:FSR")
836  && settings.flag("PartonLevel:FSRinResonances");
837  decayRHadrons = settings.flag("RHadrons:allowDecay");
838  doMomentumSpread = settings.flag("Beams:allowMomentumSpread");
839  doVertexSpread = settings.flag("Beams:allowVertexSpread");
840  abortIfVeto = settings.flag("Check:abortIfVeto");
841  checkEvent = settings.flag("Check:event");
842  checkHistory = settings.flag("Check:history");
843  nErrList = settings.mode("Check:nErrList");
844  epTolErr = settings.parm("Check:epTolErr");
845  epTolWarn = settings.parm("Check:epTolWarn");
846  mTolErr = settings.parm("Check:mTolErr");
847  mTolWarn = settings.parm("Check:mTolWarn");
848 
849  // Find out if beam is or has a resolved photon beam.
850  beamHasGamma = settings.flag("PDF:lepton2gamma");
851  gammaMode = settings.mode("Photon:ProcessType");
852  bool beamAneedResGamma = (gammaMode == 1) || (gammaMode == 2)
853  || (gammaMode == 0);
854  bool beamBneedResGamma = (gammaMode == 1) || (gammaMode == 3)
855  || (gammaMode == 0);
856  beamAisResGamma = beamAneedResGamma && idA == 22;
857  beamBisResGamma = beamBneedResGamma && idB == 22;
858  bool isChargedLeptonA = (abs(idA) == 11 || abs(idA) == 13 || abs(idA) == 15);
859  bool isChargedLeptonB = (abs(idB) == 11 || abs(idB) == 13 || abs(idB) == 15);
860  beamAhasResGamma = beamAneedResGamma && beamHasGamma && isChargedLeptonA;
861  beamBhasResGamma = beamBneedResGamma && beamHasGamma && isChargedLeptonB;
862  doVMDsideA = doSoftQCD && (beamAisResGamma || beamAhasResGamma);
863  doVMDsideB = doSoftQCD && (beamBisResGamma || beamBhasResGamma);
864 
865  // Turn off central diffraction for VMD processes.
866  if (doVMDsideA || doVMDsideB) {
867  if (doCentralDiff) {
868  info.errorMsg("Warning in Pythia::init: "
869  "Central diffractive events not implemented for gamma + p/gamma");
870  return false;
871  }
872  if (doSoftQCDall) {
873  info.errorMsg("Warning in Pythia::init: "
874  "Central diffractive events not implemented for gamma + p/gamma");
875  settings.flag("SoftQCD:all", false);
876  settings.flag("SoftQCD:elastic", true);
877  settings.flag("SoftQCD:nonDiffractive", true);
878  settings.flag("SoftQCD:singleDiffractive", true);
879  settings.flag("SoftQCD:doubleDiffractive", true);
880  }
881  if (doSoftQCDinel) {
882  info.errorMsg("Warning in Pythia::init: "
883  "Central diffractive events not implemented for gamma + p/gamma");
884  settings.flag("SoftQCD:inelastic", false);
885  settings.flag("SoftQCD:nonDiffractive", true);
886  settings.flag("SoftQCD:singleDiffractive", true);
887  settings.flag("SoftQCD:doubleDiffractive", true);
888  }
889  }
890 
891  // Initialise merging hooks.
892  if ( doMerging && (hasMergingHooks || hasOwnMergingHooks) ) {
893  mergingHooksPtr->initPtr( &settings, &info, &particleData, &partonSystems);
894  mergingHooksPtr->init();
895  }
896 
897  // Check that combinations of settings are allowed; change if not.
898  checkSettings();
899 
900  // Initialize the SM couplings (needed to initialize resonances).
901  couplingsPtr->init( settings, &rndm );
902 
903  // Initialize SLHA interface (including SLHA/BSM couplings).
904  bool useSLHAcouplings = false;
905  slhaInterface = SLHAinterface();
906  slhaInterface.setPtr( &info );
907  slhaInterface.init( settings, &rndm, couplingsPtr, &particleData,
908  useSLHAcouplings, particleDataBuffer );
909  if (useSLHAcouplings) couplingsPtr = slhaInterface.couplingsPtr;
910 
911  // Reset couplingsPtr to the correct memory address.
912  particleData.initPtr( &info, &settings, &rndm, couplingsPtr);
913  if (hasUserHooks) userHooksPtr->initPtr( &info, &settings, &particleData,
914  &rndm, &beamA, &beamB, &beamPomA, &beamPomB, couplingsPtr,
915  &partonSystems, &sigmaTot);
916 
917  // Set headers to distinguish the two event listing kinds.
918  int startColTag = settings.mode("Event:startColTag");
919  process.init("(hard process)", &particleData, startColTag);
920  event.init("(complete event)", &particleData, startColTag);
921 
922  // Final setup stage of particle data, notably resonance widths.
923  particleData.initWidths( resonancePtrs);
924 
925  // Set up R-hadrons particle data, where relevant.
926  rHadrons.init( &info, settings, &particleData, &rndm);
927 
928  // Set up and initialize setting of parton production vertices.
929  if (settings.flag("PartonVertex:setVertex")) {
930  if (partonVertexPtr == 0) {
931  partonVertexPtr = new PartonVertex();
932  useNewPartonVertex = true;
933  }
934  partonVertexPtr->initPtr( &info, &settings, &rndm);
935  partonVertexPtr->init();
936  }
937 
938  // Set up objects for timelike and spacelike showers.
939  if (timesDecPtr == 0 || timesPtr == 0) {
940  TimeShower* timesNow = new TimeShower();
941  if (timesDecPtr == 0) {
942  timesDecPtr = timesNow;
943  useNewTimesDec = true;
944  }
945  if (timesPtr == 0) {
946  timesPtr = timesNow;
947  useNewTimes = true;
948  }
949  }
950  if (spacePtr == 0) {
951  spacePtr = new SpaceShower();
952  useNewSpace = true;
953  }
954 
955  // Initialize pointers in showers.
956  timesPtr->initPtr( &info, &settings, &particleData, &rndm, couplingsPtr,
957  &partonSystems, userHooksPtr, mergingHooksPtr, partonVertexPtr);
958  timesDecPtr->initPtr( &info, &settings, &particleData, &rndm, couplingsPtr,
959  &partonSystems, userHooksPtr, mergingHooksPtr, partonVertexPtr);
960  spacePtr->initPtr( &info, &settings, &particleData, &rndm, couplingsPtr,
961  &partonSystems, userHooksPtr, mergingHooksPtr, partonVertexPtr);
962 
963  // Set up values related to beam shape.
964  if (beamShapePtr == 0) {
965  beamShapePtr = new BeamShape();
966  useNewBeamShape = true;
967  }
968  beamShapePtr->init( settings, &rndm);
969 
970  // Check that beams and beam combination can be handled.
971  if (!checkBeams()) {
972  info.errorMsg("Abort from Pythia::init: "
973  "checkBeams initialization failed");
974  return false;
975  }
976 
977  // Do not set up beam kinematics when no process level.
978  if (!doProcessLevel) boostType = 1;
979  else {
980 
981  // Set up beam kinematics.
982  if (!initKinematics()) {
983  info.errorMsg("Abort from Pythia::init: "
984  "kinematics initialization failed");
985  return false;
986  }
987 
988  // Set up pointers to PDFs.
989  if (!initPDFs()) {
990  info.errorMsg("Abort from Pythia::init: PDF initialization failed");
991  return false;
992  }
993 
994  // Set up the two beams and the common remnant system.
995  StringFlav* flavSelPtr = hadronLevel.getStringFlavPtr();
996  beamA.init( idA, pzAcm, eA, mA, &info, settings, &particleData, &rndm,
997  pdfAPtr, pdfHardAPtr, isUnresolvedA, flavSelPtr);
998  beamB.init( idB, pzBcm, eB, mB, &info, settings, &particleData, &rndm,
999  pdfBPtr, pdfHardBPtr, isUnresolvedB, flavSelPtr);
1000 
1001  // Init also unresolved PDF pointers for photon beams when needed.
1002  if ( ( beamA.isGamma() || beamAhasResGamma )
1003  && ( gammaMode == 0 || gammaMode == 3 || gammaMode == 4 ) )
1004  beamA.initUnres( pdfUnresAPtr);
1005  if ( ( beamB.isGamma() || beamBhasResGamma )
1006  && ( gammaMode == 0 || gammaMode == 2 || gammaMode == 4 ) )
1007  beamB.initUnres( pdfUnresBPtr);
1008 
1009  // Optionally set up new alternative beams for these Pomerons.
1010  if ( doDiffraction || doHardDiff ) {
1011  beamPomA.init( 990, 0.5 * eCM, 0.5 * eCM, 0., &info, settings,
1012  &particleData, &rndm, pdfPomAPtr, pdfPomAPtr, false, flavSelPtr);
1013  beamPomB.init( 990, -0.5 * eCM, 0.5 * eCM, 0., &info, settings,
1014  &particleData, &rndm, pdfPomBPtr, pdfPomBPtr, false, flavSelPtr);
1015  }
1016 
1017  // Initialise VMD beams from gammas (in leptons). Use pion PDF for VMDs.
1018  if (doVMDsideA)
1019  beamVMDA.init( 111, 0.5 * eCM, 0.5 * eCM, 0., &info, settings,
1020  &particleData, &rndm, pdfVMDAPtr, pdfVMDAPtr, false, flavSelPtr);
1021  if (doVMDsideB)
1022  beamVMDB.init( 111, 0.5 * eCM, 0.5 * eCM, 0., &info, settings,
1023  &particleData, &rndm, pdfVMDBPtr, pdfVMDBPtr, false, flavSelPtr);
1024 
1025  // Optionally set up photon beams from lepton beams if resolved photons.
1026  if (beamHasGamma && gammaMode < 4) {
1027  if ( !(beamA.isHadron()) )
1028  beamGamA.init( 22, 0.5 * eCM, 0.5 * eCM, 0., &info, settings,
1029  &particleData, &rndm, pdfGamAPtr, pdfHardGamAPtr, !beamAisResGamma,
1030  flavSelPtr);
1031  if ( !(beamB.isHadron()) )
1032  beamGamB.init( 22, 0.5 * eCM, 0.5 * eCM, 0., &info, settings,
1033  &particleData, &rndm, pdfGamBPtr, pdfHardGamBPtr, !beamBisResGamma,
1034  flavSelPtr);
1035 
1036  // Initialize also unresolved PDFs for relevant processes.
1037  if ( gammaMode == 0 || gammaMode == 3 )
1038  beamGamA.initUnres( pdfUnresGamAPtr);
1039  if ( gammaMode == 0 || gammaMode == 2 )
1040  beamGamB.initUnres( pdfUnresGamBPtr);
1041  }
1042 
1043  }
1044 
1045  // Send info/pointers to process level for initialization.
1046  if ( doProcessLevel && !processLevel.init( &info, settings, &particleData,
1047  &rndm, &beamA, &beamB, &beamGamA, &beamGamB, &beamVMDA, &beamVMDB,
1048  couplingsPtr, &sigmaTot, doLHA, &slhaInterface, userHooksPtr,
1049  sigmaPtrs, phaseSpacePtrs) ) {
1050  info.errorMsg("Abort from Pythia::init: "
1051  "processLevel initialization failed");
1052  return false;
1053  }
1054 
1055  // Initialize timelike showers already here, since needed in decays.
1056  // The pointers to the beams are needed by some external plugin showers.
1057  timesDecPtr->init( &beamA, &beamB);
1058 
1059  // Alternatively only initialize resonance decays.
1060  if ( !doProcessLevel) processLevel.initDecays( &info, settings,
1061  &particleData, &rndm, lhaUpPtr);
1062 
1063  // Send info/pointers to parton level for initialization.
1064  if ( doPartonLevel && doProcessLevel && !partonLevel.init( &info, settings,
1065  &particleData, &rndm, &beamA, &beamB, &beamPomA, &beamPomB, &beamGamA,
1066  &beamGamB, &beamVMDA, &beamVMDB, couplingsPtr, &partonSystems, &sigmaTot,
1067  timesDecPtr, timesPtr, spacePtr, &rHadrons, userHooksPtr, mergingHooksPtr,
1068  partonVertexPtr, false) ) {
1069  info.errorMsg("Abort from Pythia::init: "
1070  "partonLevel initialization failed" );
1071  return false;
1072  }
1073 
1074  // Make pointer to shower available for merging machinery.
1075  if ( doMerging && (hasMergingHooks || hasOwnMergingHooks) )
1076  mergingHooksPtr->setShowerPointer(&partonLevel);
1077 
1078  // Alternatively only initialize final-state showers in resonance decays.
1079  if ( !doProcessLevel || !doPartonLevel) partonLevel.init( &info, settings,
1080  &particleData, &rndm, 0, 0, 0, 0, 0, 0, 0, 0, couplingsPtr, &partonSystems,
1081  0, timesDecPtr, 0, 0, &rHadrons, 0, 0, partonVertexPtr, false);
1082 
1083  // Send info/pointers to parton level for trial shower initialization.
1084  if ( doMerging && !trialPartonLevel.init( &info, settings, &particleData,
1085  &rndm, &beamA, &beamB, &beamPomA, &beamPomB, &beamGamA, &beamGamB,
1086  &beamVMDA, &beamVMDB, couplingsPtr, &partonSystems, &sigmaTot,
1087  timesDecPtr, timesPtr, spacePtr, &rHadrons, userHooksPtr,
1088  mergingHooksPtr, partonVertexPtr, true) ) {
1089  info.errorMsg("Abort from Pythia::init: "
1090  "trialPartonLevel initialization failed");
1091  return false;
1092  }
1093 
1094  // Initialise the merging wrapper class.
1095  if (doMerging ) {
1096  mergingPtr->initPtr( &settings, &info, &particleData, &rndm,
1097  &beamA, &beamB, mergingHooksPtr, &trialPartonLevel, couplingsPtr );
1098  mergingPtr->init();
1099  }
1100 
1101  // Send info/pointers to hadron level for initialization.
1102  // Note: forceHadronLevel() can come, so we must always initialize.
1103  if ( !hadronLevel.init( &info, settings, &particleData, &rndm,
1104  couplingsPtr, timesDecPtr, &rHadrons, decayHandlePtr,
1105  handledParticles, userHooksPtr) ) {
1106  info.errorMsg("Abort from Pythia::init: "
1107  "hadronLevel initialization failed");
1108  return false;
1109  }
1110 
1111  // Optionally check particle data table for inconsistencies.
1112  if ( settings.flag("Check:particleData") )
1113  particleData.checkTable( settings.mode("Check:levelParticleData") );
1114 
1115  // Optionally show settings and particle data, changed or all.
1116  bool showCS = settings.flag("Init:showChangedSettings");
1117  bool showAS = settings.flag("Init:showAllSettings");
1118  bool showCPD = settings.flag("Init:showChangedParticleData");
1119  bool showCRD = settings.flag("Init:showChangedResonanceData");
1120  bool showAPD = settings.flag("Init:showAllParticleData");
1121  int show1PD = settings.mode("Init:showOneParticleData");
1122  bool showPro = settings.flag("Init:showProcesses");
1123  if (showCS) settings.listChanged();
1124  if (showAS) settings.listAll();
1125  if (show1PD > 0) particleData.list(show1PD);
1126  if (showCPD) particleData.listChanged(showCRD);
1127  if (showAPD) particleData.listAll();
1128 
1129  // Listing options for the next() routine.
1130  nCount = settings.mode("Next:numberCount");
1131  nShowLHA = settings.mode("Next:numberShowLHA");
1132  nShowInfo = settings.mode("Next:numberShowInfo");
1133  nShowProc = settings.mode("Next:numberShowProcess");
1134  nShowEvt = settings.mode("Next:numberShowEvent");
1135  showSaV = settings.flag("Next:showScaleAndVertex");
1136  showMaD = settings.flag("Next:showMothersAndDaughters");
1137 
1138  // Init colour reconnection and junction splitting.
1139  colourReconnection.init( &info, settings, &rndm, &particleData,
1140  &beamA, &beamB, &partonSystems);
1141  junctionSplitting.init(&info, settings, &rndm, &particleData);
1142 
1143  // Flags for colour reconnection.
1144  doReconnect = settings.flag("ColourReconnection:reconnect");
1145  reconnectMode = settings.mode("ColourReconnection:mode");
1146  forceHadronLevelCR = settings.flag("ColourReconnection:forceHadronLevelCR");
1147 
1148  // Succeeded.
1149  isInit = true;
1150  info.addCounter(2);
1151  if (useNewLHA && showPro) lhaUpPtr->listInit();
1152  return true;
1153 
1154 }
1155 
1156 //--------------------------------------------------------------------------
1157 
1158 // Check that combinations of settings are allowed; change if not.
1159 
1160 void Pythia::checkSettings() {
1161 
1162  // Double rescattering not allowed if ISR or FSR.
1163  if ((settings.flag("PartonLevel:ISR") || settings.flag("PartonLevel:FSR"))
1164  && settings.flag("MultipartonInteractions:allowDoubleRescatter")) {
1165  info.errorMsg("Warning in Pythia::checkSettings: "
1166  "double rescattering switched off since showering is on");
1167  settings.flag("MultipartonInteractions:allowDoubleRescatter", false);
1168  }
1169 
1170  // Collisions with direct photon(s).
1171  if ( (idA == 22 && !beamAisResGamma) || (idB == 22 && !beamBisResGamma) ) {
1172 
1173  // Turn MPIs off.
1174  if ( settings.flag("PartonLevel:MPI") ) {
1175  info.errorMsg("Warning in Pythia::checkSettings: "
1176  "MPIs turned off for collision with unresolved photon");
1177  settings.flag("PartonLevel:MPI", false);
1178  }
1179 
1180  // Check that no soft QCD processes initialized with direct photons.
1181  if ( settings.flag("SoftQCD:nonDiffractive") ) {
1182  info.errorMsg("Warning in Pythia::checkSettings: "
1183  "Soft QCD processes turned off for collision with unresolved photon");
1184  settings.flag("SoftQCD:nonDiffractive", false);
1185  }
1186 
1187  }
1188 
1189  // Lepton-lepton/hadron collisions with direct photon(s).
1190  if ( ( (abs(idA) > 10 && abs(idA) < 17)
1191  && !beamAhasResGamma && beamHasGamma)
1192  || ( (abs(idB) > 10 && abs(idB) < 17)
1193  && !beamBhasResGamma && beamHasGamma) ) {
1194 
1195  // Turn MPIs off.
1196  if ( settings.flag("PartonLevel:MPI") ) {
1197  info.errorMsg("Warning in Pythia::checkSettings: MPIs turned off for "
1198  "collision with unresolved photon");
1199  settings.flag("PartonLevel:MPI", false);
1200  }
1201 
1202  // Check that no soft QCD processes initialized with direct photons.
1203  if ( settings.flag("SoftQCD:nonDiffractive") ) {
1204  info.errorMsg("Warning in Pythia::checkSettings: "
1205  "Soft QCD processes turned off for collision with unresolved photon");
1206  settings.flag("SoftQCD:nonDiffractive", false);
1207  }
1208 
1209  }
1210 
1211 }
1212 
1213 //--------------------------------------------------------------------------
1214 
1215 // Check that beams and beam combination can be handled. Set up unresolved.
1216 
1217 bool Pythia::checkBeams() {
1218 
1219  // Absolute flavours. If not to do process level then no check needed.
1220  int idAabs = abs(idA);
1221  int idBabs = abs(idB);
1222  if (!doProcessLevel) return true;
1223 
1224  // Neutrino beams always unresolved, charged lepton ones conditionally.
1225  bool isLeptonA = (idAabs > 10 && idAabs < 17);
1226  bool isLeptonB = (idBabs > 10 && idBabs < 17);
1227  bool isUnresLep = !settings.flag("PDF:lepton");
1228  bool isGammaA = idAabs == 22;
1229  bool isGammaB = idBabs == 22;
1230  isUnresolvedA = ( isLeptonA && (idAabs%2 == 0 || isUnresLep) );
1231  isUnresolvedB = ( isLeptonB && (idBabs%2 == 0 || isUnresLep) );
1232 
1233  // Also photons may be unresolved.
1234  if ( idAabs == 22 && !beamAisResGamma ) isUnresolvedA = true;
1235  if ( idBabs == 22 && !beamBisResGamma ) isUnresolvedB = true;
1236 
1237  // If photons from beam particle, beam not unresolved.
1238  if ( beamAhasResGamma ) isUnresolvedA = false;
1239  if ( beamBhasResGamma ) isUnresolvedB = false;
1240 
1241  // Equate Dark Matter "beams" with incoming neutrinos.
1242  if (idAabs > 50 && idAabs < 61) isLeptonA = isUnresolvedA = true;
1243  if (idBabs > 50 && idBabs < 61) isLeptonB = isUnresolvedB = true;
1244 
1245  // Lepton-lepton collisions.
1246  if (isLeptonA && isLeptonB ) {
1247 
1248  // Photon-photon collision from lepton beams.
1249  if (beamHasGamma) {
1250 
1251  // Non-diffractive events only for resolved photon-photon.
1252  if ( (!beamAhasResGamma || !beamBhasResGamma)
1253  && settings.flag("SoftQCD:nonDiffractive") ) {
1254  info.errorMsg("Error in Pythia::init: Soft QCD only with resolved"
1255  " photons with lepton beams.");
1256  return false;
1257 
1258  }
1259  // Otherwise photon-photon within lepton beams OK.
1260  else return true;
1261  }
1262 
1263  // Lepton-lepton collisions OK (including neutrinos) if both (un)resolved
1264  else if (isUnresolvedA == isUnresolvedB) return true;
1265  }
1266 
1267  // MBR model only implemented for pp/ppbar/pbarp collisions.
1268  int PomFlux = settings.mode("SigmaDiffractive:PomFlux");
1269  if (PomFlux == 5) {
1270  bool ispp = (idAabs == 2212 && idBabs == 2212);
1271  bool ispbarpbar = (idA == -2212 && idB == -2212);
1272  if (ispp && !ispbarpbar) return true;
1273  info.errorMsg("Error in Pythia::init: cannot handle this beam combination"
1274  " with PomFlux == 5");
1275  return false;
1276  }
1277 
1278  // Hadron-hadron collisions OK, with Pomeron counted as hadron.
1279  bool isHadronA = (idAabs == 2212) || (idAabs == 2112) || (idA == 111)
1280  || (idAabs == 211) || (idA == 990);
1281  bool isHadronB = (idBabs == 2212) || (idBabs == 2112) || (idB == 111)
1282  || (idBabs == 211) || (idB == 990);
1283  int modeUnresolvedHadron = settings.mode("BeamRemnants:unresolvedHadron");
1284  if (isHadronA && modeUnresolvedHadron%2 == 1) isUnresolvedA = true;
1285  if (isHadronB && modeUnresolvedHadron > 1) isUnresolvedB = true;
1286  if (isHadronA && isHadronB) {
1287  // lepton2gamma flag with hadron beams may cause problems.
1288  if (beamHasGamma) {
1289  info.errorMsg("Error in Pythia::init: lepton2gamma should be off for"
1290  " hadron+hadron collision");
1291  return false;
1292  } else {
1293  return true;
1294  }
1295  }
1296 
1297  // Photon-photon collisions.
1298  if ( (idAabs == 22) && (idBabs == 22) ) {
1299 
1300  // No non-diffractive events for unresolved photon-photon.
1301  if ( ( !beamAisResGamma || !beamBisResGamma )
1302  && settings.flag("SoftQCD:nonDiffractive") ) {
1303  info.errorMsg("Error in Pythia::init: Soft QCD only with resolved"
1304  " photons.");
1305  }
1306 
1307  // lepton2gamma flag with photon beams may cause problems.
1308  if (beamHasGamma) {
1309  info.errorMsg("Error in Pythia::init: lepton2gamma should be off for"
1310  " hadron+hadron collision");
1311  return false;
1312  }
1313 
1314  // Otherwise OK.
1315  else return true;
1316  }
1317 
1318  // Gamma+hadron mode OK.
1319  if ( (isGammaA && isHadronB) || (isGammaB && isHadronA) )
1320  return true;
1321 
1322  // Lepton-hadron collisions OK for DIS processes or LHEF input,
1323  // although still primitive. Also e+p with real photons.
1324  if ( (isLeptonA && isHadronB) || (isHadronA && isLeptonB) ) {
1325  bool doDIS = settings.flag("WeakBosonExchange:all")
1326  || settings.flag("WeakBosonExchange:ff2ff(t:gmZ)")
1327  || settings.flag("WeakBosonExchange:ff2ff(t:W)")
1328  || (frameType == 4);
1329  if (doDIS || beamHasGamma ) return true;
1330  }
1331 
1332  // If no case above then failed.
1333  info.errorMsg("Error in Pythia::init: cannot handle this beam combination");
1334  return false;
1335 
1336 }
1337 
1338 //--------------------------------------------------------------------------
1339 
1340 // Calculate kinematics at initialization. Store beam four-momenta.
1341 
1342 bool Pythia::initKinematics() {
1343 
1344  // Find masses. Initial guess that we are in CM frame.
1345  mA = particleData.m0(idA);
1346  mB = particleData.m0(idB);
1347  betaZ = 0.;
1348  gammaZ = 1.;
1349 
1350  // Collinear beams not in CM frame: find CM energy.
1351  if (boostType == 2) {
1352  eA = max(eA, mA);
1353  eB = max(eB, mB);
1354  pzA = sqrt(eA*eA - mA*mA);
1355  pzB = -sqrt(eB*eB - mB*mB);
1356  pAinit = Vec4( 0., 0., pzA, eA);
1357  pBinit = Vec4( 0., 0., pzB, eB);
1358  eCM = sqrt( pow2(eA + eB) - pow2(pzA + pzB) );
1359 
1360  // Find boost to rest frame.
1361  betaZ = (pzA + pzB) / (eA + eB);
1362  gammaZ = (eA + eB) / eCM;
1363  if (abs(betaZ) < 1e-10) boostType = 1;
1364  }
1365 
1366  // Completely general beam directions: find CM energy.
1367  else if (boostType == 3) {
1368  eA = sqrt( pxA*pxA + pyA*pyA + pzA*pzA + mA*mA);
1369  eB = sqrt( pxB*pxB + pyB*pyB + pzB*pzB + mB*mB);
1370  pAinit = Vec4( pxA, pyA, pzA, eA);
1371  pBinit = Vec4( pxB, pyB, pzB, eB);
1372  eCM = (pAinit + pBinit).mCalc();
1373 
1374  // Find boost+rotation needed to move from/to CM frame.
1375  MfromCM.reset();
1376  MfromCM.fromCMframe( pAinit, pBinit);
1377  MtoCM = MfromCM;
1378  MtoCM.invert();
1379  }
1380 
1381  // Fail if CM energy below beam masses.
1382  if (eCM < mA + mB) {
1383  info.errorMsg("Error in Pythia::initKinematics: too low energy");
1384  return false;
1385  }
1386 
1387  // Set up CM-frame kinematics with beams along +-z axis.
1388  pzAcm = 0.5 * sqrtpos( (eCM + mA + mB) * (eCM - mA - mB)
1389  * (eCM - mA + mB) * (eCM + mA - mB) ) / eCM;
1390  pzBcm = -pzAcm;
1391  eA = sqrt(mA*mA + pzAcm*pzAcm);
1392  eB = sqrt(mB*mB + pzBcm*pzBcm);
1393 
1394  // If in CM frame then store beam four-vectors (else already done above).
1395  if (boostType != 2 && boostType != 3) {
1396  pAinit = Vec4( 0., 0., pzAcm, eA);
1397  pBinit = Vec4( 0., 0., pzBcm, eB);
1398  }
1399 
1400  // Store main info for access in process generation.
1401  info.setBeamA( idA, pzAcm, eA, mA);
1402  info.setBeamB( idB, pzBcm, eB, mB);
1403  info.setECM( eCM);
1404 
1405  // Must allow for generic boost+rotation when beam momentum spread.
1406  if (doMomentumSpread) boostType = 3;
1407 
1408  // Done.
1409  return true;
1410 
1411 }
1412 
1413 //--------------------------------------------------------------------------
1414 
1415 // Set up pointers to PDFs.
1416 
1417 bool Pythia::initPDFs() {
1418 
1419  // Delete any PDF's created in a previous init call.
1420  if (useNewPdfHard) {
1421  if (pdfHardAPtr != pdfAPtr) {
1422  delete pdfHardAPtr;
1423  pdfHardAPtr = 0;
1424  }
1425  if (pdfHardBPtr != pdfBPtr) {
1426  delete pdfHardBPtr;
1427  pdfHardBPtr = 0;
1428  }
1429  useNewPdfHard = false;
1430  }
1431  if (useNewPdfA) {
1432  delete pdfAPtr;
1433  useNewPdfA = false;
1434  pdfAPtr = 0;
1435  }
1436  if (useNewPdfB) {
1437  delete pdfBPtr;
1438  useNewPdfB = false;
1439  pdfBPtr = 0;
1440  }
1441  if (useNewPdfPomA) {
1442  delete pdfPomAPtr;
1443  useNewPdfPomA = false;
1444  pdfPomAPtr = 0;
1445  }
1446  if (useNewPdfPomB) {
1447  delete pdfPomBPtr;
1448  useNewPdfPomB = false;
1449  pdfPomBPtr = 0;
1450  }
1451  if (useNewPdfGamA) {
1452  delete pdfGamAPtr;
1453  useNewPdfGamA = false;
1454  pdfGamAPtr = 0;
1455  }
1456  if (useNewPdfGamB) {
1457  delete pdfGamBPtr;
1458  useNewPdfGamB = false;
1459  pdfGamBPtr = 0;
1460  }
1461  if (useNewPdfHardGamA) {
1462  delete pdfHardGamAPtr;
1463  useNewPdfHardGamA = false;
1464  pdfHardGamAPtr = 0;
1465  }
1466  if (useNewPdfHardGamB) {
1467  delete pdfHardGamBPtr;
1468  useNewPdfHardGamB = false;
1469  pdfHardGamBPtr = 0;
1470  }
1471  if (useNewPdfUnresA) {
1472  delete pdfUnresAPtr;
1473  useNewPdfUnresA = false;
1474  pdfUnresAPtr = 0;
1475  }
1476  if (useNewPdfUnresB) {
1477  delete pdfUnresBPtr;
1478  useNewPdfUnresB = false;
1479  pdfUnresBPtr = 0;
1480  }
1481  if (useNewPdfUnresGamA) {
1482  delete pdfUnresGamAPtr;
1483  useNewPdfUnresGamA = false;
1484  pdfUnresGamAPtr = 0;
1485  }
1486  if (useNewPdfUnresGamB) {
1487  delete pdfUnresGamBPtr;
1488  useNewPdfUnresGamB = false;
1489  pdfUnresGamBPtr = 0;
1490  }
1491  if (useNewPdfVMDA) {
1492  delete pdfVMDAPtr;
1493  useNewPdfVMDA = false;
1494  pdfVMDAPtr = 0;
1495  }
1496  if (useNewPdfVMDB) {
1497  delete pdfVMDBPtr;
1498  useNewPdfVMDB = false;
1499  pdfVMDBPtr = 0;
1500  }
1501 
1502  // Optionally set up photon PDF's for lepton -> gamma collisions. Done before
1503  // the main PDFs so that the gamma pointer can be used for the main PDF
1504  // (lepton). Both set also in case that only one of the photons is resolved.
1505  bool setupGammaBeams = (settings.flag("PDF:lepton2gamma")
1506  && (gammaMode < 4) );
1507  if (setupGammaBeams) {
1508  if ( (abs(idA) == 11 || abs(idA) == 13 || abs(idA) == 15)
1509  && pdfGamAPtr == 0 ) {
1510  pdfGamAPtr = getPDFPtr(22, 1, "A");
1511  if (!pdfGamAPtr->isSetup()) return false;
1512  useNewPdfGamA = true;
1513 
1514  // Set also unresolved photon beam when also unresolved photons.
1515  if (gammaMode != 1) {
1516  pdfUnresGamAPtr = getPDFPtr(22, 1, "A", false);
1517  if (!pdfUnresGamAPtr->isSetup()) return false;
1518  useNewPdfUnresGamA = true;
1519  }
1520 
1521  // Set up optional hard photon PDF pointers.
1522  if (settings.flag("PDF:useHard")) {
1523  pdfHardGamAPtr = getPDFPtr(22, 2);
1524  if (!pdfHardGamAPtr->isSetup()) return false;
1525  useNewPdfHardGamA = true;
1526  } else pdfHardGamAPtr = pdfGamAPtr;
1527  }
1528  if ( (abs(idB) == 11 || abs(idB) == 13 || abs(idB) == 15)
1529  && pdfGamBPtr == 0 ) {
1530  pdfGamBPtr = getPDFPtr(22, 1, "B");
1531  if (!pdfGamBPtr->isSetup()) return false;
1532  useNewPdfGamB = true;
1533 
1534  // Set also unresolved photon beam when also unresolved photons.
1535  if (gammaMode != 1) {
1536  pdfUnresGamBPtr = getPDFPtr(22, 1, "B", false);
1537  if (!pdfUnresGamBPtr->isSetup()) return false;
1538  useNewPdfUnresGamB = true;
1539  }
1540 
1541  // Set up optional hard photon PDF pointers.
1542  if (settings.flag("PDF:useHard")) {
1543  pdfHardGamBPtr = getPDFPtr(22, 2, "B");
1544  if (!pdfHardGamBPtr->isSetup()) return false;
1545  useNewPdfHardGamB = true;
1546  } else pdfHardGamBPtr = pdfGamBPtr;
1547  }
1548  }
1549 
1550  // Set up the PDF's, if not already done.
1551  if (pdfAPtr == 0) {
1552  pdfAPtr = getPDFPtr(idA);
1553  if (pdfAPtr == 0 || !pdfAPtr->isSetup()) {
1554  info.errorMsg("Error in Pythia::init: "
1555  "could not set up PDF for beam A");
1556  return false;
1557  }
1558  pdfHardAPtr = pdfAPtr;
1559  useNewPdfA = true;
1560  }
1561  if (pdfBPtr == 0) {
1562  pdfBPtr = getPDFPtr(idB, 1, "B");
1563  if (pdfBPtr == 0 || !pdfBPtr->isSetup()) {
1564  info.errorMsg("Error in Pythia::init: "
1565  "could not set up PDF for beam B");
1566  return false;
1567  }
1568  pdfHardBPtr = pdfBPtr;
1569  useNewPdfB = true;
1570  }
1571 
1572  // Optionally set up separate PDF's for hard process.
1573  if (settings.flag("PDF:useHard") && useNewPdfA && useNewPdfB) {
1574  pdfHardAPtr = getPDFPtr(idA, 2);
1575  if (!pdfHardAPtr->isSetup()) return false;
1576  pdfHardBPtr = getPDFPtr(idB, 2, "B");
1577  if (!pdfHardBPtr->isSetup()) return false;
1578  useNewPdfHard = true;
1579  }
1580 
1581  // Optionally use nuclear modifications for hard process PDFs.
1582  if (settings.flag("PDF:useHardNPDFA")) {
1583  int idANucleus = settings.mode("PDF:nPDFBeamA");
1584  pdfHardAPtr = getPDFPtr(idANucleus, 2, "A");
1585  if (!pdfHardAPtr->isSetup()) {
1586  info.errorMsg("Error in Pythia::init: "
1587  "could not set up nuclear PDF for beam A");
1588  return false;
1589  }
1590  useNewPdfHard = true;
1591  }
1592  if (settings.flag("PDF:useHardNPDFB")) {
1593  int idBNucleus = settings.mode("PDF:nPDFBeamB");
1594  pdfHardBPtr = getPDFPtr(idBNucleus, 2, "B");
1595  if (!pdfHardBPtr->isSetup()) {
1596  info.errorMsg("Error in Pythia::init: "
1597  "could not set up nuclear PDF for beam B");
1598  return false;
1599  }
1600  useNewPdfHard = true;
1601  }
1602 
1603  // Optionally set up additional unresolved PDFs for photon beams.
1604  if ( (idA == 22 || idB == 22) && gammaMode != 1 ) {
1605  if ( idA == 22 && pdfUnresAPtr == 0 ) {
1606  pdfUnresAPtr = getPDFPtr(idA, 1, "A", false);
1607  if (!pdfUnresAPtr->isSetup()) return false;
1608  useNewPdfUnresA = true;
1609  }
1610  if ( idB == 22 && pdfUnresBPtr == 0 ) {
1611  pdfUnresBPtr = getPDFPtr(idB, 1, "B", false);
1612  if (!pdfUnresBPtr->isSetup()) return false;
1613  useNewPdfUnresB = true;
1614  }
1615  }
1616 
1617  // Optionally set up additional unresolved PDFs for photon beam from lepton.
1618  if ( (abs(idA) == 11 || abs(idA) == 13 || abs(idA) == 15)
1619  && beamHasGamma && gammaMode != 1 ) {
1620  if ( pdfUnresAPtr == 0 ) {
1621  pdfUnresAPtr = getPDFPtr(idA, 1, "A", false);
1622  if (!pdfUnresAPtr->isSetup()) return false;
1623  useNewPdfUnresA = true;
1624  }
1625  }
1626  if ( (abs(idB) == 11 || abs(idB) == 13 || abs(idB) == 15)
1627  && beamHasGamma && gammaMode != 1 ) {
1628  if ( pdfUnresBPtr == 0 ) {
1629  pdfUnresBPtr = getPDFPtr(idB, 1, "B", false);
1630  if (!pdfUnresBPtr->isSetup()) return false;
1631  useNewPdfUnresB = true;
1632  }
1633  }
1634 
1635  // Optionally set up Pomeron PDF's for diffractive physics.
1636  if ( doDiffraction || doHardDiff) {
1637  if (pdfPomAPtr == 0) {
1638  pdfPomAPtr = getPDFPtr(990);
1639  useNewPdfPomA = true;
1640  }
1641  if (pdfPomBPtr == 0) {
1642  pdfPomBPtr = getPDFPtr(990);
1643  useNewPdfPomB = true;
1644  }
1645  }
1646 
1647  // Optionally set up Pomeron PDF's for diffractive physics.
1648  if ( doSoftQCD && (doVMDsideA || doVMDsideB)) {
1649  if (pdfVMDAPtr == 0) {
1650  pdfVMDAPtr = getPDFPtr(111);
1651  useNewPdfVMDA = true;
1652  }
1653  if (pdfVMDBPtr == 0) {
1654  pdfVMDBPtr = getPDFPtr(111);
1655  useNewPdfVMDB = true;
1656  }
1657  }
1658 
1659  // Done.
1660  return true;
1661 
1662 }
1663 
1664 //--------------------------------------------------------------------------
1665 
1666 // Main routine to generate the next event, using internal machinery.
1667 
1668 bool Pythia::next() {
1669 
1670  // Check that constructor worked.
1671  if (!isConstructed) return false;
1672 
1673  // Check if we the generation is taken over by the HeavyIons object.
1674  // Allows HeavyIons::next to call next for this Pythia object
1675  // without going into a loop.
1676  if ( doHeavyIons ) {
1677  doHeavyIons = false;
1678  bool ok = heavyIonsPtr->next();
1679  doHeavyIons = true;
1680  return ok;
1681  }
1682 
1683  // Regularly print how many events have been generated.
1684  int nPrevious = info.getCounter(3);
1685  if (nCount > 0 && nPrevious > 0 && nPrevious%nCount == 0)
1686  cout << "\n Pythia::next(): " << nPrevious
1687  << " events have been generated " << endl;
1688 
1689  // Set/reset info counters specific to each event.
1690  info.addCounter(3);
1691  for (int i = 10; i < 13; ++i) info.setCounter(i);
1692 
1693  // Simpler option when no hard process, i.e. mainly hadron level.
1694  if (!doProcessLevel) {
1695 
1696  // Optionally fetch in resonance decays from LHA interface.
1697  if (doLHA && !processLevel.nextLHAdec( event)) {
1698  if (info.atEndOfFile()) info.errorMsg("Abort from Pythia::next:"
1699  " reached end of Les Houches Events File");
1700  return false;
1701  }
1702 
1703  // Reset info and partonSystems arrays (while event record contains data).
1704  info.clear();
1705  partonSystems.clear();
1706 
1707  // Set correct energy for system.
1708  Vec4 pSum = 0.;
1709  for (int i = 1; i < event.size(); ++i)
1710  if (event[i].isFinal()) pSum += event[i].p();
1711  event[0].p( pSum );
1712  event[0].m( pSum.mCalc() );
1713 
1714  // Generate hadronization and decays.
1715  bool status = forceHadronLevel();
1716  if (status) info.addCounter(4);
1717  if (status && nPrevious < nShowEvt) event.list(showSaV, showMaD);
1718  return status;
1719  }
1720 
1721  // Reset arrays.
1722  info.clear();
1723  process.clear();
1724  event.clear();
1725  partonSystems.clear();
1726  beamA.clear();
1727  beamB.clear();
1728  beamPomA.clear();
1729  beamPomB.clear();
1730  beamGamA.clear();
1731  beamGamB.clear();
1732  beamVMDA.clear();
1733  beamVMDB.clear();
1734 
1735  // Pick current beam valence flavours (for pi0, K0S, K0L, Pomeron).
1736  beamA.newValenceContent();
1737  beamB.newValenceContent();
1738  if ( doDiffraction || doHardDiff) {
1739  beamPomA.newValenceContent();
1740  beamPomB.newValenceContent();
1741  }
1742  if (doVMDsideA) beamVMDA.newValenceContent();
1743  if (doVMDsideB) beamVMDB.newValenceContent();
1744 
1745  // Can only generate event if initialization worked.
1746  if (!isInit) {
1747  info.errorMsg("Abort from Pythia::next: "
1748  "not properly initialized so cannot generate events");
1749  return false;
1750  }
1751 
1752  // Pick beam momentum spread and beam vertex.
1753  if (doMomentumSpread || doVertexSpread) beamShapePtr->pick();
1754 
1755  // Recalculate kinematics when beam momentum spread.
1756  if (doMomentumSpread) nextKinematics();
1757 
1758  // Outer loop over hard processes; only relevant for user-set vetoes.
1759  for ( ; ; ) {
1760 
1761  info.addCounter(10);
1762  bool hasVetoed = false;
1763  bool hasVetoedDiff = false;
1764 
1765  // Provide the hard process that starts it off. Only one try.
1766  info.clear();
1767  process.clear();
1768  partonSystems.clear();
1769 
1770  // Reset the event information. Necessary if the previous event was read
1771  // from LHEF, while the current event is not read from LHEF.
1772  info.setLHEF3EventInfo();
1773 
1774  if ( !processLevel.next( process) ) {
1775  if (doLHA && info.atEndOfFile()) info.errorMsg("Abort from "
1776  "Pythia::next: reached end of Les Houches Events File");
1777  else info.errorMsg("Abort from Pythia::next: "
1778  "processLevel failed; giving up");
1779  return false;
1780  }
1781 
1782  info.addCounter(11);
1783 
1784  // Update tried and selected events immediately after next event was
1785  // generated. Note: This does not accumulate cross section.
1786  processLevel.accumulate(false);
1787 
1788  // Possibility for a user veto of the process-level event.
1789  if (doVetoProcess) {
1790  hasVetoed = userHooksPtr->doVetoProcessLevel( process);
1791  if (hasVetoed) {
1792  if (abortIfVeto) return false;
1793  continue;
1794  }
1795  }
1796 
1797  // Possibility to perform matrix element merging for this event.
1798  if (doMerging) {
1799  int veto = mergingPtr->mergeProcess( process );
1800  // Apply possible merging scale cut.
1801  if (veto == -1) {
1802  hasVetoed = true;
1803  if (abortIfVeto) return false;
1804  continue;
1805  // Exit because of vanishing no-emission probability.
1806  } else if (veto == 0) {
1807  event = process;
1808  break;
1809  }
1810 
1811  // Redo resonance decays after the merging, in case the resonance
1812  // structure has been changed because of reclusterings.
1813  if (veto == 2 && doResDec) processLevel.nextDecays( process);
1814  }
1815 
1816  // Possibility to stop the generation at this stage.
1817  if (!doPartonLevel) {
1818  boostAndVertex( true, true);
1819  processLevel.accumulate();
1820  info.addCounter(4);
1821  if (doLHA && nPrevious < nShowLHA) lhaUpPtr->listEvent();
1822  if (nPrevious < nShowInfo) info.list();
1823  if (nPrevious < nShowProc) process.list(showSaV, showMaD);
1824  return true;
1825  }
1826 
1827  // Save spare copy of process record in case of problems.
1828  Event processSave = process;
1829  int sizeMPI = info.sizeMPIarrays();
1830  info.addCounter(12);
1831  for (int i = 14; i < 19; ++i) info.setCounter(i);
1832 
1833  // Allow up to ten tries for parton- and hadron-level processing.
1834  bool physical = true;
1835  for (int iTry = 0; iTry < NTRY; ++iTry) {
1836 
1837  info.addCounter(14);
1838  physical = true;
1839  hasVetoed = false;
1840 
1841  // Restore original process record if problems.
1842  if (iTry > 0) process = processSave;
1843  if (iTry > 0) info.resizeMPIarrays( sizeMPI);
1844 
1845  // Reset event record and (extracted partons from) beam remnants.
1846  event.clear();
1847  beamA.clear();
1848  beamB.clear();
1849  beamPomA.clear();
1850  beamPomB.clear();
1851  beamGamA.clear();
1852  beamGamB.clear();
1853  beamVMDA.clear();
1854  beamVMDB.clear();
1855  partonSystems.clear();
1856 
1857  // Parton-level evolution: ISR, FSR, MPI.
1858  if ( !partonLevel.next( process, event) ) {
1859 
1860  // Abort event generation if parton level is set to abort.
1861  if (info.getAbortPartonLevel()) return false;
1862 
1863  // Skip to next hard process for failure owing to deliberate veto,
1864  // or alternatively retry for the same hard process.
1865  hasVetoed = partonLevel.hasVetoed();
1866  if (hasVetoed) {
1867  if (retryPartonLevel) {
1868  --iTry;
1869  continue;
1870  }
1871  if (abortIfVeto) return false;
1872  break;
1873  }
1874 
1875  // If hard diffractive event has been discarded retry partonLevel.
1876  hasVetoedDiff = partonLevel.hasVetoedDiff();
1877  if (hasVetoedDiff) {
1878  info.errorMsg("Warning in Pythia::next: "
1879  "discarding hard diffractive event from partonLevel; try again");
1880  break;
1881  }
1882 
1883  // Else make a new try for other failures.
1884  info.errorMsg("Error in Pythia::next: "
1885  "partonLevel failed; try again");
1886  physical = false;
1887  continue;
1888  }
1889  info.addCounter(15);
1890 
1891  // Possibility for a user veto of the parton-level event.
1892  if (doVetoPartons) {
1893  hasVetoed = userHooksPtr->doVetoPartonLevel( event);
1894  if (hasVetoed) {
1895  if (abortIfVeto) return false;
1896  break;
1897  }
1898  }
1899 
1900  // Boost to lab frame (before decays, for vertices).
1901  boostAndVertex( true, true);
1902 
1903  // Possibility to stop the generation at this stage.
1904  if (!doHadronLevel) {
1905  processLevel.accumulate();
1906  partonLevel.accumulate();
1907  // Optionally check final event for problems.
1908  if (checkEvent && !check()) {
1909  info.errorMsg("Abort from Pythia::next: "
1910  "check of event revealed problems");
1911  return false;
1912  }
1913  info.addCounter(4);
1914  if (doLHA && nPrevious < nShowLHA) lhaUpPtr->listEvent();
1915  if (nPrevious < nShowInfo) info.list();
1916  if (nPrevious < nShowProc) process.list(showSaV, showMaD);
1917  if (nPrevious < nShowEvt) event.list(showSaV, showMaD);
1918  return true;
1919  }
1920 
1921  // Hadron-level: hadronization, decays.
1922  info.addCounter(16);
1923  if ( !hadronLevel.next( event) ) {
1924  info.errorMsg("Error in Pythia::next: "
1925  "hadronLevel failed; try again");
1926  physical = false;
1927  continue;
1928  }
1929 
1930  // If R-hadrons have been formed, then (optionally) let them decay.
1931  if (decayRHadrons && rHadrons.exist() && !doRHadronDecays()) {
1932  info.errorMsg("Error in Pythia::next: "
1933  "decayRHadrons failed; try again");
1934  physical = false;
1935  continue;
1936  }
1937  info.addCounter(17);
1938 
1939  // Optionally check final event for problems.
1940  if (checkEvent && !check()) {
1941  info.errorMsg("Error in Pythia::next: "
1942  "check of event revealed problems");
1943  physical = false;
1944  continue;
1945  }
1946 
1947  // Stop parton- and hadron-level looping if you got this far.
1948  info.addCounter(18);
1949  break;
1950  }
1951 
1952  // If event vetoed then to make a new try.
1953  if (hasVetoed || hasVetoedDiff) {
1954  if (abortIfVeto) return false;
1955  continue;
1956  }
1957 
1958  // If event failed any other way (after ten tries) then give up.
1959  if (!physical) {
1960  info.errorMsg("Abort from Pythia::next: "
1961  "parton+hadronLevel failed; giving up");
1962  return false;
1963  }
1964 
1965  // Process- and parton-level statistics. Event scale.
1966  processLevel.accumulate();
1967  partonLevel.accumulate();
1968  event.scale( process.scale() );
1969 
1970  // End of outer loop over hard processes. Done with normal option.
1971  info.addCounter(13);
1972  break;
1973  }
1974 
1975  // List events.
1976  if (doLHA && nPrevious < nShowLHA) lhaUpPtr->listEvent();
1977  if (nPrevious < nShowInfo) info.list();
1978  if (nPrevious < nShowProc) process.list(showSaV,showMaD);
1979  if (nPrevious < nShowEvt) event.list(showSaV, showMaD);
1980 
1981  // Done.
1982  info.addCounter(4);
1983  return true;
1984 
1985 }
1986 
1987 //--------------------------------------------------------------------------
1988 
1989 // Generate only the hadronization/decay stage, using internal machinery.
1990 // The "event" instance should already contain a parton-level configuration.
1991 
1992 bool Pythia::forceHadronLevel(bool findJunctions) {
1993 
1994  // Can only generate event if initialization worked.
1995  if (!isInit) {
1996  info.errorMsg("Abort from Pythia::forceHadronLevel: "
1997  "not properly initialized so cannot generate events");
1998  return false;
1999  }
2000 
2001  // Check whether any junctions in system. (Normally done in ProcessLevel.)
2002  // Avoid it if there are no final-state coloured partons.
2003  if (findJunctions) {
2004  event.clearJunctions();
2005  for (int i = 0; i < event.size(); ++i)
2006  if (event[i].isFinal()
2007  && (event[i].col() != 0 || event[i].acol() != 0)) {
2008  processLevel.findJunctions( event);
2009  break;
2010  }
2011  }
2012 
2013  // Allow for CR before the hadronization.
2014  if (forceHadronLevelCR) {
2015 
2016  // Setup parton system for SK-I and SK-II colour reconnection.
2017  // Require all final state particles to have the Ws as mothers.
2018  if (reconnectMode == 3 || reconnectMode == 4) {
2019  partonSystems.clear();
2020  partonSystems.addSys();
2021  partonSystems.addSys();
2022  for (int i = 5;i < event.size();++i) {
2023  if (event[i].mother1() - 3 < 0 || event[i].mother1() - 3 > 1) {
2024  info.errorMsg("Error in Pythia::forceHadronLevel: "
2025  " Event is not setup correctly for SK-I or SK-II CR");
2026  return false;
2027  }
2028  partonSystems.addOut(event[i].mother1() - 3,i);
2029  }
2030  }
2031 
2032  // save spare copy of event in case of failure.
2033  Event spareEvent = event;
2034  bool colCorrect = false;
2035 
2036  // Allow up to ten tries for CR.
2037  for (int iTry = 0; iTry < NTRY; ++ iTry) {
2038  colourReconnection.next(event, 0);
2039  if (junctionSplitting.checkColours(event)) {
2040  colCorrect = true;
2041  break;
2042  }
2043  else event = spareEvent;
2044  }
2045 
2046  if (!colCorrect) {
2047  info.errorMsg("Error in Pythia::forceHadronLevel: "
2048  "Colour reconnection failed.");
2049  return false;
2050  }
2051  }
2052 
2053  // Save spare copy of event in case of failure.
2054  Event spareEvent = event;
2055 
2056  // Allow up to ten tries for hadron-level processing.
2057  bool physical = true;
2058  for (int iTry = 0; iTry < NTRY; ++ iTry) {
2059  physical = true;
2060 
2061  // Check whether any resonances need to be handled at process level.
2062  if (doResDec) {
2063  process = event;
2064  processLevel.nextDecays( process);
2065 
2066  // Allow for showers if decays happened at process level.
2067  if (process.size() > event.size()) {
2068  if (doFSRinRes) {
2069  partonLevel.setupShowerSys( process, event);
2070  partonLevel.resonanceShowers( process, event, false);
2071  } else event = process;
2072  }
2073  }
2074 
2075  // Hadron-level: hadronization, decays.
2076  if (hadronLevel.next( event)) break;
2077 
2078  // If failure then warn, restore original configuration and try again.
2079  info.errorMsg("Error in Pythia::forceHadronLevel: "
2080  "hadronLevel failed; try again");
2081  physical = false;
2082  event = spareEvent;
2083  }
2084 
2085  // Done for simpler option.
2086  if (!physical) {
2087  info.errorMsg("Abort from Pythia::forceHadronLevel: "
2088  "hadronLevel failed; giving up");
2089  return false;
2090  }
2091 
2092  // Optionally check final event for problems.
2093  if (checkEvent && !check()) {
2094  info.errorMsg("Abort from Pythia::forceHadronLevel: "
2095  "check of event revealed problems");
2096  return false;
2097  }
2098 
2099  // Done.
2100  return true;
2101 
2102 }
2103 
2104 //--------------------------------------------------------------------------
2105 
2106 // Recalculate kinematics for each event when beam momentum has a spread.
2107 
2108 void Pythia::nextKinematics() {
2109 
2110  // Read out momentum shift to give current beam momenta.
2111  pAnow = pAinit + beamShapePtr->deltaPA();
2112  pAnow.e( sqrt(pAnow.pAbs2() + mA * mA) );
2113  pBnow = pBinit + beamShapePtr->deltaPB();
2114  pBnow.e( sqrt(pBnow.pAbs2() + mB * mB) );
2115 
2116  // Construct CM frame kinematics.
2117  eCM = (pAnow + pBnow).mCalc();
2118  pzAcm = 0.5 * sqrtpos( (eCM + mA + mB) * (eCM - mA - mB)
2119  * (eCM - mA + mB) * (eCM + mA - mB) ) / eCM;
2120  pzBcm = -pzAcm;
2121  eA = sqrt(mA*mA + pzAcm*pzAcm);
2122  eB = sqrt(mB*mB + pzBcm*pzBcm);
2123 
2124  // Set relevant info for other classes to use.
2125  info.setBeamA( idA, pzAcm, eA, mA);
2126  info.setBeamB( idB, pzBcm, eB, mB);
2127  info.setECM( eCM);
2128  beamA.newPzE( pzAcm, eA);
2129  beamB.newPzE( pzBcm, eB);
2130 
2131  // Set boost/rotation matrices from/to CM frame.
2132  MfromCM.reset();
2133  MfromCM.fromCMframe( pAnow, pBnow);
2134  MtoCM = MfromCM;
2135  MtoCM.invert();
2136 
2137 }
2138 
2139 //--------------------------------------------------------------------------
2140 
2141 // Boost from CM frame to lab frame, or inverse. Set production vertex.
2142 
2143 void Pythia::boostAndVertex( bool toLab, bool setVertex) {
2144 
2145  // Boost process from CM frame to lab frame.
2146  if (toLab) {
2147  if (boostType == 2) process.bst(0., 0., betaZ, gammaZ);
2148  else if (boostType == 3) process.rotbst(MfromCM);
2149 
2150  // Boost nonempty event from CM frame to lab frame.
2151  if (event.size() > 0) {
2152  if (boostType == 2) event.bst(0., 0., betaZ, gammaZ);
2153  else if (boostType == 3) event.rotbst(MfromCM);
2154  }
2155 
2156  // Boost process from lab frame to CM frame.
2157  } else {
2158  if (boostType == 2) process.bst(0., 0., -betaZ, gammaZ);
2159  else if (boostType == 3) process.rotbst(MtoCM);
2160 
2161  // Boost nonempty event from lab frame to CM frame.
2162  if (event.size() > 0) {
2163  if (boostType == 2) event.bst(0., 0., -betaZ, gammaZ);
2164  else if (boostType == 3) event.rotbst(MtoCM);
2165  }
2166  }
2167 
2168  // Set production vertex; assumes particles are in lab frame and at origin.
2169  if (setVertex && doVertexSpread) {
2170  Vec4 vertex = beamShapePtr->vertex();
2171  for (int i = 0; i < process.size(); ++i) process[i].vProd( vertex);
2172  for (int i = 0; i < event.size(); ++i) event[i].vProd( vertex);
2173  }
2174 
2175 }
2176 
2177 //--------------------------------------------------------------------------
2178 
2179 // Perform R-hadron decays, either as part of normal evolution or forced.
2180 
2181 bool Pythia::doRHadronDecays( ) {
2182 
2183  // Check if R-hadrons exist to be processed.
2184  if ( !rHadrons.exist() ) return true;
2185 
2186  // Do the R-hadron decay itself.
2187  if ( !rHadrons.decay( event) ) return false;
2188 
2189  // Perform showers in resonance decay chains.
2190  if ( !partonLevel.resonanceShowers( process, event, false) ) return false;
2191 
2192  // Subsequent hadronization and decays.
2193  if ( !hadronLevel.next( event) ) return false;
2194 
2195  // Done.
2196  return true;
2197 
2198 }
2199 
2200 //--------------------------------------------------------------------------
2201 
2202 // Print statistics on event generation.
2203 
2204 void Pythia::stat() {
2205 
2206  if ( doHeavyIons ) {
2207  heavyIonsPtr->stat();
2208  return;
2209  }
2210 
2211  // Read out settings for what to include.
2212  bool showPrL = settings.flag("Stat:showProcessLevel");
2213  bool showPaL = settings.flag("Stat:showPartonLevel");
2214  bool showErr = settings.flag("Stat:showErrors");
2215  bool reset = settings.flag("Stat:reset");
2216 
2217  // Statistics on cross section and number of events.
2218  if (doProcessLevel) {
2219  if (showPrL) processLevel.statistics(false);
2220  if (reset) processLevel.resetStatistics();
2221  }
2222 
2223  // Statistics from other classes, currently multiparton interactions.
2224  if (showPaL) partonLevel.statistics(false);
2225  if (reset) partonLevel.resetStatistics();
2226 
2227  // Merging statistics.
2228  if (doMerging) mergingPtr->statistics();
2229 
2230  // Summary of which and how many warnings/errors encountered.
2231  if (showErr) info.errorStatistics();
2232  if (reset) info.errorReset();
2233 
2234 }
2235 
2236 //--------------------------------------------------------------------------
2237 
2238 // Write the Pythia banner, with symbol and version information.
2239 
2240 void Pythia::banner() {
2241 
2242  // Read in version number and last date of change.
2243  double versionNumber = settings.parm("Pythia:versionNumber");
2244  int versionDate = settings.mode("Pythia:versionDate");
2245  string month[12] = {"Jan", "Feb", "Mar", "Apr", "May", "Jun",
2246  "Jul", "Aug", "Sep", "Oct", "Nov", "Dec"};
2247 
2248  // Get date and time.
2249  time_t t = time(0);
2250  char dateNow[12];
2251  strftime(dateNow,12,"%d %b %Y",localtime(&t));
2252  char timeNow[9];
2253  strftime(timeNow,9,"%H:%M:%S",localtime(&t));
2254 
2255  cout << "\n"
2256  << " *-------------------------------------------"
2257  << "-----------------------------------------* \n"
2258  << " | "
2259  << " | \n"
2260  << " | *----------------------------------------"
2261  << "--------------------------------------* | \n"
2262  << " | | "
2263  << " | | \n"
2264  << " | | "
2265  << " | | \n"
2266  << " | | PPP Y Y TTTTT H H III A "
2267  << " Welcome to the Lund Monte Carlo! | | \n"
2268  << " | | P P Y Y T H H I A A "
2269  << " This is PYTHIA version " << fixed << setprecision(3)
2270  << setw(5) << versionNumber << " | | \n"
2271  << " | | PPP Y T HHHHH I AAAAA"
2272  << " Last date of change: " << setw(2) << versionDate%100
2273  << " " << month[ (versionDate/100)%100 - 1 ]
2274  << " " << setw(4) << versionDate/10000 << " | | \n"
2275  << " | | P Y T H H I A A"
2276  << " | | \n"
2277  << " | | P Y T H H III A A"
2278  << " Now is " << dateNow << " at " << timeNow << " | | \n"
2279  << " | | "
2280  << " | | \n"
2281  << " | | Christian Bierlich; Department of As"
2282  << "tronomy and Theoretical Physics, | | \n"
2283  << " | | Lund University, Solvegatan 14A, S"
2284  << "E-223 62 Lund, Sweden; | | \n"
2285  << " | | e-mail: christian.bierlich@thep.lu"
2286  << ".se | | \n"
2287  << " | | Nishita Desai; Laboratoire Charles C"
2288  << "oulomb (L2C), | | \n"
2289  << " | | CNRS-Universite de Montpellier, 34"
2290  << "090 Montpellier, France; | | \n"
2291  << " | | e-mail: nishita.desai@umontpellier"
2292  << ".fr | | \n"
2293  << " | | Nadine Fischer; School of Physics, "
2294  << " | | \n"
2295  << " | | Monash University, PO Box 27, 3800"
2296  << " Melbourne, Australia; | | \n"
2297  << " | | e-mail: nadine.fischer@monash.edu "
2298  << " | | \n"
2299  << " | | Ilkka Helenius; Institute for Theore"
2300  << "tical Physics, | | \n"
2301  << " | | Tuebingen University, Auf der Morge"
2302  << "nstelle 14, 72076 Tuebingen, Germany; | | \n"
2303  << " | | e-mail: ilkka.helenius@uni-tuebing"
2304  << "en.de | | \n"
2305  << " | | Philip Ilten; School of Physics and "
2306  << "Astronomy, | | \n"
2307  << " | | University of Birmingham, Birmingh"
2308  << "am, B152 2TT, UK; | | \n"
2309  << " | | e-mail: philten@cern.ch "
2310  << " | | \n"
2311  << " | | Leif Lonnblad; Department of Astrono"
2312  << "my and Theoretical Physics, | | \n"
2313  << " | | Lund University, Solvegatan 14A, S"
2314  << "E-223 62 Lund, Sweden; | | \n"
2315  << " | | e-mail: leif.lonnblad@thep.lu.se "
2316  << " | | \n"
2317  << " | | Stephen Mrenna; Computing Division, "
2318  << "Simulations Group, | | \n"
2319  << " | | Fermi National Accelerator Laborat"
2320  << "ory, MS 234, Batavia, IL 60510, USA; | | \n"
2321  << " | | e-mail: mrenna@fnal.gov "
2322  << " | | \n"
2323  << " | | Stefan Prestel; Theoretical Physics "
2324  << "Department, | | \n"
2325  << " | | Fermi National Accelerator Laborat"
2326  << "ory, MS 106, Batavia, IL 60510, USA; | | \n"
2327  << " | | e-mail: sprestel@fnal.gov "
2328  << " | | \n"
2329  << " | | Christine O. Rasmussen; Department o"
2330  << "f Astronomy and Theoretical Physics, | | \n"
2331  << " | | Lund University, Solvegatan 14A, S"
2332  << "E-223 62 Lund, Sweden; | | \n"
2333  << " | | e-mail: christine.rasmussen@thep.l"
2334  << "u.se | | \n"
2335  << " | | Torbjorn Sjostrand; Department of As"
2336  << "tronomy and Theoretical Physics, | | \n"
2337  << " | | Lund University, Solvegatan 14A, S"
2338  << "E-223 62 Lund, Sweden; | | \n"
2339  << " | | e-mail: torbjorn@thep.lu.se "
2340  << " | | \n"
2341  << " | | Peter Skands; School of Physics, "
2342  << " | | \n"
2343  << " | | Monash University, PO Box 27, 3800"
2344  << " Melbourne, Australia; | | \n"
2345  << " | | e-mail: peter.skands@monash.edu "
2346  << " | | \n"
2347  << " | | "
2348  << " | | \n"
2349  << " | | The main program reference is 'An Int"
2350  << "roduction to PYTHIA 8.2', | | \n"
2351  << " | | T. Sjostrand et al, Comput. Phys. Com"
2352  << "mun. 191 (2015) 159 | | \n"
2353  << " | | [arXiv:1410.3012 [hep-ph]] "
2354  << " | | \n"
2355  << " | | "
2356  << " | | \n"
2357  << " | | The main physics reference is the 'PY"
2358  << "THIA 6.4 Physics and Manual', | | \n"
2359  << " | | T. Sjostrand, S. Mrenna and P. Skands"
2360  << ", JHEP05 (2006) 026 [hep-ph/0603175] | | \n"
2361  << " | | "
2362  << " | | \n"
2363  << " | | An archive of program versions and do"
2364  << "cumentation is found on the web: | | \n"
2365  << " | | http://www.thep.lu.se/Pythia "
2366  << " | | \n"
2367  << " | | "
2368  << " | | \n"
2369  << " | | This program is released under the GN"
2370  << "U General Public Licence version 2. | | \n"
2371  << " | | Please respect the MCnet Guidelines f"
2372  << "or Event Generator Authors and Users. | | \n"
2373  << " | | "
2374  << " | | \n"
2375  << " | | Disclaimer: this program comes withou"
2376  << "t any guarantees. | | \n"
2377  << " | | Beware of errors and use common sense"
2378  << " when interpreting results. | | \n"
2379  << " | | "
2380  << " | | \n"
2381  << " | | Copyright (C) 2018 Torbjorn Sjostrand"
2382  << " | | \n"
2383  << " | | "
2384  << " | | \n"
2385  << " | | "
2386  << " | | \n"
2387  << " | *----------------------------------------"
2388  << "--------------------------------------* | \n"
2389  << " | "
2390  << " | \n"
2391  << " *-------------------------------------------"
2392  << "-----------------------------------------* \n" << endl;
2393 
2394 }
2395 
2396 //--------------------------------------------------------------------------
2397 
2398 // Check for lines in file that mark the beginning of new subrun.
2399 
2400 int Pythia::readSubrun(string line, bool warn) {
2401 
2402  // If empty line then done.
2403  int subrunLine = SUBRUNDEFAULT;
2404  if (line.find_first_not_of(" \n\t\v\b\r\f\a") == string::npos)
2405  return subrunLine;
2406 
2407  // If first character is not a letter, then done.
2408  string lineNow = line;
2409  int firstChar = lineNow.find_first_not_of(" \n\t\v\b\r\f\a");
2410  if (!isalpha(lineNow[firstChar])) return subrunLine;
2411 
2412  // Replace an equal sign by a blank to make parsing simpler.
2413  while (lineNow.find("=") != string::npos) {
2414  int firstEqual = lineNow.find_first_of("=");
2415  lineNow.replace(firstEqual, 1, " ");
2416  }
2417 
2418  // Get first word of a line.
2419  istringstream splitLine(lineNow);
2420  string name;
2421  splitLine >> name;
2422 
2423  // Replace two colons by one (:: -> :) to allow for such mistakes.
2424  while (name.find("::") != string::npos) {
2425  int firstColonColon = name.find_first_of("::");
2426  name.replace(firstColonColon, 2, ":");
2427  }
2428 
2429  // Convert to lowercase. If no match then done.
2430  if (toLower(name) != "main:subrun") return subrunLine;
2431 
2432  // Else find new subrun number and return it.
2433  splitLine >> subrunLine;
2434  if (!splitLine) {
2435  if (warn) cout << "\n PYTHIA Warning: Main:subrun number not"
2436  << " recognized; skip:\n " << line << endl;
2437  subrunLine = SUBRUNDEFAULT;
2438  }
2439  return subrunLine;
2440 
2441 }
2442 
2443 //--------------------------------------------------------------------------
2444 
2445 // Check for lines in file that mark the beginning or end of commented section.
2446 // Return +1 for beginning, -1 for end, 0 else.
2447 
2448 int Pythia::readCommented(string line) {
2449 
2450  // If less than two nontrivial characters on line then done.
2451  if (line.find_first_not_of(" \n\t\v\b\r\f\a") == string::npos) return 0;
2452  int firstChar = line.find_first_not_of(" \n\t\v\b\r\f\a");
2453  if (int(line.size()) < firstChar + 2) return 0;
2454 
2455  // If first two nontrivial characters are /* or */ then done.
2456  if (line.substr(firstChar, 2) == "/*") return +1;
2457  if (line.substr(firstChar, 2) == "*/") return -1;
2458 
2459  // Else done.
2460  return 0;
2461 
2462 }
2463 
2464 //--------------------------------------------------------------------------
2465 
2466 // Check that the final event makes sense: no unknown id codes;
2467 // charge and energy-momentum conserved.
2468 
2469 bool Pythia::check() {
2470 
2471  // Reset.
2472  bool physical = true;
2473  bool listVertices = false;
2474  bool listHistory = false;
2475  bool listSystems = false;
2476  bool listBeams = false;
2477  iErrId.resize(0);
2478  iErrCol.resize(0);
2479  iErrEpm.resize(0);
2480  iErrNan.resize(0);
2481  iErrNanVtx.resize(0);
2482  Vec4 pSum;
2483  double chargeSum = 0.;
2484 
2485  // Incoming beams counted with negative momentum and charge.
2486  if (doProcessLevel) {
2487  pSum = - (event[1].p() + event[2].p());
2488  chargeSum = - (event[1].charge() + event[2].charge());
2489 
2490  // If no ProcessLevel then sum final state of process record.
2491  } else if (process.size() > 0) {
2492  pSum = - process[0].p();
2493  for (int i = 0; i < process.size(); ++i)
2494  if (process[i].isFinal()) chargeSum -= process[i].charge();
2495 
2496  // If process not filled, then use outgoing primary in event.
2497  } else {
2498  pSum = - event[0].p();
2499  for (int i = 1; i < event.size(); ++i)
2500  if (event[i].statusAbs() < 10 || event[i].statusAbs() == 23)
2501  chargeSum -= event[i].charge();
2502  }
2503  double eLab = abs(pSum.e());
2504 
2505  // Loop over particles in the event.
2506  for (int i = 0; i < event.size(); ++i) {
2507 
2508  // Look for any unrecognized particle codes.
2509  int id = event[i].id();
2510  if (id == 0 || !particleData.isParticle(id)) {
2511  ostringstream errCode;
2512  errCode << ", i = " << i << ", id = " << id;
2513  info.errorMsg("Error in Pythia::check: "
2514  "unknown particle code", errCode.str());
2515  physical = false;
2516  iErrId.push_back(i);
2517 
2518  // Check that colour assignments are the expected ones.
2519  } else {
2520  int colType = event[i].colType();
2521  int col = event[i].col();
2522  int acol = event[i].acol();
2523  if ( event[i].statusAbs() / 10 == 8 ) acol = col = 0;
2524  if ( (colType == 0 && (col > 0 || acol > 0))
2525  || (colType == 1 && (col <= 0 || acol > 0))
2526  || (colType == -1 && (col > 0 || acol <= 0))
2527  || (colType == 2 && (col <= 0 || acol <= 0)) ) {
2528  ostringstream errCode;
2529  errCode << ", i = " << i << ", id = " << id << " cols = " << col
2530  << " " << acol;
2531  info.errorMsg("Error in Pythia::check: "
2532  "incorrect colours", errCode.str());
2533  physical = false;
2534  iErrCol.push_back(i);
2535  }
2536  }
2537 
2538  // Some intermediate shower partons excepted from (E, p, m) consistency.
2539  bool checkMass = event[i].statusAbs() != 49 && event[i].statusAbs() != 59;
2540 
2541  // Look for particles with mismatched or not-a-number energy/momentum/mass.
2542  if (abs(event[i].px()) >= 0. && abs(event[i].py()) >= 0.
2543  && abs(event[i].pz()) >= 0. && abs(event[i].e()) >= 0.
2544  && abs(event[i].m()) >= 0.) {
2545  double errMass = abs(event[i].mCalc() - event[i].m())
2546  / max( 1.0, event[i].e());
2547  if (checkMass && errMass > mTolErr) {
2548  info.errorMsg("Error in Pythia::check: "
2549  "unmatched particle energy/momentum/mass");
2550  physical = false;
2551  iErrEpm.push_back(i);
2552  } else if (checkMass && errMass > mTolWarn) {
2553  info.errorMsg("Warning in Pythia::check: "
2554  "not quite matched particle energy/momentum/mass");
2555  }
2556  } else {
2557  info.errorMsg("Error in Pythia::check: "
2558  "not-a-number energy/momentum/mass");
2559  physical = false;
2560  iErrNan.push_back(i);
2561  }
2562 
2563  // Look for particles with not-a-number vertex/lifetime.
2564  if (abs(event[i].xProd()) >= 0. && abs(event[i].yProd()) >= 0.
2565  && abs(event[i].zProd()) >= 0. && abs(event[i].tProd()) >= 0.
2566  && abs(event[i].tau()) >= 0.) ;
2567  else {
2568  info.errorMsg("Error in Pythia::check: "
2569  "not-a-number vertex/lifetime");
2570  physical = false;
2571  listVertices = true;
2572  iErrNanVtx.push_back(i);
2573  }
2574 
2575  // Add final-state four-momentum and charge.
2576  if (event[i].isFinal()) {
2577  pSum += event[i].p();
2578  chargeSum += event[i].charge();
2579  }
2580 
2581  // End of particle loop.
2582  }
2583 
2584  // Check energy-momentum/charge conservation.
2585  double epDev = abs(pSum.e()) + abs(pSum.px()) + abs(pSum.py())
2586  + abs(pSum.pz());
2587  if (epDev > epTolErr * eLab) {
2588  info.errorMsg("Error in Pythia::check: energy-momentum not conserved");
2589  physical = false;
2590  } else if (epDev > epTolWarn * eLab) {
2591  info.errorMsg("Warning in Pythia::check: "
2592  "energy-momentum not quite conserved");
2593  }
2594  if (abs(chargeSum) > 0.1) {
2595  info.errorMsg("Error in Pythia::check: charge not conserved");
2596  physical = false;
2597  }
2598 
2599  // Check that beams and event records agree on incoming partons.
2600  // Only meaningful for resolved beams.
2601  if (info.isResolved() && !info.hasUnresolvedBeams())
2602  for (int iSys = 0; iSys < beamA.sizeInit(); ++iSys) {
2603  int eventANw = partonSystems.getInA(iSys);
2604  int eventBNw = partonSystems.getInB(iSys);
2605  // For photon beams from leptons make sure to use correct beams.
2606  int beamANw = ( beamA.getGammaMode() == 0 || !beamHasGamma
2607  || (beamA.getGammaMode() == 2 && beamB.getGammaMode() == 2)) ?
2608  beamA[iSys].iPos() : beamGamA[iSys].iPos();
2609  int beamBNw = ( beamB.getGammaMode() == 0 || !beamHasGamma
2610  || (beamB.getGammaMode() == 2 && beamA.getGammaMode() == 2)) ?
2611  beamB[iSys].iPos() : beamGamB[iSys].iPos();
2612  if (eventANw != beamANw || eventBNw != beamBNw) {
2613  info.errorMsg("Error in Pythia::check: "
2614  "event and beams records disagree");
2615  physical = false;
2616  listSystems = true;
2617  listBeams = true;
2618  }
2619  }
2620 
2621  // Check that mother and daughter information match for each particle.
2622  vector<int> noMot;
2623  vector<int> noDau;
2624  vector< pair<int,int> > noMotDau;
2625  if (checkHistory) {
2626 
2627  // Loop through the event and check that there are beam particles.
2628  bool hasBeams = false;
2629  for (int i = 0; i < event.size(); ++i) {
2630  int status = event[i].status();
2631  if (abs(status) == 12) hasBeams = true;
2632 
2633  // Check that mother and daughter lists not empty where not expected to.
2634  vector<int> mList = event[i].motherList();
2635  vector<int> dList = event[i].daughterList();
2636  if (mList.size() == 0 && abs(status) != 11 && abs(status) != 12)
2637  noMot.push_back(i);
2638  if (dList.size() == 0 && status < 0 && status != -11)
2639  noDau.push_back(i);
2640 
2641  // Check that the particle appears in the daughters list of each mother.
2642  for (int j = 0; j < int(mList.size()); ++j) {
2643  if ( event[mList[j]].daughter1() <= i
2644  && event[mList[j]].daughter2() >= i ) continue;
2645  vector<int> dmList = event[mList[j]].daughterList();
2646  bool foundMatch = false;
2647  for (int k = 0; k < int(dmList.size()); ++k)
2648  if (dmList[k] == i) {
2649  foundMatch = true;
2650  break;
2651  }
2652  if (!hasBeams && mList.size() == 1 && mList[0] == 0) foundMatch = true;
2653  if (!foundMatch) {
2654  bool oldPair = false;
2655  for (int k = 0; k < int(noMotDau.size()); ++k)
2656  if (noMotDau[k].first == mList[j] && noMotDau[k].second == i) {
2657  oldPair = true;
2658  break;
2659  }
2660  if (!oldPair) noMotDau.push_back( make_pair( mList[j], i) );
2661  }
2662  }
2663 
2664  // Check that the particle appears in the mothers list of each daughter.
2665  for (int j = 0; j < int(dList.size()); ++j) {
2666  if ( event[dList[j]].statusAbs() > 80
2667  && event[dList[j]].statusAbs() < 90
2668  && event[dList[j]].mother1() <= i
2669  && event[dList[j]].mother2() >= i) continue;
2670  vector<int> mdList = event[dList[j]].motherList();
2671  bool foundMatch = false;
2672  for (int k = 0; k < int(mdList.size()); ++k)
2673  if (mdList[k] == i) {
2674  foundMatch = true;
2675  break;
2676  }
2677  if (!foundMatch) {
2678  bool oldPair = false;
2679  for (int k = 0; k < int(noMotDau.size()); ++k)
2680  if (noMotDau[k].first == i && noMotDau[k].second == dList[j]) {
2681  oldPair = true;
2682  break;
2683  }
2684  if (!oldPair) noMotDau.push_back( make_pair( i, dList[j]) );
2685  }
2686  }
2687  }
2688 
2689  // Warn if any errors were found.
2690  if (noMot.size() > 0 || noDau.size() > 0 || noMotDau.size() > 0) {
2691  info.errorMsg("Error in Pythia::check: "
2692  "mismatch in daughter and mother lists");
2693  physical = false;
2694  listHistory = true;
2695  }
2696  }
2697 
2698  // Done for sensible events.
2699  if (physical) return true;
2700 
2701  // Print (the first few) flawed events: local info.
2702  if (nErrEvent < nErrList) {
2703  cout << "\n PYTHIA erroneous event info: \n";
2704  if (iErrId.size() > 0) {
2705  cout << " unknown particle codes in lines ";
2706  for (int i = 0; i < int(iErrId.size()); ++i)
2707  cout << iErrId[i] << " ";
2708  cout << "\n";
2709  }
2710  if (iErrCol.size() > 0) {
2711  cout << " incorrect colour assignments in lines ";
2712  for (int i = 0; i < int(iErrCol.size()); ++i)
2713  cout << iErrCol[i] << " ";
2714  cout << "\n";
2715  }
2716  if (iErrEpm.size() > 0) {
2717  cout << " mismatch between energy/momentum/mass in lines ";
2718  for (int i = 0; i < int(iErrEpm.size()); ++i)
2719  cout << iErrEpm[i] << " ";
2720  cout << "\n";
2721  }
2722  if (iErrNan.size() > 0) {
2723  cout << " not-a-number energy/momentum/mass in lines ";
2724  for (int i = 0; i < int(iErrNan.size()); ++i)
2725  cout << iErrNan[i] << " ";
2726  cout << "\n";
2727  }
2728  if (iErrNanVtx.size() > 0) {
2729  cout << " not-a-number vertex/lifetime in lines ";
2730  for (int i = 0; i < int(iErrNanVtx.size()); ++i)
2731  cout << iErrNanVtx[i] << " ";
2732  cout << "\n";
2733  }
2734  if (epDev > epTolErr * eLab) cout << scientific << setprecision(3)
2735  << " total energy-momentum non-conservation = " << epDev << "\n";
2736  if (abs(chargeSum) > 0.1) cout << fixed << setprecision(2)
2737  << " total charge non-conservation = " << chargeSum << "\n";
2738  if (noMot.size() > 0) {
2739  cout << " missing mothers for particles ";
2740  for (int i = 0; i < int(noMot.size()); ++i) cout << noMot[i] << " ";
2741  cout << "\n";
2742  }
2743  if (noDau.size() > 0) {
2744  cout << " missing daughters for particles ";
2745  for (int i = 0; i < int(noDau.size()); ++i) cout << noDau[i] << " ";
2746  cout << "\n";
2747  }
2748  if (noMotDau.size() > 0) {
2749  cout << " inconsistent history for (mother,daughter) pairs ";
2750  for (int i = 0; i < int(noMotDau.size()); ++i)
2751  cout << "(" << noMotDau[i].first << "," << noMotDau[i].second << ") ";
2752  cout << "\n";
2753  }
2754 
2755  // Print (the first few) flawed events: standard listings.
2756  info.list();
2757  event.list(listVertices, listHistory);
2758  if (listSystems) partonSystems.list();
2759  if (listBeams) beamA.list();
2760  if (listBeams) beamB.list();
2761  }
2762 
2763  // Update error counter. Done also for flawed event.
2764  ++nErrEvent;
2765  return false;
2766 
2767 }
2768 
2769 //--------------------------------------------------------------------------
2770 
2771 // Routine to set up a PDF pointer.
2772 
2773 PDF* Pythia::getPDFPtr(int idIn, int sequence, string beam, bool resolved) {
2774 
2775  // Temporary pointer to be returned.
2776  PDF* tempPDFPtr = 0;
2777 
2778  // One option is to treat a Pomeron like a pi0.
2779  if (idIn == 990 && settings.word("PDF:PomSet") == "2") idIn = 111;
2780 
2781  // Proton beam, normal or hard choice. Also used for neutron.
2782  if (abs(idIn) == 2212 || abs(idIn) == 2112) {
2783  string pWord = settings.word("PDF:p"
2784  + string(sequence == 1 ? "" : "Hard") + "Set" + beam);
2785  if (pWord == "void" && sequence != 1 && beam == "B")
2786  pWord = settings.word("PDF:pHardSet");
2787  if (pWord == "void") pWord = settings.word("PDF:pSet");
2788  istringstream pStream(pWord);
2789  int pSet = 0;
2790  pStream >> pSet;
2791 
2792  // Use internal LHAgrid1 implementation for LHAPDF6 files.
2793  if (pSet == 0 && pWord.length() > 9
2794  && toLower(pWord).substr(0,9) == "lhagrid1:")
2795  tempPDFPtr = new LHAGrid1(idIn, pWord, xmlPath, &info);
2796 
2797  // Use sets from LHAPDF.
2798  else if (pSet == 0) tempPDFPtr = new LHAPDF(idIn, pWord, &info);
2799 
2800  // Use internal sets.
2801  else if (pSet == 1) tempPDFPtr = new GRV94L(idIn);
2802  else if (pSet == 2) tempPDFPtr = new CTEQ5L(idIn);
2803  else if (pSet <= 6)
2804  tempPDFPtr = new MSTWpdf(idIn, pSet - 2, xmlPath, &info);
2805  else if (pSet <= 12)
2806  tempPDFPtr = new CTEQ6pdf(idIn, pSet - 6, 1., xmlPath, &info);
2807  else if (pSet <= 16)
2808  tempPDFPtr = new NNPDF(idIn, pSet - 12, xmlPath, &info);
2809  else if (pSet <= 21)
2810  tempPDFPtr = new LHAGrid1(idIn, pWord, xmlPath, &info);
2811  else tempPDFPtr = 0;
2812  }
2813 
2814  // Pion beam (or, in one option, Pomeron beam).
2815  else if (abs(idIn) == 211 || idIn == 111) {
2816  string piWord = settings.word("PDF:piSet" + beam);
2817  if (piWord == "void" && beam == "B") piWord = settings.word("PDF:piSet");
2818  istringstream piStream(piWord);
2819  int piSet = 0;
2820  piStream >> piSet;
2821 
2822  // If VMD process then scale PDF accordingly:
2823  // f_a^VMD = alphaEM * (1/f_rho^2 + 1/f_omega^2 + 1/f_phi^2) * f_a^pi0.
2824  double rescale = (doVMDsideA || doVMDsideB) ? 0.00402068 : 1.;
2825 
2826  // Use internal LHAgrid1 implementation for LHAPDF6 files.
2827  if (piSet == 0 && piWord.length() > 9
2828  && toLower(piWord).substr(0,9) == "lhagrid1:")
2829  tempPDFPtr = new LHAGrid1(idIn, piWord, xmlPath, &info);
2830 
2831  // Use sets from LHAPDF.
2832  else if (piSet == 0) tempPDFPtr = new LHAPDF(idIn, piWord, &info);
2833 
2834  // Use internal set.
2835  else if (piSet == 1) tempPDFPtr = new GRVpiL(idIn, rescale);
2836  else tempPDFPtr = 0;
2837  }
2838 
2839  // Pomeron beam, if not treated like a pi0 beam.
2840  else if (idIn == 990) {
2841  string pomWord = settings.word("PDF:PomSet");
2842  double rescale = settings.parm("PDF:PomRescale");
2843  istringstream pomStream(pomWord);
2844  int pomSet = 0;
2845  pomStream >> pomSet;
2846 
2847  // Use internal LHAgrid1 implementation for LHAPDF6 files.
2848  if (pomSet == 0 && pomWord.length() > 9
2849  && toLower(pomWord).substr(0,9) == "lhagrid1:")
2850  tempPDFPtr = new LHAGrid1(idIn, pomWord, xmlPath, &info);
2851 
2852  // Use sets from LHAPDF.
2853  else if (pomSet == 0) tempPDFPtr = new LHAPDF(idIn, pomWord, &info);
2854 
2855  // A generic Q2-independent parametrization.
2856  else if (pomSet == 1) {
2857  double gluonA = settings.parm("PDF:PomGluonA");
2858  double gluonB = settings.parm("PDF:PomGluonB");
2859  double quarkA = settings.parm("PDF:PomQuarkA");
2860  double quarkB = settings.parm("PDF:PomQuarkB");
2861  double quarkFrac = settings.parm("PDF:PomQuarkFrac");
2862  double strangeSupp = settings.parm("PDF:PomStrangeSupp");
2863  tempPDFPtr = new PomFix( 990, gluonA, gluonB, quarkA, quarkB,
2864  quarkFrac, strangeSupp);
2865  }
2866 
2867  // The H1 Q2-dependent parametrizations. Initialization requires files.
2868  else if (pomSet == 3 || pomSet == 4)
2869  tempPDFPtr = new PomH1FitAB( 990, pomSet - 2, rescale, xmlPath, &info);
2870  else if (pomSet == 5)
2871  tempPDFPtr = new PomH1Jets( 990, 1, rescale, xmlPath, &info);
2872  else if (pomSet == 6)
2873  tempPDFPtr = new PomH1FitAB( 990, 3, rescale, xmlPath, &info);
2874  // The parametrizations of Alvero, Collins, Terron and Whitmore.
2875  else if (pomSet > 6 && pomSet < 11) {
2876  tempPDFPtr = new CTEQ6pdf( 990, pomSet + 4, rescale, xmlPath, &info);
2877  info.errorMsg("Warning: Pomeron flux parameters forced for ACTW PDFs");
2878  settings.mode("SigmaDiffractive:PomFlux", 4);
2879  double pomFluxEps = (pomSet == 10) ? 0.19 : 0.14;
2880  settings.parm("SigmaDiffractive:PomFluxEpsilon", pomFluxEps);
2881  settings.parm("SigmaDiffractive:PomFluxAlphaPrime", 0.25);
2882  }
2883  else if (pomSet == 11 )
2884  tempPDFPtr = new PomHISASD(990, getPDFPtr(2212), settings, &info);
2885  else if (pomSet == 12 || pomSet == 13)
2886  tempPDFPtr = new LHAGrid1(idIn, "1" + pomWord, xmlPath, &info);
2887  else tempPDFPtr = 0;
2888  }
2889 
2890  // Set up nuclear PDFs.
2891  else if (idIn > 100000000) {
2892 
2893  // Which nPDF set to use.
2894  int nPDFSet = (beam == "B") ? settings.mode("PDF:nPDFSetB")
2895  : settings.mode("PDF:nPDFSetA");
2896 
2897  // Temporary pointer for storing proton PDF pointer.
2898  PDF* tempProtonPDFPtr = (beam == "B") ? pdfHardBPtr : pdfHardAPtr;
2899  if (nPDFSet == 0) tempPDFPtr = new Isospin(idIn, tempProtonPDFPtr);
2900  else if (nPDFSet == 1 || nPDFSet == 2) tempPDFPtr = new EPS09(idIn,
2901  nPDFSet, 1, xmlPath, tempProtonPDFPtr, &info);
2902  else if (nPDFSet == 3) tempPDFPtr = new EPPS16(idIn, 1, xmlPath,
2903  tempProtonPDFPtr, &info);
2904  else tempPDFPtr = 0;
2905  }
2906 
2907  // Photon beam, either point-like (unresolved) or resolved.
2908  else if (abs(idIn) == 22) {
2909 
2910  // For unresolved beam use the point-like PDF.
2911  if (!resolved) {
2912  tempPDFPtr = new GammaPoint(idIn);
2913  } else {
2914  int gammaSet = settings.mode("PDF:GammaSet");
2915 
2916  // Point-like beam if unresolved photons.
2917  bool beamAisPoint = ( !beamAisResGamma && !beamAhasResGamma );
2918  bool beamBisPoint = ( !beamBisResGamma && !beamBhasResGamma );
2919  bool beamIsPoint = ( beamAisPoint && !(beam == "B") )
2920  || ( beamBisPoint && (beam == "B") );
2921 
2922  // Use different PDFs for hard process.
2923  if ( sequence == 2) {
2924 
2925  // Find the name or number of the hard PDF set.
2926  string gmWord = settings.word("PDF:GammaHardSet");
2927  int gmSet = 0;
2928  if (gmWord == "void") gmSet = settings.mode("PDF:GammaSet");
2929  else {
2930  istringstream gmStream(gmWord);
2931  gmStream >> gmSet;
2932  }
2933 
2934  // Use sets from LHAPDF. Only available for hard processes.
2935  if (gmSet == 0 && !beamIsPoint) {
2936  tempPDFPtr = new LHAPDF(idIn, gmWord, &info);
2937  return tempPDFPtr;
2938  }
2939 
2940  // Or set up an internal set (though currently only one).
2941  gammaSet = gmSet;
2942  }
2943 
2944  // Set up the PDF.
2945  if (beamIsPoint) tempPDFPtr = new GammaPoint(idIn);
2946  else if (gammaSet == 1) tempPDFPtr = new CJKL(idIn, &rndm);
2947  else tempPDFPtr = 0;
2948  }
2949  }
2950 
2951  // Lepton beam: neutrino, resolved charged lepton or unresolved ditto.
2952  // Also photon inside lepton PDFs.
2953  else if (abs(idIn) > 10 && abs(idIn) < 17) {
2954  if (abs(idIn)%2 == 0) tempPDFPtr = new NeutrinoPoint(idIn);
2955 
2956  // Set up resolved photon inside lepton for beam A.
2957  if ( beamAhasResGamma && !(beam == "B") && resolved ) {
2958 
2959  // Find the pre-set photon PDF, hard or normal.
2960  PDF* tempGammaPDFPtr = 0;
2961  if ( sequence == 2) tempGammaPDFPtr = pdfHardGamAPtr;
2962  else tempGammaPDFPtr = pdfGamAPtr;
2963 
2964  // Get the mass of lepton and maximum virtuality of the photon.
2965  double m2beam = pow2(particleData.m0(idIn));
2966  double Q2maxGamma = settings.parm("Photon:Q2max");
2967 
2968  // Initialize the gamma-inside-lepton PDFs with internal photon flux.
2969  if (settings.mode("PDF:lepton2gammaSet") == 1) {
2970  tempPDFPtr = new Lepton2gamma(idIn, m2beam, Q2maxGamma,
2971  tempGammaPDFPtr, &info, &rndm);
2972 
2973  // Initialize the gamma-inside-lepton PDFs with external photon flux.
2974  // Requires that the pointer to the flux set.
2975  } else if ( settings.mode("PDF:lepton2gammaSet") == 2 ) {
2976  PDF* tempGammaFluxPtr = pdfGamFluxAPtr;
2977  if ( tempGammaFluxPtr != 0) tempPDFPtr = new EPAexternal(idIn, m2beam,
2978  tempGammaFluxPtr, tempGammaPDFPtr, &settings, &info, &rndm );
2979  else {
2980  tempPDFPtr = 0;
2981  info.errorMsg("Error in Pythia::getPDFPtr: "
2982  "No external photon flux provided with PDF:lepton2gammaSet == 2");
2983  }
2984  } else tempPDFPtr = 0;
2985 
2986  // Set up resolved photon inside lepton for beam B.
2987  } else if ( beamBhasResGamma && (beam == "B") && resolved ) {
2988 
2989  // Find the pre-set photon PDF, hard or normal.
2990  PDF* tempGammaPDFPtr = 0;
2991  if ( sequence == 2) tempGammaPDFPtr = pdfHardGamBPtr;
2992  else tempGammaPDFPtr = pdfGamBPtr;
2993 
2994  // Get the mass of lepton and maximum virtuality of the photon.
2995  double m2beam = pow2(particleData.m0(idIn));
2996  double Q2maxGamma = settings.parm("Photon:Q2max");
2997 
2998  // Initialize the gamma-inside-lepton PDFs with internal photon flux.
2999  if (settings.mode("PDF:lepton2gammaSet") == 1) {
3000  tempPDFPtr = new Lepton2gamma(idIn, m2beam, Q2maxGamma,
3001  tempGammaPDFPtr, &info, &rndm);
3002 
3003  // Initialize the gamma-inside-lepton PDFs with external photon flux.
3004  } else if ( settings.mode("PDF:lepton2gammaSet") == 2 ) {
3005  PDF* tempGammaFluxPtr = pdfGamFluxBPtr;
3006  if ( tempGammaFluxPtr != 0) tempPDFPtr = new EPAexternal(idIn, m2beam,
3007  tempGammaFluxPtr, tempGammaPDFPtr, &settings, &info, &rndm );
3008  else {
3009  tempPDFPtr = 0;
3010  info.errorMsg("Error in Pythia::getPDFPtr: "
3011  "No external photon flux provided with PDF:lepton2gammaSet == 2");
3012  }
3013  } else tempPDFPtr = 0;
3014 
3015  // Usual lepton PDFs.
3016  } else if (settings.flag("PDF:lepton")) {
3017  double m2beam = pow2(particleData.m0(idIn));
3018  double Q2maxGamma = settings.parm("Photon:Q2max");
3019  if (settings.mode("PDF:lepton2gammaSet") == 1 ) {
3020  tempPDFPtr = new Lepton(idIn, Q2maxGamma, &info, &rndm);
3021 
3022  // External photon flux for direct-photon processes.
3023  } else if (settings.mode("PDF:lepton2gammaSet") == 2 ) {
3024  PDF* tempGammaPDFPtr = 0;
3025  PDF* tempGammaFluxPtr = (beam == "B") ?
3026  pdfGamFluxBPtr : pdfGamFluxAPtr;
3027  if ( tempGammaFluxPtr != 0) tempPDFPtr = new EPAexternal(idIn, m2beam,
3028  tempGammaFluxPtr, tempGammaPDFPtr, &settings, &info, &rndm );
3029  else {
3030  tempPDFPtr = 0;
3031  info.errorMsg("Error in Pythia::getPDFPtr: "
3032  "No external photon flux provided with PDF:lepton2gammaSet == 2");
3033  }
3034  } else tempPDFPtr = 0;
3035  }
3036  else tempPDFPtr = new LeptonPoint(idIn);
3037 
3038  // Dark matter beam set up as pointlike lepton.
3039  } else if (abs(idIn) > 50 && abs(idIn) < 60) {
3040  tempPDFPtr = new LeptonPoint(idIn);
3041  }
3042 
3043  // Optionally allow extrapolation beyond x and Q2 limits.
3044  if (tempPDFPtr)
3045  tempPDFPtr->setExtrapolate( settings.flag("PDF:extrapolate") );
3046 
3047  // Done.
3048  return tempPDFPtr;
3049 }
3050 
3051 //==========================================================================
3052 
3053 } // end namespace Pythia8
bool setHIUserHooksPtr(HIUserHooks *userHooksPtrIn)
Possibility to pass in pointer for special heavy ion user hooks.
Definition: HeavyIons.h:69
static void addSpecialSettings(Settings &settings)
Definition: HeavyIons.cc:24
Definition: beam.h:43
Definition: AgUStep.h:26
virtual bool init()=0
virtual bool next()=0
static bool isHeavyIon(Settings &settings)
Definition: HeavyIons.cc:258