PYTHIA  8.312
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
UserHooks Class Reference

UserHooks is base class for user access to program execution. More...

#include <UserHooks.h>

Inheritance diagram for UserHooks:
PhysicsBase AlpgenHooks amcnlo_unitarised_interface HeavyIons::InfoGrabber JetMatching MBReconUserHooks MergeResScaleHook PowhegHooks ResonanceDecayFilterHook SetLHEDecayProductHook SuppressSmallPT TopReconUserHooks UserHooksVector VinciaDiagnostics VinciaEWVetoHook

Public Member Functions

virtual ~UserHooks ()
 Destructor.
 
virtual bool initAfterBeams ()
 Initialisation after beams have been set by Pythia::init().
 
virtual bool canModifySigma ()
 Possibility to modify cross section of process.
 
virtual double multiplySigmaBy (const SigmaProcess *sigmaProcessPtr, const PhaseSpace *phaseSpacePtr, bool inEvent)
 Multiplicative factor modifying the cross section of a hard process. More...
 
virtual bool canBiasSelection ()
 Possibility to bias selection of events, compensated by a weight.
 
virtual double biasSelectionBy (const SigmaProcess *sigmaProcessPtr, const PhaseSpace *phaseSpacePtr, bool inEvent)
 Multiplicative factor in the phase space selection of a hard process. More...
 
virtual double biasedSelectionWeight ()
 Event weight to compensate for selection weight above.
 
virtual bool canVetoProcessLevel ()
 Possibility to veto event after process-level selection.
 
virtual bool doVetoProcessLevel (Event &)
 
virtual bool canSetLowEnergySigma (int, int) const
 
virtual double doSetLowEnergySigma (int, int, double, double, double) const
 
virtual bool canVetoResonanceDecays ()
 Possibility to veto resonance decay chain.
 
virtual bool doVetoResonanceDecays (Event &)
 
virtual bool canVetoPT ()
 
virtual double scaleVetoPT ()
 Transverse-momentum scale for veto test.
 
virtual bool doVetoPT (int, const Event &)
 
virtual bool canVetoStep ()
 
virtual int numberVetoStep ()
 Up to how many ISR + FSR steps of hardest interaction should be checked.
 
virtual bool doVetoStep (int, int, int, const Event &)
 
virtual bool canVetoMPIStep ()
 
virtual int numberVetoMPIStep ()
 Up to how many MPI steps should be checked.
 
virtual bool doVetoMPIStep (int, const Event &)
 
virtual bool canVetoPartonLevelEarly ()
 
virtual bool doVetoPartonLevelEarly (const Event &)
 
virtual bool retryPartonLevel ()
 
virtual bool canVetoPartonLevel ()
 Possibility to veto event after parton-level selection.
 
virtual bool doVetoPartonLevel (const Event &)
 
virtual bool canSetResonanceScale ()
 Possibility to set initial scale in TimeShower for resonance decay.
 
virtual double scaleResonance (int, const Event &)
 
virtual bool canVetoISREmission ()
 Possibility to veto an emission in the ISR machinery.
 
virtual bool doVetoISREmission (int, const Event &, int)
 
virtual bool canVetoFSREmission ()
 Possibility to veto an emission in the FSR machinery.
 
virtual bool doVetoFSREmission (int, const Event &, int, bool=false)
 
virtual bool canVetoMPIEmission ()
 Possibility to veto an MPI.
 
virtual bool doVetoMPIEmission (int, const Event &)
 
virtual bool canReconnectResonanceSystems ()
 Possibility to reconnect colours from resonance decay systems.
 
virtual bool doReconnectResonanceSystems (int, Event &)
 
virtual bool canChangeFragPar ()
 Can change fragmentation parameters.
 
virtual void setStringEnds (const StringEnd *, const StringEnd *, vector< int >)
 
virtual bool doChangeFragPar (StringFlav *, StringZ *, StringPT *, int, double, vector< int >, const StringEnd *)
 
virtual bool doVetoFragmentation (Particle, const StringEnd *)
 
virtual bool doVetoFragmentation (Particle, Particle, const StringEnd *, const StringEnd *)
 
virtual bool canVetoAfterHadronization ()
 
virtual bool doVetoAfterHadronization (const Event &)
 Do the actual veto after hadronization.
 
