herwig is hosted by Hepforge, IPPP Durham
Herwig 7.3.0
Herwig::SplittingFunction Class Referenceabstract

This is an abstract class which defines the common interface for all \(1\to2\) splitting functions, for both initial-state and final-state radiation. More...

#include <SplittingFunction.h>

Inheritance diagram for Herwig::SplittingFunction:

Public Member Functions

 SplittingFunction ()
 The default constructor.
 
ShowerInteraction interactionType () const
 Methods to return the interaction type and order for the splitting function.
 
ColourStructure colourStructure () const
 Return the colour structure.
 
double colourFactor (const IdList &ids) const
 Return the colour factor.
 
virtual bool accept (const IdList &ids) const =0
 Purely virtual method which should determine whether this splitting function can be used for a given set of particles.
 
virtual bool checkColours (const IdList &ids) const
 Method to check the colours are correct.
 
virtual double P (const double z, const Energy2 t, const IdList &ids, const bool mass, const RhoDMatrix &rho) const =0
 Methods to return the splitting function.
 
virtual double overestimateP (const double z, const IdList &ids) const =0
 Purely virtual method which should return an overestimate of the splitting function, \(P_{\rm over}\) such that the result \(P_{\rm over}\geq P\).
 
virtual double ratioP (const double z, const Energy2 t, const IdList &ids, const bool mass, const RhoDMatrix &rho) const =0
 Purely virtual method which should return the ratio of the splitting function to the overestimate, i.e.
 
virtual double integOverP (const double z, const IdList &ids, unsigned int PDFfactor=0) const =0
 Purely virtual method which should return the indefinite integral of the overestimated splitting function, \(P_{\rm over}\).
 
virtual double invIntegOverP (const double r, const IdList &ids, unsigned int PDFfactor=0) const =0
 Purely virtual method which should return the inverse of the indefinite integral of the overestimated splitting function, \(P_{\rm over}\) which is used to generate the value of \(z\).
 
virtual void colourConnection (tShowerParticlePtr parent, tShowerParticlePtr first, tShowerParticlePtr second, ShowerPartnerType partnerType, const bool back) const
 Purely virtual method which should make the proper colour connection between the emitting parent and the branching products.
 
virtual vector< pair< int, Complex > > generatePhiForward (const double z, const Energy2 t, const IdList &ids, const RhoDMatrix &)=0
 Method to calculate the azimuthal angle for forward evolution.
 
virtual vector< pair< int, Complex > > generatePhiBackward (const double z, const Energy2 t, const IdList &ids, const RhoDMatrix &)=0
 Method to calculate the azimuthal angle for backward evolution.
 
virtual DecayMEPtr matrixElement (const double z, const Energy2 t, const IdList &ids, const double phi, bool timeLike)=0
 Calculate the matrix element for the splitting.
 
bool angularOrdered () const
 Whether or not the interaction is angular ordered.
 
bool pTScale () const
 Scale choice.
 
void evaluateFinalStateScales (ShowerPartnerType type, Energy scale, double z, tShowerParticlePtr parent, tShowerParticlePtr first, tShowerParticlePtr second)
 Functions to state scales after branching happens.
 
void evaluateInitialStateScales (ShowerPartnerType type, Energy scale, double z, tShowerParticlePtr parent, tShowerParticlePtr first, tShowerParticlePtr second)
 Sort out scales for initial-state emission.
 
void evaluateDecayScales (ShowerPartnerType type, Energy scale, double z, tShowerParticlePtr parent, tShowerParticlePtr first, tShowerParticlePtr second)
 Sort out scales for decay emission.
 
- Public Member Functions inherited from ThePEG::Interfaced
virtual bool defaultInit ()
 
PPtr getParticle (PID) const
 
PDPtr getParticleData (PID) const
 
bool used () const
 
void useMe () const
 
tEGPtr generator () const
 
void persistentOutput (PersistentOStream &os) const
 
void persistentInput (PersistentIStream &is, int version)
 
PPtr getParticle (PID) const
 
PDPtr getParticleData (PID) const
 
bool used () const
 
void useMe () const
 
tEGPtr generator () const
 
- Public Member Functions inherited from ThePEG::InterfacedBase
string fullName () const
 
