PYTHIA  8.311
Classes | Public Member Functions | Public Attributes | Friends | List of all members
MergingHooks Class Reference

MergingHooks is base class for user input to the merging procedure. More...

#include <MergingHooks.h>

Inheritance diagram for MergingHooks:
PhysicsBase DireMergingHooks VinciaMergingHooks

Classes

struct  IndividualWeights
 Struct to save individual weights. More...
 

Public Member Functions

 MergingHooks ()
 Constructor.
 
virtual ~MergingHooks ()
 Destructor. More...
 
virtual double tmsDefinition (const Event &event)
 Function encoding the functional definition of the merging scale.
 
virtual double dampenIfFailCuts (const Event &inEvent)
 
virtual bool canCutOnRecState ()
 
virtual bool doCutOnRecState (const Event &event)
 
virtual bool canVetoTrialEmission ()
 Function to allow not counting a trial emission.
 
virtual bool doVetoTrialEmission (const Event &, const Event &)
 Function to check if trial emission should be rejected.
 
virtual double hardProcessME (const Event &inEvent)
 Function to calculate the hard process matrix element. More...
 
virtual void init ()
 Initialise MergingHooks class. More...
 
double tms ()
 Function returning the value of the merging scale. More...
 
double tmsCut ()
 
void tms (double tmsIn)
 
double dRijMS ()
 
double pTiMS ()
 
double QijMS ()
 
int nMaxJets ()
 Function returning the value of the maximal number of merged jets.
 
int nMaxJetsNLO ()
 
string getProcessString ()
 Function to return hard process string.
 
int nHardOutPartons ()
 Function to return the number of outgoing partons in the core process.
 
int nHardOutLeptons ()
 Function to return the number of outgoing leptons in the core process.
 
int nHardOutBosons ()
 
int nHardInPartons ()
 
int nHardInLeptons ()
 Function to return the number of incoming leptons in the core process.
 
int nResInCurrent ()
 
bool doUserMerging ()
 Function to determine if user defined merging should be applied.
 
bool doMGMerging ()
 Function to determine if automated MG/ME merging should be applied.
 
bool doKTMerging ()
 Function to determine if KT merging should be applied.
 
bool doPTLundMerging ()
 Function to determine if PTLund merging should be applied.
 
bool doCutBasedMerging ()
 Function to determine if cut based merging should be applied.
 
bool doCKKWLMerging ()
 
bool doUMEPSTree ()
 
bool doUMEPSSubt ()
 
bool doUMEPSMerging ()
 
bool doNL3Tree ()
 
bool doNL3Loop ()
 
bool doNL3Subt ()
 
bool doNL3Merging ()
 
bool doUNLOPSTree ()
 
bool doUNLOPSLoop ()
 
bool doUNLOPSSubt ()
 
bool doUNLOPSSubtNLO ()
 
bool doUNLOPSMerging ()
 
bool doRuntimeAMCATNLOInterface ()
 Function to determine if we have a runtime interface to aMC.
 
int nRecluster ()
 Return the number clustering steps that have actually been done.
 
int nRequested ()
 Return number of requested additional jets on top of the Born process.
 
bool isFirstEmission (const Event &event)
 Output functions to analyse/prepare event for merging. More...
 
bool hasEffectiveG2EW ()
 Function to allow effective gg -> EW boson couplings.
 
bool allowEffectiveVertex (vector< int > in, vector< int > out)
 Function to allow effective gg -> EW boson couplings.
 
Event bareEvent (const Event &inputEventIn, bool storeInputEvent)
 Return event stripped from decay products. More...
 
bool reattachResonanceDecays (Event &process)
 Write event with decay products attached to argument. More...
 
bool isInHard (int iPos, const Event &event)
 Check if particle at position iPos is part of the hard sub-system. More...
 
virtual int getNumberOfClusteringSteps (const Event &event, bool resetNjetMax=false)
 Function to return the number of clustering steps for the current event. More...
 
void orderHistories (bool doOrderHistoriesIn)
 Functions to steer construction of histories. More...
 
void allowCutOnRecState (bool doCutOnRecStateIn)
 
void doWeakClustering (bool doWeakClusteringIn)
 Function to allow final state clusterings of weak bosons.
 
bool checkAgainstCut (const Particle &particle)
 Functions used as default merging scales. More...
 
virtual double tmsNow (const Event &event)
 
double rhoms (const Event &event, bool withColour)
 Find the minimal Lund pT between coloured partons in the event. More...
 
double kTms (const Event &event)
 Function to calculate the minimal kT in the event. More...
 
double cutbasedms (const Event &event)
 
void doIgnoreEmissions (bool doIgnoreIn)
 Functions to steer shower evolution (public to allow for PS plugin) More...
 
virtual bool canVetoEmission ()
 Function to allow not counting a trial emission.
 
virtual bool doVetoEmission (const Event &)
 Function to check if emission should be rejected. More...
 
virtual bool usesVincia ()
 
virtual bool useShowerPlugin ()
 
bool includeWGTinXSEC ()
 
int nHardNow ()
 Functions to retrieve event veto information.
 
double tmsHardNow ()
 
int nJetsNow ()
 
double tmsNow ()
 