virtual bool canSetImpactParameter () const
 Can set the overall impact parameter for the MPI treatment.
 
virtual double doSetImpactParameter ()
 Set the overall impact parameter for the MPI treatment.
 
virtual bool onEndHadronLevel (HadronLevel &, Event &)
 Custom processing at the end of HadronLevel::next.
 
- Public Member Functions inherited from PhysicsBase
void initInfoPtr (Info &infoPtrIn)
 This function is called from above for physics objects used in a run. More...
 
virtual ~PhysicsBase ()
 Empty virtual destructor.
 
bool flag (string key) const
 Shorthand to read settings values.
 
int mode (string key) const
 
double parm (string key) const
 
string word (string key) const
 
vector< bool > fvec (string key) const
 
vector< int > mvec (string key) const
 
vector< double > pvec (string key) const
 
vector< string > wvec (string key) const
 

Protected Member Functions

 UserHooks ()
 Constructor.
 
virtual void onInitInfoPtr () override
 After initInfoPtr, initialize workEvent. More...
 
void omitResonanceDecays (const Event &process, bool finalOnly=false)
 omitResonanceDecays omits resonance decay chains from process record. More...
 
void subEvent (const Event &event, bool isHardest=true)
 subEvent extracts currently resolved partons in the hard process. More...
 
- Protected Member Functions inherited from PhysicsBase
 PhysicsBase ()
 Default constructor.
 
virtual void onBeginEvent ()
 This function is called in the very beginning of each Pythia::next call.
 
virtual void onEndEvent (Status)
 
virtual void onStat ()
 This function is called from the Pythia::stat() call.
 
void registerSubObject (PhysicsBase &pb)
 Register a sub object that should have its information in sync with this.
 

Protected Attributes

Event workEvent = {}
 Have one event object around as work area.
 
double selBias = 1.
 User-imposed selection bias.
 
double enhancedEventWeight = {}
 Bookkept quantities for boosted event weights.
 
double pTEnhanced = {}
 
double wtEnhanced = {}
 
- Protected Attributes inherited from PhysicsBase
InfoinfoPtr = {}
 
SettingssettingsPtr = {}
 Pointer to the settings database.
 
ParticleDataparticleDataPtr = {}
 Pointer to the particle data table.
 
LoggerloggerPtr = {}
 Pointer to logger.
 
HadronWidthshadronWidthsPtr = {}
 Pointer to the hadron widths data table.
 
RndmrndmPtr = {}
 Pointer to the random number generator.
 
CoupSMcoupSMPtr = {}
 Pointers to SM and SUSY couplings.
 
CoupSUSYcoupSUSYPtr = {}
 
BeamSetupbeamSetupPtr = {}
 
BeamParticlebeamAPtr = {}
 
BeamParticlebeamBPtr = {}
 
BeamParticlebeamPomAPtr = {}
 
BeamParticlebeamPomBPtr = {}
 
BeamParticlebeamGamAPtr = {}
 
BeamParticlebeamGamBPtr = {}
 
BeamParticlebeamVMDAPtr = {}
 
BeamParticlebeamVMDBPtr = {}
 
PartonSystemspartonSystemsPtr = {}
 Pointer to information on subcollision parton locations.
 
SigmaTotalsigmaTotPtr = {}
 Pointers to the total/elastic/diffractive cross sections.
 
SigmaCombinedsigmaCmbPtr = {}
 
set< PhysicsBase * > subObjects
 
UserHooksPtr userHooksPtr
 

Additional Inherited Members

- Public Types inherited from PhysicsBase
enum  Status {
  INCOMPLETE = -1, COMPLETE = 0, CONSTRUCTOR_FAILED, INIT_FAILED,
  LHEF_END, LOWENERGY_FAILED, PROCESSLEVEL_FAILED, PROCESSLEVEL_USERVETO,
  MERGING_FAILED, PARTONLEVEL_FAILED, PARTONLEVEL_USERVETO, HADRONLEVEL_FAILED,
  CHECK_FAILED, OTHER_UNPHYSICAL, HEAVYION_FAILED, HADRONLEVEL_USERVETO
}
 Enumerate the different status codes the event generation can have.
 

Detailed Description

UserHooks is base class for user access to program execution.