string name () const
 
string path () const
 
string comment () const
 
void setup (istream &is)
 
void update ()
 
void init ()
 
virtual bool preInitialize () const
 
void initrun ()
 
void finish ()
 
void touch ()
 
void reset ()
 
void clear ()
 
InitState state () const
 
bool locked () const
 
bool touched () const
 
virtual IBPtr fullclone () const
 
void persistentOutput (PersistentOStream &os) const
 
void persistentInput (PersistentIStream &is, int version)
 
virtual void debugme () const
 
void update ()
 
void init ()
 
virtual bool preInitialize () const
 
void initrun ()
 
void finish ()
 
void touch ()
 
void reset ()
 
void clear ()
 
InitState state () const
 
bool locked () const
 
bool touched () const
 
virtual IBPtr fullclone () const
 
- Public Member Functions inherited from ThePEG::Base
void debug () const
 
virtual void debugme () const
 
- Public Member Functions inherited from ThePEG::Pointer::ReferenceCounted
CounterType referenceCount () const
 
- Public Member Functions inherited from ThePEG::Named
 Named (const string &newName=string())
 
 Named (const Named &)=default
 
const string & name () const
 
bool operator== (const Named &other) const
 
bool operator< (const Named &other) const
 

Standard Interfaced functions.

ShowerInteraction _interactionType
 The interaction type for the splitting function.
 
ColourStructure _colourStructure
 The colour structure.
 
double _colourFactor
 The colour factor.
 
bool angularOrdered_
 Whether or not this interaction is angular-ordered.
 
unsigned int scaleChoice_
 The choice of scale.
 
bool strictAO_
 Enforce strict AO.
 
virtual void doinit ()
 Initialize this object after the setup phase before saving an EventGenerator to disk.
 
void colourFactor (double in)
 Set the colour factor.
 
SplittingFunctionoperator= (const SplittingFunction &)=delete
 The assignment operator is private and must never be called.
 

Functions used by the persistent I/O system.

void persistentOutput (PersistentOStream &os) const
 Function used to write out object persistently.
 
void persistentInput (PersistentIStream &is, int version)
 Function used to read in object persistently.
 
static void Init ()
 The standard Init function used to initialize the interfaces.
 

Additional Inherited Members

- Public Types inherited from ThePEG::InterfacedBase
enum  InitState
 
- Public Types inherited from ThePEG::Pointer::ReferenceCounted
typedef unsigned int CounterType
 
- Static Public Member Functions inherited from ThePEG::Interfaced
static void Init ()
 
- Static Public Member Functions inherited from ThePEG::InterfacedBase
static void Init ()
 
- Static Public Member Functions inherited from ThePEG::Base
static void Init ()
 
- Public Attributes inherited from ThePEG::InterfacedBase
 initializing
 
 uninitialized
 
 initialized
 
 runready
 
- Public Attributes inherited from ThePEG::Pointer::ReferenceCounted
const unsigned long uniqueId
 
- Protected Member Functions inherited from ThePEG::Interfaced
void reporeg (IBPtr object, string name) const
 
bool setDefaultReference (PtrT &ptr, string classname, string objectname)
 
 Interfaced (const string &newName)
 
 Interfaced (const Interfaced &i)
 
void setGenerator (tEGPtr generator)
 
- Protected Member Functions inherited from ThePEG::InterfacedBase
virtual void readSetup (istream &is)
 
virtual void doupdate ()
 
virtual void doinit ()
 
virtual void doinitrun ()
 
virtual void dofinish ()
 
virtual IVector getReferences ()
 
virtual void rebind (const TranslationMap &)
 
virtual IBPtr clone () const=0
 
 InterfacedBase (string newName)
 
 InterfacedBase (const InterfacedBase &i)
 
virtual void readSetup (istream &is)
 
virtual void doupdate ()
 
virtual void doinit ()
 
virtual void doinitrun ()
 
virtual void dofinish ()
 
virtual IVector getReferences ()
 
virtual void rebind (const TranslationMap &)
 
- Protected Member Functions inherited from ThePEG::Pointer::ReferenceCounted
 ReferenceCounted (const ReferenceCounted &)
 
ReferenceCountedoperator= (const ReferenceCounted &)
 