void setHardProcessPtr (HardProcess *hardProcIn)
 
int nMuRVar ()
 Functions related to renormalization scale variations.
 
void printIndividualWeights ()
 Function to print individual merging weight components, for debugging.
 
void setShowerPointer (PartonLevel *psIn)
 
void storeHardProcessCandidates (const Event &event)
 Generic setup functions. More...
 
void setLHEInputFile (string lheFile)
 
AlphaStrongAlphaS_FSR ()
 Functions for output of class members. More...
 
AlphaStrongAlphaS_ISR ()
 
AlphaEMAlphaEM_FSR ()
 
AlphaEMAlphaEM_ISR ()
 
bool includeMassive ()
 
bool enforceStrongOrdering ()
 Prefer strongly ordered histories.
 
bool orderInRapidity ()
 Prefer histories ordered in rapidity and evolution pT.
 
bool pickByFull ()
 Pick history probabilistically by full (correct) splitting probabilities.
 
bool pickByPoPT2 ()
 Pick history probabilistically, with easier form of probabilities.
 
bool includeRedundant ()
 Include redundant terms (e.g. PDF ratios) in the splitting probabilities.
 
bool pickBySumPT ()
 Pick by winner-takes-it-all, with minimum sum of scalar evolution pT.
 
int unorderedScalePrescip ()
 
int unorderedASscalePrescip ()
 
int unorderedPDFscalePrescip ()
 
int incompleteScalePrescip ()
 
bool allowColourShuffling ()
 Allow swapping one colour index while reclustering.
 
bool resetHardQRen ()
 Allow use of dynamical renormalisation scale of the core 2-> 2 process.
 
bool resetHardQFac ()
 Allow use of dynamical factorisation scale of the core 2-> 2 process.
 
double scaleSeparationFactor ()
 
double nonJoinedNorm ()
 
double fsrInRecNorm ()
 
double herwigAcollFSR ()
 
double herwigAcollISR ()
 
double pT0ISR ()
 ISR regularisation scale.
 
double pTcut ()
 Shower cut-off scale.
 
void muMI (double mu)
 MI starting scale.
 
double muMI ()
 
double kFactor (int njet=0)
 Full k-Factor for current event.
 
double k1Factor (int njet=0)
 O()-term of the k-Factor for current event.
 
bool orderHistories ()
 
bool allowCutOnRecState ()
 
bool doWeakClustering ()
 INTERNAL Hooks to allow clustering W bosons.
 
bool doSQCDClustering ()
 INTERNAL Hooks to allow clustering clustering of gluons to squarks.
 
double muF ()
 Store / get first scale in PDF's that Pythia should have used.
 
double muR ()
 
double muFinME ()
 Store / get factorisation scale used in matrix element calculation. More...
 
double muRinME ()
 
void doIgnoreStep (bool doIgnoreIn)
 Functions to steer merging code. More...
 
virtual bool canVetoStep ()
 Function to allow event veto.
 
void storeWeights (vector< double > weight)
 Stored weights in case veto needs to be revoked.
 
virtual bool doVetoStep (const Event &process, const Event &event, bool doResonance=false)
 Function to check event veto. More...
 
virtual bool setShowerStartingScales (bool isTrial, bool doMergeFirstEmm, double &pTscaleIn, const Event &event, double &pTmaxFSRIn, bool &limitPTmaxFSRin, double &pTmaxISRIn, bool &limitPTmaxISRin, double &pTmaxMPIIn, bool &limitPTmaxMPIin)
 Set starting scales. More...
 
void setShowerStoppingScale (double scale=0.)
 
double getShowerStoppingScale ()
 
void nMinMPI (int nMinMPIIn)
 
int nMinMPI ()
 
double kTdurham (const Particle &RadAfterBranch, const Particle &EmtAfterBranch, int Type, double D)
 Functions for internal merging scale definions. More...
 
double rhoPythia (const Event &event, int rad, int emt, int rec, int ShowerType)
 Function to compute "pythia pT separation" from Particle input. More...
 
int findColour (int col, int iExclude1, int iExclude2, const Event &event, int type, bool isHardIn)
 
double deltaRij (Vec4 jet1, Vec4 jet2)
 Function to compute Delta R separation from 4-vector input. More...
 
double getWeightNLO (int i=0)
 Functions for weight management. More...
 
vector< double > getWeightCKKWL ()
 Return CKKW-L weight.
 
vector< double > getWeightFIRST ()
 Return O() weight.
 
void setWeightCKKWL (vector< double > weightIn)
 Set CKKW-L weight.
 
void setWeightFIRST (vector< double > weightIn)
 Set O() weight.
 
vector< double > getSudakovWeight ()
 
vector< double > getCouplingWeight ()
 Function to return coupling weight.
 
void setEventVetoInfo (int nJetNowIn, double tmsNowIn)
 Functions and members to store the event veto information. More...
 
void setHardProcessInfo (int nHardNowIn, double tmsHardNowIn)
 Set the hard process information.
 
void addVetoInMainShower ()
 
int getNumberVetoedInMainShower ()
 
- 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
 

Public Attributes

bool useShowerPluginSave
 Functions used as clusterings / probabilities.
 
