PYTHIA  8.311
Public Member Functions | Public Attributes | Static Public Attributes | List of all members
StringEnd Class Reference

#include <StringFragmentation.h>

Public Member Functions

 StringEnd ()
 Constructor.
 
void init (ParticleData *particleDataPtrIn, StringFlav *flavSelPtrIn, StringPT *pTSelPtrIn, StringZ *zSelPtrIn, Settings &settings)
 Save pointers.
 
void setUp (bool fromPosIn, int iEndIn, int idOldIn, int iMaxIn, double pxIn, double pyIn, double GammaIn, double xPosIn, double xNegIn, int colIn)
 Set up initial endpoint values from input. More...
 
void newHadron (double kappaRatio, bool forbidPopcornNow=false, bool allowPop=true, double strangeFac=0., double probQQmod=1.)
 Fragment off one hadron from the string system, in flavour and pT. More...
 
void pearlHadron (StringSystem &system, int idPearlIn, Vec4 pPearlIn)
 Creation of pearl hadron. More...
 
Vec4 kinematicsHadron (StringSystem &system, StringVertex &newVertex, bool useInputZ=false, double zHadIn=0., bool pearlIn=false, Vec4 pPearlIn={0., 0., 0., 0.})
 
Vec4 kinematicsHadronTmp (StringSystem system, Vec4 pRem, double phi, double mult)
 
void update ()
 Update string end information after a hadron has been removed.
 
void storePrev ()
 Store old string end information.
 
void updateToPrev ()
 Update string end information to previous string break.
 

Public Attributes

ParticleDataparticleDataPtr
 Pointer to the particle data table.
 
StringFlavflavSelPtr
 Pointers to classes for flavour, pT and z generation.
 
StringPTpTSelPtr
 
StringZzSelPtr
 
StringFlav flavSelNow
 Local copy of flavSelPtr for modified flavour selection.
 
bool fromPos
 Data members.
 
bool thermalModel
 
bool mT2suppression
 
bool closePacking
 
int iEnd
 
int iMax
 
int idHad
 
int iPosOld
 
int iNegOld
 
int iPosNew
 
int iNegNew
 
int hadSoFar
 
int colOld
 
int colNew
 
double pxOld
 
double pyOld
 
double pxNew
 
double pyNew
 
double pxHad
 
double pyHad
 
double mHad
 
double mT2Had
 
double zHad
 
double GammaOld
 
double GammaNew
 
double xPosOld
 
double xPosNew
 
double xPosHad
 
double xNegOld
 
double xNegNew
 
double xNegHad
 
double aLund
 
double bLund
 
int iPosOldPrev
 
int iNegOldPrev
 
int colOldPrev
 
double pxOldPrev
 
double pyOldPrev
 
double GammaOldPrev
 
double xPosOldPrev
 
double xNegOldPrev
 
FlavContainer flavOld
 
FlavContainer flavNew
 
FlavContainer flavOldPrev
 
Vec4 pHad
 
Vec4 pSoFar
 

Static Public Attributes

static const double TINY = 1e-6
 Constants: could only be changed in the code itself. More...
 
static const double PT2SAME = 0.01
 Assume two (eX, eY) regions are related if pT2 differs by less.
 
static const double MEANMMIN = 0.2
 
static const double MEANM = 0.8
 
static const double MEANPT = 0.4
 

Detailed Description

The StringEnd class contains the information related to one of the current endpoints of the string system. Only to be used inside StringFragmentation, so no private members.

Member Function Documentation

Vec4 kinematicsHadron ( StringSystem system,
StringVertex newVertex,
bool  useInputZ = false,
double  zHadIn = 0.,
bool  pearlIn = false,
Vec4  pPearlIn = { 0., 0., 0., 0.} 
)

Fragment off one hadron from the string system, in momentum space, by taking steps either from positive or from negative end.

Fragment off one hadron from the string system, in momentum space, by taking steps from positive end.

Ensure pearl momentum is not used if pearlIn = false.

Pick fragmentation step z and calculate new Gamma.

Set up references that are direction-neutral; ...Dir for direction of iteration and ...Inv for its inverse.

Start search for new breakup in the old region.

Each step corresponds to trying a new string region.

Referance to current string region.

Now begin special section for rapid processing of low region.

A first step within a low region is easy.

Translate into x coordinates.

Store breakup vertex information from the fragmentation process.

Find and return four-momentum of the produced particle.

A first step out of a low region also OK, if there are more regions. Negative energy signals failure, i.e. in last region.

Momentum taken by stepping out of region. Continue to next region.

Else, for first step, take into account starting pT.

Now begin normal treatment of nontrivial regions. Set up four-vectors in a region not visited before.

