PYTHIA  8.317
UserHooks.h
1 // UserHooks.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2026 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  // Can set the enhancement for the MPI treatment.
241  virtual bool canSetEnhanceB() const { return false; }
242 
243  // Set the enhancement for the MPI treatment.
244  virtual double doSetEnhanceB() { return 0.0; }
245 
246  // Custom processing at the end of HadronLevel::next.
247  virtual bool onEndHadronLevel(HadronLevel&, Event&) { return true; }
248 
249 protected:
250 
251  // Constructor.
253 
254  // After initInfoPtr, initialize workEvent.
255  virtual void onInitInfoPtr() override {
256  // Set smart pointer to null, in order to avoid circular dependency.
257  userHooksPtr = nullptr;
258  workEvent.init("(work event)", particleDataPtr);
259  }
260 
261  // omitResonanceDecays omits resonance decay chains from process record.
262  void omitResonanceDecays(const Event& process, bool finalOnly = false);
263 
264  // subEvent extracts currently resolved partons in the hard process.
265  void subEvent(const Event& event, bool isHardest = true);
266 
267  // Have one event object around as work area.
269 
270  // User-imposed selection bias.
271  double selBias = 1.;
272 
273  // Bookkept quantities for boosted event weights.
274  double enhancedEventWeight = {}, pTEnhanced = {}, wtEnhanced = {};
275 
276 };
277 
278 //==========================================================================
279 
280 // SuppressSmallPT is a derived class for user access to program execution.
281 // It is a simple example, illustrating how to suppress the cross section
282 // of 2 -> 2 processes by a factor pT^4 / (pT0^2 + pT^2)^2, with pT0 input,
283 // and also modify alpha_strong scale similarly.
284 
285 class SuppressSmallPT : public UserHooks {
286 
287 public:
288 
289  // Constructor.
290  SuppressSmallPT( double pT0timesMPIIn = 1., int numberAlphaSIn = 0,
291  bool useSameAlphaSasMPIIn = true) : pT20(0.) {isInit = false;
292  pT0timesMPI = pT0timesMPIIn; numberAlphaS = numberAlphaSIn;
293  useSameAlphaSasMPI = useSameAlphaSasMPIIn;}
294 
295  // Possibility to modify cross section of process.
296  virtual bool canModifySigma() {return true;}
297 
298  // Multiplicative factor modifying the cross section of a hard process.
299  // Usage: inEvent is true for event generation, false for initialization.
300  virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
301  const PhaseSpace* phaseSpacePtr, bool );
302 
303 private:
304 
305  // Save input properties and the squared pT0 scale.
306  bool isInit, useSameAlphaSasMPI;
307  int numberAlphaS;
308  double pT0timesMPI, pT20;
309 
310  // Alpha_strong calculation.
311  AlphaStrong alphaS;
312 
313 };
314 
315 //==========================================================================
316 
317 // UserHooksVector implements a vector of UserHooks and is itself a UserHooks.
318 
319 class UserHooksVector: public UserHooks {
320 
321 public:
322 
323  // The default constructor is private, and should only be used
324  // internally in Pythia.
326  friend class Pythia;
327 
328  // Destructor.
329  virtual ~UserHooksVector() {}
330 
331  // Initialisation after beams have been set by Pythia::init().
332  // Check that there are no (obvious) clashes.
333  virtual bool initAfterBeams() {
334  int nCanSetResonanceScale = 0;
335  int nCanChangeFragPar = 0;
336  int nCanSetImpactParameter = 0;
337  int nCanSetEnhanceB = 0;
338  for ( int i = 0, N = hooks.size(); i < N; ++i ) {
339  registerSubObject(*hooks[i]);
340  if ( !hooks[i]->initAfterBeams() ) return false;
341  if (hooks[i]->canSetResonanceScale()) ++nCanSetResonanceScale;
342  if (hooks[i]->canChangeFragPar()) ++nCanChangeFragPar;
343  if (hooks[i]->canSetImpactParameter()) ++nCanSetImpactParameter;
344  if (hooks[i]->canSetEnhanceB()) ++nCanSetEnhanceB;
345  }
346  if (nCanSetResonanceScale > 1) {
347  loggerPtr->ERROR_MSG(
348  "multiple UserHooks with canSetResonanceScale() not allowed");
349  return false;
350  }
351  if (nCanChangeFragPar > 1) {
352  loggerPtr->ERROR_MSG(
353  "multiple UserHooks with canChangeFragPar() not allowed");
354  return false;
355  }
356  if (nCanSetImpactParameter > 1) {
357  loggerPtr->ERROR_MSG(
358  "multiple UserHooks with canSetImpactParameter() not allowed");
359  return false;
360  }
361  if (nCanSetEnhanceB > 1) {
362  loggerPtr->ERROR_MSG(
363  "multiple UserHooks with canSetEnhanceB() not allowed");
364  return false;
365  }
366  return true;
367  }
368 
369  // Possibility to modify cross section of process.
370  virtual bool canModifySigma() {
371  for ( int i = 0, N = hooks.size(); i < N; ++i )
372  if ( hooks[i]->canModifySigma() ) return true;
373  return false;
374  }
375 
376  // Multiplicative factor modifying the cross section of a hard process.
377  virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
378  const PhaseSpace* phaseSpacePtr, bool inEvent) {
379  double f = 1.0;
380  for ( int i = 0, N = hooks.size(); i < N; ++i )
381  if ( hooks[i]->canModifySigma() )
382  f *= hooks[i]->multiplySigmaBy(sigmaProcessPtr, phaseSpacePtr,inEvent);
383  return f;
384  }
385 
386  // Possibility to bias selection of events, compensated by a weight.
387  virtual bool canBiasSelection() {
388  for ( int i = 0, N = hooks.size(); i < N; ++i )
389  if ( hooks[i]->canBiasSelection() ) return true;
390  return false;
391  }
392 
393  // Multiplicative factor in the phase space selection of a hard process.
394  virtual double biasSelectionBy(const SigmaProcess* sigmaProcessPtr,
395  const PhaseSpace* phaseSpacePtr, bool inEvent) {
396  double f = 1.0;
397  for ( int i = 0, N = hooks.size(); i < N; ++i )
398  if ( hooks[i]->canBiasSelection() )
399  f *= hooks[i]->biasSelectionBy(sigmaProcessPtr, phaseSpacePtr,
400  inEvent);
401  return f;
402  }
403 
404  // Event weight to compensate for selection weight above.
405  virtual double biasedSelectionWeight() {
406  double f = 1.0;
407  for ( int i = 0, N = hooks.size(); i < N; ++i )
408  if ( hooks[i]->canBiasSelection() )
409  f *= hooks[i]->biasedSelectionWeight();
410  return f;
411  }
412 
413  // Possibility to veto event after process-level selection.
414  virtual bool canVetoProcessLevel() {
415  for ( int i = 0, N = hooks.size(); i < N; ++i )
416  if ( hooks[i]->canVetoProcessLevel() ) return true;
417  return false;
418  }
419 
420  // Decide whether to veto current process or not, based on process record.
421  // Usage: doVetoProcessLevel( process).
422  virtual bool doVetoProcessLevel(Event& e) {
423  for ( int i = 0, N = hooks.size(); i < N; ++i )
424  if ( hooks[i]->canVetoProcessLevel() &&
425  hooks[i]->doVetoProcessLevel(e) ) return true;
426  return false;
427  }
428 
429  // Possibility to veto resonance decay chain.
430  virtual bool canVetoResonanceDecays() {
431  for ( int i = 0, N = hooks.size(); i < N; ++i )
432  if ( hooks[i]->canVetoResonanceDecays() ) return true;
433  return false;
434  }
435 
436  // Decide whether to veto current resonance decay chain or not, based on
437  // process record. Usage: doVetoProcessLevel( process).
438  virtual bool doVetoResonanceDecays(Event& e) {
439  for ( int i = 0, N = hooks.size(); i < N; ++i )
440  if ( hooks[i]->canVetoResonanceDecays() &&
441  hooks[i]->doVetoResonanceDecays(e) ) return true;
442  return false;
443  }
444 
445  // Possibility to veto MPI + ISR + FSR evolution and kill event,
446  // making decision at a fixed pT scale. Useful for MLM-style matching.
447  virtual bool canVetoPT() {
448  for ( int i = 0, N = hooks.size(); i < N; ++i )
449  if ( hooks[i]->canVetoPT() ) return true;
450  return false;
451  }
452 
453  // Transverse-momentum scale for veto test.
454  virtual double scaleVetoPT() {
455  double s = 0.0;
456  for ( int i = 0, N = hooks.size(); i < N; ++i )
457  if ( hooks[i]->canVetoPT() ) s = max(s, hooks[i]->scaleVetoPT());
458  return s;
459  }
460 
461  // Decide whether to veto current event or not, based on event record.
462  // Usage: doVetoPT( iPos, event), where iPos = 0: no emissions so far;
463  // iPos = 1/2/3 joint evolution, latest step was MPI/ISR/FSR;
464  // iPos = 4: FSR only afterwards; iPos = 5: FSR in resonance decay.
465  virtual bool doVetoPT( int iPos, const Event& e) {
466  for ( int i = 0, N = hooks.size(); i < N; ++i )
467  if ( hooks[i]->canVetoPT() && hooks[i]->doVetoPT(iPos, e) ) return true;
468  return false;
469  }
470 
471  // Possibility to veto MPI + ISR + FSR evolution and kill event,
472  // making decision after fixed number of ISR or FSR steps.
473  virtual bool canVetoStep() {
474  for ( int i = 0, N = hooks.size(); i < N; ++i )
475  if ( hooks[i]->canVetoStep() ) return true;
476  return false;
477  }
478 
479  // Up to how many ISR + FSR steps of hardest interaction should be checked.
480  virtual int numberVetoStep() {
481  int n = 1;
482  for ( int i = 0, N = hooks.size(); i < N; ++i )
483  if ( hooks[i]->canVetoStep() ) n = max(n, hooks[i]->numberVetoStep());
484  return n;
485  }
486 
487  // Decide whether to veto current event or not, based on event record.
488  // Usage: doVetoStep( iPos, nISR, nFSR, event), where iPos as above,
489  // nISR and nFSR number of emissions so far for hard interaction only.
490  virtual bool doVetoStep( int iPos, int nISR, int nFSR, const Event& e) {
491  for ( int i = 0, N = hooks.size(); i < N; ++i )
492  if ( hooks[i]->canVetoStep()
493  && hooks[i]->doVetoStep(iPos, nISR, nFSR, e) ) return true;
494  return false;
495  }
496 
497  // Possibility to veto MPI + ISR + FSR evolution and kill event,
498  // making decision after fixed number of MPI steps.
499  virtual bool canVetoMPIStep() {
500  for ( int i = 0, N = hooks.size(); i < N; ++i )
501  if ( hooks[i]->canVetoMPIStep() ) return true;
502  return false;
503  }
504 
505  // Up to how many MPI steps should be checked.
506  virtual int numberVetoMPIStep() {
507  int n = 1;
508  for ( int i = 0, N = hooks.size(); i < N; ++i )
509  if ( hooks[i]->canVetoMPIStep() )
510  n = max(n, hooks[i]->numberVetoMPIStep());
511  return n;
512  }
513 
514  // Decide whether to veto current event or not, based on event record.
515  // Usage: doVetoMPIStep( nMPI, event), where nMPI is number of MPI's so far.
516  virtual bool doVetoMPIStep( int nMPI, const Event& e) {
517  for ( int i = 0, N = hooks.size(); i < N; ++i )
518  if ( hooks[i]->canVetoMPIStep() && hooks[i]->doVetoMPIStep(nMPI, e) )
519  return true;
520  return false;
521  }
522 
523  // Possibility to veto event after ISR + FSR + MPI in parton level,
524  // but before beam remnants and resonance decays.
525  virtual bool canVetoPartonLevelEarly() {
526  for ( int i = 0, N = hooks.size(); i < N; ++i )
527  if ( hooks[i]->canVetoPartonLevelEarly() ) return true;
528  return false;
529  }
530 
531  // Decide whether to veto current partons or not, based on event record.
532  // Usage: doVetoPartonLevelEarly( event).
533  virtual bool doVetoPartonLevelEarly( const Event& e) {
534  for ( int i = 0, N = hooks.size(); i < N; ++i )
535  if ( hooks[i]->canVetoPartonLevelEarly()
536  && hooks[i]->doVetoPartonLevelEarly(e) ) return true;
537  return false;
538  }
539 
540  // Retry same ProcessLevel with a new PartonLevel after a veto in
541  // doVetoPT, doVetoStep, doVetoMPIStep or doVetoPartonLevelEarly
542  // if you overload this method to return true.
543  virtual bool retryPartonLevel() {
544  for ( int i = 0, N = hooks.size(); i < N; ++i )
545  if ( hooks[i]->retryPartonLevel() ) return true;
546  return false;
547  }
548 
549  // Possibility to veto event after parton-level selection.
550  virtual bool canVetoPartonLevel() {
551  for ( int i = 0, N = hooks.size(); i < N; ++i )
552  if ( hooks[i]->canVetoPartonLevel() ) return true;
553  return false;
554  }
555 
556  // Decide whether to veto current partons or not, based on event record.
557  // Usage: doVetoPartonLevel( event).
558  virtual bool doVetoPartonLevel( const Event& e) {
559  for ( int i = 0, N = hooks.size(); i < N; ++i )
560  if ( hooks[i]->canVetoPartonLevel()
561  && hooks[i]->doVetoPartonLevel(e) ) return true;
562  return false;
563  }
564 
565  // Possibility to set initial scale in TimeShower for resonance decay.
566  virtual bool canSetResonanceScale() {
567  for ( int i = 0, N = hooks.size(); i < N; ++i )
568  if ( hooks[i]->canSetResonanceScale() ) return true;
569  return false;
570  }
571 
572  // Initial scale for TimeShower evolution.
573  // Usage: scaleResonance( iRes, event), where iRes is location
574  // of decaying resonance in the event record.
575  virtual double scaleResonance( int iRes, const Event& e) {
576  double s = 0.0;
577  for ( int i = 0, N = hooks.size(); i < N; ++i )
578  if ( hooks[i]->canSetResonanceScale() )
579  s = max(s, hooks[i]->scaleResonance(iRes, e));
580  return s;
581  }
582 
583  // Possibility to veto an emission in the ISR machinery.
584  virtual bool canVetoISREmission() {
585  for ( int i = 0, N = hooks.size(); i < N; ++i )
586  if ( hooks[i]->canVetoISREmission() ) return true;
587  return false;
588  }
589 
590  // Decide whether to veto current emission or not, based on event record.
591  // Usage: doVetoISREmission( sizeOld, event, iSys) where sizeOld is size
592  // of event record before current emission-to-be-scrutinized was added,
593  // and iSys is the system of the radiation (according to PartonSystems).
594  virtual bool doVetoISREmission( int sizeOld, const Event& e, int iSys) {
595  for ( int i = 0, N = hooks.size(); i < N; ++i )
596  if ( hooks[i]->canVetoISREmission()
597  && hooks[i]->doVetoISREmission(sizeOld, e, iSys) ) return true;
598  return false;
599  }
600 
601  // Possibility to veto an emission in the FSR machinery.
602  virtual bool canVetoFSREmission() {
603  for ( int i = 0, N = hooks.size(); i < N; ++i )
604  if ( hooks[i]->canVetoFSREmission() ) return true;
605  return false;
606  }
607 
608  // Decide whether to veto current emission or not, based on event record.
609  // Usage: doVetoFSREmission( sizeOld, event, iSys, inResonance) where
610  // sizeOld is size of event record before current emission-to-be-scrutinized
611  // was added, iSys is the system of the radiation (according to
612  // PartonSystems), and inResonance is true if the emission takes place in a
613  // resonance decay.
614  virtual bool doVetoFSREmission(int sizeOld, const Event& e,
615  int iSys, bool inResonance = false ) {
616  for ( int i = 0, N = hooks.size(); i < N; ++i )
617  if ( hooks[i]->canVetoFSREmission()
618  && hooks[i]->doVetoFSREmission(sizeOld, e, iSys, inResonance) )
619  return true;
620  return false;
621  }
622 
623  // Possibility to veto an MPI.
624  virtual bool canVetoMPIEmission() {
625  for ( int i = 0, N = hooks.size(); i < N; ++i )
626  if ( hooks[i]->canVetoMPIEmission() ) return true;
627  return false;
628  }
629 
630  // Decide whether to veto an MPI based on event record.
631  // Usage: doVetoMPIEmission( sizeOld, event) where sizeOld
632  // is size of event record before the current MPI.
633  virtual bool doVetoMPIEmission( int sizeOld, const Event & e) {
634  for ( int i = 0, N = hooks.size(); i < N; ++i )
635  if ( hooks[i]->canVetoMPIEmission()
636  && hooks[i]->doVetoMPIEmission(sizeOld, e) )
637  return true;
638  return false;
639  }
640 
641  // Possibility to reconnect colours from resonance decay systems.
643  for ( int i = 0, N = hooks.size(); i < N; ++i )
644  if ( hooks[i]->canReconnectResonanceSystems() ) return true;
645  return false;
646  }
647 
648  // Do reconnect colours from resonance decay systems.
649  // Usage: doVetoFSREmission( oldSizeEvt, event)
650  // where oldSizeEvent is the event size before resonance decays.
651  // Should normally return true, while false means serious failure.
652  // Value of PartonLevel:earlyResDec determines where method is called.
653  virtual bool doReconnectResonanceSystems( int j, Event & e) {
654  for ( int i = 0, N = hooks.size(); i < N; ++i )
655  if ( hooks[i]->canReconnectResonanceSystems()
656  && hooks[i]->doReconnectResonanceSystems(j, e) ) return true;
657  return false;
658  }
659 
660  // Can change fragmentation parameters.
661  virtual bool canChangeFragPar() {
662  for ( int i = 0, N = hooks.size(); i < N; ++i )
663  if ( hooks[i]->canChangeFragPar() ) return true;
664  return false;
665  }
666 
667  // Set initial ends of a string to be fragmented. This is done once
668  // for each string. Note that the second string end may be zero in
669  // case we are hadronising a string piece leading to a junction.
670  virtual void setStringEnds( const StringEnd* pos, const StringEnd* neg,
671  vector<int> iPart) {
672  for ( int i = 0, N = hooks.size(); i < N; ++i )
673  hooks[i]->setStringEnds( pos, neg, iPart);
674  }
675 
676  // Do change fragmentation parameters.
677  // Input: flavPtr, zPtr, pTPtr, idEnd, m2Had, iParton and posEnd (or
678  // negEnd).
679  virtual bool doChangeFragPar( StringFlav* sfIn, StringZ* zIn,
680  StringPT* ptIn, int idIn, double mIn, vector<int> parIn,
681  const StringEnd* endIn) {
682  for ( int i = 0, N = hooks.size(); i < N; ++i )
683  if ( hooks[i]->canChangeFragPar()
684  && hooks[i]->doChangeFragPar(sfIn, zIn, ptIn, idIn,
685  mIn, parIn, endIn) ) return true;
686  return false;}
687 
688  // Can veto hadrons in the fragmentation.
689  virtual bool canVetoFragmentation() {
690  for ( int i = 0, N = hooks.size(); i < N; ++i )
691  if ( hooks[i]->canVetoFragmentation() ) return true;
692  return false;}
693 
694  // Do a veto on a hadron just before it is added to the final state.
695  virtual bool doVetoFragmentation(Particle p, const StringEnd* nowEnd) {
696  for ( int i = 0, N = hooks.size(); i < N; ++i )
697  if ( hooks[i]->canChangeFragPar()
698  && hooks[i]->doVetoFragmentation(p, nowEnd) ) return true;
699  return false;
700  }
701 
702  virtual bool doVetoFragmentation(Particle p1, Particle p2,
703  const StringEnd* e1, const StringEnd* e2) {
704  for ( int i = 0, N = hooks.size(); i < N; ++i )
705  if ( hooks[i]->canChangeFragPar()
706  && hooks[i]->doVetoFragmentation(p1, p2, e1, e2) ) return true;
707  return false;
708  }
709 
710  virtual bool doVetoFinalTwo(Particle p1, Particle p2,
711  const StringEnd* e1, const StringEnd* e2) {
712  for ( int i = 0, N = hooks.size(); i < N; ++i )
713  if ( hooks[i]->canChangeFragPar()
714  && hooks[i]->doVetoFinalTwo(p1, p2, e1, e2) ) return true;
715  return false;
716  }
717 
718  // Possibility to veto an event after hadronization based on event
719  // contents. Works as an early trigger to avoid running the time
720  // consuming rescattering process on uninteresting events.
721  virtual bool canVetoAfterHadronization() {
722  for ( int i = 0, N = hooks.size(); i < N; ++i )
723  if ( hooks[i]->canVetoAfterHadronization() ) return true;
724  return false;
725  }
726 
727  // Do the actual veto after hadronization.
728  virtual bool doVetoAfterHadronization(const Event& e) {
729  for ( int i = 0, N = hooks.size(); i < N; ++i )
730  if ( hooks[i]->canVetoAfterHadronization()
731  && hooks[i]->doVetoAfterHadronization(e) ) return true;
732  return false;
733  }
734 
735  // Can set the overall impact parameter for the MPI treatment.
736  virtual bool canSetImpactParameter() const {
737  for ( int i = 0, N = hooks.size(); i < N; ++i )
738  if ( hooks[i]->canSetImpactParameter() ) return true;
739  return false;
740  }
741 
742  // Set the overall impact parameter for the MPI treatment.
743  virtual double doSetImpactParameter() {
744  for ( int i = 0, N = hooks.size(); i < N; ++i )
745  if ( hooks[i]->canSetImpactParameter() )
746  return hooks[i]->doSetImpactParameter();
747  return 0.0;
748  }
749 
750  // Can set the enhancement for the MPI treatment.
751  virtual bool canSetEnhanceB() const {
752  for ( int i = 0, N = hooks.size(); i < N; ++i )
753  if ( hooks[i]->canSetEnhanceB() ) return true;
754  return false;
755  }
756 
757  // Set the enhancement for the MPI treatment.
758  virtual double doSetEnhanceB() {
759  for ( int i = 0, N = hooks.size(); i < N; ++i )
760  if ( hooks[i]->canSetEnhanceB() )
761  return hooks[i]->doSetEnhanceB();
762  return 0.0;
763  }
764 
765  // The vector of user hooks.
766  vector< shared_ptr<UserHooks> > hooks = {};
767 
768 };
769 
770 //==========================================================================
771 
772 } // end namespace Pythia8
773 
774 #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:480
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:710
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:447
virtual double scaleVetoPT()
Transverse-momentum scale for veto test.
Definition: UserHooks.h:454
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:325
Definition: PhysicsBase.h:26
virtual bool canModifySigma()
Possibility to modify cross section of process.
Definition: UserHooks.h:296
void registerSubObject(PhysicsBase &pb)
Register a sub object that should have its information in sync with this.
Definition: PhysicsBase.cc:57
virtual double doSetImpactParameter()
Set the overall impact parameter for the MPI treatment.
Definition: UserHooks.h:743
virtual bool canSetEnhanceB() const
Can set the enhancement for the MPI treatment.
Definition: UserHooks.h:241
The Event class holds all info on the generated event.
Definition: Event.h:408
UserHooks()
Constructor.
Definition: UserHooks.h:252
virtual bool doVetoISREmission(int sizeOld, const Event &e, int iSys)
Definition: UserHooks.h:594
virtual bool retryPartonLevel()
Definition: UserHooks.h:129
virtual bool onEndHadronLevel(HadronLevel &, Event &)
Custom processing at the end of HadronLevel::next.
Definition: UserHooks.h:247
SigmaProcess is the base class for cross section calculations.
Definition: SigmaProcess.h:86
virtual bool doVetoPT(int iPos, const Event &e)
Definition: UserHooks.h:465
virtual bool doVetoPartonLevel(const Event &e)
Definition: UserHooks.h:558
virtual bool doVetoMPIStep(int nMPI, const Event &e)
Definition: UserHooks.h:516
virtual double biasedSelectionWeight()
Event weight to compensate for selection weight above.
Definition: UserHooks.h:405
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:414
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:300
virtual bool canVetoMPIEmission()
Possibility to veto an MPI.
Definition: UserHooks.h:624
virtual bool doVetoStep(int iPos, int nISR, int nFSR, const Event &e)
Definition: UserHooks.h:490
virtual bool canVetoStep()
Definition: UserHooks.h:473
virtual bool canVetoAfterHadronization()
Definition: UserHooks.h:721
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:670
virtual bool doVetoFSREmission(int sizeOld, const Event &e, int iSys, bool inResonance=false)
Definition: UserHooks.h:614
virtual bool canSetEnhanceB() const
Can set the enhancement for the MPI treatment.
Definition: UserHooks.h:751
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:679
virtual bool doVetoResonanceDecays(Event &e)
Definition: UserHooks.h:438
Definition: UserHooks.h:285
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:88
virtual bool doVetoAfterHadronization(const Event &e)
Do the actual veto after hadronization.
Definition: UserHooks.h:728
virtual bool initAfterBeams()
Definition: UserHooks.h:333
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:219
virtual bool canVetoResonanceDecays()
Possibility to veto resonance decay chain.
Definition: UserHooks.h:430
virtual double scaleResonance(int, const Event &)
Definition: UserHooks.h:144
virtual bool canVetoPartonLevel()
Possibility to veto event after parton-level selection.
Definition: UserHooks.h:550
Definition: HadronLevel.h:44
virtual bool canVetoMPIEmission()
Possibility to veto an MPI.
Definition: UserHooks.h:168
UserHooksPtr userHooksPtr
Definition: PhysicsBase.h:128
virtual void onInitInfoPtr() override
After initInfoPtr, initialize workEvent.
Definition: UserHooks.h:255
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:422
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:566
virtual bool doVetoPartonLevelEarly(const Event &e)
Definition: UserHooks.h:533
virtual bool canChangeFragPar()
Can change fragmentation parameters.
Definition: UserHooks.h:186
Definition: StringFragmentation.h:23
double selBias
User-imposed selection bias.
Definition: UserHooks.h:271
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:370
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:695
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:506
SuppressSmallPT(double pT0timesMPIIn=1., int numberAlphaSIn=0, bool useSameAlphaSasMPIIn=true)
Constructor.
Definition: UserHooks.h:290
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:387
virtual bool doReconnectResonanceSystems(int j, Event &e)
Definition: UserHooks.h:653
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:329
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:633
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:661
double enhancedEventWeight
Bookkept quantities for boosted event weights.
Definition: UserHooks.h:274
virtual bool canVetoMPIStep()
Definition: UserHooks.h:499
virtual bool canVetoFSREmission()
Possibility to veto an emission in the FSR machinery.
Definition: UserHooks.h:602
Logger * loggerPtr
Pointer to logger.
Definition: PhysicsBase.h:91
virtual bool doVetoISREmission(int, const Event &, int)
Definition: UserHooks.h:153
virtual bool retryPartonLevel()
Definition: UserHooks.h:543
virtual bool canSetImpactParameter() const
Can set the overall impact parameter for the MPI treatment.
Definition: UserHooks.h:736
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:72
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:268
UserHooksVector implements a vector of UserHooks and is itself a UserHooks.
Definition: UserHooks.h:319
The StringFlav class is used to select quark and hadron flavours.
Definition: FragmentationFlavZpT.h:76
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:394
virtual double doSetEnhanceB()
Set the enhancement for the MPI treatment.
Definition: UserHooks.h:758
virtual double multiplySigmaBy(const SigmaProcess *sigmaProcessPtr, const PhaseSpace *phaseSpacePtr, bool inEvent)
Multiplicative factor modifying the cross section of a hard process.
Definition: UserHooks.h:377
virtual bool canReconnectResonanceSystems()
Possibility to reconnect colours from resonance decay systems.
Definition: UserHooks.h:642
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 double doSetEnhanceB()
Set the enhancement for the MPI treatment.
Definition: UserHooks.h:244
virtual bool doVetoPartonLevel(const Event &)
Definition: UserHooks.h:136
virtual bool canVetoFragmentation()
Can veto hadrons in the fragmentation.
Definition: UserHooks.h:689
virtual bool canVetoPartonLevelEarly()
Definition: UserHooks.h:525
virtual bool doVetoFragmentation(Particle p1, Particle p2, const StringEnd *e1, const StringEnd *e2)
Definition: UserHooks.h:702
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:584
virtual double scaleResonance(int iRes, const Event &e)
Definition: UserHooks.h:575
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