bool useOwnHardProcess
 The members, switches etc. More...
 
HardProcesshardProcess
 
PartonLevelshowers
 
AlphaStrong AlphaS_FSRSave
 AlphaS objects for alphaS reweighting use.
 
AlphaStrong AlphaS_ISRSave
 
AlphaEM AlphaEM_FSRSave
 
AlphaEM AlphaEM_ISRSave
 
string lheInputFile
 Saved path to LHE file for more automated merging.
 
bool doUserMergingSave
 Flags for merging procedure definition.
 
bool doMGMergingSave
 
bool doKTMergingSave
 
bool doPTLundMergingSave
 
bool doCutBasedMergingSave
 
bool includeMassiveSave
 
bool enforceStrongOrderingSave
 
bool orderInRapiditySave
 
bool pickByFullPSave
 
bool pickByPoPT2Save
 
bool includeRedundantSave
 
bool pickBySumPTSave
 
bool allowColourShufflingSave
 
bool resetHardQRenSave
 
bool resetHardQFacSave
 
int unorderedScalePrescipSave
 
int unorderedASscalePrescipSave
 
int unorderedPDFscalePrescipSave
 
int incompleteScalePrescipSave
 
int ktTypeSave
 
int nReclusterSave
 
int nQuarksMergeSave
 
int nRequestedSave
 
double scaleSeparationFactorSave
 
double nonJoinedNormSave
 
double fsrInRecNormSave
 
double herwigAcollFSRSave
 
double herwigAcollISRSave
 
double pT0ISRSave
 
double pTcutSave
 
double pTminISRSave
 
double pTminFSRSave
 
bool doNL3TreeSave
 
bool doNL3LoopSave
 
bool doNL3SubtSave
 
bool doUNLOPSTreeSave
 
bool doUNLOPSLoopSave
 
bool doUNLOPSSubtSave
 
bool doUNLOPSSubtNLOSave
 
bool doUMEPSTreeSave
 
bool doUMEPSSubtSave
 
bool doEstimateXSection
 Flag to only do phase space cut, rejecting events below the tms cut.
 
bool doRuntimeAMCATNLOInterfaceSave
 Flag for runtime aMC interface. Needed for aMC-Delta.
 
bool applyVeto
 
Event inputEvent
 Save input event in case decay products need to be detached.
 
vector< pair< int, int > > resonances
 
bool doRemoveDecayProducts
 
double muMISave
 Starting scale for attaching MPI.
 
double kFactor0jSave
 Precalculated K-factors.
 
double kFactor1jSave
 
double kFactor2jSave
 
double tmsValueSave
 Saved members.
 
double tmsValueNow
 
double DparameterSave
 
int nJetMaxSave
 
int nJetMaxNLOSave
 
string processSave
 
string processNow
 
vector< double > tmsListSave
 
bool doOrderHistoriesSave
 
bool doCutOnRecStateSave
 
bool doWeakClusteringSave
 INTERNAL Hooks to allow clustering W bosons.
 
bool doSQCDClusteringSave
 
double muFSave
 Store / get first scale in PDF's that Pythia should have used.
 
double muRSave
 
double muFinMESave
 Store / get factorisation scale used in matrix element calculation.
 
double muRinMESave
 
bool doIgnoreEmissionsSave
 Flag to indicate trial shower usage.
 
bool doIgnoreStepSave
 Flag to indicate if events should be vetoed.
 
double pTsave
 Stored weights in case veot needs to be revoked.
 
vector< double > weightCKKWL1Save
 
vector< double > weightCKKWL2Save
 
int nMinMPISave
 
vector< double > weightCKKWLSave
 Save CKKW-L weight / O() weight.
 
vector< double > weightFIRSTSave
 
IndividualWeights individualWeights
 
bool doVariations
 Flag to indicate whether renormalization scale variations are performed.
 
vector< double > muRVarFactors
 Vector of variation factors applied to renormalization scale.
 
int nWgts
 Number of weights, nominal + variations.
 
int nJetMaxLocal
 Local copies of nJetMax inputs, if recalculation is necessary.
 
int nJetMaxNLOLocal
 
bool hasJetMaxLocal
 
bool includeWGTinXSECSave
 
int nHardNowSave
 
int nJetNowSave
 
double tmsHardNowSave
 
double tmsNowSave
 
double stopScaleSave
 
int nVetoedInMainShower
 Statistics.
 

Friends

class History
 Make History class friend to allow access to advanced switches.
 
class Pythia
 Make Pythia class friend.
 
class PartonLevel
 Make PartonLevel class friend.
 
class SpaceShower
 Make SpaceShower class friend.
 
class TimeShower
 Make TimeShower class friend.
 
class Merging
 Make Merging class friend.
 

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.
 
- Protected Member Functions inherited from PhysicsBase
 PhysicsBase ()
 Default constructor.
 
virtual void onInitInfoPtr ()
 
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 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
 

Detailed Description

MergingHooks is base class for user input to the merging procedure.

Constructor & Destructor Documentation

~MergingHooks ( )
virtual

Destructor.

The MergingHooks class.

Destructor

Member Function Documentation

void allowCutOnRecState ( bool  doCutOnRecStateIn)
inline