Member Function Documentation

double biasSelectionBy ( const SigmaProcess sigmaProcessPtr,
const PhaseSpace phaseSpacePtr,
bool  inEvent 
)
virtual

Multiplicative factor in the phase space selection of a hard process.

biasSelectionBy allows the user to introduce a multiplicative factor that modifies the cross section of a hard process. The event is assigned a wegith that is the inverse of the selection bias, such that the cross section is unchanged. Since it is called from before the event record is generated in full, the normal analysis does not work. The code here provides a rather extensive summary of which methods actually do work. It is a convenient starting point for writing your own derived routine.

Process code, necessary when some to be treated differently. int code = sigmaProcessPtr->code();

Final multiplicity, i.e. whether 2 -> 1 or 2 -> 2. int nFinal = sigmaProcessPtr->nFinal();

Incoming x1 and x2 to the hard collision, and factorization scale. double x1 = phaseSpacePtr->x1(); double x2 = phaseSpacePtr->x2(); double Q2Fac = sigmaProcessPtr->Q2Fac();

Renormalization scale and assumed alpha_strong and alpha_EM. double Q2Ren = sigmaProcessPtr->Q2Ren(); double alphaS = sigmaProcessPtr->alphaSRen(); double alphaEM = sigmaProcessPtr->alphaEMRen();

Subprocess mass-square. double sHat = phaseSpacePtr->sHat();

Now methods only relevant for 2 -> 2. if (nFinal == 2) {

Mandelstam variables and hard-process pT. double tHat = phaseSpacePtr->tHat(); double uHat = phaseSpacePtr->uHat(); double pTHat = phaseSpacePtr->pTHat();

Masses of the final-state particles. (Here 0 for light quarks.) double m3 = sigmaProcessPtr->m(3); double m4 = sigmaProcessPtr->m(4);

}

Insert here your calculation of the selection bias. Here illustrated by a weighting up of events at high pT. selBias = pow4(phaseSpacePtr->pTHat());

Return the selBias weight. Warning: if you use another variable than selBias the compensating weight will not be set correctly. return selBias;

Dummy statement to avoid compiler warnings.

Reimplemented in UserHooksVector.

virtual bool canSetLowEnergySigma ( int  ,
int   
) const
inlinevirtual

Possibility to set low energy cross sections. Usage: canSetLowEnergySigma(idA, idB).

virtual bool canVetoAfterHadronization ( )
inlinevirtual

Possibility to veto an event after hadronization based on event contents. Works as an early trigger to avoid running the time consuming rescattering process on uninteresting events.

Reimplemented in UserHooksVector.

virtual bool canVetoMPIStep ( )
inlinevirtual

Possibility to veto MPI + ISR + FSR evolution and kill event, making decision after fixed number of MPI steps.

Reimplemented in UserHooksVector, and PowhegHooks.

virtual bool canVetoPartonLevelEarly ( )
inlinevirtual

Possibility to veto event after ISR + FSR + MPI in parton level, but before beam remnants and resonance decays.

Reimplemented in UserHooksVector, JetMatching, JetMatchingMadgraphInputAlpgen, and JetMatchingAlpgenInputAlpgen.

virtual bool canVetoPT ( )
inlinevirtual

Possibility to veto MPI + ISR + FSR evolution and kill event, making decision at a fixed pT scale. Useful for MLM-style matching.

Reimplemented in UserHooksVector.

virtual bool canVetoStep ( )
inlinevirtual

Possibility to veto MPI + ISR + FSR evolution and kill event, making decision after fixed number of ISR or FSR steps.

Reimplemented in UserHooksVector, JetMatchingMadgraph, and JetMatching.

virtual bool doChangeFragPar ( StringFlav ,
StringZ ,
StringPT ,
int  ,
double  ,
vector< int >  ,
const StringEnd  
)
inlinevirtual

Do change fragmentation parameters. Input: flavPtr, zPtr, pTPtr, idEnd, m2Had, iParton and posEnd (or negEnd).

Reimplemented in UserHooksVector.

virtual bool doReconnectResonanceSystems ( int  ,
Event  
)
inlinevirtual

