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) 2012 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the Pythia class.
7 
8 #include "Pythia.h"
9 
10 // Access time information.
11 #include <ctime>
12 
13 // Allow string and character manipulation.
14 #include <cctype>
15 
16 namespace Pythia8 {
17 
18 //==========================================================================
19 
20 // The Pythia class.
21 
22 //--------------------------------------------------------------------------
23 
24 // Constants: could be changed here if desired, but normally should not.
25 // These are of technical nature, as described for each.
26 
27 // Maximum number of tries to produce parton level from given input.
28 const int Pythia::NTRY = 10;
29 
30 // Negative integer to denote that no subrun has been set.
31 const int Pythia::SUBRUNDEFAULT = -999;
32 
33 //--------------------------------------------------------------------------
34 
35 // Constructor.
36 
37 Pythia::Pythia(string xmlDir) {
38 
39  // Initial values for pointers to PDF's.
40  useNewPdfA = false;
41  useNewPdfB = false;
42  useNewPdfHard = false;
43  useNewPdfPomA = false;
44  useNewPdfPomB = false;
45  pdfAPtr = 0;
46  pdfBPtr = 0;
47  pdfHardAPtr = 0;
48  pdfHardBPtr = 0;
49  pdfPomAPtr = 0;
50  pdfPomBPtr = 0;
51 
52  // Initial values for pointers to Les Houches Event objects.
53  doLHA = false;
54  useNewLHA = false;
55  lhaUpPtr = 0;
56 
57  //Initial value for couplings pointer
58  couplingsPtr = &couplings;
59 
60  // Initial value for pointer to external decay handler.
61  decayHandlePtr = 0;
62 
63  // Initial value for pointer to user hooks.
64  userHooksPtr = 0;
65 
66  // Initial value for pointer to merging hooks.
67  doMerging = false;
68  hasMergingHooks = false;
69  hasOwnMergingHooks = false;
70  mergingHooksPtr = 0;
71 
72  // Initial value for pointer to beam shape.
73  useNewBeamShape = false;
74  beamShapePtr = 0;
75 
76  // Initial values for pointers to timelike and spacelike showers.
77  useNewTimes = false;
78  useNewSpace = false;
79  timesDecPtr = 0;
80  timesPtr = 0;
81  spacePtr = 0;
82 
83  // Find path to data files, i.e. xmldoc directory location.
84  // Environment variable takes precedence, else use constructor input.
85  xmlPath = "";
86  const char* PYTHIA8DATA = "PYTHIA8DATA";
87  char* envPath = getenv(PYTHIA8DATA);
88  if (envPath != 0 && *envPath != '\0') {
89  int i = 0;
90  while (*(envPath+i) != '\0') xmlPath += *(envPath+(i++));
91  }
92  else xmlPath = xmlDir;
93  if (xmlPath[ xmlPath.length() - 1 ] != '/') xmlPath += "/";
94 
95  // Read in files with all flags, modes, parms and words.
96  settings.initPtr( &info);
97  string initFile = xmlPath + "Index.xml";
98  isConstructed = settings.init( initFile);
99  if (!isConstructed) {
100  info.errorMsg("Abort from Pythia::Pythia: settings unavailable");
101  return;
102  }
103 
104  // Read in files with all particle data.
105  particleData.initPtr( &info, &settings, &rndm, couplingsPtr);
106  string dataFile = xmlPath + "ParticleData.xml";
107  isConstructed = particleData.init( dataFile);
108  if (!isConstructed) {
109  info.errorMsg("Abort from Pythia::Pythia: particle data unavailable");
110  return;
111  }
112 
113  // Write the Pythia banner to output.
114  banner();
115 
116  // Not initialized until at the end of the init() call.
117  isInit = false;
118  info.addCounter(0);
119 
120 }
121 
122 //--------------------------------------------------------------------------
123 
124 // Destructor.
125 
126 Pythia::~Pythia() {
127 
128  // Delete the PDF's created with new.
129  if (useNewPdfHard && pdfHardAPtr != pdfAPtr) delete pdfHardAPtr;
130  if (useNewPdfHard && pdfHardBPtr != pdfBPtr) delete pdfHardBPtr;
131  if (useNewPdfA) delete pdfAPtr;
132  if (useNewPdfB) delete pdfBPtr;
133  if (useNewPdfPomA) delete pdfPomAPtr;
134  if (useNewPdfPomB) delete pdfPomBPtr;
135 
136  // Delete the Les Houches object created with new.
137  if (useNewLHA) delete lhaUpPtr;
138 
139  // Delete the MergingHooks object created with new.
140  if (hasOwnMergingHooks) delete mergingHooksPtr;
141 
142  // Delete the BeamShape object created with new.
143  if (useNewBeamShape) delete beamShapePtr;
144 
145  // Delete the timelike and spacelike showers created with new.
146  if (useNewTimes) delete timesPtr;
147  if (useNewSpace) delete spacePtr;
148 
149 }
150 
151 //--------------------------------------------------------------------------
152 
153 // Read in one update for a setting or particle data from a single line.
154 
155 bool Pythia::readString(string line, bool warn) {
156 
157  // Check that constructor worked.
158  if (!isConstructed) return false;
159 
160  // If empty line then done.
161  if (line.find_first_not_of(" \n\t\v\b\r\f\a") == string::npos) return true;
162 
163  // If first character is not a letter/digit, then taken to be a comment.
164  int firstChar = line.find_first_not_of(" \n\t\v\b\r\f\a");
165  if (!isalnum(line[firstChar])) return true;
166 
167  // Send on particle data to the ParticleData database.
168  if (isdigit(line[firstChar]))
169  return particleData.readString(line, warn);
170 
171  // Everything else sent on to Settings.
172  return settings.readString(line, warn);
173 
174 }
175 
176 //--------------------------------------------------------------------------
177 
178 // Read in updates for settings or particle data from user-defined file.
179 
180 bool Pythia::readFile(string fileName, bool warn, int subrun) {
181 
182  // Check that constructor worked.
183  if (!isConstructed) return false;
184 
185  // Open file for reading.
186  const char* cstring = fileName.c_str();
187  ifstream is(cstring);
188  if (!is.good()) {
189  info.errorMsg("Error in Pythia::readFile: did not find file", fileName);
190  return false;
191  }
192 
193  // Hand over real work to next method.
194  return readFile( is, warn, subrun);
195 
196 }
197 
198 //--------------------------------------------------------------------------
199 
200 // Read in updates for settings or particle data
201 // from user-defined stream (or file).
202 
203 bool Pythia::readFile(istream& is, bool warn, int subrun) {
204 
205  // Check that constructor worked.
206  if (!isConstructed) return false;
207 
208  // Read in one line at a time.
209  string line;
210  bool accepted = true;
211  int subrunNow = SUBRUNDEFAULT;
212  while ( getline(is, line) ) {
213 
214  // Check whether entered new subrun.
215  int subrunLine = readSubrun( line, warn);
216  if (subrunLine >= 0) subrunNow = subrunLine;
217 
218  // Process the line if in correct subrun.
219  if ( (subrunNow == subrun || subrunNow == SUBRUNDEFAULT)
220  && !readString( line, warn) ) accepted = false;
221 
222  // Reached end of input file.
223  };
224  return accepted;
225 
226 }
227 
228 //--------------------------------------------------------------------------
229 
230 // Routine to pass in pointers to PDF's. Usage optional.
231 
232 bool Pythia::setPDFPtr( PDF* pdfAPtrIn, PDF* pdfBPtrIn, PDF* pdfHardAPtrIn,
233  PDF* pdfHardBPtrIn, PDF* pdfPomAPtrIn, PDF* pdfPomBPtrIn) {
234 
235  // Delete any PDF's created in a previous init call.
236  if (useNewPdfHard && pdfHardAPtr != pdfAPtr) delete pdfHardAPtr;
237  if (useNewPdfHard && pdfHardBPtr != pdfBPtr) delete pdfHardBPtr;
238  if (useNewPdfA) delete pdfAPtr;
239  if (useNewPdfB) delete pdfBPtr;
240  if (useNewPdfPomA) delete pdfPomAPtr;
241  if (useNewPdfPomB) delete pdfPomBPtr;
242 
243  // Reset pointers to be empty.
244  useNewPdfA = false;
245  useNewPdfB = false;
246  useNewPdfHard = false;
247  useNewPdfPomA = false;
248  useNewPdfPomB = false;
249  pdfAPtr = 0;
250  pdfBPtr = 0;
251  pdfHardAPtr = 0;
252  pdfHardBPtr = 0;
253  pdfPomAPtr = 0;
254  pdfPomBPtr = 0;
255 
256  // Switch off external PDF's by zero as input.
257  if (pdfAPtrIn == 0 && pdfBPtrIn == 0) return true;
258 
259  // The two PDF objects cannot be one and the same.
260  if (pdfAPtrIn == pdfBPtrIn) return false;
261 
262  // Save pointers.
263  pdfAPtr = pdfAPtrIn;
264  pdfBPtr = pdfBPtrIn;
265 
266  // By default same pointers for hard-process PDF's.
267  pdfHardAPtr = pdfAPtrIn;
268  pdfHardBPtr = pdfBPtrIn;
269 
270  // Optionally allow separate pointers for hard process.
271  if (pdfHardAPtrIn != 0 && pdfHardBPtrIn != 0) {
272  if (pdfHardAPtrIn == pdfHardBPtrIn) return false;
273  pdfHardAPtr = pdfHardAPtrIn;
274  pdfHardBPtr = pdfHardBPtrIn;
275  }
276 
277  // Optionally allow pointers for Pomerons in the proton.
278  if (pdfPomAPtrIn != 0 && pdfPomBPtrIn != 0) {
279  if (pdfPomAPtrIn == pdfPomBPtrIn) return false;
280  pdfPomAPtr = pdfPomAPtrIn;
281  pdfPomBPtr = pdfPomBPtrIn;
282  }
283 
284  // Done.
285  return true;
286 }
287 
288 //--------------------------------------------------------------------------
289 
290 // Routine to initialize with the variable values of the Beams kind.
291 
292 bool Pythia::init() {
293 
294  // Check that constructor worked.
295  isInit = false;
296  if (!isConstructed) {
297  info.errorMsg("Abort from Pythia::init: constructur "
298  "initialization failed");
299  return false;
300  }
301 
302  // Begin initialization. Find which frame type to use.
303  info.addCounter(1);
304  frameType = mode("Beams:frameType");
305 
306  // Initialization with internal processes: read in and set values.
307  if (frameType < 4 ) {
308  doLHA = false;
309  boostType = frameType;
310  idA = mode("Beams:idA");
311  idB = mode("Beams:idB");
312  eCM = parm("Beams:eCM");
313  eA = parm("Beams:eA");
314  eB = parm("Beams:eB");
315  pxA = parm("Beams:pxA");
316  pyA = parm("Beams:pyA");
317  pzA = parm("Beams:pzA");
318  pxB = parm("Beams:pxB");
319  pyB = parm("Beams:pyB");
320  pzB = parm("Beams:pzB");
321 
322  // Initialization with a Les Houches Event File or an LHAup object.
323  } else {
324  doLHA = true;
325  boostType = 2;
326  string lhef = word("Beams:LHEF");
327  bool skipInit = flag("Beams:newLHEFsameInit");
328 
329  // For file input: delete any old LHAup object and create new one.
330  if (frameType == 4) {
331  if (useNewLHA) delete lhaUpPtr;
332  const char* cstring = lhef.c_str();
333  lhaUpPtr = new LHAupLHEF(cstring);
334  useNewLHA = true;
335 
336  // Check that file was properly opened.
337  if (!lhaUpPtr->fileFound()) {
338  info.errorMsg("Abort from Pythia::init: "
339  "Les Houches Event File not found");
340  return false;
341  }
342 
343  // For object input: at least check that not null pointer.
344  } else {
345  if (lhaUpPtr == 0) {
346  info.errorMsg("Abort from Pythia::init: "
347  "LHAup object not found");
348  return false;
349  }
350  }
351 
352  // Send in pointer to info. Store or replace LHA pointer in other classes.
353  lhaUpPtr->setPtr( &info);
354  processLevel.setLHAPtr( lhaUpPtr);
355 
356  // If second time around, only with new file, then simplify.
357  if (skipInit) {
358  isInit = true;
359  info.addCounter(2);
360  return true;
361  }
362 
363  // Set up values related to user hooks.
364  if (frameType == 4) {
365  doUserMerging = settings.flag("Merging:doUserMerging");
366  doMGMerging = settings.flag("Merging:doMGMerging");
367  doKTMerging = settings.flag("Merging:doKTMerging");
368  doMerging = doUserMerging || doMGMerging || doKTMerging;
369  if (!doUserMerging){
370  if (mergingHooksPtr) delete mergingHooksPtr;
371  mergingHooksPtr = new MergingHooks();
372  hasOwnMergingHooks = true;
373  }
374  hasMergingHooks = (mergingHooksPtr > 0);
375  if (hasMergingHooks) mergingHooksPtr->setLHEInputFile( lhef);
376  }
377 
378  // Set LHAinit information (in some external program).
379  if (!lhaUpPtr->setInit()) {
380  info.errorMsg("Abort from Pythia::init: "
381  "Les Houches initialization failed");
382  return false;
383  }
384 
385  // Extract beams from values set in an LHAinit object.
386  idA = lhaUpPtr->idBeamA();
387  idB = lhaUpPtr->idBeamB();
388  eA = lhaUpPtr->eBeamA();
389  eB = lhaUpPtr->eBeamB();
390  }
391 
392  // Back to common initialization. Reset error counters.
393  nErrEvent = 0;
394  info.errorReset();
395  info.setTooLowPTmin(false);
396  info.setSigma( 0, 0, 0, 0., 0., 0.);
397 
398  // Initialize data members extracted from database.
399  doProcessLevel = settings.flag("ProcessLevel:all");
400  doPartonLevel = settings.flag("PartonLevel:all");
401  doHadronLevel = settings.flag("HadronLevel:all");
402  doDiffraction = settings.flag("SoftQCD:all")
403  || settings.flag("SoftQCD:singleDiffractive")
404  || settings.flag("SoftQCD:doubleDiffractive");
405  decayRHadrons = settings.flag("RHadrons:allowDecay");
406  doMomentumSpread = settings.flag("Beams:allowMomentumSpread");
407  doVertexSpread = settings.flag("Beams:allowVertexSpread");
408  abortIfVeto = settings.flag("Check:abortIfVeto");
409  checkEvent = settings.flag("Check:event");
410  nErrList = settings.mode("Check:nErrList");
411  epTolErr = settings.parm("Check:epTolErr");
412  epTolWarn = settings.parm("Check:epTolWarn");
413 
414  // Initialise merging hooks.
415  if (hasMergingHooks && doMerging)
416  mergingHooksPtr->init( settings, &info, &particleData );
417 
418  // Initialize the random number generator.
419  if ( settings.flag("Random:setSeed") )
420  rndm.init( settings.mode("Random:seed") );
421 
422  // Check that combinations of settings are allowed; change if not.
423  checkSettings();
424 
425  // Initialize couplings (needed to initialize resonances).
426  // Check if SUSY couplings need to be read in
427  if( !initSLHA()) info.errorMsg("Error in Pythia::init: "
428  "Could not read SLHA file");
429  if (couplingsPtr->isSUSY) {
430  // Initialize the SM and SUSY.
431  coupSUSY.init( settings, &rndm);
432  coupSUSY.initSUSY(&slha, &settings, &particleData);
433  couplingsPtr = (Couplings *) &coupSUSY;
434  } else {
435  // Initialize the SM couplings.
436  couplingsPtr->init( settings, &rndm);
437  }
438 
439  // Reset couplingsPtr to the correct place.
440  particleData.initPtr( &info, &settings, &rndm, couplingsPtr);
441 
442  // Set headers to distinguish the two event listing kinds.
443  int startColTag = settings.mode("Event:startColTag");
444  process.init("(hard process)", &particleData, startColTag);
445  event.init("(complete event)", &particleData, startColTag);
446 
447  // Final setup stage of particle data, notably resonance widths.
448  particleData.initWidths( resonancePtrs);
449 
450  // Set up R-hadrons particle data, where relevant.
451  rHadrons.init( &info, settings, &particleData, &rndm);
452 
453  // Set up values related to user hooks.
454  hasUserHooks = (userHooksPtr > 0);
455  doVetoProcess = (hasUserHooks)
456  ? userHooksPtr->canVetoProcessLevel() : false;
457  doVetoPartons = (hasUserHooks)
458  ? userHooksPtr->canVetoPartonLevel() : false;
459  if (hasUserHooks) userHooksPtr->initPtr( &info, &settings, &particleData,
460  &rndm, &beamA, &beamB, &beamPomA, &beamPomB, couplingsPtr, &partonSystems,
461  &sigmaTot);
462 
463  // Set up objects for timelike and spacelike showers.
464  if (timesDecPtr == 0 || timesPtr == 0) {
465  TimeShower* timesNow = new TimeShower();
466  if (timesDecPtr == 0) timesDecPtr = timesNow;
467  if (timesPtr == 0) timesPtr = timesNow;
468  useNewTimes = true;
469  }
470  if (spacePtr == 0) {
471  spacePtr = new SpaceShower();
472  useNewSpace = true;
473  }
474 
475  // Initialize showers, especially for simple showers in decays.
476  timesPtr->initPtr( &info, &settings, &particleData, &rndm, couplingsPtr,
477  &partonSystems, userHooksPtr);
478  timesDecPtr->initPtr( &info, &settings, &particleData, &rndm, couplingsPtr,
479  &partonSystems, userHooksPtr);
480  spacePtr->initPtr( &info, &settings, &particleData, &rndm,
481  &partonSystems, userHooksPtr);
482  timesDecPtr->init( 0, 0);
483 
484  // Set up values related to beam shape.
485  if (beamShapePtr == 0) {
486  beamShapePtr = new BeamShape();
487  useNewBeamShape = true;
488  }
489  beamShapePtr->init( settings, &rndm);
490 
491  // Check that beams and beam combination can be handled.
492  if (!checkBeams()) {
493  info.errorMsg("Abort from Pythia::init: "
494  "checkBeams initialization failed");
495  return false;
496  }
497 
498  // Do not set up beam kinematics when no process level.
499  if (!doProcessLevel) boostType = 1;
500  else {
501 
502  // Set up beam kinematics.
503  if (!initKinematics()) {
504  info.errorMsg("Abort from Pythia::init: "
505  "kinematics initialization failed");
506  return false;
507  }
508 
509  // Set up pointers to PDFs.
510  if (!initPDFs()) {
511  info.errorMsg("Abort from Pythia::init: PDF initialization failed");
512  return false;
513  }
514 
515  // Set up the two beams and the common remnant system.
516  StringFlav* flavSelPtr = hadronLevel.getStringFlavPtr();
517  beamA.init( idA, pzAcm, eA, mA, &info, settings, &particleData, &rndm,
518  pdfAPtr, pdfHardAPtr, isUnresolvedA, flavSelPtr);
519  beamB.init( idB, pzBcm, eB, mB, &info, settings, &particleData, &rndm,
520  pdfBPtr, pdfHardBPtr, isUnresolvedB, flavSelPtr);
521 
522  // Optionally set up new alternative beams for these Pomerons.
523  if ( doDiffraction) {
524  beamPomA.init( 990, 0.5 * eCM, 0.5 * eCM, 0., &info, settings,
525  &particleData, &rndm, pdfPomAPtr, pdfPomAPtr, false, flavSelPtr);
526  beamPomB.init( 990, -0.5 * eCM, 0.5 * eCM, 0., &info, settings,
527  &particleData, &rndm, pdfPomBPtr, pdfPomBPtr, false, flavSelPtr);
528  }
529  }
530 
531  // Send info/pointers to process level for initialization.
532  if ( doProcessLevel && !processLevel.init( &info, settings, &particleData,
533  &rndm, &beamA, &beamB, couplingsPtr, &sigmaTot, doLHA, &slha, userHooksPtr,
534  sigmaPtrs) ) {
535  info.errorMsg("Abort from Pythia::init: "
536  "processLevel initialization failed");
537  return false;
538  }
539 
540  // Send info/pointers to parton level for initialization.
541  if ( doPartonLevel && doProcessLevel && !partonLevel.init( &info, settings,
542  &particleData, &rndm, &beamA, &beamB, &beamPomA, &beamPomB, couplingsPtr,
543  &partonSystems, &sigmaTot, timesDecPtr, timesPtr, spacePtr, &rHadrons,
544  userHooksPtr, mergingHooksPtr, false) ) {
545  info.errorMsg("Abort from Pythia::init: "
546  "partonLevel initialization failed" );
547  return false;
548  }
549 
550  // Send info/pointers to parton level for trial shower initialization.
551  if ( doMerging && !trialPartonLevel.init( &info, settings, &particleData,
552  &rndm, &beamA, &beamB, &beamPomA, &beamPomB, couplingsPtr, &partonSystems,
553  &sigmaTot, timesDecPtr, timesPtr, spacePtr, &rHadrons, NULL,
554  mergingHooksPtr, true) ) {
555  info.errorMsg("Abort from Pythia::init: "
556  "trialPartonLevel initialization failed");
557  return false;
558  }
559 
560  // Send info/pointers to hadron level for initialization.
561  // Note: forceHadronLevel() can come, so we must always initialize.
562  if ( !hadronLevel.init( &info, settings, &particleData, &rndm,
563  couplingsPtr, timesDecPtr, &rHadrons, decayHandlePtr,
564  handledParticles) ) {
565  info.errorMsg("Abort from Pythia::init: "
566  "hadronLevel initialization failed");
567  return false;
568  }
569 
570  // Optionally check particle data table for inconsistencies.
571  if ( settings.flag("Check:particleData") )
572  particleData.checkTable( settings.mode("Check:levelParticleData") );
573 
574  // Optionally show settings and particle data, changed or all.
575  bool showCS = settings.flag("Init:showChangedSettings");
576  bool showAS = settings.flag("Init:showAllSettings");
577  bool showCPD = settings.flag("Init:showChangedParticleData");
578  bool showCRD = settings.flag("Init:showChangedResonanceData");
579  bool showAPD = settings.flag("Init:showAllParticleData");
580  int show1PD = settings.mode("Init:showOneParticleData");
581  bool showPro = settings.flag("Init:showProcesses");
582  if (showCS) settings.listChanged();
583  if (showAS) settings.listAll();
584  if (show1PD > 0) particleData.list(show1PD);
585  if (showCPD) particleData.listChanged(showCRD);
586  if (showAPD) particleData.listAll();
587 
588  // Listing options for the next() routine.
589  nCount = settings.mode("Next:numberCount");
590  nShowLHA = settings.mode("Next:numberShowLHA");
591  nShowInfo = settings.mode("Next:numberShowInfo");
592  nShowProc = settings.mode("Next:numberShowProcess");
593  nShowEvt = settings.mode("Next:numberShowEvent");
594  showSaV = settings.flag("Next:showScaleAndVertex");
595  showMaD = settings.flag("Next:showMothersAndDaughters");
596 
597  // Succeeded.
598  isInit = true;
599  info.addCounter(2);
600  if (useNewLHA && showPro) lhaUpPtr->listInit();
601  return true;
602 
603 }
604 
605 //--------------------------------------------------------------------------
606 
607 // Routine to initialize with CM energy.
608 
609 bool Pythia::init( int idAin, int idBin, double eCMin) {
610 
611  // Set arguments in Settings database.
612  settings.mode("Beams:idA", idAin);
613  settings.mode("Beams:idB", idBin);
614  settings.mode("Beams:frameType", 1);
615  settings.parm("Beams:eCM", eCMin);
616 
617  // Send on to common initialization.
618  return init();
619 
620 }
621 
622 //--------------------------------------------------------------------------
623 
624 // Routine to initialize with two collinear beams, energies specified.
625 
626 bool Pythia::init( int idAin, int idBin, double eAin, double eBin) {
627 
628  // Set arguments in Settings database.
629  settings.mode("Beams:idA", idAin);
630  settings.mode("Beams:idB", idBin);
631  settings.mode("Beams:frameType", 2);
632  settings.parm("Beams:eA", eAin);
633  settings.parm("Beams:eB", eBin);
634 
635  // Send on to common initialization.
636  return init();
637 
638 }
639 
640 //--------------------------------------------------------------------------
641 
642 // Routine to initialize with two beams specified by three-momenta.
643 
644 bool Pythia::init( int idAin, int idBin, double pxAin, double pyAin,
645  double pzAin, double pxBin, double pyBin, double pzBin) {
646 
647  // Set arguments in Settings database.
648  settings.mode("Beams:idA", idAin);
649  settings.mode("Beams:idB", idBin);
650  settings.mode("Beams:frameType", 3);
651  settings.parm("Beams:pxA", pxAin);
652  settings.parm("Beams:pyA", pyAin);
653  settings.parm("Beams:pzA", pzAin);
654  settings.parm("Beams:pxB", pxBin);
655  settings.parm("Beams:pyB", pyBin);
656  settings.parm("Beams:pzB", pzBin);
657 
658  // Send on to common initialization.
659  return init();
660 
661 }
662 
663 //--------------------------------------------------------------------------
664 
665 // Routine to initialize when all info is given in a Les Houches Event File.
666 
667 bool Pythia::init( string LesHouchesEventFile, bool skipInit) {
668 
669  // Set arguments in Settings database.
670  settings.mode("Beams:frameType", 4);
671  settings.word("Beams:LHEF", LesHouchesEventFile);
672  settings.flag("Beams:newLHEFsameInit", skipInit);
673 
674  // Send on to common initialization.
675  return init();
676 
677 }
678 
679 //--------------------------------------------------------------------------
680 
681 // Routine to initialize when beam info is given in an LHAup object.
682 
683 bool Pythia::init( LHAup* lhaUpPtrIn) {
684 
685  // Store LHAup object and set arguments in Settings database.
686  setLHAupPtr( lhaUpPtrIn);
687  settings.mode("Beams:frameType", 5);
688 
689  // Send on to common initialization.
690  return init();
691 
692 }
693 
694 //--------------------------------------------------------------------------
695 
696 // Check that combinations of settings are allowed; change if not.
697 
698 void Pythia::checkSettings() {
699 
700  // Double rescattering not allowed if ISR or FSR.
701  if ((settings.flag("PartonLevel:ISR") || settings.flag("PartonLevel:FSR"))
702  && settings.flag("MultipartonInteractions:allowDoubleRescatter")) {
703  info.errorMsg("Warning in Pythia::checkSettings: "
704  "double rescattering switched off since showering is on");
705  settings.flag("MultipartonInteractions:allowDoubleRescatter", false);
706  }
707 
708 }
709 
710 //--------------------------------------------------------------------------
711 
712 // Check that beams and beam combination can be handled. Set up unresolved.
713 
714 bool Pythia::checkBeams() {
715 
716  // Absolute flavours. If not to do process level then no check needed.
717  int idAabs = abs(idA);
718  int idBabs = abs(idB);
719  if (!doProcessLevel) return true;
720 
721  // Neutrino beams always unresolved, charged lepton ones conditionally.
722  bool isLeptonA = (idAabs > 10 && idAabs < 17);
723  bool isLeptonB = (idBabs > 10 && idBabs < 17);
724  bool isUnresLep = !settings.flag("PDF:lepton");
725  isUnresolvedA = isLeptonA && (idAabs%2 == 0 || isUnresLep);
726  isUnresolvedB = isLeptonB && (idBabs%2 == 0 || isUnresLep);
727 
728  // Lepton-lepton collisions OK (including neutrinos) if both (un)resolved.
729  if (isLeptonA && isLeptonB && isUnresolvedA == isUnresolvedB) return true;
730 
731  // Hadron-hadron collisions OK, with Pomeron counted as hadron.
732  bool isHadronA = (idAabs == 2212) || (idA == 111) || (idAabs == 211)
733  || (idA == 990);
734  bool isHadronB = (idBabs == 2212) || (idA == 111)|| (idBabs == 211)
735  || (idB == 990);
736  if (isHadronA && isHadronB) return true;
737 
738  // If no case above then failed.
739  info.errorMsg("Error in Pythia::init: cannot handle this beam combination");
740  return false;
741 
742 }
743 
744 //--------------------------------------------------------------------------
745 
746 // Calculate kinematics at initialization. Store beam four-momenta.
747 
748 bool Pythia::initKinematics() {
749 
750  // Find masses. Initial guess that we are in CM frame.
751  mA = particleData.m0(idA);
752  mB = particleData.m0(idB);
753  betaZ = 0.;
754  gammaZ = 1.;
755 
756  // Collinear beams not in CM frame: find CM energy.
757  if (boostType == 2) {
758  eA = max(eA, mA);
759  eB = max(eB, mB);
760  pzA = sqrt(eA*eA - mA*mA);
761  pzB = -sqrt(eB*eB - mB*mB);
762  pAinit = Vec4( 0., 0., pzA, eA);
763  pBinit = Vec4( 0., 0., pzB, eB);
764  eCM = sqrt( pow2(eA + eB) - pow2(pzA + pzB) );
765 
766  // Find boost to rest frame.
767  betaZ = (pzA + pzB) / (eA + eB);
768  gammaZ = (eA + eB) / eCM;
769  if (abs(betaZ) < 1e-10) boostType = 1;
770  }
771 
772  // Completely general beam directions: find CM energy.
773  else if (boostType == 3) {
774  eA = sqrt( pxA*pxA + pyA*pyA + pzA*pzA + mA*mA);
775  eB = sqrt( pxB*pxB + pyB*pyB + pzB*pzB + mB*mB);
776  pAinit = Vec4( pxA, pyA, pzA, eA);
777  pBinit = Vec4( pxB, pyB, pzB, eB);
778  eCM = (pAinit + pBinit).mCalc();
779 
780  // Find boost+rotation needed to move from/to CM frame.
781  MfromCM.reset();
782  MfromCM.fromCMframe( pAinit, pBinit);
783  MtoCM = MfromCM;
784  MtoCM.invert();
785  }
786 
787  // Fail if CM energy below beam masses.
788  if (eCM < mA + mB) {
789  info.errorMsg("Error in Pythia::initKinematics: too low energy");
790  return false;
791  }
792 
793  // Set up CM-frame kinematics with beams along +-z axis.
794  pzAcm = 0.5 * sqrtpos( (eCM + mA + mB) * (eCM - mA - mB)
795  * (eCM - mA + mB) * (eCM + mA - mB) ) / eCM;
796  pzBcm = -pzAcm;
797  eA = sqrt(mA*mA + pzAcm*pzAcm);
798  eB = sqrt(mB*mB + pzBcm*pzBcm);
799 
800  // If in CM frame then store beam four-vectors (else already done above).
801  if (boostType != 2 && boostType != 3) {
802  pAinit = Vec4( 0., 0., pzAcm, eA);
803  pBinit = Vec4( 0., 0., pzBcm, eB);
804  }
805 
806  // Store main info for access in process generation.
807  info.setBeamA( idA, pzAcm, eA, mA);
808  info.setBeamB( idB, pzBcm, eB, mB);
809  info.setECM( eCM);
810 
811  // Must allow for generic boost+rotation when beam momentum spread.
812  if (doMomentumSpread) boostType = 3;
813 
814  // Done.
815  return true;
816 
817 }
818 
819 //--------------------------------------------------------------------------
820 
821 // Set up pointers to PDFs.
822 
823 bool Pythia::initPDFs() {
824 
825  // Delete any PDF's created in a previous init call.
826  if (useNewPdfHard) {
827  if (pdfHardAPtr != pdfAPtr) {
828  delete pdfHardAPtr;
829  pdfHardAPtr = 0;
830  }
831  if (pdfHardBPtr != pdfBPtr) {
832  delete pdfHardBPtr;
833  pdfHardBPtr = 0;
834  }
835  useNewPdfHard = false;
836  }
837  if (useNewPdfA) {
838  delete pdfAPtr;
839  useNewPdfA = false;
840  pdfAPtr = 0;
841  }
842  if (useNewPdfB) {
843  delete pdfBPtr;
844  useNewPdfB = false;
845  pdfBPtr = 0;
846  }
847  if (useNewPdfPomA) {
848  delete pdfPomAPtr;
849  useNewPdfPomA = false;
850  pdfPomAPtr = 0;
851  }
852  if (useNewPdfPomB) {
853  delete pdfPomBPtr;
854  useNewPdfPomB = false;
855  pdfPomBPtr = 0;
856  }
857 
858  // Set up the PDF's, if not already done.
859  if (pdfAPtr == 0) {
860  pdfAPtr = getPDFPtr(idA);
861  if (pdfAPtr == 0 || !pdfAPtr->isSetup()) {
862  info.errorMsg("Error in Pythia::init: "
863  "could not set up PDF for beam A");
864  return false;
865  }
866  pdfHardAPtr = pdfAPtr;
867  useNewPdfA = true;
868  }
869  if (pdfBPtr == 0) {
870  pdfBPtr = getPDFPtr(idB);
871  if (pdfBPtr == 0 || !pdfBPtr->isSetup()) {
872  info.errorMsg("Error in Pythia::init: "
873  "could not set up PDF for beam B");
874  return false;
875  }
876  pdfHardBPtr = pdfBPtr;
877  useNewPdfB = true;
878  }
879 
880  // Optionally set up separate PDF's for hard process.
881  if (settings.flag("PDF:useHard") && useNewPdfA && useNewPdfB) {
882  pdfHardAPtr = getPDFPtr(idA, 2);
883  if (!pdfHardAPtr->isSetup()) return false;
884  pdfHardBPtr = getPDFPtr(idB, 2);
885  if (!pdfHardBPtr->isSetup()) return false;
886  useNewPdfHard = true;
887  }
888 
889  // Optionally set up Pomeron PDF's for diffractive physics.
890  if ( doDiffraction) {
891  if (pdfPomAPtr == 0) {
892  pdfPomAPtr = getPDFPtr(990);
893  useNewPdfPomA = true;
894  }
895  if (pdfPomBPtr == 0) {
896  pdfPomBPtr = getPDFPtr(990);
897  useNewPdfPomB = true;
898  }
899  }
900 
901  // Done.
902  return true;
903 
904 }
905 
906 //--------------------------------------------------------------------------
907 
908 // Main routine to generate the next event, using internal machinery.
909 
910 bool Pythia::next() {
911 
912  // Regularly print how many events have been generated.
913  int nPrevious = info.getCounter(3);
914  if (nCount > 0 && nPrevious > 0 && nPrevious%nCount == 0)
915  cout << "\n Pythia::next(): " << nPrevious
916  << " events have been generated " << endl;
917 
918  // Set/reset info counters specific to each event.
919  info.addCounter(3);
920  for (int i = 10; i < 13; ++i) info.setCounter(i);
921 
922  // Simpler option when only HadronLevel to be generated.
923  if (!doProcessLevel) {
924 
925  // Reset info array (while event record contains data).
926  info.clear();
927 
928  // Set correct energy for system.
929  Vec4 pSum = 0.;
930  for (int i = 1; i < event.size(); ++i)
931  if (event[i].isFinal()) pSum += event[i].p();
932  event[0].p( pSum );
933  event[0].m( pSum.mCalc() );
934 
935  // Generate hadronization and decays.
936  bool status = forceHadronLevel();
937  if (status) info.addCounter(4);
938  if (status && nPrevious < nShowEvt) event.list(showSaV, showMaD);
939  return status;
940  }
941 
942  // Reset arrays.
943  info.clear();
944  process.clear();
945  event.clear();
946  partonSystems.clear();
947  beamA.clear();
948  beamB.clear();
949  beamPomA.clear();
950  beamPomB.clear();
951 
952  // Can only generate event if initialization worked.
953  if (!isInit) {
954  info.errorMsg("Abort from Pythia::next: "
955  "not properly initialized so cannot generate events");
956  return false;
957  }
958 
959  // Pick beam momentum spread and beam vertex.
960  if (doMomentumSpread || doVertexSpread) beamShapePtr->pick();
961 
962  // Recalculate kinematics when beam momentum spread.
963  if (doMomentumSpread) nextKinematics();
964 
965  // Outer loop over hard processes; only relevant for user-set vetoes.
966  for ( ; ; ) {
967  info.addCounter(10);
968  bool hasVetoed = false;
969 
970  // Provide the hard process that starts it off. Only one try.
971  info.clear();
972  process.clear();
973 
974  if ( !processLevel.next( process) ) {
975  if (doLHA && info.atEndOfFile()) info.errorMsg("Abort from "
976  "Pythia::next: reached end of Les Houches Events File");
977  else info.errorMsg("Abort from Pythia::next: "
978  "processLevel failed; giving up");
979  return false;
980  }
981 
982  info.addCounter(11);
983 
984  // Possiblility to perform CKKWL merging on this event
985  if(doMerging) {
986  hasVetoed = mergeProcess();
987  if (hasVetoed) {
988  if (abortIfVeto) return false;
989  continue;
990  }
991  }
992 
993  // Possibility for a user veto of the process-level event.
994  if (doVetoProcess) {
995  hasVetoed = userHooksPtr->doVetoProcessLevel( process);
996  if (hasVetoed) {
997  if (abortIfVeto) return false;
998  continue;
999  }
1000  }
1001 
1002  // Possibility to stop the generation at this stage.
1003  if (!doPartonLevel) {
1004  boostAndVertex( true, true);
1005  processLevel.accumulate();
1006  info.addCounter(4);
1007  if (doLHA && nPrevious < nShowLHA) lhaUpPtr->listEvent();
1008  if (nPrevious < nShowInfo) info.list();
1009  if (nPrevious < nShowProc) process.list(showSaV, showMaD);
1010  return true;
1011  }
1012 
1013  // Save spare copy of process record in case of problems.
1014  Event processSave = process;
1015  info.addCounter(12);
1016  for (int i = 14; i < 19; ++i) info.setCounter(i);
1017 
1018  // Allow up to ten tries for parton- and hadron-level processing.
1019  bool physical = true;
1020  for (int iTry = 0; iTry < NTRY; ++ iTry) {
1021  info.addCounter(14);
1022  physical = true;
1023 
1024  // Restore original process record if problems.
1025  if (iTry > 0) process = processSave;
1026 
1027  // Reset event record and (extracted partons from) beam remnants.
1028  event.clear();
1029  beamA.clear();
1030  beamB.clear();
1031  beamPomA.clear();
1032  beamPomB.clear();
1033  partonSystems.clear();
1034 
1035  // Parton-level evolution: ISR, FSR, MPI.
1036  if ( !partonLevel.next( process, event) ) {
1037  // Skip to next hard process for failure owing to deliberate veto.
1038  hasVetoed = partonLevel.hasVetoed();
1039  if (hasVetoed) {
1040  if (abortIfVeto) return false;
1041  break;
1042  }
1043  // Else make a new try for other failures.
1044  info.errorMsg("Error in Pythia::next: "
1045  "partonLevel failed; try again");
1046  physical = false;
1047  continue;
1048  }
1049  info.addCounter(15);
1050 
1051  // Possibility for a user veto of the parton-level event.
1052  if (doVetoPartons) {
1053  hasVetoed = userHooksPtr->doVetoPartonLevel( event);
1054  if (hasVetoed) {
1055  if (abortIfVeto) return false;
1056  break;
1057  }
1058  }
1059 
1060  // Boost to lab frame (before decays, for vertices).
1061  boostAndVertex( true, true);
1062 
1063  // Possibility to stop the generation at this stage.
1064  if (!doHadronLevel) {
1065  processLevel.accumulate();
1066  partonLevel.accumulate();
1067 
1068  // Optionally check final event for problems.
1069  if (checkEvent && !check()) {
1070  info.errorMsg("Abort from Pythia::next: "
1071  "check of event revealed problems");
1072  return false;
1073  }
1074  info.addCounter(4);
1075  if (doLHA && nPrevious < nShowLHA) lhaUpPtr->listEvent();
1076  if (nPrevious < nShowInfo) info.list();
1077  if (nPrevious < nShowProc) process.list(showSaV, showMaD);
1078  if (nPrevious < nShowEvt) event.list(showSaV, showMaD);
1079  return true;
1080  }
1081 
1082  // Hadron-level: hadronization, decays.
1083  info.addCounter(16);
1084  if ( !hadronLevel.next( event) ) {
1085  info.errorMsg("Error in Pythia::next: "
1086  "hadronLevel failed; try again");
1087  physical = false;
1088  continue;
1089  }
1090 
1091  // If R-hadrons have been formed, then (optionally) let them decay.
1092  if (decayRHadrons && rHadrons.exist() && !doRHadronDecays()) {
1093  info.errorMsg("Error in Pythia::next: "
1094  "decayRHadrons failed; try again");
1095  physical = false;
1096  continue;
1097  }
1098  info.addCounter(17);
1099 
1100  // Optionally check final event for problems.
1101  if (checkEvent && !check()) {
1102  info.errorMsg("Error in Pythia::next: "
1103  "check of event revealed problems");
1104  physical = false;
1105  continue;
1106  }
1107 
1108  // Stop parton- and hadron-level looping if you got this far.
1109  info.addCounter(18);
1110  break;
1111  }
1112 
1113  // If event vetoed then to make a new try.
1114  if (hasVetoed) {
1115  if (abortIfVeto) return false;
1116  continue;
1117  }
1118 
1119  // If event failed any other way (after ten tries) then give up.
1120  if (!physical) {
1121  info.errorMsg("Abort from Pythia::next: "
1122  "parton+hadronLevel failed; giving up");
1123  return false;
1124  }
1125 
1126  // Process- and parton-level statistics. Event scale.
1127  processLevel.accumulate();
1128  partonLevel.accumulate();
1129  event.scale( process.scale() );
1130 
1131  // End of outer loop over hard processes. Done with normal option.
1132  info.addCounter(13);
1133  break;
1134  }
1135 
1136  // List events.
1137  if (doLHA && nPrevious < nShowLHA) lhaUpPtr->listEvent();
1138  if (nPrevious < nShowInfo) info.list();
1139  if (nPrevious < nShowProc) process.list(showSaV,showMaD);
1140  if (nPrevious < nShowEvt) event.list(showSaV, showMaD);
1141 
1142  // Done.
1143  info.addCounter(4);
1144  return true;
1145 
1146 }
1147 
1148 //--------------------------------------------------------------------------
1149 
1150 // Generate only the hadronization/decay stage, using internal machinery.
1151 // The "event" instance should already contain a parton-level configuration.
1152 
1153 bool Pythia::forceHadronLevel(bool findJunctions) {
1154 
1155  // Can only generate event if initialization worked.
1156  if (!isInit) {
1157  info.errorMsg("Abort from Pythia::forceHadronLevel: "
1158  "not properly initialized so cannot generate events");
1159  return false;
1160  }
1161 
1162  // Check whether any junctions in system. (Normally done in ProcessLevel.)
1163  if (findJunctions) {
1164  event.clearJunctions();
1165  processLevel.findJunctions( event);
1166  }
1167 
1168  // Save spare copy of event in case of failure.
1169  Event spareEvent = event;
1170 
1171  // Allow up to ten tries for hadron-level processing.
1172  bool physical = true;
1173  for (int iTry = 0; iTry < NTRY; ++ iTry) {
1174  physical = true;
1175 
1176  // Hadron-level: hadronization, decays.
1177  if (hadronLevel.next( event)) break;
1178 
1179  // If failure then warn, restore original configuration and try again.
1180  info.errorMsg("Error in Pythia::forceHadronLevel: "
1181  "hadronLevel failed; try again");
1182  physical = false;
1183  event = spareEvent;
1184  }
1185 
1186  // Done for simpler option.
1187  if (!physical) {
1188  info.errorMsg("Abort from Pythia::forceHadronLevel: "
1189  "hadronLevel failed; giving up");
1190  return false;
1191  }
1192 
1193  // Optionally check final event for problems.
1194  if (checkEvent && !check()) {
1195  info.errorMsg("Abort from Pythia::forceHadronLevel: "
1196  "check of event revealed problems");
1197  return false;
1198  }
1199 
1200  // Done.
1201  return true;
1202 
1203 }
1204 
1205 //--------------------------------------------------------------------------
1206 
1207 // Recalculate kinematics for each event when beam momentum has a spread.
1208 
1209 void Pythia::nextKinematics() {
1210 
1211  // Read out momentum shift to give current beam momenta.
1212  pAnow = pAinit + beamShapePtr->deltaPA();
1213  pAnow.e( sqrt(pAnow.pAbs2() + mA * mA) );
1214  pBnow = pBinit + beamShapePtr->deltaPB();
1215  pBnow.e( sqrt(pBnow.pAbs2() + mB * mB) );
1216 
1217  // Construct CM frame kinematics.
1218  eCM = (pAnow + pBnow).mCalc();
1219  pzAcm = 0.5 * sqrtpos( (eCM + mA + mB) * (eCM - mA - mB)
1220  * (eCM - mA + mB) * (eCM + mA - mB) ) / eCM;
1221  pzBcm = -pzAcm;
1222  eA = sqrt(mA*mA + pzAcm*pzAcm);
1223  eB = sqrt(mB*mB + pzBcm*pzBcm);
1224 
1225  // Set relevant info for other classes to use.
1226  info.setBeamA( idA, pzAcm, eA, mA);
1227  info.setBeamB( idB, pzBcm, eB, mB);
1228  info.setECM( eCM);
1229  beamA.newPzE( pzAcm, eA);
1230  beamB.newPzE( pzBcm, eB);
1231 
1232  // Set boost/rotation matrices from/to CM frame.
1233  MfromCM.reset();
1234  MfromCM.fromCMframe( pAnow, pBnow);
1235  MtoCM = MfromCM;
1236  MtoCM.invert();
1237 
1238 }
1239 
1240 //--------------------------------------------------------------------------
1241 
1242 // Boost from CM frame to lab frame, or inverse. Set production vertex.
1243 
1244 void Pythia::boostAndVertex( bool toLab, bool setVertex) {
1245 
1246  // Boost process from CM frame to lab frame.
1247  if (toLab) {
1248  if (boostType == 2) process.bst(0., 0., betaZ, gammaZ);
1249  else if (boostType == 3) process.rotbst(MfromCM);
1250 
1251  // Boost nonempty event from CM frame to lab frame.
1252  if (event.size() > 0) {
1253  if (boostType == 2) event.bst(0., 0., betaZ, gammaZ);
1254  else if (boostType == 3) event.rotbst(MfromCM);
1255  }
1256 
1257  // Boost process from lab frame to CM frame.
1258  } else {
1259  if (boostType == 2) process.bst(0., 0., -betaZ, gammaZ);
1260  else if (boostType == 3) process.rotbst(MtoCM);
1261 
1262  // Boost nonempty event from lab frame to CM frame.
1263  if (event.size() > 0) {
1264  if (boostType == 2) event.bst(0., 0., -betaZ, gammaZ);
1265  else if (boostType == 3) event.rotbst(MtoCM);
1266  }
1267  }
1268 
1269  // Set production vertex; assumes particles are in lab frame and at origin.
1270  if (setVertex && doVertexSpread) {
1271  Vec4 vertex = beamShapePtr->vertex();
1272  for (int i = 0; i < process.size(); ++i) process[i].vProd( vertex);
1273  for (int i = 0; i < event.size(); ++i) event[i].vProd( vertex);
1274  }
1275 
1276 }
1277 
1278 //--------------------------------------------------------------------------
1279 
1280 // Perform R-hadron decays, either as part of normal evolution or forced.
1281 
1282 bool Pythia::doRHadronDecays( ) {
1283 
1284  // Check if R-hadrons exist to be processed.
1285  if ( !rHadrons.exist() ) return true;
1286 
1287  // Do the R-hadron decay itself.
1288  if ( !rHadrons.decay( event) ) return false;
1289 
1290  // Perform showers in resonance decay chains.
1291  if ( !partonLevel.resonanceShowers( process, event, false) ) return false;
1292 
1293  // Subsequent hadronization and decays.
1294  if ( !hadronLevel.next( event) ) return false;
1295 
1296  // Done.
1297  return true;
1298 
1299 }
1300 
1301 //--------------------------------------------------------------------------
1302 
1303 // Print statistics on event generation.
1304 
1305 void Pythia::stat() {
1306 
1307  // Read out settings for what to include.
1308  bool showPrL = settings.flag("Stat:showProcessLevel");
1309  bool showPaL = settings.flag("Stat:showPartonLevel");
1310  bool showErr = settings.flag("Stat:showErrors");
1311  bool reset = settings.flag("Stat:reset");
1312 
1313  // Statistics on cross section and number of events.
1314  if (doProcessLevel) {
1315  if (showPrL) processLevel.statistics(false);
1316  if (reset) processLevel.resetStatistics();
1317  }
1318 
1319  // Statistics from other classes, currently multiparton interactions.
1320  if (showPaL) partonLevel.statistics(false);
1321  if (reset) partonLevel.resetStatistics();
1322 
1323  // Summary of which and how many warnings/errors encountered.
1324  if (showErr) info.errorStatistics();
1325  if (reset) info.errorReset();
1326 
1327 }
1328 
1329 //--------------------------------------------------------------------------
1330 
1331 // Print statistics on event generation - deprecated version.
1332 
1333 void Pythia::statistics(bool all, bool reset) {
1334 
1335  // Statistics on cross section and number of events.
1336  if (doProcessLevel) processLevel.statistics(reset);
1337 
1338  // Statistics from other classes, e.g. multiparton interactions.
1339  if (all) partonLevel.statistics(reset);
1340 
1341  // Summary of which and how many warnings/errors encountered.
1342  info.errorStatistics();
1343  if (reset) info.errorReset();
1344 
1345 }
1346 
1347 //--------------------------------------------------------------------------
1348 
1349 // Write the Pythia banner, with symbol and version information.
1350 
1351 void Pythia::banner(ostream& os) {
1352 
1353  // Read in version number and last date of change.
1354  double versionNumber = settings.parm("Pythia:versionNumber");
1355  int versionDate = settings.mode("Pythia:versionDate");
1356  string month[12] = {"Jan", "Feb", "Mar", "Apr", "May", "Jun",
1357  "Jul", "Aug", "Sep", "Oct", "Nov", "Dec"};
1358 
1359  // Get date and time.
1360  time_t t = time(0);
1361  char dateNow[12];
1362  strftime(dateNow,12,"%d %b %Y",localtime(&t));
1363  char timeNow[9];
1364  strftime(timeNow,9,"%H:%M:%S",localtime(&t));
1365 
1366  os << "\n"
1367  << " *-------------------------------------------"
1368  << "-----------------------------------------* \n"
1369  << " | "
1370  << " | \n"
1371  << " | *----------------------------------------"
1372  << "--------------------------------------* | \n"
1373  << " | | "
1374  << " | | \n"
1375  << " | | "
1376  << " | | \n"
1377  << " | | PPP Y Y TTTTT H H III A "
1378  << " Welcome to the Lund Monte Carlo! | | \n"
1379  << " | | P P Y Y T H H I A A "
1380  << " This is PYTHIA version " << fixed << setprecision(3)
1381  << setw(5) << versionNumber << " | | \n"
1382  << " | | PPP Y T HHHHH I AAAAA"
1383  << " Last date of change: " << setw(2) << versionDate%100
1384  << " " << month[ (versionDate/100)%100 - 1 ]
1385  << " " << setw(4) << versionDate/10000 << " | | \n"
1386  << " | | P Y T H H I A A"
1387  << " | | \n"
1388  << " | | P Y T H H III A A"
1389  << " Now is " << dateNow << " at " << timeNow << " | | \n"
1390  << " | | "
1391  << " | | \n"
1392  << " | | Torbjorn Sjostrand; Department of As"
1393  << "tronomy and Theoretical Physics, | | \n"
1394  << " | | Lund University, Solvegatan 14A, S"
1395  << "E-223 62 Lund, Sweden; | | \n"
1396  << " | | phone: + 46 - 46 - 222 48 16; e-ma"
1397  << "il: torbjorn@thep.lu.se | | \n"
1398  << " | | Stefan Ask; Department of Physics, U"
1399  << "niversity of Cambridge, | | \n"
1400  << " | | Cavendish Laboratory, JJ Thomson A"
1401  << "ve., Cambridge CB3 0HE, UK; | | \n"
1402  << " | | phone: + 41 - 22 - 767 6707; e-mai"
1403  << "l: Stefan.Ask@cern.ch | | \n"
1404  << " | | Richard Corke; Department of Astrono"
1405  << "my and Theoretical Physics, | | \n"
1406  << " | | Lund University, Solvegatan 14A, S"
1407  << "E-223 62 Lund, Sweden; | | \n"
1408  << " | | e-mail: richard.corke@thep.lu.se "
1409  << " | | \n"
1410  << " | | Stephen Mrenna; Computing Division, "
1411  << "Simulations Group, | | \n"
1412  << " | | Fermi National Accelerator Laborat"
1413  << "ory, MS 234, Batavia, IL 60510, USA; | | \n"
1414  << " | | phone: + 1 - 630 - 840 - 2556; e-m"
1415  << "ail: mrenna@fnal.gov | | \n"
1416  << " | | Stefan Prestel; Department of Astron"
1417  << "omy and Theoretical Physics, | | \n"
1418  << " | | Lund University, Solvegatan 14A, S"
1419  << "E-223 62 Lund, Sweden; | | \n"
1420  << " | | phone: + 46 - 46 - 222 06 64; e-ma"
1421  << "il: stefan.prestel@thep.lu.se | | \n"
1422  << " | | Peter Skands; Theoretical Physics, C"
1423  << "ERN, CH-1211 Geneva 23, Switzerland; | | \n"
1424  << " | | phone: + 41 - 22 - 767 2447; e-mai"
1425  << "l: peter.skands@cern.ch | | \n"
1426  << " | | "
1427  << " | | \n"
1428  << " | | The main program reference is the 'Br"
1429  << "ief Introduction to PYTHIA 8.1', | | \n"
1430  << " | | T. Sjostrand, S. Mrenna and P. Skands"
1431  << ", Comput. Phys. Comm. 178 (2008) 85 | | \n"
1432  << " | | [arXiv:0710.3820] "
1433  << " | | \n"
1434  << " | | "
1435  << " | | \n"
1436  << " | | The main physics reference is the 'PY"
1437  << "THIA 6.4 Physics and Manual', | | \n"
1438  << " | | T. Sjostrand, S. Mrenna and P. Skands"
1439  << ", JHEP05 (2006) 026 [hep-ph/0603175]. | | \n"
1440  << " | | "
1441  << " | | \n"
1442  << " | | An archive of program versions and do"
1443  << "cumentation is found on the web: | | \n"
1444  << " | | http://www.thep.lu.se/~torbjorn/Pythi"
1445  << "a.html | | \n"
1446  << " | | "
1447  << " | | \n"
1448  << " | | This program is released under the GN"
1449  << "U General Public Licence version 2. | | \n"
1450  << " | | Please respect the MCnet Guidelines f"
1451  << "or Event Generator Authors and Users. | | \n"
1452  << " | | "
1453  << " | | \n"
1454  << " | | Disclaimer: this program comes withou"
1455  << "t any guarantees. | | \n"
1456  << " | | Beware of errors and use common sense"
1457  << " when interpreting results. | | \n"
1458  << " | | "
1459  << " | | \n"
1460  << " | | Copyright (C) 2012 Torbjorn Sjostrand"
1461  << " | | \n"
1462  << " | | "
1463  << " | | \n"
1464  << " | | "
1465  << " | | \n"
1466  << " | *----------------------------------------"
1467  << "--------------------------------------* | \n"
1468  << " | "
1469  << " | \n"
1470  << " *-------------------------------------------"
1471  << "-----------------------------------------* \n" << endl;
1472 
1473 }
1474 
1475 //--------------------------------------------------------------------------
1476 
1477 // Check for lines in file that mark the beginning of new subrun.
1478 
1479 int Pythia::readSubrun(string line, bool warn, ostream& os) {
1480 
1481  // If empty line then done.
1482  int subrunLine = SUBRUNDEFAULT;
1483  if (line.find_first_not_of(" \n\t\v\b\r\f\a") == string::npos)
1484  return subrunLine;
1485 
1486  // If first character is not a letter, then done.
1487  string lineNow = line;
1488  int firstChar = lineNow.find_first_not_of(" \n\t\v\b\r\f\a");
1489  if (!isalpha(lineNow[firstChar])) return subrunLine;
1490 
1491  // Replace an equal sign by a blank to make parsing simpler.
1492  while (lineNow.find("=") != string::npos) {
1493  int firstEqual = lineNow.find_first_of("=");
1494  lineNow.replace(firstEqual, 1, " ");
1495  }
1496 
1497  // Get first word of a line.
1498  istringstream splitLine(lineNow);
1499  string name;
1500  splitLine >> name;
1501 
1502  // Replace two colons by one (:: -> :) to allow for such mistakes.
1503  while (name.find("::") != string::npos) {
1504  int firstColonColon = name.find_first_of("::");
1505  name.replace(firstColonColon, 2, ":");
1506  }
1507 
1508  // Convert to lowercase.
1509  for (int i = 0; i < int(name.length()); ++i)
1510  name[i] = std::tolower(name[i]);
1511 
1512  // If no match then done.
1513  if (name != "main:subrun") return subrunLine;
1514 
1515  // Else find new subrun number and return it.
1516  splitLine >> subrunLine;
1517  if (!splitLine) {
1518  if (warn) os << "\n PYTHIA Warning: Main:subrun number not"
1519  << " recognized; skip:\n " << line << endl;
1520  subrunLine = SUBRUNDEFAULT;
1521  }
1522  return subrunLine;
1523 
1524 }
1525 
1526 //--------------------------------------------------------------------------
1527 
1528 // Check that the final event makes sense: no unknown id codes;
1529 // charge and energy-momentum conserved.
1530 
1531 bool Pythia::check(ostream& os) {
1532 
1533  // Reset.
1534  bool physical = true;
1535  bool listVertices = false;
1536  bool listSystems = false;
1537  bool listBeams = false;
1538  iErrId.resize(0);
1539  iErrCol.resize(0);
1540  iErrNan.resize(0);
1541  iErrNanVtx.resize(0);
1542  Vec4 pSum;
1543  double chargeSum = 0.;
1544 
1545  // Incoming beams counted with negative momentum and charge.
1546  if (doProcessLevel) {
1547  pSum = - (event[1].p() + event[2].p());
1548  chargeSum = - (event[1].charge() + event[2].charge());
1549 
1550  // If no ProcessLevel then sum final state of process record.
1551  } else if (process.size() > 0) {
1552  pSum = - process[0].p();
1553  for (int i = 0; i < process.size(); ++i)
1554  if (process[i].isFinal()) chargeSum -= process[i].charge();
1555 
1556  // If process not filled, then use outgoing primary in event.
1557  } else {
1558  pSum = - event[0].p();
1559  for (int i = 1; i < event.size(); ++i)
1560  if (event[i].statusAbs() < 10 || event[i].statusAbs() == 23)
1561  chargeSum -= event[i].charge();
1562  }
1563  double eLab = abs(pSum.e());
1564 
1565  // Loop over particles in the event.
1566  for (int i = 0; i < event.size(); ++i) {
1567 
1568  // Look for any unrecognized particle codes.
1569  int id = event[i].id();
1570  if (id == 0 || !particleData.isParticle(id)) {
1571  ostringstream errCode;
1572  errCode << ", id = " << id;
1573  info.errorMsg("Error in Pythia::check: "
1574  "unknown particle code", errCode.str());
1575  physical = false;
1576  iErrId.push_back(i);
1577 
1578  // Check that colour assignments are the expected ones.
1579  } else {
1580  int colType = event[i].colType();
1581  int col = event[i].col();
1582  int acol = event[i].acol();
1583  if ( (colType == 0 && (col > 0 || acol > 0))
1584  || (colType == 1 && (col <= 0 || acol > 0))
1585  || (colType == -1 && (col > 0 || acol <= 0))
1586  || (colType == 2 && (col <= 0 || acol <= 0)) ) {
1587  ostringstream errCode;
1588  errCode << ", id = " << id << " cols = " << col << " " << acol;
1589  info.errorMsg("Error in Pythia::check: "
1590  "incorrect colours", errCode.str());
1591  physical = false;
1592  iErrCol.push_back(i);
1593  }
1594  }
1595 
1596  // Look for particles with not-a-number energy/momentum/mass.
1597  if (abs(event[i].px()) >= 0. && abs(event[i].py()) >= 0.
1598  && abs(event[i].pz()) >= 0. && abs(event[i].e()) >= 0.
1599  && abs(event[i].m()) >= 0.) ;
1600  else {
1601  info.errorMsg("Error in Pythia::check: "
1602  "not-a-number energy/momentum/mass");
1603  physical = false;
1604  iErrNan.push_back(i);
1605  }
1606 
1607  // Look for particles with not-a-number vertex/lifetime.
1608  if (abs(event[i].xProd()) >= 0. && abs(event[i].yProd()) >= 0.
1609  && abs(event[i].zProd()) >= 0. && abs(event[i].tProd()) >= 0.
1610  && abs(event[i].tau()) >= 0.) ;
1611  else {
1612  info.errorMsg("Error in Pythia::check: "
1613  "not-a-number vertex/lifetime");
1614  physical = false;
1615  listVertices = true;
1616  iErrNanVtx.push_back(i);
1617  }
1618 
1619  // Add final-state four-momentum and charge.
1620  if (event[i].isFinal()) {
1621  pSum += event[i].p();
1622  chargeSum += event[i].charge();
1623  }
1624 
1625  // End of particle loop.
1626  }
1627 
1628  // Check energy-momentum/charge conservation.
1629  double epDev = abs(pSum.e()) + abs(pSum.px()) + abs(pSum.py())
1630  + abs(pSum.pz());
1631  if (epDev > epTolErr * eLab) {
1632  info.errorMsg("Error in Pythia::check: energy-momentum not conserved");
1633  physical = false;
1634  } else if (epDev > epTolWarn * eLab) {
1635  info.errorMsg("Warning in Pythia::check: "
1636  "energy-momentum not quite conserved");
1637  }
1638  if (abs(chargeSum) > 0.1) {
1639  info.errorMsg("Error in Pythia::check: charge not conserved");
1640  physical = false;
1641  }
1642 
1643  // Check that beams and event records agree on incoming partons.
1644  // Only meaningful for resolved beams.
1645  if (info.isResolved())
1646  for (int iSys = 0; iSys < beamA.sizeInit(); ++iSys) {
1647  int eventANw = partonSystems.getInA(iSys);
1648  int eventBNw = partonSystems.getInB(iSys);
1649  int beamANw = beamA[iSys].iPos();
1650  int beamBNw = beamB[iSys].iPos();
1651  if (eventANw != beamANw || eventBNw != beamBNw) {
1652  info.errorMsg("Error in Pythia::check: "
1653  "event and beams records disagree");
1654  physical = false;
1655  listSystems = true;
1656  listBeams = true;
1657  }
1658  }
1659 
1660  // Done for sensible events.
1661  if (physical) return true;
1662 
1663  // Print (the first few) flawed events.
1664  if (nErrEvent < nErrList) {
1665  os << " PYTHIA erroneous event info: \n";
1666  if (iErrId.size() > 0) {
1667  os << " unknown particle codes in lines ";
1668  for (int i = 0; i < int(iErrId.size()); ++i)
1669  os << iErrId[i] << " ";
1670  os << "\n";
1671  }
1672  if (iErrCol.size() > 0) {
1673  os << " incorrect colour assignments in lines ";
1674  for (int i = 0; i < int(iErrCol.size()); ++i)
1675  os << iErrCol[i] << " ";
1676  os << "\n";
1677  }
1678  if (iErrNan.size() > 0) {
1679  os << " not-a-number energy/momentum/mass in lines ";
1680  for (int i = 0; i < int(iErrNan.size()); ++i)
1681  os << iErrNan[i] << " ";
1682  os << "\n";
1683  }
1684  if (iErrNanVtx.size() > 0) {
1685  os << " not-a-number vertex/lifetime in lines ";
1686  for (int i = 0; i < int(iErrNanVtx.size()); ++i)
1687  os << iErrNanVtx[i] << " ";
1688  os << "\n";
1689  }
1690  if (epDev > epTolErr * eLab) os << scientific << setprecision(3)
1691  << " total energy-momentum non-conservation = " << epDev << "\n";
1692  if (abs(chargeSum) > 0.1) os << fixed << setprecision(2)
1693  << " total charge non-conservation = " << chargeSum << "\n";
1694  info.list();
1695  event.list(listVertices);
1696  if (listSystems) partonSystems.list();
1697  if (listBeams) beamA.list();
1698  if (listBeams) beamB.list();
1699  }
1700 
1701  // Update error counter. Done also for flawed event.
1702  ++nErrEvent;
1703  return false;
1704 
1705 }
1706 
1707 //--------------------------------------------------------------------------
1708 
1709 // Routine to set up a PDF pointer.
1710 
1711 PDF* Pythia::getPDFPtr(int idIn, int sequence) {
1712 
1713  // Temporary pointer to be returned.
1714  PDF* tempPDFPtr = 0;
1715 
1716  // One option is to treat a Pomeron like a pi0.
1717  if (idIn == 990 && settings.mode("PDF:PomSet") == 2) idIn = 111;
1718 
1719  // Proton beam, normal choice.
1720  if (abs(idIn) == 2212 && sequence == 1) {
1721  int pSet = settings.mode("PDF:pSet");
1722  bool useLHAPDF = settings.flag("PDF:useLHAPDF");
1723 
1724  // Use internal sets.
1725  if (!useLHAPDF) {
1726  if (pSet == 1) tempPDFPtr = new GRV94L(idIn);
1727  else if (pSet == 2) tempPDFPtr = new CTEQ5L(idIn);
1728  else if (pSet <= 6)
1729  tempPDFPtr = new MSTWpdf(idIn, pSet - 2, xmlPath, &info);
1730  else tempPDFPtr = new CTEQ6pdf(idIn, pSet - 6, xmlPath, &info);
1731  }
1732 
1733  // Use sets from LHAPDF.
1734  else {
1735  string LHAPDFset = settings.word("PDF:LHAPDFset");
1736  int LHAPDFmember = settings.mode("PDF:LHAPDFmember");
1737  tempPDFPtr = new LHAPDF(idIn, LHAPDFset, LHAPDFmember, 1, &info);
1738 
1739  // Optionally allow extrapolation beyond x and Q2 limits.
1740  tempPDFPtr->setExtrapolate( settings.flag("PDF:extrapolateLHAPDF") );
1741  }
1742  }
1743 
1744  // Proton beam, special choice for the hard process.
1745  else if (abs(idIn) == 2212) {
1746  int pSet = settings.mode("PDF:pHardSet");
1747  bool useLHAPDF = settings.flag("PDF:useHardLHAPDF");
1748 
1749  // Use internal sets.
1750  if (!useLHAPDF) {
1751  if (pSet == 1) tempPDFPtr = new GRV94L(idIn);
1752  else if (pSet == 2) tempPDFPtr = new CTEQ5L(idIn);
1753  else if (pSet <= 6)
1754  tempPDFPtr = new MSTWpdf(idIn, pSet - 2, xmlPath, &info);
1755  else tempPDFPtr = new CTEQ6pdf(idIn, pSet - 6, xmlPath, &info);
1756  }
1757 
1758  // Use sets from LHAPDF.
1759  else {
1760  string LHAPDFset = settings.word("PDF:hardLHAPDFset");
1761  int LHAPDFmember = settings.mode("PDF:hardLHAPDFmember");
1762  tempPDFPtr = new LHAPDF(idIn, LHAPDFset, LHAPDFmember, 2, &info);
1763 
1764  // Optionally allow extrapolation beyond x and Q2 limits.
1765  tempPDFPtr->setExtrapolate( settings.flag("PDF:extrapolateLHAPDF") );
1766  }
1767  }
1768 
1769  // Pion beam (or, in one option, Pomeron beam).
1770  else if (abs(idIn) == 211 || idIn == 111) {
1771  bool useLHAPDF = settings.flag("PDF:piUseLHAPDF");
1772 
1773  // Use internal sets.
1774  if (!useLHAPDF) {
1775  tempPDFPtr = new GRVpiL(idIn);
1776  }
1777 
1778  // Use sets from LHAPDF.
1779  else {
1780  string LHAPDFset = settings.word("PDF:piLHAPDFset");
1781  int LHAPDFmember = settings.mode("PDF:piLHAPDFmember");
1782  tempPDFPtr = new LHAPDF(idIn, LHAPDFset, LHAPDFmember, 1, &info);
1783 
1784  // Optionally allow extrapolation beyond x and Q2 limits.
1785  tempPDFPtr->setExtrapolate( settings.flag("PDF:extrapolateLHAPDF") );
1786  }
1787  }
1788 
1789  // Pomeron beam, if not treated like a pi0 beam.
1790  else if (idIn == 990) {
1791  int pomSet = settings.mode("PDF:PomSet");
1792  double rescale = settings.parm("PDF:PomRescale");
1793 
1794  // A generic Q2-independent parametrization.
1795  if (pomSet == 1) {
1796  double gluonA = settings.parm("PDF:PomGluonA");
1797  double gluonB = settings.parm("PDF:PomGluonB");
1798  double quarkA = settings.parm("PDF:PomQuarkA");
1799  double quarkB = settings.parm("PDF:PomQuarkB");
1800  double quarkFrac = settings.parm("PDF:PomQuarkFrac");
1801  double strangeSupp = settings.parm("PDF:PomStrangeSupp");
1802  tempPDFPtr = new PomFix( 990, gluonA, gluonB, quarkA, quarkB,
1803  quarkFrac, strangeSupp);
1804  }
1805 
1806  // The H1 Q2-dependent parametrizations. Initialization requires files.
1807  else if (pomSet == 3 || pomSet == 4)
1808  tempPDFPtr = new PomH1FitAB( 990, pomSet - 2, rescale, xmlPath, &info);
1809  else if (pomSet == 5)
1810  tempPDFPtr = new PomH1Jets( 990, rescale, xmlPath, &info);
1811  else if (pomSet == 6)
1812  tempPDFPtr = new PomH1FitAB( 990, 3, rescale, xmlPath, &info);
1813  }
1814 
1815  // Lepton beam: neutrino, resolved charged lepton or unresolved ditto.
1816  else if (abs(idIn) > 10 && abs(idIn) < 17) {
1817  if (abs(idIn)%2 == 0) tempPDFPtr = new NeutrinoPoint(idIn);
1818  else if (settings.flag("PDF:lepton")) tempPDFPtr = new Lepton(idIn);
1819  else tempPDFPtr = new LeptonPoint(idIn);
1820  }
1821 
1822  // Failure for unrecognized particle.
1823  else tempPDFPtr = 0;
1824 
1825  // Done.
1826  return tempPDFPtr;
1827 }
1828 
1829 //--------------------------------------------------------------------------
1830 
1831 // Initialize SUSY Les Houches Accord data.
1832 
1833 bool Pythia::initSLHA() {
1834 
1835  // Initial and settings values.
1836  int ifailLHE = 1;
1837  int ifailSpc = 1;
1838  int readFrom = settings.mode("SLHA:readFrom");
1839  string lhefFile = settings.word("Beams:LHEF");
1840  string slhaFile = settings.word("SLHA:file");
1841  int verboseSLHA = settings.mode("SLHA:verbose");
1842  bool slhaUseDec = settings.flag("SLHA:useDecayTable");
1843 
1844  // No SUSY by default
1845  couplingsPtr->isSUSY = false;
1846 
1847  // Option with no SLHA read-in at all.
1848  if (readFrom == 0) return true;
1849 
1850  // First check LHEF header (if reading from LHEF)
1851  if (readFrom == 1 && lhefFile != "void") {
1852  ifailLHE = slha.readFile(lhefFile, verboseSLHA, slhaUseDec );
1853  }
1854 
1855  // If LHEF read successful, everything needed should already be ready
1856  if (ifailLHE == 0) {
1857  ifailSpc = 0;
1858  couplingsPtr->isSUSY = true;
1859  // If no LHEF file or no SLHA info in header, read from SLHA:file
1860  } else {
1861  lhefFile = "void";
1862  if ( settings.word("SLHA:file") == "none"
1863  || settings.word("SLHA:file") == "void"
1864  || settings.word("SLHA:file") == ""
1865  || settings.word("SLHA:file") == " ") return true;
1866  ifailSpc = slha.readFile(slhaFile,verboseSLHA, slhaUseDec);
1867  }
1868 
1869  // In case of problems, print error and fail init.
1870  if (ifailSpc != 0) {
1871  info.errorMsg("Error in Pythia::initSLHA: "
1872  "problem reading SLHA file", slhaFile);
1873  return false;
1874  } else {
1875  couplingsPtr->isSUSY = true;
1876  }
1877 
1878  // Check spectrum for consistency. Switch off SUSY if necessary.
1879  ifailSpc = slha.checkSpectrum();
1880 
1881  // ifail >= 1 : no MODSEL found -> don't switch on SUSY
1882  if (ifailSpc == 1) {
1883  // no SUSY, but MASS ok
1884  couplingsPtr->isSUSY = false;
1885  info.errorMsg("Warning in Pythia::initSLHA: "
1886  "No MODSEL found, keeping internal SUSY switched off");
1887  } else if (ifailSpc >= 2) {
1888  // no SUSY, but problems
1889  info.errorMsg("Warning in Pythia::initSLHA: "
1890  "Problem with SLHA MASS or QNUMBERS.");
1891  couplingsPtr->isSUSY = false;
1892  }
1893  // ifail = 0 : MODSEL found, spectrum OK
1894  else if (ifailSpc == 0) {
1895  // Print spectrum. Done.
1896  slha.printSpectrum(0);
1897  }
1898  else if (ifailSpc < 0) {
1899  info.errorMsg("Warning in Pythia::initSLHA: "
1900  "Problem with SLHA spectrum.",
1901  "\n Only using masses and switching off SUSY.");
1902  settings.flag("SUSY:all", false);
1903  couplingsPtr->isSUSY = false;
1904  slha.printSpectrum(ifailSpc);
1905  }
1906 
1907  // Import qnumbers
1908  if ( (ifailSpc == 1 || ifailSpc == 0) && slha.qnumbers.size() > 0) {
1909  for (int iQnum=0; iQnum < int(slha.qnumbers.size()); iQnum++) {
1910  // Always use positive id codes
1911  int id = abs(slha.qnumbers[iQnum](0));
1912  ostringstream idCode;
1913  idCode << id;
1914  if (particleData.isParticle(id)) {
1915  info.errorMsg("Warning in Pythia::initSLHA: "
1916  "ignoring QNUMBERS", "for id = "+idCode.str()
1917  +" (already exists)", true);
1918  } else {
1919  int qEM3 = slha.qnumbers[iQnum](1);
1920  int nSpins = slha.qnumbers[iQnum](2);
1921  int colRep = slha.qnumbers[iQnum](3);
1922  int hasAnti = slha.qnumbers[iQnum](4);
1923  // Translate colRep to PYTHIA colType
1924  int colType = 0;
1925  if (colRep == 3) colType = 1;
1926  else if (colRep == -3) colType = -1;
1927  else if (colRep == 8) colType = 2;
1928  else if (colRep == 6) colType = 3;
1929  else if (colRep == -6) colType = -3;
1930  // Default name: PDG code
1931  string name, antiName;
1932  ostringstream idStream;
1933  idStream<<id;
1934  name = idStream.str();
1935  antiName = "-"+name;
1936  if (iQnum < int(slha.qnumbersName.size())) {
1937  name = slha.qnumbersName[iQnum];
1938  if (iQnum < int(slha.qnumbersAntiName.size()))
1939  antiName = slha.qnumbersAntiName[iQnum];
1940  if (antiName == "") {
1941  if (name.find("+") != string::npos) {
1942  antiName = name;
1943  antiName.replace(antiName.find("+"),1,"-");
1944  } else if (name.find("-") != string::npos) {
1945  antiName = name;
1946  antiName.replace(antiName.find("+"),1,"-");
1947  } else {
1948  antiName = name+"bar";
1949  }
1950  }
1951  }
1952  if ( hasAnti == 0) {
1953  antiName = "";
1954  particleData.addParticle(id, name, nSpins, qEM3, colType);
1955  } else {
1956  particleData.addParticle(id, name, antiName, nSpins, qEM3, colType);
1957  }
1958  }
1959  }
1960  }
1961 
1962  // Import mass spectrum.
1963  bool keepSM = settings.flag("SLHA:keepSM");
1964  double minMassSM = settings.parm("SLHA:minMassSM");
1965  double massMargin = settings.parm("SLHA:minDecayDeltaM");
1966  if (ifailSpc == 1 || ifailSpc == 0) {
1967 
1968  // Loop through to update particle data.
1969  int id = slha.mass.first();
1970  for (int i = 1; i <= slha.mass.size() ; i++) {
1971  double mass = abs(slha.mass(id));
1972 
1973  // Ignore masses for known SM particles or particles with
1974  // default masses < minMassSM; overwrite masses for rest.
1975  if (keepSM && (id < 25 || (id > 80 && id < 1000000))) ;
1976  else if (id < 1000000 && particleData.m0(id) < minMassSM) {
1977  ostringstream idCode;
1978  idCode << id;
1979  info.errorMsg("Warning in Pythia::initSLHA: "
1980  "ignoring MASS entry", "for id = "+idCode.str()
1981  +" (m0 < SLHA:minMassSM)", true);
1982  }
1983  else particleData.m0(id,mass);
1984  id = slha.mass.next();
1985  };
1986 
1987  }
1988 
1989  // Update decay data.
1990  for (int iTable=0; iTable < int(slha.decays.size()); iTable++) {
1991 
1992  // Pointer to this SLHA table
1993  LHdecayTable* slhaTable=&(slha.decays[iTable]);
1994 
1995  // Extract ID and create pointer to corresponding particle data object
1996  int idRes = slhaTable->getId();
1997  ParticleDataEntry* particlePtr
1998  = particleData.particleDataEntryPtr(idRes);
1999 
2000  // Ignore decay channels for known SM particles or particles with
2001  // default masses < minMassSM; overwrite masses for rest.
2002  if (keepSM && (idRes < 25 || (idRes > 80 && idRes < 1000000))) continue;
2003  else if (idRes < 1000000 && particleData.m0(idRes) < minMassSM) {
2004  ostringstream idCode;
2005  idCode << idRes;
2006  info.errorMsg("Warning in Pythia::initSLHA: "
2007  "ignoring DECAY table", "for id = " + idCode.str()
2008  + " (m0 < SLHA:minMassSM)", true);
2009  continue;
2010  }
2011 
2012  // Extract and store total width (absolute value, neg -> switch off)
2013  double widRes = abs(slhaTable->getWidth());
2014  particlePtr->setMWidth(widRes);
2015 
2016  // Reset decay table of the particle. Allow decays, treat as resonance.
2017  if (slhaTable->size() > 0) {
2018  particlePtr->clearChannels();
2019  particleData.mayDecay(idRes,true);
2020  particleData.isResonance(idRes,true);
2021  }
2022 
2023  // Reset to stable if width <= 0.0
2024  if (slhaTable->getWidth() <= 0.0) particleData.mayDecay(idRes,false);
2025 
2026  // Set initial minimum mass.
2027  double brWTsum = 0.;
2028  double massWTsum = 0.;
2029 
2030  // Loop over SLHA channels, import into Pythia, treating channels
2031  // with negative branching fractions as having the equivalent positive
2032  // branching fraction, but being switched off for this run
2033  for (int iChannel=0 ; iChannel<slhaTable->size(); iChannel++) {
2034  LHdecayChannel slhaChannel = slhaTable->getChannel(iChannel);
2035  double brat = slhaChannel.getBrat();
2036  vector<int> idDa = slhaChannel.getIdDa();
2037  if (idDa.size() >= 9) {
2038  info.errorMsg("Error in Pythia::initSLHA: "
2039  "max number of decay products is 8.");
2040  } else if (idDa.size() <= 1) {
2041  info.errorMsg("Error in Pythia::initSLHA: "
2042  "min number of decay products is 2.");
2043  }
2044  else {
2045  int onMode = 1;
2046  if (brat < 0.0) onMode = 0;
2047 
2048  // Check phase space, including margin
2049  double massSum = massMargin;
2050  for (int jDa=0; jDa<int(idDa.size()); ++jDa)
2051  massSum += particleData.m0( idDa[jDa] );
2052  if (onMode == 1 && brat > 0.0 && massSum > particleData.m0(idRes) ) {
2053  // String containing decay name
2054  ostringstream errCode;
2055  errCode << idRes <<" ->";
2056  for (int jDa=0; jDa<int(idDa.size()); ++jDa) errCode<<" "<<idDa[jDa];
2057  info.errorMsg("Warning in Pythia::initSLHA: switching off decay",
2058  errCode.str() + " (mRes - mDa < minDecayDeltaM)"
2059  "\n (Note: cross sections will be scaled by remaining"
2060  " open branching fractions!)" , true);
2061  onMode=0;
2062  }
2063 
2064  // Branching-ratio-weighted average mass in decay.
2065  brWTsum += abs(brat);
2066  massWTsum += abs(brat) * massSum;
2067 
2068  // Add channel
2069  int id0 = idDa[0];
2070  int id1 = idDa[1];
2071  int id2 = (idDa.size() >= 3) ? idDa[2] : 0;
2072  int id3 = (idDa.size() >= 4) ? idDa[3] : 0;
2073  int id4 = (idDa.size() >= 5) ? idDa[4] : 0;
2074  int id5 = (idDa.size() >= 6) ? idDa[5] : 0;
2075  int id6 = (idDa.size() >= 7) ? idDa[6] : 0;
2076  int id7 = (idDa.size() >= 8) ? idDa[7] : 0;
2077  particlePtr->addChannel(onMode,abs(brat),101,
2078  id0,id1,id2,id3,id4,id5,id6,id7);
2079 
2080  }
2081  }
2082 
2083  // Set minimal mass, but always below nominal one.
2084  if (slhaTable->size() > 0) {
2085  double massAvg = massWTsum / brWTsum;
2086  double massMin = min( massAvg, particlePtr->m0()) - massMargin;
2087  particlePtr->setMMin(massMin);
2088  }
2089  }
2090 
2091  return true;
2092 
2093 }
2094 
2095 //--------------------------------------------------------------------------
2096 
2097 // Function to perform CKKWL merging on the input event.
2098 
2099 bool Pythia::mergeProcess() {
2100 
2101  // Reset weight of the event.
2102  mergingHooksPtr->setWeight(1.);
2103  info.setMergingWeight(1.);
2104  double wgt = 1.0;
2105 
2106  // Store candidates for the splitting V -> qqbar'.
2107  mergingHooksPtr->storeHardProcessCandidates( process);
2108 
2109  // Calculate number of clustering steps.
2110  int nSteps = mergingHooksPtr->getNumberOfClusteringSteps( process);
2111 
2112  // Save number of hard final state quarks / leptons.
2113  int nQuarks = mergingHooksPtr->nHardOutPartons();
2114  int nLeptons = mergingHooksPtr->nHardOutLeptons();
2115  bool isHadronic = (mergingHooksPtr->nHardInLeptons() == 0);
2116 
2117  // Set QCD 2->2 starting scale different from arbirtrary scale in LHEF!
2118  // --> Set to pT of partons.
2119  double pT = 0.;
2120  for( int i=0; i < process.size(); ++i)
2121  if(process[i].isFinal() && process[i].colType() != 0) {
2122  pT = process[i].pT();
2123  break;
2124  }
2125 
2126  // Count number of top quarks to distinguish process with stable top pair
2127  // from pure QCD di-jets.
2128  int nTops = 0;
2129  for( int i=0; i < process.size(); ++i)
2130  if(process[i].isFinal() && process[i].idAbs() == 6) ++nTops;
2131 
2132  // Declare pdfWeight and alpha_s weight.
2133  double RN = rndm.flat();
2134  process.scale(0.0);
2135 
2136  // Generate all histories.
2137  History FullHistory( nSteps, 0.0, process, Clustering(), mergingHooksPtr,
2138  beamA, beamB, &particleData, &info, true, true, 1.0, 0);
2139 
2140  // Perform reweighting with Sudakov factors, save alpha_s ratios and
2141  // PDF ratio weights.
2142  wgt = FullHistory.weightTREE( &trialPartonLevel,
2143  mergingHooksPtr->AlphaS_FSR(), mergingHooksPtr->AlphaS_ISR(), RN);
2144 
2145  // Event with production scales set for further (trial) showering
2146  // and starting conditions for the shower.
2147  FullHistory.getStartingConditions( RN, process );
2148 
2149  // Allow to dampen histories in which the lowest multiplicity reclustered
2150  // state does not pass the lowest multiplicity cut of the matrix element.
2151  double dampWeight = mergingHooksPtr->dampenIfFailCuts(
2152  FullHistory.lowestMultProc(RN) );
2153 
2154  // Save the weight of the event for histogramming. Only change the
2155  // event weight after trial shower on the matrix element
2156  // multiplicity event (= in doVetoStep).
2157  wgt *= dampWeight;
2158 
2159  // For pure QCD dijet events (only!), set the process scale to the
2160  // transverse momentum of the outgoing partons.
2161  if ( isHadronic && nQuarks == 2 && nLeptons == 0 && nTops == 0
2162  && nSteps == 0 ) process.scale(pT);
2163 
2164  // Save the weight of the event for histogramming.
2165  mergingHooksPtr->setWeight(wgt);
2166  info.setMergingWeight(wgt);
2167 
2168  // If the weight of the event is zero, do not continue evolution.
2169  if(wgt == 0.) return true;
2170 
2171  // Done
2172  return false;
2173 
2174 }
2175 
2176 //==========================================================================
2177 
2178 } // end namespace Pythia8
Definition: AgUStep.h:26