- Protected Member Functions inherited from ThePEG::Named
const Namedoperator= (const Named &other)
 
const string & name (const string &newName)
 
- Static Protected Member Functions inherited from ThePEG::Interfaced
static void registerRepository (IBPtr)
 
static void registerRepository (IBPtr, string newName)
 

Detailed Description

This is an abstract class which defines the common interface for all \(1\to2\) splitting functions, for both initial-state and final-state radiation.

The SplittingFunction class contains a number of purely virtual members which must be implemented in the inheriting classes. The class also stores the interaction type of the spltting function.

The inheriting classes need to specific the splitting function \(P(z,2p_j\cdot p_k)\), in terms of the energy fraction \(z\) and the evolution scale. In order to allow the splitting functions to be used with different choices of evolution functions the scale is given by

\[2p_j\cdot p_k=(p_j+p_k)^2-m_{jk}^2=Q^2-(p_j+p_k)^2=z(1-z)\tilde{q}^2= \frac{p_T^2}{z(1-z)}-m_{jk}^2+\frac{m_j^2}{z}+\frac{m_k^2}{1-z},\]

where \(Q^2\) is the virtuality of the branching particle, $p_T$ is the relative transverse momentum of the branching products and \(\tilde{q}^2\) is the angular variable described in hep-ph/0310083.

In addition an overestimate of the splitting function, \(P_{\rm over}(z)\) which only depends upon \(z\), the integral, inverse of the integral for this overestimate and ratio of the true splitting function to the overestimate must be provided as they are necessary for the veto alogrithm used to implement the evolution.

See also
The interfaces defined for SplittingFunction.

Definition at line 68 of file SplittingFunction.h.

Constructor & Destructor Documentation

◆ SplittingFunction()

Herwig::SplittingFunction::SplittingFunction ( )
inline

The default constructor.


Parameters
bAll splitting functions must have an interaction order

Definition at line 76 of file SplittingFunction.h.

Member Function Documentation

◆ accept()

virtual bool Herwig::SplittingFunction::accept ( const IdList ids) const
pure virtual

Purely virtual method which should determine whether this splitting function can be used for a given set of particles.

Parameters
idsThe PDG codes for the particles in the splitting.

Implemented in Herwig::HalfHalfOneEWSplitFn, Herwig::HalfHalfOneSplitFn, Herwig::HalfHalfZeroEWSplitFn, Herwig::HalfOneHalfSplitFn, Herwig::OneHalfHalfSplitFn, Herwig::OneOneOneEWSplitFn, Herwig::OneOneOneMassiveSplitFn, Herwig::OneOneOneQEDSplitFn, Herwig::OneOneOneSplitFn, Herwig::OneOneZeroEWSplitFn, and Herwig::ZeroZeroOneSplitFn.

◆ angularOrdered()

bool Herwig::SplittingFunction::angularOrdered ( ) const
inline

Whether or not the interaction is angular ordered.

Definition at line 262 of file SplittingFunction.h.

References angularOrdered_.

◆ colourConnection()

virtual void Herwig::SplittingFunction::colourConnection ( tShowerParticlePtr  parent,
tShowerParticlePtr  first,
tShowerParticlePtr  second,
ShowerPartnerType  partnerType,
const bool  back 
) const
virtual

Purely virtual method which should make the proper colour connection between the emitting parent and the branching products.

Parameters
parentThe parent for the branching
firstThe first branching product
secondThe second branching product
partnerTypeThe type of evolution partner
backWhether this is foward or backward evolution.

◆ colourFactor() [1/2]

double Herwig::SplittingFunction::colourFactor ( const IdList ids) const
inline

Return the colour factor.

Definition at line 100 of file SplittingFunction.h.

References _colourFactor, _colourStructure, and ThePEG::sqr().

◆ colourFactor() [2/2]

void Herwig::SplittingFunction::colourFactor ( double  in)
inlineprotected

Set the colour factor.

Definition at line 345 of file SplittingFunction.h.

References _colourFactor.

◆ colourStructure()

ColourStructure Herwig::SplittingFunction::colourStructure ( ) const
inline

Return the colour structure.

Definition at line 95 of file SplittingFunction.h.