Do reconnect colours from resonance decay systems. Usage: doVetoFSREmission( oldSizeEvt, event) where oldSizeEvent is the event size before resonance decays. Should normally return true, while false means serious failure. Value of PartonLevel:earlyResDec determines where method is called.

Reimplemented in TopReconUserHooks, UserHooksVector, and MBReconUserHooks.

virtual double doSetLowEnergySigma ( int  ,
int  ,
double  ,
double  ,
double   
) const
inlinevirtual

Set low energy cross sections for ids where canSetLowEnergySigma is true. Usage: doSetLowEnergySigma(idA, idB, eCM, mA, mB).

virtual bool doVetoFragmentation ( Particle  ,
const StringEnd  
)
inlinevirtual

Do a veto on a hadron just before it is added to the final state. The StringEnd from which the the hadron was produced is included for information.

Reimplemented in UserHooksVector.

virtual bool doVetoFragmentation ( Particle  ,
Particle  ,
const StringEnd ,
const StringEnd  
)
inlinevirtual

Do a veto on a hadron just before it is added to the final state (final two hadron case).

Reimplemented in UserHooksVector.

virtual bool doVetoFSREmission ( int  ,
const Event ,
int  ,
bool  = false 
)
inlinevirtual

Decide whether to veto current emission or not, based on event record. Usage: doVetoFSREmission( sizeOld, event, iSys, inResonance) where sizeOld is size of event record before current emission-to-be-scrutinized was added, iSys is the system of the radiation (according to PartonSystems), and inResonance is true if the emission takes place in a resonance decay.

Reimplemented in VinciaEWVetoHook, PowhegHooks, and UserHooksVector.

virtual bool doVetoISREmission ( int  ,
const Event ,
int   
)
inlinevirtual

Decide whether to veto current emission or not, based on event record. Usage: doVetoISREmission( sizeOld, event, iSys) where sizeOld is size of event record before current emission-to-be-scrutinized was added, and iSys is the system of the radiation (according to PartonSystems).

Reimplemented in VinciaEWVetoHook, UserHooksVector, and PowhegHooks.

virtual bool doVetoMPIEmission ( int  ,
const Event  
)
inlinevirtual

Decide whether to veto an MPI based on event record. Usage: doVetoMPIEmission( sizeOld, event) where sizeOld is size of event record before the current MPI.

Reimplemented in PowhegHooks, and UserHooksVector.

virtual bool doVetoMPIStep ( int  ,
const Event  
)
inlinevirtual

Decide whether to veto current event or not, based on event record. Usage: doVetoMPIStep( nMPI, event), where nMPI is number of MPI's so far.

Reimplemented in UserHooksVector, and PowhegHooks.

virtual bool doVetoPartonLevel ( const Event )
inlinevirtual

Decide whether to veto current partons or not, based on event record. Usage: doVetoPartonLevel( event).

Reimplemented in UserHooksVector.

virtual bool doVetoPartonLevelEarly ( const Event )
inlinevirtual

Decide whether to veto current partons or not, based on event record. Usage: doVetoPartonLevelEarly( event).

Reimplemented in UserHooksVector, JetMatching, JetMatchingMadgraphInputAlpgen, and JetMatchingAlpgenInputAlpgen.

virtual bool doVetoProcessLevel ( Event )
inlinevirtual

Decide whether to veto current process or not, based on process record. Usage: doVetoProcessLevel( process).

Reimplemented in UserHooksVector, JetMatchingMadgraph, JetMatching, JetMatchingMadgraphInputAlpgen, JetMatchingAlpgenInputAlpgen, amcnlo_unitarised_interface, and SetLHEDecayProductHook.

virtual bool doVetoPT ( int  ,
const Event  
)
inlinevirtual

Decide whether to veto current event or not, based on event record. Usage: doVetoPT( iPos, event), where iPos = 0: no emissions so far; iPos = 1/2/3 joint evolution, latest step was MPI/ISR/FSR; iPos = 4: FSR only afterwards; iPos = 5: FSR in resonance decay.

Reimplemented in UserHooksVector.

virtual bool doVetoResonanceDecays ( Event )
inlinevirtual

Decide whether to veto current resonance decay chain or not, based on process record. Usage: doVetoProcessLevel( process).

Reimplemented in UserHooksVector, and ResonanceDecayFilterHook.

