Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-09-18 09:25:40

0001 // UserHooks.h is a part of the PYTHIA event generator.
0002 // Copyright (C) 2025 Torbjorn Sjostrand.
0003 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
0004 // Please respect the MCnet Guidelines, see GUIDELINES for details.
0005 
0006 // Header file to allow user access to program at different stages.
0007 // UserHooks: almost empty base class, with user to write the rela code.
0008 // MyUserHooks: derived class, only intended as an example.
0009 
0010 #ifndef Pythia8_UserHooks_H
0011 #define Pythia8_UserHooks_H
0012 
0013 #include "Pythia8/Event.h"
0014 #include "Pythia8/PartonSystems.h"
0015 #include "Pythia8/PhysicsBase.h"
0016 #include "Pythia8/PythiaStdlib.h"
0017 #include "Pythia8/SigmaProcess.h"
0018 
0019 namespace Pythia8 {
0020 
0021 //==========================================================================
0022 
0023 // Forward references.
0024 class PhaseSpace;
0025 class StringEnd;
0026 class HadronLevel;
0027 
0028 //==========================================================================
0029 
0030 // UserHooks is base class for user access to program execution.
0031 
0032 class UserHooks : public PhysicsBase {
0033 
0034 public:
0035 
0036   // Destructor.
0037   virtual ~UserHooks() {}
0038 
0039   // Initialisation after beams have been set by Pythia::init().
0040   virtual bool initAfterBeams() { return true; }
0041 
0042   // Possibility to modify cross section of process.
0043   virtual bool canModifySigma() {return false;}
0044 
0045   // Multiplicative factor modifying the cross section of a hard process.
0046   virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
0047     const PhaseSpace* phaseSpacePtr, bool inEvent);
0048 
0049   // Possibility to bias selection of events, compensated by a weight.
0050   virtual bool canBiasSelection() {return false;}
0051 
0052   // Multiplicative factor in the phase space selection of a hard process.
0053   virtual double biasSelectionBy(const SigmaProcess* sigmaProcessPtr,
0054     const PhaseSpace* phaseSpacePtr, bool inEvent);
0055 
0056   // Event weight to compensate for selection weight above.
0057   virtual double biasedSelectionWeight() {return 1./selBias;}
0058 
0059   // Possibility to veto event after process-level selection.
0060   virtual bool canVetoProcessLevel() {return false;}
0061 
0062   // Decide whether to veto current process or not, based on process record.
0063   // Usage: doVetoProcessLevel( process).
0064   virtual bool doVetoProcessLevel(Event& ) {return false;}
0065 
0066   // Possibility to set low energy cross sections.
0067   // Usage: canSetLowEnergySigma(idA, idB).
0068   virtual bool canSetLowEnergySigma(int, int) const {return false;}
0069 
0070   // Set low energy cross sections for ids where canSetLowEnergySigma is true.
0071   // Usage: doSetLowEnergySigma(idA, idB, eCM, mA, mB).
0072   virtual double doSetLowEnergySigma(int, int, double, double, double) const {
0073     return 0.;}
0074 
0075   // Possibility to veto resonance decay chain.
0076   virtual bool canVetoResonanceDecays() {return false;}
0077 
0078   // Decide whether to veto current resonance decay chain or not, based on
0079   // process record. Usage: doVetoProcessLevel( process).
0080   virtual bool doVetoResonanceDecays(Event& ) {return false;}
0081 
0082   // Possibility to veto MPI + ISR + FSR evolution and kill event,
0083   // making decision at a fixed pT scale. Useful for MLM-style matching.
0084   virtual bool canVetoPT() {return false;}
0085 
0086   // Transverse-momentum scale for veto test.
0087   virtual double scaleVetoPT() {return 0.;}
0088 
0089   // Decide whether to veto current event or not, based on event record.
0090   // Usage: doVetoPT( iPos, event), where iPos = 0: no emissions so far;
0091   // iPos = 1/2/3 joint evolution, latest step was MPI/ISR/FSR;
0092   // iPos = 4: FSR only afterwards; iPos = 5: FSR in resonance decay.
0093   virtual bool doVetoPT( int , const Event& ) {return false;}
0094 
0095   // Possibility to veto MPI + ISR + FSR evolution and kill event,
0096   // making decision after fixed number of ISR or FSR steps.
0097   virtual bool canVetoStep() {return false;}
0098 
0099   // Up to how many ISR + FSR steps of hardest interaction should be checked.
0100   virtual int numberVetoStep() {return 1;}
0101 
0102   // Decide whether to veto current event or not, based on event record.
0103   // Usage: doVetoStep( iPos, nISR, nFSR, event), where iPos as above,
0104   // nISR and nFSR number of emissions so far for hard interaction only.
0105   virtual bool doVetoStep( int , int , int , const Event& ) {return false;}
0106 
0107   // Possibility to veto MPI + ISR + FSR evolution and kill event,
0108   // making decision after fixed number of MPI steps.
0109   virtual bool canVetoMPIStep() {return false;}
0110 
0111   // Up to how many MPI steps should be checked.
0112   virtual int numberVetoMPIStep() {return 1;}
0113 
0114   // Decide whether to veto current event or not, based on event record.
0115   // Usage: doVetoMPIStep( nMPI, event), where nMPI is number of MPI's so far.
0116   virtual bool doVetoMPIStep( int , const Event& ) {return false;}
0117 
0118   // Possibility to veto event after ISR + FSR + MPI in parton level,
0119   // but before beam remnants and resonance decays.
0120   virtual bool canVetoPartonLevelEarly() {return false;}
0121 
0122   // Decide whether to veto current partons or not, based on event record.
0123   // Usage: doVetoPartonLevelEarly( event).
0124   virtual bool doVetoPartonLevelEarly( const Event& ) {return false;}
0125 
0126   // Retry same ProcessLevel with a new PartonLevel after a veto in
0127   // doVetoPT, doVetoStep, doVetoMPIStep or doVetoPartonLevelEarly
0128   // if you overload this method to return true.
0129   virtual bool retryPartonLevel() {return false;}
0130 
0131   // Possibility to veto event after parton-level selection.
0132   virtual bool canVetoPartonLevel() {return false;}
0133 
0134   // Decide whether to veto current partons or not, based on event record.
0135   // Usage: doVetoPartonLevel( event).
0136   virtual bool doVetoPartonLevel( const Event& ) {return false;}
0137 
0138   // Possibility to set initial scale in TimeShower for resonance decay.
0139   virtual bool canSetResonanceScale() {return false;}
0140 
0141   // Initial scale for TimeShower evolution.
0142   // Usage: scaleResonance( iRes, event), where iRes is location
0143   // of decaying resonance in the event record.
0144   virtual double scaleResonance( int, const Event& ) {return 0.;}
0145 
0146   // Possibility to veto an emission in the ISR machinery.
0147   virtual bool canVetoISREmission() {return false;}
0148 
0149   // Decide whether to veto current emission or not, based on event record.
0150   // Usage: doVetoISREmission( sizeOld, event, iSys) where sizeOld is size
0151   // of event record before current emission-to-be-scrutinized was added,
0152   // and iSys is the system of the radiation (according to PartonSystems).
0153   virtual bool doVetoISREmission( int, const Event&, int ) {return false;}
0154 
0155   // Possibility to veto an emission in the FSR machinery.
0156   virtual bool canVetoFSREmission() {return false;}
0157 
0158   // Decide whether to veto current emission or not, based on event record.
0159   // Usage: doVetoFSREmission( sizeOld, event, iSys, inResonance) where
0160   // sizeOld is size of event record before current emission-to-be-scrutinized
0161   // was added, iSys is the system of the radiation (according to
0162   // PartonSystems), and inResonance is true if the emission takes place in a
0163   // resonance decay.
0164   virtual bool doVetoFSREmission( int, const Event&, int, bool = false )
0165       {return false;}
0166 
0167   // Possibility to veto an MPI.
0168   virtual bool canVetoMPIEmission() { return false; }
0169 
0170   // Decide whether to veto an MPI based on event record.
0171   // Usage: doVetoMPIEmission( sizeOld, event) where sizeOld
0172   // is size of event record before the current MPI.
0173   virtual bool doVetoMPIEmission( int, const Event &) { return false; }
0174 
0175   // Possibility to reconnect colours from resonance decay systems.
0176   virtual bool canReconnectResonanceSystems() { return false; }
0177 
0178   // Do reconnect colours from resonance decay systems.
0179   // Usage: doVetoFSREmission( oldSizeEvt, event)
0180   // where oldSizeEvent is the event size before resonance decays.
0181   // Should normally return true, while false means serious failure.
0182   // Value of PartonLevel:earlyResDec determines where method is called.
0183   virtual bool doReconnectResonanceSystems( int, Event &) {return true;}
0184 
0185   // Can change fragmentation parameters.
0186   virtual bool canChangeFragPar() { return false;}
0187 
0188   // Set initial ends of a string to be fragmented. This is done once
0189   // for each string. Note that the second string end may be zero in case
0190   // we are hadronising a string piece leading to a junction.
0191   virtual void setStringEnds( const StringEnd*, const StringEnd*,
0192     vector<int>) {}
0193 
0194   // Do change fragmentation parameters.
0195   // Input: flavPtr, zPtr, pTPtr, idEnd, m2Had, iParton and posEnd (or
0196   // negEnd).
0197   virtual bool doChangeFragPar( StringFlav*, StringZ*, StringPT*, int,
0198     double, vector<int>, const StringEnd* ) { return false;}
0199 
0200   // Can veto hadrons in the fragmentation.
0201   virtual bool canVetoFragmentation() { return false;}
0202 
0203   // Do a veto on a hadron just before it is added to the final state.
0204   // The StringEnd from which the the hadron was produced is included
0205   // for information.
0206   virtual bool doVetoFragmentation( Particle, const StringEnd *)
0207     { return false;}
0208 
0209   // Do a veto on a hadron just before it is added to the final state
0210   // (final two hadron case).  Note that a veto from this function
0211   // will refragment the whole string. If only the regeneration of
0212   // the final two is desired, the doVetoFinalTwo function should be
0213   // used instead.
0214   virtual bool doVetoFragmentation(Particle, Particle,
0215     const StringEnd*, const StringEnd* ) { return false;}
0216 
0217   // Do a veto on a hadron just before it is added to the final state
0218   // (final two hadron case). Note that a veto from this function will
0219   // regenerate the final two. If the refragmentation of the whole
0220   // string is desired, the doVetoFragmentation function should be
0221   // used instead.
0222   virtual bool doVetoFinalTwo(Particle, Particle,
0223     const StringEnd*, const StringEnd* ) { return false;}
0224 
0225   // Possibility to veto an event after hadronization based
0226   // on event contents. Works as an early trigger to avoid
0227   // running the time consuming rescattering process on
0228   // uninteresting events.
0229   virtual bool canVetoAfterHadronization() {return false;}
0230 
0231   // Do the actual veto after hadronization.
0232   virtual bool doVetoAfterHadronization(const Event& ) {return false;}
0233 
0234   // Can set the overall impact parameter for the MPI treatment.
0235   virtual bool canSetImpactParameter() const { return false; }
0236 
0237   // Set the overall impact parameter for the MPI treatment.
0238   virtual double doSetImpactParameter() { return 0.0; }
0239 
0240   // Custom processing at the end of HadronLevel::next.
0241   virtual bool onEndHadronLevel(HadronLevel&, Event&) { return true; }
0242 
0243 protected:
0244 
0245   // Constructor.
0246   UserHooks() {}
0247 
0248   // After initInfoPtr, initialize workEvent.
0249   virtual void onInitInfoPtr() override {
0250     // Set smart pointer to null, in order to avoid circular dependency.
0251     userHooksPtr = nullptr;
0252     workEvent.init("(work event)", particleDataPtr);
0253   }
0254 
0255   // omitResonanceDecays omits resonance decay chains from process record.
0256   void omitResonanceDecays(const Event& process, bool finalOnly = false);
0257 
0258   // subEvent extracts currently resolved partons in the hard process.
0259   void subEvent(const Event& event, bool isHardest = true);
0260 
0261   // Have one event object around as work area.
0262   Event workEvent = {};
0263 
0264   // User-imposed selection bias.
0265   double selBias = 1.;
0266 
0267   // Bookkept quantities for boosted event weights.
0268   double enhancedEventWeight = {}, pTEnhanced = {}, wtEnhanced = {};
0269 
0270 };
0271 
0272 //==========================================================================
0273 
0274 // SuppressSmallPT is a derived class for user access to program execution.
0275 // It is a simple example, illustrating how to suppress the cross section
0276 // of 2 -> 2 processes by a factor pT^4 / (pT0^2 + pT^2)^2, with pT0 input,
0277 // and also modify alpha_strong scale similarly.
0278 
0279 class SuppressSmallPT : public UserHooks {
0280 
0281 public:
0282 
0283   // Constructor.
0284   SuppressSmallPT( double pT0timesMPIIn = 1., int numberAlphaSIn = 0,
0285     bool useSameAlphaSasMPIIn = true) : pT20(0.) {isInit = false;
0286     pT0timesMPI = pT0timesMPIIn; numberAlphaS = numberAlphaSIn;
0287     useSameAlphaSasMPI = useSameAlphaSasMPIIn;}
0288 
0289   // Possibility to modify cross section of process.
0290   virtual bool canModifySigma() {return true;}
0291 
0292   // Multiplicative factor modifying the cross section of a hard process.
0293   // Usage: inEvent is true for event generation, false for initialization.
0294   virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
0295     const PhaseSpace* phaseSpacePtr, bool );
0296 
0297 private:
0298 
0299   // Save input properties and the squared pT0 scale.
0300   bool   isInit, useSameAlphaSasMPI;
0301   int    numberAlphaS;
0302   double pT0timesMPI, pT20;
0303 
0304   // Alpha_strong calculation.
0305   AlphaStrong alphaS;
0306 
0307 };
0308 
0309 //==========================================================================
0310 
0311 // UserHooksVector implements a vector of UserHooks and is itself a UserHooks.
0312 
0313 class UserHooksVector: public UserHooks {
0314 
0315 public:
0316 
0317   // The default constructor is private, and should only be used
0318   // internally in Pythia.
0319   UserHooksVector() {};
0320   friend class Pythia;
0321 
0322   // Destructor.
0323   virtual ~UserHooksVector() {}
0324 
0325   // Initialisation after beams have been set by Pythia::init().
0326   // Check that there are no (obvious) clashes.
0327   virtual bool initAfterBeams() {
0328     int nCanSetResonanceScale  = 0;
0329     int nCanChangeFragPar      = 0;
0330     int nCanSetImpactParameter = 0;
0331     for ( int i = 0, N = hooks.size(); i < N; ++i ) {
0332       registerSubObject(*hooks[i]);
0333       if ( !hooks[i]->initAfterBeams() ) return false;
0334       if (hooks[i]->canSetResonanceScale()) ++nCanSetResonanceScale;
0335       if (hooks[i]->canChangeFragPar()) ++nCanChangeFragPar;
0336       if (hooks[i]->canSetImpactParameter()) ++nCanSetImpactParameter;
0337     }
0338     if (nCanSetResonanceScale > 1) {
0339       loggerPtr->ERROR_MSG(
0340         "multiple UserHooks with canSetResonanceScale() not allowed");
0341       return false;
0342     }
0343     if (nCanChangeFragPar > 1) {
0344       loggerPtr->ERROR_MSG(
0345         "multiple UserHooks with canChangeFragPar() not allowed");
0346       return false;
0347     }
0348     if (nCanSetImpactParameter > 1) {
0349       loggerPtr->ERROR_MSG(
0350         "multiple UserHooks with canSetImpactParameter() not allowed");
0351       return false;
0352     }
0353     return true;
0354   }
0355 
0356   // Possibility to modify cross section of process.
0357   virtual bool canModifySigma() {
0358     for ( int i = 0, N = hooks.size(); i < N; ++i )
0359       if ( hooks[i]->canModifySigma() ) return true;
0360     return false;
0361   }
0362 
0363   // Multiplicative factor modifying the cross section of a hard process.
0364   virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
0365     const PhaseSpace* phaseSpacePtr, bool inEvent) {
0366     double f = 1.0;
0367     for ( int i = 0, N = hooks.size(); i < N; ++i )
0368       if ( hooks[i]->canModifySigma() )
0369         f *= hooks[i]->multiplySigmaBy(sigmaProcessPtr, phaseSpacePtr,inEvent);
0370     return f;
0371   }
0372 
0373   // Possibility to bias selection of events, compensated by a weight.
0374   virtual bool canBiasSelection() {
0375     for ( int i = 0, N = hooks.size(); i < N; ++i )
0376       if ( hooks[i]->canBiasSelection() ) return true;
0377     return false;
0378   }
0379 
0380   // Multiplicative factor in the phase space selection of a hard process.
0381   virtual double biasSelectionBy(const SigmaProcess* sigmaProcessPtr,
0382     const PhaseSpace* phaseSpacePtr, bool inEvent) {
0383     double f = 1.0;
0384     for ( int i = 0, N = hooks.size(); i < N; ++i )
0385       if ( hooks[i]->canBiasSelection() )
0386         f *= hooks[i]->biasSelectionBy(sigmaProcessPtr, phaseSpacePtr,
0387              inEvent);
0388     return f;
0389   }
0390 
0391   // Event weight to compensate for selection weight above.
0392   virtual double biasedSelectionWeight() {
0393     double f = 1.0;
0394     for ( int i = 0, N = hooks.size(); i < N; ++i )
0395       if ( hooks[i]->canBiasSelection() )
0396         f *= hooks[i]->biasedSelectionWeight();
0397     return f;
0398   }
0399 
0400   // Possibility to veto event after process-level selection.
0401   virtual bool canVetoProcessLevel() {
0402     for ( int i = 0, N = hooks.size(); i < N; ++i )
0403       if ( hooks[i]->canVetoProcessLevel() ) return true;
0404     return false;
0405   }
0406 
0407   // Decide whether to veto current process or not, based on process record.
0408   // Usage: doVetoProcessLevel( process).
0409   virtual bool doVetoProcessLevel(Event& e) {
0410     for ( int i = 0, N = hooks.size(); i < N; ++i )
0411       if ( hooks[i]->canVetoProcessLevel() &&
0412            hooks[i]->doVetoProcessLevel(e) ) return true;
0413     return false;
0414   }
0415 
0416   // Possibility to veto resonance decay chain.
0417   virtual bool canVetoResonanceDecays() {
0418     for ( int i = 0, N = hooks.size(); i < N; ++i )
0419       if ( hooks[i]->canVetoResonanceDecays() ) return true;
0420     return false;
0421   }
0422 
0423   // Decide whether to veto current resonance decay chain or not, based on
0424   // process record. Usage: doVetoProcessLevel( process).
0425   virtual bool doVetoResonanceDecays(Event& e) {
0426     for ( int i = 0, N = hooks.size(); i < N; ++i )
0427       if ( hooks[i]->canVetoResonanceDecays() &&
0428            hooks[i]->doVetoResonanceDecays(e) ) return true;
0429     return false;
0430   }
0431 
0432   // Possibility to veto MPI + ISR + FSR evolution and kill event,
0433   // making decision at a fixed pT scale. Useful for MLM-style matching.
0434   virtual bool canVetoPT() {
0435     for ( int i = 0, N = hooks.size(); i < N; ++i )
0436       if ( hooks[i]->canVetoPT() ) return true;
0437     return false;
0438   }
0439 
0440   // Transverse-momentum scale for veto test.
0441   virtual double scaleVetoPT() {
0442     double s = 0.0;
0443     for ( int i = 0, N = hooks.size(); i < N; ++i )
0444       if ( hooks[i]->canVetoPT() ) s = max(s, hooks[i]->scaleVetoPT());
0445     return s;
0446   }
0447 
0448   // Decide whether to veto current event or not, based on event record.
0449   // Usage: doVetoPT( iPos, event), where iPos = 0: no emissions so far;
0450   // iPos = 1/2/3 joint evolution, latest step was MPI/ISR/FSR;
0451   // iPos = 4: FSR only afterwards; iPos = 5: FSR in resonance decay.
0452   virtual bool doVetoPT( int iPos, const Event& e) {
0453     for ( int i = 0, N = hooks.size(); i < N; ++i )
0454       if ( hooks[i]->canVetoPT() && hooks[i]->doVetoPT(iPos, e) ) return true;
0455     return false;
0456   }
0457 
0458   // Possibility to veto MPI + ISR + FSR evolution and kill event,
0459   // making decision after fixed number of ISR or FSR steps.
0460   virtual bool canVetoStep() {
0461     for ( int i = 0, N = hooks.size(); i < N; ++i )
0462       if ( hooks[i]->canVetoStep() ) return true;
0463     return false;
0464   }
0465 
0466   // Up to how many ISR + FSR steps of hardest interaction should be checked.
0467   virtual int numberVetoStep() {
0468     int n = 1;
0469     for ( int i = 0, N = hooks.size(); i < N; ++i )
0470       if ( hooks[i]->canVetoStep() ) n = max(n, hooks[i]->numberVetoStep());
0471     return n;
0472   }
0473 
0474   // Decide whether to veto current event or not, based on event record.
0475   // Usage: doVetoStep( iPos, nISR, nFSR, event), where iPos as above,
0476   // nISR and nFSR number of emissions so far for hard interaction only.
0477   virtual bool doVetoStep( int iPos, int nISR, int nFSR, const Event& e) {
0478     for ( int i = 0, N = hooks.size(); i < N; ++i )
0479       if ( hooks[i]->canVetoStep()
0480         && hooks[i]->doVetoStep(iPos, nISR, nFSR, e) ) return true;
0481     return false;
0482   }
0483 
0484   // Possibility to veto MPI + ISR + FSR evolution and kill event,
0485   // making decision after fixed number of MPI steps.
0486   virtual bool canVetoMPIStep() {
0487     for ( int i = 0, N = hooks.size(); i < N; ++i )
0488       if ( hooks[i]->canVetoMPIStep() ) return true;
0489     return false;
0490   }
0491 
0492   // Up to how many MPI steps should be checked.
0493   virtual int numberVetoMPIStep() {
0494     int n = 1;
0495     for ( int i = 0, N = hooks.size(); i < N; ++i )
0496       if ( hooks[i]->canVetoMPIStep() )
0497         n = max(n, hooks[i]->numberVetoMPIStep());
0498     return n;
0499   }
0500 
0501   // Decide whether to veto current event or not, based on event record.
0502   // Usage: doVetoMPIStep( nMPI, event), where nMPI is number of MPI's so far.
0503   virtual bool doVetoMPIStep( int nMPI, const Event& e) {
0504     for ( int i = 0, N = hooks.size(); i < N; ++i )
0505       if ( hooks[i]->canVetoMPIStep() && hooks[i]->doVetoMPIStep(nMPI, e) )
0506         return true;
0507     return false;
0508   }
0509 
0510   // Possibility to veto event after ISR + FSR + MPI in parton level,
0511   // but before beam remnants and resonance decays.
0512   virtual bool canVetoPartonLevelEarly() {
0513     for ( int i = 0, N = hooks.size(); i < N; ++i )
0514       if ( hooks[i]->canVetoPartonLevelEarly() ) return true;
0515     return false;
0516   }
0517 
0518   // Decide whether to veto current partons or not, based on event record.
0519   // Usage: doVetoPartonLevelEarly( event).
0520   virtual bool doVetoPartonLevelEarly( const Event& e) {
0521     for ( int i = 0, N = hooks.size(); i < N; ++i )
0522       if ( hooks[i]->canVetoPartonLevelEarly()
0523         && hooks[i]->doVetoPartonLevelEarly(e) ) return true;
0524     return false;
0525   }
0526 
0527   // Retry same ProcessLevel with a new PartonLevel after a veto in
0528   // doVetoPT, doVetoStep, doVetoMPIStep or doVetoPartonLevelEarly
0529   // if you overload this method to return true.
0530   virtual bool retryPartonLevel() {
0531     for ( int i = 0, N = hooks.size(); i < N; ++i )
0532       if ( hooks[i]->retryPartonLevel() ) return true;
0533     return false;
0534   }
0535 
0536   // Possibility to veto event after parton-level selection.
0537   virtual bool canVetoPartonLevel() {
0538     for ( int i = 0, N = hooks.size(); i < N; ++i )
0539       if ( hooks[i]->canVetoPartonLevel() ) return true;
0540     return false;
0541   }
0542 
0543   // Decide whether to veto current partons or not, based on event record.
0544   // Usage: doVetoPartonLevel( event).
0545   virtual bool doVetoPartonLevel( const Event& e) {
0546     for ( int i = 0, N = hooks.size(); i < N; ++i )
0547       if ( hooks[i]->canVetoPartonLevel()
0548         && hooks[i]->doVetoPartonLevel(e) ) return true;
0549    return false;
0550   }
0551 
0552   // Possibility to set initial scale in TimeShower for resonance decay.
0553   virtual bool canSetResonanceScale() {
0554     for ( int i = 0, N = hooks.size(); i < N; ++i )
0555       if ( hooks[i]->canSetResonanceScale() ) return true;
0556     return false;
0557   }
0558 
0559   // Initial scale for TimeShower evolution.
0560   // Usage: scaleResonance( iRes, event), where iRes is location
0561   // of decaying resonance in the event record.
0562   virtual double scaleResonance( int iRes, const Event& e) {
0563     double s = 0.0;
0564     for ( int i = 0, N = hooks.size(); i < N; ++i )
0565       if ( hooks[i]->canSetResonanceScale() )
0566         s = max(s, hooks[i]->scaleResonance(iRes, e));
0567     return s;
0568   }
0569 
0570   // Possibility to veto an emission in the ISR machinery.
0571   virtual bool canVetoISREmission() {
0572     for ( int i = 0, N = hooks.size(); i < N; ++i )
0573       if ( hooks[i]->canVetoISREmission() ) return true;
0574     return false;
0575   }
0576 
0577   // Decide whether to veto current emission or not, based on event record.
0578   // Usage: doVetoISREmission( sizeOld, event, iSys) where sizeOld is size
0579   // of event record before current emission-to-be-scrutinized was added,
0580   // and iSys is the system of the radiation (according to PartonSystems).
0581   virtual bool doVetoISREmission( int sizeOld, const Event& e, int iSys) {
0582     for ( int i = 0, N = hooks.size(); i < N; ++i )
0583       if ( hooks[i]->canVetoISREmission()
0584         && hooks[i]->doVetoISREmission(sizeOld, e, iSys) ) return true;
0585     return false;
0586   }
0587 
0588   // Possibility to veto an emission in the FSR machinery.
0589   virtual bool canVetoFSREmission() {
0590     for ( int i = 0, N = hooks.size(); i < N; ++i )
0591       if ( hooks[i]->canVetoFSREmission() ) return true;
0592     return false;
0593   }
0594 
0595   // Decide whether to veto current emission or not, based on event record.
0596   // Usage: doVetoFSREmission( sizeOld, event, iSys, inResonance) where
0597   // sizeOld is size of event record before current emission-to-be-scrutinized
0598   // was added, iSys is the system of the radiation (according to
0599   // PartonSystems), and inResonance is true if the emission takes place in a
0600   // resonance decay.
0601   virtual bool doVetoFSREmission(int sizeOld, const Event& e,
0602     int iSys, bool inResonance = false ) {
0603     for ( int i = 0, N = hooks.size(); i < N; ++i )
0604       if ( hooks[i]->canVetoFSREmission()
0605         && hooks[i]->doVetoFSREmission(sizeOld, e, iSys, inResonance) )
0606         return true;
0607     return false;
0608   }
0609 
0610   // Possibility to veto an MPI.
0611   virtual bool canVetoMPIEmission() {
0612     for ( int i = 0, N = hooks.size(); i < N; ++i )
0613       if ( hooks[i]->canVetoMPIEmission() ) return true;
0614     return false;
0615   }
0616 
0617   // Decide whether to veto an MPI based on event record.
0618   // Usage: doVetoMPIEmission( sizeOld, event) where sizeOld
0619   // is size of event record before the current MPI.
0620   virtual bool doVetoMPIEmission( int sizeOld, const Event & e) {
0621     for ( int i = 0, N = hooks.size(); i < N; ++i )
0622       if ( hooks[i]->canVetoMPIEmission()
0623         && hooks[i]->doVetoMPIEmission(sizeOld, e) )
0624         return true;
0625     return false;
0626   }
0627 
0628   // Possibility to reconnect colours from resonance decay systems.
0629   virtual bool canReconnectResonanceSystems() {
0630     for ( int i = 0, N = hooks.size(); i < N; ++i )
0631       if ( hooks[i]->canReconnectResonanceSystems() ) return true;
0632     return false;
0633   }
0634 
0635   // Do reconnect colours from resonance decay systems.
0636   // Usage: doVetoFSREmission( oldSizeEvt, event)
0637   // where oldSizeEvent is the event size before resonance decays.
0638   // Should normally return true, while false means serious failure.
0639   // Value of PartonLevel:earlyResDec determines where method is called.
0640   virtual bool doReconnectResonanceSystems( int j, Event & e) {
0641     for ( int i = 0, N = hooks.size(); i < N; ++i )
0642       if ( hooks[i]->canReconnectResonanceSystems()
0643         && hooks[i]->doReconnectResonanceSystems(j, e) ) return true;
0644     return false;
0645   }
0646 
0647   // Can change fragmentation parameters.
0648   virtual bool canChangeFragPar() {
0649     for ( int i = 0, N = hooks.size(); i < N; ++i )
0650       if ( hooks[i]->canChangeFragPar() ) return true;
0651     return false;
0652   }
0653 
0654   // Set initial ends of a string to be fragmented. This is done once
0655   // for each string. Note that the second string end may be zero in
0656   // case we are hadronising a string piece leading to a junction.
0657   virtual void setStringEnds( const StringEnd* pos, const StringEnd* neg,
0658     vector<int> iPart) {
0659     for ( int i = 0, N = hooks.size(); i < N; ++i )
0660       hooks[i]->setStringEnds( pos, neg, iPart);
0661   }
0662 
0663   // Do change fragmentation parameters.
0664   // Input: flavPtr, zPtr, pTPtr, idEnd, m2Had, iParton and posEnd (or
0665   // negEnd).
0666   virtual bool doChangeFragPar( StringFlav* sfIn, StringZ* zIn,
0667     StringPT* ptIn, int idIn, double mIn, vector<int> parIn,
0668     const StringEnd* endIn) {
0669     for ( int i = 0, N = hooks.size(); i < N; ++i )
0670       if ( hooks[i]->canChangeFragPar()
0671            && hooks[i]->doChangeFragPar(sfIn, zIn, ptIn, idIn,
0672                                      mIn, parIn, endIn) ) return true;
0673     return false;}
0674 
0675   // Can veto hadrons in the fragmentation.
0676   virtual bool canVetoFragmentation() {
0677     for ( int i = 0, N = hooks.size(); i < N; ++i )
0678       if ( hooks[i]->canVetoFragmentation() ) return true;
0679     return false;}
0680 
0681   // Do a veto on a hadron just before it is added to the final state.
0682   virtual bool doVetoFragmentation(Particle p, const StringEnd* nowEnd) {
0683     for ( int i = 0, N = hooks.size(); i < N; ++i )
0684       if ( hooks[i]->canChangeFragPar()
0685         && hooks[i]->doVetoFragmentation(p, nowEnd) ) return true;
0686     return false;
0687   }
0688 
0689   virtual bool doVetoFragmentation(Particle p1, Particle p2,
0690     const StringEnd* e1, const StringEnd* e2) {
0691     for ( int i = 0, N = hooks.size(); i < N; ++i )
0692       if ( hooks[i]->canChangeFragPar()
0693         && hooks[i]->doVetoFragmentation(p1, p2, e1, e2) ) return true;
0694     return false;
0695   }
0696 
0697   virtual bool doVetoFinalTwo(Particle p1, Particle p2,
0698     const StringEnd* e1, const StringEnd* e2) {
0699     for ( int i = 0, N = hooks.size(); i < N; ++i )
0700       if ( hooks[i]->canChangeFragPar()
0701         && hooks[i]->doVetoFinalTwo(p1, p2, e1, e2) ) return true;
0702     return false;
0703   }
0704 
0705   // Possibility to veto an event after hadronization based on event
0706   // contents. Works as an early trigger to avoid running the time
0707   // consuming rescattering process on uninteresting events.
0708   virtual bool canVetoAfterHadronization() {
0709     for ( int i = 0, N = hooks.size(); i < N; ++i )
0710       if ( hooks[i]->canVetoAfterHadronization() ) return true;
0711     return false;
0712   }
0713 
0714   // Do the actual veto after hadronization.
0715   virtual bool doVetoAfterHadronization(const Event& e) {
0716     for ( int i = 0, N = hooks.size(); i < N; ++i )
0717       if ( hooks[i]->canVetoAfterHadronization()
0718         && hooks[i]->doVetoAfterHadronization(e) ) return true;
0719     return false;
0720   }
0721 
0722   // Can set the overall impact parameter for the MPI treatment.
0723   virtual bool canSetImpactParameter() const {
0724     for ( int i = 0, N = hooks.size(); i < N; ++i )
0725       if ( hooks[i]->canSetImpactParameter() ) return true;
0726     return false;
0727   }
0728 
0729   // Set the overall impact parameter for the MPI treatment.
0730   virtual double doSetImpactParameter() {
0731     for ( int i = 0, N = hooks.size(); i < N; ++i )
0732       if ( hooks[i]->canSetImpactParameter() )
0733         return hooks[i]->doSetImpactParameter();
0734     return 0.0;
0735   }
0736 
0737   // The vector of user hooks.
0738   vector< shared_ptr<UserHooks> > hooks = {};
0739 
0740 };
0741 
0742 //==========================================================================
0743 
0744 } // end namespace Pythia8
0745 
0746 #endif // Pythia8_UserHooks_H