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