virtual bool doVetoStep ( int  ,
int  ,
int  ,
const Event  
)
inlinevirtual

Decide whether to veto current event or not, based on event record. Usage: doVetoStep( iPos, nISR, nFSR, event), where iPos as above, nISR and nFSR number of emissions so far for hard interaction only.

Reimplemented in UserHooksVector, JetMatchingMadgraph, and JetMatching.

double multiplySigmaBy ( const SigmaProcess sigmaProcessPtr,
const PhaseSpace phaseSpacePtr,
bool  inEvent 
)
virtual

Multiplicative factor modifying the cross section of a hard process.

The UserHooks class.

multiplySigmaBy allows the user to introduce a multiplicative factor that modifies the cross section of a hard process. Since it is called from before the event record is generated in full, the normal analysis does not work. The code here provides a rather extensive summary of which methods actually do work. It is a convenient starting point for writing your own derived routine.

Process code, necessary when some to be treated differently. int code = sigmaProcessPtr->code();

Final multiplicity, i.e. whether 2 -> 1 or 2 -> 2. int nFinal = sigmaProcessPtr->nFinal();

Incoming x1 and x2 to the hard collision, and factorization scale. double x1 = phaseSpacePtr->x1(); double x2 = phaseSpacePtr->x2(); double Q2Fac = sigmaProcessPtr->Q2Fac();

Renormalization scale and assumed alpha_strong and alpha_EM. double Q2Ren = sigmaProcessPtr->Q2Ren(); double alphaS = sigmaProcessPtr->alphaSRen(); double alphaEM = sigmaProcessPtr->alphaEMRen();

Subprocess mass-square. double sHat = phaseSpacePtr->sHat();

Now methods only relevant for 2 -> 2. if (nFinal == 2) {

Mandelstam variables and hard-process pT. double tHat = phaseSpacePtr->tHat(); double uHat = phaseSpacePtr->uHat(); double pTHat = phaseSpacePtr->pTHat();

Masses of the final-state particles. (Here 0 for light quarks.) double m3 = sigmaProcessPtr->m(3); double m4 = sigmaProcessPtr->m(4);

}

Dummy statement to avoid compiler warnings.

Reimplemented in UserHooksVector, and SuppressSmallPT.

void omitResonanceDecays ( const Event process,
bool  finalOnly = false 
)
protected

omitResonanceDecays omits resonance decay chains from process record.

Reset work event to be empty

Loop through all partons. Beam particles should be copied.

Daughters of beams should normally be copied.

Granddaughters of beams should normally be copied and are final.

Optionally non-final are not copied.

Do copying and modify status/daughters of final.

When final only : no mothers; position in full event as daughters.

virtual void onInitInfoPtr ( )
inlineoverrideprotectedvirtual

After initInfoPtr, initialize workEvent.

Set smart pointer to null, in order to avoid circular dependency

Reimplemented from PhysicsBase.

virtual bool retryPartonLevel ( )
inlinevirtual

Retry same ProcessLevel with a new PartonLevel after a veto in doVetoPT, doVetoStep, doVetoMPIStep or doVetoPartonLevelEarly if you overload this method to return true.

Reimplemented in UserHooksVector.

virtual double scaleResonance ( int  ,
const Event  
)
inlinevirtual

Initial scale for TimeShower evolution. Usage: scaleResonance( iRes, event), where iRes is location of decaying resonance in the event record.

Reimplemented in MergeResScaleHook, and UserHooksVector.

virtual void setStringEnds ( const StringEnd ,
const StringEnd ,
vector< int >   
)
inlinevirtual

Set initial ends of a string to be fragmented. This is done once for each string. Note that the second string end may be zero in case we are hadronising a string piece leading to a junction.

Reimplemented in UserHooksVector.

void subEvent ( const Event event,
bool  isHardest = true 
)
protected

subEvent extracts currently resolved partons in the hard process.

Reset work event to be empty.

At the PartonLevel final partons are bookkept by subsystem.

Find which subsystem to study.

Loop through all the final partons of the given subsystem.

Copy partons to work event.

No mothers. Position in full event as daughters.

At the ProcessLevel no subsystems have been defined.

Loop through all partons, and copy all final ones.

No mothers. Position in full event as daughters.


The documentation for this class was generated from the following files: