Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 10:06:33

0001 // UserHooks.h is a part of the PYTHIA event generator.
0002 // Copyright (C) 2024 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   // Do a veto on a hadron just before it is added to the final state.
0201   // The StringEnd from which the the hadron was produced is included
0202   // for information.
0203   virtual bool doVetoFragmentation( Particle, const StringEnd *)
0204     { return false;}
0205 
0206   // Do a veto on a hadron just before it is added to the final state
0207   // (final two hadron case).
0208   virtual bool doVetoFragmentation(Particle, Particle,
0209     const StringEnd*, const StringEnd* ) { return false;}
0210 
0211   // Possibility to veto an event after hadronization based
0212   // on event contents. Works as an early trigger to avoid
0213   // running the time consuming rescattering process on
0214   // uninteresting events.
0215   virtual bool canVetoAfterHadronization() {return false;}
0216 
0217   // Do the actual veto after hadronization.
0218   virtual bool doVetoAfterHadronization(const Event& ) {return false;}
0219 
0220   // Can set the overall impact parameter for the MPI treatment.
0221   virtual bool canSetImpactParameter() const { return false; }
0222 
0223   // Set the overall impact parameter for the MPI treatment.
0224   virtual double doSetImpactParameter() { return 0.0; }
0225 
0226   // Custom processing at the end of HadronLevel::next.
0227   virtual bool onEndHadronLevel(HadronLevel&, Event&) { return true; }
0228 
0229 protected:
0230 
0231   // Constructor.
0232   UserHooks() {}
0233 
0234   // After initInfoPtr, initialize workEvent
0235   virtual void onInitInfoPtr() override {
0236     // Set smart pointer to null, in order to avoid circular dependency
0237     userHooksPtr = nullptr;
0238     workEvent.init("(work event)", particleDataPtr);
0239   }
0240 
0241   // omitResonanceDecays omits resonance decay chains from process record.
0242   void omitResonanceDecays(const Event& process, bool finalOnly = false);
0243 
0244   // subEvent extracts currently resolved partons in the hard process.
0245   void subEvent(const Event& event, bool isHardest = true);
0246 
0247   // Have one event object around as work area.
0248   Event workEvent = {};
0249 
0250   // User-imposed selection bias.
0251   double selBias = 1.;
0252 
0253   // Bookkept quantities for boosted event weights.
0254   double enhancedEventWeight = {}, pTEnhanced = {}, wtEnhanced = {};
0255 
0256 };
0257 
0258 //==========================================================================
0259 
0260 // SuppressSmallPT is a derived class for user access to program execution.
0261 // It is a simple example, illustrating how to suppress the cross section
0262 // of 2 -> 2 processes by a factor pT^4 / (pT0^2 + pT^2)^2, with pT0 input,
0263 // and also modify alpha_strong scale similarly.
0264 
0265 class SuppressSmallPT : public UserHooks {
0266 
0267 public:
0268 
0269   // Constructor.
0270   SuppressSmallPT( double pT0timesMPIIn = 1., int numberAlphaSIn = 0,
0271     bool useSameAlphaSasMPIIn = true) : pT20(0.) {isInit = false;
0272     pT0timesMPI = pT0timesMPIIn; numberAlphaS = numberAlphaSIn;
0273     useSameAlphaSasMPI = useSameAlphaSasMPIIn;}
0274 
0275   // Possibility to modify cross section of process.
0276   virtual bool canModifySigma() {return true;}
0277 
0278   // Multiplicative factor modifying the cross section of a hard process.
0279   // Usage: inEvent is true for event generation, false for initialization.
0280   virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
0281     const PhaseSpace* phaseSpacePtr, bool );
0282 
0283 private:
0284 
0285   // Save input properties and the squared pT0 scale.
0286   bool   isInit, useSameAlphaSasMPI;
0287   int    numberAlphaS;
0288   double pT0timesMPI, pT20;
0289 
0290   // Alpha_strong calculation.
0291   AlphaStrong alphaS;
0292 
0293 };
0294 
0295 //==========================================================================
0296 
0297 // UserHooksVector implements a vector of UserHooks and is itself a UserHooks.
0298 
0299 class UserHooksVector: public UserHooks {
0300 
0301 public:
0302 
0303   // The default constructor is private, and should only be used
0304   // internally in Pythia.
0305   UserHooksVector() {};
0306   friend class Pythia;
0307 
0308   // Destructor.
0309   virtual ~UserHooksVector() {}
0310 
0311   // Initialisation after beams have been set by Pythia::init().
0312   // Check that there are no (obvious) clashes.
0313   virtual bool initAfterBeams() {
0314     int nCanSetResonanceScale  = 0;
0315     int nCanChangeFragPar      = 0;
0316     int nCanSetImpactParameter = 0;
0317     for ( int i = 0, N = hooks.size(); i < N; ++i ) {
0318       registerSubObject(*hooks[i]);
0319       if ( !hooks[i]->initAfterBeams() ) return false;
0320       if (hooks[i]->canSetResonanceScale()) ++nCanSetResonanceScale;
0321       if (hooks[i]->canChangeFragPar()) ++nCanChangeFragPar;
0322       if (hooks[i]->canSetImpactParameter()) ++nCanSetImpactParameter;
0323     }
0324     if (nCanSetResonanceScale > 1) {
0325       loggerPtr->ERROR_MSG(
0326         "multiple UserHooks with canSetResonanceScale() not allowed");
0327       return false;
0328     }
0329     if (nCanChangeFragPar > 1) {
0330       loggerPtr->ERROR_MSG(
0331         "multiple UserHooks with canChangeFragPar() not allowed");
0332       return false;
0333     }
0334     if (nCanSetImpactParameter > 1) {
0335       loggerPtr->ERROR_MSG(
0336         "multiple UserHooks with canSetImpactParameter() not allowed");
0337       return false;
0338     }
0339     return true;
0340   }
0341 
0342   // Possibility to modify cross section of process.
0343   virtual bool canModifySigma() {
0344     for ( int i = 0, N = hooks.size(); i < N; ++i )
0345       if ( hooks[i]->canModifySigma() ) return true;
0346     return false;
0347   }
0348 
0349   // Multiplicative factor modifying the cross section of a hard process.
0350   virtual double multiplySigmaBy(const SigmaProcess* sigmaProcessPtr,
0351     const PhaseSpace* phaseSpacePtr, bool inEvent) {
0352     double f = 1.0;
0353     for ( int i = 0, N = hooks.size(); i < N; ++i )
0354       if ( hooks[i]->canModifySigma() )
0355         f *= hooks[i]->multiplySigmaBy(sigmaProcessPtr, phaseSpacePtr,inEvent);
0356     return f;
0357   }
0358 
0359   // Possibility to bias selection of events, compensated by a weight.
0360   virtual bool canBiasSelection() {
0361     for ( int i = 0, N = hooks.size(); i < N; ++i )
0362       if ( hooks[i]->canBiasSelection() ) return true;
0363     return false;
0364   }
0365 
0366   // Multiplicative factor in the phase space selection of a hard process.
0367   virtual double biasSelectionBy(const SigmaProcess* sigmaProcessPtr,
0368     const PhaseSpace* phaseSpacePtr, bool inEvent) {
0369     double f = 1.0;
0370     for ( int i = 0, N = hooks.size(); i < N; ++i )
0371       if ( hooks[i]->canBiasSelection() )
0372         f *= hooks[i]->biasSelectionBy(sigmaProcessPtr, phaseSpacePtr,
0373              inEvent);
0374     return f;
0375   }
0376 
0377   // Event weight to compensate for selection weight above.
0378   virtual double biasedSelectionWeight() {
0379     double f = 1.0;
0380     for ( int i = 0, N = hooks.size(); i < N; ++i )
0381       if ( hooks[i]->canBiasSelection() )
0382         f *= hooks[i]->biasedSelectionWeight();
0383     return f;
0384   }
0385 
0386   // Possibility to veto event after process-level selection.
0387   virtual bool canVetoProcessLevel() {
0388     for ( int i = 0, N = hooks.size(); i < N; ++i )
0389       if ( hooks[i]->canVetoProcessLevel() ) return true;
0390     return false;
0391   }
0392 
0393   // Decide whether to veto current process or not, based on process record.
0394   // Usage: doVetoProcessLevel( process).
0395   virtual bool doVetoProcessLevel(Event& e) {
0396     for ( int i = 0, N = hooks.size(); i < N; ++i )
0397       if ( hooks[i]->canVetoProcessLevel() &&
0398            hooks[i]->doVetoProcessLevel(e) ) return true;
0399     return false;
0400   }
0401 
0402   // Possibility to veto resonance decay chain.
0403   virtual bool canVetoResonanceDecays() {
0404     for ( int i = 0, N = hooks.size(); i < N; ++i )
0405       if ( hooks[i]->canVetoResonanceDecays() ) return true;
0406     return false;
0407   }
0408 
0409   // Decide whether to veto current resonance decay chain or not, based on
0410   // process record. Usage: doVetoProcessLevel( process).
0411   virtual bool doVetoResonanceDecays(Event& e) {
0412     for ( int i = 0, N = hooks.size(); i < N; ++i )
0413       if ( hooks[i]->canVetoResonanceDecays() &&
0414            hooks[i]->doVetoResonanceDecays(e) ) return true;
0415     return false;
0416   }
0417 
0418   // Possibility to veto MPI + ISR + FSR evolution and kill event,
0419   // making decision at a fixed pT scale. Useful for MLM-style matching.
0420   virtual bool canVetoPT() {
0421     for ( int i = 0, N = hooks.size(); i < N; ++i )
0422       if ( hooks[i]->canVetoPT() ) return true;
0423     return false;
0424   }
0425 
0426   // Transverse-momentum scale for veto test.
0427   virtual double scaleVetoPT() {
0428     double s = 0.0;
0429     for ( int i = 0, N = hooks.size(); i < N; ++i )
0430       if ( hooks[i]->canVetoPT() ) s = max(s, hooks[i]->scaleVetoPT());
0431     return s;
0432   }
0433 
0434   // Decide whether to veto current event or not, based on event record.
0435   // Usage: doVetoPT( iPos, event), where iPos = 0: no emissions so far;
0436   // iPos = 1/2/3 joint evolution, latest step was MPI/ISR/FSR;
0437   // iPos = 4: FSR only afterwards; iPos = 5: FSR in resonance decay.
0438   virtual bool doVetoPT( int iPos, const Event& e) {
0439     for ( int i = 0, N = hooks.size(); i < N; ++i )
0440       if ( hooks[i]->canVetoPT() && hooks[i]->doVetoPT(iPos, e) ) return true;
0441     return false;
0442   }
0443 
0444   // Possibility to veto MPI + ISR + FSR evolution and kill event,
0445   // making decision after fixed number of ISR or FSR steps.
0446   virtual bool canVetoStep() {
0447     for ( int i = 0, N = hooks.size(); i < N; ++i )
0448       if ( hooks[i]->canVetoStep() ) return true;
0449     return false;
0450   }
0451 
0452   // Up to how many ISR + FSR steps of hardest interaction should be checked.
0453   virtual int numberVetoStep() {
0454     int n = 1;
0455     for ( int i = 0, N = hooks.size(); i < N; ++i )
0456       if ( hooks[i]->canVetoStep() ) n = max(n, hooks[i]->numberVetoStep());
0457     return n;
0458   }
0459 
0460   // Decide whether to veto current event or not, based on event record.
0461   // Usage: doVetoStep( iPos, nISR, nFSR, event), where iPos as above,
0462   // nISR and nFSR number of emissions so far for hard interaction only.
0463   virtual bool doVetoStep( int iPos, int nISR, int nFSR, const Event& e) {
0464     for ( int i = 0, N = hooks.size(); i < N; ++i )
0465       if ( hooks[i]->canVetoStep()
0466         && hooks[i]->doVetoStep(iPos, nISR, nFSR, e) ) return true;
0467     return false;
0468   }
0469 
0470   // Possibility to veto MPI + ISR + FSR evolution and kill event,
0471   // making decision after fixed number of MPI steps.
0472   virtual bool canVetoMPIStep() {
0473     for ( int i = 0, N = hooks.size(); i < N; ++i )
0474       if ( hooks[i]->canVetoMPIStep() ) return true;
0475     return false;
0476   }
0477 
0478   // Up to how many MPI steps should be checked.
0479   virtual int numberVetoMPIStep() {
0480     int n = 1;
0481     for ( int i = 0, N = hooks.size(); i < N; ++i )
0482       if ( hooks[i]->canVetoMPIStep() )
0483         n = max(n, hooks[i]->numberVetoMPIStep());
0484     return n;
0485   }
0486 
0487   // Decide whether to veto current event or not, based on event record.
0488   // Usage: doVetoMPIStep( nMPI, event), where nMPI is number of MPI's so far.
0489   virtual bool doVetoMPIStep( int nMPI, const Event& e) {
0490     for ( int i = 0, N = hooks.size(); i < N; ++i )
0491       if ( hooks[i]->canVetoMPIStep() && hooks[i]->doVetoMPIStep(nMPI, e) )
0492         return true;
0493     return false;
0494   }
0495 
0496   // Possibility to veto event after ISR + FSR + MPI in parton level,
0497   // but before beam remnants and resonance decays.
0498   virtual bool canVetoPartonLevelEarly() {
0499     for ( int i = 0, N = hooks.size(); i < N; ++i )
0500       if ( hooks[i]->canVetoPartonLevelEarly() ) return true;
0501     return false;
0502   }
0503 
0504   // Decide whether to veto current partons or not, based on event record.
0505   // Usage: doVetoPartonLevelEarly( event).
0506   virtual bool doVetoPartonLevelEarly( const Event& e) {
0507     for ( int i = 0, N = hooks.size(); i < N; ++i )
0508       if ( hooks[i]->canVetoPartonLevelEarly()
0509         && hooks[i]->doVetoPartonLevelEarly(e) ) return true;
0510     return false;
0511   }
0512 
0513   // Retry same ProcessLevel with a new PartonLevel after a veto in
0514   // doVetoPT, doVetoStep, doVetoMPIStep or doVetoPartonLevelEarly
0515   // if you overload this method to return true.
0516   virtual bool retryPartonLevel() {
0517     for ( int i = 0, N = hooks.size(); i < N; ++i )
0518       if ( hooks[i]->retryPartonLevel() ) return true;
0519     return false;
0520   }
0521 
0522   // Possibility to veto event after parton-level selection.
0523   virtual bool canVetoPartonLevel() {
0524     for ( int i = 0, N = hooks.size(); i < N; ++i )
0525       if ( hooks[i]->canVetoPartonLevel() ) return true;
0526     return false;
0527   }
0528 
0529   // Decide whether to veto current partons or not, based on event record.
0530   // Usage: doVetoPartonLevel( event).
0531   virtual bool doVetoPartonLevel( const Event& e) {
0532     for ( int i = 0, N = hooks.size(); i < N; ++i )
0533       if ( hooks[i]->canVetoPartonLevel()
0534         && hooks[i]->doVetoPartonLevel(e) ) return true;
0535    return false;
0536   }
0537 
0538   // Possibility to set initial scale in TimeShower for resonance decay.
0539   virtual bool canSetResonanceScale() {
0540     for ( int i = 0, N = hooks.size(); i < N; ++i )
0541       if ( hooks[i]->canSetResonanceScale() ) return true;
0542     return false;
0543   }
0544 
0545   // Initial scale for TimeShower evolution.
0546   // Usage: scaleResonance( iRes, event), where iRes is location
0547   // of decaying resonance in the event record.
0548   virtual double scaleResonance( int iRes, const Event& e) {
0549     double s = 0.0;
0550     for ( int i = 0, N = hooks.size(); i < N; ++i )
0551       if ( hooks[i]->canSetResonanceScale() )
0552         s = max(s, hooks[i]->scaleResonance(iRes, e));
0553     return s;
0554   }
0555 
0556   // Possibility to veto an emission in the ISR machinery.
0557   virtual bool canVetoISREmission() {
0558     for ( int i = 0, N = hooks.size(); i < N; ++i )
0559       if ( hooks[i]->canVetoISREmission() ) return true;
0560     return false;
0561   }
0562 
0563   // Decide whether to veto current emission or not, based on event record.
0564   // Usage: doVetoISREmission( sizeOld, event, iSys) where sizeOld is size
0565   // of event record before current emission-to-be-scrutinized was added,
0566   // and iSys is the system of the radiation (according to PartonSystems).
0567   virtual bool doVetoISREmission( int sizeOld, const Event& e, int iSys) {
0568     for ( int i = 0, N = hooks.size(); i < N; ++i )
0569       if ( hooks[i]->canVetoISREmission()
0570         && hooks[i]->doVetoISREmission(sizeOld, e, iSys) ) return true;
0571     return false;
0572   }
0573 
0574   // Possibility to veto an emission in the FSR machinery.
0575   virtual bool canVetoFSREmission() {
0576     for ( int i = 0, N = hooks.size(); i < N; ++i )
0577       if ( hooks[i]->canVetoFSREmission() ) return true;
0578     return false;
0579   }
0580 
0581   // Decide whether to veto current emission or not, based on event record.
0582   // Usage: doVetoFSREmission( sizeOld, event, iSys, inResonance) where
0583   // sizeOld is size of event record before current emission-to-be-scrutinized
0584   // was added, iSys is the system of the radiation (according to
0585   // PartonSystems), and inResonance is true if the emission takes place in a
0586   // resonance decay.
0587   virtual bool doVetoFSREmission(int sizeOld, const Event& e,
0588     int iSys, bool inResonance = false ) {
0589     for ( int i = 0, N = hooks.size(); i < N; ++i )
0590       if ( hooks[i]->canVetoFSREmission()
0591         && hooks[i]->doVetoFSREmission(sizeOld, e, iSys, inResonance) )
0592         return true;
0593     return false;
0594   }
0595 
0596   // Possibility to veto an MPI.
0597   virtual bool canVetoMPIEmission() {
0598     for ( int i = 0, N = hooks.size(); i < N; ++i )
0599       if ( hooks[i]->canVetoMPIEmission() ) return true;
0600     return false;
0601   }
0602 
0603   // Decide whether to veto an MPI based on event record.
0604   // Usage: doVetoMPIEmission( sizeOld, event) where sizeOld
0605   // is size of event record before the current MPI.
0606   virtual bool doVetoMPIEmission( int sizeOld, const Event & e) {
0607     for ( int i = 0, N = hooks.size(); i < N; ++i )
0608       if ( hooks[i]->canVetoMPIEmission()
0609         && hooks[i]->doVetoMPIEmission(sizeOld, e) )
0610         return true;
0611     return false;
0612   }
0613 
0614   // Possibility to reconnect colours from resonance decay systems.
0615   virtual bool canReconnectResonanceSystems() {
0616     for ( int i = 0, N = hooks.size(); i < N; ++i )
0617       if ( hooks[i]->canReconnectResonanceSystems() ) return true;
0618     return false;
0619   }
0620 
0621   // Do reconnect colours from resonance decay systems.
0622   // Usage: doVetoFSREmission( oldSizeEvt, event)
0623   // where oldSizeEvent is the event size before resonance decays.
0624   // Should normally return true, while false means serious failure.
0625   // Value of PartonLevel:earlyResDec determines where method is called.
0626   virtual bool doReconnectResonanceSystems( int j, Event & e) {
0627     for ( int i = 0, N = hooks.size(); i < N; ++i )
0628       if ( hooks[i]->canReconnectResonanceSystems()
0629         && hooks[i]->doReconnectResonanceSystems(j, e) ) return true;
0630     return false;
0631   }
0632 
0633   // Can change fragmentation parameters.
0634   virtual bool canChangeFragPar() {
0635     for ( int i = 0, N = hooks.size(); i < N; ++i )
0636       if ( hooks[i]->canChangeFragPar() ) return true;
0637     return false;
0638   }
0639 
0640   // Set initial ends of a string to be fragmented. This is done once
0641   // for each string. Note that the second string end may be zero in
0642   // case we are hadronising a string piece leading to a junction.
0643   virtual void setStringEnds( const StringEnd* pos, const StringEnd* neg,
0644     vector<int> iPart) {
0645     for ( int i = 0, N = hooks.size(); i < N; ++i )
0646       hooks[i]->setStringEnds( pos, neg, iPart);
0647   }
0648 
0649   // Do change fragmentation parameters.
0650   // Input: flavPtr, zPtr, pTPtr, idEnd, m2Had, iParton and posEnd (or
0651   // negEnd).
0652   virtual bool doChangeFragPar( StringFlav* sfIn, StringZ* zIn,
0653     StringPT* ptIn, int idIn, double mIn, vector<int> parIn,
0654     const StringEnd* endIn) {
0655     for ( int i = 0, N = hooks.size(); i < N; ++i )
0656       if ( hooks[i]->canChangeFragPar()
0657            && hooks[i]->doChangeFragPar(sfIn, zIn, ptIn, idIn,
0658                                      mIn, parIn, endIn) ) return true;
0659     return false;}
0660 
0661   // Do a veto on a hadron just before it is added to the final state.
0662   virtual bool doVetoFragmentation(Particle p, const StringEnd* nowEnd) {
0663     for ( int i = 0, N = hooks.size(); i < N; ++i )
0664       if ( hooks[i]->canChangeFragPar()
0665         && hooks[i]->doVetoFragmentation(p, nowEnd) ) return true;
0666     return false;
0667   }
0668 
0669   virtual bool doVetoFragmentation(Particle p1, Particle p2,
0670     const StringEnd* e1, const StringEnd* e2) {
0671     for ( int i = 0, N = hooks.size(); i < N; ++i )
0672       if ( hooks[i]->canChangeFragPar()
0673         && hooks[i]->doVetoFragmentation(p1, p2, e1, e2) ) return true;
0674     return false;
0675   }
0676 
0677   // Possibility to veto an event after hadronization based on event
0678   // contents. Works as an early trigger to avoid running the time
0679   // consuming rescattering process on uninteresting events.
0680   virtual bool canVetoAfterHadronization() {
0681     for ( int i = 0, N = hooks.size(); i < N; ++i )
0682       if ( hooks[i]->canVetoAfterHadronization() ) return true;
0683     return false;
0684   }
0685 
0686   // Do the actual veto after hadronization.
0687   virtual bool doVetoAfterHadronization(const Event& e) {
0688     for ( int i = 0, N = hooks.size(); i < N; ++i )
0689       if ( hooks[i]->canVetoAfterHadronization()
0690         && hooks[i]->doVetoAfterHadronization(e) ) return true;
0691     return false;
0692   }
0693 
0694   // Can set the overall impact parameter for the MPI treatment.
0695   virtual bool canSetImpactParameter() const {
0696     for ( int i = 0, N = hooks.size(); i < N; ++i )
0697       if ( hooks[i]->canSetImpactParameter() ) return true;
0698     return false;
0699   }
0700 
0701   // Set the overall impact parameter for the MPI treatment.
0702   virtual double doSetImpactParameter() {
0703     for ( int i = 0, N = hooks.size(); i < N; ++i )
0704       if ( hooks[i]->canSetImpactParameter() )
0705         return hooks[i]->doSetImpactParameter();
0706     return 0.0;
0707   }
0708 
0709   // The vector of user hooks.
0710   vector< shared_ptr<UserHooks> > hooks = {};
0711 
0712 };
0713 
0714 //==========================================================================
0715 
0716 } // end namespace Pythia8
0717 
0718 #endif // Pythia8_UserHooks_H