PYTHIA  8.313
UserHooks.h
1 // UserHooks.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 to allow user access to program at different stages.
7 // UserHooks: almost empty base class, with user to write the rela code.
8 // MyUserHooks: derived class, only intended as an example.
9 
10 #ifndef Pythia8_UserHooks_H
11 #define Pythia8_UserHooks_H
12 
13 #include "Pythia8/Event.h"
14 #include "Pythia8/PartonSystems.h"
15 #include "Pythia8/PhysicsBase.h"
16 #include "Pythia8/PythiaStdlib.h"
17 #include "Pythia8/SigmaProcess.h"
18 
19 namespace Pythia8 {
20 
21 //==========================================================================
22 
23 // Forward references.
24 class PhaseSpace;
25 class StringEnd;
26 class HadronLevel;
27 
28 //==========================================================================
29 
30 // UserHooks is base class for user access to program execution.
31 
32 class UserHooks : public PhysicsBase {
33 
34 public:
35 
36  // Destructor.
37  virtual ~UserHooks() {}
38 
39  // Initialisation after beams have been set by Pythia::init().
40  virtual bool initAfterBeams() { return true; }
41 
42  // Possibility to modify cross section of process.
43  virtual bool canModifySigma() {return false;}
44 
45  // Multiplicative factor modifying the cross section of a hard process.
46  virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
47  const PhaseSpace* phaseSpacePtr, bool inEvent);
48 
49  // Possibility to bias selection of events, compensated by a weight.
50  virtual bool canBiasSelection() {return false;}
51 
52  // Multiplicative factor in the phase space selection of a hard process.
53  virtual double biasSelectionBy(const SigmaProcess* sigmaProcessPtr,
54  const PhaseSpace* phaseSpacePtr, bool inEvent);
55 
56  // Event weight to compensate for selection weight above.
57  virtual double biasedSelectionWeight() {return 1./selBias;}
58 
59  // Possibility to veto event after process-level selection.
60  virtual bool canVetoProcessLevel() {return false;}
61 
62  // Decide whether to veto current process or not, based on process record.
63  // Usage: doVetoProcessLevel( process).
64  virtual bool doVetoProcessLevel(Event& ) {return false;}
65 
66  // Possibility to set low energy cross sections.
67  // Usage: canSetLowEnergySigma(idA, idB).
68  virtual bool canSetLowEnergySigma(int, int) const {return false;}
69 
70  // Set low energy cross sections for ids where canSetLowEnergySigma is true.
71  // Usage: doSetLowEnergySigma(idA, idB, eCM, mA, mB).
72  virtual double doSetLowEnergySigma(int, int, double, double, double) const {
73  return 0.;}
74 
75  // Possibility to veto resonance decay chain.
76  virtual bool canVetoResonanceDecays() {return false;}
77 
78  // Decide whether to veto current resonance decay chain or not, based on
79  // process record. Usage: doVetoProcessLevel( process).
80  virtual bool doVetoResonanceDecays(Event& ) {return false;}
81 
82  // Possibility to veto MPI + ISR + FSR evolution and kill event,
83  // making decision at a fixed pT scale. Useful for MLM-style matching.
84  virtual bool canVetoPT() {return false;}
85 
86  // Transverse-momentum scale for veto test.
87  virtual double scaleVetoPT() {return 0.;}
88 
89  // Decide whether to veto current event or not, based on event record.
90  // Usage: doVetoPT( iPos, event), where iPos = 0: no emissions so far;
91  // iPos = 1/2/3 joint evolution, latest step was MPI/ISR/FSR;
92  // iPos = 4: FSR only afterwards; iPos = 5: FSR in resonance decay.
93  virtual bool doVetoPT( int , const Event& ) {return false;}
94 
95  // Possibility to veto MPI + ISR + FSR evolution and kill event,
96  // making decision after fixed number of ISR or FSR steps.
97  virtual bool canVetoStep() {return false;}
98 
99  // Up to how many ISR + FSR steps of hardest interaction should be checked.
100  virtual int numberVetoStep() {return 1;}
101 
102  // Decide whether to veto current event or not, based on event record.
103  // Usage: doVetoStep( iPos, nISR, nFSR, event), where iPos as above,
104  // nISR and nFSR number of emissions so far for hard interaction only.
105  virtual bool doVetoStep( int , int , int , const Event& ) {return false;}
106 
107  // Possibility to veto MPI + ISR + FSR evolution and kill event,
108  // making decision after fixed number of MPI steps.
109  virtual bool canVetoMPIStep() {return false;}
110 
111  // Up to how many MPI steps should be checked.
112  virtual int numberVetoMPIStep() {return 1;}
113 
114  // Decide whether to veto current event or not, based on event record.
115  // Usage: doVetoMPIStep( nMPI, event), where nMPI is number of MPI's so far.
116  virtual bool doVetoMPIStep( int , const Event& ) {return false;}
117 
118  // Possibility to veto event after ISR + FSR + MPI in parton level,
119  // but before beam remnants and resonance decays.
120  virtual bool canVetoPartonLevelEarly() {return false;}
121 
122  // Decide whether to veto current partons or not, based on event record.
123  // Usage: doVetoPartonLevelEarly( event).
124  virtual bool doVetoPartonLevelEarly( const Event& ) {return false;}
125 
126  // Retry same ProcessLevel with a new PartonLevel after a veto in
127  // doVetoPT, doVetoStep, doVetoMPIStep or doVetoPartonLevelEarly
128  // if you overload this method to return true.
129  virtual bool retryPartonLevel() {return false;}
130 
131  // Possibility to veto event after parton-level selection.
132  virtual bool canVetoPartonLevel() {return false;}
133 
134  // Decide whether to veto current partons or not, based on event record.
135  // Usage: doVetoPartonLevel( event).
136  virtual bool doVetoPartonLevel( const Event& ) {return false;}
137 
138  // Possibility to set initial scale in TimeShower for resonance decay.
139  virtual bool canSetResonanceScale() {return false;}
140 
141  // Initial scale for TimeShower evolution.
142  // Usage: scaleResonance( iRes, event), where iRes is location
143  // of decaying resonance in the event record.
144  virtual double scaleResonance( int, const Event& ) {return 0.;}
145 
146  // Possibility to veto an emission in the ISR machinery.
147  virtual bool canVetoISREmission() {return false;}
148 
149  // Decide whether to veto current emission or not, based on event record.
150  // Usage: doVetoISREmission( sizeOld, event, iSys) where sizeOld is size
151  // of event record before current emission-to-be-scrutinized was added,
152  // and iSys is the system of the radiation (according to PartonSystems).
153  virtual bool doVetoISREmission( int, const Event&, int ) {return false;}
154 
155  // Possibility to veto an emission in the FSR machinery.
156  virtual bool canVetoFSREmission() {return false;}
157 
158  // Decide whether to veto current emission or not, based on event record.
159  // Usage: doVetoFSREmission( sizeOld, event, iSys, inResonance) where
160  // sizeOld is size of event record before current emission-to-be-scrutinized
161  // was added, iSys is the system of the radiation (according to
162  // PartonSystems), and inResonance is true if the emission takes place in a
163  // resonance decay.
164  virtual bool doVetoFSREmission( int, const Event&, int, bool = false )
165  {return false;}
166 
167  // Possibility to veto an MPI.
168  virtual bool canVetoMPIEmission() { return false; }
169 
170  // Decide whether to veto an MPI based on event record.
171  // Usage: doVetoMPIEmission( sizeOld, event) where sizeOld
172  // is size of event record before the current MPI.
173  virtual bool doVetoMPIEmission( int, const Event &) { return false; }
174 
175  // Possibility to reconnect colours from resonance decay systems.
176  virtual bool canReconnectResonanceSystems() { return false; }
177 
178  // Do reconnect colours from resonance decay systems.
179  // Usage: doVetoFSREmission( oldSizeEvt, event)
180  // where oldSizeEvent is the event size before resonance decays.
181  // Should normally return true, while false means serious failure.
182  // Value of PartonLevel:earlyResDec determines where method is called.
183  virtual bool doReconnectResonanceSystems( int, Event &) {return true;}
184 
185  // Can change fragmentation parameters.
186  virtual bool canChangeFragPar() { return false;}
187 
188  // Set initial ends of a string to be fragmented. This is done once
189  // for each string. Note that the second string end may be zero in case
190  // we are hadronising a string piece leading to a junction.
191  virtual void setStringEnds( const StringEnd*, const StringEnd*,
192  vector<int>) {}
193 
194  // Do change fragmentation parameters.
195  // Input: flavPtr, zPtr, pTPtr, idEnd, m2Had, iParton and posEnd (or
196  // negEnd).
197  virtual bool doChangeFragPar( StringFlav*, StringZ*, StringPT*, int,
198  double, vector<int>, const StringEnd* ) { return false;}
199 
200  // Can veto hadrons in the fragmentation.
201  virtual bool canVetoFragmentation() { return false;}
202 
203  // Do a veto on a hadron just before it is added to the final state.
204  // The StringEnd from which the the hadron was produced is included
205  // for information.
206  virtual bool doVetoFragmentation( Particle, const StringEnd *)
207  { return false;}
208 
209  // Do a veto on a hadron just before it is added to the final state
210  // (final two hadron case). Note that a veto from this function
211  // will refragment the whole string. If only the regeneration of
212  // the final two is desired, the doVetoFinalTwo function should be
213  // used instead.
215  const StringEnd*, const StringEnd* ) { return false;}
216 
217  // Do a veto on a hadron just before it is added to the final state
218  // (final two hadron case). Note that a veto from this function will
219  // regenerate the final two. If the refragmentation of the whole
220  // string is desired, the doVetoFragmentation function should be
221  // used instead.
223  const StringEnd*, const StringEnd* ) { return false;}
224 
225  // Possibility to veto an event after hadronization based
226  // on event contents. Works as an early trigger to avoid
227  // running the time consuming rescattering process on
228  // uninteresting events.
229  virtual bool canVetoAfterHadronization() {return false;}
230 
231  // Do the actual veto after hadronization.
232  virtual bool doVetoAfterHadronization(const Event& ) {return false;}
233 
234  // Can set the overall impact parameter for the MPI treatment.
235  virtual bool canSetImpactParameter() const { return false; }
236 
237  // Set the overall impact parameter for the MPI treatment.
238  virtual double doSetImpactParameter() { return 0.0; }
239 
240  // Custom processing at the end of HadronLevel::next.
241  virtual bool onEndHadronLevel(HadronLevel&, Event&) { return true; }
242 
243 protected:
244 
245  // Constructor.
247 
248  // After initInfoPtr, initialize workEvent
249  virtual void onInitInfoPtr() override {
250  // Set smart pointer to null, in order to avoid circular dependency
251  userHooksPtr = nullptr;
252  workEvent.init("(work event)", particleDataPtr);
253  }
254 
255  // omitResonanceDecays omits resonance decay chains from process record.
256  void omitResonanceDecays(const Event& process, bool finalOnly = false);
257 
258  // subEvent extracts currently resolved partons in the hard process.
259  void subEvent(const Event& event, bool isHardest = true);
260 
261  // Have one event object around as work area.
263 
264  // User-imposed selection bias.
265  double selBias = 1.;
266 
267  // Bookkept quantities for boosted event weights.
268  double enhancedEventWeight = {}, pTEnhanced = {}, wtEnhanced = {};
269 
270 };
271 
272 //==========================================================================
273 
274 // SuppressSmallPT is a derived class for user access to program execution.
275 // It is a simple example, illustrating how to suppress the cross section
276 // of 2 -> 2 processes by a factor pT^4 / (pT0^2 + pT^2)^2, with pT0 input,
277 // and also modify alpha_strong scale similarly.
278 
279 class SuppressSmallPT : public UserHooks {
280 
281 public:
282 
283  // Constructor.
284  SuppressSmallPT( double pT0timesMPIIn = 1., int numberAlphaSIn = 0,
285  bool useSameAlphaSasMPIIn = true) : pT20(0.) {isInit = false;
286  pT0timesMPI = pT0timesMPIIn; numberAlphaS = numberAlphaSIn;
287  useSameAlphaSasMPI = useSameAlphaSasMPIIn;}
288 
289  // Possibility to modify cross section of process.
290  virtual bool canModifySigma() {return true;}
291 
292  // Multiplicative factor modifying the cross section of a hard process.
293  // Usage: inEvent is true for event generation, false for initialization.
294  virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
295  const PhaseSpace* phaseSpacePtr, bool );
296 
297 private:
298 
299  // Save input properties and the squared pT0 scale.
300  bool isInit, useSameAlphaSasMPI;
301  int numberAlphaS;
302  double pT0timesMPI, pT20;
303 
304  // Alpha_strong calculation.
305  AlphaStrong alphaS;
306 
307 };
308 
309 //==========================================================================
310 
311 // UserHooksVector implements a vector of UserHooks and is itself a UserHooks.
312 
313 class UserHooksVector: public UserHooks {
314 
315 public:
316 
317  // The default constructor is private, and should only be used
318  // internally in Pythia.
320  friend class Pythia;
321 
322  // Destructor.
323  virtual ~UserHooksVector() {}
324 
325  // Initialisation after beams have been set by Pythia::init().
326  // Check that there are no (obvious) clashes.
327  virtual bool initAfterBeams() {
328  int nCanSetResonanceScale = 0;
329  int nCanChangeFragPar = 0;
330  int nCanSetImpactParameter = 0;
331  for ( int i = 0, N = hooks.size(); i < N; ++i ) {
332  registerSubObject(*hooks[i]);
333  if ( !hooks[i]->initAfterBeams() ) return false;
334  if (hooks[i]->canSetResonanceScale()) ++nCanSetResonanceScale;
335  if (hooks[i]->canChangeFragPar()) ++nCanChangeFragPar;
336  if (hooks[i]->canSetImpactParameter()) ++nCanSetImpactParameter;
337  }
338  if (nCanSetResonanceScale > 1) {
339  loggerPtr->ERROR_MSG(
340  "multiple UserHooks with canSetResonanceScale() not allowed");
341  return false;
342  }
343  if (nCanChangeFragPar > 1) {
344  loggerPtr->ERROR_MSG(
345  "multiple UserHooks with canChangeFragPar() not allowed");
346  return false;
347  }
348  if (nCanSetImpactParameter > 1) {
349  loggerPtr->ERROR_MSG(
350  "multiple UserHooks with canSetImpactParameter() not allowed");
351  return false;
352  }
353  return true;
354  }
355 
356  // Possibility to modify cross section of process.
357  virtual bool canModifySigma() {
358  for ( int i = 0, N = hooks.size(); i < N; ++i )
359  if ( hooks[i]->canModifySigma() ) return true;
360  return false;
361  }
362 
363  // Multiplicative factor modifying the cross section of a hard process.
364  virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
365  const PhaseSpace* phaseSpacePtr, bool inEvent) {
366  double f = 1.0;
367  for ( int i = 0, N = hooks.size(); i < N; ++i )
368  if ( hooks[i]->canModifySigma() )
369  f *= hooks[i]->multiplySigmaBy(sigmaProcessPtr, phaseSpacePtr,inEvent);
370  return f;
371  }
372 
373  // Possibility to bias selection of events, compensated by a weight.
374  virtual bool canBiasSelection() {
375  for ( int i = 0, N = hooks.size(); i < N; ++i )
376  if ( hooks[i]->canBiasSelection() ) return true;
377  return false;
378  }
379 
380  // Multiplicative factor in the phase space selection of a hard process.
381  virtual double biasSelectionBy(const SigmaProcess* sigmaProcessPtr,
382  const PhaseSpace* phaseSpacePtr, bool inEvent) {
383  double f = 1.0;
384  for ( int i = 0, N = hooks.size(); i < N; ++i )
385  if ( hooks[i]->canBiasSelection() )
386  f *= hooks[i]->biasSelectionBy(sigmaProcessPtr, phaseSpacePtr,
387  inEvent);
388  return f;
389  }
390 
391  // Event weight to compensate for selection weight above.
392  virtual double biasedSelectionWeight() {
393  double f = 1.0;
394  for ( int i = 0, N = hooks.size(); i < N; ++i )
395  if ( hooks[i]->canBiasSelection() )
396  f *= hooks[i]->biasedSelectionWeight();
397  return f;
398  }
399 
400  // Possibility to veto event after process-level selection.
401  virtual bool canVetoProcessLevel() {
402  for ( int i = 0, N = hooks.size(); i < N; ++i )
403  if ( hooks[i]->canVetoProcessLevel() ) return true;
404  return false;
405  }
406 
407  // Decide whether to veto current process or not, based on process record.
408  // Usage: doVetoProcessLevel( process).
409  virtual bool doVetoProcessLevel(Event& e) {
410  for ( int i = 0, N = hooks.size(); i < N; ++i )
411  if ( hooks[i]->canVetoProcessLevel() &&
412  hooks[i]->doVetoProcessLevel(e) ) return true;
413  return false;
414  }
415 
416  // Possibility to veto resonance decay chain.
417  virtual bool canVetoResonanceDecays() {
418  for ( int i = 0, N = hooks.size(); i < N; ++i )
419  if ( hooks[i]->canVetoResonanceDecays() ) return true;
420  return false;
421  }
422 
423  // Decide whether to veto current resonance decay chain or not, based on
424  // process record. Usage: doVetoProcessLevel( process).
425  virtual bool doVetoResonanceDecays(Event& e) {
426  for ( int i = 0, N = hooks.size(); i < N; ++i )
427  if ( hooks[i]->canVetoResonanceDecays() &&
428  hooks[i]->doVetoResonanceDecays(e) ) return true;
429  return false;
430  }
431 
432  // Possibility to veto MPI + ISR + FSR evolution and kill event,
433  // making decision at a fixed pT scale. Useful for MLM-style matching.
434  virtual bool canVetoPT() {
435  for ( int i = 0, N = hooks.size(); i < N; ++i )
436  if ( hooks[i]->canVetoPT() ) return true;
437  return false;
438  }
439 
440  // Transverse-momentum scale for veto test.
441  virtual double scaleVetoPT() {
442  double s = 0.0;
443  for ( int i = 0, N = hooks.size(); i < N; ++i )
444  if ( hooks[i]->canVetoPT() ) s = max(s, hooks[i]->scaleVetoPT());
445  return s;
446  }
447 
448  // Decide whether to veto current event or not, based on event record.
449  // Usage: doVetoPT( iPos, event), where iPos = 0: no emissions so far;
450  // iPos = 1/2/3 joint evolution, latest step was MPI/ISR/FSR;
451  // iPos = 4: FSR only afterwards; iPos = 5: FSR in resonance decay.
452  virtual bool doVetoPT( int iPos, const Event& e) {
453  for ( int i = 0, N = hooks.size(); i < N; ++i )
454  if ( hooks[i]->canVetoPT() && hooks[i]->doVetoPT(iPos, e) ) return true;
455  return false;
456  }
457 
458  // Possibility to veto MPI + ISR + FSR evolution and kill event,
459  // making decision after fixed number of ISR or FSR steps.
460  virtual bool canVetoStep() {
461  for ( int i = 0, N = hooks.size(); i < N; ++i )
462  if ( hooks[i]->canVetoStep() ) return true;
463  return false;
464  }
465 
466  // Up to how many ISR + FSR steps of hardest interaction should be checked.
467  virtual int numberVetoStep() {
468  int n = 1;
469  for ( int i = 0, N = hooks.size(); i < N; ++i )
470  if ( hooks[i]->canVetoStep() ) n = max(n, hooks[i]->numberVetoStep());
471  return n;
472  }
473 
474  // Decide whether to veto current event or not, based on event record.
475  // Usage: doVetoStep( iPos, nISR, nFSR, event), where iPos as above,
476  // nISR and nFSR number of emissions so far for hard interaction only.
477  virtual bool doVetoStep( int iPos, int nISR, int nFSR, const Event& e) {
478  for ( int i = 0, N = hooks.size(); i < N; ++i )
479  if ( hooks[i]->canVetoStep()
480  && hooks[i]->doVetoStep(iPos, nISR, nFSR, e) ) return true;
481  return false;
482  }
483 
484  // Possibility to veto MPI + ISR + FSR evolution and kill event,
485  // making decision after fixed number of MPI steps.
486  virtual bool canVetoMPIStep() {
487  for ( int i = 0, N = hooks.size(); i < N; ++i )
488  if ( hooks[i]->canVetoMPIStep() ) return true;
489  return false;
490  }
491 
492  // Up to how many MPI steps should be checked.
493  virtual int numberVetoMPIStep() {
494  int n = 1;
495  for ( int i = 0, N = hooks.size(); i < N; ++i )
496  if ( hooks[i]->canVetoMPIStep() )
497  n = max(n, hooks[i]->numberVetoMPIStep());
498  return n;
499  }
500 
501  // Decide whether to veto current event or not, based on event record.
502  // Usage: doVetoMPIStep( nMPI, event), where nMPI is number of MPI's so far.
503  virtual bool doVetoMPIStep( int nMPI, const Event& e) {
504  for ( int i = 0, N = hooks.size(); i < N; ++i )
505  if ( hooks[i]->canVetoMPIStep() && hooks[i]->doVetoMPIStep(nMPI, e) )
506  return true;
507  return false;
508  }
509 
510  // Possibility to veto event after ISR + FSR + MPI in parton level,
511  // but before beam remnants and resonance decays.
512  virtual bool canVetoPartonLevelEarly() {
513  for ( int i = 0, N = hooks.size(); i < N; ++i )
514  if ( hooks[i]->canVetoPartonLevelEarly() ) return true;
515  return false;
516  }
517 
518  // Decide whether to veto current partons or not, based on event record.
519  // Usage: doVetoPartonLevelEarly( event).
520  virtual bool doVetoPartonLevelEarly( const Event& e) {
521  for ( int i = 0, N = hooks.size(); i < N; ++i )
522  if ( hooks[i]->canVetoPartonLevelEarly()
523  && hooks[i]->doVetoPartonLevelEarly(e) ) return true;
524  return false;
525  }
526 
527  // Retry same ProcessLevel with a new PartonLevel after a veto in
528  // doVetoPT, doVetoStep, doVetoMPIStep or doVetoPartonLevelEarly
529  // if you overload this method to return true.
530  virtual bool retryPartonLevel() {
531  for ( int i = 0, N = hooks.size(); i < N; ++i )
532  if ( hooks[i]->retryPartonLevel() ) return true;
533  return false;
534  }
535 
536  // Possibility to veto event after parton-level selection.
537  virtual bool canVetoPartonLevel() {
538  for ( int i = 0, N = hooks.size(); i < N; ++i )
539  if ( hooks[i]->canVetoPartonLevel() ) return true;
540  return false;
541  }
542 
543  // Decide whether to veto current partons or not, based on event record.
544  // Usage: doVetoPartonLevel( event).
545  virtual bool doVetoPartonLevel( const Event& e) {
546  for ( int i = 0, N = hooks.size(); i < N; ++i )
547  if ( hooks[i]->canVetoPartonLevel()
548  && hooks[i]->doVetoPartonLevel(e) ) return true;
549  return false;
550  }
551 
552  // Possibility to set initial scale in TimeShower for resonance decay.
553  virtual bool canSetResonanceScale() {
554  for ( int i = 0, N = hooks.size(); i < N; ++i )
555  if ( hooks[i]->canSetResonanceScale() ) return true;
556  return false;
557  }
558 
559  // Initial scale for TimeShower evolution.
560  // Usage: scaleResonance( iRes, event), where iRes is location
561  // of decaying resonance in the event record.
562  virtual double scaleResonance( int iRes, const Event& e) {
563  double s = 0.0;
564  for ( int i = 0, N = hooks.size(); i < N; ++i )
565  if ( hooks[i]->canSetResonanceScale() )
566  s = max(s, hooks[i]->scaleResonance(iRes, e));
567  return s;
568  }
569 
570  // Possibility to veto an emission in the ISR machinery.
571  virtual bool canVetoISREmission() {
572  for ( int i = 0, N = hooks.size(); i < N; ++i )
573  if ( hooks[i]->canVetoISREmission() ) return true;
574  return false;
575  }
576 
577  // Decide whether to veto current emission or not, based on event record.
578  // Usage: doVetoISREmission( sizeOld, event, iSys) where sizeOld is size
579  // of event record before current emission-to-be-scrutinized was added,
580  // and iSys is the system of the radiation (according to PartonSystems).
581  virtual bool doVetoISREmission( int sizeOld, const Event& e, int iSys) {
582  for ( int i = 0, N = hooks.size(); i < N; ++i )
583  if ( hooks[i]->canVetoISREmission()
584  && hooks[i]->doVetoISREmission(sizeOld, e, iSys) ) return true;
585  return false;
586  }
587 
588  // Possibility to veto an emission in the FSR machinery.
589  virtual bool canVetoFSREmission() {
590  for ( int i = 0, N = hooks.size(); i < N; ++i )
591  if ( hooks[i]->canVetoFSREmission() ) return true;
592  return false;
593  }
594 
595  // Decide whether to veto current emission or not, based on event record.
596  // Usage: doVetoFSREmission( sizeOld, event, iSys, inResonance) where
597  // sizeOld is size of event record before current emission-to-be-scrutinized
598  // was added, iSys is the system of the radiation (according to
599  // PartonSystems), and inResonance is true if the emission takes place in a
600  // resonance decay.
601  virtual bool doVetoFSREmission(int sizeOld, const Event& e,
602  int iSys, bool inResonance = false ) {
603  for ( int i = 0, N = hooks.size(); i < N; ++i )
604  if ( hooks[i]->canVetoFSREmission()
605  && hooks[i]->doVetoFSREmission(sizeOld, e, iSys, inResonance) )
606  return true;
607  return false;
608  }
609 
610  // Possibility to veto an MPI.
611  virtual bool canVetoMPIEmission() {
612  for ( int i = 0, N = hooks.size(); i < N; ++i )
613  if ( hooks[i]->canVetoMPIEmission() ) return true;
614  return false;
615  }
616 
617  // Decide whether to veto an MPI based on event record.
618  // Usage: doVetoMPIEmission( sizeOld, event) where sizeOld
619  // is size of event record before the current MPI.
620  virtual bool doVetoMPIEmission( int sizeOld, const Event & e) {
621  for ( int i = 0, N = hooks.size(); i < N; ++i )
622  if ( hooks[i]->canVetoMPIEmission()
623  && hooks[i]->doVetoMPIEmission(sizeOld, e) )
624  return true;
625  return false;
626  }
627 
628  // Possibility to reconnect colours from resonance decay systems.
630  for ( int i = 0, N = hooks.size(); i < N; ++i )
631  if ( hooks[i]->canReconnectResonanceSystems() ) return true;
632  return false;
633  }
634 
635  // Do reconnect colours from resonance decay systems.
636  // Usage: doVetoFSREmission( oldSizeEvt, event)
637  // where oldSizeEvent is the event size before resonance decays.
638  // Should normally return true, while false means serious failure.
639  // Value of PartonLevel:earlyResDec determines where method is called.
640  virtual bool doReconnectResonanceSystems( int j, Event & e) {
641  for ( int i = 0, N = hooks.size(); i < N; ++i )
642  if ( hooks[i]->canReconnectResonanceSystems()
643  && hooks[i]->doReconnectResonanceSystems(j, e) ) return true;
644  return false;
645  }
646 
647  // Can change fragmentation parameters.
648  virtual bool canChangeFragPar() {
649  for ( int i = 0, N = hooks.size(); i < N; ++i )
650  if ( hooks[i]->canChangeFragPar() ) return true;
651  return false;
652  }
653 
654  // Set initial ends of a string to be fragmented. This is done once
655  // for each string. Note that the second string end may be zero in
656  // case we are hadronising a string piece leading to a junction.
657  virtual void setStringEnds( const StringEnd* pos, const StringEnd* neg,
658  vector<int> iPart) {
659  for ( int i = 0, N = hooks.size(); i < N; ++i )
660  hooks[i]->setStringEnds( pos, neg, iPart);
661  }
662 
663  // Do change fragmentation parameters.
664  // Input: flavPtr, zPtr, pTPtr, idEnd, m2Had, iParton and posEnd (or
665  // negEnd).
666  virtual bool doChangeFragPar( StringFlav* sfIn, StringZ* zIn,
667  StringPT* ptIn, int idIn, double mIn, vector<int> parIn,
668  const StringEnd* endIn) {
669  for ( int i = 0, N = hooks.size(); i < N; ++i )
670  if ( hooks[i]->canChangeFragPar()
671  && hooks[i]->doChangeFragPar(sfIn, zIn, ptIn, idIn,
672  mIn, parIn, endIn) ) return true;
673  return false;}
674 
675  // Can veto hadrons in the fragmentation.
676  virtual bool canVetoFragmentation() {
677  for ( int i = 0, N = hooks.size(); i < N; ++i )
678  if ( hooks[i]->canVetoFragmentation() ) return true;
679  return false;}
680 
681  // Do a veto on a hadron just before it is added to the final state.
682  virtual bool doVetoFragmentation(Particle p, const StringEnd* nowEnd) {
683  for ( int i = 0, N = hooks.size(); i < N; ++i )
684  if ( hooks[i]->canChangeFragPar()
685  && hooks[i]->doVetoFragmentation(p, nowEnd) ) return true;
686  return false;
687  }
688 
689  virtual bool doVetoFragmentation(Particle p1, Particle p2,
690  const StringEnd* e1, const StringEnd* e2) {
691  for ( int i = 0, N = hooks.size(); i < N; ++i )
692  if ( hooks[i]->canChangeFragPar()
693  && hooks[i]->doVetoFragmentation(p1, p2, e1, e2) ) return true;
694  return false;
695  }
696 
697  virtual bool doVetoFinalTwo(Particle p1, Particle p2,
698  const StringEnd* e1, const StringEnd* e2) {
699  for ( int i = 0, N = hooks.size(); i < N; ++i )
700  if ( hooks[i]->canChangeFragPar()
701  && hooks[i]->doVetoFinalTwo(p1, p2, e1, e2) ) return true;
702  return false;
703  }
704 
705  // Possibility to veto an event after hadronization based on event
706  // contents. Works as an early trigger to avoid running the time
707  // consuming rescattering process on uninteresting events.
708  virtual bool canVetoAfterHadronization() {
709  for ( int i = 0, N = hooks.size(); i < N; ++i )
710  if ( hooks[i]->canVetoAfterHadronization() ) return true;
711  return false;
712  }
713 
714  // Do the actual veto after hadronization.
715  virtual bool doVetoAfterHadronization(const Event& e) {
716  for ( int i = 0, N = hooks.size(); i < N; ++i )
717  if ( hooks[i]->canVetoAfterHadronization()
718  && hooks[i]->doVetoAfterHadronization(e) ) return true;
719  return false;
720  }
721 
722  // Can set the overall impact parameter for the MPI treatment.
723  virtual bool canSetImpactParameter() const {
724  for ( int i = 0, N = hooks.size(); i < N; ++i )
725  if ( hooks[i]->canSetImpactParameter() ) return true;
726  return false;
727  }
728 
729  // Set the overall impact parameter for the MPI treatment.
730  virtual double doSetImpactParameter() {
731  for ( int i = 0, N = hooks.size(); i < N; ++i )
732  if ( hooks[i]->canSetImpactParameter() )
733  return hooks[i]->doSetImpactParameter();
734  return 0.0;
735  }
736 
737  // The vector of user hooks.
738  vector< shared_ptr<UserHooks> > hooks = {};
739 
740 };
741 
742 //==========================================================================
743 
744 } // end namespace Pythia8
745 
746 #endif // Pythia8_UserHooks_H
virtual bool canBiasSelection()
Possibility to bias selection of events, compensated by a weight.
Definition: UserHooks.h:50
virtual int numberVetoStep()
Up to how many ISR + FSR steps of hardest interaction should be checked.
Definition: UserHooks.h:467
virtual double scaleVetoPT()
Transverse-momentum scale for veto test.
Definition: UserHooks.h:87
virtual bool doVetoFinalTwo(Particle p1, Particle p2, const StringEnd *e1, const StringEnd *e2)
Definition: UserHooks.h:697
virtual int numberVetoStep()
Up to how many ISR + FSR steps of hardest interaction should be checked.
Definition: UserHooks.h:100
virtual bool canVetoPT()
Definition: UserHooks.h:434
virtual double scaleVetoPT()
Transverse-momentum scale for veto test.
Definition: UserHooks.h:441
virtual bool canVetoPartonLevelEarly()
Definition: UserHooks.h:120
virtual double multiplySigmaBy(const SigmaProcess *sigmaProcessPtr, const PhaseSpace *phaseSpacePtr, bool inEvent)
Multiplicative factor modifying the cross section of a hard process.
Definition: UserHooks.cc:27
Definition: StandardModel.h:23
UserHooksVector()
Definition: UserHooks.h:319
Definition: PhysicsBase.h:27
virtual bool canModifySigma()
Possibility to modify cross section of process.
Definition: UserHooks.h:290
void registerSubObject(PhysicsBase &pb)
Register a sub object that should have its information in sync with this.
Definition: PhysicsBase.cc:56
virtual double doSetImpactParameter()
Set the overall impact parameter for the MPI treatment.
Definition: UserHooks.h:730
The Event class holds all info on the generated event.
Definition: Event.h:408
UserHooks()
Constructor.
Definition: UserHooks.h:246
virtual bool doVetoISREmission(int sizeOld, const Event &e, int iSys)
Definition: UserHooks.h:581
virtual bool retryPartonLevel()
Definition: UserHooks.h:129
virtual bool onEndHadronLevel(HadronLevel &, Event &)
Custom processing at the end of HadronLevel::next.
Definition: UserHooks.h:241
SigmaProcess is the base class for cross section calculations.
Definition: SigmaProcess.h:86
virtual bool doVetoPT(int iPos, const Event &e)
Definition: UserHooks.h:452
virtual bool doVetoPartonLevel(const Event &e)
Definition: UserHooks.h:545
virtual bool doVetoMPIStep(int nMPI, const Event &e)
Definition: UserHooks.h:503
virtual double biasedSelectionWeight()
Event weight to compensate for selection weight above.
Definition: UserHooks.h:392
virtual bool doChangeFragPar(StringFlav *, StringZ *, StringPT *, int, double, vector< int >, const StringEnd *)
Definition: UserHooks.h:197
virtual bool canVetoProcessLevel()
Possibility to veto event after process-level selection.
Definition: UserHooks.h:401
virtual void setStringEnds(const StringEnd *, const StringEnd *, vector< int >)
Definition: UserHooks.h:191
The StringPT class is used to select select transverse momenta.
Definition: FragmentationFlavZpT.h:326
virtual bool canVetoMPIEmission()
Possibility to veto an MPI.
Definition: UserHooks.h:611
virtual bool doVetoStep(int iPos, int nISR, int nFSR, const Event &e)
Definition: UserHooks.h:477
virtual bool canVetoStep()
Definition: UserHooks.h:460
virtual bool canVetoAfterHadronization()
Definition: UserHooks.h:708
virtual bool canSetResonanceScale()
Possibility to set initial scale in TimeShower for resonance decay.
Definition: UserHooks.h:139
virtual bool doVetoMPIStep(int, const Event &)
Definition: UserHooks.h:116
virtual bool doVetoFSREmission(int, const Event &, int, bool=false)
Definition: UserHooks.h:164
virtual bool doVetoPartonLevelEarly(const Event &)
Definition: UserHooks.h:124
virtual ~UserHooks()
Destructor.
Definition: UserHooks.h:37
virtual int numberVetoMPIStep()
Up to how many MPI steps should be checked.
Definition: UserHooks.h:112
virtual void setStringEnds(const StringEnd *pos, const StringEnd *neg, vector< int > iPart)
Definition: UserHooks.h:657
virtual bool doVetoFSREmission(int sizeOld, const Event &e, int iSys, bool inResonance=false)
Definition: UserHooks.h:601
virtual bool doVetoResonanceDecays(Event &)
Definition: UserHooks.h:80
virtual bool doChangeFragPar(StringFlav *sfIn, StringZ *zIn, StringPT *ptIn, int idIn, double mIn, vector< int > parIn, const StringEnd *endIn)
Definition: UserHooks.h:666
virtual bool doVetoResonanceDecays(Event &e)
Definition: UserHooks.h:425
Definition: UserHooks.h:279
virtual double biasedSelectionWeight()
Event weight to compensate for selection weight above.
Definition: UserHooks.h:57
virtual bool doVetoPT(int, const Event &)
Definition: UserHooks.h:93
ParticleData * particleDataPtr
Pointer to the particle data table.
Definition: PhysicsBase.h:84
virtual bool doVetoAfterHadronization(const Event &e)
Do the actual veto after hadronization.
Definition: UserHooks.h:715
virtual bool initAfterBeams()
Definition: UserHooks.h:327
virtual bool canVetoISREmission()
Possibility to veto an emission in the ISR machinery.
Definition: UserHooks.h:147
The StringZ class is used to sample the fragmentation function f(z).
Definition: FragmentationFlavZpT.h:265
virtual bool canVetoResonanceDecays()
Possibility to veto resonance decay chain.
Definition: UserHooks.h:417
virtual double scaleResonance(int, const Event &)
Definition: UserHooks.h:144
virtual bool canVetoPartonLevel()
Possibility to veto event after parton-level selection.
Definition: UserHooks.h:537
Definition: HadronLevel.h:44
virtual bool canVetoMPIEmission()
Possibility to veto an MPI.
Definition: UserHooks.h:168
UserHooksPtr userHooksPtr
Definition: PhysicsBase.h:124
virtual void onInitInfoPtr() override
After initInfoPtr, initialize workEvent.
Definition: UserHooks.h:249
virtual bool doVetoFragmentation(Particle, const StringEnd *)
Definition: UserHooks.h:206
virtual bool canVetoResonanceDecays()
Possibility to veto resonance decay chain.
Definition: UserHooks.h:76
virtual bool doVetoProcessLevel(Event &e)
Definition: UserHooks.h:409
void omitResonanceDecays(const Event &process, bool finalOnly=false)
omitResonanceDecays omits resonance decay chains from process record.
Definition: UserHooks.cc:132
virtual bool canSetResonanceScale()
Possibility to set initial scale in TimeShower for resonance decay.
Definition: UserHooks.h:553
virtual bool doVetoPartonLevelEarly(const Event &e)
Definition: UserHooks.h:520
virtual bool canChangeFragPar()
Can change fragmentation parameters.
Definition: UserHooks.h:186
Definition: StringFragmentation.h:23
double selBias
User-imposed selection bias.
Definition: UserHooks.h:265
virtual bool doVetoFinalTwo(Particle, Particle, const StringEnd *, const StringEnd *)
Definition: UserHooks.h:222
virtual bool canVetoFragmentation()
Can veto hadrons in the fragmentation.
Definition: UserHooks.h:201
virtual bool canModifySigma()
Possibility to modify cross section of process.
Definition: UserHooks.h:357
virtual bool canModifySigma()
Possibility to modify cross section of process.
Definition: UserHooks.h:43
virtual bool doVetoFragmentation(Particle p, const StringEnd *nowEnd)
Do a veto on a hadron just before it is added to the final state.
Definition: UserHooks.h:682
virtual bool canReconnectResonanceSystems()
Possibility to reconnect colours from resonance decay systems.
Definition: UserHooks.h:176
virtual int numberVetoMPIStep()
Up to how many MPI steps should be checked.
Definition: UserHooks.h:493
SuppressSmallPT(double pT0timesMPIIn=1., int numberAlphaSIn=0, bool useSameAlphaSasMPIIn=true)
Constructor.
Definition: UserHooks.h:284
virtual bool doVetoFragmentation(Particle, Particle, const StringEnd *, const StringEnd *)
Definition: UserHooks.h:214
virtual bool canBiasSelection()
Possibility to bias selection of events, compensated by a weight.
Definition: UserHooks.h:374
virtual bool doReconnectResonanceSystems(int j, Event &e)
Definition: UserHooks.h:640
virtual double biasSelectionBy(const SigmaProcess *sigmaProcessPtr, const PhaseSpace *phaseSpacePtr, bool inEvent)
Multiplicative factor in the phase space selection of a hard process.
Definition: UserHooks.cc:79
virtual bool canSetLowEnergySigma(int, int) const
Definition: UserHooks.h:68
virtual ~UserHooksVector()
Destructor.
Definition: UserHooks.h:323
virtual bool canVetoProcessLevel()
Possibility to veto event after process-level selection.
Definition: UserHooks.h:60
virtual bool doVetoMPIEmission(int sizeOld, const Event &e)
Definition: UserHooks.h:620
virtual bool canSetImpactParameter() const
Can set the overall impact parameter for the MPI treatment.
Definition: UserHooks.h:235
Definition: Event.h:32
virtual bool canVetoStep()
Definition: UserHooks.h:97
virtual bool canVetoPartonLevel()
Possibility to veto event after parton-level selection.
Definition: UserHooks.h:132
virtual bool doReconnectResonanceSystems(int, Event &)
Definition: UserHooks.h:183
virtual bool canChangeFragPar()
Can change fragmentation parameters.
Definition: UserHooks.h:648
double enhancedEventWeight
Bookkept quantities for boosted event weights.
Definition: UserHooks.h:268
virtual bool canVetoMPIStep()
Definition: UserHooks.h:486
virtual bool canVetoFSREmission()
Possibility to veto an emission in the FSR machinery.
Definition: UserHooks.h:589
Logger * loggerPtr
Pointer to logger.
Definition: PhysicsBase.h:87
virtual bool doVetoISREmission(int, const Event &, int)
Definition: UserHooks.h:153
virtual bool retryPartonLevel()
Definition: UserHooks.h:530
virtual bool canSetImpactParameter() const
Can set the overall impact parameter for the MPI treatment.
Definition: UserHooks.h:723
virtual bool canVetoPT()
Definition: UserHooks.h:84
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The Pythia class contains the top-level routines to generate an event.
Definition: Pythia.h:71
virtual bool doVetoMPIEmission(int, const Event &)
Definition: UserHooks.h:173
virtual bool canVetoAfterHadronization()
Definition: UserHooks.h:229
Event workEvent
Have one event object around as work area.
Definition: UserHooks.h:262
UserHooksVector implements a vector of UserHooks and is itself a UserHooks.
Definition: UserHooks.h:313
The StringFlav class is used to select quark and hadron flavours.
Definition: FragmentationFlavZpT.h:84
UserHooks is base class for user access to program execution.
Definition: UserHooks.h:32
void subEvent(const Event &event, bool isHardest=true)
subEvent extracts currently resolved partons in the hard process.
Definition: UserHooks.cc:182
Definition: PhaseSpace.h:42
virtual bool canVetoFSREmission()
Possibility to veto an emission in the FSR machinery.
Definition: UserHooks.h:156
void init(string headerIn="", ParticleData *particleDataPtrIn=0, int startColTagIn=100)
Initialize header for event listing, particle data table, and colour.
Definition: Event.h:422
virtual bool initAfterBeams()
Initialisation after beams have been set by Pythia::init().
Definition: UserHooks.h:40
virtual double biasSelectionBy(const SigmaProcess *sigmaProcessPtr, const PhaseSpace *phaseSpacePtr, bool inEvent)
Multiplicative factor in the phase space selection of a hard process.
Definition: UserHooks.h:381
virtual double multiplySigmaBy(const SigmaProcess *sigmaProcessPtr, const PhaseSpace *phaseSpacePtr, bool inEvent)
Multiplicative factor modifying the cross section of a hard process.
Definition: UserHooks.h:364
virtual bool canReconnectResonanceSystems()
Possibility to reconnect colours from resonance decay systems.
Definition: UserHooks.h:629
virtual bool doVetoProcessLevel(Event &)
Definition: UserHooks.h:64
virtual bool canVetoMPIStep()
Definition: UserHooks.h:109
virtual bool doVetoStep(int, int, int, const Event &)
Definition: UserHooks.h:105
virtual bool doVetoPartonLevel(const Event &)
Definition: UserHooks.h:136
virtual bool canVetoFragmentation()
Can veto hadrons in the fragmentation.
Definition: UserHooks.h:676
virtual bool canVetoPartonLevelEarly()
Definition: UserHooks.h:512
virtual bool doVetoFragmentation(Particle p1, Particle p2, const StringEnd *e1, const StringEnd *e2)
Definition: UserHooks.h:689
virtual bool doVetoAfterHadronization(const Event &)
Do the actual veto after hadronization.
Definition: UserHooks.h:232
virtual bool canVetoISREmission()
Possibility to veto an emission in the ISR machinery.
Definition: UserHooks.h:571
virtual double scaleResonance(int iRes, const Event &e)
Definition: UserHooks.h:562
virtual double doSetLowEnergySigma(int, int, double, double, double) const
Definition: UserHooks.h:72
virtual double doSetImpactParameter()
Set the overall impact parameter for the MPI treatment.
Definition: UserHooks.h:238