PYTHIA  8.311
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
SplitOnia Class Referenceabstract

#include <SplittingsOnia.h>

Inheritance diagram for SplitOnia:
Split2g2QQbar1S01g Split2g2QQbar3PJ1g Split2g2QQbar3S11gg Split2g2QQbarX8 Split2Q2QQbar1P11Q Split2Q2QQbar1S01Q Split2Q2QQbar3PJ1Q Split2Q2QQbar3PJ8Q Split2Q2QQbar3S11Q Split2QQbarXq82QQbarX8g

Public Member Functions

 SplitOnia (int idAIn, int idBIn, int idCIn, double ldmeIn, Info *infoPtrIn, AlphaStrong *alphaSPtrIn)
 Constructor that needs to be called by the subclass.
 
virtual ~SplitOnia ()=default
 Virtual destructor.
 
virtual double overestimate (const TimeDipoleEnd &dip, double pT2Min, bool enh)
 Calculate the splitting overestimate. More...
 
virtual double weight (const TimeDipoleEnd &dip) const =0
 
virtual double generateZ (const TimeDipoleEnd &)
 
bool isActive (const TimeDipoleEnd &dip, int id, double m)
 
bool isOctet ()
 Check if the emission is a colour octet state.
 
bool isEnhanced ()
 Check if the emission is enhanced.
 
void setEnhance (const unordered_map< string, double > &enhanceFSRs)
 Set the enhancement factor.
 
double enhanceWeight ()
 Return the enhancement weight.
 
virtual string enhanceName () const
 Return the name to be used if this emission is enhanced.
 
virtual vector< int > addEmitted (Event &, int, int, int, int, vector< TimeDipoleEnd > &)
 Possibility to add more than one emitted particle.
 
virtual void updateDipole (TimeDipoleEnd &dip)
 Update the dipole end with the generated values.
 
bool updateBranchVars (const TimeDipoleEnd *dip, Event &event, int &idRadIn, int &idEmtIn, int &colRadIn, int &acolRadIn, int &colEmtIn, int &acolEmtIn, int &appendEmtIn, double &pTorigIn, double &pTcorrIn, double &pzRadPlusEmtIn, double &pzRadIn, double &pzEmtIn, double &mRadIn, double &m2RadIn, double &mEmtIn)
 Update the passed branching variables with the internal values. More...
 

Protected Member Functions

virtual bool kinematics (const TimeDipoleEnd *dip, Event &event)
 Update the internal kinematics. More...
 
virtual void overestimate (const TimeDipoleEnd &, double)
 Set the z limits and return the difference.
 
virtual double integrateZ () const
 
void setOctetID (int state, double mSplit, Info *infoPtr)
 Set the colour octet ID and ensure in particle database. More...
 
double alphaScale (double m2, double pT2, double s) const
 Return alpha_s at the requested scale.
 

Protected Attributes

int idA
 
int idB
 
int idC
 
double mA
 
double mB
 
double mC
 
double m2A
 
double m2B
 
double m2C
 
double enhance {1}
 Enhancement, max value of alpha_S, and long-distance matrix element.
 
double ldme {-1}
 
double cFac {0}
 The constant and variable overestimate prefactors.
 
double oFac {0}
 
double zMin {0}
 The z limits used in the overestimated z-integral and generated z value.
 
double zMax {1}
 
double zGen {0}
 
int idRad {0}
 IDs, colors, and number of emitters to append.
 
int idEmt {0}
 
int colRad {0}
 
int acolRad {0}
 
int colEmt {0}
 
int acolEmt {0}
 
int appendEmt {1}
 
double pTorig {0}
 Transverse and longitudinal momentum, and masses.
 
double pTcorr {0}
 
double pzRadPlusEmt {0}
 
double pzRad {0}
 
double pzEmt {0}
 
double mRad {0}
 
double m2Rad {0}
 
double mEmt {0}
 
int alphaMode {1}
 Mode to evaluate the final alpha_s scale.
 
