PYTHIA  8.313
SimpleTimeShower.h
1 // SimpleTimeShower.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2025 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Header file for the original simple timelike final-state showers.
7 // TimeDipoleEnd: data on a radiating dipole end in FSR.
8 // SimpleTimeShower: handles the showering description.
9 
10 #ifndef Pythia8_SimpleTimeShower_H
11 #define Pythia8_SimpleTimeShower_H
12 
13 #include "Pythia8/TimeShower.h"
14 #include "Pythia8/SimpleWeakShowerMEs.h"
15 
16 namespace Pythia8 {
17 
18 // Forward declaration of SplitOnia;
19 class SplitOnia;
20 
21 //==========================================================================
22 
23 // Data on radiating dipole ends; only used inside SimpleTimeShower class.
24 
26 
27 public:
28 
29  // Constructors.
30  TimeDipoleEnd() : iRadiator(-1), iRecoiler(-1), pTmax(0.), colType(0),
31  chgType(0), gamType(0), weakType(0), isrType(0), system(0), systemRec(0),
32  MEtype(0), iMEpartner(-1), weakPol(0), oniumType(0), isHiddenValley(false),
33  colvType(0), MEmix(0.), MEorder(true), MEsplit(true), MEgluinoRec(false),
34  isFlexible(false), hasJunction(false), flavour(), iAunt(),
35  mRad(), m2Rad(), mRec(), m2Rec(), mDip(), m2Dip(), m2DipCorr(), pT2(),
36  m2(), z(), mFlavour(), asymPol(), flexFactor(), pAccept() {}
37  TimeDipoleEnd(int iRadiatorIn, int iRecoilerIn, double pTmaxIn = 0.,
38  int colIn = 0, int chgIn = 0, int gamIn = 0, int weakTypeIn = 0,
39  int isrIn = 0, int systemIn = 0, int MEtypeIn = 0, int iMEpartnerIn = -1,
40  int weakPolIn = 0, int oniumTypeIn = 0,
41  bool isHiddenValleyIn = false, int colvTypeIn = 0, double MEmixIn = 0.,
42  bool MEorderIn = true, bool MEsplitIn = true, bool MEgluinoRecIn = false,
43  bool isFlexibleIn = false, bool hasJunctionIn = false) :
44  iRadiator(iRadiatorIn), iRecoiler(iRecoilerIn), pTmax(pTmaxIn),
45  colType(colIn), chgType(chgIn), gamType(gamIn), weakType(weakTypeIn),
46  isrType(isrIn), system(systemIn), systemRec(systemIn), MEtype(MEtypeIn),
47  iMEpartner(iMEpartnerIn), weakPol(weakPolIn), oniumType(oniumTypeIn),
48  isHiddenValley(isHiddenValleyIn), colvType(colvTypeIn), MEmix(MEmixIn),
49  MEorder (MEorderIn), MEsplit(MEsplitIn), MEgluinoRec(MEgluinoRecIn),
50  isFlexible(isFlexibleIn), hasJunction(hasJunctionIn),
51  flavour(), iAunt(), mRad(), m2Rad(), mRec(), m2Rec(), mDip(), m2Dip(),
52  m2DipCorr(), pT2(), m2(), z(), mFlavour(), asymPol(), flexFactor(),
53  pAccept() {}
54 
55  // Basic properties related to dipole and matrix element corrections.
56  int iRadiator, iRecoiler;
57  double pTmax;
58  int colType, chgType, gamType, weakType, isrType, system, systemRec,
59  MEtype, iMEpartner, weakPol, oniumType;
60  bool isHiddenValley;
61  int colvType;
62  double MEmix;
63  bool MEorder, MEsplit, MEgluinoRec, isFlexible;
64  bool hasJunction;
65 
66 
67  // Properties specific to current trial emission.
68  int flavour, iAunt;
69  double mRad, m2Rad, mRec, m2Rec, mDip, m2Dip, m2DipCorr,
70  pT2, m2, z, mFlavour, asymPol, flexFactor, pAccept;
71  double m2A{0}, m2B{0}, m2C{0}, m2gg{-1};
72 
73  // Pointer to an onium emission object, if present.
74  shared_ptr<SplitOnia> emissionPtr{};
75 
76 };
77 
78 //==========================================================================
79 
80 // The SimpleTimeShower class does timelike showers.
81 
82 class SimpleTimeShower : public TimeShower {
83 
84 public:
85 
86  // Constructor.
87  SimpleTimeShower() : hasWeaklyRadiated(), iSysSel(), pTmaxFudge(),
88  pTLastBranch(), doQCDshower(), doQEDshowerByQ(), doQEDshowerByL(),
89  doQEDshowerByOther(), doQEDshowerByGamma(), doWeakShower(),
90  doMEcorrections(), doMEextended(), doMEafterFirst(), doPhiPolAsym(),
91  doPhiPolAsymHard(), doInterleave(), doInterleaveResDec(),
92  allowBeamRecoil(), dampenBeamRecoil(), useFixedFacScale(),
93  allowRescatter(), canVetoEmission(), doHVshower(), brokenHVsym(),
94  setLambdaHV(), globalRecoil(), useLocalRecoilNow(), doSecondHard(),
95  hasUserHooks(), singleWeakEmission(), alphaSuseCMW(), vetoWeakJets(),
96  allowMPIdipole(), weakExternal(), recoilDeadCone(), doDipoleRecoil(),
97  doPartonVertex(), pTmaxMatch(), pTdampMatch(), alphaSorder(),
98  alphaSnfmax(), nGluonToQuark(), weightGluonToQuark(), recoilStrategyRF(),
99  alphaEMorder(), nGammaToQuark(), nGammaToLepton(), nCHV(), nFlavHV(),
100  idHV(), alphaHVorder(), nMaxGlobalRecoil(), weakMode(), pTdampFudge(),
101  mc(), mb(), m2c(), m2b(), renormMultFac(), factorMultFac(),
102  fixedFacScale2(), alphaSvalue(), alphaS2pi(), Lambda3flav(),
103  Lambda4flav(), Lambda5flav(), Lambda3flav2(), Lambda4flav2(),
104  Lambda5flav2(), scaleGluonToQuark(), extraGluonToQuark(), pTcolCutMin(),
105  pTcolCut(), pT2colCut(), pTchgQCut(), pT2chgQCut(), pTchgLCut(),
106  pT2chgLCut(), pTweakCut(), pT2weakCut(), mMaxGamma(), m2MaxGamma(),
107  mZ(), gammaZ(), thetaWRat(),
108  mW(), gammaW(), CFHV(), alphaHVfix(), alphaHVref(), LambdaHV(),
109  pThvCut(), pT2hvCut(), mHV(), pTmaxFudgeMPI(), weakEnhancement(),
110  vetoWeakDeltaR2(), twoHard(), dopTlimit1(), dopTlimit2(), dopTdamp(),
111  pT2damp(), kRad(), kEmt(), pdfScale2(), doTrialNow(), canEnhanceEmission(),
112  canEnhanceTrial(), canEnhanceET(), doUncertaintiesNow(), dipSel(),
113  iDipSel(), nHard(), nFinalBorn(), nMaxGlobalBranch(), nGlobal(),
114  globalRecoilMode(), limitMUQ(), weakHardSize() { beamOffset = 0;
115  pdfMode = 0; useSystems = true; }
116 
117  // Destructor.
118  virtual ~SimpleTimeShower() override {}
119 
120  // Initialize alphaStrong and related pTmin parameters.
121  virtual void init( BeamParticle* beamAPtrIn = 0,
122  BeamParticle* beamBPtrIn = 0) override;
123 
124  // Find whether to limit maximum scale of emissions, and whether to dampen.
125  virtual bool limitPTmax( Event& event, double Q2Fac = 0.,
126  double Q2Ren = 0.) override;
127 
128  // Top-level routine to do a full time-like shower in resonance decay.
129  virtual int shower( int iBeg, int iEnd, Event& event, double pTmax,
130  int nBranchMax = 0) override;
131 
132  // Top-level routine for QED radiation in hadronic decay to two leptons.
133  virtual int showerQED( int i1, int i2, Event& event, double pTmax = -1.)
134  override;
135 
136  // Prepare process-level event for shower + interleaved resonance decays.
137  // Usage: prepareProcess( process, event, iPos).
138  // iPos provides mapping from process to event entries (before showering).
139  virtual void prepareProcess( Event& process, Event&, vector<int>&) override;
140 
141  // Global recoil: reset counters and store locations of outgoing partons.
142  virtual void prepareGlobal( Event& event) override;
143 
144  // Prepare system for evolution after each new interaction; identify ME.
145  virtual void prepare( int iSys, Event& event, bool limitPTmaxIn = true)
146  override;
147 
148  // Update dipole list after a multiparton interactions rescattering.
149  virtual void rescatterUpdate( int iSys, Event& event) override;
150 
151  // Update dipole list after each ISR emission.
152  virtual void update( int iSys, Event& event, bool hasWeakRad = false)
153  override;
154 
155  // Select next pT for FSR in downwards evolution.
156  virtual double pTnext( Event& event, double pTbegAll, double pTendAll,
157  bool isFirstTrial = false, bool doTrialIn = false) override;
158 
159  // Select next pT for interleaved resonance decays.
160  virtual double pTnextResDec() override;
161 
162  // ME corrections and kinematics that may give failure.
163  virtual bool branch( Event& event, bool isInterleaved = false) override;
164 
165  // Do a resonance decay + resonance shower (including any nested decays).
166  // May be called recursively for nested decays.
167  // Usage: resonanceShower( process, event, iPos, pTmerge), where iPos
168  // maps process to event entries, and pTmerge is the scale at which this
169  // system should be merged into its parent system.
170  virtual bool resonanceShower(Event& process, Event& event,
171  vector<int>& iPos, double qRestart) override;
172 
173  // Print dipole list; for debug mainly.
174  virtual void list() const override;
175 
176  // Initialize data members for calculation of uncertainty bands.
177  virtual bool initUncertainties() override;
178 
179  // Initialize data members for application of enhancements.
180  virtual bool initEnhancements() override;
181 
182  // Tell whether FSR has done a weak emission.
183  virtual bool getHasWeaklyRadiated() override {return hasWeaklyRadiated;}
184 
185  // Tell which system was the last processed one.
186  virtual int system() const override {return iSysSel;}
187 
188  // Potential enhancement factor of pTmax scale for hardest emission.
189  virtual double enhancePTmax() override {return pTmaxFudge;}
190 
191  // Provide the pT scale of the last branching in the above shower.
192  virtual double pTLastInShower() override {return pTLastBranch;}
193 
194  // Functions to directly extract the probability of no emission between two
195  // scales. These functions are not used in the Pythia core code, but can be
196  // used by external programs to interface with the shower directly.
197  double noEmissionProbability( double pTbegAll, double pTendAll, double m2dip,
198  int id, int type, double s = -1., double x = -1.) override;
199  double pTnext( vector<TimeDipoleEnd> dipEnds, Event event, double pTbegAll,
200  double pTendAll, double m2dip, int id, int type, double s = -1.,
201  double x = -1.);
202  int pdfMode;
203  bool useSystems;
204 
205 private:
206 
207  // Constants: could only be changed in the code itself.
208  static const double MCMIN, MBMIN, SIMPLIFYROOT, XMARGIN, XMARGINCOMB,
209  WTRATIOMAX,TINYPDF, LARGEM2, THRESHM2, LAMBDA3MARGIN1ORD,
210  LAMBDA3MARGIN2ORD,WEAKPSWEIGHT, WG2QEXTRA, REJECTFACTOR,
211  PROBLIMIT;
212  static const int NLOOPMAX;
213  // Rescatter: try to fix up recoil between systems
214  static const bool FIXRESCATTER, VETONEGENERGY;
215  static const double MAXVIRTUALITYFRACTION, MAXNEGENERGYFRACTION;
216 
217  // Store properties to be returned by methods.
218  bool hasWeaklyRadiated;
219  int iSysSel;
220  double pTmaxFudge, pTLastBranch;
221 
222  // Initialization data, normally only set once.
223  bool doQCDshower, doQEDshowerByQ, doQEDshowerByL, doQEDshowerByOther,
224  doQEDshowerByGamma, doWeakShower, doMEcorrections, doMEextended,
225  doMEafterFirst, doPhiPolAsym, doPhiPolAsymHard, doInterleave,
226  doInterleaveResDec, allowBeamRecoil, dampenBeamRecoil,
227  useFixedFacScale, allowRescatter, canVetoEmission, doHVshower,
228  brokenHVsym, setLambdaHV, globalRecoil, useLocalRecoilNow,
229  doSecondHard, hasUserHooks, singleWeakEmission, alphaSuseCMW,
230  vetoWeakJets, allowMPIdipole, weakExternal, recoilDeadCone,
231  doDipoleRecoil, doPartonVertex;
232  int pdfModeSave;
233  int pTmaxMatch, pTdampMatch, alphaSorder, alphaSnfmax, nGluonToQuark,
234  weightGluonToQuark, recoilStrategyRF, alphaEMorder, nGammaToQuark,
235  nGammaToLepton, nCHV, nFlavHV, idHV, alphaHVorder, nMaxGlobalRecoil,
236  weakMode;
237  double pTdampFudge, mc, mb, m2c, m2b, renormMultFac, factorMultFac,
238  fixedFacScale2, alphaSvalue, alphaS2pi, Lambda3flav, Lambda4flav,
239  Lambda5flav, Lambda3flav2, Lambda4flav2, Lambda5flav2,
240  scaleGluonToQuark, extraGluonToQuark, pTcolCutMin, pTcolCut,
241  pT2colCut, pTchgQCut, pT2chgQCut, pTchgLCut, pT2chgLCut,
242  pTweakCut, pT2weakCut, mMaxGamma, m2MaxGamma,
243  mZ, gammaZ, thetaWRat, mW, gammaW, CFHV,
244  alphaHVfix, alphaHVref, LambdaHV, pThvCut, pT2hvCut, mHV,
245  pTmaxFudgeMPI, weakEnhancement, vetoWeakDeltaR2;
246 
247  // alphaStrong, alphaEM and alpha_HV calculations.
248  AlphaStrong alphaS;
249  AlphaEM alphaEM;
250  AlphaSUN alphaHV;
251 
252  // Weak matrix elements used for corrections both of ISR and FSR.
253  SimpleWeakShowerMEs simpleWeakShowerMEs;
254 
255  // Some current values.
256  bool twoHard, dopTlimit1, dopTlimit2, dopTdamp;
257  double pT2damp, kRad, kEmt, pdfScale2;
258 
259  // Bookkeeping of enhanced actual or trial emissions (see EPJC (2013) 73).
260  bool doTrialNow, canEnhanceEmission, canEnhanceTrial, canEnhanceET,
261  doUncertaintiesNow;
262  string splittingNameNow, splittingNameSel;
263  map< double, pair<string,double> > enhanceFactors;
264  void storeEnhanceFactor(double pT2, string name, double enhanceFactorIn)
265  { enhanceFactors.insert(make_pair(pT2,make_pair(name,enhanceFactorIn)));}
266 
267  // All dipole ends and a pointer to the selected hardest dipole end.
268  vector<TimeDipoleEnd> dipEnd;
269  TimeDipoleEnd* dipSel;
270  int iDipSel;
271 
272  // Setup a dipole end, either QCD, QED/photon, weak or Hidden Valley one.
273  void setupQCDdip( int iSys, int i, int colTag, int colSign, Event& event,
274  bool isOctetOnium = false, bool limitPTmaxIn = true);
275  void setupQEDdip( int iSys, int i, int chgType, int gamType, Event& event,
276  bool limitPTmaxIn = true);
277  void setupWeakdip( int iSys, int i,int weakType, Event& event,
278  bool limitPTmaxIn = true);
279  // Special setup for weak dipoles if already specified in info ptr.
280  void setupWeakdipExternal(Event& event, bool limitPTmaxIn = true);
281  void setupHVdip( int iSys, int i, int colvType, Event& event,
282  bool limitPTmaxIn = true);
283 
284  // Apply ME corrections for a specific subsystem.
285  bool applyMECorrections(const Event& event, TimeDipoleEnd* dipBranch,
286  int iSysBranch);
287 
288  // Special setup for onium.
289  void regenerateOniumDipoles(Event & event);
290 
291  // Evolve a QCD dipole end.
292  void pT2nextQCD( double pT2begDip, double pT2sel, TimeDipoleEnd& dip,
293  Event& event);
294 
295  // Evolve a QED dipole end, either charged or photon.
296  void pT2nextQED( double pT2begDip, double pT2sel, TimeDipoleEnd& dip,
297  Event& event);
298 
299  // Evolve a weak dipole end.
300  void pT2nextWeak( double pT2begDip, double pT2sel, TimeDipoleEnd& dip,
301  Event& event);
302 
303  // Evolve a Hidden Valley dipole end.
304  void pT2nextHV( double pT2begDip, double pT2sel, TimeDipoleEnd& dip,
305  Event& );
306 
307  // Evolve an onium dipole end.
308  void pT2nextOnium( double pT2begDip, double pT2sel, TimeDipoleEnd& dip,
309  Event& );
310 
311  // Find kind of QCD ME correction.
312  void findMEtype( Event& event, TimeDipoleEnd& dip);
313 
314  // Find type of particle; used by findMEtype.
315  int findMEparticle( int id, bool isHiddenColour = false);
316 
317  // Find mixture of V and A in gamma/Z: energy- and flavour-dependent.
318  double gammaZmix( Event& event, int iRes, int iDau1, int iDau2);
319 
320  // Set up to calculate QCD ME correction with calcMEcorr.
321  double findMEcorr(TimeDipoleEnd* dip, Particle& rad, Particle& partner,
322  Particle& emt, bool cutEdge = true);
323 
324  // Set up to calculate weak ME corrections for t-channel processes.
325  double findMEcorrWeak(TimeDipoleEnd* dip, Vec4 rad, Vec4 rec,
326  Vec4 emt, Vec4 p3, Vec4 p4, Vec4 radBef, Vec4 recBef);
327 
328  // Calculate value of QCD ME correction.
329  double calcMEcorr( int kind, int combiIn, double mixIn, double x1,
330  double x2, double r1, double r2, double r3 = 0., bool cutEdge = true);
331 
332  // Find coefficient of azimuthal asymmetry from gluon polarization.
333  void findAsymPol( Event& event, TimeDipoleEnd* dip);
334 
335  // Compute scale for interleaved resonance decays
336  double calcPTresDec(Particle& res);
337 
338  // Rescatter: propagate dipole recoil to internal lines connecting systems.
339  bool rescatterPropagateRecoil( Event& event, Vec4& pNew);
340 
341  // Properties stored for (some) global recoil schemes.
342  // Vectors of event indices defining the hard process.
343  vector<int> hardPartons;
344  // Number of partons in current hard event, number of partons in Born-type
345  // hard event (to distinguish between S and H), maximally allowed number of
346  // global recoil branchings.
347  int nHard, nFinalBorn, nMaxGlobalBranch;
348  // Number of proposed splittings in hard scattering systems.
349  map<int,int> nProposed;
350  // Number of splittings with global recoil (currently only 1).
351  int nGlobal, globalRecoilMode;
352  // Switch to constrain recoiling system.
353  bool limitMUQ;
354 
355  // Calculate uncertainty-band weights for accepted/rejected trial branching.
356  void calcUncertainties( bool accept, double pAccept,
357  double enhance, TimeDipoleEnd* dip, Particle* radPtr,
358  Particle* emtPtr, Particle* recPtr );
359 
360  // 2 -> 2 information needed for the external weak setup.
361  vector<Vec4> weakMomenta;
362  vector<int> weak2to2lines;
363  int weakHardSize;
364 
365  // Settings and member variables for interleaved resonance decays.
366  bool doFSRinResonances{};
367  int resDecScaleChoice{-1}, iHardResDecSav{}, nRecurseResDec{};
368  vector<int> idResDecSav;
369  vector<double> pTresDecSav;
370 
371  // Settings and containers for control of MECs.
372  bool skipFirstMECinHardProc;
373  vector<int> skipFirstMECinResDecIDs{};
374 
375  // Onium emissions.
376  bool doOniumShower{false};
377  vector<SplitOniaPtr> oniumEmissions;
378  set<double> oniumThresholds;
379 
380 };
381 
382 //==========================================================================
383 
384 } // end namespace Pythia8
385 
386 #endif // Pythia8_SimpleTimeShower_H
SimpleTimeShower()
Constructor.
Definition: SimpleTimeShower.h:87
Definition: StandardModel.h:23
virtual double enhancePTmax() override
Potential enhancement factor of pTmax scale for hardest emission.
Definition: SimpleTimeShower.h:189
The Event class holds all info on the generated event.
Definition: Event.h:408
The TimeShower class does timelike showers.
Definition: TimeShower.h:33
Definition: BeamParticle.h:133
TimeDipoleEnd()
Constructors.
Definition: SimpleTimeShower.h:30
int flavour
Properties specific to current trial emission.
Definition: SimpleTimeShower.h:68
Definition: StandardModel.h:106
shared_ptr< SplitOnia > emissionPtr
Pointer to an onium emission object, if present.
Definition: SimpleTimeShower.h:74
Definition: StandardModel.h:226
virtual bool getHasWeaklyRadiated() override
Tell whether FSR has done a weak emission.
Definition: SimpleTimeShower.h:183
Data on radiating dipole ends; only used inside SimpleTimeShower class.
Definition: SimpleTimeShower.h:25
virtual double pTLastInShower() override
Provide the pT scale of the last branching in the above shower.
Definition: SimpleTimeShower.h:192
Definition: Event.h:32
virtual ~SimpleTimeShower() override
Destructor.
Definition: SimpleTimeShower.h:118
Definition: SimpleWeakShowerMEs.h:26
virtual int system() const override
Tell which system was the last processed one.
Definition: SimpleTimeShower.h:186
int iRadiator
Basic properties related to dipole and matrix element corrections.
Definition: SimpleTimeShower.h:56
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The SimpleTimeShower class does timelike showers.
Definition: SimpleTimeShower.h:82
Definition: Basics.h:32