StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
ProcessLevel.cc
1 // ProcessLevel.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 ProcessLevel class.
7 
8 #include "Pythia8/ProcessLevel.h"
9 
10 namespace Pythia8 {
11 
12 //==========================================================================
13 
14 // The ProcessLevel class.
15 
16 //--------------------------------------------------------------------------
17 
18 // Constants: could be changed here if desired, but normally should not.
19 // These are of technical nature, as described for each.
20 
21 // Allow a few failures in final construction of events.
22 const int ProcessLevel::MAXLOOP = 5;
23 
24 //--------------------------------------------------------------------------
25 
26 // Destructor.
27 
28 ProcessLevel::~ProcessLevel() {
29 
30  // Run through list of first hard processes and delete them.
31  for (int i = 0; i < int(containerPtrs.size()); ++i)
32  delete containerPtrs[i];
33 
34  // Run through list of second hard processes and delete them.
35  for (int i = 0; i < int(container2Ptrs.size()); ++i)
36  delete container2Ptrs[i];
37 
38 }
39 
40 //--------------------------------------------------------------------------
41 
42 // Main routine to initialize generation process.
43 
44 bool ProcessLevel::init( Info* infoPtrIn, Settings& settings,
45  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn,
46  BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
47  BeamParticle* beamGamAPtrIn, BeamParticle* beamGamBPtrIn,
48  BeamParticle* beamVMDAPtrIn, BeamParticle* beamVMDBPtrIn,
49  Couplings* couplingsPtrIn, SigmaTotal* sigmaTotPtrIn, bool doLHA,
50  SLHAinterface* slhaInterfacePtrIn, UserHooks* userHooksPtrIn,
51  vector<SigmaProcess*>& sigmaPtrs, vector<PhaseSpace*>& phaseSpacePtrs) {
52 
53  // Store input pointers for future use.
54  infoPtr = infoPtrIn;
55  particleDataPtr = particleDataPtrIn;
56  rndmPtr = rndmPtrIn;
57  beamAPtr = beamAPtrIn;
58  beamBPtr = beamBPtrIn;
59  couplingsPtr = couplingsPtrIn;
60  sigmaTotPtr = sigmaTotPtrIn;
61  userHooksPtr = userHooksPtrIn;
62  slhaInterfacePtr = slhaInterfacePtrIn;
63 
64  // Check whether photon inside lepton and save the mode.
65  beamHasGamma = settings.flag("PDF:lepton2gamma");
66  gammaMode = settings.mode("Photon:ProcessType");
67  bool beamA2gamma = beamAPtr->isLepton() && beamHasGamma;
68  bool beamB2gamma = beamBPtr->isLepton() && beamHasGamma;
69 
70  // initialize gammaKinematics when relevant.
71  if (beamHasGamma)
72  gammaKin.init(infoPtr, &settings, rndmPtr, beamAPtr, beamBPtr);
73 
74  // Initialize variables related to photon-inside-lepton.
75  beamGamAPtr = beamGamAPtrIn;
76  beamGamBPtr = beamGamBPtrIn;
77 
78  // Initialize variables related to VMD-inside-photon.
79  beamVMDAPtr = beamVMDAPtrIn;
80  beamVMDBPtr = beamVMDBPtrIn;
81 
82  // Send on some input pointers.
83  resonanceDecays.init( infoPtr, particleDataPtr, rndmPtr);
84 
85  // Set up SigmaTotal. Store sigma_nondiffractive for future use.
86  sigmaTotPtr->init( infoPtr, settings, particleDataPtr, rndmPtr);
87  int idA = infoPtr->idA();
88  int idB = infoPtr->idB();
89  double eCM = infoPtr->eCM();
90  if (beamHasGamma) {
91  int idAin = beamA2gamma ? 22 : idA;
92  int idBin = beamB2gamma ? 22 : idB;
93  sigmaTotPtr->calc( idAin, idBin, eCM);
94  }
95  else sigmaTotPtr->calc( idA, idB, eCM);
96  sigmaND = sigmaTotPtr->sigmaND();
97 
98  // Options to allow second hard interaction and resonance decays.
99  doSecondHard = settings.flag("SecondHard:generate");
100  doSameCuts = settings.flag("PhaseSpace:sameForSecond");
101  doResDecays = settings.flag("ProcessLevel:resonanceDecays");
102  startColTag = settings.mode("Event:startColTag");
103 
104  // Check whether ISR or MPI applied. Affects processes with photon beams.
105  doISR = ( settings.flag("PartonLevel:ISR")
106  && settings.flag("PartonLevel:all") );
107  doMPI = ( settings.flag("PartonLevel:MPI")
108  && settings.flag("PartonLevel:all") );
109 
110  // Second interaction not to be combined with biased phase space.
111  if (doSecondHard && userHooksPtr != 0
112  && userHooksPtr->canBiasSelection()) {
113  infoPtr->errorMsg("Error in ProcessLevel::init: "
114  "cannot combine second interaction with biased phase space");
115  return false;
116  }
117 
118  // Mass and pT cuts for two hard processes.
119  mHatMin1 = settings.parm("PhaseSpace:mHatMin");
120  mHatMax1 = settings.parm("PhaseSpace:mHatMax");
121  if (mHatMax1 < mHatMin1) mHatMax1 = eCM;
122  pTHatMin1 = settings.parm("PhaseSpace:pTHatMin");
123  pTHatMax1 = settings.parm("PhaseSpace:pTHatMax");
124  if (pTHatMax1 < pTHatMin1) pTHatMax1 = eCM;
125  mHatMin2 = settings.parm("PhaseSpace:mHatMinSecond");
126  mHatMax2 = settings.parm("PhaseSpace:mHatMaxSecond");
127  if (mHatMax2 < mHatMin2) mHatMax2 = eCM;
128  pTHatMin2 = settings.parm("PhaseSpace:pTHatMinSecond");
129  pTHatMax2 = settings.parm("PhaseSpace:pTHatMaxSecond");
130  if (pTHatMax2 < pTHatMin2) pTHatMax2 = eCM;
131 
132  // Check whether mass and pT ranges agree or overlap.
133  cutsAgree = doSameCuts;
134  if (mHatMin2 == mHatMin1 && mHatMax2 == mHatMax1 && pTHatMin2 == pTHatMin1
135  && pTHatMax2 == pTHatMax1) cutsAgree = true;
136  cutsOverlap = cutsAgree;
137  if (!cutsAgree) {
138  bool mHatOverlap = (max( mHatMin1, mHatMin2)
139  < min( mHatMax1, mHatMax2));
140  bool pTHatOverlap = (max( pTHatMin1, pTHatMin2)
141  < min( pTHatMax1, pTHatMax2));
142  if (mHatOverlap && pTHatOverlap) cutsOverlap = true;
143  }
144 
145  // Set up containers for all the internal hard processes.
146  SetupContainers setupContainers;
147  setupContainers.init(containerPtrs, infoPtr, settings, particleDataPtr,
148  couplingsPtr);
149 
150  // Append containers for external hard processes, if any.
151  if (sigmaPtrs.size() > 0) {
152  for (int iSig = 0; iSig < int(sigmaPtrs.size()); ++iSig)
153  containerPtrs.push_back( new ProcessContainer(sigmaPtrs[iSig],
154  true, phaseSpacePtrs[iSig]) );
155  }
156 
157  // Append single container for Les Houches processes, if any.
158  if (doLHA) {
159  SigmaProcess* sigmaPtr = new SigmaLHAProcess();
160  containerPtrs.push_back( new ProcessContainer(sigmaPtr) );
161 
162  // Store location of this container, and send in LHA pointer.
163  iLHACont = containerPtrs.size() - 1;
164  containerPtrs[iLHACont]->setLHAPtr(lhaUpPtr);
165  }
166 
167  // If no processes found then refuse to do anything.
168  if ( int(containerPtrs.size()) == 0) {
169  infoPtr->errorMsg("Error in ProcessLevel::init: "
170  "no process switched on");
171  return false;
172  }
173 
174  // Check whether pT-based weighting in 2 -> 2 is requested.
175  if (settings.flag("PhaseSpace:bias2Selection")) {
176  bool bias2Sel = false;
177  if (sigmaPtrs.size() == 0 && !doLHA && !doSecondHard) {
178  bias2Sel = true;
179  for (int i = 0; i < int(containerPtrs.size()); ++i) {
180  if (containerPtrs[i]->nFinal() != 2) bias2Sel = false;
181  int code = containerPtrs[i]->code();
182  if (code > 100 && code < 110) bias2Sel = false;
183  }
184  }
185  if (!bias2Sel) {
186  infoPtr->errorMsg("Error in ProcessLevel::init: "
187  "requested event weighting not possible");
188  return false;
189  }
190  }
191 
192  // Check that SUSY couplings were indeed initialized where necessary.
193  bool hasSUSY = false;
194  for (int i = 0; i < int(containerPtrs.size()); ++i)
195  if (containerPtrs[i]->isSUSY()) hasSUSY = true;
196 
197  // If SUSY processes requested but no SUSY couplings present
198  if (hasSUSY && !couplingsPtr->isSUSY) {
199  infoPtr->errorMsg("Error in ProcessLevel::init: "
200  "SUSY process switched on but no SUSY couplings found");
201  return false;
202  }
203 
204  // Fill SLHA blocks SMINPUTS and MASS from PYTHIA SM parameter values.
205  slhaInterfacePtr->pythia2slha(particleDataPtr);
206 
207  // Initialize each process.
208  int numberOn = 0;
209  for (int i = 0; i < int(containerPtrs.size()); ++i)
210  if (containerPtrs[i]->init(true, infoPtr, settings, particleDataPtr,
211  rndmPtr, beamAPtr, beamBPtr, couplingsPtr, sigmaTotPtr,
212  &resonanceDecays, slhaInterfacePtr, userHooksPtr, &gammaKin))
213  ++numberOn;
214 
215  // Sum maxima for Monte Carlo choice.
216  sigmaMaxSum = 0.;
217  for (int i = 0; i < int(containerPtrs.size()); ++i)
218  sigmaMaxSum += containerPtrs[i]->sigmaMax();
219 
220  // Option to pick a second hard interaction: repeat as above.
221  int number2On = 0;
222  if (doSecondHard) {
223  setupContainers.init2(container2Ptrs, settings);
224  if ( int(container2Ptrs.size()) == 0) {
225  infoPtr->errorMsg("Error in ProcessLevel::init: "
226  "no second hard process switched on");
227  return false;
228  }
229  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2)
230  if (container2Ptrs[i2]->init(false, infoPtr, settings, particleDataPtr,
231  rndmPtr, beamAPtr, beamBPtr, couplingsPtr, sigmaTotPtr,
232  &resonanceDecays, slhaInterfacePtr, userHooksPtr, &gammaKin))
233  ++number2On;
234 
235  sigma2MaxSum = 0.;
236  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2)
237  sigma2MaxSum += container2Ptrs[i2]->sigmaMax();
238  }
239 
240  // Printout during initialization is optional.
241  if (settings.flag("Init:showProcesses")) {
242 
243  // Construct string with incoming beams and for cm energy.
244  string collision = "We collide " + particleDataPtr->name(idA)
245  + " with " + particleDataPtr->name(idB) + " at a CM energy of ";
246  string pad( max( 0, 51 - int(collision.length())), ' ');
247 
248  // Print initialization information: header.
249  cout << "\n *------- PYTHIA Process Initialization ---------"
250  << "-----------------*\n"
251  << " | "
252  << " |\n"
253  << " | " << collision << scientific << setprecision(3)
254  << setw(9) << eCM << " GeV" << pad << " |\n"
255  << " | "
256  << " |\n"
257  << " |---------------------------------------------------"
258  << "---------------|\n"
259  << " | "
260  << " | |\n"
261  << " | Subprocess Code"
262  << " | Estimated |\n"
263  << " | "
264  << " | max (mb) |\n"
265  << " | "
266  << " | |\n"
267  << " |---------------------------------------------------"
268  << "---------------|\n"
269  << " | "
270  << " | |\n";
271 
272  // Loop over existing processes: print individual process info.
273  map<int, double> sigmaMaxM;
274  map<int, string> nameM;
275  for (int i = 0; i < int(containerPtrs.size()); ++i) {
276  int code = containerPtrs[i]->code();
277  nameM[code] = containerPtrs[i]->name();
278  sigmaMaxM[code] = containerPtrs[i]->sigmaMax() > sigmaMaxM[code] ?
279  containerPtrs[i]->sigmaMax() : sigmaMaxM[code];
280  }
281  for (map<int, string>::iterator i = nameM.begin(); i != nameM.end(); ++i)
282  cout << " | " << left << setw(45) << i->second
283  << right << setw(5) << i->first << " | "
284  << scientific << setprecision(3) << setw(11)
285  << sigmaMaxM[i->first] << " |\n";
286 
287  // Loop over second hard processes, if any, and repeat as above.
288  if (doSecondHard) {
289  cout << " | "
290  << " | |\n"
291  << " |---------------------------------------------------"
292  <<"---------------|\n"
293  << " | "
294  <<" | |\n";
295  sigmaMaxM.clear();
296  nameM.clear();
297  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2) {
298  int code = container2Ptrs[i2]->code();
299  nameM[code] = container2Ptrs[i2]->name();
300  sigmaMaxM[code] = container2Ptrs[i2]->sigmaMax() > sigmaMaxM[code] ?
301  container2Ptrs[i2]->sigmaMax() : sigmaMaxM[code];
302  }
303  for (map<int, string>::iterator i2 = nameM.begin(); i2 != nameM.end();
304  ++i2)
305  cout << " | " << left << setw(45) << i2->second
306  << right << setw(5) << i2->first << " | "
307  << scientific << setprecision(3) << setw(11)
308  << sigmaMaxM[i2->first] << " |\n";
309  }
310 
311  // Listing finished.
312  cout << " | "
313  << " |\n"
314  << " *------- End PYTHIA Process Initialization ----------"
315  <<"-------------*" << endl;
316  }
317 
318  // If sum of maxima vanishes then refuse to do anything.
319  if ( numberOn == 0 || sigmaMaxSum <= 0.) {
320  infoPtr->errorMsg("Error in ProcessLevel::init: "
321  "all processes have vanishing cross sections");
322  return false;
323  }
324  if ( doSecondHard && (number2On == 0 || sigma2MaxSum <= 0.) ) {
325  infoPtr->errorMsg("Error in ProcessLevel::init: "
326  "all second hard processes have vanishing cross sections");
327  return false;
328  }
329 
330  // If two hard processes then check whether some (but not all) agree.
331  allHardSame = true;
332  noneHardSame = true;
333  if (doSecondHard) {
334  bool foundMatch = false;
335 
336  // Check for each first process if matched in second.
337  for (int i = 0; i < int(containerPtrs.size()); ++i) {
338  foundMatch = false;
339  if (cutsOverlap)
340  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2)
341  if (container2Ptrs[i2]->code() == containerPtrs[i]->code())
342  foundMatch = true;
343  containerPtrs[i]->isSame( foundMatch );
344  if (!foundMatch) allHardSame = false;
345  if ( foundMatch) noneHardSame = false;
346  }
347 
348  // Check for each second process if matched in first.
349  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2) {
350  foundMatch = false;
351  if (cutsOverlap)
352  for (int i = 0; i < int(containerPtrs.size()); ++i)
353  if (containerPtrs[i]->code() == container2Ptrs[i2]->code())
354  foundMatch = true;
355  container2Ptrs[i2]->isSame( foundMatch );
356  if (!foundMatch) allHardSame = false;
357  if ( foundMatch) noneHardSame = false;
358  }
359  }
360 
361  // Concluding classification, including cuts.
362  if (!cutsAgree) allHardSame = false;
363  someHardSame = !allHardSame && !noneHardSame;
364 
365  // Done.
366  return true;
367 }
368 
369 //--------------------------------------------------------------------------
370 
371 // Main routine to generate the hard process.
372 
373 bool ProcessLevel::next( Event& process) {
374 
375  // Generate the next event with two or one hard interactions.
376  bool physical = (doSecondHard) ? nextTwo( process) : nextOne( process);
377 
378  // Check that colour assignments make sense.
379  if (physical) physical = checkColours( process);
380 
381  // Done.
382  return physical;
383 }
384 
385 //--------------------------------------------------------------------------
386 
387 // Generate (= read in) LHA input of resonance decay only.
388 
389 bool ProcessLevel::nextLHAdec( Event& process) {
390 
391  // Read resonance decays from LHA interface.
392  infoPtr->setEndOfFile(false);
393  if (!lhaUpPtr->setEvent()) {
394  infoPtr->setEndOfFile(true);
395  return false;
396  }
397 
398  // Store LHA output in standard event record format.
399  containerLHAdec.constructDecays( process);
400 
401  // Done.
402  return true;
403 }
404 
405 //--------------------------------------------------------------------------
406 
407 // Accumulate and update statistics (after possible user veto).
408 
409 void ProcessLevel::accumulate( bool doAccumulate) {
410 
411  // Increase number of accepted events.
412  if (doAccumulate) containerPtrs[iContainer]->accumulate();
413 
414  // Provide current generated cross section estimate.
415  long nTrySum = 0;
416  long nSelSum = 0;
417  long nAccSum = 0;
418  double sigmaSum = 0.;
419  double delta2Sum = 0.;
420  double sigSelSum = 0.;
421  double weightSum = 0.;
422  int codeNow;
423  long nTryNow, nSelNow, nAccNow;
424  double sigmaNow, deltaNow, sigSelNow, weightNow;
425  map<int, bool> duplicate;
426  for (int i = 0; i < int(containerPtrs.size()); ++i)
427  if (containerPtrs[i]->sigmaMax() != 0.) {
428  codeNow = containerPtrs[i]->code();
429  nTryNow = containerPtrs[i]->nTried();
430  nSelNow = containerPtrs[i]->nSelected();
431  nAccNow = containerPtrs[i]->nAccepted();
432  sigmaNow = containerPtrs[i]->sigmaMC(doAccumulate);
433  deltaNow = containerPtrs[i]->deltaMC(doAccumulate);
434  sigSelNow = containerPtrs[i]->sigmaSelMC(doAccumulate);
435  weightNow = containerPtrs[i]->weightSum();
436  nTrySum += nTryNow;
437  nSelSum += nSelNow;
438  nAccSum += nAccNow;
439  sigmaSum += sigmaNow;
440  delta2Sum += pow2(deltaNow);
441  sigSelSum += sigSelNow;
442  weightSum += weightNow;
443  if (!doSecondHard) {
444  if (!duplicate[codeNow])
445  infoPtr->setSigma( codeNow, containerPtrs[i]->name(),
446  nTryNow, nSelNow, nAccNow, sigmaNow, deltaNow, weightNow);
447  else
448  infoPtr->addSigma( codeNow, nTryNow, nSelNow, nAccNow, sigmaNow,
449  deltaNow);
450  duplicate[codeNow] = true;
451  }
452  }
453 
454  // Normally only one hard interaction. Then store info and done.
455  if (!doSecondHard) {
456  double deltaSum = sqrtpos(delta2Sum);
457  infoPtr->setSigma( 0, "sum", nTrySum, nSelSum, nAccSum, sigmaSum, deltaSum,
458  weightSum);
459  return;
460  }
461 
462  // Increase counter for a second hard interaction.
463  if (doAccumulate) container2Ptrs[i2Container]->accumulate();
464 
465  // Cross section estimate for second hard process.
466  double sigma2Sum = 0.;
467  double sig2SelSum = 0.;
468  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2)
469  if (container2Ptrs[i2]->sigmaMax() != 0.) {
470  nTrySum += container2Ptrs[i2]->nTried();
471  if (doAccumulate) sigma2Sum += container2Ptrs[i2]->sigmaMC();
472  if (doAccumulate) sig2SelSum += container2Ptrs[i2]->sigmaSelMC();
473  }
474 
475  // Average impact-parameter factor.
476  double impactFac = max( 1., infoPtr->enhanceMPIavg() );
477 
478  // Cross section estimate for combination of first and second process.
479  // Combine two possible ways and take average.
480  double sigmaComb = 0.5 * (sigmaSum * sig2SelSum + sigSelSum * sigma2Sum);
481  sigmaComb *= impactFac / sigmaND;
482  if (allHardSame) sigmaComb *= 0.5;
483  double deltaComb = (nAccSum == 0) ? 0. : sqrtpos(2. / nAccSum) * sigmaComb;
484 
485  // Store info and done.
486  infoPtr->setSigma( 0, "sum", nTrySum, nSelSum, nAccSum, sigmaComb, deltaComb,
487  weightSum);
488 
489 }
490 
491 //--------------------------------------------------------------------------
492 
493 // Print statistics on cross sections and number of events.
494 
495 void ProcessLevel::statistics(bool reset) {
496 
497  // Special processing if two hard interactions selected.
498  if (doSecondHard) {
499  statistics2(reset);
500  return;
501  }
502 
503  // Header.
504  cout << "\n *------- PYTHIA Event and Cross Section Statistics ------"
505  << "-------------------------------------------------------*\n"
506  << " | "
507  << " |\n"
508  << " | Subprocess Code | "
509  << " Number of events | sigma +- delta |\n"
510  << " | | "
511  << "Tried Selected Accepted | (estimated) (mb) |\n"
512  << " | | "
513  << " | |\n"
514  << " |------------------------------------------------------------"
515  << "-----------------------------------------------------|\n"
516  << " | | "
517  << " | |\n";
518 
519  // Reset sum counters.
520  long nTrySum = 0;
521  long nSelSum = 0;
522  long nAccSum = 0;
523  double sigmaSum = 0.;
524  double delta2Sum = 0.;
525 
526  // Reset process maps.
527  map<int, string> nameM;
528  map<int, long> nTryM, nSelM, nAccM;
529  map<int, double> sigmaM, delta2M;
530  vector<ProcessContainer*> lheContainerPtrs;
531 
532  // Loop over existing processes.
533  for (int i = 0; i < int(containerPtrs.size()); ++i)
534  if (containerPtrs[i]->sigmaMax() != 0.) {
535 
536  // Read info for process. Sum counters.
537  nTrySum += containerPtrs[i]->nTried();
538  nSelSum += containerPtrs[i]->nSelected();
539  nAccSum += containerPtrs[i]->nAccepted();
540  sigmaSum += containerPtrs[i]->sigmaMC();
541  delta2Sum += pow2(containerPtrs[i]->deltaMC());
542 
543  // Skip Les Houches containers.
544  if (containerPtrs[i]->code() == 9999) {
545  lheContainerPtrs.push_back(containerPtrs[i]);
546  continue;
547  }
548 
549  // Internal process info.
550  int code = containerPtrs[i]->code();
551  nameM[code] = containerPtrs[i]->name();
552  nTryM[code] += containerPtrs[i]->nTried();
553  nSelM[code] += containerPtrs[i]->nSelected();
554  nAccM[code] += containerPtrs[i]->nAccepted();
555  sigmaM[code] += containerPtrs[i]->sigmaMC();
556  delta2M[code]+= pow2(containerPtrs[i]->deltaMC());
557  }
558 
559  // Print internal process info.
560  for (map<int, string>::iterator i = nameM.begin(); i != nameM.end(); ++i) {
561  int code = i->first;
562  cout << " | " << left << setw(45) << i->second
563  << right << setw(5) << code << " | "
564  << setw(11) << nTryM[code] << " " << setw(10) << nSelM[code] << " "
565  << setw(10) << nAccM[code] << " | " << scientific << setprecision(3)
566  << setw(11) << sigmaM[code]
567  << setw(11) << sqrtpos(delta2M[code]) << " |\n";
568  }
569 
570  // Print Les Houches process info.
571  for (int i = 0; i < int(lheContainerPtrs.size()); ++i) {
572  ProcessContainer *ptr = lheContainerPtrs[i];
573  cout << " | " << left << setw(45) << ptr->name()
574  << right << setw(5) << ptr->code() << " | "
575  << setw(11) << ptr->nTried() << " " << setw(10) << ptr->nSelected()
576  << " " << setw(10) << ptr->nAccepted() << " | " << scientific
577  << setprecision(3) << setw(11) << ptr->sigmaMC() << setw(11)
578  << ptr->deltaMC() << " |\n";
579 
580  // Print subdivision by user code for Les Houches process.
581  for (int j = 0; j < ptr->codeLHASize(); ++j)
582  cout << " | ... whereof user classification code " << setw(10)
583  << ptr->subCodeLHA(j) << " | " << setw(11) << ptr->nTriedLHA(j)
584  << " " << setw(10) << ptr->nSelectedLHA(j) << " " << setw(10)
585  << ptr->nAcceptedLHA(j) << " | | \n";
586  }
587 
588  // Print summed process info.
589  cout << " | | "
590  << " | |\n"
591  << " | " << left << setw(50) << "sum" << right << " | " << setw(11)
592  << nTrySum << " " << setw(10) << nSelSum << " " << setw(10)
593  << nAccSum << " | " << scientific << setprecision(3) << setw(11)
594  << sigmaSum << setw(11) << sqrtpos(delta2Sum) << " |\n";
595 
596  // Listing finished.
597  cout << " | "
598  << " |\n"
599  << " *------- End PYTHIA Event and Cross Section Statistics -----"
600  << "-----------------------------------------------------*" << endl;
601 
602  // Optionally reset statistics contants.
603  if (reset) resetStatistics();
604 
605 }
606 
607 //--------------------------------------------------------------------------
608 
609 // Reset statistics on cross sections and number of events.
610 
611 void ProcessLevel::resetStatistics() {
612 
613  for (int i = 0; i < int(containerPtrs.size()); ++i)
614  containerPtrs[i]->reset();
615  if (doSecondHard)
616  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2)
617  container2Ptrs[i2]->reset();
618 
619 }
620 
621 //--------------------------------------------------------------------------
622 
623 // Generate the next event with one interaction.
624 
625 bool ProcessLevel::nextOne( Event& process) {
626 
627  // Update CM energy for phase space selection.
628  double eCM = infoPtr->eCM();
629  for (int i = 0; i < int(containerPtrs.size()); ++i)
630  containerPtrs[i]->newECM(eCM);
631 
632  // Outer loop in case of rare failures.
633  bool physical = true;
634  for (int loop = 0; loop < MAXLOOP; ++loop) {
635  if (!physical) process.clear();
636  physical = true;
637 
638  // Loop over tries until trial event succeeds.
639  for ( ; ; ) {
640 
641  // Pick one of the subprocesses.
642  double sigmaMaxNow = sigmaMaxSum * rndmPtr->flat();
643  int iMax = containerPtrs.size() - 1;
644  iContainer = -1;
645  do sigmaMaxNow -= containerPtrs[++iContainer]->sigmaMax();
646  while (sigmaMaxNow > 0. && iContainer < iMax);
647 
648  // Do a trial event of this subprocess; accept or not.
649  if (containerPtrs[iContainer]->trialProcess()) break;
650 
651  // Check for end-of-file condition for Les Houches events.
652  if (infoPtr->atEndOfFile()) return false;
653  }
654 
655  // Update sum of maxima if current maximum violated.
656  if (containerPtrs[iContainer]->newSigmaMax()) {
657  sigmaMaxSum = 0.;
658  for (int i = 0; i < int(containerPtrs.size()); ++i)
659  sigmaMaxSum += containerPtrs[i]->sigmaMax();
660  }
661 
662  // Construct kinematics of acceptable process.
663  containerPtrs[iContainer]->constructState();
664  if ( !containerPtrs[iContainer]->constructProcess( process) )
665  physical = false;
666 
667  // For photon beams from leptons copy the state to additional photon beams.
668  if (beamHasGamma) {
669  beamGamAPtr->setGammaMode(beamAPtr->getGammaMode());
670  beamGamBPtr->setGammaMode(beamBPtr->getGammaMode());
671  }
672 
673  // Do all resonance decays.
674  if ( physical && doResDecays
675  && !containerPtrs[iContainer]->decayResonances( process) )
676  physical = false;
677 
678  // Retry process for unphysical states.
679  for (int i =1; i < process.size(); ++i)
680  if (process[i].e() < 0.) {
681  infoPtr->errorMsg("Error in ProcessLevel::nextOne: "
682  "Constructed particle with negative energy.");
683  physical = false;
684  }
685 
686  // Add any junctions to the process event record list.
687  if (physical) findJunctions( process);
688 
689  // Check that enough room for beam remnants in the photon beams and
690  // set the valence content for photon beams.
691  // Do not check for softQCD processes since no initiators yet.
692  if ( ( ( ( beamAPtr->isGamma() && !beamAPtr->isUnresolved() )
693  || ( beamBPtr->isGamma() && !beamBPtr->isUnresolved() ) )
694  || ( beamAPtr->hasResGamma() || beamBPtr->hasResGamma() ) )
695  && !containerPtrs[iContainer]->isSoftQCD() ) {
696  if ( !roomForRemnants() ) {
697  physical = false;
698  continue;
699  }
700  }
701 
702  // Outer loop should normally work first time around.
703  if (physical) break;
704  }
705 
706  // Update information in VMD beams according to chosen state.
707  if (infoPtr->isVMDstateA()) {
708  beamVMDAPtr->setGammaMode(beamAPtr->getGammaMode());
709  beamVMDAPtr->setVMDstate(true, infoPtr->idVMDA(), infoPtr->mVMDA(),
710  infoPtr->scaleVMDA(), true);
711  }
712  if (infoPtr->isVMDstateB()) {
713  beamVMDBPtr->setGammaMode(beamBPtr->getGammaMode());
714  beamVMDBPtr->setVMDstate(true, infoPtr->idVMDB(), infoPtr->mVMDB(),
715  infoPtr->scaleVMDB(), true);
716  }
717 
718  // Done.
719  return physical;
720 }
721 
722 //--------------------------------------------------------------------------
723 
724 // Generate the next event with two hard interactions.
725 
726 bool ProcessLevel::nextTwo( Event& process) {
727 
728  // Update CM energy for phase space selection.
729  double eCM = infoPtr->eCM();
730  for (int i = 0; i < int(containerPtrs.size()); ++i)
731  containerPtrs[i]->newECM(eCM);
732  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2)
733  container2Ptrs[i2]->newECM(eCM);
734 
735  // Outer loop in case of rare failures.
736  bool physical = true;
737  for (int loop = 0; loop < MAXLOOP; ++loop) {
738  if (!physical) process.clear();
739  physical = true;
740 
741  // Loop over both hard processes to find consistent common kinematics.
742  for ( ; ; ) {
743 
744  // Loop internally over tries for hardest process until succeeds.
745  for ( ; ; ) {
746 
747  // Pick one of the subprocesses.
748  double sigmaMaxNow = sigmaMaxSum * rndmPtr->flat();
749  int iMax = containerPtrs.size() - 1;
750  iContainer = -1;
751  do sigmaMaxNow -= containerPtrs[++iContainer]->sigmaMax();
752  while (sigmaMaxNow > 0. && iContainer < iMax);
753 
754  // Do a trial event of this subprocess; accept or not.
755  if (containerPtrs[iContainer]->trialProcess()) break;
756 
757  // Check for end-of-file condition for Les Houches events.
758  if (infoPtr->atEndOfFile()) return false;
759  }
760 
761  // Update sum of maxima if current maximum violated.
762  if (containerPtrs[iContainer]->newSigmaMax()) {
763  sigmaMaxSum = 0.;
764  for (int i = 0; i < int(containerPtrs.size()); ++i)
765  sigmaMaxSum += containerPtrs[i]->sigmaMax();
766  }
767 
768  // Loop internally over tries for second hardest process until succeeds.
769  for ( ; ; ) {
770 
771  // Pick one of the subprocesses.
772  double sigma2MaxNow = sigma2MaxSum * rndmPtr->flat();
773  int i2Max = container2Ptrs.size() - 1;
774  i2Container = -1;
775  do sigma2MaxNow -= container2Ptrs[++i2Container]->sigmaMax();
776  while (sigma2MaxNow > 0. && i2Container < i2Max);
777 
778  // Do a trial event of this subprocess; accept or not.
779  if (container2Ptrs[i2Container]->trialProcess()) break;
780  }
781 
782  // Update sum of maxima if current maximum violated.
783  if (container2Ptrs[i2Container]->newSigmaMax()) {
784  sigma2MaxSum = 0.;
785  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2)
786  sigma2MaxSum += container2Ptrs[i2]->sigmaMax();
787  }
788 
789  // Pick incoming flavours (etc), needed for PDF reweighting.
790  containerPtrs[iContainer]->constructState();
791  container2Ptrs[i2Container]->constructState();
792 
793  // Check whether common set of x values is kinematically possible.
794  double xA1 = containerPtrs[iContainer]->x1();
795  double xB1 = containerPtrs[iContainer]->x2();
796  double xA2 = container2Ptrs[i2Container]->x1();
797  double xB2 = container2Ptrs[i2Container]->x2();
798  if (xA1 + xA2 >= 1. || xB1 + xB2 >= 1.) continue;
799 
800  // Reset beam contents. Naive parton densities for second interaction.
801  // (Subsequent procedure could be symmetrized, but would be overkill.)
802  beamAPtr->clear();
803  beamBPtr->clear();
804  int idA1 = containerPtrs[iContainer]->id1();
805  int idB1 = containerPtrs[iContainer]->id2();
806  int idA2 = container2Ptrs[i2Container]->id1();
807  int idB2 = container2Ptrs[i2Container]->id2();
808  double Q2Fac1 = containerPtrs[iContainer]->Q2Fac();
809  double Q2Fac2 = container2Ptrs[i2Container]->Q2Fac();
810  double pdfA2Raw = beamAPtr->xf( idA2, xA2,Q2Fac2);
811  double pdfB2Raw = beamBPtr->xf( idB2, xB2,Q2Fac2);
812 
813  // Remove partons in first interaction from beams.
814  beamAPtr->append( 3, idA1, xA1);
815  beamAPtr->xfISR( 0, idA1, xA1, Q2Fac1);
816  beamAPtr->pickValSeaComp();
817  beamBPtr->append( 4, idB1, xB1);
818  beamBPtr->xfISR( 0, idB1, xB1, Q2Fac1);
819  beamBPtr->pickValSeaComp();
820 
821  // Reevaluate pdf's for second interaction and weight by reduction.
822  double pdfA2Mod = beamAPtr->xfMPI( idA2, xA2,Q2Fac2);
823  double pdfB2Mod = beamBPtr->xfMPI( idB2, xB2,Q2Fac2);
824  double wtPdfMod = (pdfA2Mod * pdfB2Mod) / (pdfA2Raw * pdfB2Raw);
825  if (wtPdfMod < rndmPtr->flat()) continue;
826 
827  // Reduce by a factor of 2 for identical processes when others not,
828  // and when in same phase space region.
829  bool toLoop = false;
830  if ( someHardSame && containerPtrs[iContainer]->isSame()
831  && container2Ptrs[i2Container]->isSame()) {
832  if (cutsAgree) {
833  if (rndmPtr->flat() > 0.5) toLoop = true;
834  } else {
835  double mHat1 = containerPtrs[iContainer]->mHat();
836  double pTHat1 = containerPtrs[iContainer]->pTHat();
837  double mHat2 = container2Ptrs[i2Container]->mHat();
838  double pTHat2 = container2Ptrs[i2Container]->pTHat();
839  if (mHat1 > mHatMin2 && mHat1 < mHatMax2
840  && pTHat1 > pTHatMin2 && pTHat1 < pTHatMax2
841  && mHat2 > mHatMin1 && mHat2 < mHatMax1
842  && pTHat2 > pTHatMin1 && pTHat2 < pTHatMax1
843  && rndmPtr->flat() > 0.5) toLoop = true;
844  }
845  }
846  if (toLoop) continue;
847 
848  // If come this far then acceptable event.
849  break;
850  }
851 
852  // Construct kinematics of acceptable processes.
853  Event process2;
854  process2.init( "(second hard)", particleDataPtr, startColTag);
855  process2.initColTag();
856  if ( !containerPtrs[iContainer]->constructProcess( process) )
857  physical = false;
858  if (physical && !container2Ptrs[i2Container]->constructProcess( process2,
859  false) ) physical = false;
860 
861  // Do all resonance decays.
862  if ( physical && doResDecays
863  && !containerPtrs[iContainer]->decayResonances( process) )
864  physical = false;
865  if ( physical && doResDecays
866  && !container2Ptrs[i2Container]->decayResonances( process2) )
867  physical = false;
868 
869  // Append second hard interaction to normal process object.
870  if (physical) combineProcessRecords( process, process2);
871 
872  // Add any junctions to the process event record list.
873  if (physical) findJunctions( process);
874 
875  // Outer loop should normally work first time around.
876  if (physical) break;
877  }
878 
879  // Done.
880  return physical;
881 }
882 
883 //--------------------------------------------------------------------------
884 
885 // Check that enough room for beam remnants is left and fix the valence
886 // content for photon beams.
887 
888 bool ProcessLevel::roomForRemnants() {
889 
890  // Set the beamParticle pointers to photon beams.
891  bool beamAhasResGamma = beamAPtr->hasResGamma();
892  bool beamBhasResGamma = beamBPtr->hasResGamma();
893  bool beamHasResGamma = beamAhasResGamma || beamBhasResGamma;
894  BeamParticle* tmpBeamAPtr = beamAhasResGamma ? beamGamAPtr : beamAPtr;
895  BeamParticle* tmpBeamBPtr = beamBhasResGamma ? beamGamBPtr : beamBPtr;
896 
897  // Check whether photons are unresolved.
898  bool resGammaA = !(beamGamAPtr->isUnresolved());
899  bool resGammaB = !(beamGamBPtr->isUnresolved());
900 
901  // Get the x_gamma values from beam particle.
902  double xGamma1 = beamAPtr->xGamma();
903  double xGamma2 = beamBPtr->xGamma();
904 
905  // Clear the previous choice.
906  tmpBeamAPtr->posVal(-1);
907  tmpBeamBPtr->posVal(-1);
908 
909  // Store the relevant information.
910  int id1 = containerPtrs[iContainer]->id1();
911  int id2 = containerPtrs[iContainer]->id2();
912  double x1 = containerPtrs[iContainer]->x1();
913  double x2 = containerPtrs[iContainer]->x2();
914  double Q2 = containerPtrs[iContainer]->Q2Fac();
915 
916  // Invariant mass of gamma-gamma system.
917  double eCMgmgm = infoPtr->eCM();
918 
919  // If photon(s) from beam particle(s) use the derived gmgm CM energy.
920  if (beamHasResGamma) eCMgmgm = infoPtr->eCMsub();
921 
922  // Calculate the mT left for the beams after hard interaction.
923  double mTRem = 0;
924 
925  // Direct-resolved processes with photons from lepton beams.
926  if ( ( (resGammaA && !resGammaB) || (!resGammaA && resGammaB) )
927  && beamHasResGamma ) {
928  double wTot = infoPtr->eCMsub();
929  double w2scat = infoPtr->sHatNew();
930  mTRem = wTot - sqrt(w2scat);
931 
932  // Direct-resolved processes with photons beams.
933  } else if ( tmpBeamAPtr->isUnresolved() && !tmpBeamBPtr->isUnresolved() ) {
934  mTRem = eCMgmgm*( 1. - sqrt( x2) );
935  } else if ( !tmpBeamAPtr->isUnresolved() && tmpBeamBPtr->isUnresolved() ) {
936  mTRem = eCMgmgm*( 1. - sqrt( x1) );
937 
938  // Resolved-resolved processes.
939  } else {
940 
941  // Photons from lepton beams.
942  if ( beamAhasResGamma && beamBhasResGamma) {
943 
944  // Rescale the x_gamma values according to the new invariant mass
945  // of the gamma+gamma system and rescale x values.
946  double sCM = infoPtr->s();
947  x1 /= xGamma1 * pow2(eCMgmgm) / ( xGamma1 * xGamma2 * sCM);
948  x2 /= xGamma2 * pow2(eCMgmgm) / ( xGamma1 * xGamma2 * sCM);
949  }
950 
951  // Invariant mass left with resolved photons.
952  mTRem = eCMgmgm * sqrt( (1. - x1)*(1. - x2) );
953  }
954 
955  // Initial values.
956  double m1 = 0, m2 = 0;
957  bool physical = false;
958 
959  // If no ISR or MPI decide the valence content according to the hard process.
960  if ( !doISR && !doMPI ) {
961 
962  // For physical process a few tries for the valence content should suffice.
963  int nTry = 0, maxTry = 4;
964 
965  while (!physical) {
966 
967  // Check whether the hard parton is a valence quark and if not use the
968  // parametrization for the valence flavor ratios.
969  bool init1Val = tmpBeamAPtr->isGamma() ?
970  tmpBeamAPtr->gammaInitiatorIsVal(0, id1, x1, Q2) : false;
971  bool init2Val = tmpBeamBPtr->isGamma() ?
972  tmpBeamBPtr->gammaInitiatorIsVal(0, id2, x2, Q2) : false;
973 
974  // If no ISR is generated must leave room for beam remnants.
975  // Calculate the required amount of energy for three different cases:
976  // Hard parton is a valence quark: Remnant mass = m_val.
977  // Hard parton is a gluon: Remnant mass = 2*m_val.
978  // Hard parton is not valence quark: Remnant mass = 2*m_val + m_hard.
979  // Unresolved photon: Remnant mass = 0.
980  if ( tmpBeamAPtr->isGamma() ) {
981  if ( !resGammaA) {
982  m1 = 0;
983  } else if (init1Val) {
984  m1 = particleDataPtr->m0(id1);
985  } else {
986  m1 = 2*( particleDataPtr->m0( tmpBeamAPtr->getGammaValFlavour() ) );
987  if (id1 != 21) m1 += particleDataPtr->m0(id1);
988  }
989 
990  // For hadrons start with hadron mass.
991  } else if ( tmpBeamAPtr->isHadron() ) {
992  m1 = particleDataPtr->m0(tmpBeamAPtr->id());
993 
994  // If a valence flavor, remove the mass from remnants, else add.
995  int valSign1 = (tmpBeamAPtr->nValence(id1) > 0) ? -1 : 1;
996  m1 = m1 + valSign1 * particleDataPtr->m0(id1);
997  }
998 
999  // Photons.
1000  if ( tmpBeamBPtr->isGamma() ) {
1001  if ( !resGammaB) {
1002  m2 = 0;
1003  } else if (init2Val) {
1004  m2 = particleDataPtr->m0(id2);
1005  } else {
1006  m2 = 2*( particleDataPtr->m0( tmpBeamBPtr->getGammaValFlavour() ) );
1007  if (id2 != 21) m2 += particleDataPtr->m0(id2);
1008  }
1009 
1010  // For hadrons start with hadron mass.
1011  } else if ( tmpBeamBPtr->isHadron() ) {
1012  m2 = particleDataPtr->m0(tmpBeamBPtr->id());
1013 
1014  // If a valence flavor, remove the mass from remnants, else add.
1015  int valSign2 = (tmpBeamBPtr->nValence(id2) > 0) ? -1 : 1;
1016  m2 = m2 + valSign2 * particleDataPtr->m0(id2);
1017  }
1018 
1019  // Check whether room for remnants.
1020  physical = ( (m1 + m2) < mTRem );
1021 
1022  // Check that maximum number of trials not exceeded.
1023  ++nTry;
1024  if (nTry == maxTry) {
1025  if ( abs(id1) == 5 || abs(id2) == 5 ) {
1026  infoPtr->errorMsg("Warning in ProcessLevel::roomForRemnants: "
1027  "No room for bottom quarks in beam remnants");
1028  } else if ( abs(id1) == 4 || abs(id2) == 4 ) {
1029  infoPtr->errorMsg("Warning in ProcessLevel::roomForRemnants: "
1030  "No room for charm quarks in beam remnants");
1031  } else {
1032  infoPtr->errorMsg("Warning in ProcessLevel::roomForRemnants: "
1033  "No room for light quarks in beam remnants");
1034  }
1035  break;
1036  }
1037  }
1038 
1039  // With ISR or MPI the initiator list can change so reject only
1040  // process which will surely fail.
1041  } else {
1042 
1043  // Do not allow processes that ISR cannot turn into physical one.
1044  int idLight = 2;
1045 
1046  // Side A with a photon or hadron.
1047  if ( tmpBeamAPtr->isGamma() ) {
1048  m1 = (id1 == 21) ? 2*( particleDataPtr->m0( idLight ) ) :
1049  particleDataPtr->m0( id1 );
1050  } else if ( tmpBeamAPtr->isHadron() ) {
1051  m1 = particleDataPtr->m0( tmpBeamAPtr->id() );
1052  int valSign1 = (tmpBeamAPtr->nValence(id1) > 0) ? -1 : 1;
1053  m1 = m1 + valSign1 * particleDataPtr->m0(id1);
1054  }
1055 
1056  // Side B with a photon or hadron.
1057  if ( tmpBeamBPtr->isGamma() ) {
1058  m2 = (id2 == 21) ? 2*( particleDataPtr->m0( idLight ) ) :
1059  particleDataPtr->m0( id2 );
1060  } else if ( tmpBeamBPtr->isHadron() ) {
1061  m2 = particleDataPtr->m0( tmpBeamBPtr->id() );
1062  int valSign2 = (tmpBeamBPtr->nValence(id2) > 0) ? -1 : 1;
1063  m2 = m2 + valSign2 * particleDataPtr->m0(id2);
1064  }
1065 
1066  // No remnants needed for direct photons.
1067  if ( !resGammaA && !(tmpBeamAPtr->isHadron()) ) m1 = 0;
1068  if ( !resGammaB && !(tmpBeamBPtr->isHadron()) ) m2 = 0;
1069 
1070  // Check whether room for remnants.
1071  physical = ( (m1 + m2) < mTRem );
1072  if (physical == false) {
1073  if ( abs(id1) == 5 || abs(id2) == 5 ) {
1074  infoPtr->errorMsg("Warning in ProcessLevel::roomForRemnants: "
1075  "No room for bottom quarks in beam remnants");
1076  } else if ( abs(id1) == 4 || abs(id2) == 4 ) {
1077  infoPtr->errorMsg("Warning in ProcessLevel::roomForRemnants: "
1078  "No room for charm quarks in beam remnants");
1079  } else {
1080  infoPtr->errorMsg("Warning in ProcessLevel::roomForRemnants: "
1081  "No room for light quarks in beam remnants");
1082  }
1083  }
1084  }
1085 
1086  return physical;
1087 }
1088 
1089 //--------------------------------------------------------------------------
1090 
1091 // Append second hard interaction to normal process object.
1092 // Complication: all resonance decay chains must be put at the end.
1093 
1094 void ProcessLevel::combineProcessRecords( Event& process, Event& process2) {
1095 
1096  // Find first event record size, excluding resonances.
1097  int nSize = process.size();
1098  int nHard = 5;
1099  while (nHard < nSize && process[nHard].mother1() == 3) ++nHard;
1100 
1101  // Save resonance products temporarily elsewhere.
1102  vector<Particle> resProd;
1103  if (nSize > nHard) {
1104  for (int i = nHard; i < nSize; ++i) resProd.push_back( process[i] );
1105  process.popBack(nSize - nHard);
1106  }
1107 
1108  // Find second event record size, excluding resonances.
1109  int nSize2 = process2.size();
1110  int nHard2 = 5;
1111  while (nHard2 < nSize2 && process2[nHard2].mother1() == 3) ++nHard2;
1112 
1113  // Find amount of necessary position and colour offset for second process.
1114  int addPos = nHard - 3;
1115  int addCol = process.lastColTag() - startColTag;
1116 
1117  // Loop over all particles (except beams) from second process.
1118  for (int i = 3; i < nSize2; ++i) {
1119 
1120  // Offset mother and daughter pointers and colour tags of particle.
1121  process2[i].offsetHistory( 2, addPos, 2, addPos);
1122  process2[i].offsetCol( addCol);
1123 
1124  // Append hard-process particles from process2 to process.
1125  if (i < nHard2) process.append( process2[i] );
1126  }
1127 
1128  // Reinsert resonance decay chains of first hard process.
1129  int addPos2 = nHard2 - 3;
1130  if (nHard < nSize) {
1131 
1132  // Offset daughter pointers of unmoved mothers.
1133  for (int i = 5; i < nHard; ++i)
1134  process[i].offsetHistory( 0, 0, nHard - 1, addPos2);
1135 
1136  // Modify history of resonance products when restoring.
1137  for (int i = 0; i < int(resProd.size()); ++i) {
1138  resProd[i].offsetHistory( nHard - 1, addPos2, nHard - 1, addPos2);
1139  process.append( resProd[i] );
1140  }
1141  }
1142 
1143  // Insert resonance decay chains of second hard process.
1144  if (nHard2 < nSize2) {
1145  int nHard3 = nHard + nHard2 - 3;
1146  int addPos3 = nSize - nHard;
1147 
1148  // Offset daughter pointers of second-process mothers.
1149  for (int i = nHard + 2; i < nHard3; ++i)
1150  process[i].offsetHistory( 0, 0, nHard3 - 1, addPos3);
1151 
1152  // Modify history of second-process resonance products and insert.
1153  for (int i = nHard2; i < nSize2; ++i) {
1154  process2[i].offsetHistory( nHard3 - 1, addPos3, nHard3 - 1, addPos3);
1155  process.append( process2[i] );
1156  }
1157  }
1158 
1159  // Store PDF scale for second interaction.
1160  process.scaleSecond( process2.scale() );
1161 
1162 }
1163 
1164 //--------------------------------------------------------------------------
1165 
1166 // Add any junctions to the process event record list.
1167 // Also check that do not doublebook if called repeatedly.
1168 
1169 void ProcessLevel::findJunctions( Event& junEvent) {
1170 
1171  // Check all hard vertices for BNV
1172  for (int i = 1; i<junEvent.size(); i++) {
1173 
1174  // Ignore colorless particles and stages before hard-scattering
1175  // final state.
1176  if (abs(junEvent[i].status()) <= 21 || junEvent[i].colType() == 0)
1177  continue;
1178  vector<int> motherList = junEvent[i].motherList();
1179  int iMot1 = motherList[0];
1180  vector<int> sisterList = junEvent[iMot1].daughterList();
1181 
1182  // Check baryon number of vertex.
1183  int barSum = 0;
1184  map<int,int> colVertex, acolVertex;
1185 
1186  // Loop over mothers (enter with crossed colors and negative sign).
1187  for (unsigned int indx = 0; indx < motherList.size(); indx++) {
1188  int iMot = motherList[indx];
1189  if ( abs(junEvent[iMot].colType()) == 1 )
1190  barSum -= junEvent[iMot].colType();
1191  else if ( abs(junEvent[iMot].colType()) == 3)
1192  barSum -= 2*junEvent[iMot].colType()/3;
1193  int col = junEvent[iMot].acol();
1194  int acol = junEvent[iMot].col();
1195 
1196  // If unmatched (so far), add end. Else erase matching parton.
1197  if (col > 0) {
1198  if (acolVertex.find(col) == acolVertex.end() ) colVertex[col] = iMot;
1199  else acolVertex.erase(col);
1200  } else if (col < 0) {
1201  if (colVertex.find(-col) == colVertex.end() ) acolVertex[-col] = iMot;
1202  else colVertex.erase(-col);
1203  }
1204  if (acol > 0) {
1205  if (colVertex.find(acol) == colVertex.end()) acolVertex[acol] = iMot;
1206  else colVertex.erase(acol);
1207  } else if (acol < 0) {
1208  if (acolVertex.find(-acol) == acolVertex.end())
1209  colVertex[-acol] = iMot;
1210  else acolVertex.erase(-acol);
1211  }
1212  }
1213 
1214  // Loop over sisters.
1215  for (unsigned int indx = 0; indx < sisterList.size(); indx++) {
1216  int iDau = sisterList[indx];
1217  if ( abs(junEvent[iDau].colType()) == 1 )
1218  barSum += junEvent[iDau].colType();
1219  else if ( abs(junEvent[iDau].colType()) == 3)
1220  barSum += 2*junEvent[iDau].colType()/3;
1221  int col = junEvent[iDau].col();
1222  int acol = junEvent[iDau].acol();
1223 
1224  // If unmatched (so far), add end. Else erase matching parton.
1225  if (col > 0) {
1226  if (acolVertex.find(col) == acolVertex.end() ) colVertex[col] = iDau;
1227  else acolVertex.erase(col);
1228  } else if (col < 0) {
1229  if (colVertex.find(-col) == colVertex.end() ) acolVertex[-col] = iDau;
1230  else colVertex.erase(-col);
1231  }
1232  if (acol > 0) {
1233  if (colVertex.find(acol) == colVertex.end()) acolVertex[acol] = iDau;
1234  else colVertex.erase(acol);
1235  } else if (acol < 0) {
1236  if (acolVertex.find(-acol) == acolVertex.end())
1237  colVertex[-acol] = iDau;
1238  else acolVertex.erase(-acol);
1239  }
1240 
1241  }
1242 
1243  // Skip if baryon number conserved in this vertex.
1244  if (barSum == 0) continue;
1245 
1246  // Check and skip any junctions that have already been added.
1247  for (int iJun = 0; iJun < junEvent.sizeJunction(); ++iJun) {
1248  // Remove the tags corresponding to each of the 3 existing junction legs.
1249  for (int j = 0; j < 3; ++j) {
1250  int colNow = junEvent.colJunction(iJun, j);
1251  if (junEvent.kindJunction(iJun) % 2 == 1) colVertex.erase(colNow);
1252  else acolVertex.erase(colNow);
1253  }
1254  }
1255 
1256  // Skip if no junction colors remain.
1257  if (colVertex.size() == 0 && acolVertex.size() == 0) continue;
1258 
1259  // If baryon number violated, is B = +1 or -1 (larger values not handled).
1260  int kindJun = 0;
1261  if (colVertex.size() == 3 && acolVertex.size() == 0) kindJun = 1;
1262  else if (colVertex.size() == 0 && acolVertex.size() == 3) kindJun = 2;
1263  else {
1264  infoPtr->errorMsg("Error in ProcessLevel::findJunctions: "
1265  "N(unmatched (anti)colour tags) != 3");
1266  return;
1267  }
1268 
1269  // From now on, use colJun as shorthand for colVertex or acolVertex.
1270  map<int,int> colJun = (kindJun == 1) ? colVertex : acolVertex;
1271 
1272  // Order so incoming tags appear first in colVec, outgoing tags last.
1273  vector<int> colVec;
1274  for (map<int,int>::iterator it = colJun.begin();
1275  it != colJun.end(); it++) {
1276  int col = it->first;
1277  int iCol = it->second;
1278  // Step across final-state gluons (if they come from ISR => kindJun += 2)
1279  int iColNow = iCol;
1280  int colNow = col;
1281  int nLoop = 0;
1282  while (junEvent[iColNow].isFinal() && junEvent[iColNow].id() == 21) {
1283  colNow = (kindJun % 2 == 1) ? junEvent[iColNow].acol()
1284  : junEvent[iColNow].col();
1285  ++nLoop;
1286  for (int j=0; j<(int)junEvent.size(); ++j) {
1287  // Check for matching initial-state (anti)colour
1288  if ( !junEvent[j].isFinal() ) {
1289  if ( (kindJun%2 == 1 && junEvent[j].acol() == colNow)
1290  || (kindJun%2 == 0 && junEvent[j].col() == colNow) ) {
1291  iColNow = j;
1292  break;
1293  }
1294  }
1295  // Step across final-state gluon
1296  else if ( (kindJun%2 == 1 && junEvent[j].col() == colNow)
1297  || (kindJun%2 == 0 && junEvent[j].acol() == colNow) ) {
1298  iColNow = j;
1299  break;
1300  }
1301  }
1302  // Check for infinite loop
1303  if (nLoop > (int)junEvent.size()) {
1304  infoPtr->errorMsg("Error in ProcessLevel::findJunctions: "
1305  "failed to trace across final-state gluons");
1306  iColNow = iCol;
1307  break;
1308  }
1309  }
1310  for (unsigned int indx = 0; indx < motherList.size(); indx++) {
1311  if (iColNow == motherList[indx]) {
1312  kindJun += 2;
1313  colVec.insert(colVec.begin(),col);
1314  }
1315  }
1316  if (colVec.size() == 0 || colVec[0] != col) colVec.push_back(col);
1317  }
1318 
1319  // Add junction with these tags.
1320  junEvent.appendJunction( kindJun, colVec[0], colVec[1], colVec[2]);
1321 
1322  }
1323 
1324  // Check if any junction colour lines appear both as incoming and outgoing
1325  // E.g. MadGraph writes out 501 + 502 -> -503 -> 501 + 502. Repaint such
1326  // cases so that the outgoing tags are different from the incoming ones.
1327  bool foundMatch = true;
1328  while (foundMatch) {
1329  foundMatch = false;
1330  for (int iJun=0; iJun<junEvent.sizeJunction(); ++iJun) {
1331  int kindJunA = junEvent.kindJunction(iJun);
1332  for (int jJun=iJun+1; jJun<junEvent.sizeJunction(); ++jJun) {
1333  int kindJunB = junEvent.kindJunction(jJun);
1334  // Only consider junction-antijunction combinations
1335  if ( kindJunA % 2 == kindJunB % 2 ) continue;
1336  // Check if all tags same
1337  int nMatch = 0;
1338  for (int iLeg=0; iLeg<3; ++iLeg)
1339  for (int jLeg=0; jLeg<3; ++jLeg)
1340  if (junEvent.colJunction(iJun, iLeg) ==
1341  junEvent.colJunction(jJun, jLeg)) ++nMatch;
1342  if (nMatch == 3) {
1343  foundMatch = true;
1344  // Decide which junction to repaint the final-state legs of
1345  // (If both are types 3-4, arbitrarily decide to repaint iJun)
1346  int kJun = 0;
1347  if (kindJunA >= 5 || kindJunA <= 2) kJun = jJun;
1348  else kJun = iJun;
1349  int kindJun = junEvent.kindJunction(kJun);
1350  int col = junEvent.colJunction(kJun,0);
1351  // Find the corresponding decay vertex: repaint daughters recursively
1352  for (int i=0; i<junEvent.size(); ++i) {
1353  // Find a resonance with the right colour
1354  if ( kindJun % 2 == 0 && junEvent[i].col() != col ) continue;
1355  else if ( kindJun % 2 == 1 && junEvent[i].acol() != col ) continue;
1356  else if ( junEvent[i].status() != -22 ) continue;
1357  // Check if colour is conserved in decay
1358  int iDau1 = junEvent[i].daughter1();
1359  int iDau2 = junEvent[i].daughter2();
1360  bool isBNV = true;
1361  for (int iDau = iDau1; iDau <= iDau2; ++iDau)
1362  if ( (kindJun % 2 == 0 && junEvent[iDau].col() == col)
1363  || (kindJun % 2 == 1 && junEvent[iDau].acol() == col) )
1364  isBNV = false;
1365  if ( !isBNV ) continue;
1366  vector<int> daughters = junEvent[i].daughterListRecursive();
1367  int lastColTag = junEvent.lastColTag();
1368  for (int iLeg = 1; iLeg <= 2; ++iLeg) {
1369  // Encode new colour tag so last digit remains the same
1370  // (That way, new CR type models would still allow reconnection)
1371  int colOld = junEvent.colJunction(kJun, iLeg);
1372  int colNew = (lastColTag/10) * 10 + 10 + colOld % 10 ;
1373  // Count up used colour tags until we reach colNew
1374  while (junEvent.lastColTag() < colNew) junEvent.nextColTag();
1375  junEvent.colJunction(kJun,iLeg,colNew);
1376  for (int jDau = 0; jDau < (int)daughters.size(); ++jDau) {
1377  int iDau = daughters[jDau];
1378  if ( kindJun % 2 == 1 && junEvent[iDau].col() == colOld)
1379  junEvent[iDau].col(colNew);
1380  else if ( kindJun % 2 == 0 && junEvent[iDau].acol() == colOld)
1381  junEvent[iDau].acol(colNew);
1382  }
1383  }
1384  // Done (we found the right BNV vertex and acted recursively)
1385  break;
1386  }
1387  }
1388  }
1389  }
1390  }
1391 }
1392 //--------------------------------------------------------------------------
1393 
1394 // Check that colours match up.
1395 
1396 bool ProcessLevel::checkColours( Event& process) {
1397 
1398  // Variables and arrays for common usage.
1399  bool physical = true;
1400  bool match;
1401  int colType, col, acol, iPos, iNow, iNowA;
1402  vector<int> colTags, colPos, acolPos;
1403 
1404  // Check that each particle has the kind of colours expected of it.
1405  for (int i = 0; i < process.size(); ++i) {
1406  colType = process[i].colType();
1407  col = process[i].col();
1408  acol = process[i].acol();
1409  if (colType == 0 && (col != 0 || acol != 0)) physical = false;
1410  else if (colType == 1 && (col <= 0 || acol != 0)) physical = false;
1411  else if (colType == -1 && (col != 0 || acol <= 0)) physical = false;
1412  else if (colType == 2 && (col <= 0 || acol <= 0)) physical = false;
1413  // Colour-sextet assignments (colType == 3 for sextet, -3 for antisextet).
1414  // Sextet (two colours) represented by (colour, negative anticolour) tags.
1415  // Antisextet (two anticolours) by (negative colour, anticolour) tags.
1416  else if (colType == 3 && (col <= 0 || acol >= 0)) physical = false;
1417  else if (colType == -3 && (col >= 0 || acol <= 0)) physical = false;
1418  // All other cases
1419  else if (colType == -2 || colType < -3 || colType > 3) physical = false;
1420 
1421  // Add to the list of colour tags.
1422  if (col > 0) {
1423  match = false;
1424  for (int ic = 0; ic < int(colTags.size()) ; ++ic)
1425  if (col == colTags[ic]) match = true;
1426  if (!match) colTags.push_back(col);
1427  } else if (acol > 0) {
1428  match = false;
1429  for (int ic = 0; ic < int(colTags.size()) ; ++ic)
1430  if (acol == colTags[ic]) match = true;
1431  if (!match) colTags.push_back(acol);
1432  }
1433  // Colour sextets : map negative colour -> anticolour and vice versa
1434  if (col < 0) {
1435  match = false;
1436  for (int ic = 0; ic < int(colTags.size()) ; ++ic)
1437  if (-col == colTags[ic]) match = true;
1438  if (!match) colTags.push_back(-col);
1439  } else if (acol < 0) {
1440  match = false;
1441  for (int ic = 0; ic < int(colTags.size()) ; ++ic)
1442  if (-acol == colTags[ic]) match = true;
1443  if (!match) colTags.push_back(-acol);
1444  }
1445  }
1446 
1447  // Warn and give up if particles did not have the expected colours.
1448  if (!physical) {
1449  infoPtr->errorMsg("Error in ProcessLevel::checkColours: "
1450  "incorrect colour assignment");
1451  return false;
1452  }
1453 
1454  // Remove (anti)colours coming from an (anti)junction.
1455  for (int iJun = 0; iJun < process.sizeJunction(); ++iJun) {
1456  for (int j = 0; j < 3; ++j) {
1457  int colJun = process.colJunction(iJun, j);
1458  for (int ic = 0; ic < int(colTags.size()) ; ++ic)
1459  if (colJun == colTags[ic]) {
1460  colTags[ic] = colTags[colTags.size() - 1];
1461  colTags.pop_back();
1462  break;
1463  }
1464  }
1465  }
1466 
1467  // Loop through all colour tags and find their positions (by sign).
1468  for (int ic = 0; ic < int(colTags.size()); ++ic) {
1469  col = colTags[ic];
1470  colPos.resize(0);
1471  acolPos.resize(0);
1472  for (int i = 0; i < process.size(); ++i) {
1473  if (process[i].col() == col || process[i].acol() == -col)
1474  colPos.push_back(i);
1475  if (process[i].acol() == col || process[i].col() == -col)
1476  acolPos.push_back(i);
1477  }
1478 
1479  // Trace colours back through decays; remove daughters.
1480  while (colPos.size() > 1) {
1481  iPos = colPos.size() - 1;
1482  iNow = colPos[iPos];
1483  if ( process[iNow].mother1() == colPos[iPos - 1]
1484  && process[iNow].mother2() == 0) colPos.pop_back();
1485  else break;
1486  }
1487  while (acolPos.size() > 1) {
1488  iPos = acolPos.size() - 1;
1489  iNow = acolPos[iPos];
1490  if ( process[iNow].mother1() == acolPos[iPos - 1]
1491  && process[iNow].mother2() == 0) acolPos.pop_back();
1492  else break;
1493  }
1494 
1495  // Now colour should exist in only 2 copies.
1496  if (colPos.size() + acolPos.size() != 2) physical = false;
1497 
1498  // If both colours or both anticolours then one mother of the other.
1499  else if (colPos.size() == 2) {
1500  iNow = colPos[1];
1501  if ( process[iNow].mother1() != colPos[0]
1502  && process[iNow].mother2() != colPos[0] ) physical = false;
1503  }
1504  else if (acolPos.size() == 2) {
1505  iNowA = acolPos[1];
1506  if ( process[iNowA].mother1() != acolPos[0]
1507  && process[iNowA].mother2() != acolPos[0] ) physical = false;
1508  }
1509 
1510  // If one of each then should have same mother(s), or point to beams.
1511  else {
1512  iNow = colPos[0];
1513  iNowA = acolPos[0];
1514  if ( process[iNow].status() == -21 && process[iNowA].status() == -21 );
1515  else if ( (process[iNow].mother1() != process[iNowA].mother1())
1516  || (process[iNow].mother2() != process[iNowA].mother2()) )
1517  physical = false;
1518  }
1519 
1520  }
1521 
1522  // Error message if problem found. Done.
1523  if (!physical) infoPtr->errorMsg("Error in ProcessLevel::checkColours: "
1524  "unphysical colour flow");
1525  return physical;
1526 
1527 }
1528 
1529 //--------------------------------------------------------------------------
1530 
1531 // Print statistics when two hard processes allowed.
1532 
1533 void ProcessLevel::statistics2(bool reset) {
1534 
1535  // Average impact-parameter factor.
1536  double impactFac = max( 1., infoPtr->enhanceMPIavg() );
1537 
1538  // Derive scaling factor to be applied to first set of processes.
1539  double sigma2SelSum = 0.;
1540  int n2SelSum = 0;
1541  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2) {
1542  sigma2SelSum += container2Ptrs[i2]->sigmaSelMC();
1543  n2SelSum += container2Ptrs[i2]->nSelected();
1544  }
1545  double factor1 = impactFac * sigma2SelSum / sigmaND;
1546  double rel1Err = sqrt(1. / max(1, n2SelSum));
1547  if (allHardSame) factor1 *= 0.5;
1548 
1549  // Derive scaling factor to be applied to second set of processes.
1550  double sigma1SelSum = 0.;
1551  int n1SelSum = 0;
1552  for (int i = 0; i < int(containerPtrs.size()); ++i) {
1553  sigma1SelSum += containerPtrs[i]->sigmaSelMC();
1554  n1SelSum += containerPtrs[i]->nSelected();
1555  }
1556  double factor2 = impactFac * sigma1SelSum / sigmaND;
1557  if (allHardSame) factor2 *= 0.5;
1558  double rel2Err = sqrt(1. / max(1, n1SelSum));
1559 
1560  // Header.
1561  cout << "\n *------- PYTHIA Event and Cross Section Statistics ------"
1562  << "--------------------------------------------------*\n"
1563  << " | "
1564  << " |\n"
1565  << " | Subprocess Code | "
1566  << "Number of events | sigma +- delta |\n"
1567  << " | | Tried"
1568  << " Selected Accepted | (estimated) (mb) |\n"
1569  << " | | "
1570  << " | |\n"
1571  << " |------------------------------------------------------------"
1572  << "------------------------------------------------|\n"
1573  << " | | "
1574  << " | |\n"
1575  << " | First hard process: | "
1576  << " | |\n"
1577  << " | | "
1578  << " | |\n";
1579 
1580  // Reset sum counters.
1581  long nTrySum = 0;
1582  long nSelSum = 0;
1583  long nAccSum = 0;
1584  double sigmaSum = 0.;
1585  double delta2Sum = 0.;
1586 
1587  // Reset process maps.
1588  map<int, string> nameM;
1589  map<int, long> nTryM, nSelM, nAccM;
1590  map<int, double> sigmaM, delta2M;
1591 
1592  // Loop over existing first processes.
1593  for (int i = 0; i < int(containerPtrs.size()); ++i)
1594  if (containerPtrs[i]->sigmaMax() != 0.) {
1595 
1596  // Read info for process. Sum counters.
1597  int code = containerPtrs[i]->code();
1598  nTrySum += containerPtrs[i]->nTried();
1599  nSelSum += containerPtrs[i]->nSelected();
1600  nAccSum += containerPtrs[i]->nAccepted();
1601  sigmaSum += containerPtrs[i]->sigmaMC() * factor1;
1602  delta2Sum += pow2(containerPtrs[i]->deltaMC() * factor1);
1603  nameM[code] = containerPtrs[i]->name();
1604  nTryM[code] += containerPtrs[i]->nTried();
1605  nSelM[code] += containerPtrs[i]->nSelected();
1606  nAccM[code] += containerPtrs[i]->nAccepted();
1607  sigmaM[code] += containerPtrs[i]->sigmaMC() * factor1;
1608  delta2M[code] += pow2(containerPtrs[i]->deltaMC() * factor1);
1609  delta2M[code] += pow2(containerPtrs[i]->sigmaMC() * factor1 * rel1Err);
1610  }
1611 
1612  // Print first process info.
1613  for (map<int, string>::iterator i = nameM.begin(); i != nameM.end(); ++i) {
1614  int code = i->first;
1615  cout << " | " << left << setw(40) << i->second
1616  << right << setw(5) << code << " | "
1617  << setw(11) << nTryM[code] << " " << setw(10) << nSelM[code] << " "
1618  << setw(10) << nAccM[code] << " | " << scientific << setprecision(3)
1619  << setw(11) << sigmaM[code]
1620  << setw(11) << sqrtpos(delta2M[code]) << " |\n";
1621  }
1622 
1623  // Print summed info for first processes.
1624  delta2Sum += pow2( sigmaSum * rel1Err );
1625  cout << " | | "
1626  << " | |\n"
1627  << " | " << left << setw(45) << "sum" << right << " | " << setw(11)
1628  << nTrySum << " " << setw(10) << nSelSum << " " << setw(10)
1629  << nAccSum << " | " << scientific << setprecision(3) << setw(11)
1630  << sigmaSum << setw(11) << sqrtpos(delta2Sum) << " |\n";
1631 
1632  // Separation lines to second hard processes.
1633  cout << " | | "
1634  << " | |\n"
1635  << " |------------------------------------------------------------"
1636  << "------------------------------------------------|\n"
1637  << " | | "
1638  << " | |\n"
1639  << " | Second hard process: | "
1640  << " | |\n"
1641  << " | | "
1642  << " | |\n";
1643 
1644  // Reset sum counters.
1645  nTrySum = 0;
1646  nSelSum = 0;
1647  nAccSum = 0;
1648  sigmaSum = 0.;
1649  delta2Sum = 0.;
1650 
1651  // Reset process maps.
1652  nameM.clear();
1653  nTryM.clear();
1654  nSelM.clear();
1655  nAccM.clear();
1656  sigmaM.clear();
1657  delta2M.clear();
1658 
1659  // Loop over existing second processes.
1660  for (int i2 = 0; i2 < int(container2Ptrs.size()); ++i2)
1661  if (container2Ptrs[i2]->sigmaMax() != 0.) {
1662 
1663  // Read info for process. Sum counters.
1664  int code = container2Ptrs[i2]->code();
1665  nTrySum += container2Ptrs[i2]->nTried();
1666  nSelSum += container2Ptrs[i2]->nSelected();
1667  nAccSum += container2Ptrs[i2]->nAccepted();
1668  sigmaSum += container2Ptrs[i2]->sigmaMC() * factor2;
1669  delta2Sum += pow2(container2Ptrs[i2]->deltaMC() * factor2);
1670  nameM[code] = container2Ptrs[i2]->name();
1671  nTryM[code] += container2Ptrs[i2]->nTried();
1672  nSelM[code] += container2Ptrs[i2]->nSelected();
1673  nAccM[code] += container2Ptrs[i2]->nAccepted();
1674  sigmaM[code] += container2Ptrs[i2]->sigmaMC() * factor2;
1675  delta2M[code] += pow2(container2Ptrs[i2]->deltaMC() * factor2);
1676  delta2M[code] += pow2(container2Ptrs[i2]->sigmaMC() * factor2 * rel2Err);
1677  }
1678 
1679  // Print second process info.
1680  for (map<int, string>::iterator i2 = nameM.begin(); i2 != nameM.end();
1681  ++i2) {
1682  int code = i2->first;
1683  cout << " | " << left << setw(40) << i2->second
1684  << right << setw(5) << code << " | "
1685  << setw(11) << nTryM[code] << " " << setw(10) << nSelM[code] << " "
1686  << setw(10) << nAccM[code] << " | " << scientific << setprecision(3)
1687  << setw(11) << sigmaM[code]
1688  << setw(11) << sqrtpos(delta2M[code]) << " |\n";
1689  }
1690 
1691  // Print summed info for second processes.
1692  delta2Sum += pow2( sigmaSum * rel2Err );
1693  cout << " | | "
1694  << " | |\n"
1695  << " | " << left << setw(45) << "sum" << right << " | " << setw(11)
1696  << nTrySum << " " << setw(10) << nSelSum << " " << setw(10)
1697  << nAccSum << " | " << scientific << setprecision(3) << setw(11)
1698  << sigmaSum << setw(11) << sqrtpos(delta2Sum) << " |\n";
1699 
1700  // Print information on how the two processes were combined.
1701  cout << " | | "
1702  << " | |\n"
1703  << " |------------------------------------------------------------"
1704  << "------------------------------------------------|\n"
1705  << " | "
1706  << " |\n"
1707  << " | Uncombined cross sections for the two event sets were "
1708  << setw(10) << sigma1SelSum << " and " << sigma2SelSum << " mb, "
1709  << "respectively, combined |\n"
1710  << " | using a sigma(nonDiffractive) of " << setw(10) << sigmaND
1711  << " mb and an impact-parameter enhancement factor of "
1712  << setw(10) << impactFac << ". |\n";
1713 
1714  // Listing finished.
1715  cout << " | "
1716  << " |\n"
1717  << " *------- End PYTHIA Event and Cross Section Statistics -----"
1718  << "------------------------------------------------*" << endl;
1719 
1720  // Optionally reset statistics contants.
1721  if (reset) resetStatistics();
1722 
1723 }
1724 
1725 //==========================================================================
1726 
1727 } // end namespace Pythia8
Definition: AgUStep.h:26