LoggerloggerPtr {}
 Logger pointer.
 
AlphaStrongalphaSPtr {}
 The alphaS object in current use.
 
RndmrndmPtr {}
 The random number generator in current use.
 

Detailed Description

The SplitOnia encodes the generation of a particular emission of an onium state from a parton in a dipole.

The basic generation in the shower assumes a function of the form dP(pT2, z) = alphaS/2pi dpT2/pT2 dz F(pT2, z) in an A -> B + C splitting, where B takes a fraction z of the momentum of A.

A subclass of this base clase should encode not only the function F(pT2, z), but also an overestimate and the z-integral of this overestimate:

G(pT2, z) >= F(pT2, z) and OIZ = ^1 G(pT2, z) dz

The ratio F/G for a generated point should then be used as a weight for the phase space point generated. The subclass should also be able to generate the z according to the overestimate and possible other internal variables that have been integrated out.

Member Function Documentation

virtual double generateZ ( const TimeDipoleEnd )
inlinevirtual

Generate internal degrees of freedom (in particluar z). Default behaviour is to generate z from a flat distribution.

Reimplemented in Split2QQbarXg82QQbarX8g, Split2QQbarXq82QQbarX8g, and Split2g2QQbar3S11gg.

virtual double integrateZ ( ) const
inlineprotectedvirtual

Integrate the distribution used to generate z. Default behaviour is to integrate a flat distribution.

bool isActive ( const TimeDipoleEnd dip,
int  id,
double  m 
)
inline

Check if this emission is available for the current dipole end with the given ID of the radiator.

bool kinematics ( const TimeDipoleEnd dip,
Event event 
)
protectedvirtual

Update the internal kinematics.

Set the kinematics. Suitable for 1 -> 2 onium emissions.

Reimplemented in Split2g2QQbar3S11gg.

double overestimate ( const TimeDipoleEnd dip,
double  pT2Min,
bool  enh 
)
virtual

Calculate the splitting overestimate.

Implementation of the SplitOnia class.

Set up the the overestimate. This calculates the splitting independent factor.

Reimplemented in Split2g2QQbarX8.

void setOctetID ( int  state,
double  mSplit,
Info infoPtr 
)
protected

Set the colour octet ID and ensure in particle database.

Set the colour octet ID and ensure in particle database. Here, the state corresponds to 0 is 3S1, 1 is 1S0, and 2 is 3PJ.

Determine quark composition and quantum numbers of the physical state.

Set the name of the physical and octet state.

Ensure the dummy particle for the colour-octet state is valid.

bool updateBranchVars ( const TimeDipoleEnd dip,
Event event,
int &  idRadIn,
int &  idEmtIn,
int &  colRadIn,
int &  acolRadIn,
int &  colEmtIn,
int &  acolEmtIn,
int &  appendEmtIn,
double &  pTorigIn,
double &  pTcorrIn,
double &  pzRadPlusEmtIn,
double &  pzRadIn,
double &  pzEmtIn,
double &  mRadIn,
double &  m2RadIn,
double &  mEmtIn 
)

Update the passed branching variables with the internal values.

Set branch variables. Call class specific version first, then set passed variables.

virtual double weight ( const TimeDipoleEnd dip) const
pure virtual

Return the weight of the splitting function and Jacobian at the generated point divided by the overestimate.

Implemented in Split2QQbarXg82QQbarX8g, Split2QQbarXq82QQbarX8g, Split2g2QQbarX8, Split2Q2QQbar3PJ8Q, Split2g2QQbar3PJ1g, Split2Q2QQbar3PJ1Q, Split2Q2QQbar1P11Q, Split2g2QQbar3S11gg, Split2Q2QQbar3S11Q, Split2g2QQbar1S01g, and Split2Q2QQbar1S01Q.

Member Data Documentation

int idA
protected

The ID and mass (squared) of the radiating particle (A), the scattered radiated particle (B, taking a fraction z of A's momentum), and the emitted particle (C).


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