References _colourStructure.

◆ doinit()

virtual void Herwig::SplittingFunction::doinit ( )
protectedvirtual

Initialize this object after the setup phase before saving an EventGenerator to disk.

Exceptions
InitExceptionif object could not be initialized properly.

Reimplemented from ThePEG::InterfacedBase.

Reimplemented in Herwig::CMWHalfHalfOneSplitFn, Herwig::CMWOneOneOneSplitFn, Herwig::HalfHalfOneEWSplitFn, Herwig::HalfHalfZeroEWSplitFn, Herwig::OneOneOneEWSplitFn, Herwig::OneOneOneQEDSplitFn, and Herwig::OneOneZeroEWSplitFn.

◆ evaluateFinalStateScales()

void Herwig::SplittingFunction::evaluateFinalStateScales ( ShowerPartnerType  type,
Energy  scale,
double  z,
tShowerParticlePtr  parent,
tShowerParticlePtr  first,
tShowerParticlePtr  second 
)

Functions to state scales after branching happens.

Sort out scales for final-state emission

◆ generatePhiBackward()

virtual vector< pair< int, Complex > > Herwig::SplittingFunction::generatePhiBackward ( const double  z,
const Energy2  t,
const IdList ids,
const RhoDMatrix  
)
pure virtual

Method to calculate the azimuthal angle for backward evolution.

Parameters
zThe energy fraction
tThe scale \(t=2p_j\cdot p_k\).
idsThe PDG codes for the particles in the splitting.
Returns
The weight

Implemented in Herwig::HalfHalfOneEWSplitFn, Herwig::HalfHalfOneSplitFn, Herwig::HalfHalfZeroEWSplitFn, Herwig::HalfOneHalfSplitFn, Herwig::OneHalfHalfSplitFn, Herwig::OneOneOneEWSplitFn, Herwig::OneOneOneMassiveSplitFn, Herwig::OneOneOneQEDSplitFn, Herwig::OneOneOneSplitFn, Herwig::OneOneZeroEWSplitFn, and Herwig::ZeroZeroOneSplitFn.

◆ generatePhiForward()

virtual vector< pair< int, Complex > > Herwig::SplittingFunction::generatePhiForward ( const double  z,
const Energy2  t,
const IdList ids,
const RhoDMatrix  
)
pure virtual

Method to calculate the azimuthal angle for forward evolution.

Parameters
zThe energy fraction
tThe scale \(t=2p_j\cdot p_k\).
idsThe PDG codes for the particles in the splitting.
Theazimuthal angle, \(\phi\).
Returns
The weight

Implemented in Herwig::HalfHalfOneEWSplitFn, Herwig::HalfHalfOneSplitFn, Herwig::HalfHalfZeroEWSplitFn, Herwig::HalfOneHalfSplitFn, Herwig::OneHalfHalfSplitFn, Herwig::OneOneOneEWSplitFn, Herwig::OneOneOneMassiveSplitFn, Herwig::OneOneOneQEDSplitFn, Herwig::OneOneOneSplitFn, Herwig::OneOneZeroEWSplitFn, and Herwig::ZeroZeroOneSplitFn.

◆ Init()

static void Herwig::SplittingFunction::Init ( )
static

The standard Init function used to initialize the interfaces.

Called exactly once for each class by the class description system before the main function starts or when this class is dynamically loaded.

◆ integOverP()

virtual double Herwig::SplittingFunction::integOverP ( const double  z,
const IdList ids,
unsigned int  PDFfactor = 0 
) const
pure virtual

Purely virtual method which should return the indefinite integral of the overestimated splitting function, \(P_{\rm over}\).

Parameters
zThe energy fraction.
idsThe PDG codes for the particles in the splitting.
PDFfactorWhich additional factor to include for the PDF 0 is no additional factor, 1 is \(1/z\), 2 is \(1/(1-z)\) and 3 is \(1/z/(1-z)\)

Implemented in Herwig::HalfHalfOneEWSplitFn, Herwig::HalfHalfOneSplitFn, Herwig::HalfHalfZeroEWSplitFn, Herwig::HalfOneHalfSplitFn, Herwig::OneHalfHalfSplitFn, Herwig::OneOneOneEWSplitFn, Herwig::OneOneOneMassiveSplitFn, Herwig::OneOneOneQEDSplitFn, Herwig::OneOneOneSplitFn, Herwig::OneOneZeroEWSplitFn, and Herwig::ZeroZeroOneSplitFn.

