StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
PartonLevel.cc
1 // PartonLevel.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 // Hard diffraction added by Christine Rasmussen.
6 
7 // Function definitions (not found in the header) for the PartonLevel class.
8 
9 #include "Pythia8/PartonLevel.h"
10 
11 namespace Pythia8 {
12 
13 //==========================================================================
14 
15 // The PartonLevel class.
16 
17 //--------------------------------------------------------------------------
18 
19 // Constants: could be changed here if desired, but normally should not.
20 // These are of technical nature, as described for each.
21 
22 // Maximum number of tries to produce parton level from given input.
23 const int PartonLevel::NTRY = 10;
24 
25 //--------------------------------------------------------------------------
26 
27 // Main routine to initialize the parton-level generation process.
28 
29 bool PartonLevel::init( Info* infoPtrIn, Settings& settings,
30  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn,
31  BeamParticle* beamAPtrIn, BeamParticle* beamBPtrIn,
32  BeamParticle* beamPomAPtrIn, BeamParticle* beamPomBPtrIn,
33  BeamParticle* beamGamAPtrIn, BeamParticle* beamGamBPtrIn,
34  BeamParticle* beamVMDAPtrIn, BeamParticle* beamVMDBPtrIn,
35  Couplings* couplingsPtrIn, PartonSystems* partonSystemsPtrIn,
36  SigmaTotal* sigmaTotPtr, TimeShower* timesDecPtrIn, TimeShower* timesPtrIn,
37  SpaceShower* spacePtrIn, RHadrons* rHadronsPtrIn, UserHooks* userHooksPtrIn,
38  MergingHooks* mergingHooksPtrIn, PartonVertex* partonVertexPtrIn,
39  bool useAsTrial ) {
40 
41  // Store input pointers and modes for future use.
42  infoPtr = infoPtrIn;
43  particleDataPtr = particleDataPtrIn;
44  rndmPtr = rndmPtrIn;
45  beamAPtr = beamAPtrIn;
46  beamBPtr = beamBPtrIn;
47  beamHadAPtr = beamAPtr;
48  beamHadBPtr = beamBPtr;
49  beamPomAPtr = beamPomAPtrIn;
50  beamPomBPtr = beamPomBPtrIn;
51  beamGamAPtr = beamGamAPtrIn;
52  beamGamBPtr = beamGamBPtrIn;
53  beamVMDAPtr = beamVMDAPtrIn;
54  beamVMDBPtr = beamVMDBPtrIn;
55  couplingsPtr = couplingsPtrIn;
56  partonSystemsPtr = partonSystemsPtrIn;
57  timesDecPtr = timesDecPtrIn;
58  timesPtr = timesPtrIn;
59  spacePtr = spacePtrIn;
60  rHadronsPtr = rHadronsPtrIn;
61  userHooksPtr = userHooksPtrIn;
62  mergingHooksPtr = mergingHooksPtrIn;
63  partonVertexPtr = partonVertexPtrIn;
64 
65  // Min bias and diffraction processes need special treatment.
66  bool doSQ = settings.flag("SoftQCD:all")
67  || settings.flag("SoftQCD:inelastic");
68  bool doND = settings.flag("SoftQCD:nonDiffractive");
69  bool doSD = settings.flag("SoftQCD:singleDiffractive");
70  bool doDD = settings.flag("SoftQCD:doubleDiffractive");
71  bool doCD = settings.flag("SoftQCD:centralDiffractive");
72  doNonDiff = doSQ || doND;
73  doDiffraction = doSQ || doSD || doDD || doCD;
74  doHardDiff = settings.flag("Diffraction:doHard");
75  sampleTypeDiff = (doHardDiff) ? settings.mode("Diffraction:sampleType")
76  : 0;
77 
78  // Separate low-mass (unresolved) and high-mass (perturbative) diffraction.
79  mMinDiff = settings.parm("Diffraction:mMinPert");
80  mWidthDiff = settings.parm("Diffraction:mWidthPert");
81  pMaxDiff = settings.parm("Diffraction:probMaxPert");
82  if (mMinDiff > infoPtr->eCM()) doDiffraction = false;
83 
84  // Set whether photon inside lepton. Mode updated event-by-event.
85  gammaMode = settings.mode("Photon:ProcessType");
86  gammaModeEvent = 0;
87  beamHasGamma = settings.flag("PDF:lepton2gamma")
88  && (beamAPtr != 0) && (beamBPtr != 0);
89  hasGammaA = false;
90  hasGammaB = false;
91  beamAisGamma = (beamAPtr != 0) ? beamAPtr->isGamma() : false;
92  beamBisGamma = (beamBPtr != 0) ? beamBPtr->isGamma() : false;
93  beamAhasGamma = (beamHasGamma && beamAPtr->isLepton());
94  beamBhasGamma = (beamHasGamma && beamBPtr->isLepton());
95  beamAhasResGamma = (beamAPtr != 0) ? beamAPtr->hasResGamma() : false;
96  beamBhasResGamma = (beamBPtr != 0) ? beamBPtr->hasResGamma() : false;
97  beamHasResGamma = (gammaMode < 4) && beamHasGamma;
98  isGammaHadronDir = false;
99  bool isAgamBhad = (beamAisGamma || beamAhasGamma) && beamBPtr->isHadron();
100  bool isBgamAhad = (beamBisGamma || beamBhasGamma) && beamAPtr->isHadron();
101  bool isAgamBgam = (beamAisGamma || beamAhasGamma)
102  && (beamBisGamma || beamBhasGamma);
103  bool onlyDirGamma = false;
104  if ( gammaMode == 4 || ( gammaMode == 3 && isAgamBhad )
105  || ( gammaMode == 2 && isBgamAhad ) || ( gammaMode > 1 && isAgamBgam) )
106  onlyDirGamma = true;
107 
108  // Show the copies of beam photon if found in ISR.
109  showUnresGamma = settings.flag("Photon:showUnres");
110 
111  // Need MPI initialization for soft QCD processes, even if only first MPI.
112  // But no need to initialize MPI if never going to use it.
113  doMPI = settings.flag("PartonLevel:MPI");
114  doMPIMB = doMPI;
115  doMPISDA = doMPI;
116  doMPISDB = doMPI;
117  doMPICD = doMPI;
118  doMPIinit = doMPI;
119  doMPIgmgm = doMPI;
120  if (doNonDiff || doDiffraction) doMPIinit = true;
121  if (!settings.flag("PartonLevel:all")) doMPIinit = false;
122 
123  // Nature of MPI matching also used here for one case.
124  pTmaxMatchMPI = settings.mode("MultipartonInteractions:pTmaxMatch");
125 
126  // Initialise trial shower switch.
127  doTrial = useAsTrial;
128  // Merging initialization.
129  bool hasMergingHooks = (mergingHooksPtr != 0);
130  canRemoveEvent = !doTrial && hasMergingHooks
131  && ( mergingHooksPtr->doCKKWLMerging() || mergingHooksPtr->doNL3Merging());
132  canRemoveEmission = !doTrial && hasMergingHooks
133  && ( mergingHooksPtr->doUMEPSMerging() || mergingHooksPtr->doNL3Merging()
134  || mergingHooksPtr->doUNLOPSMerging() );
135  nTrialEmissions = 1;
136  pTLastBranch = 0.0;
137  typeLastBranch = 0;
138 
139  // Flags for showers: ISR and FSR.
140  doISR = settings.flag("PartonLevel:ISR");
141  bool FSR = settings.flag("PartonLevel:FSR");
142  bool FSRinProcess = settings.flag("PartonLevel:FSRinProcess");
143  bool interleaveFSR = settings.flag("TimeShower:interleave");
144  doFSRduringProcess = FSR && FSRinProcess && interleaveFSR;
145  doFSRafterProcess = FSR && FSRinProcess && !interleaveFSR;
146  doFSRinResonances = FSR && settings.flag("PartonLevel:FSRinResonances");
147 
148  // Flags for colour reconnection.
149  doReconnect = settings.flag("ColourReconnection:reconnect");
150  reconnectMode = settings.mode("ColourReconnection:mode");
151  forceResonanceCR = settings.flag("ColourReconnection:forceResonance");
152 
153  // Some other flags.
154  doRemnants = settings.flag("PartonLevel:Remnants");
155  doSecondHard = settings.flag("SecondHard:generate");
156  earlyResDec = settings.flag("PartonLevel:earlyResDec");
157 
158  // Allow R-hadron formation.
159  allowRH = settings.flag("RHadrons:allow");
160 
161  // Possibility to allow user veto during evolution.
162  canVetoPT = (userHooksPtr != 0)
163  ? userHooksPtr->canVetoPT() : false;
164  pTvetoPT = (canVetoPT)
165  ? userHooksPtr->scaleVetoPT() : -1.;
166  canVetoStep = (userHooksPtr != 0)
167  ? userHooksPtr->canVetoStep() : false;
168  nVetoStep = (canVetoStep)
169  ? userHooksPtr->numberVetoStep() : -1;
170  canVetoMPIStep = (userHooksPtr != 0)
171  ? userHooksPtr->canVetoMPIStep() : false;
172  nVetoMPIStep = (canVetoMPIStep)
173  ? userHooksPtr->numberVetoMPIStep() : -1;
174  canVetoEarly = (userHooksPtr != 0)
175  ? userHooksPtr->canVetoPartonLevelEarly() : false;
176 
177  // Settings for vetoing of QCD emission for Drell-Yan weak boson production.
178  vetoWeakJets = settings.flag("WeakShower:vetoQCDjets");
179  vetoWeakDeltaR2 = pow2(settings.parm("WeakShower:vetoWeakDeltaR"));
180 
181  // Possibility to set maximal shower scale in resonance decays.
182  canSetScale = (userHooksPtr != 0)
183  ? userHooksPtr->canSetResonanceScale() : false;
184 
185  // Possibility to reconnect specifically for resonance decays.
186  canReconResSys = (userHooksPtr != 0)
187  ? userHooksPtr->canReconnectResonanceSystems() : false;
188 
189  // Done with initialization only for FSR in resonance decays.
190  if (beamAPtr == 0 || beamBPtr == 0) return true;
191 
192  // Make sure that photons are in resolved mode when mixing with unresolved
193  // before initializing MPIs.
194  if ( (beamAisGamma || beamAhasGamma) && gammaMode == 0) {
195  beamAPtr->setGammaMode(1);
196  if (beamAhasGamma) beamGamAPtr->setGammaMode(1);
197  }
198  if ( (beamBisGamma || beamBhasGamma) && gammaMode == 0) {
199  beamBPtr->setGammaMode(1);
200  if (beamBhasGamma) beamGamBPtr->setGammaMode(1);
201  }
202 
203  // If only direct photons do not initialize diffractive MPI systems.
204  if ( ( (gammaMode == 3) &&
205  ( (beamAPtr->isGamma() || beamAhasGamma) && beamBPtr->isHadron() ) )
206  || ( (gammaMode == 2) &&
207  ( (beamBPtr->isGamma() || beamBhasGamma) && beamAPtr->isHadron() ) )
208  || ( (gammaMode > 1) &&
209  ( (beamAPtr->isGamma() || beamAhasGamma)
210  && (beamBPtr->isGamma() || beamBhasGamma) ) ) )
211  onlyDirGamma = true;
212 
213  // Flag if lepton beams, and if non-resolved ones. May change main flags.
214  hasTwoLeptonBeams = beamAPtr->isLepton() && beamBPtr->isLepton();
215  hasOneLeptonBeam = (beamAPtr->isLepton() || beamBPtr->isLepton())
216  && !hasTwoLeptonBeams;
217  hasPointLeptons = (hasOneLeptonBeam || hasTwoLeptonBeams)
218  && (beamAPtr->isUnresolved() || beamBPtr->isUnresolved());
219  if ( (hasOneLeptonBeam || hasTwoLeptonBeams) && !beamHasResGamma ) {
220  doMPIMB = false;
221  doMPISDA = false;
222  doMPISDB = false;
223  doMPICD = false;
224  doMPIinit = false;
225  doMPIgmgm = false;
226  }
227  if (hasTwoLeptonBeams && hasPointLeptons) {
228  doISR = false;
229  doRemnants = false;
230  }
231 
232  // For ND events in lepton->gamma events no need to initialize MPIs for l+l-.
233  doNDgamma = false;
234  if (beamHasResGamma) doMPIinit = false;
235  if (beamHasResGamma && doND) doNDgamma = true;
236 
237  // Set info and initialize the respective program elements.
238  timesPtr->init( beamAPtr, beamBPtr);
239  if (doISR) spacePtr->init( beamAPtr, beamBPtr);
240 
241  doMPIMB = multiMB.init( doMPIinit, 0, infoPtr, settings, particleDataPtr,
242  rndmPtr, beamAPtr, beamBPtr, couplingsPtr,
243  partonSystemsPtr, sigmaTotPtr, userHooksPtr, partonVertexPtr);
244 
245  // Initialize MPIs for diffractive system, possibly photon beam from
246  // lepton, possibly VMD from photon.
247  if (doSD || doDD || doSQ || (doHardDiff && beamBPtr->getGammaMode() < 2)) {
248  BeamParticle* tmpBeamA = (beamAhasGamma) ? beamGamAPtr : beamAPtr;
249  if (infoPtr->isVMDstateA()) tmpBeamA = beamVMDAPtr;
250  doMPISDA = multiSDA.init( !onlyDirGamma, 1, infoPtr, settings,
251  particleDataPtr, rndmPtr, tmpBeamA, beamPomBPtr, couplingsPtr,
252  partonSystemsPtr, sigmaTotPtr, userHooksPtr, partonVertexPtr,
253  (beamAisGamma || beamAhasGamma) );
254  }
255  if (doSD || doDD || doSQ || (doHardDiff && beamAPtr->getGammaMode() < 2)) {
256  BeamParticle* tmpBeamB = (beamBhasGamma) ? beamGamBPtr : beamBPtr;
257  if (infoPtr->isVMDstateB()) tmpBeamB = beamVMDBPtr;
258  doMPISDB = multiSDB.init( !onlyDirGamma, 2, infoPtr, settings,
259  particleDataPtr, rndmPtr, beamPomAPtr, tmpBeamB, couplingsPtr,
260  partonSystemsPtr, sigmaTotPtr, userHooksPtr, partonVertexPtr,
261  (beamBisGamma || beamBhasGamma) );
262  }
263  if (doCD || doSQ) doMPICD = multiCD.init( doMPIinit, 3, infoPtr, settings,
264  particleDataPtr, rndmPtr, beamPomAPtr, beamPomBPtr, couplingsPtr,
265  partonSystemsPtr, sigmaTotPtr, userHooksPtr, partonVertexPtr);
266  if (!remnants.init( infoPtr, settings, rndmPtr, beamAPtr, beamBPtr,
267  partonSystemsPtr, partonVertexPtr, particleDataPtr, &colourReconnection))
268  return false;
269  resonanceDecays.init( infoPtr, particleDataPtr, rndmPtr);
270  colourReconnection.init( infoPtr, settings, rndmPtr, particleDataPtr,
271  beamAPtr, beamBPtr, partonSystemsPtr);
272  junctionSplitting.init(infoPtr, settings, rndmPtr, particleDataPtr);
273 
274  // Initialize hard diffraction with possible photon beams from leptons.
275  if (doHardDiff && (gammaMode !=4) )
276  hardDiffraction.init(infoPtr, settings, rndmPtr, beamAhasGamma ?
277  beamGamAPtr : beamAPtr, beamBhasGamma ? beamGamBPtr : beamBPtr,
278  beamPomAPtr, beamPomBPtr, sigmaTotPtr);
279 
280  // Initialize an MPI instance for photons from leptons.
281  if ( beamHasResGamma && (doMPI || doNDgamma) ) {
282  doMPIinit = true;
283  // Lepton-hadron.
284  if (beamAPtr->isLepton() && beamBPtr->isHadron() ) {
285  doMPIgmgm = multiGmGm.init( doMPIinit, 0, infoPtr, settings,
286  particleDataPtr, rndmPtr, beamGamAPtr, beamBPtr, couplingsPtr,
287  partonSystemsPtr, sigmaTotPtr, userHooksPtr, partonVertexPtr, true);
288  // Hadron-lepton.
289  } else if (beamBPtr->isLepton() && beamAPtr->isHadron() ) {
290  doMPIgmgm = multiGmGm.init( doMPIinit, 0, infoPtr, settings,
291  particleDataPtr, rndmPtr, beamAPtr, beamGamBPtr, couplingsPtr,
292  partonSystemsPtr, sigmaTotPtr, userHooksPtr, partonVertexPtr, true);
293  // Lepton-lepton.
294  } else {
295  doMPIgmgm = multiGmGm.init( doMPIinit, 0, infoPtr, settings,
296  particleDataPtr, rndmPtr, beamGamAPtr, beamGamBPtr, couplingsPtr,
297  partonSystemsPtr, sigmaTotPtr, userHooksPtr, partonVertexPtr, true);
298  }
299  doMPIMB = doMPIgmgm;
300  }
301 
302  // Succeeded, or not.
303  multiPtr = &multiMB;
304  if (doMPIinit && !doMPIMB) return false;
305  if (doMPIinit && (doSD || doDD || doSQ) && (!doMPISDA || !doMPISDB))
306  return false;
307  if (doMPIinit && (doCD || doSQ) && !doMPICD) return false;
308  if (!doMPIMB || !doMPISDA || !doMPISDB || !doMPICD) doMPI = false;
309  return true;
310 
311 }
312 
313 //--------------------------------------------------------------------------
314 
315 // Function to reset PartonLevel object for trial shower usage.
316 
317 void PartonLevel::resetTrial() {
318 
319  // Clear input pointers.
320  partonSystemsPtr->clear();
321  beamAPtr->clear();
322  beamBPtr->clear();
323  beamHadAPtr->clear();
324  beamHadBPtr->clear();
325  beamPomAPtr->clear();
326  beamPomBPtr->clear();
327  beamGamAPtr->clear();
328  beamGamBPtr->clear();
329  beamVMDAPtr->clear();
330  beamVMDBPtr->clear();
331 
332  // Clear last branching return values.
333  pTLastBranch = 0.0;
334  typeLastBranch = 0;
335 
336 }
337 
338 //--------------------------------------------------------------------------
339 
340 // Main routine to do the parton-level evolution.
341 
342 bool PartonLevel::next( Event& process, Event& event) {
343 
344  // Current event classification.
345  isResolved = infoPtr->isResolved();
346  isResolvedA = isResolved;
347  isResolvedB = isResolved;
348  isResolvedC = isResolved;
349  isDiffA = infoPtr->isDiffractiveA();
350  isDiffB = infoPtr->isDiffractiveB();
351  isDiffC = infoPtr->isDiffractiveC();
352  isDiff = isDiffA || isDiffB || isDiffC;
353  isCentralDiff = isDiffC;
354  isDoubleDiff = isDiffA && isDiffB;
355  isSingleDiff = isDiff && !isDoubleDiff && !isCentralDiff;
356  isNonDiff = infoPtr->isNonDiffractive();
357 
358  // Default values for what is to come with diffraction.
359  isHardDiffA = false;
360  isHardDiffB = false;
361  isHardDiff = false;
362  doDiffVeto = false;
363  // Mark hard diffractive events to handle CR correctly.
364  bool doDiffCR = false;
365  // The setup of the diffractive events can come after the first evolution.
366  int nHardDiffLoop = 1;
367  // Flag to check whether hard diffraction system is set up for the process.
368  hardDiffSet = false;
369  // Offset for the initiator position when photons from leptons.
370  gammaOffset = 0;
371 
372  // Parton-level vetoes for matching and merging.
373  doVeto = false;
374  infoPtr->setAbortPartonLevel(false);
375 
376  // Update photon state according to beams set in processContainer.
377  beamAhasResGamma = beamAPtr->hasResGamma();
378  beamBhasResGamma = beamBPtr->hasResGamma();
379  beamHasResGamma = beamAhasResGamma || beamBhasResGamma;
380 
381  // Save if photoproduction from either side.
382  hasGammaA = (beamAPtr != 0) ? beamAPtr->getGammaMode() > 0 : false;
383  hasGammaB = (beamBPtr != 0) ? beamBPtr->getGammaMode() > 0 : false;
384 
385  // Save current photon mode when mixing processes.
386  gammaModeEvent = gammaMode;
387  if ( hasGammaA || hasGammaB ) {
388  if (beamAPtr->getGammaMode() < 2 && beamBPtr->getGammaMode() < 2)
389  gammaModeEvent = 1;
390  if (beamAPtr->getGammaMode() < 2 && beamBPtr->getGammaMode() == 2)
391  gammaModeEvent = 2;
392  if (beamAPtr->getGammaMode() == 2 && beamBPtr->getGammaMode() < 2)
393  gammaModeEvent = 3;
394  if (beamAPtr->getGammaMode() == 2 && beamBPtr->getGammaMode() == 2)
395  gammaModeEvent = 4;
396  }
397 
398  // Check if direct-photon + hadron to set hard system correctly.
399  isGammaHadronDir = !(hasGammaA || hasGammaB) ? false :
400  ( (beamAPtr->getGammaMode() == 2) && (beamBPtr->getGammaMode() == 0) )
401  || ( (beamAPtr->getGammaMode() == 0) && (beamBPtr->getGammaMode() == 2) );
402 
403  // Set up gamma+gamma/p subcollision. May fail due to extreme kinematics.
404  if (beamHasGamma) {
405  if ( !setupResolvedLeptonGamma( process) ) return false;
406  }
407 
408  // Prepare for a potential hard diffractive event.
409  if (doHardDiff) {
410 
411  // Preliminary decision based on diffractive-to-inclusive PDF ratio.
412  // If Pomeron taken from side A(=1), then B is the diffractive system.
413  // If Pomeron taken from side B(=2), then A is the diffractive system.
414  // Do not check for diffraction if direct photons.
415  // Rescale the x values when photon from lepton and calculate PDFs
416  // wrt. photon beam.
417  if ( beamBPtr->isHadron() || ((beamBisGamma || beamBhasGamma)
418  && beamBPtr->getGammaMode() == 1)) {
419  double xGammaB = beamBhasGamma ? beamHadBPtr->xGamma() : 1.;
420  double xBrescaled = infoPtr->x2pdf() / xGammaB;
421  double pdfB = beamBhasGamma ? beamGamBPtr->xf(infoPtr->id2pdf(),
422  xBrescaled, infoPtr->Q2Fac()) : infoPtr->pdf2();
423  isHardDiffA = hardDiffraction.isDiffractive(2, infoPtr->id2pdf(),
424  xBrescaled, infoPtr->Q2Fac(), pdfB);
425  }
426  if ( beamAPtr->isHadron() || ((beamAisGamma || beamAhasGamma)
427  && beamAPtr->getGammaMode() == 1)) {
428  double xGammaA = beamAhasGamma ? beamHadAPtr->xGamma() : 1.;
429  double xArescaled = infoPtr->x1pdf() / xGammaA;
430  double pdfA = beamAhasGamma ? beamGamAPtr->xf(infoPtr->id1pdf(),
431  xArescaled, infoPtr->Q2Fac()) : infoPtr->pdf1();
432  isHardDiffB = hardDiffraction.isDiffractive(1, infoPtr->id1pdf(),
433  xArescaled, infoPtr->Q2Fac(), pdfA);
434  }
435 
436  // No hard double diffraction yet, so randomly choose one of the sides.
437  if (isHardDiffA && isHardDiffB) {
438  if (rndmPtr->flat() < 0.5) isHardDiffA = false;
439  else isHardDiffB = false;
440  }
441  isHardDiff = isHardDiffA || isHardDiffB;
442 
443  // Save diffractive values.
444  double xPomA = (isHardDiffB) ? hardDiffraction.getXPomeronA() : 0.;
445  double xPomB = (isHardDiffA) ? hardDiffraction.getXPomeronB() : 0.;
446  double tPomA = (isHardDiffB) ? hardDiffraction.getTPomeronA() : 0.;
447  double tPomB = (isHardDiffA) ? hardDiffraction.getTPomeronB() : 0.;
448  infoPtr->setHardDiff( false, false, isHardDiffA, isHardDiffB,
449  xPomA, xPomB, tPomA, tPomB);
450 
451  // Discard all nondiffractive events if only diffractive sample is wanted.
452  if (!isHardDiff && sampleTypeDiff > 2) {
453  doDiffVeto = true;
454  // Reset beam pointers and set event back to original frame.
455  if (beamHasGamma) leaveResolvedLeptonGamma( process, event, false);
456  return false;
457  }
458 
459  if (isHardDiff) {
460  // Set up the diffractive system if run without MPI veto.
461  if (sampleTypeDiff%2 == 1) setupHardDiff( process);
462  // Allow for second loop if run with MPI veto.
463  else nHardDiffLoop = 2;
464  }
465  }
466 
467  // nHardLoop counts how many hard-scattering subsystems are to be processed.
468  // Almost always 1, but elastic and low-mass diffraction gives 0, while
469  // double diffraction can give up to 2. Not to be confused with SecondHard.
470  int nHardLoop = 1;
471  if (!isResolved) nHardLoop = (isDiff) ? decideResolvedDiff( process) : 0;
472 
473  // Handle unresolved subsystems. Done if no resolved ones.
474  sizeProcess = 0;
475  sizeEvent = 0;
476  if (!isResolvedA || !isResolvedB || !isResolvedC) {
477  bool physical = setupUnresolvedSys( process, event);
478  if (!physical || nHardLoop == 0) return physical;
479  sizeProcess = process.size();
480  sizeEvent = event.size();
481  }
482 
483  // Number of actual branchings.
484  int nBranch = 0;
485  // Number of desired branchings, negative value means no restriction.
486  int nBranchMax = (doTrial) ? nTrialEmissions : -1;
487 
488  // Store merging weight.
489  bool hasMergingHooks = (mergingHooksPtr != 0);
490  if ( hasMergingHooks && canRemoveEvent )
491  mergingHooksPtr->storeWeights(infoPtr->getWeightCKKWL());
492 
493  // Reset event weight coming from enhanced branchings.
494  if (userHooksPtr != 0) userHooksPtr->setEnhancedEventWeight(1.);
495 
496  // Loop to set up diffractive system if run with MPI veto.
497  for (int iHardDiffLoop = 1; iHardDiffLoop <= nHardDiffLoop;
498  ++iHardDiffLoop) {
499 
500  // Big outer loop to handle up to two systems (in double diffraction),
501  // but normally one. (Not indented in following, but end clearly marked.)
502  for (int iHardLoop = 1; iHardLoop <= nHardLoop; ++iHardLoop) {
503  infoPtr->setCounter(20, iHardLoop);
504  infoPtr->setCounter(21);
505 
506  // Classification of diffractive system: 1 = A, 2 = B, 3 = central.
507  iDS = 0;
508  if (isDiffA || isDiffB) iDS = (iHardLoop == 2 || !isResolvedA) ? 2 : 1;
509  if (isDiffC) iDS = 3;
510 
511  // Process and event records can be out of step for diffraction.
512  if (iHardLoop == 2) {
513  sizeProcess = process.size();
514  sizeEvent = event.size();
515  partonSystemsPtr->clear();
516  if (event.lastColTag() > process.lastColTag())
517  process.initColTag(event.lastColTag());
518  }
519 
520  // If you need to restore then do not throw existing diffractive system.
521  if (isDiff) {
522  event.saveSize();
523  event.saveJunctionSize();
524 
525  // Allow special treatment of diffractive systems.
526  setupResolvedDiff( process);
527  }
528 
529  // Prepare to do multiparton interactions; at new mass for diffraction.
530  if (doMPIinit) multiPtr->reset();
531 
532  // Special case if nondiffractive: do hardest interaction.
533  if (isNonDiff || isDiff) {
534  multiPtr->pTfirst();
535  multiPtr->setupFirstSys( process);
536  }
537 
538  // Allow up to ten tries; failure possible for beam remnants.
539  // Main cause: inconsistent colour flow at the end of the day.
540  bool physical = true;
541  int nRad = 0;
542  for (int iTry = 0; iTry < NTRY; ++ iTry) {
543  infoPtr->addCounter(21);
544  for (int i = 22; i < 32; ++i) infoPtr->setCounter(i);
545 
546  // Reset flag, counters and max scales.
547  physical = true;
548  nMPI = (doSecondHard) ? 2 : 1;
549  nISR = 0;
550  nFSRinProc = 0;
551  nFSRinRes = 0;
552  nISRhard = 0;
553  nFSRhard = 0;
554  pTsaveMPI = 0.;
555  pTsaveISR = 0.;
556  pTsaveFSR = 0.;
557 
558  // Reset nMPI and nISR for showers.
559  infoPtr->setPartEvolved(nMPI, nISR);
560 
561  // Reset parameters related to valence content and remnants of photon
562  // beam if ISR or MPI generated.
563  if (beamAPtr->isGamma() && ( doMPI || doISR ) ) beamAPtr->resetGamma();
564  if (beamBPtr->isGamma() && ( doMPI || doISR ) ) beamBPtr->resetGamma();
565 
566  // Identify hard interaction system for showers.
567  setupHardSys( process, event);
568 
569  // Optionally check for a veto after the hardest interaction.
570  if (canVetoMPIStep) {
571  doVeto = userHooksPtr->doVetoMPIStep( 1, event);
572  // Abort event if vetoed.
573  if (doVeto) {
574  if (isDiff) leaveResolvedDiff( iHardLoop, process, event);
575  if (beamHasResGamma) leaveResolvedLeptonGamma( process, event, false);
576  return false;
577  }
578  }
579 
580  // Check matching of process scale to maximum ISR/FSR/MPI scales.
581  double Q2Fac = infoPtr->Q2Fac();
582  double Q2Ren = infoPtr->Q2Ren();
583  bool limitPTmaxISR = (doISR)
584  ? spacePtr->limitPTmax( event, Q2Fac, Q2Ren) : false;
585  bool limitPTmaxFSR = (doFSRduringProcess)
586  ? timesPtr->limitPTmax( event, Q2Fac, Q2Ren) : false;
587  bool limitPTmaxMPI = (doMPI) ? multiPtr->limitPTmax( event) : false;
588 
589  // Global recoil: reset counters and store locations of outgoing partons.
590  timesPtr->prepareGlobal( event);
591  bool isFirstTrial = true;
592 
593  // Set hard scale, maximum for showers and multiparton interactions.
594  double pTscaleRad = process.scale();
595  double pTscaleMPI = (doMPI && pTmaxMatchMPI == 3)
596  ? multiPtr->scaleLimitPT() : pTscaleRad;
597  if (doSecondHard) {
598  pTscaleRad = max( pTscaleRad, process.scaleSecond() );
599  pTscaleMPI = min( pTscaleMPI, process.scaleSecond() );
600  }
601  double pTmaxMPI = (limitPTmaxMPI) ? pTscaleMPI : infoPtr->eCM();
602  double pTmaxISR = (limitPTmaxISR) ? spacePtr->enhancePTmax() * pTscaleRad
603  : infoPtr->eCM();
604  double pTmaxFSR = (limitPTmaxFSR) ? timesPtr->enhancePTmax() * pTscaleRad
605  : infoPtr->eCM();
606 
607  // Store the starting scale to use it for valence selection for gamma
608  // beam. In case of MPIs use the pT scale of the last one.
609  beamAPtr->pTMPI( process.scale() );
610  beamBPtr->pTMPI( process.scale() );
611 
612  // Potentially reset starting scales for matrix element merging.
613  if ( hasMergingHooks && (doTrial || canRemoveEvent || canRemoveEmission) )
614  mergingHooksPtr->setShowerStartingScales( doTrial,
615  (canRemoveEvent || canRemoveEmission), pTscaleRad, process, pTmaxFSR,
616  limitPTmaxFSR, pTmaxISR, limitPTmaxISR, pTmaxMPI, limitPTmaxMPI );
617  double pTmax = max( pTmaxMPI, max( pTmaxISR, pTmaxFSR) );
618  pTsaveMPI = pTmaxMPI;
619  pTsaveISR = pTmaxISR;
620  pTsaveFSR = pTmaxFSR;
621 
622  // Prepare the classes to begin the generation.
623  if (doMPI) multiPtr->prepare( event, pTmaxMPI, (iHardDiffLoop == 2) );
624  if (doISR) spacePtr->prepare( 0, event, limitPTmaxISR);
625  if (doFSRduringProcess) timesPtr->prepare( 0, event, limitPTmaxFSR);
626  if (doSecondHard && doISR) spacePtr->prepare( 1, event, limitPTmaxISR);
627  if (doSecondHard && doFSRduringProcess) timesPtr->prepare( 1, event,
628  limitPTmaxFSR);
629 
630  // Impact parameter has now been chosen, but not usefully for diffraction.
631  // (Never chosen for low-mass diffraction, twice for double diffraction.)
632  if (!isDiff) infoPtr->setImpact( multiPtr->bMPI(), multiPtr->enhanceMPI(),
633  multiPtr->enhanceMPIavg(), true, (iHardDiffLoop == 2) );
634 
635  // Set up initial veto scale.
636  doVeto = false;
637  double pTveto = pTvetoPT;
638  typeLatest = 0;
639 
640  // Begin evolution down in pT from hard pT scale.
641  do {
642  infoPtr->addCounter(22);
643  typeVetoStep = 0;
644  nRad = nISR + nFSRinProc;
645 
646  // Check whether the beam photon has unresolved during the evolution.
647  bool unresolvedGammaA = (beamAPtr->isGamma()
648  && !(beamAPtr->resolvedGamma()) );
649  bool unresolvedGammaB = (beamBPtr->isGamma()
650  && !(beamBPtr->resolvedGamma()) );
651  bool unresolvedGamma = (unresolvedGammaA || unresolvedGammaB)
652  || gammaModeEvent == 4;
653 
654  // Find next pT value for FSR, MPI and ISR.
655  // Order calls to minimize time expenditure.
656  double pTgen = 0.;
657 
658  // Potentially increase shower stopping scale for trial showers, to
659  // avoid accumulating low-pT emissions (and weights thereof)
660  if ( hasMergingHooks && doTrial)
661  pTgen = max( pTgen, mergingHooksPtr->getShowerStoppingScale() );
662 
663  double pTtimes = (doFSRduringProcess)
664  ? timesPtr->pTnext( event, pTmaxFSR, pTgen, isFirstTrial, doTrial)
665  : -1.;
666  pTgen = max( pTgen, pTtimes);
667  // No MPIs for unresolved photons.
668  double pTmulti = (doMPI && !unresolvedGamma)
669  ? multiPtr->pTnext( pTmaxMPI, pTgen, event) : -1.;
670  pTgen = max( pTgen, pTmulti);
671  double pTspace = (doISR)
672  ? spacePtr->pTnext( event, pTmaxISR, pTgen, nRad, doTrial) : -1.;
673  double pTnow = max( pTtimes, max( pTmulti, pTspace));
674 
675  // Update information.
676  infoPtr->setPTnow( pTnow);
677  isFirstTrial = false;
678 
679  // Allow a user veto. Only do it once, so remember to change pTveto.
680  if (pTveto > 0. && pTveto > pTnow) {
681  pTveto = -1.;
682  doVeto = userHooksPtr->doVetoPT( typeLatest, event);
683  // Abort event if vetoed.
684  if (doVeto) {
685  if (isDiff) leaveResolvedDiff( iHardLoop, process, event);
686  if (beamHasResGamma)
687  leaveResolvedLeptonGamma( process, event, false);
688  return false;
689  }
690  }
691 
692  // Do a multiparton interaction (if allowed).
693  if (pTmulti > 0. && pTmulti > pTspace && pTmulti > pTtimes) {
694  infoPtr->addCounter(23);
695  if (multiPtr->scatter( event)) {
696  typeLatest = 1;
697  ++nMPI;
698  if (canVetoMPIStep && nMPI <= nVetoMPIStep) typeVetoStep = 1;
699 
700  // Break for hard diffraction with MPI veto.
701  if (isHardDiff && sampleTypeDiff == 4 && iHardDiffLoop == 1) {
702  infoPtr->setHardDiff( false, false, false, false, 0., 0., 0., 0.);
703  doDiffVeto = true;
704  // Reset beam pointers and set event back to original frame.
705  if (beamHasResGamma)
706  leaveResolvedLeptonGamma( process, event, false);
707  return false;
708  }
709 
710  // Update ISR and FSR dipoles.
711  if (doISR) spacePtr->prepare( nMPI - 1, event);
712  if (doFSRduringProcess) timesPtr->prepare( nMPI - 1, event);
713  nBranch++;
714  pTLastBranch = pTmulti;
715  typeLastBranch = 1;
716  }
717 
718  // Set maximal scales for next pT to pick.
719  pTmaxMPI = pTmulti;
720  pTmaxISR = min(pTmulti, pTmaxISR);
721  pTmaxFSR = min(pTmulti, pTmaxFSR);
722  pTmax = pTmulti;
723  }
724 
725  // Do an initial-state emission (if allowed).
726  else if (pTspace > 0. && pTspace > pTtimes) {
727  infoPtr->addCounter(24);
728 
729  // If MPIs, construct the gamma->qqbar branching in beamRemnants.
730  if (spacePtr->branch( event)
731  && ( !(nMPI > 1 && spacePtr->wasGamma2qqbar()) ) ) {
732  typeLatest = 2;
733  iSysNow = spacePtr->system();
734  ++nISR;
735  if (iSysNow == 0) ++nISRhard;
736  if (canVetoStep && iSysNow == 0 && nISRhard <= nVetoStep)
737  typeVetoStep = 2;
738 
739  // Update FSR dipoles.
740  if (doFSRduringProcess) timesPtr->update( iSysNow, event,
741  spacePtr->getHasWeaklyRadiated());
742  nBranch++;
743  pTLastBranch = pTspace;
744  typeLastBranch = 2;
745 
746  // Rescatter: it is possible for kinematics to fail, in which
747  // case we need to restart the parton level processing.
748  } else if (spacePtr->doRestart()) {
749  physical = false;
750  break;
751  }
752 
753  // Set maximal scales for next pT to pick.
754  pTmaxMPI = min( min(pTspace,pTmaxISR), pTmaxMPI);
755  pTmaxISR = min(pTspace,pTmaxISR);
756  pTmaxFSR = min( min(pTspace,pTmaxISR), pTmaxFSR);
757  pTmax = pTspace;
758  }
759 
760  // Do a final-state emission (if allowed).
761  else if (pTtimes > 0.) {
762  infoPtr->addCounter(25);
763  if (timesPtr->branch( event, true)) {
764  typeLatest = 3;
765  iSysNow = timesPtr->system();
766  ++nFSRinProc;
767  if (iSysNow == 0) ++nFSRhard;
768  if (canVetoStep && iSysNow == 0 && nFSRhard <= nVetoStep)
769  typeVetoStep = 3;
770 
771  // Update ISR dipoles.
772  if (doISR) spacePtr->update( iSysNow, event,
773  timesPtr->getHasWeaklyRadiated());
774  nBranch++;
775  pTLastBranch = pTtimes;
776  typeLastBranch = 3;
777 
778  }
779 
780  // Set maximal scales for next pT to pick.
781  pTmaxMPI = min( min(pTtimes,pTmaxFSR), pTmaxMPI);
782  pTmaxISR = min( min(pTtimes,pTmaxFSR), pTmaxISR);
783  pTmaxFSR = min(pTtimes, pTmaxFSR);
784  pTmax = pTtimes;
785  }
786 
787  // If no pT scales above zero then nothing to be done.
788  else pTmax = 0.;
789 
790  // Check for double counting for Drell-Yan weak production.
791  // Only look at the second emission.
792  if ( (infoPtr->code() == 221 || infoPtr->code() == 222) &&
793  nISRhard + nFSRhard == 2 && vetoWeakJets) {
794  int id1 = event[partonSystemsPtr->getOut(0,0)].id();
795  int id2 = event[partonSystemsPtr->getOut(0,1)].id();
796  int id3 = event[partonSystemsPtr->getOut(0,2)].id();
797  Vec4 p1 = event[partonSystemsPtr->getOut(0,0)].p();
798  Vec4 p2 = event[partonSystemsPtr->getOut(0,1)].p();
799  Vec4 p3 = event[partonSystemsPtr->getOut(0,2)].p();
800 
801  // Make sure id1 is weak boson, and check that there
802  // only is a single weak boson and no photons.
803  bool doubleCountEvent = true;
804  if (abs(id1) == 24 || abs(id1) == 23) {
805  if (abs(id2) > 21 || abs(id3) > 21)
806  doubleCountEvent = false;
807  } else if (abs(id2) == 24 || abs(id2) == 23) {
808  swap(id1,id2);
809  swap(p1,p2);
810  if (abs(id3) > 21)
811  doubleCountEvent = false;
812  } else if ( abs(id3) == 24 || abs(id3) == 23) {
813  swap(id1,id3);
814  swap(p1,p3);
815  }
816 
817  if (doubleCountEvent) {
818  double d = p1.pT2();
819  bool cut = true;
820  if (p2.pT2() < d) {d = p2.pT2(); cut = false;}
821  if (p3.pT2() < d) {d = p3.pT2(); cut = false;}
822 
823  // Check for angle between weak boson and quarks.
824  // (require final state particle to be a fermion)
825  if (abs(id2) < 20) {
826  double dij = min(p1.pT2(),p2.pT2())
827  * pow2(RRapPhi(p1,p2)) / vetoWeakDeltaR2;
828  if (dij < d) {
829  d = dij;
830  cut = true;
831  }
832  }
833 
834  if (abs(id3) < 20) {
835  double dij = min(p1.pT2(),p3.pT2())
836  * pow2(RRapPhi(p1,p3)) / vetoWeakDeltaR2;
837  if (dij < d) {
838  d = dij;
839  cut = true;
840  }
841  }
842 
843  // Check for angle between recoiler and radiator,
844  // if it is a quark anti-quark pair
845  // or if the recoiler is a gluon.
846  if (abs(id2) == 21 || abs(id3) == 21 || id2 == - id3) {
847  double dij = min(p2.pT2(),p3.pT2())
848  * pow2(RRapPhi(p2,p3)) / vetoWeakDeltaR2;
849  if (dij < d) {
850  d = dij;
851  cut = false;
852  }
853  }
854 
855  // Veto event if it does not belong to Drell-Yan production.
856  if (cut) return false;
857  }
858  }
859 
860  // Optionally check for a veto after the first few interactions,
861  // or after the first few emissions, ISR or FSR, in the hardest system.
862  if (typeVetoStep == 1) {
863  doVeto = userHooksPtr->doVetoMPIStep( nMPI, event);
864  } else if (typeVetoStep > 1) {
865  doVeto = userHooksPtr->doVetoStep( typeVetoStep, nISRhard,
866  nFSRhard, event);
867  }
868 
869  // Abort event if vetoed.
870  if (doVeto) {
871  if (isDiff) leaveResolvedDiff( iHardLoop, process, event);
872  if (beamHasResGamma) leaveResolvedLeptonGamma( process, event);
873  return false;
874  }
875 
876  // Keep on evolving until nothing is left to be done.
877  if (typeLatest > 0 && typeLatest < 4)
878  infoPtr->addCounter(25 + typeLatest);
879  if (!isDiff) infoPtr->setPartEvolved( nMPI, nISR);
880 
881  // Handle potential merging veto.
882  if ( canRemoveEvent && nISRhard + nFSRhard == 1 ) {
883  // Simply check, and possibly reset weights.
884  mergingHooksPtr->doVetoStep( process, event );
885  }
886 
887  // End loop evolution down in pT from hard pT scale.
888  } while (pTmax > 0. && (nBranchMax <= 0 || nBranch < nBranchMax) );
889 
890  // Do all final-state emissions if not already considered above.
891  if (doFSRafterProcess && (nBranchMax <= 0 || nBranch < nBranchMax) ) {
892 
893  // Find largest scale for final partons.
894  pTmax = 0.;
895  for (int i = 0; i < event.size(); ++i)
896  if (event[i].isFinal() && event[i].scale() > pTmax)
897  pTmax = event[i].scale();
898  pTsaveFSR = pTmax;
899 
900  // Prepare all subsystems for evolution.
901  for (int iSys = 0; iSys < partonSystemsPtr->sizeSys(); ++iSys)
902  timesPtr->prepare( iSys, event);
903 
904  // Set up initial veto scale.
905  doVeto = false;
906  pTveto = pTvetoPT;
907 
908  // Begin evolution down in pT from hard pT scale.
909  do {
910  infoPtr->addCounter(29);
911  typeVetoStep = 0;
912  double pTtimes = timesPtr->pTnext( event, pTmax, 0.);
913  infoPtr->setPTnow( pTtimes);
914 
915  // Allow a user veto. Only do it once, so remember to change pTveto.
916  if (pTveto > 0. && pTveto > pTtimes) {
917  pTveto = -1.;
918  doVeto = userHooksPtr->doVetoPT( 4, event);
919  // Abort event if vetoed.
920  if (doVeto) {
921  if (isDiff) leaveResolvedDiff( iHardLoop, process, event);
922  if (beamHasResGamma) leaveResolvedLeptonGamma( process, event);
923  return false;
924  }
925  }
926 
927  // Do a final-state emission (if allowed).
928  if (pTtimes > 0.) {
929  infoPtr->addCounter(30);
930  if (timesPtr->branch( event, true)) {
931  iSysNow = timesPtr->system();
932  ++nFSRinProc;
933  if (iSysNow == 0) ++nFSRhard;
934  if (canVetoStep && iSysNow == 0 && nFSRhard <= nVetoStep)
935  typeVetoStep = 4;
936 
937  nBranch++;
938  pTLastBranch = pTtimes;
939  typeLastBranch = 4;
940 
941  }
942  pTmax = pTtimes;
943  }
944 
945  // If no pT scales above zero then nothing to be done.
946  else pTmax = 0.;
947 
948  // Optionally check for a veto after the first few emissions.
949  if (typeVetoStep > 0) {
950  doVeto = userHooksPtr->doVetoStep( typeVetoStep, nISRhard,
951  nFSRhard, event);
952  // Abort event if vetoed.
953  if (doVeto) {
954  if (isDiff) leaveResolvedDiff( iHardLoop, process, event);
955  if (beamHasResGamma) leaveResolvedLeptonGamma( process, event);
956  return false;
957  }
958  }
959 
960  // Handle potential merging veto.
961  if ( canRemoveEvent && nISRhard + nFSRhard == 1 ) {
962  // Simply check, and possibly reset weights.
963  mergingHooksPtr->doVetoStep( process, event );
964  }
965 
966  // Keep on evolving until nothing is left to be done.
967  infoPtr->addCounter(31);
968 
969  } while (pTmax > 0. && (nBranchMax <= 0 || nBranch < nBranchMax) );
970  }
971 
972  // Handle veto after ISR + FSR + MPI, but before beam remnants
973  // and resonance decays, e.g. for MLM matching.
974  if (canVetoEarly && userHooksPtr->doVetoPartonLevelEarly( event)) {
975  doVeto = true;
976  if (isDiff) leaveResolvedDiff( iHardLoop, process, event);
977  if (beamHasResGamma) leaveResolvedLeptonGamma( process, event);
978  return false;
979  }
980 
981  // Perform showers in resonance decay chains before beams & reconnection.
982  if (earlyResDec) {
983  int oldSizeEvt = event.size();
984  int oldSizeSys = partonSystemsPtr->sizeSys();
985  if (nBranchMax <= 0 || nBranch < nBranchMax)
986  doVeto = !resonanceShowers( process, event, true);
987  // Abort event if vetoed.
988  if (doVeto) return false;
989 
990  // Reassign new decay products to original system.
991  for (int iSys = oldSizeSys; iSys < partonSystemsPtr->sizeSys(); ++iSys)
992  for (int iOut = 0; iOut < partonSystemsPtr->sizeOut(iSys); ++iOut)
993  partonSystemsPtr->addOut(0, partonSystemsPtr->getOut( iSys, iOut) );
994  partonSystemsPtr->setSizeSys( oldSizeSys);
995 
996  // Perform decays and showers of W and Z emitted in shower.
997  // To do:check if W/Z emission is on in ISR or FSR??
998  if (!wzDecayShowers( event)) return false;
999 
1000  // User hook to reconnect colours specifically in resonance decays.
1001  if (canReconResSys && !userHooksPtr->doReconnectResonanceSystems(
1002  oldSizeEvt, event)) return false;
1003  }
1004 
1005  // Find the first particle in the current diffractive system.
1006  int iFirst = 0;
1007  if (isDiff) {
1008  doDiffCR = isDiff;
1009  iFirst = (iHardLoop == 1) ? 5 + sizeEvent - sizeProcess : sizeEvent;
1010  if (isDiffC) iFirst = 6 + sizeEvent - sizeProcess;
1011  }
1012 
1013  // Change the first particle for hard diffraction.
1014  if (infoPtr->hasPomPsystem()) {
1015  doDiffCR = true;
1016  iFirst = 5;
1017  }
1018 
1019  // Add beam remnants, including primordial kT kick and colour tracing.
1020  if (!doTrial && physical && doRemnants
1021  && (!beamHasGamma || gammaModeEvent != 4)
1022  && !remnants.add( event, iFirst, doDiffCR)) physical = false;
1023 
1024  // If no problems then done.
1025  if (physical) break;
1026 
1027  // Else restore and loop, but do not throw existing diffractive system.
1028  if (!isDiff) event.clear();
1029  else {
1030  event.restoreSize();
1031  event.restoreJunctionSize();
1032  }
1033  beamAPtr->clear();
1034  beamBPtr->clear();
1035  partonSystemsPtr->clear();
1036 
1037  // Restore also the lepton beams if include photons.
1038  if (beamAhasResGamma) beamHadAPtr->clear();
1039  if (beamBhasResGamma) beamHadBPtr->clear();
1040  if (infoPtr->isVMDstateA()) beamVMDAPtr->clear();
1041  if (infoPtr->isVMDstateB()) beamVMDBPtr->clear();
1042 
1043  // End loop over ten tries. Restore from diffraction. Hopefully it worked.
1044  }
1045  if (isDiff) leaveResolvedDiff( iHardLoop, process, event);
1046 
1047  if (!physical) {
1048  // Leave hard diffractive system properly if beam remnant failed.
1049  if (infoPtr->hasPomPsystem()) leaveHardDiff( process, event);
1050  // Leave also photon from lepton framework.
1051  if (beamHasGamma) leaveResolvedLeptonGamma( process, event, false);
1052  return false;
1053  }
1054 
1055  // End big outer loop to handle two systems in double diffraction.
1056  }
1057 
1058  // If no additional MPI has been found then set up the diffractive
1059  // system the first time around.
1060  if (isHardDiff && sampleTypeDiff%2 == 0 && iHardDiffLoop == 1 && nMPI == 1) {
1061  event.clear();
1062  beamAPtr->clear();
1063  beamBPtr->clear();
1064  partonSystemsPtr->clear();
1065  setupHardDiff( process);
1066  continue;
1067  }
1068 
1069  // Do colour reconnection for non-diffractive events before resonance decays.
1070  if (doReconnect && !doDiffCR && reconnectMode > 0) {
1071  Event eventSave = event;
1072  bool colCorrect = false;
1073  for (int i = 0; i < 10; ++i) {
1074  colourReconnection.next(event, 0);
1075  if (junctionSplitting.checkColours(event)) {
1076  colCorrect = true;
1077  break;
1078  }
1079  else event = eventSave;
1080  }
1081  if (!colCorrect) {
1082  infoPtr->errorMsg("Error in PartonLevel::next: "
1083  "Colour reconnection failed.");
1084  return false;
1085  }
1086  }
1087 
1088  // Perform showers in resonance decay chains after beams & reconnection.
1089  int oldSizeEvt = event.size();
1090  if (!earlyResDec) {
1091  if (nBranchMax <= 0 || nBranch < nBranchMax)
1092  doVeto = !resonanceShowers( process, event, true);
1093  // Abort event if vetoed.
1094  if (doVeto) return false;
1095 
1096  // Perform decays and showers of W and Z emitted in shower.
1097  // To do:check if W/Z emission is on in ISR or FSR??
1098  if (!wzDecayShowers( event)) return false;
1099 
1100  // User hook to reconnect colours specifically in resonance decays.
1101  if (canReconResSys && !userHooksPtr->doReconnectResonanceSystems(
1102  oldSizeEvt, event)) return false;
1103  }
1104 
1105  // Store event properties. Not available for diffraction.
1106  if (!isDiff) infoPtr->setEvolution( pTsaveMPI, pTsaveISR, pTsaveFSR,
1107  nMPI, nISR, nFSRinProc, nFSRinRes);
1108  if (isDiff) {
1109  multiPtr->setEmpty();
1110  infoPtr->setImpact( multiPtr->bMPI(), multiPtr->enhanceMPI(),
1111  multiPtr->enhanceMPIavg(), false);
1112  }
1113 
1114  // Do colour reconnection for resonance decays.
1115  if (!earlyResDec && forceResonanceCR && doReconnect &&
1116  !doDiffCR && reconnectMode != 0) {
1117  Event eventSave = event;
1118  bool colCorrect = false;
1119  for (int i = 0; i < 10; ++i) {
1120  colourReconnection.next(event, oldSizeEvt);
1121  if (junctionSplitting.checkColours(event)) {
1122  colCorrect = true;
1123  break;
1124  }
1125  else event = eventSave;
1126  }
1127  if (!colCorrect) {
1128  infoPtr->errorMsg("Error in PartonLevel::next: "
1129  "Colour reconnection failed.");
1130  return false;
1131  }
1132  }
1133 
1134  // Leave diffractive events.
1135  if (isHardDiff) {
1136 
1137  // If inclusive sample wanted for MPI veto and nMPI > 1
1138  // then event is non-diffractive and we can break the loop.
1139  if (sampleTypeDiff == 2 && iHardDiffLoop == 1 && nMPI > 1) {
1140  infoPtr->setHardDiff( false, false, false, false, 0., 0., 0., 0.);
1141  break;
1142  }
1143 
1144  // Leave diffractive system properly.
1145  if (infoPtr->hasPomPsystem()) leaveHardDiff( process, event);
1146  }
1147 
1148  // End big outer loop to handle the setup of the diffractive system.
1149  }
1150 
1151  // If beam photon unresolved during evolution remove the copies of the
1152  // beam particle from the event record.
1153  if ( ( beamAPtr->isGamma() || beamBPtr->isGamma() )
1154  && ( !beamAPtr->resolvedGamma() || !beamBPtr->resolvedGamma() ) ) {
1155  if (!showUnresGamma && (gammaModeEvent != 4) ) cleanEventFromGamma( event);
1156  }
1157 
1158  // After parton level generation, add scattered leptons, restore the event.
1159  if (beamHasGamma) leaveResolvedLeptonGamma( process, event, true);
1160 
1161  // Done.
1162  return true;
1163 
1164 }
1165 
1166 //--------------------------------------------------------------------------
1167 
1168 // Decide which diffractive subsystems are resolved (= perturbative).
1169 
1170 int PartonLevel::decideResolvedDiff( Event& process) {
1171 
1172  // Loop over two systems.
1173  int nHighMass = 0;
1174  int iDSmin = (isDiffC) ? 3 : 1;
1175  int iDSmax = (isDiffC) ? 3 : 2;
1176  for (int iDSnow = iDSmin; iDSnow <= iDSmax; ++iDSnow) {
1177  int iDiffMot = iDSnow + 2;
1178 
1179  // Only high-mass diffractive systems should be resolved.
1180  double mDiff = process[iDiffMot].m();
1181  bool isHighMass = ( mDiff > mMinDiff && rndmPtr->flat()
1182  < pMaxDiff * ( 1. - exp( -(mDiff - mMinDiff) / mWidthDiff ) ) );
1183 
1184  // Set outcome and done.
1185  if (isHighMass) ++nHighMass;
1186  if (iDSnow == 1) isResolvedA = isHighMass;
1187  if (iDSnow == 2) isResolvedB = isHighMass;
1188  if (iDSnow == 3) isResolvedC = isHighMass;
1189  }
1190  return nHighMass;
1191 
1192 }
1193 
1194 //--------------------------------------------------------------------------
1195 
1196 // Set up an unresolved process, i.e. elastic or diffractive.
1197 
1198 bool PartonLevel::setupUnresolvedSys( Event& process, Event& event) {
1199 
1200  // No hard scale in event.
1201  process.scale( 0.);
1202 
1203  // Copy particles from process to event.
1204  for (int i = 0; i < process.size(); ++ i) event.append( process[i]);
1205 
1206  // Loop to find diffractively excited beams.
1207  for (iDS = 1; iDS < 4; ++iDS)
1208  if ( (iDS == 1 && isDiffA && !isResolvedA)
1209  || (iDS == 2 && isDiffB && !isResolvedB)
1210  || (iDS == 3 && isDiffC && !isResolvedC) ) {
1211  int iBeam = iDS + 2;
1212 
1213  // Diffractive mass. Boost and rotation from diffractive system
1214  // rest frame, aligned along z axis, to event cm frame.
1215  double mDiff = process[iBeam].m();
1216  double m2Diff = mDiff * mDiff;
1217  Vec4 pDiffA = (iDS == 1) ? process[1].p()
1218  : process[1].p() - process[3].p();
1219  Vec4 pDiffB = (iDS == 2) ? process[2].p()
1220  : process[2].p() - process[4].p();
1221  RotBstMatrix MtoCM;
1222  MtoCM.fromCMframe( pDiffA, pDiffB);
1223 
1224  // Beam Particle used for flavour content kicked out by Pomeron.
1225  // Randomize for central diffraction; misses closed gluon loop case.
1226  bool beamSideA = (iDS == 1 || (iDS == 3 && rndmPtr->flat() < 0.5));
1227  BeamParticle* beamPtr = (beamSideA) ? beamAPtr : beamBPtr;
1228  if (iDS == 3) beamPtr = (beamSideA) ? beamPomAPtr : beamPomBPtr;
1229 
1230  // Pick quark or gluon kicked out and flavour subdivision.
1231  beamPtr->newValenceContent();
1232  bool gluonIsKicked = beamPtr->pickGluon(mDiff);
1233  int id1 = beamPtr->pickValence();
1234  int id2 = beamPtr->pickRemnant();
1235 
1236  // Find flavour masses. Scale them down if too big.
1237  double m1 = particleDataPtr->constituentMass(id1);
1238  double m2 = particleDataPtr->constituentMass(id2);
1239  if (m1 + m2 > 0.5 * mDiff) {
1240  double reduce = 0.5 * mDiff / (m1 + m2);
1241  m1 *= reduce;
1242  m2 *= reduce;
1243  }
1244 
1245  // If quark is kicked out, then trivial kinematics in rest frame.
1246  if (!gluonIsKicked) {
1247  double pAbs = sqrt( pow2(m2Diff - m1*m1 - m2*m2)
1248  - pow2(2. * m1 * m2) ) / (2. * mDiff);
1249  if (!beamSideA) pAbs = -pAbs;
1250  double e1 = (m2Diff + m1*m1 - m2*m2) / (2. * mDiff);
1251  double e2 = (m2Diff + m2*m2 - m1*m1) / (2. * mDiff);
1252  Vec4 p1( 0., 0., -pAbs, e1);
1253  Vec4 p2( 0., 0., pAbs, e2);
1254 
1255  // Boost and rotate to event cm frame.
1256  p1.rotbst( MtoCM);
1257  p2.rotbst( MtoCM);
1258 
1259  // Set colours.
1260  int col1, acol1, col2, acol2;
1261  if (particleDataPtr->colType(id1) == 1) {
1262  col1 = event.nextColTag();
1263  acol1 = 0;
1264  col2 = 0;
1265  acol2 = col1;
1266  } else {
1267  col1 = 0;
1268  acol1 = event.nextColTag();
1269  col2 = acol1;
1270  acol2 = 0;
1271  }
1272  // Update process colours to stay in step.
1273  process.nextColTag();
1274 
1275  // Store partons of diffractive system and mark system decayed.
1276  int iDauBeg = event.append( id1, 24, iBeam, 0, 0, 0, col1, acol1,
1277  p1, m1);
1278  int iDauEnd = event.append( id2, 63, iBeam, 0, 0, 0, col2, acol2,
1279  p2, m2);
1280  event[iBeam].statusNeg();
1281  event[iBeam].daughters(iDauBeg, iDauEnd);
1282 
1283  // If gluon is kicked out: share momentum between two remnants.
1284  } else {
1285  double m2Sys, zSys, pxSys, pySys, mTS1, mTS2;
1286  zSys = beamPtr->zShare(mDiff, m1, m2);
1287 
1288  // Provide relative pT kick in remnant. Construct (transverse) masses.
1289  pxSys = beamPtr->pxShare();
1290  pySys = beamPtr->pyShare();
1291  mTS1 = m1*m1 + pxSys*pxSys + pySys*pySys;
1292  mTS2 = m2*m2 + pxSys*pxSys + pySys*pySys;
1293  m2Sys = mTS1 / zSys + mTS2 / (1. - zSys);
1294 
1295  // Momentum of kicked-out massless gluon in diffractive rest frame.
1296  double pAbs = (m2Diff - m2Sys) / (2. * mDiff);
1297  double pLRem = (beamSideA) ? pAbs : -pAbs;
1298  Vec4 pG( 0., 0., -pLRem, pAbs);
1299  Vec4 pRem(0., 0., pLRem, mDiff - pAbs);
1300 
1301  // Momenta of the two beam remnant flavours. (Lightcone p+ = m_diff!)
1302  double e1 = 0.5 * (zSys * mDiff + mTS1 / (zSys * mDiff));
1303  double pL1 = 0.5 * (zSys * mDiff - mTS1 / (zSys * mDiff));
1304  if (!beamSideA) pL1 = -pL1;
1305  Vec4 p1(pxSys, pySys, pL1, e1);
1306  Vec4 p2 = pRem - p1;
1307 
1308  // Boost and rotate to event cm frame. Improve precision.
1309  pG.rotbst( MtoCM);
1310  p1.rotbst( MtoCM);
1311  p2.rotbst( MtoCM);
1312  pG.e( pG.pAbs());
1313 
1314  // Set colours.
1315  int colG, acolG, col1, acol1, col2, acol2;
1316  if (particleDataPtr->colType(id1) == 1) {
1317  col1 = event.nextColTag();
1318  acol1 = 0;
1319  colG = event.nextColTag();
1320  acolG = col1;
1321  col2 = 0;
1322  acol2 = colG;
1323  } else {
1324  col1 = 0;
1325  acol1 = event.nextColTag();
1326  colG = acol1;
1327  acolG = event.nextColTag();
1328  col2 = acolG;
1329  acol2 = 0;
1330  }
1331  // Update process colours to stay in step.
1332  process.nextColTag();
1333  process.nextColTag();
1334 
1335  // Store partons of diffractive system and mark system decayed.
1336  int iDauBeg = event.append( 21, 24, iBeam, 0, 0, 0, colG, acolG, pG, 0.);
1337  event.append( id1, 63, iBeam, 0, 0, 0, col1, acol1, p1, m1);
1338  int iDauEnd = event.append( id2, 63, iBeam, 0, 0, 0, col2, acol2,
1339  p2, m2);
1340  event[iBeam].statusNeg();
1341  event[iBeam].daughters(iDauBeg, iDauEnd);
1342  }
1343 
1344  // End loop over beams. Done.
1345  }
1346  return true;
1347 
1348 }
1349 
1350 //--------------------------------------------------------------------------
1351 
1352 // Set up the hard process(es), excluding subsequent resonance decays.
1353 
1354 void PartonLevel::setupHardSys( Event& process, Event& event) {
1355 
1356  // Incoming partons to hard process are stored in slots 3 and 4.
1357  int inS = 0;
1358  int inP = 3;
1359  int inM = 4;
1360 
1361  // Mother and last entry of diffractive system. Offset.
1362  int iDiffMot = iDS + 2;
1363  int iDiffDau = process.size() - 1;
1364  int nOffset = sizeEvent - sizeProcess;
1365 
1366  // Corrected information for hard diffraction.
1367  if (infoPtr->hasPomPsystem()) {
1368  iDiffMot = (isHardDiffB) ? 4 : 3;
1369  inS = iDiffMot;
1370  inP = 7;
1371  inM = 8;
1372  }
1373 
1374  // Resolved diffraction means more entries.
1375  if (isDiff) {
1376  inS = iDiffMot;
1377  inP = iDiffDau - 3;
1378  inM = iDiffDau - 2;
1379 
1380  // Diffractively excited particle described as Pomeron-hadron beams.
1381  event[inS].statusNeg();
1382  event[inS].daughters( inP - 2 + nOffset, inM - 2 + nOffset);
1383  }
1384 
1385  // If photons inside leptons more entries in event.
1386  // Further offset if hard diffraction set up.
1387  int nGammaOffset = 0;
1388  int nGammaDiffOffset = 0;
1389  if ( beamHasResGamma || (beamHasGamma && isGammaHadronDir) ) {
1390  nGammaOffset = 2;
1391  inP += nGammaOffset;
1392  inM += nGammaOffset;
1393  if ( (isHardDiffA || isHardDiffB) && hardDiffSet ) nGammaDiffOffset = 4;
1394  }
1395 
1396  // If two hard interactions then find where second begins.
1397  int iBeginSecond = process.size();
1398  if (doSecondHard) {
1399  iBeginSecond = 5;
1400  while (process[iBeginSecond].status() != -21) ++iBeginSecond;
1401  }
1402 
1403  // If incoming partons are massive then recalculate to put them massless.
1404  if (process[inP].m() != 0. || process[inM].m() != 0.) {
1405  double pPos = process[inP].pPos() + process[inM].pPos();
1406  double pNeg = process[inP].pNeg() + process[inM].pNeg();
1407  process[inP].pz( 0.5 * pPos);
1408  process[inP].e( 0.5 * pPos);
1409  process[inP].m( 0.);
1410  process[inM].pz(-0.5 * pNeg);
1411  process[inM].e( 0.5 * pNeg);
1412  process[inM].m( 0.);
1413  }
1414 
1415  // Add incoming hard-scattering partons to list in beam remnants.
1416  double x1 = process[inP].pPos() / process[inS].m();
1417  double x2 = process[inM].pNeg() / process[inS].m();
1418 
1419  // If photon inside gamma calculate x with respect to photon beams.
1420  if ( beamHasResGamma || (beamHasGamma && isGammaHadronDir) ) {
1421  if ( hasGammaA)
1422  beamHadAPtr->append( 3, process[3].id(), beamHadAPtr->xGamma() );
1423  if ( hasGammaB)
1424  beamHadBPtr->append( 4, process[4].id(), beamHadBPtr->xGamma() );
1425  x1 = process[inP].pPos() / ( process[3 + nGammaDiffOffset].p()
1426  + process[4 + nGammaDiffOffset].p() ).mCalc();
1427  x2 = process[inM].pNeg() / ( process[3 + nGammaDiffOffset].p()
1428  + process[4 + nGammaDiffOffset].p() ).mCalc();
1429  }
1430  beamAPtr->append( inP + nOffset, process[inP].id(), x1);
1431  beamBPtr->append( inM + nOffset, process[inM].id(), x2);
1432 
1433  // Scale. Find whether incoming partons are valence or sea. Store.
1434  // When an x-dependent matter profile is used with nonDiffractive,
1435  // trial interactions mean that the valence/sea choice has already
1436  // been made and should be restored here.
1437  double scale = process.scale();
1438  int vsc1, vsc2;
1439  beamAPtr->xfISR( 0, process[inP].id(), x1, scale*scale);
1440  if (isNonDiff && (vsc1 = multiPtr->getVSC1()) != 0)
1441  (*beamAPtr)[0].companion(vsc1);
1442  else if (beamAPtr->isGamma()) vsc1 = beamAPtr->gammaValSeaComp(0);
1443  else vsc1 = beamAPtr->pickValSeaComp();
1444  beamBPtr->xfISR( 0, process[inM].id(), x2, scale*scale);
1445  if (isNonDiff && (vsc2 = multiPtr->getVSC2()) != 0)
1446  (*beamBPtr)[0].companion(vsc2);
1447  else if (beamBPtr->isGamma()) vsc2 = beamBPtr->gammaValSeaComp(0);
1448  else vsc2 = beamBPtr->pickValSeaComp();
1449  bool isVal1 = (vsc1 == -3);
1450  bool isVal2 = (vsc2 == -3);
1451  infoPtr->setValence( isVal1, isVal2);
1452 
1453  // Initialize info needed for subsequent sequential decays + showers.
1454  nHardDone = sizeProcess;
1455  iPosBefShow.resize( process.size() );
1456  fill( iPosBefShow.begin(), iPosBefShow.end(), 0);
1457 
1458  // Add the beam and hard subprocess partons to the event record.
1459  for (int i = sizeProcess; i < iBeginSecond; ++i) {
1460  if (process[i].mother1() > inM) break;
1461  int j = event.append(process[i]);
1462  iPosBefShow[i] = i;
1463 
1464  // Offset history if process and event not in step.
1465  if (nOffset != 0) {
1466  int iOrd = i - iBeginSecond + 7;
1467  if (iOrd == 1 || iOrd == 2)
1468  event[j].offsetHistory( 0, 0, 0, nOffset);
1469  else if (iOrd == 3 || iOrd == 4)
1470  event[j].offsetHistory( 0, nOffset, 0, nOffset);
1471  else if (iOrd == 5 || iOrd == 6)
1472  event[j].offsetHistory( 0, nOffset, 0, 0);
1473  }
1474 
1475  // Currently outgoing ones should not count as decayed.
1476  if (event[j].status() == -22) {
1477  event[j].statusPos();
1478  event[j].daughters(0, 0);
1479  }
1480 
1481  // Complete task of copying hard subsystem into event record.
1482  ++nHardDone;
1483  }
1484 
1485  // Store participating partons as first set in list of all systems.
1486  partonSystemsPtr->addSys();
1487  partonSystemsPtr->setInA(0, inP + nOffset);
1488  partonSystemsPtr->setInB(0, inM + nOffset);
1489  for (int i = inM + 1; i < nHardDone; ++i)
1490  partonSystemsPtr->addOut(0, i + nOffset);
1491  partonSystemsPtr->setSHat( 0,
1492  (event[inP + nOffset].p() + event[inM + nOffset].p()).m2Calc() );
1493  partonSystemsPtr->setPTHat( 0, scale);
1494 
1495  // Identify second hard process where applicable.
1496  // Since internally generated incoming partons are guaranteed massless.
1497  if (doSecondHard) {
1498  int inP2 = iBeginSecond;
1499  int inM2 = iBeginSecond + 1;
1500 
1501  // Add incoming hard-scattering partons to list in beam remnants.
1502  // Not valid if not in rest frame??
1503  x1 = process[inP2].pPos() / process[0].e();
1504  beamAPtr->append( inP2, process[inP2].id(), x1);
1505  x2 = process[inM2].pNeg() / process[0].e();
1506  beamBPtr->append( inM2, process[inM2].id(), x2);
1507 
1508  // Find whether incoming partons are valence or sea.
1509  scale = process.scaleSecond();
1510  beamAPtr->xfISR( 1, process[inP2].id(), x1, scale*scale);
1511  beamAPtr->pickValSeaComp();
1512  beamBPtr->xfISR( 1, process[inM2].id(), x2, scale*scale);
1513  beamBPtr->pickValSeaComp();
1514 
1515  // Add the beam and hard subprocess partons to the event record.
1516  for (int i = inP2; i < process.size(); ++ i) {
1517  int mother = process[i].mother1();
1518  if ( (mother > 2 && mother < inP2) || mother > inM2 ) break;
1519  event.append(process[i]);
1520  iPosBefShow[i] = i;
1521 
1522  // Currently outgoing ones should not count as decayed.
1523  if (event[i].status() == -22) {
1524  event[i].statusPos();
1525  event[i].daughters(0, 0);
1526  }
1527 
1528  // Complete task of copying hard subsystem into event record.
1529  ++nHardDone;
1530  }
1531 
1532  // Store participating partons as second set in list of all systems.
1533  partonSystemsPtr->addSys();
1534  partonSystemsPtr->setInA(1, inP2);
1535  partonSystemsPtr->setInB(1, inM2);
1536  for (int i = inM2 + 1; i < nHardDone; ++i)
1537  partonSystemsPtr->addOut(1, i);
1538  partonSystemsPtr->setSHat( 1,
1539  (event[inP2].p() + event[inM2].p()).m2Calc() );
1540  partonSystemsPtr->setPTHat( 1, scale);
1541 
1542  // End code for second hard process.
1543  }
1544 
1545  // Update event colour tag to maximum in whole process.
1546  int maxColTag = 0;
1547  for (int i = 0; i < process.size(); ++ i) {
1548  if (process[i].col() > maxColTag) maxColTag = process[i].col();
1549  if (process[i].acol() > maxColTag) maxColTag = process[i].acol();
1550  }
1551  event.initColTag(maxColTag);
1552 
1553  // Copy junctions from process to event.
1554  for (int iJun = 0; iJun < process.sizeJunction(); ++iJun) {
1555  // Resonance decay products may not have been copied from process to
1556  // event yet. If so, do not add junctions associated with decays yet.
1557  int kindJunction = process.kindJunction(iJun);
1558  bool doCopy = true;
1559  // For junction types <= 4, check if final-state legs were copied.
1560  if (kindJunction <= 4) {
1561  int iLegF1 = (kindJunction - 1) / 2;
1562  for (int iLeg = iLegF1; iLeg <= 2; ++iLeg) {
1563  bool colFound = false;
1564  for (int i = inM + 1; i < event.size(); ++i) {
1565  int col = (kindJunction % 2 == 1) ? event[i].col() : event[i].acol();
1566  if (col == process.colJunction(iJun,iLeg)) colFound = true;
1567  }
1568  if (!colFound) doCopy = false;
1569  }
1570  }
1571  if (doCopy) {
1572  event.appendJunction( process.getJunction(iJun));
1573  }
1574  }
1575 
1576  // Done.
1577 }
1578 
1579 //--------------------------------------------------------------------------
1580 
1581 // Set up the event for subsequent resonance decays and showers.
1582 
1583 void PartonLevel::setupShowerSys( Event& process, Event& event) {
1584 
1585  // Reset event record to only contain line 0.
1586  event.clear();
1587  event.append( process[0]);
1588 
1589  // Initialize info needed for subsequent sequential decays + showers.
1590  nHardDone = 1;
1591  iPosBefShow.resize( process.size());
1592  fill( iPosBefShow.begin(), iPosBefShow.end(), 0);
1593 
1594  // Add the hard subprocess partons to the event record.
1595  for (int i = 1; i < process.size(); ++i) {
1596  if (process[i].mother1() > 0) break;
1597  int j = event.append(process[i]);
1598  iPosBefShow[i] = i;
1599 
1600  // Currently outgoing ones should not count as decayed.
1601  if (event[j].status() == -22) {
1602  event[j].statusPos();
1603  event[j].daughters(0, 0);
1604  }
1605 
1606  // Complete task of copying hard subsystem into event record.
1607  ++nHardDone;
1608  }
1609 
1610  // Store participating partons as first set in list of all systems.
1611  partonSystemsPtr->clear();
1612  partonSystemsPtr->addSys();
1613  for (int i = 1; i < nHardDone; ++i) partonSystemsPtr->addOut(0, i);
1614  partonSystemsPtr->setSHat( 0, pow2(process[0].m()) );
1615  partonSystemsPtr->setPTHat( 0, 0.5 * process[0].m() );
1616 
1617  // Copy junctions from process to event.
1618  for (int iJun = 0; iJun < process.sizeJunction(); ++iJun) {
1619  // Resonance decay products may not have been copied from process to
1620  // event yet. If so, do not add junctions associated with decays yet.
1621  int kindJunction = process.kindJunction(iJun);
1622  bool doCopy = true;
1623  // For junction types <= 4, check if final-state legs were copied.
1624  if (kindJunction <= 4) {
1625  int iLegF1 = (kindJunction - 1) / 2;
1626  for (int iLeg = iLegF1; iLeg <= 2; ++iLeg) {
1627  bool colFound = false;
1628  for (int i = 1; i < event.size(); ++i) {
1629  int col = (kindJunction % 2 == 1) ? event[i].col() : event[i].acol();
1630  if (col == process.colJunction(iJun,iLeg)) colFound = true;
1631  }
1632  if (!colFound) doCopy = false;
1633  }
1634  }
1635  if (doCopy) {
1636  event.appendJunction( process.getJunction(iJun));
1637  }
1638  }
1639 
1640  // Done.
1641 }
1642 
1643 //--------------------------------------------------------------------------
1644 
1645 // Resolved diffraction: replace full event with diffractive subsystem.
1646 
1647 void PartonLevel::setupResolvedDiff( Event& process) {
1648 
1649  // Mother and last entry of diffractive system.
1650  int iDiffMot = iDS + 2;
1651  int iDiffDau = process.size() - 1;
1652 
1653  // Diffractively excited particle to be replaced by Pomeron-hadron beams
1654  // (or Pomeron-Pomeron beams for central diffraction).
1655  process[iDiffMot].statusNeg();
1656  process[iDiffMot].daughters( iDiffDau + 1, iDiffDau + 2);
1657 
1658  // Diffractive system mass.
1659  double mDiff = process[iDiffMot].m();
1660  double m2Diff = mDiff * mDiff;
1661 
1662  // Set up Pomeron-proton or Pomeron-Pomeron system as if it were
1663  // the complete collision. Set Pomeron "beam particle" massless.
1664  int idDiffA = (iDS == 1) ? process[1].id() : 990;
1665  int idDiffB = (iDS == 2) ? process[2].id() : 990;
1666  double mDiffA = (iDS == 1) ? process[1].m() : 0.;
1667  double mDiffB = (iDS == 2) ? process[2].m() : 0.;
1668  if (idDiffA == 22 && infoPtr->isVMDstateA()) {
1669  idDiffA = (iDS == 1) ? infoPtr->idVMDA() : 990;
1670  mDiffA = (iDS == 1) ? infoPtr->mVMDA() : 0.;
1671  }
1672  if (idDiffB == 22 && infoPtr->isVMDstateB()) {
1673  idDiffB = (iDS == 2) ? infoPtr->idVMDB() : 990;
1674  mDiffB = (iDS == 2) ? infoPtr->mVMDB() : 0.;
1675  }
1676  double m2DiffA = mDiffA * mDiffA;
1677  double m2DiffB = mDiffB * mDiffB;
1678  double eDiffA = 0.5 * (m2Diff + m2DiffA - m2DiffB) / mDiff;
1679  double eDiffB = 0.5 * (m2Diff + m2DiffB - m2DiffA) / mDiff;
1680  double pzDiff = 0.5 * sqrtpos( pow2(m2Diff - m2DiffA - m2DiffB)
1681  - 4. * m2DiffA * m2DiffB ) / mDiff;
1682  process.append( idDiffA, 13, iDiffMot, 0, 0, 0, 0, 0,
1683  0., 0., pzDiff, eDiffA, mDiffA);
1684  process.append( idDiffB, 13, iDiffMot, 0, 0, 0, 0, 0,
1685  0., 0., -pzDiff, eDiffB, mDiffB);
1686 
1687  // Reassign beam pointers to refer to subsystem effective beams.
1688  // If VMD state in gamma, then use VMD instead of gamma as beam pointer.
1689  beamAPtr = (iDS == 1) ? beamHadAPtr : beamPomAPtr;
1690  beamBPtr = (iDS == 2) ? beamHadBPtr : beamPomBPtr;
1691  if (infoPtr->isVMDstateA())
1692  beamAPtr = (iDS == 1) ? beamVMDAPtr : beamPomAPtr;
1693  if (infoPtr->isVMDstateB())
1694  beamBPtr = (iDS == 2) ? beamVMDBPtr : beamPomBPtr;
1695 
1696  // Pretend that the diffractive system is the whole collision.
1697  eCMsave = infoPtr->eCM();
1698  infoPtr->setECM( mDiff);
1699  beamAPtr->newPzE( pzDiff, eDiffA);
1700  beamBPtr->newPzE( -pzDiff, eDiffB);
1701 
1702  // Keep track of pomeron momentum fraction.
1703  if ( beamAPtr->id() == 990 )
1704  beamAPtr->xPom(pow2(mDiff/eCMsave));
1705  if ( beamBPtr->id() == 990 )
1706  beamBPtr->xPom(pow2(mDiff/eCMsave));
1707 
1708  // Beams not found in normal slots 1 and 2.
1709  int beamOffset = (sizeEvent > 0) ? sizeEvent - 1 : 4;
1710  int beamRemnOffset = iDS;
1711  if (beamAPtr->isGamma() || beamBPtr->isGamma()) beamRemnOffset = 4;
1712 
1713  // Reassign beam pointers in other classes.
1714  timesPtr->reassignBeamPtrs( beamAPtr, beamBPtr, beamOffset);
1715  timesDecPtr->reassignBeamPtrs( beamAPtr, beamBPtr, beamOffset);
1716  spacePtr->reassignBeamPtrs( beamAPtr, beamBPtr, beamOffset);
1717  remnants.reassignBeamPtrs( beamAPtr, beamBPtr, beamRemnOffset);
1718  colourReconnection.reassignBeamPtrs( beamAPtr, beamBPtr);
1719 
1720 
1721  // Reassign multiparton interactions pointer to right object.
1722  if (iDS == 1) multiPtr = &multiSDA;
1723  else if (iDS == 2) multiPtr = &multiSDB;
1724  else multiPtr = &multiCD;
1725 
1726 }
1727 
1728 //--------------------------------------------------------------------------
1729 
1730 // Resolved diffraction: restore to original behaviour.
1731 
1732 void PartonLevel::leaveResolvedDiff( int iHardLoop, Event& process,
1733  Event& event) {
1734 
1735  // Reconstruct boost and rotation to event cm frame.
1736  Vec4 pDiffA = (iDS == 1) ? process[1].p()
1737  : process[1].p() - process[3].p();
1738  Vec4 pDiffB = (iDS == 2) ? process[2].p()
1739  : process[2].p() - process[4].p();
1740  RotBstMatrix MtoCM;
1741  MtoCM.fromCMframe( pDiffA, pDiffB);
1742 
1743  // Perform rotation and boost on diffractive system.
1744  for (int i = sizeProcess; i < process.size(); ++i)
1745  process[i].rotbst( MtoCM);
1746  int iFirst = (iHardLoop == 1) ? 5 + sizeEvent - sizeProcess : sizeEvent;
1747  if (isDiffC) iFirst = 6 + sizeEvent - sizeProcess;
1748  for (int i = iFirst; i < event.size(); ++i)
1749  event[i].rotbst( MtoCM);
1750 
1751  // Restore cm energy.
1752  infoPtr->setECM( eCMsave);
1753  beamAPtr->newPzE( event[1].pz(), event[1].e());
1754  beamBPtr->newPzE( event[2].pz(), event[2].e());
1755  // Keeping track of pomeron momentum fraction.
1756  beamAPtr->xPom();
1757  beamBPtr->xPom();
1758 
1759  // Restore beam pointers to incoming hadrons.
1760  beamAPtr = beamHadAPtr;
1761  beamBPtr = beamHadBPtr;
1762 
1763  // Reassign beam pointers in other classes.
1764  timesPtr->reassignBeamPtrs( beamAPtr, beamBPtr, 0);
1765  timesDecPtr->reassignBeamPtrs( beamAPtr, beamBPtr, 0);
1766  spacePtr->reassignBeamPtrs( beamAPtr, beamBPtr, 0);
1767  remnants.reassignBeamPtrs( beamAPtr, beamBPtr, 0);
1768  colourReconnection.reassignBeamPtrs( beamAPtr, beamBPtr);
1769 
1770  // Restore multiparton interactions pointer to default object.
1771  multiPtr = &multiMB;
1772 
1773 }
1774 
1775 //--------------------------------------------------------------------------
1776 
1777 // Set up special handling of hard diffraction.
1778 
1779 void PartonLevel::setupHardDiff( Event& process) {
1780 
1781  // Create a temporary event record holding the info of the hard process.
1782  Event tmpProcess = process;
1783  process.clear();
1784  process.scale(tmpProcess.scale());
1785 
1786  // Add the first three entries: system + incoming beams.
1787  for (int iEntry = 0; iEntry < 3; ++iEntry)
1788  process.append( tmpProcess[iEntry]);
1789 
1790  // Get system info and calculate diffractive system mass.
1791  double eCM = infoPtr->eCM();
1792  double sNow = eCM * eCM;
1793  double xPom = (isHardDiffB) ? infoPtr->xPomeronA()
1794  : infoPtr->xPomeronB();
1795  double phiPom = 2. * M_PI * rndmPtr->flat();
1796  double thetaPom = (isHardDiffB) ? hardDiffraction.getThetaPomeronA()
1797  : hardDiffraction.getThetaPomeronB();
1798  double m2Diff = xPom * sNow;
1799  double mDiff = sqrt(m2Diff);
1800 
1801  // Add possible photon beams from leptons.
1802  if (beamAhasGamma || beamBhasGamma) {
1803  process.append( tmpProcess[3]);
1804  process.append( tmpProcess[4]);
1805  }
1806 
1807  // Particle masses and ids.
1808  // If isHardDiffB (or A) then m1 + m2 -> m1 + mDiff (or mDiff + m2).
1809  // Change an elastically scattered gamma to a rho, since a Pomeron only
1810  // can be taken from a VMD particle, not from the gamma proper.
1811  int id1 = process[1 + gammaOffset].id();
1812  int id2 = process[2 + gammaOffset].id();
1813  int idEl = (isHardDiffB) ? id1 : id2;
1814  if (idEl == 22) idEl = 113;
1815  int idX = (isHardDiffB) ? id2 : id1;
1816  double m1 = process[1 + gammaOffset].m();
1817  double m2 = process[2 + gammaOffset].m();
1818  double mEl = (isHardDiffB) ? m1 : m2;
1819 
1820  // Evaluate momenta of outgoing particles, initially along beam axis.
1821  // Special handling for rho with too high a Breit-Wigner-selected mass,
1822  // where kinematics failures leads to retries.
1823  double m3, m4, s3, s4, lambda34;
1824  if (idEl == 113) {
1825  int nTry = 0;
1826  do {
1827  mEl = particleDataPtr->mSel(113);
1828  m3 = (isHardDiffB) ? mEl : mDiff;
1829  m4 = (isHardDiffA) ? mEl : mDiff;
1830  s3 = pow2(m3);
1831  s4 = pow2(m4);
1832  lambda34 = sqrtpos( pow2( sNow - s3 - s4) - 4. * s3 * s4 );
1833  ++nTry;
1834  } while (lambda34 <= 0. && nTry < 10);
1835  } else {
1836  m3 = (isHardDiffB) ? mEl : mDiff;
1837  m4 = (isHardDiffA) ? mEl : mDiff;
1838  s3 = pow2(m3);
1839  s4 = pow2(m4);
1840  lambda34 = sqrtpos( pow2( sNow - s3 - s4) - 4. * s3 * s4 );
1841  }
1842  double pAbs = 0.5 * lambda34 / eCM;
1843  Vec4 p3 = Vec4( 0., 0., pAbs, 0.5 * (sNow + s3 - s4) / eCM);
1844  Vec4 p4 = Vec4( 0., 0., -pAbs, 0.5 * (sNow + s4 - s3) / eCM);
1845 
1846  // Take copies for later longitudinal boost; then rotate outgoing beams.
1847  Vec4 pD3 = p3;
1848  Vec4 pD4 = p4;
1849  p3.rot( thetaPom, phiPom);
1850  p4.rot( thetaPom, phiPom);
1851 
1852  // Append intermediate states to the event record.
1853  int status3 = (isHardDiffB) ? 14 : 15;
1854  int status4 = (isHardDiffA) ? 14 : 15;
1855  int sign = (idX > 0) ? 1 : -1;
1856  int idDiff = (idX == 22) ? 9900020 : sign * 9902210;
1857  int id3 = (isHardDiffB) ? idEl : idDiff;
1858  int id4 = (isHardDiffA) ? idEl : idDiff;
1859  process.append( id3, status3, 1 + gammaOffset, 0, 0, 0, 0, 0, p3, m3);
1860  process.append( id4, status4, 2 + gammaOffset, 0, 0, 0, 0, 0, p4, m4);
1861 
1862  // Correct event record history accordingly.
1863  process[1 + gammaOffset].daughters(3 + gammaOffset, 0);
1864  process[2 + gammaOffset].daughters(4 + gammaOffset, 0);
1865  int iDiffMot = (isHardDiffB) ? 4 + gammaOffset: 3 + gammaOffset;
1866  int iDiffRad = process.size() - 1;
1867  process[iDiffMot].statusNeg();
1868  process[iDiffMot].daughters( iDiffRad + 1, iDiffRad + 2);
1869 
1870  // Set up Pomeron-particle system as if it were the complete collision.
1871  // Set Pomeron "beam particle" massless.
1872  int idDiffA = (isHardDiffB) ? 990 : process[1 + gammaOffset].id();
1873  int idDiffB = (isHardDiffA) ? 990 : process[2 + gammaOffset].id();
1874  double mDiffA = (isHardDiffB) ? 0. : process[1 + gammaOffset].m();
1875  double mDiffB = (isHardDiffA) ? 0. : process[2 + gammaOffset].m();
1876  double m2DiffA = mDiffA * mDiffA;
1877  double m2DiffB = mDiffB * mDiffB;
1878  double eDiffA = 0.5 * (m2Diff + m2DiffA - m2DiffB) / mDiff;
1879  double eDiffB = 0.5 * (m2Diff + m2DiffB - m2DiffA) / mDiff;
1880  double pzDiff = 0.5 * sqrtpos( pow2(m2Diff - m2DiffA - m2DiffB)
1881  - 4. * m2DiffA * m2DiffB ) / mDiff;
1882  process.append( idDiffA, 13, iDiffMot, 0, 0, 0, 0, 0,
1883  0., 0., pzDiff, eDiffA, mDiffA);
1884  process.append( idDiffB, 13, iDiffMot, 0, 0, 0, 0, 0,
1885  0., 0., -pzDiff, eDiffB, mDiffB);
1886 
1887  // Append hard process.
1888  vector<int> hardParton;
1889  for (int iHard = 3 + gammaOffset; iHard < tmpProcess.size(); ++iHard)
1890  hardParton.push_back( process.append(tmpProcess[iHard]) );
1891 
1892  // Boost the hard partons in z-direction (from pp to Pp system).
1893  Vec4 pDiffA = (isHardDiffA) ? process[1 + gammaOffset].p()
1894  : process[1 + gammaOffset].p() - pD3;
1895  Vec4 pDiffB = (isHardDiffB) ? process[2 + gammaOffset].p()
1896  : process[2 + gammaOffset].p() - pD4;
1897  RotBstMatrix MtoCM;
1898  MtoCM.toCMframe( pDiffA, pDiffB);
1899  for (unsigned int i = 0; i < hardParton.size(); ++i)
1900  process[hardParton[i]].rotbst(MtoCM);
1901 
1902  // Change mothers and daughters after appending hard process.
1903  for (unsigned int j = 0; j < hardParton.size(); ++j) {
1904  int mother1 = (tmpProcess[j + 3 + gammaOffset].mother1() == 0)
1905  ? 0 : tmpProcess[j + 3 + gammaOffset].mother1() + 4;
1906  int mother2 = (tmpProcess[j + 3 + gammaOffset].mother2() == 0)
1907  ? 0 : tmpProcess[j + 3 + gammaOffset].mother2() + 4;
1908  int daughter1 = (tmpProcess[j + 3 + gammaOffset].daughter1() == 0)
1909  ? 0 : tmpProcess[j + 3 + gammaOffset].daughter1() + 4;
1910  int daughter2 = (tmpProcess[j + 3 + gammaOffset].daughter2() == 0)
1911  ? 0 : tmpProcess[j + 3 + gammaOffset].daughter2() + 4;
1912  process[hardParton[j]].mothers( mother1,mother2);
1913  process[hardParton[j]].daughters( daughter1, daughter2);
1914  }
1915 
1916  // Search for pomeron and particle with status codes 13 (beam-inside-beam).
1917  int iPomeron = 0;
1918  int iPrtcl = 0;
1919  for (int i = 0; i < process.size(); ++i) {
1920  if (process[i].id() == 990 && process[i].status() == 13) iPomeron = i;
1921  if (process[i].idAbs() == idX && process[i].status() == 13) iPrtcl = i;
1922  }
1923 
1924  if (isHardDiffB) {
1925  process[iPomeron].daughters(hardParton[0], 0);
1926  process[iPrtcl].daughters(hardParton[1],0);
1927  process[hardParton[0]].mothers(iPomeron,0);
1928  process[hardParton[1]].mothers(iPrtcl, 0);
1929  } else {
1930  process[iPomeron].daughters(hardParton[1], 0);
1931  process[iPrtcl].daughters(hardParton[0],0);
1932  process[hardParton[1]].mothers(iPomeron,0);
1933  process[hardParton[0]].mothers(iPrtcl, 0);
1934  }
1935 
1936  // Negate status of Pomeron and proton
1937  process[iPomeron].statusNeg();
1938  process[iPrtcl].statusNeg();
1939 
1940  // Change state of system to unresolved to avoid aborting from Pythia.
1941  infoPtr->setHasUnresolvedBeams( true);
1942 
1943  // Reassign beam pointers to refer to subsystem effective beams.
1944  // If gamma-in-lepton change to beamGamPtr.
1945  beamAPtr = (isHardDiffB) ? beamPomAPtr
1946  : (beamAhasGamma ? beamGamAPtr : beamHadAPtr);
1947  beamBPtr = (isHardDiffA) ? beamPomBPtr
1948  : (beamBhasGamma ? beamGamBPtr : beamHadBPtr);
1949 
1950  // Pretend that the diffractive system is the whole collision.
1951  eCMsave = infoPtr->eCM();
1952  infoPtr->setECM( mDiff);
1953  beamAPtr->newPzE( pzDiff, eDiffA);
1954  beamBPtr->newPzE( -pzDiff, eDiffB);
1955 
1956  // Beams not found in normal slots 1 and 2.
1957  int beamOffset = 4 + gammaOffset;
1958  int beamRemnOffset = (isHardDiffB) ? 2 : 1;
1959  if (beamAPtr->isGamma() || beamBPtr->isGamma())
1960  beamRemnOffset = 4 + gammaOffset;
1961 
1962  // Reassign beam pointers in other classes.
1963  timesPtr->reassignBeamPtrs( beamAPtr, beamBPtr, beamOffset);
1964  timesDecPtr->reassignBeamPtrs( beamAPtr, beamBPtr, beamOffset);
1965  spacePtr->reassignBeamPtrs( beamAPtr, beamBPtr, beamOffset);
1966  remnants.reassignBeamPtrs( beamAPtr, beamBPtr, beamRemnOffset);
1967  colourReconnection.reassignBeamPtrs( beamAPtr, beamBPtr);
1968 
1969  // Reassign multiparton interactions pointer to right object.
1970  if (isHardDiffA) multiPtr = &multiSDA;
1971  else if (isHardDiffB) multiPtr = &multiSDB;
1972 
1973  // Set the beam offset for MPIs.
1974  multiPtr->setBeamOffset(beamOffset);
1975 
1976  // Hard diffractive system set.
1977  hardDiffSet = true;
1978 
1979  // Done.
1980  infoPtr->setHasPomPsystem( true);
1981 
1982 }
1983 
1984 //--------------------------------------------------------------------------
1985 
1986 // Leave special handling of hard diffraction.
1987 
1988 void PartonLevel::leaveHardDiff( Event& process, Event& event) {
1989 
1990  // Reconstruct boost and rotation to event cm frame.
1991  Vec4 pDiffA = (isHardDiffA) ? process[1 + gammaOffset].p()
1992  : process[1 + gammaOffset].p() - process[3 + gammaOffset].p();
1993  Vec4 pDiffB = (isHardDiffB) ? process[2 + gammaOffset].p()
1994  : process[2 + gammaOffset].p() - process[4 + gammaOffset].p();
1995  RotBstMatrix MtoCM;
1996  MtoCM.fromCMframe( pDiffA, pDiffB);
1997 
1998  // Perform rotation and boost on diffractive system.
1999  for (int i = 5 + gammaOffset; i < process.size(); ++i)
2000  process[i].rotbst( MtoCM);
2001  for (int i = 5 + gammaOffset; i < event.size(); ++i)
2002  event[i].rotbst( MtoCM);
2003 
2004  // Clear diffractive info.
2005  isHardDiffA = isHardDiffB = isHardDiff = false;
2006 
2007  // Restore cm energy.
2008  infoPtr->setECM( eCMsave);
2009  beamAPtr->newPzE( event[1 + gammaOffset].pz(), event[1 + gammaOffset].e());
2010  beamBPtr->newPzE( event[2 + gammaOffset].pz(), event[2 + gammaOffset].e());
2011 
2012  // Restore beam pointers to incoming hadrons.
2013  // If gamma-in-lepton need to change to beamGamPtr.
2014  beamAPtr = beamAhasGamma ? beamGamAPtr : beamHadAPtr;
2015  beamBPtr = beamBhasGamma ? beamGamBPtr : beamHadBPtr;
2016 
2017  // Reassign beam pointers in other classes.
2018  timesPtr->reassignBeamPtrs( beamAPtr, beamBPtr, 0);
2019  timesDecPtr->reassignBeamPtrs( beamAPtr, beamBPtr, 0);
2020  spacePtr->reassignBeamPtrs( beamAPtr, beamBPtr, 0);
2021  remnants.reassignBeamPtrs( beamAPtr, beamBPtr, 0);
2022  colourReconnection.reassignBeamPtrs( beamAPtr, beamBPtr);
2023 
2024  // Reset the beam offset to normal.
2025  multiPtr->setBeamOffset(0);
2026 
2027  // Restore multiparton interactions pointer to default object.
2028  multiPtr = &multiMB;
2029 
2030 }
2031 
2032 //--------------------------------------------------------------------------
2033 
2034 // Replace full event with photon-photon/hadron subsystem.
2035 // Use the sampled kinematics to construct the momenta.
2036 
2037 bool PartonLevel::setupResolvedLeptonGamma( Event& process) {
2038 
2039  // Save the collision energy of the lepton system.
2040  eCMsaveGamma = infoPtr->eCM();
2041 
2042  // Beams not found in normal slots 1 and 2 but 2 step forward.
2043  int beamOffset = 2;
2044  int iBeamA = 3;
2045  int iBeamB = 4;
2046  gammaOffset = beamOffset;
2047 
2048  // Retrieve the information set on GammaKinematics. Use sHat for 2->1.
2049  double mGmGm = (infoPtr->nFinal() > 1) || (gammaModeEvent != 4)
2050  ? infoPtr->eCMsub() : sqrt( infoPtr->sHat());
2051  double m2GmGm = pow2(mGmGm);
2052 
2053  // Massless photons here, virtualities added after parton level evolution.
2054  double m2Gamma1 = hasGammaA ? 0. : pow2(beamAPtr->m());
2055  double m2Gamma2 = hasGammaB ? 0. : pow2(beamBPtr->m());
2056 
2057  // Derive the new momenta in the CM frame of the gamma-gamma system.
2058  double eGamA = 0.5 * (m2GmGm + m2Gamma1 - m2Gamma2) / mGmGm;
2059  double eGamB = 0.5 * (m2GmGm + m2Gamma2 - m2Gamma1) / mGmGm;
2060  double pzGam = 0.5 * sqrtpos( pow2(m2GmGm - m2Gamma1 - m2Gamma2)
2061  - 4. * m2Gamma1 * m2Gamma2 ) / mGmGm;
2062  Vec4 pGammaANew( 0, 0, pzGam, eGamA);
2063  Vec4 pGammaBNew( 0, 0, -pzGam, eGamB);
2064 
2065  // Set the beam momenta to new rest frame of gamma-gamma.
2066  beamGamAPtr->newPzE( pzGam, eGamA);
2067  beamGamBPtr->newPzE( -pzGam, eGamB);
2068 
2069  // Save the original photon momenta.
2070  Vec4 pGammaA = process[iBeamA].p();
2071  Vec4 pGammaB = process[iBeamB].p();
2072 
2073  // Boost the process to gamma-gamma rest frame (no kT added yet).
2074  RotBstMatrix MtoGammaGamma;
2075  MtoGammaGamma.toCMframe( pGammaA, pGammaB);
2076  process.rotbst(MtoGammaGamma);
2077 
2078  // Set momenta of photons according to new m2GmGm.
2079  process[iBeamA].p(pGammaANew);
2080  process[iBeamB].p(pGammaBNew);
2081 
2082  // If another beam not a photon set the mass as well.
2083  if ( !hasGammaA && !(beamBPtr->getGammaMode() == 2) )
2084  process[iBeamA].m( sqrt(m2Gamma1) );
2085  if ( !hasGammaB && !(beamAPtr->getGammaMode() == 2) )
2086  process[iBeamB].m( sqrt(m2Gamma2) );
2087 
2088  // Done for direct-direct processes since no need to reassign beams.
2089  if ( gammaModeEvent == 4) return true;
2090 
2091  // Reassign beam pointers to refer to subsystem effective beams.
2092  if ( hasGammaA) beamAPtr = beamGamAPtr;
2093  else beamAPtr->newPzE( pzGam, eGamA);
2094  if ( hasGammaB) beamBPtr = beamGamBPtr;
2095  else beamBPtr->newPzE( -pzGam, eGamB);
2096 
2097  // Change state of system to unresolved to avoid aborting from Pythia.
2098  if ( (beamAhasResGamma && (!beamBhasResGamma && hasGammaB))
2099  || ( (!beamAhasResGamma && hasGammaA) && beamBhasResGamma) )
2100  infoPtr->setHasUnresolvedBeams( true);
2101 
2102  // Pretend that the gamma-gamma system is the whole collision.
2103  infoPtr->setECM( mGmGm);
2104 
2105  // Reassign beam pointers in other classes.
2106  timesPtr->reassignBeamPtrs( beamAPtr, beamBPtr, beamOffset);
2107  timesDecPtr->reassignBeamPtrs( beamAPtr, beamBPtr, beamOffset);
2108  spacePtr->reassignBeamPtrs( beamAPtr, beamBPtr, beamOffset);
2109  remnants.reassignBeamPtrs( beamAPtr, beamBPtr, beamOffset);
2110  colourReconnection.reassignBeamPtrs( beamAPtr, beamBPtr);
2111 
2112  // Set the MPI to point the gamma-gamma system.
2113  multiPtr = &multiGmGm;
2114  multiPtr->setBeamOffset(beamOffset);
2115 
2116  // Done.
2117  return true;
2118 
2119 }
2120 
2121 //--------------------------------------------------------------------------
2122 
2123 // Move back to CM-frame of colliding leptons. Add also the scatterd lepton
2124 // if remnants are constructed.
2125 
2126 void PartonLevel::leaveResolvedLeptonGamma( Event& process, Event& event,
2127  bool physical) {
2128 
2129  // Find the momenta of incoming leptons and photons.
2130  Vec4 pLeptonA = process[1].p();
2131  Vec4 pLeptonB = process[2].p();
2132 
2133  // Momenta for scattered leptons to be set later according to final state.
2134  Vec4 pLepton1scat;
2135  Vec4 pLepton2scat;
2136 
2137  // Skip the boost of the extra intermediate photon in hard diffraction.
2138  int iSkipHardDiff = -1;
2139 
2140  // Reconstruct boost and rotation to event cm frame.
2141  RotBstMatrix MtoLeptonLepton;
2142  MtoLeptonLepton.toCMframe( pLeptonA, pLeptonB);
2143 
2144  // Boost event to cm frame.
2145  process.rotbst(MtoLeptonLepton);
2146  event.rotbst(MtoLeptonLepton);
2147 
2148  // Restore cm energy.
2149  infoPtr->setECM( eCMsaveGamma);
2150 
2151  // Restore beam pointers to incoming leptons.
2152  if ( hasGammaA) beamAPtr = beamHadAPtr;
2153  if ( hasGammaB) beamBPtr = beamHadBPtr;
2154 
2155  // Get the masses of beam particles.
2156  double m2BeamA = pow2( beamAPtr->m());
2157  double m2BeamB = pow2( beamBPtr->m());
2158 
2159  // Get the original collision energy and derive the lepton energies in CMS.
2160  double sCM = infoPtr->s();
2161  double eCM2A = 0.25 * pow2(sCM + m2BeamA - m2BeamB) / sCM;
2162  double eCM2B = 0.25 * pow2(sCM - m2BeamA + m2BeamB) / sCM;
2163 
2164  // Find the current momenta of photons.
2165  Vec4 pGamma1Orig = process[3].p();
2166  Vec4 pGamma2Orig = process[4].p();
2167  Vec4 pGamma1 = pGamma1Orig;
2168  Vec4 pGamma2 = pGamma2Orig;
2169  double mGamma1 = sqrt(m2BeamA);
2170  double mGamma2 = sqrt(m2BeamB);
2171 
2172  // Separate treatment for 2 -> 1 processes with two direct gammas to preserve
2173  // sampled invariant mass. Corresponds to the usual beam remnant handling.
2174  if ( (infoPtr->nFinal() == 1) && (gammaModeEvent == 4)
2175  && hasGammaA && hasGammaB ) {
2176 
2177  // Save the virtualities and transverse mometum of the photons.
2178  mGamma1 = -sqrt( beamAPtr->Q2Gamma());
2179  mGamma2 = -sqrt( beamBPtr->Q2Gamma());
2180  double kTxA = beamAPtr->gammaKTx();
2181  double kTyA = beamAPtr->gammaKTy();
2182  double kTxB = beamBPtr->gammaKTx();
2183  double kTyB = beamBPtr->gammaKTy();
2184 
2185  // Derive transverse mass and new sHat.
2186  double mT2A = pow2( beamAPtr->gammaKT() ) - pow2(mGamma1);
2187  double mT2B = pow2( beamBPtr->gammaKT() ) - pow2(mGamma2);
2188  double sHatBef = infoPtr->sHat();
2189  double sHatAft = infoPtr->sHat() + pow2(kTxA + kTxB) + pow2(kTyA + kTyB);
2190  double lambda = pow2(sHatAft) + pow2(mT2A) + pow2(mT2B)
2191  - 2. * ( sHatAft * (mT2A + mT2B) + mT2A * mT2B );
2192  double kz2New = 0.25 * lambda / sHatAft;
2193 
2194  // Photon momenta in CM frame.
2195  Vec4 pGamma1New(kTxA, kTyA, sqrt(kz2New), sqrt(kz2New + mT2A) );
2196  Vec4 pGamma2New(kTxB, kTyB, -sqrt(kz2New), sqrt(kz2New + mT2B) );
2197 
2198  // Boost and rotate the new photon momenta.
2199  RotBstMatrix MfromGmGmOrig;
2200  MfromGmGmOrig.toCMframe( pGamma1Orig, pGamma2Orig);
2201  MfromGmGmOrig.invert();
2202  pGamma1New.rotbst(MfromGmGmOrig);
2203  pGamma2New.rotbst(MfromGmGmOrig);
2204 
2205  // Set the new photon momenta.
2206  pGamma1 = pGamma1New;
2207  pGamma2 = pGamma2New;
2208 
2209  // Modify the mass and momenta in event record.
2210  event[3].p( pGamma1);
2211  event[3].m( mGamma1);
2212  event[4].p( pGamma2);
2213  event[4].m( mGamma2);
2214 
2215  // Light-cone momentum removed by the photons.
2216  double rescale = sqrt(sHatAft / sHatBef);
2217  double wPosBef = beamAPtr->xGamma() * infoPtr->eCM();
2218  double wNegBef = beamBPtr->xGamma() * infoPtr->eCM();
2219 
2220  // Transverse mass of the leptons and invariant mass left for remnants.
2221  double w2BeamA = pow2(beamAPtr->gammaKT()) + m2BeamA;
2222  double w2BeamB = pow2(beamBPtr->gammaKT()) + m2BeamB;
2223  double wPosRem = infoPtr->eCM() - rescale * wPosBef;
2224  double wNegRem = infoPtr->eCM() - rescale * wNegBef;
2225  double w2Rem = wPosRem * wNegRem;
2226  double xLepA = 1. - beamAPtr->xGamma();
2227  double xLepB = 1. - beamBPtr->xGamma();
2228 
2229  // Rescaling factors for leptons.
2230  double lambdaRoot = sqrtpos( pow2(w2Rem - w2BeamA - w2BeamB)
2231  - 4. * w2BeamA * w2BeamB );
2232  double rescaleA = (w2Rem + w2BeamA - w2BeamB + lambdaRoot)
2233  / (2. * w2Rem * xLepA);
2234  double rescaleB = (w2Rem + w2BeamB - w2BeamA + lambdaRoot)
2235  / (2. * w2Rem * xLepB);
2236 
2237  // Momenta of the scattered leptons.
2238  double pPosA = rescaleA * xLepA * wPosRem;
2239  double pNegA = w2BeamA / pPosA;
2240  double eLepA = 0.5 * (pPosA + pNegA);
2241  double pzLepA = 0.5 * (pPosA - pNegA);
2242  double pNegB = rescaleB * xLepB * wNegRem;
2243  double pPosB = w2BeamB / pNegB;
2244  double eLepB = 0.5 * (pPosB + pNegB);
2245  double pzLepB = 0.5 * (pPosB - pNegB);
2246  Vec4 pLeptonANew( -kTxA, -kTyA, pzLepA, eLepA);
2247  Vec4 pLeptonBNew( -kTxB, -kTyB, pzLepB, eLepB);
2248 
2249  // Save the momenta of scattered leptons.
2250  pLepton1scat = pLeptonANew;
2251  pLepton2scat = pLeptonBNew;
2252 
2253  // Save the scattering angles of the leptons.
2254  double theta1 = pLepton1scat.theta();
2255  double theta2 = M_PI - pLepton2scat.theta();
2256  infoPtr->setTheta1(theta1);
2257  infoPtr->setTheta2(theta2);
2258  infoPtr->setECMsub(sqrt(sHatBef));
2259  infoPtr->setsHatNew(sHatBef);
2260 
2261  // Otherwise derive the kinematics according to sampled virtualities.
2262  } else {
2263 
2264  // Calculate the new momentum for virtual photon if present for side A.
2265  if ( hasGammaA ) {
2266 
2267  // Get the x_gamma and virtuality and derive mass.
2268  double xGamma1 = beamAPtr->xGamma();
2269  double Q2gamma1 = beamAPtr->Q2Gamma();
2270  mGamma1 = -sqrt( Q2gamma1);
2271  beamGamAPtr->newM( mGamma1);
2272 
2273  // Derive the kinematics with virtuality and kT.
2274  double eGamma1 = xGamma1 * sqrt( eCM2A);
2275  double kz1 = (eCM2A * xGamma1 + 0.5 * Q2gamma1)
2276  / sqrt(eCM2A - m2BeamA);
2277 
2278  // Set the new momemtum and mass.
2279  pGamma1 = Vec4( beamAPtr->gammaKTx(), beamAPtr->gammaKTy(),
2280  kz1, eGamma1);
2281  event[3].p( pGamma1);
2282  event[3].m( mGamma1);
2283 
2284  // For hard diffraction set the new kinematics also for the copy.
2285  if (infoPtr->isHardDiffractiveA()) {
2286  event[7].p( pGamma1);
2287  event[7].m( mGamma1);
2288  iSkipHardDiff = 7;
2289  }
2290  }
2291 
2292  // Calculate the new momentum for virtual photon if present for side B.
2293  if ( hasGammaB ) {
2294 
2295  // Get the x_gamma and virtuality and derive mass.
2296  double xGamma2 = beamBPtr->xGamma();
2297  double Q2gamma2 = beamBPtr->Q2Gamma();
2298  mGamma2 = -sqrt( Q2gamma2);
2299  beamGamBPtr->newM( mGamma2);
2300 
2301  // Derive the kinematics with virtuality and kT.
2302  double eGamma2 = xGamma2 * sqrt( eCM2B);
2303  double kz2 = (eCM2B * xGamma2 + 0.5 * Q2gamma2)
2304  / sqrt(eCM2B - m2BeamB);
2305 
2306  // Save the 4-momentum of photons with sampled kT.
2307  pGamma2 = Vec4( beamBPtr->gammaKTx(), beamBPtr->gammaKTy(),
2308  -kz2, eGamma2);
2309  event[4].p( pGamma2);
2310  event[4].m( mGamma2);
2311 
2312  // For hard diffraction set the new kinematics also for the copy.
2313  if (infoPtr->isHardDiffractiveB()) {
2314  event[8].p( pGamma2);
2315  event[8].m( mGamma2);
2316  iSkipHardDiff = 8;
2317  }
2318  }
2319 
2320  // Find momenta for scattered lepton.
2321  Vec4 pLepton1 = process[1].p();
2322  Vec4 pLepton2 = process[2].p();
2323  pLepton1scat = pLepton1 - pGamma1;
2324  pLepton2scat = pLepton2 - pGamma2;
2325  }
2326 
2327  // Find the boost from rest frame of collinear photons to rest frame of
2328  // photons with kT.
2329  RotBstMatrix MfromGmGm;
2330  MfromGmGm.toCMframe( pGamma1, pGamma2);
2331  MfromGmGm.fromCMframe( pGamma1Orig, pGamma2Orig);
2332  MfromGmGm.invert();
2333 
2334  // Copy the momentum and mass of the unresolved photon for direct-resolved
2335  // processes to have correct virtualities in the event record.
2336  int iSkipGamma = -1;
2337  if ( gammaModeEvent == 3 ) {
2338  iSkipGamma = 5;
2339  event[iSkipGamma].m( mGamma1);
2340  event[iSkipGamma].p( pGamma1);
2341  } else if ( gammaModeEvent == 2 ) {
2342  iSkipGamma = 6;
2343  event[iSkipGamma].m( mGamma2);
2344  event[iSkipGamma].p( pGamma2);
2345  }
2346 
2347  // Boost scattered system to frame where photon beam has non-zero kT.
2348  for (int i = 5; i < event.size(); ++i) {
2349  if ( (i != iSkipGamma) && (i != iSkipHardDiff) )
2350  event[i].rotbst( MfromGmGm);
2351  }
2352 
2353  // Add the scattered leptons if remnants are constructed and event allowed.
2354  if ( doRemnants && physical) {
2355 
2356  // Add scattered leptons and fix the daughter codes.
2357  if ( hasGammaA) {
2358  int iPosLepton1 = event.append( beamAPtr->id(), 63, 1, 0, 0, 0, 0, 0,
2359  pLepton1scat, beamAPtr->m());
2360  event[1].daughter2( event[1].daughter1());
2361  event[1].daughter1( iPosLepton1);
2362  }
2363  if ( hasGammaB) {
2364  int iPosLepton2 = event.append( beamBPtr->id(), 63, 2, 0, 0, 0, 0, 0,
2365  pLepton2scat, beamBPtr->m());
2366  event[2].daughter2( event[2].daughter1());
2367  event[2].daughter1( iPosLepton2);
2368  }
2369  }
2370 
2371  // Done for direct-direct processes.
2372  if ( gammaModeEvent == 4) return;
2373 
2374  // Reassign beam pointers in other classes.
2375  timesPtr->reassignBeamPtrs( beamAPtr, beamBPtr, 0);
2376  timesDecPtr->reassignBeamPtrs( beamAPtr, beamBPtr, 0);
2377  spacePtr->reassignBeamPtrs( beamAPtr, beamBPtr, 0);
2378  remnants.reassignBeamPtrs( beamAPtr, beamBPtr, 0);
2379  colourReconnection.reassignBeamPtrs( beamAPtr, beamBPtr);
2380 
2381  // Set the MPI pointer back to the original collisions.
2382  multiPtr = &multiMB;
2383  multiPtr->setBeamOffset(0);
2384 
2385 }
2386 
2387 //--------------------------------------------------------------------------
2388 
2389 // Remove the copies of the beam photon from the event record.
2390 
2391 void PartonLevel::cleanEventFromGamma( Event& event) {
2392 
2393  // Offset to normal beam position when photons emitted from a lepton beam.
2394  int beamOffset = 2;
2395  // If hard diffraction more offset for sub-system beams.
2396  if (infoPtr->isHardDiffractiveA() || infoPtr->isHardDiffractiveB())
2397  beamOffset += 4;
2398  int iPosBeam1 = 1 + beamOffset;
2399  int iPosBeam2 = 2 + beamOffset;
2400 
2401  // Go through the event record from the end and find the copies.
2402  int iPosGamma1 = 0;
2403  int iPosGamma2 = 0;
2404  for (int i = event.size() - 1; i > 0; --i) {
2405  if ( (event[i].id() == 22 && event[i].mother1() == iPosBeam1)
2406  && beamAhasResGamma ) iPosGamma1 = i;
2407  if ( (event[i].id() == 22 && event[i].mother1() == iPosBeam2)
2408  && beamBhasResGamma ) iPosGamma2 = i;
2409  }
2410 
2411  // Check how many unresolved photons are present in the event.
2412  int nGamma = 0;
2413  if (iPosGamma1 > 0) ++nGamma;
2414  if (iPosGamma2 > 0) ++nGamma;
2415 
2416  // Exit if no copies found.
2417  if ( nGamma == 0 ) return;
2418 
2419  // Loop over two beams if required.
2420  for (int i = 0; i < nGamma; ++i) {
2421 
2422  // Set the positions to match the beam.
2423  int iPosGamma = (iPosGamma1 > 0 && i == 0) ? iPosGamma1 : iPosGamma2;
2424  int iPosBeam = (iPosGamma1 > 0 && i == 0) ? iPosBeam1 : iPosBeam2;
2425 
2426  // Go through the history of the beam photon.
2427  while ( iPosGamma > iPosBeam ) {
2428  int iDaughter1 = event[iPosGamma].daughter1();
2429  int iDaughter2 = event[iPosGamma].daughter2();
2430  int iMother1 = event[iPosGamma].mother1();
2431  int iMother2 = event[iPosGamma].mother2();
2432 
2433  // If equal daughters photon just a carbon copy.
2434  if ( iDaughter1 == iDaughter2 ) {
2435  event[iDaughter1].mothers( iMother1, iMother2 );
2436  event.remove( iPosGamma, iPosGamma);
2437  iPosGamma = iDaughter1;
2438 
2439  // If non-equal daughters the photon from ISR branching.
2440  } else {
2441  event[iMother1].daughters( iDaughter1, iDaughter2 );
2442  event[iDaughter1].mother1( iMother1 );
2443  event[iDaughter2].mother1( iMother1 );
2444  event.remove( iPosGamma, iPosGamma);
2445  iPosGamma = iMother1;
2446  }
2447 
2448  // If both beams unresolved fix the position of the latter one.
2449  if ( (i == 0 && nGamma > 1) && iPosGamma2 > iPosGamma ) --iPosGamma2;
2450 
2451  }
2452  }
2453 }
2454 
2455 //--------------------------------------------------------------------------
2456 
2457 // Handle showers in successive resonance decays.
2458 
2459 bool PartonLevel::resonanceShowers( Event& process, Event& event,
2460  bool skipForR) {
2461 
2462  // Prepare to start over from beginning for R-hadron decays.
2463  if (allowRH) {
2464  if (skipForR) {
2465  nHardDoneRHad = nHardDone;
2466  inRHadDecay.resize(0);
2467  for (int i = 0; i < process.size(); ++i)
2468  inRHadDecay.push_back( false);
2469  } else nHardDone = nHardDoneRHad;
2470  }
2471 
2472  // Isolate next system to be processed, if anything remains.
2473  int nRes = 0;
2474  int nFSRres = 0;
2475  // Number of desired branchings, negative value means no restriction.
2476  int nBranchMax = (doTrial) ? nTrialEmissions : -1;
2477  // Vector to tell which junctions have already been copied
2478  vector<int> iJunCopied;
2479 
2480  while (nHardDone < process.size()) {
2481  ++nRes;
2482  int iBegin = nHardDone;
2483 
2484  // In first call (skipForR = true) skip over daughters
2485  // of resonances that should form R-hadrons
2486  if (allowRH) {
2487  if (skipForR) {
2488  bool comesFromR = false;
2489  int iTraceUp = process[iBegin].mother1();
2490  do {
2491  if ( rHadronsPtr->givesRHadron(process[iTraceUp].id()) )
2492  comesFromR = true;
2493  iTraceUp = process[iTraceUp].mother1();
2494  } while (iTraceUp > 0 && !comesFromR);
2495  if (comesFromR) {
2496  inRHadDecay[iBegin] = true;
2497  ++nHardDone;
2498  continue;
2499  }
2500 
2501  // In optional second call (skipForR = false) process decay chains
2502  // inside R-hadrons.
2503  } else if (!inRHadDecay[iBegin]) {
2504  ++nHardDone;
2505  continue;
2506  }
2507  }
2508 
2509  // Mother in hard process and in complete event (after shower).
2510  int iHardMother = process[iBegin].mother1();
2511  Particle& hardMother = process[iHardMother];
2512  int iBefMother = iPosBefShow[iHardMother];
2513  int iAftMother = event[iBefMother].iBotCopyId();
2514  // Possibly trace across intermediate R-hadron state.
2515  if (allowRH) {
2516  int iTraceRHadron = rHadronsPtr->trace( iAftMother);
2517  if (iTraceRHadron > 0) iAftMother = iTraceRHadron;
2518  }
2519  Particle& aftMother = event[iAftMother];
2520 
2521  // From now on mother counts as decayed.
2522  aftMother.statusNeg();
2523 
2524  // Mother can have been moved by showering (in any of previous steps),
2525  // so prepare to update colour and momentum information for system.
2526  int colBef = hardMother.col();
2527  int acolBef = hardMother.acol();
2528  int colAft = aftMother.col();
2529  int acolAft = aftMother.acol();
2530  RotBstMatrix M;
2531  M.bst( hardMother.p(), aftMother.p());
2532 
2533  // New colour reconnection can not handle late resonance decay
2534  // of coloured particles so abort event.
2535  if ( (colBef != 0 || acolBef != 0) && doReconnect && reconnectMode == 1
2536  && forceResonanceCR && !earlyResDec) {
2537  infoPtr->errorMsg("Abort from PartonLevel::resonanceShower: "
2538  "new CR can't handle separate CR for coloured resonance decays");
2539  infoPtr->setAbortPartonLevel(true);
2540  return false;
2541  }
2542 
2543  // Extract next partons from hard event into normal event record.
2544  vector<bool> doCopyJun;
2545  for (int i = iBegin; i < process.size(); ++i) {
2546  if (process[i].mother1() != iHardMother) break;
2547  int iNow = event.append( process[i] );
2548  iPosBefShow[i] = iNow;
2549  Particle& now = event.back();
2550 
2551  // Currently outgoing ones should not count as decayed.
2552  if (now.status() == -22) {
2553  now.statusPos();
2554  now.daughters(0, 0);
2555  }
2556 
2557  // Update daughter and mother information.
2558  if (i == iBegin) event[iAftMother].daughter1( iNow);
2559  else event[iAftMother].daughter2( iNow);
2560  now.mother1(iAftMother);
2561 
2562  // Check if this parton came from a BNV (junction) decay in hard event.
2563  for (int iJun = 0; iJun < process.sizeJunction(); ++iJun) {
2564  if (iJun >= int(doCopyJun.size())) doCopyJun.push_back(false);
2565  // Skip if we already decided we're going to copy this junction
2566  if (doCopyJun[iJun]) continue;
2567  // Only consider junctions that can appear in decays.
2568  int kindJunction = process.kindJunction(iJun);
2569  if (kindJunction <= 2) {
2570  // Junction Kinds 1 and 2: all legs in final state
2571  int nMatch = 0;
2572  // Loop over resonance-decay daughters
2573  int iDau1 = hardMother.daughter1();
2574  int iDau2 = hardMother.daughter2();
2575  // Must have at least 3 decay products
2576  if (iDau1 == 0 || iDau2 - iDau1 < 2) continue;
2577  for (int iDau=iDau1; iDau<=iDau2; ++iDau) {
2578  int colDau = (kindJunction == 1 ? process[iDau].col()
2579  : process[iDau].acol());
2580  for (int iLeg = 0; iLeg <= 2; ++iLeg)
2581  if ( process.colJunction(iJun,iLeg) == colDau ) nMatch += 1;
2582  }
2583  // If three legs match
2584  if (nMatch == 3) doCopyJun[iJun] = true;
2585  } else if (kindJunction <= 4) {
2586  // Junction Kinds 3 and 4: copy if initial-state leg matches
2587  // this resonance.
2588  int col = (kindJunction == 3 ? hardMother.acol() : hardMother.col());
2589  if ( process.colJunction(iJun,0) == col ) doCopyJun[iJun] = true;
2590  }
2591  // Extra safety: Check if this junction has already been copied
2592  // (e.g., in setupHardSys). If so, do not copy again.
2593  for (int kJun = 0; kJun < event.sizeJunction(); ++kJun) {
2594  int nMatch = 0;
2595  for (int iLeg = 0; iLeg <= 2; ++iLeg)
2596  if (event.colJunction(kJun,iLeg) == process.colJunction(iJun,iLeg))
2597  ++nMatch;
2598  if (nMatch == 3) doCopyJun[iJun] = false;
2599  }
2600  }
2601 
2602  // Update colour and momentum information.
2603  if (now.col() == colBef) now.col( colAft);
2604  if (now.acol() == acolBef) now.acol( acolAft);
2605  // Sextet mothers have additional (negative) tag
2606  if (now.col() == -acolBef) now.col( -acolAft);
2607  if (now.acol() == -colBef) now.acol( -colAft);
2608  now.rotbst( M);
2609 
2610  // Update vertex information.
2611  if (now.hasVertex()) now.vProd( event[iAftMother].vDec() );
2612 
2613  // Complete task of copying next subsystem into event record.
2614  ++nHardDone;
2615  }
2616  int iEnd = nHardDone - 1;
2617 
2618  // Copy down junctions from hard event into normal event record.
2619  for (int iJun = 0; iJun < int(doCopyJun.size()); ++iJun) {
2620  // Check if this junction was already copied
2621  for (int jJun = 0; jJun < int(iJunCopied.size()); ++jJun)
2622  if (iJunCopied[jJun] == iJun) doCopyJun[iJun] = false;
2623  // Skip if not doing anything
2624  if (!doCopyJun[iJun]) continue;
2625  // Check for changed colors and update as necessary.
2626  Junction junCopy = process.getJunction(iJun);
2627  for (int iLeg = 0; iLeg <= 2; ++iLeg) {
2628  int colLeg = junCopy.col(iLeg);
2629  if (colLeg == colBef) junCopy.col(iLeg, colAft);
2630  if (colLeg == acolBef) junCopy.col(iLeg, acolAft);
2631  }
2632  event.appendJunction(junCopy);
2633  // Mark junction as copied (to avoid later recopying)
2634  iJunCopied.push_back(iJun);
2635  }
2636 
2637  // Reset pT of last branching
2638  pTLastBranch = 0.0;
2639 
2640  // Add new system, automatically with two empty beam slots.
2641  int iSys = partonSystemsPtr->addSys();
2642  partonSystemsPtr->setSHat(iSys, pow2(hardMother.m()) );
2643  partonSystemsPtr->setPTHat(iSys, 0.5 * hardMother.m() );
2644 
2645  // Loop over allowed range to find all final-state particles.
2646  for (int i = iPosBefShow[iBegin]; i <= iPosBefShow[iEnd]; ++i)
2647  if (event[i].isFinal()) partonSystemsPtr->addOut( iSys, i);
2648 
2649  // Do parton showers inside subsystem: maximum scale by mother mass.
2650  if (doFSRinResonances) {
2651  double pTmax = 0.5 * hardMother.m();
2652  if (canSetScale) pTmax
2653  = userHooksPtr->scaleResonance( iAftMother, event);
2654  nFSRhard = 0;
2655 
2656  // Set correct scale for trial showers.
2657  if (doTrial) pTmax = process.scale();
2658 
2659  // Set correct scale for showers off multi-parton events when
2660  // merging e+e- -> V -> jets.
2661  int iMother1 = hardMother.mother1();
2662  int iMother2 = hardMother.mother2();
2663  if ( canRemoveEvent
2664  && event[iMother1].colType() == 0 && event[iMother2].colType() == 0)
2665  pTmax = process.scale();
2666 
2667  // Let prepare routine do the setup.
2668  timesDecPtr->prepare( iSys, event);
2669 
2670  // Number of actual branchings
2671  int nBranch = 0;
2672 
2673  // Set up initial veto scale.
2674  doVeto = false;
2675  double pTveto = pTvetoPT;
2676 
2677  // Begin evolution down in pT from hard pT scale.
2678  do {
2679  typeVetoStep = 0;
2680  double pTtimes = timesDecPtr->pTnext( event, pTmax, 0.);
2681 
2682  // Allow a user veto. Only do it once, so remember to change pTveto.
2683  if (pTveto > 0. && pTveto > pTtimes) {
2684  pTveto = -1.;
2685  doVeto = userHooksPtr->doVetoPT( 5, event);
2686  // Abort event if vetoed.
2687  if (doVeto) return false;
2688  }
2689 
2690  // Do a final-state emission (if allowed).
2691  if (pTtimes > 0.) {
2692  if (timesDecPtr->branch( event)) {
2693  ++nFSRres;
2694  ++nFSRhard;
2695  if (canVetoStep && nFSRhard <= nVetoStep) typeVetoStep = 5;
2696 
2697  nBranch++;
2698  pTLastBranch = pTtimes;
2699  typeLastBranch = 5;
2700 
2701  }
2702  pTmax = pTtimes;
2703  }
2704 
2705  // If no pT scales above zero then nothing to be done.
2706  else pTmax = 0.;
2707 
2708  // Optionally check for a veto after the first few emissions.
2709  if (typeVetoStep > 0) {
2710  doVeto = userHooksPtr->doVetoStep( typeVetoStep, 0, nFSRhard,
2711  event);
2712  // Abort event if vetoed.
2713  if (doVeto) return false;
2714  }
2715 
2716  // Handle potential merging veto.
2717  if ( canRemoveEvent && nFSRhard == 1 ) {
2718  // Simply check, and possibly reset weights.
2719  mergingHooksPtr->doVetoStep( process, event, true );
2720  }
2721 
2722  // Keep on evolving until nothing is left to be done.
2723  } while (pTmax > 0. && (nBranchMax <= 0 || nBranch < nBranchMax) );
2724 
2725  }
2726 
2727  // No more systems to be processed. Set total number of emissions.
2728  }
2729  if (skipForR) nFSRinRes = nFSRres;
2730  return true;
2731 
2732 }
2733 
2734 //--------------------------------------------------------------------------
2735 
2736 // Perform decays and showers of W and Z emitted in shower.
2737 
2738 bool PartonLevel::wzDecayShowers( Event& event) {
2739 
2740  // Identify W/Z produced by a parton shower.
2741  for (int iWZ = 0; iWZ < event.size(); ++iWZ)
2742  if (event[iWZ].isFinal()
2743  && (event[iWZ].id() == 23 || event[iWZ].idAbs() == 24) ) {
2744 
2745  // Do nothing if particle should not be decayed.
2746  if ( event[iWZ].canDecay() && event[iWZ].mayDecay()
2747  && event[iWZ].isResonance() ) ;
2748  else continue;
2749 
2750  int iWZtop = event[iWZ].iTopCopy();
2751  int typeWZ = 0;
2752  if (event[iWZtop].statusAbs() == 56) typeWZ = 1;
2753  if (event[iWZtop].statusAbs() == 47) typeWZ = 2;
2754  int sizeSave = event.size();
2755 
2756  // Map id_Z = 23 -> 93 and id_W = 24 -> 94, for separate decay settings.
2757  // Let W/Z resonance decay. Restore correct identity and status codes.
2758  if (typeWZ > 0) {
2759  int idSave = event[iWZ].id();
2760  event[iWZ].id( (idSave > 0) ? idSave + 70 : idSave - 70);
2761  int statusSave = event[iWZ].status();
2762  resonanceDecays.next( event, iWZ);
2763  event[iWZ].id( idSave);
2764  if (event.size() - sizeSave != 2) return true;
2765  event[iWZ].status( -statusSave);
2766  }
2767 
2768  // FSR decays.
2769  if (typeWZ == 1) {
2770 
2771  // Identify fermion after W/Z emission.
2772  vector<int> iSisters = event[iWZtop].sisterList();
2773  if (iSisters.size() != 1) {
2774  infoPtr->errorMsg("Error in PartonLevel::wzDecayShowers: "
2775  "Not able to find a single sister particle");
2776  return false;
2777  }
2778  int iEmitter = iSisters[0];
2779 
2780  // Boosts to study decay in W/Z rest frame.
2781  RotBstMatrix MtoNew, MtoRest, MtoCM;
2782  MtoNew.bst( event[iWZtop].p(), event[iWZ].p());
2783  MtoRest.bstback( event[iWZ].p());
2784  MtoCM.bst( event[iWZ].p());
2785 
2786  // Emitter and recoiler in W/Z rest frame.
2787  Vec4 pEmitter = event[iEmitter].p();
2788  pEmitter.rotbst( MtoNew);
2789  pEmitter.rotbst( MtoRest);
2790  if (event[iWZtop + 1].statusAbs() != 52) {
2791  infoPtr->errorMsg("Error in PartonLevel::wzDecayShowers: "
2792  "Found wrong recoiler");
2793  return false;
2794  }
2795  Vec4 pRecoiler = event[iWZtop + 1].p();
2796  pRecoiler.rotbst( MtoNew);
2797  pRecoiler.rotbst( MtoRest);
2798  Vec4 pWZRest = event[iWZ].p();
2799  pWZRest.rotbst( MtoRest);
2800 
2801  // Always choose p4 as the particle and p5 as the anti-particle.
2802  Vec4 p4 = pEmitter;
2803  Vec4 p5 = pRecoiler;
2804  if (event[iEmitter].id() < 0) swap( p4, p5);
2805 
2806  // Decay daughters in W/Z rest frame.
2807  // Always choose pDec1 as the particle and p2Dec as the anti-particle.
2808  Vec4 pDec1 = event[sizeSave].p();
2809  Vec4 pDec2 = event[sizeSave + 1].p();
2810  if (event[sizeSave].id() < 0) swap( pDec1, pDec2);
2811  pDec1.rotbst( MtoRest);
2812  pDec2.rotbst( MtoRest);
2813 
2814  // Couplings.
2815  double li2, ri2, lf2, rf2;
2816  // Z couplings: make use of initial fermion polarization if set.
2817  if (event[iWZ].id() == 23) {
2818  li2 = pow2(couplingsPtr->lf( event[iEmitter].idAbs() ));
2819  ri2 = pow2(couplingsPtr->rf( event[iEmitter].idAbs() ));
2820  lf2 = pow2(couplingsPtr->lf( event[sizeSave].idAbs() ));
2821  rf2 = pow2(couplingsPtr->rf( event[sizeSave].idAbs() ));
2822  if ( abs( event[iEmitter].pol() + 1.) < 0.1) ri2 = 0.;
2823  if ( abs( event[iEmitter].pol() - 1.) < 0.1) li2 = 0.;
2824  // W couplings.
2825  } else {
2826  li2 = 1.;
2827  ri2 = 0.;
2828  lf2 = 1.;
2829  rf2 = 0.;
2830  }
2831 
2832  // Different needed kinematic variables.
2833  double sWZER = (p4 + pWZRest + p5).m2Calc();
2834  double x1 = 2. * p4 * (p4 + pWZRest + p5) / sWZER;
2835  double x2 = 2. * p5 * (p4 + pWZRest + p5) / sWZER;
2836  double x1s = x1 * x1;
2837  double x2s = x2 * x2;
2838  double m2Sel = pWZRest.m2Calc();
2839  double rWZER = m2Sel / sWZER;
2840 
2841  // Calculate constants needed in correction.
2842  double con[9];
2843  con[0] = 2. * m2Sel * (1.-x1) * ((x2s+1.-x1-x2) - rWZER * (1.-x2))
2844  * (li2 * lf2 + ri2 * rf2);
2845  con[1] = 2. * m2Sel * (1.-x2) * ((x1s+1.-x1-x2) - rWZER * (1.-x1))
2846  * (li2 * rf2 + ri2 * lf2);
2847  con[2] = 2. * m2Sel * (1.-x1) * ((x2s+1.-x1-x2) - rWZER * (1.-x2))
2848  * (li2 * rf2 + ri2 * lf2);
2849  con[3] = 2. * m2Sel * (1.-x2) * ((x1s+1.-x1-x2) - rWZER * (1.-x1))
2850  * (li2 * lf2 + ri2 * rf2);
2851  con[4] = m2Sel * sWZER * (1.-x1) * (1.-x2) * ((x1+x2-1.) + rWZER)
2852  * (li2 + ri2) * (lf2 + rf2);
2853  con[5] = -4. * (1.-x1) * (1.-x2) * (li2 + ri2) * (lf2 + rf2);
2854  con[6] = -4. * (1.-x1) * (1.-x2) * (li2 + ri2) * (lf2 + rf2);
2855  con[7] = 4. * (1.-x1) * ((1.-x1) - rWZER * (1.-x2))
2856  * (li2 + ri2) * (lf2 + rf2);
2857  con[8] = 4. * (1.-x2) * ((1.-x2) - rWZER * (1.-x1))
2858  * (li2 + ri2) * (lf2 + rf2);
2859 
2860  // Find maximum value: try pDec1 and pDec2 = -pDec1 along +-x, +-y, +-z.
2861  double wtMax = 0.;
2862  double pAbs12 = pDec1.pAbs();
2863  for (int j = 0; j < 6; ++j) {
2864  Vec4 pDec1Test( 0., 0., 0., pDec1.e());
2865  Vec4 pDec2Test( 0., 0., 0., pDec2.e());
2866  if (j == 0) { pDec1Test.px( pAbs12); pDec1Test.px( -pAbs12);}
2867  else if (j == 1) { pDec1Test.px( -pAbs12); pDec1Test.px( pAbs12);}
2868  else if (j == 2) { pDec1Test.py( pAbs12); pDec1Test.py( -pAbs12);}
2869  else if (j == 3) { pDec1Test.py( -pAbs12); pDec1Test.py( pAbs12);}
2870  else if (j == 4) { pDec1Test.pz( pAbs12); pDec1Test.pz( -pAbs12);}
2871  else if (j == 5) { pDec1Test.pz( -pAbs12); pDec1Test.pz( pAbs12);}
2872 
2873  // Evaluate matrix element and compare with current maximum.
2874  double p2p4Test = p4 * pDec1Test;
2875  double p3p4Test = p4 * pDec2Test;
2876  double p2p5Test = p5 * pDec1Test;
2877  double p3p5Test = p5 * pDec2Test;
2878  double testValues[9] = { p2p4Test, p2p5Test, p3p4Test, p3p5Test, 1.,
2879  p2p5Test * p3p4Test, p2p4Test * p3p5Test, p2p4Test * p3p4Test,
2880  p2p5Test * p3p5Test};
2881  double wtTest = 0.;
2882  for (int i = 0; i < 9; ++i) wtTest += con[i] * testValues[i];
2883  if (wtTest > wtMax) wtMax = wtTest;
2884  }
2885 
2886  // Multiply by four to ensure maximum is an overestimate.
2887  wtMax *= 4.;
2888 
2889  // Iterate with new angles until weighting succeeds.
2890  int nRot = -1;
2891  double wt = 0.;
2892  do {
2893  ++nRot;
2894  if (nRot > 0) {
2895  RotBstMatrix MrndmRot;
2896  MrndmRot.rot( acos(2. * rndmPtr->flat() - 1.),
2897  2. * M_PI * rndmPtr->flat());
2898  pDec1.rotbst(MrndmRot);
2899  pDec2.rotbst(MrndmRot);
2900  }
2901 
2902  // p2 is decay product, p3 is anti decay product,
2903  // p4 is dipole particle, p5 is dipole anti particle.
2904  // So far assumed that we always have qQ-dipole.
2905  double p2p4 = p4 * pDec1;
2906  double p3p4 = p4 * pDec2;
2907  double p2p5 = p5 * pDec1;
2908  double p3p5 = p5 * pDec2;
2909 
2910  // Calculate weight and compare with maximum weight.
2911  double wtValues[9] = { p2p4, p2p5, p3p4, p3p5, 1., p2p5 * p3p4,
2912  p2p4 * p3p5, p2p4 * p3p4, p2p5 * p3p5};
2913  wt = 0.;
2914  for (int i = 0; i < 9; ++i) wt += con[i] * wtValues[i];
2915  if (wt > wtMax || wt < 0.) {
2916  infoPtr->errorMsg("Error in PartonLevel::wzDecayShowers: "
2917  "wt bigger than wtMax or less than zero");
2918  return false;
2919  }
2920  } while (wt < wtMax * rndmPtr->flat());
2921 
2922  // If momenta rotated then store new ones.
2923  if (nRot > 0) {
2924  pDec1.rotbst( MtoCM);
2925  pDec2.rotbst( MtoCM);
2926  if (event[sizeSave].id() > 0) {
2927  event[sizeSave].p( pDec1);
2928  event[sizeSave + 1].p( pDec2);
2929  }
2930  else {
2931  event[sizeSave].p( pDec2);
2932  event[sizeSave + 1].p( pDec1);
2933  }
2934  }
2935  }
2936 
2937  // ISR decays.
2938  if (typeWZ == 2) {
2939 
2940  // Identify mother of W/Z emission.
2941  double iMother = event[iWZtop].mother1();
2942 
2943  // Boosts to study decay in W/Z rest frame.
2944  RotBstMatrix MtoNew, MtoRest, MtoCM;
2945  MtoNew.bst( event[iWZtop].p(), event[iWZ].p());
2946  MtoRest.bstback( event[iWZ].p());
2947  MtoCM.bst( event[iWZ].p());
2948 
2949  // Find recoiler.
2950  int iRecoiler;
2951  if (event[iMother + 1].statusAbs() == 42) iRecoiler = iMother + 1;
2952  else if (event[iMother - 1].statusAbs() == 42) iRecoiler = iMother - 1;
2953  else {
2954  infoPtr->errorMsg("Error in PartonLevel::wzDecayShowers: "
2955  "Found wrong recoiler");
2956  return false;
2957  }
2958 
2959  // Emitter and recoiler in W/Z rest frame.
2960  Vec4 pMother = event[iMother].p();
2961  pMother.rotbst( MtoNew);
2962  pMother.rotbst( MtoRest);
2963  Vec4 pRecoiler = event[iRecoiler].p();
2964  pRecoiler.rotbst( MtoNew);
2965  pRecoiler.rotbst( MtoRest);
2966  Vec4 pWZRest = event[iWZ].p();
2967  pWZRest.rotbst( MtoRest);
2968 
2969  // Always choose p1 as the particle and p2 as the anti-particle.
2970  // If no anti-particles just continue.
2971  Vec4 p1 = pMother;
2972  Vec4 p2 = pRecoiler;
2973  if (event[iMother].id() < 0) swap( p1, p2);
2974 
2975  // Decay daughters in W/Z rest frame.
2976  // Always choose pDec1 as the particle and p2Dec as the anti-particle.
2977  Vec4 pDec1 = event[sizeSave].p();
2978  Vec4 pDec2 = event[sizeSave + 1].p();
2979  if (event[sizeSave].id() < 0) swap( pDec1, pDec2);
2980  pDec1.rotbst( MtoRest);
2981  pDec2.rotbst( MtoRest);
2982 
2983  // Couplings.
2984  double li2, ri2, lf2, rf2;
2985  // Z couplings: make use of initial fermion polarization if set.
2986  if (event[iWZ].id() == 23) {
2987  li2 = pow2(couplingsPtr->lf( event[iMother].idAbs() ));
2988  ri2 = pow2(couplingsPtr->rf( event[iMother].idAbs() ));
2989  lf2 = pow2(couplingsPtr->lf( event[sizeSave].idAbs() ));
2990  rf2 = pow2(couplingsPtr->rf( event[sizeSave].idAbs() ));
2991  if ( abs( event[iMother].pol() + 1.) < 0.1) ri2 = 0.;
2992  if ( abs( event[iMother].pol() - 1.) < 0.1) li2 = 0.;
2993  // W couplings.
2994  } else {
2995  li2 = 1.;
2996  ri2 = 0.;
2997  lf2 = 1.;
2998  rf2 = 0.;
2999  }
3000 
3001  // Different needed kinematic variables.
3002  double s = (p1 + p2).m2Calc();
3003  double t = (p1-pWZRest).m2Calc();
3004  double u = (p2-pWZRest).m2Calc();
3005  double m2Sel = pWZRest.m2Calc();
3006  double m2Final = t + u + s - m2Sel;
3007 
3008  // Calculate constants needed in correction.
3009  double con[9];
3010  con[0] = 2. * m2Sel * (u * (1. - m2Final / t) + s)
3011  * (li2 * rf2 + ri2 * lf2);
3012  con[1] = 2. * m2Sel * (u * (1. - m2Final / t) + s)
3013  * (li2 * lf2 + ri2 * rf2);
3014  con[2] = 2. * m2Sel * (t * (1. - m2Final / u) + s)
3015  * (li2 * lf2 + ri2 * rf2);
3016  con[3] = 2. * m2Sel * (t * (1. - m2Final / u) + s)
3017  * (li2 * rf2 + ri2 * lf2);
3018  con[4] = m2Sel * m2Final * s * (li2 + ri2) * (lf2 + rf2);
3019  con[5] = -4. * m2Final * (li2 + ri2) * (lf2 + rf2);
3020  con[6] = -4. * m2Final * (li2 + ri2) * (lf2 + rf2);
3021  con[7] = 4. * (m2Final * u / t - m2Sel) * (li2 + ri2) * (lf2 + rf2);
3022  con[8] = 4. * (m2Final * t / u - m2Sel) * (li2 + ri2) * (lf2 + rf2);
3023 
3024  // Find maximum value: try pDec1 and pDec2 = -pDec1 along +-x, +-y, +-z.
3025  double wtMax = 0.;
3026  double pAbs12 = pDec1.pAbs();
3027  for (int j = 0; j < 6; ++j) {
3028  Vec4 pDec1Test( 0., 0., 0., pDec1.e());
3029  Vec4 pDec2Test( 0., 0., 0., pDec2.e());
3030  if (j == 0) { pDec1Test.px( pAbs12); pDec1Test.px( -pAbs12);}
3031  else if (j == 1) { pDec1Test.px( -pAbs12); pDec1Test.px( pAbs12);}
3032  else if (j == 2) { pDec1Test.py( pAbs12); pDec1Test.py( -pAbs12);}
3033  else if (j == 3) { pDec1Test.py( -pAbs12); pDec1Test.py( pAbs12);}
3034  else if (j == 4) { pDec1Test.pz( pAbs12); pDec1Test.pz( -pAbs12);}
3035  else if (j == 5) { pDec1Test.pz( -pAbs12); pDec1Test.pz( pAbs12);}
3036 
3037  // Evaluate matrix element and compare with current maximum.
3038  double p1p4Test = p1 * pDec1Test;
3039  double p1p5Test = p1 * pDec2Test;
3040  double p2p4Test = p2 * pDec1Test;
3041  double p2p5Test = p2 * pDec2Test;
3042  double testValues[9] = { p1p4Test, p1p5Test, p2p4Test, p2p5Test, 1.,
3043  p1p4Test * p2p5Test, p1p5Test * p2p4Test, p1p4Test * p1p5Test,
3044  p2p4Test * p2p5Test};
3045  double wtTest = 0.;
3046  for (int i = 0; i < 9; ++i) wtTest += con[i] * testValues[i];
3047  if (wtTest > wtMax) wtMax = wtTest;
3048  }
3049 
3050  // Multiply by four to ensure maximum is an overestimate.
3051  wtMax *= 4.;
3052 
3053  // Iterate with new angles until weighting succeeds.
3054  int nRot = -1;
3055  double wt = 0.;
3056  do {
3057  ++nRot;
3058  if (nRot > 0) {
3059  RotBstMatrix MrndmRot;
3060  MrndmRot.rot( acos(2. * rndmPtr->flat() - 1.),
3061  2. * M_PI * rndmPtr->flat());
3062  pDec1.rotbst(MrndmRot);
3063  pDec2.rotbst(MrndmRot);
3064  }
3065 
3066  // p2 is decay product, p3 is anti decay product,
3067  // p4 is dipole particle, p5 is dipole anti particle.
3068  // So far assumed that we always have qQ-dipole.
3069  double p1p4 = p1 * pDec1;
3070  double p1p5 = p1 * pDec2;
3071  double p2p4 = p2 * pDec1;
3072  double p2p5 = p2 * pDec2;
3073 
3074  // Calculate weight and compare with maximum weight.
3075  double wtValues[9] = { p1p4, p1p5, p2p4, p2p5, 1., p1p4 * p2p5,
3076  p1p5 * p2p4, p1p4 * p1p5, p2p4 * p2p5};
3077  wt = 0.;
3078  for (int i = 0; i < 9; ++i) wt += con[i] * wtValues[i];
3079  if (wt > wtMax || wt < 0.) {
3080  infoPtr->errorMsg("Error in PartonLevel::wzDecayShowers: "
3081  "wt bigger than wtMax or less than zero");
3082  return false;
3083  }
3084  } while (wt < wtMax * rndmPtr->flat());
3085 
3086  // If momenta rotated then store new ones.
3087  if (nRot > 0) {
3088  pDec1.rotbst( MtoCM);
3089  pDec2.rotbst( MtoCM);
3090  if (event[sizeSave].id() > 0) {
3091  event[sizeSave].p( pDec1);
3092  event[sizeSave + 1].p( pDec2);
3093  }
3094  else {
3095  event[sizeSave].p( pDec2);
3096  event[sizeSave + 1].p( pDec1);
3097  }
3098  }
3099  }
3100 
3101  // Add new system, automatically with two empty beam slots.
3102  if (typeWZ > 0) {
3103  // Maximum shower scale set by mother mass.
3104  double pTmax = 0.5 * event[iWZ].m();
3105  int iSys = partonSystemsPtr->addSys();
3106  partonSystemsPtr->setSHat(iSys, pow2(event[iWZ].m()) );
3107  partonSystemsPtr->setPTHat(iSys, pTmax );
3108  for (int i = sizeSave; i < event.size(); ++i)
3109  partonSystemsPtr->addOut( iSys, i);
3110 
3111  // Do parton showers inside subsystem.
3112  if (doFSRinResonances) {
3113 
3114  // Let prepare routine do the setup.
3115  timesDecPtr->prepare( iSys, event);
3116 
3117  // Begin evolution down in pT from hard pT scale.
3118  do {
3119  double pTtimes = timesDecPtr->pTnext( event, pTmax, 0.);
3120 
3121  // Do a final-state emission (if allowed).
3122  if (pTtimes > 0.) {
3123  timesDecPtr->branch( event);
3124  pTmax = pTtimes;
3125  }
3126 
3127  // If no pT scales above zero then nothing to be done.
3128  else pTmax = 0.;
3129 
3130  // Keep on evolving until nothing is left to be done.
3131  } while (pTmax > 0.);
3132  }
3133  }
3134 
3135  // End loop over event to find W/Z gauge bosons.
3136  }
3137 
3138  // Done.
3139  return true;
3140 
3141 }
3142 
3143 //==========================================================================
3144 
3145 } // end namespace Pythia8
Definition: AgUStep.h:26