Function to force cut on reconstructed states internally, as needed for gg -> Higgs to ensure that e.g. uu~ -> Higgs is not constructed.

bool allowCutOnRecState ( )
inline

INTERNAL Hooks to disallow states in the construction of all histories, e.g. because jets are below the merging scale, of to avoid the construction of uu~ -> Higgs histories.

AlphaStrong* AlphaS_FSR ( )
inline

Functions for output of class members.

Return AlphaStrong objects

Event bareEvent ( const Event inputEventIn,
bool  storeInputEvent 
)

Return event stripped from decay products.

Find and detach decay products.

If desired, store input event.

Now remove decay products.

Add the beam and initial partons to the event record.

Add the intermediate particles to the event record.

Add remaining outgoing particles to the event record.

Update event colour tag to maximum in whole process.

Copy junctions from process to newProcess.

Remember scale

Done

virtual bool canCutOnRecState ( )
inlinevirtual

Hooks to disallow states in the construction of all histories, e.g. because jets are below the merging scale or fail the matrix element cuts Function to allow interference in the construction of histories

Reimplemented in VinciaMergingHooks.

bool checkAgainstCut ( const Particle particle)

Functions used as default merging scales.

Function to check if the input particle is a light jet, i.e. should be checked against the merging scale defintion.

Function to check if the properties of the input particle should be checked against the cut-based merging scale defintion.

Do not check uncoloured particles.

By default, use u-, d-, c-, s- and b-quarks.

Done

double cutbasedms ( const Event event)

Find the if the event passes the Delta R_{ij}, pT_{i} and Q_{ij} cuts on the matrix element, as needed for cut-based merging scale definition

Find the if the event passes the Delta R_{ij}, pT_{i} and Q_{ij} cuts on the matrix element, as needed for cut-based merging scale definition.

Only check first emission.

Save allowed final state particles

Declare overall veto

Declare vetoes

Declare cuts used in matrix element

Declare minimum values

Check matrix element cuts

Save pT value

Check two-parton cuts

Save delta R value

Save delta R value

Done with cut evaluation

Check if all partons are in the matrix element region

In the matrix element calculation, we have rejected the events if any of the cuts had not been fulfilled, i.e. we are in the matrix element domain if all cuts are fulfilled. We veto any emission in the ME region. Disregard the two-parton cuts if only one parton in the final state

Veto if the combination of cuts would be in the ME region

If event is above merging scale, veto

Else, do nothing

virtual double dampenIfFailCuts ( const Event inEvent)
inlinevirtual

Function to dampen weights calculated from histories with lowest multiplicity reclustered events that do not pass the ME cuts

Dummy statement to avoid compiler warnings

Reimplemented in VinciaMergingHooks.

double deltaRij ( Vec4  jet1,
Vec4  jet2 
)

Function to compute Delta R separation from 4-vector input.

Declare return variable

Get delta_eta and cosh(Delta_eta) for hadronic collisions

Get delta_phi and cos(Delta_phi) for hadronic collisions

Calculate kT durham like fastjet

Return kT

virtual bool doCutOnRecState ( const Event event)
inlinevirtual

Function to check reclustered state while generating all possible histories Function implementing check of reclustered events while constructing all possible histories

Dummy statement to avoid compiler warnings.

Count number of final state partons.

For gg -> h, allow only histories with gluons in initial state

Reimplemented in VinciaMergingHooks.

void doIgnoreEmissions ( bool  doIgnoreIn)
inline

Functions to steer shower evolution (public to allow for PS plugin)

Flag to indicate trial shower usage.

void doIgnoreStep ( bool  doIgnoreIn)
inline

Functions to steer merging code.

Flag to indicate if events should be vetoed.

bool doNL3Tree ( )
inline

Functions to determine if and which part of NL3 merging should be applied

bool doUMEPSTree ( )
inline

Functions to determine if and which part of UMEPS merging should be applied

bool doUNLOPSTree ( )
inline

Functions to determine if and which part of UNLOPS merging should be applied

bool doVetoEmission ( const Event event)
virtual

Function to check if emission should be rejected.

Do nothing in trial showers, or after first step.

Do nothing in CKKW-L

For NLO merging, count and veto emissions above the merging scale

Get number of clustering steps

Get merging scale in current event

Get maximal number of additional jets

Always remove emissions above the merging scale for samples containing reclusterings!

Check veto condition

Do not veto if state already includes MPI.

When performing NL3 merging of tree-level events, reset the CKKWL weight.

If the emission is allowed, do not check any further emissions

Done

Reimplemented in VinciaMergingHooks, and DireMergingHooks.

bool doVetoStep ( const Event process,
const Event event,
bool  doResonance = false 
)
virtual

Function to check event veto.

Function to check if emission should be rejected.

Do nothing in trial showers, or after first step.

Do nothing in UMEPS or UNLOPS

Get number of clustering steps. If necessary, remove resonance decay products first.

Get maximal number of additional jets.

Get merging scale in current event.

Do not print zero-weight events. For non-resonant showers, simply check veto. If event should indeed be vetoed, save the current pT and weights in case the veto needs to be revoked.

Store pT to check if veto needs to be revoked later.