◆ interactionType()

ShowerInteraction Herwig::SplittingFunction::interactionType ( ) const
inline

Methods to return the interaction type and order for the splitting function.

Return the type of the interaction

Definition at line 90 of file SplittingFunction.h.

References _interactionType.

◆ invIntegOverP()

virtual double Herwig::SplittingFunction::invIntegOverP ( const double  r,
const IdList ids,
unsigned int  PDFfactor = 0 
) const
pure virtual

Purely virtual method which should return the inverse of the indefinite integral of the overestimated splitting function, \(P_{\rm over}\) which is used to generate the value of \(z\).

Parameters
rValue of the splitting function to be inverted
idsThe PDG codes for the particles in the splitting.
PDFfactorWhich additional factor to include for the PDF 0 is no additional factor, 1 is \(1/z\), 2 is \(1/(1-z)\) and 3 is \(1/z/(1-z)\)

Implemented in Herwig::HalfHalfOneEWSplitFn, Herwig::HalfHalfOneSplitFn, Herwig::HalfHalfZeroEWSplitFn, Herwig::HalfOneHalfSplitFn, Herwig::OneHalfHalfSplitFn, Herwig::OneOneOneEWSplitFn, Herwig::OneOneOneMassiveSplitFn, Herwig::OneOneOneQEDSplitFn, Herwig::OneOneOneSplitFn, Herwig::OneOneZeroEWSplitFn, and Herwig::ZeroZeroOneSplitFn.

◆ matrixElement()

virtual DecayMEPtr Herwig::SplittingFunction::matrixElement ( const double  z,
const Energy2  t,
const IdList ids,
const double  phi,
bool  timeLike 
)
pure virtual

Calculate the matrix element for the splitting.

Parameters
zThe energy fraction
tThe scale \(t=2p_j\cdot p_k\).
idsThe PDG codes for the particles in the splitting.
phiThe azimuthal angle, \(\phi\).
timeLikeWhether timelike or spacelike, affects inclusive of mass terms

Implemented in Herwig::HalfHalfOneEWSplitFn, Herwig::HalfHalfOneSplitFn, Herwig::HalfHalfZeroEWSplitFn, Herwig::HalfOneHalfSplitFn, Herwig::OneHalfHalfSplitFn, Herwig::OneOneOneEWSplitFn, Herwig::OneOneOneMassiveSplitFn, Herwig::OneOneOneQEDSplitFn, Herwig::OneOneOneSplitFn, Herwig::OneOneZeroEWSplitFn, and Herwig::ZeroZeroOneSplitFn.

◆ operator=()

SplittingFunction & Herwig::SplittingFunction::operator= ( const SplittingFunction )
privatedelete

The assignment operator is private and must never be called.

In fact, it should not even be implemented.

◆ overestimateP()

virtual double Herwig::SplittingFunction::overestimateP ( const double  z,
const IdList ids 
) const
pure virtual

Purely virtual method which should return an overestimate of the splitting function, \(P_{\rm over}\) such that the result \(P_{\rm over}\geq P\).

This function should be simple enough that it does not depend on the evolution scale.

Parameters
zThe energy fraction.
idsThe PDG codes for the particles in the splitting.

Implemented in Herwig::HalfHalfOneEWSplitFn, Herwig::HalfHalfOneSplitFn, Herwig::HalfHalfZeroEWSplitFn, Herwig::HalfOneHalfSplitFn, Herwig::OneHalfHalfSplitFn, Herwig::OneOneOneEWSplitFn, Herwig::OneOneOneMassiveSplitFn, Herwig::OneOneOneQEDSplitFn, Herwig::OneOneOneSplitFn, Herwig::OneOneZeroEWSplitFn, and Herwig::ZeroZeroOneSplitFn.

◆ P()

virtual double Herwig::SplittingFunction::P ( const double  z,
const Energy2  t,
const IdList ids,
const bool  mass,
const RhoDMatrix rho 
) const
pure virtual