If new region is vanishingly small, continue immediately to next. Negative energy signals failure to do this, i.e. moved too low.

Reexpress pTNew w.r.t. base vectors in new region, if possible. Recall minus sign from normalization e_x^2 = e_y^2 = -1.

Four-momentum taken so far, including new pT.

Derive coefficients for m2 expression. cM2 * x+ + cM3 * x- + cM4 * x+ * x- = m^2 - cM1;

Derive coefficients for Gamma expression. cGam2 * x+ + cGam3 * x- + cGam4 * x+ * x- = Gamma_new - cGam1;

Solve (m2, Gamma) equation system => r2 * x-^2 + r1 * x- + r0 = 0.

Define position of new trial vertex.

Step up to new region if new x- > 1.

Step down to new region if new x+ < 0.

Store breakup vertex information from the fragmentation process.

Else we have found the correct region, and can set the new colour index and return four-momentum.

End of "infinite" loop of stepping to new region.

Vec4 kinematicsHadronTmp ( StringSystem  system,
Vec4  pRem,
double  phi,
double  mult 
)

Generate momentum for some possible next hadron, based on mean values to get an estimate for rapidity and pT.

Now estimate the energy the next hadron will take.

Modify Gamma value in case of earlier fails.

Special case of first hadron.

Negative energy signals failure.

First make a copy of all variables to not overwrite anything.

Set up references that are direction-neutral; ...Dir for direction of iteration and ...Inv for its inverse.

Start search for new breakup in the old region.

Each step corresponds to trying a new string region.

Referance to current string region.

Now begin special section for rapid processing of low region.

A first step within a low region is easy. Make sure we use this region in case it's the last one.

Translate into x coordinates.

Find and return four-momentum of the produced particle.

A first step out of a low region also OK, if there are more regions. Negative energy signals failure, i.e. in last region.

Should be covered by the above check.

Momentum taken by stepping out of region. Continue to next region.

Else, for first step, take into account starting pT.

Now begin normal treatment of nontrivial regions. Set up four-vectors in a region not visited before.

If new region is vanishingly small, continue immediately to next. Negative energy signals failure to do this, i.e. moved too low.

Reexpress pTNew w.r.t. base vectors in new region, if possible. Recall minus sign from normalization e_x^2 = e_y^2 = -1.

Four-momentum taken so far, including new pT.

Derive coefficients for m2 expression. cM2 * x+ + cM3 * x- + cM4 * x+ * x- = m^2 - cM1;

Derive coefficients for Gamma expression. cGam2 * x+ + cGam3 * x- + cGam4 * x+ * x- = Gamma_new - cGam1;

Solve (m2, Gamma) equation system => r2 * x-^2 + r1 * x- + r0 = 0.

Define position of new trial vertex.

Step up to new region if new x- > 1.

Step down to new region if new x+ < 0.

Else we have found the correct region, and can return four-momentum.

End of "infinite" loop of stepping to new region.

void newHadron ( double  kappaRatio,
bool  forbidPopcornNow = false,
bool  allowPop = true,
double  strangeFac = 0.,
double  probQQmod = 1. 
)

Fragment off one hadron from the string system, in flavour and pT.

Avoid contradiction in allowPop and forbidPopcornNow.

In case we are using the thermal model or Gaussian with mT2 suppression we have to pick the pT first.

Pick its transverse momentum.

Pick new flavour and form a new hadron. For forbidPopcornNow == true it must be a baryon.

Get its mass and thereby define its transverse mass.

In case of the Gaussian without mT2 suppression we pick the new flavour first to make the width flavour dependent.

Pick new flavour and form a new hadron. For forbidPopcornNow == true it must be a baryon.

Reinitialise probabilities if close-packing.

Pick its transverse momentum.

Pick its mass and thereby define its transverse mass.

void pearlHadron ( StringSystem system,
int  idPearlIn,
Vec4  pPearlIn 
)

Creation of pearl hadron.

Make (anti)baryon from stepping over pearl (anti)quark.

Find transverse momentum contribution of pearl in current region.

Define baryon including pearl.

void setUp ( bool  fromPosIn,
int  iEndIn,
int  idOldIn,
int  iMaxIn,
double  pxIn,
double  pyIn,
double  GammaIn,
double  xPosIn,
double  xNegIn,
int  colIn 
)

Set up initial endpoint values from input.

Simple transcription from input.

Member Data Documentation

const double MEANMMIN = 0.2
static

Fictitious typical mass and pT used to look ahead for approximately where the next hadron will be produced, to quantify string close-packing there.

const double TINY = 1e-6
static

Constants: could only be changed in the code itself.

The StringEnd class.

Constants: could be changed here if desired, but normally should not. Avoid unphysical solutions to equation system.


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