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