Store veto inputs to perform veto at a later stage.

Set weight to zero if event should be vetoed.

Save weight before veto, in case veto needs to be revoked.

Reset stored weights.

Check merging veto condition.

Set weight to zero if event should be vetoed.

Save weight before veto, in case veto needs to be revoked.

Reset stored weights.

Done

For resonant showers, check if any previous veto should be revoked. This means we treat showers off resonance decay products identical to MPI: If a hard resonance emission has been produced, the event should have been kept. Following this reasoning, it might be necessary to revoke any previous veto.

Initialise switch to revoke vetoing.

Nothing to check if pTsave was not stored, i.e. no emission to possibly veto was recorded.

For hadronic resonance decays at hadron colliders, do not veto events with a hard emission of the resonance decay products, since such states are not included in the matrix element

Check how many resonance decay systems are allowed

Find systems really containing only emission off a resonance decay

Resonance decay systems are considered last, thus at the end of the list

Check the members of the resonance decay systems

Save the (three) members of the resonance decay system

Now find emitted gluon or gamma

Per system, only one emission

Revoke previous veto of last emission if a splitting of the resonance produced a harder parton, i.e. we are inside the PS region

Do nothing (i.e. allow other first emission veto) for soft splitting

Done with one system

Done with all systems

Check veto condition

Check veto condition.

Set stored weights to zero.

Now allow veto.

If the emission is allowed, do not check any further emissions

Done

Done

Reimplemented in VinciaMergingHooks, and DireMergingHooks.

double dRijMS ( )
inline

Function returning the value of the Delta R_{ij} cut for cut based merging scale definition.

int findColour ( int  col,
int  iExclude1,
int  iExclude2,
const Event event,
int  type,
bool  isHardIn 
)

Function to find a colour (anticolour) index in the input event, used to find colour-connected recoilers

Function to find a colour (anticolour) index in the input event. Helper for rhoms IN int col : Colour tag to be investigated int iExclude1 : Identifier of first particle to be excluded from search int iExclude2 : Identifier of second particle to be excluded from search Event event : event to be searched for colour tag int type : Tag to define if col should be counted as colour (type = 1) [->find anti-colour index contracted with col] anticolour (type = 2) [->find colour index contracted with col] OUT int : Position of particle in event record contraced with col [0 if col is free tag]

Search event record for matching colour & anticolour

Check outgoing

Check incoming

Search event record for matching colour & anticolour

Check outgoing from ISR

Check outgoing from FSR

Check outgoing from FSR

first initial

second initial

if no matching colour / anticolour has been found, return false

double fsrInRecNorm ( )
inline

Absolute normalization of splitting probability for final state splittings with initial state recoiler

int getNumberOfClusteringSteps ( const Event event,
bool  resetNjetMax = false 
)
virtual

Function to return the number of clustering steps for the current event.

Count the number of final state partons

Count the number of final state leptons

Add neutralinos to number of leptons

Add sleptons to number of leptons

Count the number of final state electroweak bosons

Save sum of all final state particles

Return the difference to the core process outgoing particles

For inclusive handling, the number of reclustering steps can be different within a single sample.

Final particle counters

Set steps for QCD or QCD+QED events: Need at least two massless particles at lowest multiplicity.

Set steps for events containing heavy bosons. Need at least one massive particle at lowest multiplicity.

dynamical handling of steps

Return the difference to the core process outgoing particles

Reimplemented in VinciaMergingHooks, and DireMergingHooks.

vector<double> getSudakovWeight ( )
inline

Function to return Sudakov weight as calculated before, also include MPI weight. Only call after regular weight functions, since it is calculated there.

double getWeightNLO ( int  i = 0)
inline

Functions for weight management.

Function to get the CKKW-L weight for the current event

virtual double hardProcessME ( const Event inEvent)
inlinevirtual

Function to calculate the hard process matrix element.

Dummy statement to avoid compiler warnings.

Reimplemented in VinciaMergingHooks.

double herwigAcollFSR ( )
inline