Methods to return the splitting function.

Purely virtual method which should return the exact value of the splitting function, \(P\) evaluated in terms of the energy fraction, \(z\), and the evolution scale \(\tilde{q}^2\).

Parameters
zThe energy fraction.
tThe scale \(t=2p_j\cdot p_k\).
idsThe PDG codes for the particles in the splitting.
massWhether or not to include the mass dependent terms
rhoThe spin density matrix

Implemented in Herwig::ZeroZeroOneSplitFn, Herwig::CMWHalfHalfOneSplitFn, Herwig::CMWOneOneOneSplitFn, Herwig::HalfHalfOneEWSplitFn, Herwig::HalfHalfOneSplitFn, Herwig::HalfHalfZeroEWSplitFn, Herwig::HalfOneHalfSplitFn, Herwig::OneHalfHalfSplitFn, Herwig::OneOneOneEWSplitFn, Herwig::OneOneOneMassiveSplitFn, Herwig::OneOneOneQEDSplitFn, Herwig::OneOneOneSplitFn, and Herwig::OneOneZeroEWSplitFn.

◆ persistentInput()

void Herwig::SplittingFunction::persistentInput ( PersistentIStream is,
int  version 
)

Function used to read in object persistently.

Parameters
isthe persistent input stream read from.
versionthe version number of the object when written.

◆ persistentOutput()

void Herwig::SplittingFunction::persistentOutput ( PersistentOStream os) const

Function used to write out object persistently.

Parameters
osthe persistent output stream written to.

◆ pTScale()

bool Herwig::SplittingFunction::pTScale ( ) const
inline

Scale choice.

Definition at line 267 of file SplittingFunction.h.

References angularOrdered_, and scaleChoice_.

◆ ratioP()

virtual double Herwig::SplittingFunction::ratioP ( const double  z,
const Energy2  t,
const IdList ids,
const bool  mass,
const RhoDMatrix rho 
) const
pure virtual

Purely virtual method which should return the ratio of the splitting function to the overestimate, i.e.

\(P(z,\tilde{q}^2)/P_{\rm over}(z)\).

Parameters
zThe energy fraction.
tThe scale \(t=2p_j\cdot p_k\).
idsThe PDG codes for the particles in the splitting.
massWhether or not to include the mass dependent terms
rhoThe spin density matrix

Implemented in Herwig::ZeroZeroOneSplitFn, Herwig::CMWHalfHalfOneSplitFn, Herwig::CMWOneOneOneSplitFn, Herwig::HalfHalfOneEWSplitFn, Herwig::HalfHalfOneSplitFn, Herwig::HalfHalfZeroEWSplitFn, Herwig::HalfOneHalfSplitFn, Herwig::OneHalfHalfSplitFn, Herwig::OneOneOneEWSplitFn, Herwig::OneOneOneMassiveSplitFn, Herwig::OneOneOneQEDSplitFn, Herwig::OneOneOneSplitFn, and Herwig::OneOneZeroEWSplitFn.

Member Data Documentation

◆ _colourFactor

double Herwig::SplittingFunction::_colourFactor
private

The colour factor.

Definition at line 370 of file SplittingFunction.h.

Referenced by colourFactor().

◆ _colourStructure

ColourStructure Herwig::SplittingFunction::_colourStructure
private

The colour structure.

Definition at line 365 of file SplittingFunction.h.

Referenced by colourFactor(), and colourStructure().

◆ _interactionType

ShowerInteraction Herwig::SplittingFunction::_interactionType
private

The interaction type for the splitting function.

Definition at line 360 of file SplittingFunction.h.

Referenced by interactionType().

◆ angularOrdered_

bool Herwig::SplittingFunction::angularOrdered_
private

Whether or not this interaction is angular-ordered.

Definition at line 375 of file SplittingFunction.h.

Referenced by angularOrdered(), and pTScale().

◆ scaleChoice_

unsigned int Herwig::SplittingFunction::scaleChoice_
private

The choice of scale.

Definition at line 380 of file SplittingFunction.h.

Referenced by pTScale().

◆ strictAO_

bool Herwig::SplittingFunction::strictAO_
private

Enforce strict AO.

Definition at line 385 of file SplittingFunction.h.


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