Factor multiplying scalar evolution pT for FSR splitting, when picking history by minimum scalar pT (see Jonathan Tully's thesis)

double herwigAcollISR ( )
inline

Factor multiplying scalar evolution pT for ISR splitting, when picking history by minimum scalar pT (see Jonathan Tully's thesis)

bool includeMassive ( )
inline

Functions to return advanced merging switches Include masses in definition of evolution pT and splitting kernels

bool includeWGTinXSEC ( )
inline

Functions to retrieve if merging weight should countin the internal cross section and the event weight.

int incompleteScalePrescip ( )
inline

Prescription for starting scale of incomplete histories 0: use factorization scale 1: use sHat 2: use s

void init ( )
virtual

Initialise MergingHooks class.

Functions for internal use inside Pythia source code Initialize.

Save pointers

Initialise AlphaS objects for reweighting

Initialise AlphaEM objects for reweighting

Initialise merging switches

Initialise automated MadGraph kT merging

Initialise kT merging

Initialise evolution-pT merging

Initialise {ij}, pT_i Q_{ij} merging

Initialise exact definition of kT

Initialise NL3 switches.

Initialise UNLOPS switches.

Initialise UMEPS switches

Flag to only do phase space cut.

Flag to check if we have an aMC runtime interface.

Flag to check if merging weight should directly be included in the cross section.

Flag to check if CKKW-L event veto should be applied.

Get core process from user input

If the process string is "guess", temporarily set it to something safe for initialization.

Clear hard process

Initialise input event.

Initialise the hard process

Remove whitespace from process string

Parameters for reconstruction of evolution scales

Parameters for choosing history probabilistically

Parameters for scale choices

Parameter for allowing swapping of one colour index while reclustering

Parameters to allow setting hard process scales to default (dynamical) Pythia values.

Parameters for choosing history by sum(|pT|)

Information on the shower cut-off scale

Information on renormalization scale variations

Initialise CKKWL weight

Initialize merging weights in weight container

Initialise merging scale

Save merging scale on maximal number of jets

Read merging scale (defined in kT) from input parameter.

Read merging scale (defined in kT) from LHE file.

Save list of cuts defining the merging scale.

Write tms cut values to list of cut values, ordered by DeltaR_{ij}, pT_{i}, Q_{ij}.

Read additional settings for NLO merging methods.

Internal Pythia cross section should not include NLO merging weights.

Check if external shower plugin should be used.

Write banner.

Reimplemented in VinciaMergingHooks, and DireMergingHooks.

bool isFirstEmission ( const Event event)

Output functions to analyse/prepare event for merging.

Function to check if event contains an emission not present in the hard process.

If the beam remnant treatment or hadronisation has already started, do no veto.

Count particle types

Return highest value if the event does not contain any final state coloured particles.

Use MergingHooks functions to get information on the hard process.

The state is already in the PS region if the number of leptons had been increased by QED splittings.

If the mumber of photons if larger than in the hard process, QED radiation has pushed the state into the PS region.

Done

bool isInHard ( int  iPos,
const Event event 
)

Check if particle at position iPos is part of the hard sub-system.

MPI not part of hard process

Beam remnants and hadronisation not part of hard process

Still MPI: Check that the particle is not due to radiation off MPIs. First get all intermediate MPI partons in the state.

Disregard any parton iPos that has MPI ancestors.

Disallow other systems. Get sub-system of particles for iPos

First do simple check if the system is sensible. Might not be the case when constructing a process record by hand, yet still keeping the partonSystems. Remember to not check problematic systems.

Check all partons belonging to the same system as iPos. If any is produced in MPI or has MPI ancestors, the whole system is not the hard subprocess, i.e. iPos is not in the hard subprocess.

MPI not part of hard process

Disregard any parton iPos that has MPI ancestors.

Beam remnants and hadronisation not part of hard process

Check if any ancestor contains the hard incoming partons as daughters. Begin loop to trace upwards from the daughter.

If out of range then failed to find match.

If positive match then done.

If unique mother then keep on moving up the chain.

Done

double kTdurham ( const Particle RadAfterBranch,
const Particle EmtAfterBranch,
int  Type,
double  D 
)

Functions for internal merging scale definions.

Function to compute durham y separation from Particle input.

Function to calculate the kT separation between two particles

Declare return variable

Save 4-momenta of final state particles

Get angle between jets for e+e- collisions, make sure that -1 <= cos(theta) <= 1

Calculate kt durham separation between jets for e+e- collisions

Get delta_y for hadronic collisions: Get mT of first jet

Get mT of second jet

Get rapidity of first jet

Get rapidity of second jet

Get delta_phi for hadronic collisions

Calculate kT durham like fastjet, but with rapidity instead of pseudo-rapidity

Get mT of first jet

Get mT of second jet

Get pseudo-rapidity of first jet

Get pseudo-rapidity of second jet

Get delta_phi and cos(Delta_phi) for hadronic collisions

Calculate kT durham like fastjet

Get cosh(Delta_eta) for hadronic collisions

Get delta_phi and cos(Delta_phi) for hadronic collisions

Calculate kT durham separation "SHERPA-like"

Return kT

double kTms ( const Event event)

Function to calculate the minimal kT in the event.

Function to return the minimal kT in the event. If doKTMerging = true, this function will be used as a merging scale definition.

Only check first emission.

Find all electroweak decayed bosons in the state.

Declare final parton vectors

Search in Event record for final state partons. Exclude decay products of ew resonance.

Except for e+e- -> jets, do not check radiation in resonance decays.

Find minimal Durham kT in event, using own function: Check definition of separation

Find minimal kT

Compute separation to the beam axis for hadronic collisions

Compute separation to other final state jets

Keep the minimal Durham separation

Return minimal Durham kT

double muFinME ( )
inline

Store / get factorisation scale used in matrix element calculation.

Start with checking the event attribute called "muf".

Check the scales tag of the event.

double muRinME ( )
inline

Start with checking the event attribute called "mur2".

Check the scales tag of the event.

int nHardInPartons ( )
inline

Function to return the number of incoming partons (hadrons) in the core process.

int nHardOutBosons ( )
inline

Function to return the number of outgoing electroweak bosons in the core process.

int nMaxJetsNLO ( )
inline

Function returning the value of the maximal number of merged jets, for which NLO corrections are available.

double nonJoinedNorm ( )
inline

Absolute normalization of splitting probability for non-joined evolution.

int nResInCurrent ( )
inline

Function to report the number of resonace decays in the 2->2 sub-process of the current state.

void orderHistories ( bool  doOrderHistoriesIn)
inline

Functions to steer construction of histories.

Function to force preferred picking of ordered histories. By default, unordered histories will only be considered if no ordered histories were found.

bool orderHistories ( )
inline

Function to return if construction of histories is biased towards ordered histories.

double pTiMS ( )
inline

Function returning the value of the pT_{i} cut for cut based merging scale definition.

double QijMS ( )
inline

Function returning the value of the pT_{i} cut for cut based merging scale definition.

bool reattachResonanceDecays ( Event process)

Write event with decay products attached to argument.

Write event with decay products attached to argument. Only possible if an input event with decay producs had been stored before.

Now reattach the decay products.

Vector of resonances for which the decay products were already attached.

Reset daughters and status of intermediate particles.

Get momenta in case of reclustering.

Find current mother copy (after clustering).

Only attempt if decays of this resonance were not attached.

Resonance can have been moved by clustering, so prepare to update colour and momentum information for system.

Attach resonance decay products.

Update colour and momentum information.

Update vertex information.

Update mothers.

Loop through event and attach remaining decays

Done if no daughters exist.

Attach daughters.

Update colour and momentum information.

Update vertex information.

Update mothers.

Modify mother status and daughters.

End of loop over all entries.

End loop over process entries.

End loop over resonances.

Update event colour tag to maximum in whole process.

Done.

double rhoms ( const Event event,
bool  withColour 
)

Find the minimal Lund pT between coloured partons in the event.

Find the minimal Lund pT between coloured partons in the input event. If doPTLundMerging = true, this function will be used as a merging scale definition.

Only check first emission.

Find all electroweak decayed bosons in the state.

Declare final parton vectors

Search inEvent record for final state partons. Exclude decay products of ew resonance.

Except for e+e- -> jets, do not check radiation in resonance decays.

Include photons into the tms definition for "weak+QCD merging".

Include Z-bosons into the tms definition for "weak+QCD merging".

Include W-bosons into the tms definition for "weak+QCD merging".

Get index of first incoming

Get index of second incoming

If no incoming of the cascade are found, try incoming

Find current incoming partons

For pure QCD set the cut to the pT of the dijet system.

Find minimal pythia pt in event

Compute pythia ISR separation i-jet and first incoming

Compute pythia ISR separation i-jet and second incoming

Compute pythia FSR separation between two jets, with knowledge of colour connections

Find recoiler in event

Check in final state

Check in initial state

Compute pythia FSR separation between two jets, without any knowledge of colour connections

Allow any parton as recoiler

Only check splittings allowed by flavour, e.g. only q -> qg and g -> qqbar

Compute pythia FSR separation between two jets, with initial recoiler without any knowledge of colour connections

Allow both initial partons as recoiler

Check with first initial as recoiler

Check with second initial as recoiler

Reset minimal y separation

double rhoPythia ( const Event event,
int  rad,
int  emt,
int  rec,
int  ShowerType 
)

Function to compute "pythia pT separation" from Particle input.

Function to compute "pythia pT separation" from Particle input, as a helper for rhoms(...).

Use external shower for merging. Ask showers for evolution variable.

Note: If massive particles are involved, this definition slightly differs from History:pTLund(), as we need to ensure consistency with aMC (!). In the latter, no masses are available at the point where the merging scale value is calculated, and thus masses are set by hand there, and consequently here.

Save type: 1 = FSR pT definition, else ISR definition

Set masses (as used in MG5).

Calculate virtuality of splitting

Splitting not possible for negative virtuality.

Construct 2->3 variables for FSR

Try to reconstruct flavour of radiator before emission.

gluon radiation: idBef = idAft

final state gluon splitting: idBef = 21

final state quark -> gluon conversion

initial state gluon splitting: idBef = -idEmt

initial state gluon -> quark conversion

W-boson radiation

ee flavour sensitive cut

Store masses both after and prior to emission.

Final state splitting not possible for negative "dipole mass".

Rescaling of recoiler for FSR with initial state recoiler.

Final-initial splitting not possible for negative rescaling.

Final state splitting not possible for ill-defined 3-body-variables.

Construct momenta of dipole before/after splitting for ISR

Prepare for more complicated z definition for massive splittings.

Calculate z of splitting, different for FSR and ISR

Splitting not possible for ill-defined energy sharing.

Separation of splitting, different for FSR and ISR

pT^2 = separation*virtuality

Check for threshold in ISR, only relevant for c and b. Use pT2 = (1 - z) * (Qsq + m^2).

Kinematically impossible splittings should not be included in the pT definition!

Return pT

double scaleSeparationFactor ( )
inline

Factor by which two scales should differ to be classified strongly ordered.

void setEventVetoInfo ( int  nJetNowIn,
double  tmsNowIn 
)
inline

Functions and members to store the event veto information.

Set CKKWL veto information.

void setLHEInputFile ( string  lheFile)
inline

Function to set the path to the LHE file, so that more automated merging can be used. Remove "_1.lhe" suffix from LHE file name. This is done so that the HarsProcess class can access both the +0 and +1 LHE files to find both the merging scale and the core process string Store.

bool setShowerStartingScales ( bool  isTrial,
bool  doMergeFirstEmm,
double &  pTscaleIn,
const Event event,
double &  pTmaxFSRIn,
bool &  limitPTmaxFSRIn,
double &  pTmaxISRIn,
bool &  limitPTmaxISRIn,
double &  pTmaxMPIIn,
bool &  limitPTmaxMPIIn 
)
virtual

Set starting scales.

Function to set the correct starting scales of the shower. Note: 2 -> 2 QCD systems can be produced by MPI. Hence, there is an overlap between MPI and "hard" 2 -> 2 QCD systems which needs to be removed by no-MPI probabilities. This means that for any "hard" 2 -> 2 QCD system, multiparton interactions should start at the maximal scale of multiple interactions. The same argument holds for any "hard" process that overlaps with MPI.

Local copies of power/wimpy shower booleans and scales.

Merging of EW+QCD showers with matrix elements: Ensure that

  1. any event with more than one final state particle will be showered from the reconstructed transverse momentum of the last emission, even if the factorisation scale is low.
  2. the shower starting scale for events with no emission is given by the (user-defined) choice.

Check if the process only contains two outgoing partons. If so, then this process could also have been produced by MPI. Thus, the MPI starting scale would need to be set accordingly to correctly attach a "no-MPI-probability" to multi-jet events. ("Hard" MPI are included by not restricting MPI when showering the lowest-multiplicity sample.)

SET THE STARTING SCALES FOR TRIAL SHOWERS.

Reset shower and MPI scales.

Reset to minimal scale for wimpy showers. Keep scales for EW+QCD merging.

For EW+QCD merging, apply wimpy shower only to 2->1 processes.

For pure QCD set the PS starting scales to the pT of the dijet system.

If necessary, set the MPI starting scale to the collider energy.

Reset phase space limitation flags

If necessary, set the MPI starting scale to the collider energy.

SET THE STARTING SCALES FOR REGULAR SHOWERS.

Remember if this is a "regular" shower off a reclustered event.

Reset shower and MPI scales.

Reset to minimal scale for wimpy showers. Keep scales for EW+QCD merging.

For EW+QCD merging, apply wimpy shower only to 2->1 processes.

For pure QCD set the PS starting scales to the pT of the dijet system.

If necessary, set the MPI starting scale to the collider energy.

For reclustered events, no-MPI-probability between "pTmaxMPI" and "scale" already included in the event weight.

Reset power/wimpy shower switches iand scales if necessary.

Done

Reimplemented in VinciaMergingHooks, and DireMergingHooks.

void storeHardProcessCandidates ( const Event event)
inline

Generic setup functions.

Function storing candidates for the hard process in the current event Needed in order not to cluster members of the core process

double tms ( )
inline

Function returning the value of the merging scale.

else return tmsValueSave;

double tmsNow ( const Event event)
virtual

Function to return the value of the merging scale function in the current event.

Get merging scale in current event.

Use KT/Durham merging scale definition.

Use Lund PT merging scale definition.

Use DeltaR_{ij}, pT_i, Q_{ij} combination merging scale definition.

Use NLO merging (Lund PT) merging scale definition.

Use NLO merging (Lund PT) merging scale definition.

Use UMEPS (Lund PT) merging scale definition.

Use user-defined merging scale.

Return merging scale value. Done

Reimplemented in VinciaMergingHooks, and DireMergingHooks.

int unorderedASscalePrescip ( )
inline

Prescription for combined scale used in alpha_s for unordered emissions 0 : use combined emission scale in alpha_s weight for both (!) splittings 1 : use original reclustered scales of each emission in alpha_s weight

int unorderedPDFscalePrescip ( )
inline

Prescription for combined scale used in PDF ratios for unordered emissions 0 : use combined emission scale in PDFs for both (!) splittings 1 : use original reclustered scales of each emission in PDF ratiost

int unorderedScalePrescip ( )
inline

Prescription for combined scale of unordered emissions 0 : use larger scale 1 : use smaller scale

Member Data Documentation

bool applyVeto

Flag for postponing event vetos. If false, events are not vetoed based on the merging scale in CKKW-L merging, but have to be vetoed by the user in the main program.

bool doCutOnRecStateSave

INTERNAL Hooks to disallow states in the construction of all histories, e.g. because jets are below the merging scale, of to avoid the construction of uu~ -> Higgs histories.

bool doOrderHistoriesSave

INTERNAL Hooks to allow construction of all histories, including un-ordered ones

bool includeWGTinXSECSave

Event veto and hard process information, if veto should not applied be directly, but is up to the user.

double stopScaleSave

Set shower stopping scale. Necessary to e.g. avoid accumulation of incorrect (low-pT) shower weights through trial showering.

vector<double> tmsListSave

List of cut values to used to define a merging scale. Ordering: 0: DeltaR_{jet_i,jet_j,min} 1: p_{T,jet_i,min} 2: Q_{jet_i,jet_j,min}

bool useOwnHardProcess

The members, switches etc.

Helper class doing all the core process book-keeping


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