PYTHIA
8.312
|
The DireTimes class does timelike showers. More...
#include <DireTimes.h>
Public Member Functions | |
DireTimes () | |
Constructor. | |
DireTimes (MergingHooksPtr mergingHooksPtrIn, PartonVertexPtr) | |
virtual | ~DireTimes () |
Destructor. | |
virtual void | init (BeamParticle *beamAPtrIn=nullptr, BeamParticle *beamBPtrIn=nullptr) |
Initialize alphaStrong and related pTmin parameters. More... | |
bool | initSplits () |
void | reinitPtr (Info *infoPtrIn, MergingHooksPtr mergingHooksPtrIn, DireSplittingLibrary *splittingsPtrIn, DireInfo *direInfoPtrIn) |
void | initVariations () |
Initialize bookkeeping of shower variations. More... | |
void | clear () |
Reset parton shower. More... | |
void | setWeightContainerPtr (DireWeightContainer *weightsIn) |
virtual bool | limitPTmax (Event &event, double Q2Fac=0., double Q2Ren=0.) |
Find whether to limit maximum scale of emissions, and whether to dampen. More... | |
virtual double | enhancePTmax () |
Potential enhancement factor of pTmax scale for hardest emission. | |
virtual int | shower (int iBeg, int iEnd, Event &event, double pTmax, int nBranchMax=0) |
Top-level routine to do a full time-like shower in resonance decay. More... | |
virtual int | showerQED (int i1, int i2, Event &event, double pTmax) |
Top-level routine for QED radiation in hadronic decay to two leptons. More... | |
virtual void | prepareGlobal (Event &) |
Global recoil: reset counters and store locations of outgoing partons. More... | |
virtual void | prepare (int iSys, Event &event, bool limitPTmaxIn=true) |
Prepare system for evolution after each new interaction; identify ME. More... | |
void | finalize (Event &event) |
Finalize event after evolution. More... | |
virtual void | rescatterUpdate (int iSys, Event &event) |
Update dipole list after a multiparton interactions rescattering. More... | |
virtual void | update (int iSys, Event &event, bool=false) |
Update dipole list after each ISR emission. More... | |
void | updateAfterFF (int iSysSelNow, int iSysSelRec, Event &event, int iRadBef, int iRecBef, int iRad, int iEmt, int iRec, int flavour, int colType, double pTsel) |
Update dipole list after final-final splitting. More... | |
void | updateAfterFI (int iSysSelNow, int iSysSelRec, Event &event, int iRadBef, int iRecBef, int iRad, int iEmt, int iRec, int flavour, int colType, double pTsel, double xNew) |
Update dipole list after final-final splitting. More... | |
virtual double | pTnext (Event &event, double pTbegAll, double pTendAll, bool=false, bool=false) |
double | newPoint (const Event &event) |
virtual bool | branch (Event &event, bool=false) |
Setup branching kinematics. More... | |
bool | branch_FF (Event &event, bool=false, DireSplitInfo *split=nullptr) |
bool | branch_FI (Event &event, bool=false, DireSplitInfo *split=nullptr) |
pair< Vec4, Vec4 > | decayWithOnshellRec (double zCS, double yCS, double phi, double m2Rec, double m2RadAft, double m2EmtAft, Vec4 pRadBef, Vec4 pRecBef) |
pair< Vec4, Vec4 > | decayWithOffshellRec (double zCS, double yCS, double phi, double m2RadBef, double m2RadAft, double m2EmtAft, Vec4 pRadBef, Vec4 pRecBef) |
bool | getHasWeaklyRadiated () |
Tell whether FSR has done a weak emission. | |
int | system () const |
Tell which system was the last processed one. | |
virtual Event | clustered (const Event &state, int iRad, int iEmt, int iRec, string name) |
Setup clustering kinematics. | |
pair< Event, pair< int, int > > | clustered_internal (const Event &state, int iRad, int iEmt, int iRec, string name) |
bool | cluster_FF (const Event &state, int iRad, int iEmt, int iRec, int idRadBef, Particle &radBef, Particle &recBef) |
bool | cluster_FI (const Event &state, int iRad, int iEmt, int iRec, int idRadBef, Particle &radBef, Particle &recBef) |
double | pT2Times (const Particle &rad, const Particle &emt, const Particle &rec) |
From Pythia version 8.215 onwards no longer virtual. | |
double | pT2_FF (const Particle &rad, const Particle &emt, const Particle &rec) |
double | pT2_FI (const Particle &rad, const Particle &emt, const Particle &rec) |
double | zTimes (const Particle &rad, const Particle &emt, const Particle &rec) |
From Pythia version 8.215 onwards no longer virtual. | |
double | z_FF (const Particle &rad, const Particle &emt, const Particle &rec) |
double | z_FI (const Particle &rad, const Particle &emt, const Particle &rec) |
double | z_FF_fromVec (const Vec4 &rad, const Vec4 &emt, const Vec4 &rec) |
double | m2dipTimes (const Particle &rad, const Particle &emt, const Particle &rec) |
double | m2dip_FF (const Particle &rad, const Particle &emt, const Particle &rec) |
double | m2dip_FI (const Particle &rad, const Particle &emt, const Particle &rec) |
virtual map< string, double > | getStateVariables (const Event &state, int rad, int emt, int rec, string name) |
virtual bool | isTimelike (const Event &state, int iRad, int, int, string) |
virtual vector< string > | getSplittingName (const Event &state, int iRad, int iEmt, int) |
virtual double | getSplittingProb (const Event &state, int iRad, int iEmt, int iRec, string name) |
Compute splitting probability. More... | |
virtual bool | allowedSplitting (const Event &state, int iRad, int iEmt) |
virtual vector< int > | getRecoilers (const Event &state, int iRad, int iEmt, string name) |
virtual double | getCoupling (double mu2Ren, string name) |
bool | isSymmetric (string name, const Particle *rad, const Particle *emt) |
int | FindParticle (const Particle &particle, const Event &event, bool checkStatus=true) |
virtual void | list () const |
Print dipole list; for debug mainly. More... | |
Event | makeHardEvent (int iSys, const Event &state, bool isProcess=false) |
bool | validMomentum (const Vec4 &p, int id, int status) |
Check that particle has sensible momentum. More... | |
bool | validEvent (const Event &state, bool isProcess=false, int iSysCheck=-1) |
Check colour/flavour correctness of state. More... | |
bool | validMotherDaughter (const Event &state) |
Check that mother-daughter-relations are correctly set. More... | |
int | FindCol (int col, vector< int > iExclude, const Event &event, int type, int iSys=-1) |
Find index of colour partner for input colour. More... | |
BeamParticle * | getBeamA () |
Pointers to the two incoming beams. | |
BeamParticle * | getBeamB () |
double | alphasNow (double pT2, double renormMultFacNow=1., int iSys=0) |
double | alphaemNow (double pT2, double renormMultFacNow=1., int iSys=0) |
Function to calculate the correct alphaEM/2*Pi value. More... | |
bool | isInit () |
double | m2Max (int iDip, const Event &state) |
Function to calculate the absolute phase-sace boundary for emissions. | |
Public Member Functions inherited from TimeShower | |
TimeShower ()=default | |
Constructor. | |
virtual | ~TimeShower () |
Destructor. | |
void | initPtrs (MergingHooksPtr mergingHooksPtrIn, PartonVertexPtr partonVertexPtrIn, WeightContainer *weightContainerPtrIn) |
void | reassignBeamPtrs (BeamParticle *beamAPtrIn, BeamParticle *beamBPtrIn, int beamOffsetIn=0) |
New beams possible for handling of hard diffraction. (Not virtual.) | |
virtual int | showerQEDafterRemnants (Event &) |
virtual void | prepareProcess (Event &, Event &, vector< int > &) |
virtual double | pTnextResDec () |
Select next pT for interleaved resonance decays. | |
virtual bool | resonanceShower (Event &, Event &, vector< int > &, double=0.) |
virtual bool | initUncertainties () |
Initialize data members for calculation of uncertainty bands. | |
virtual bool | initEnhancements () |
Initialize data members for application of enhancements. | |
virtual double | pTLastInShower () |
Provide the pT scale of the last branching in the above shower. | |
virtual double | enhanceFactor (const string &name) |
virtual double | noEmissionProbability (double, double, double, int, int, double, double) |
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 | dryrun |
DireWeightContainer * | weights |
DireInfo * | direInfoPtr |
ProcessLevel | processLevel |
unordered_map< string, DireSplitting * > | splits |
vector< int > | bornColors |
Public Attributes inherited from TimeShower | |
MergingHooksPtr | mergingHooksPtr {} |
Pointer to MergingHooks object for NLO merging. | |
WeightContainer * | weightContainerPtr {} |
Protected Attributes | |
int | iSysSel |
Store properties to be returned by methods. | |
double | pTmaxFudge |
double | pTLastBranch |
Protected Attributes inherited from TimeShower | |
int | beamOffset {} |
Beam location offset in event. | |
PartonVertexPtr | partonVertexPtr {} |
Pointer to assign space-time vertices during parton evolution. | |
bool | doUncertainties {} |
Store uncertainty variations relevant to TimeShower. | |
bool | uVarMuSoftCorr {} |
bool | uVarMPIshowers {} |
bool | noResVariations {} |
bool | noProcVariations {} |
int | nUncertaintyVariations {} |
int | nVarQCD {} |
int | uVarNflavQ {} |
double | dASmax {} |
double | cNSpTmin {} |
double | uVarpTmin2 {} |
double | overFactor {} |
double | overFactorEnhance {} |
map< int, double > | varG2GGmuRfac |
map< int, double > | varQ2QGmuRfac |
map< int, double > | varG2QQmuRfac |
map< int, double > | varX2XGmuRfac |
map< int, double > | varG2GGcNS |
map< int, double > | varQ2QGcNS |
map< int, double > | varG2QQcNS |
map< int, double > | varX2XGcNS |
map< int, double > * | varPDFplus |
map< int, double > * | varPDFminus |
map< int, double > * | varPDFmember |
unordered_map< string, double > | enhanceFSR |
Protected Attributes inherited from PhysicsBase | |
Info * | infoPtr = {} |
Settings * | settingsPtr = {} |
Pointer to the settings database. | |
ParticleData * | particleDataPtr = {} |
Pointer to the particle data table. | |
Logger * | loggerPtr = {} |
Pointer to logger. | |
HadronWidths * | hadronWidthsPtr = {} |
Pointer to the hadron widths data table. | |
Rndm * | rndmPtr = {} |
Pointer to the random number generator. | |
CoupSM * | coupSMPtr = {} |
Pointers to SM and SUSY couplings. | |
CoupSUSY * | coupSUSYPtr = {} |
BeamSetup * | beamSetupPtr = {} |
BeamParticle * | beamAPtr = {} |
BeamParticle * | beamBPtr = {} |
BeamParticle * | beamPomAPtr = {} |
BeamParticle * | beamPomBPtr = {} |
BeamParticle * | beamGamAPtr = {} |
BeamParticle * | beamGamBPtr = {} |
BeamParticle * | beamVMDAPtr = {} |
BeamParticle * | beamVMDBPtr = {} |
PartonSystems * | partonSystemsPtr = {} |
Pointer to information on subcollision parton locations. | |
SigmaTotal * | sigmaTotPtr = {} |
Pointers to the total/elastic/diffractive cross sections. | |
SigmaCombined * | sigmaCmbPtr = {} |
set< PhysicsBase * > | subObjects |
UserHooksPtr | userHooksPtr |
Friends | |
class | DireSplitting |
class | DireSpace |
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. | |
The DireTimes class does timelike showers.
|
virtual |
Only consider final-state emissions.
Gluon emission is allowed.
Gluon emission is allowed.
Gluon branching to quarks is allowed.
Photon emission. Photon emission from quarks.
Photon emission from charged leptons.
Z-boson emission. Z-boson emission from quarks.
Z-boson emission from charged leptons.
Photon splitting. Photon branching to quarks is allowed.
Photon branching to leptons is allowed.
W-boson splitting. W-boson branching to quarks is allowed.
H-boson splitting. H-boson branching to photons is allowed.
Reimplemented from TimeShower.
double alphaemNow | ( | double | pT2, |
double | renormMultFacNow = 1. , |
||
int | iSys = 0 |
||
) |
Function to calculate the correct alphaEM/2*Pi value.
Function to calculate the correct alphaem/2*Pi value, including renormalisation scale variations + threshold matching.
Get alphaEM(k*pT^2) and subtractions.
Done.
double alphasNow | ( | double | pT2, |
double | renormMultFacNow = 1. , |
||
int | iSys = 0 |
||
) |
Function to calculate the correct alphaS/2*Pi value, including renormalisation scale variations + threshold matching.
Get beam for PDF alphaS, if necessary.
Get alphaS(k*pT^2) and subtractions.
Get kernel order.
Use simple kernels for showering secondary scatterings.
Now find the necessary thresholds so that alphaS can be matched correctly.
Done.
|
virtual |
Setup branching kinematics.
This function is a wrapper for setting up the branching kinematics.
Done.
Reimplemented from TimeShower.
bool branch_FF | ( | Event & | event, |
bool | trial = false , |
||
DireSplitInfo * | split = nullptr |
||
) |
ME corrections and kinematics that may give failure. Notation: radBef, recBef = radiator, recoiler before emission, rad, rec, emt = radiator, recoiler, emitted efter emission. (rad, emt distinguished by colour flow for g -> q qbar.)
Check if the first emission should be studied for removal.
Find initial radiator and recoiler particles in dipole branching.
Find their momenta, with special sum for global recoil.
Get splitting variables.
Allow splitting kernel to overwrite phase space variables.
Get flavour of splitting.
Name of the splitting.
Default flavours and colour tags for new particles in dipole branching.
New colour tag required for gluon emission.
New flavours for g -> q qbar; split colours.
Now reset if splitting information is available.
Get particle masses. Radiator before splitting.
Radiator after splitting.
Recoiler.
Emission.
Adjust the dipole kinematical mass to accomodate masses after branching.
Q2 already correctly stored in splitInfo.
Calculate CS variables.
Allow splitting kernel to overwrite phase space variables.
Auxiliary angle.
Second angle for 1->3 splitting.
Allow splitting kernel to overwrite phase space variables.
Get dipole 4-momentum.
1->3 splittings generated in CS variables directly.
Calculate derived variables.
Not possible to construct kinematics if kT2 < 0.0
NaN kT2 can happen for a 1->3 splitting in which the g->QQ~ produces massive quarks Q.
Now construct the new recoiler momentum in the lab frame.
Construct left-over dipole momentum by momentum conservation.
Set up transverse momentum vector by using two perpendicular four-vectors.
Construct new radiator momentum.
Contruct the emission momentum by momentum conservation.
Ensure that radiator is on mass-shell
Ensure that emission is on mass-shell
Ensure that recoiler is on mass-shell
Swap emitted and radiator properties for first part of 1->3 splitting (q -> "massive gluon" + q)
For emitted color singlet, redefine the colors of the "massive gluon".
Define new particles from dipole branching.
Exempt off-shell radiator from Pythia momentum checks.
Default to stored color info for intermediate step in 1->3 branching.
Special checks to set weak particles status equal to 56. This is needed for decaying the particles. Also set polarisation.
Save properties to be restored in case of user-hook veto of emission.
Shower may occur at a displaced vertex.
Put new particles into the event record. Mark original dipole partons as branched and set daughters/mothers.
Store flavour again, in case dipSel gets removed.
Check user veto for 1->2 branchings.
Check momenta.
Apply ME correction if necessary.
Try to find incoming particle in other systems, i.e. if the current system arose from a resonance decay.
&& pT2 > pT2minMECs && checkSIJ(event,1.)) {
Finally update the list of all partons in all systems.
Update dipoles and beams.
Heavy particle 1->2 decay for "second step" in 1->3 splitting.
Check momenta.
Swap emitted and radiator indices.
Flavours already fixed by 1->3 kernel.
Colour tags for new particles in branching.
Already correctly read id and colors from SplitInfo object.
Get particle masses.
Radiator after splitting.
Recoiler.
Emission.
Construct FF dipole momentum.
Calculate CS variables.
Calculate derived variables.
Construct left-over dipole momentum by momentum conservation.
Set up transverse momentum vector by using two perpendicular 4-vectors.
Construct new radiator momentum.
Contruct the emission momentum by momentum conservation.
Recoiler unchanged.
Define new particles from dipole branching.
Check momenta.
Check invariants.
Update bookkeeping
Update dipoles and beams.
Shower may occur at a displaced vertex.
Put new particles into the event record. Mark original dipole partons as branched and set daughters/mothers.
Update dipoles and beams.
Allow veto of branching. If so restore event record to before emission.
This case is identical to the case where the probability to accept the emission was indeed zero all along. In this case, neither acceptProbability nor rejectProbability would have been filled. Thus, remove the relevant entries from the weight container!
Clear accept/reject weights.
Store positions of new particles.
Set shower weight.
Store positions of new soft particles.
For 1->3 splitting, do not tag emissions as "soft".
Clear accept/reject weights.
Done.
bool branch_FI | ( | Event & | event, |
bool | trial = false , |
||
DireSplitInfo * | split = nullptr |
||
) |
ME corrections and kinematics that may give failure. Notation: radBef, recBef = radiator, recoiler before emission, rad, rec, emt = radiator, recoiler, emitted efter emission. (rad, emt distinguished by colour flow for g -> q qbar.)
Check if the first emission should be studied for removal.
Find initial radiator and recoiler particles in dipole branching.
Find their momenta, with special sum for global recoil.
Get splitting variables.
Allow splitting kernel to overwrite phase space variables.
Calculate CS variables.
Get flavour of splitting.
Store flavour again, in case dipSel gets removed or flavour gets reset.
Name of the splitting.
Default flavours and colour tags for new particles in dipole branching.
New colour tag required for gluon emission.
New flavours for g -> q qbar; split colours.
Now reset if splitting information is available.
Get particle masses. Radiator before splitting.
Radiator after splitting.
Emission.
Initial state recoiler always assumed massless.
Recoiler mass
Recalculate the kinematicaly available dipole mass.
Calculate CS variables.
Allow splitting kernel to overwrite phase space variables.
Auxiliary angle.
Second angle for 1->3 splitting.
Allow splitting kernel to overwrite phase space variables.
Get dipole 4-momentum.
1->3 splittings generated in CS variables directly.
Calculate derived variables.
Not possible to construct kinematics if kT2 < 0.0
NaN kT2 can happen for a 1->3 splitting in which the g->QQ~ produces massive quarks Q.
Now construct the new recoiler momentum in the lab frame.
Construct left-over dipole momentum by momentum conservation.
Set up transverse momentum vector by using two perpendicular four-vectors.
Construct new radiator momentum.
Contruct the emission momentum by momentum conservation.
Construct FF dipole momentum.
Calculate derived variables.
Ensure that radiator is on mass-shell
Ensure that emission is on mass-shell
Ensure that recoiler is on mass-shell
New: Return if the x-value for the incoming recoiler is nonsense.
Swap emitted and radiator properties for first part of 1->3 splitting (q -> "massive gluon" + q)
For emitted color singlet, redefine the colors of the "massive gluon".
Define new particles from dipole branching.
Exempt off-shell radiator from Pythia momentum checks.
Special checks to set weak particles status equal to 56. This is needed for decaying the particles. Also set polarisation.
Default to stored color info for intermediate step in 1->3 branching.
Save properties to be restored in case of user-hook veto of emission.
Shower may occur at a displaced vertex.
Put new particles into the event record. Mark original dipole partons as branched and set daughters/mothers.
Check momenta.
Check that beam still has leftover momentum.
Restore old beams.
Apply ME correction if necessary.
&& pT2 > pT2minMECs && checkSIJ(event,1.)) {
Temporarily update parton systems.
Undo update of parton systems.
Just update dipoles and beams.
Heavy particle 1->2 decay for "second step" in 1->3 splitting.
Check momenta.
Swap emitted and radiator indices.
Flavours already fixed by 1->3 kernel.
Colour tags for new particles in branching.
Already correctly read id and colors from SplitInfo object.
Get particle masses.
Construct FF dipole momentum.
Calculate CS variables.
Calculate derived variables.
Not possible to construct kinematics if kT2 < 0.0
NaN kT2 can happen for a 1->3 splitting in which the g->QQ~ produces massive quarks Q.
Update dipoles and beams.
Check that beam still has leftover momentum.
Restore old beams.
Boost the transverse momentum vector into the lab frame. Construct left-over dipole momentum by momentum conservation.
Set up transverse momentum vector by using two perpendicular four-vectors.
Construct new radiator momentum.
Contruct the emission momentum by momentum conservation.
Recoiler unchanged.
Check invariants.
Check momenta.
Check that beam still has leftover momentum.
Restore old beams.
Update bookkeeping
Define new particles from dipole branching.
colRad, acolRad, pRad, sqrt(m2r), pTsel);
colEmt, acolEmt, pEmt, sqrt(m2e), pTsel);
Shower may occur at a displaced vertex.
Put new particles into the event record. Mark original dipole partons as branched and set daughters/mothers.
Update dipoles and beams.
Temporarily set the daughters in the beams to zero, to allow mother-daughter relation checks.
Check if mother-daughter relations are correctly set. Check only possible if no MPI are present.
Restore correct daughters in the beams.
Allow veto of branching. If so restore event record to before emission.
This case is identical to the case where the probability to accept the emission was indeed zero all along. In this case, neither acceptProbability nor rejectProbability would have been filled. Thus, remove the relevant entries from the weight container!
Clear accept/reject weights.
Store positions of new particles.
Set shower weight.
Store positions of new soft particles.
Clear accept/reject weights.
Done.
void clear | ( | ) |
Reset parton shower.
Clear weighted shower book-keeping.
bool cluster_FF | ( | const Event & | state, |
int | iRad, | ||
int | iEmt, | ||
int | iRec, | ||
int | idRadBef, | ||
Particle & | radBef, | ||
Particle & | recBef | ||
) |
Calculate CS variables.
Get particle masses.
Set resonance mass to virtuality.
Get dipole 4-momentum.
Upate type if this is a massive splitting.
Check phase space constraints.
Calculate derived variables.
Now construct the new recoiler momentum in the lab frame.
Get momentum of radiator my momentum conservation.
Done
bool cluster_FI | ( | const Event & | state, |
int | iRad, | ||
int | iEmt, | ||
int | iRec, | ||
int | idRadBef, | ||
Particle & | radBef, | ||
Particle & | recBef | ||
) |
Calculate CS variables.
Get particle masses.
Set resonance mass to virtuality.
Get dipole 4-momentum.
Get momentum of radiator my momentum conservation.
Ensure that radiator is on mass-shell
Ensure that recoiler is on mass-shell
Done
Check phase space constraints.
Calculate CS variables.
Upate type if this is a massive splitting.
Calculate derived variables.
Now construct the new recoiler momentum in the lab frame.
Get momentum of radiator my momentum conservation.
Ensure that radiator is on mass-shell
Ensure that recoiler is on mass-shell
Done
pair< Event, pair< int, int > > clustered_internal | ( | const Event & | state, |
int | iRad, | ||
int | iEmt, | ||
int | iRec, | ||
string | name | ||
) |
Flags for type of radiation
Construct the clustered event
Copy all unchanged particles to NewEvent
Copy all the junctions one by one
Find an appropriate scale for the hard process
Initialise scales for new event
Set properties of radiator/recoiler after the clustering Recoiler properties will be unchanged
Find flavour of radiator before splitting
Put dummy values for colours
Reset status if the reclustered radiator is a resonance.
Put mass for radiator and recoiler
Construct momenta and colours of clustered particles.
Put some dummy production scales for RecBefore, RadBefore
Append new recoiler and find new radiator colour
Assign the correct colour to re-clustered radiator.
Build the clustered event
Clustering might not be possible due to phase space constraints.
Copy system and incoming beam particles to outState
Copy all the junctions one by one
Initialise scales for new event
Save position of radiator in new event record
Append first incoming particle
Find second incoming in input event
Append second incoming particle
Find second incoming in input event
Append new recoiler if not done already
Append new radiator if not done already
Append intermediate particle (careful not to append reclustered recoiler)
Append final state particles, resonances first
Then start appending partons
Then partons (not reclustered recoiler)
Then the rest
Find intermediate and respective daughters
Find daughters in output state
If both daughters found, done Else put first final particle as first daughter and last final particle as second daughter
Set daughters and mothers
Force outgoing particles to be part of hard process.
Find range of final state partons
Update event properties
Almost there... If an intermediate coloured parton exists which was directly colour connected to the radiator before the splitting, and the radiator before and after the splitting had only one colour, problems will arise since the colour of the radiator will be changed, whereas the intermediate parton still has the old colour. In effect, this means that when setting up a event for trial showering, one colour will be free. Hence, check for an intermediate coloured triplet resonance has been colour-connected to the "old" radiator. Find resonance
Find resonance connected to initial colour
Find resonance connected to initial anticolour
Find resonance connected to final state colour
Find resonance connected to final state anticolour
Now find this resonance in the reclustered state
Find reclustered radiator colours
Find resonance radiator colours
Check if any of the reclustered radiators colours match the resonance
If a resonance has been found, but no colours match, change the colour of the resonance
If a resonance has been found, but no colours match, and the position of the resonance in the event record has been changed, update the radiator mother
Force outgoing particles to be part of hard process.
Force incoming to be part of the hard process.
Now check event.
Check if the state is valid. If not, return empty state.
Done
pair< Vec4, Vec4 > decayWithOffshellRec | ( | double | zCS, |
double | yCS, | ||
double | phi, | ||
double | m2RadBef, | ||
double | m2RadAft, | ||
double | m2EmtAft, | ||
Vec4 | pRadBef, | ||
Vec4 | pRecBef | ||
) |
Calculate derived variables.
Not possible to construct kinematics if kT2 < 0.0
Construct left-over dipole momentum by momentum conservation.
Set up transverse momentum vector by using two perpendicular four-vectors.
Construct new radiator momentum.
Contruct the emission momentum by momentum conservation.
Store and return.
pair< Vec4, Vec4 > decayWithOnshellRec | ( | double | zCS, |
double | yCS, | ||
double | phi, | ||
double | m2Rec, | ||
double | m2RadAft, | ||
double | m2EmtAft, | ||
Vec4 | pRadBef, | ||
Vec4 | pRecBef | ||
) |
Construct FF dipole momentum.
Calculate derived variables.
Construct left-over dipole momentum by momentum conservation.
Set up transverse momentum vector by using two perpendicular 4-vectors.
Construct new radiator momentum.
Contruct the emission momentum by momentum conservation.
Recoiler unchanged.
Store and return.
void finalize | ( | Event & | event | ) |
Finalize event after evolution.
Perform resonance decays if some shower-produced resonances are still left-over after evolution terminates.
Construct the decay products.
Attach decay products to event.
Add new system, automatically with two empty beam slots.
Loop over allowed range to find all final-state particles.
Clear resonances (even if not decayed)
Currently, do nothing to not break other things.
int FindCol | ( | int | col, |
vector< int > | iExclude, | ||
const Event & | event, | ||
int | type, | ||
int | iSys = -1 |
||
) |
Find index of colour partner for input colour.
Unset if the incoming particles are flagged as outgoing. Instead, try to resort to information stored in 0th event entry.
Search event record for matching colour & anticolour
Skip if this index is excluded.
Search event record for matching colour & anticolour
Skip if this index is excluded.
Check incoming
if no matching colour / anticolour has been found, return false
Auxiliary function to return the position of a particle. Should go int Event class eventually!
Function to in the input event find a particle with quantum numbers matching those of the input particle IN Particle : Particle to be searched for Event : Event to be searched in OUT int : > 0 : Position of matching particle in event < 0 : No match in event
|
virtual |
List of recoilers.
Reimplemented from TimeShower.
|
inlinevirtual |
From Pythia version 8.215 onwards. Return a string identifier of a splitting. Usage: getSplittingName( const Event& event, int iRad, int iEmt, int iRec)
Reimplemented from TimeShower.
|
virtual |
Compute splitting probability.
From Pythia version 8.215 onwards. Return the splitting probability. Usage: getSplittingProb( const Event& event, int iRad, int iEmt, int iRec)
From Pythia version 8.215 onwards.
Get kernel order.
Do nothing if kernel says so, e.g. to avoid infinite loops if the kernel uses the History class.
Disallow below cut-off.
Upate type if this is a massive splitting.
Recalculate the kinematicaly available dipole mass.
Check phase space constraints.
Disallow below cut-off. if ( pT2cut(state[iEmt].id()) > pT2) return 0.;
Calculate splitting probability.
Get phi angle.
Get splitting probability.
Get complete kernel.
Reset again.
Multiply with 1/pT^2. Note: No additional Jacobian factors, since for our choice of variables, we always have Jacobian_{mass to CS} * Jacobian_{CS to DIRE} * Propagator = 1/pT2
Make splitting probability positive if ME corrections are available.
Note: The additional factor 1/xCS for rescaling the initial flux is NOT included, so that we can apply PDF ratios [x1 f(x1)] / [x0 f(x0) ] later.
Reimplemented from TimeShower.
|
virtual |
From Pythia version 8.218 onwards. Return the evolution variable. Usage: getStateVariables( const Event& event, int iRad, int iEmt, int iRec, string name) Important note:
From Pythia version 8.218 onwards. Return the evolution variable and splitting information. More comments in the header.
Kinematical variables.
Book-keeping for particle before emission.
Variables defining the PS starting scales.
In this case, insert only dummy information except for PDF scale.
Find the shower starting scale.
Loop through final state of system to find possible dipole ends.
Find dipole end formed by colour index.
Find dipole end formed by anticolour index.
Now find non-QCD dipoles and/or update the existing dipoles.
Get x for both beams.
Store invariant masses of all dipole ends.
Reimplemented from TimeShower.
|
virtual |
Initialize alphaStrong and related pTmin parameters.
Initialize alphaStrong, alphaEM and related pTmin parameters.
Colour factors.
Alternatively only initialize resonance decays.
Store input pointers for future use.
Main flags.
Matching in pT of hard interaction or MPI to shower evolution.
Charm and bottom mass thresholds.
Parameters of scale choices (inherited from Pythia).
Parameters of alphaStrong generation.
Set flavour thresholds by default Pythia masses, unless zero.
Initialize alphaStrong generation.
Lambda for 5, 4 and 3 flavours.
Parameters of QCD evolution. Warn if pTmin must be raised.
Parameters of Pythia QED evolution.
Parameters of alphaEM generation.
Initialize alphaEM generation.
Parameters of QED evolution, sums of charges, as necessary to pick flavor.
Allow massive incoming particles. Currently not supported by Pythia. useMassiveBeams = settingsPtr->flag("Beams:massiveLeptonBeams");
Z0 and W+- properties needed for gamma/Z0 mixing and weak showers.
Mode for higher-order kernels.
Create maps of accept/reject weights
Number of MPI, in case MPI forces intervention in shower weights.
Set splitting library, if already exists.
May have to fix up recoils related to rescattering.
Possibility of two predetermined hard emissions in event.
Possibility to allow user veto of emission step.
Set initial value, just in case.
Done.
Reimplemented from TimeShower.
void initVariations | ( | ) |
Initialize bookkeeping of shower variations.
Create maps of accept/reject weights
Done.
|
inlinevirtual |
From Pythia version 8.215 onwards. Check if attempted clustering is handled by timelike shower Usage: isTimelike( const Event& event, int iRad, int iEmt, int iRec, string name)
Reimplemented from TimeShower.
|
virtual |
Find whether to limit maximum scale of emissions, and whether to dampen.
Find whether to limit maximum scale of emissions. Also allow for dampening at factorization or renormalization scale.
Find whether to limit pT. Begin by user-set cases.
Always restrict SoftQCD processes.
Look if any quark (u, d, s, c, b), gluon or photon in final state. Also count number of heavy coloured particles, like top.
Dampening at factorization or renormalization scale; only for hardest.
Done.
Reimplemented from TimeShower.
|
virtual |
Print dipole list; for debug mainly.
Print the list of dipoles.
Header.
Loop over dipole list and print it.
Done.
Reimplemented from TimeShower.
Try to find incoming particle in other systems, i.e. if the current system arose from a resonance decay.
Attach the first incoming particle.
Append one time (as beams)
Append second time as incoming particles.
Careful when builing the sub-events: A particle that is currently intermediate in one system could be the progenitor of another system, i.e. when resonance decays are present. In this case, the intermediate particle in the current system should be final.
Set daughters of initial particles.
double newPoint | ( | const Event & | event | ) |
Begin loop over all possible radiating dipole ends.
Dipole properties.
Find maximum mass scale for dipole.
Reset emission properties.
Reset properties of 1->3 splittings.
Store start/end scales.
Do evolution if it makes sense.
Unable to produce splitting.
Return nonvanishing value if found pT bigger than already found.
|
virtual |
Prepare system for evolution after each new interaction; identify ME.
Prepare system for evolution; identify ME.
Calculate remainder shower weight after MPI.
Clear accept/reject weights.
Store number of MPI, in case a subsequent MPI forces calculation and reset of shower weights.
Reset dipole-ends list for first interaction and for resonance decays.
Set splitting library.
No dipoles for 2 -> 1 processes.
Loop through final state of system to find possible decays.
Setup decay dipoles.
In case of DPS overwrite limitPTmaxIn by saved value.
Loop through final state of system to find possible dipole ends.
Find dipole end formed by colour index.
Find dipole end formed by anticolour index.
Now find non-QCD dipoles and/or update the existing dipoles.
End loop over system final state. Have now found the dipole ends.
Setup decay dipoles.
Loop through dipole ends and set matrix element correction flag.
Now update masses and allowed emissions. (Not necessary here, since ensured in setupQCDdip etc.)
Update dipole list after a multiparton interactions rescattering.
Counter of proposed emissions.
Clear weighted shower book-keeping.
Reimplemented from TimeShower.
|
virtual |
Global recoil: reset counters and store locations of outgoing partons.
Global recoil not used, but abuse function to reset some generic things.
Initialize weight container.
Clear event-by-event diagnostic messages.
Clear accept/reject weights.
Now also attempt to reset ISR weights!
Reimplemented from TimeShower.
return sij*saj / (sai+saj) * (sij+saj+sai) / (sai+saj) ;
|
virtual |
Select next pT in downwards evolution. Wrapper function inherited from Pythia.
Select next pT in downwards evolution of the existing dipoles. Classical Sudakov veto algorithm to produce next set of state variables.
Begin loop over all possible radiating dipole ends.
Remember if this is a trial emission.
Limit final state multiplicity. For debugging only
Dipole properties.
Find maximum mass scale for dipole.
Reset emission properties.
Reset properties of 1->3 splittings.
Do not try splitting if the corrected dipole mass is negative.
Evolution.
Store start/end scales.
Do evolution if it makes sense.
Update if found pT larger than current maximum.
Update the number of proposed timelike emissions.
Insert additional weights.
Timelike shower of resonances: Finalize in case resonances where produced.
Return nonvanishing value if found pT bigger than already found.
Reimplemented from TimeShower.
|
inline |
Initialize various pointers. (Separated from rest of init since not virtual.)
|
virtual |
Update dipole list after a multiparton interactions rescattering.
Loop over two incoming partons in system; find their rescattering mother. (iOut is outgoing from old system = incoming iIn of rescattering system.)
Loop over all dipoles.
Kill dipoles where rescattered parton is radiator.
No matrix element for dipoles between scatterings.
Update dipoles where outgoing rescattered parton is recoiler.
Colour dipole: recoil in final state, initial state or new.
This line in case mother is a rescattered parton.
If above options failed, then create new dipole.
Anticolour dipole: recoil in final state, initial state or new.
This line in case mother is a rescattered parton.
If above options failed, then create new dipole.
End of loop over dipoles and two incoming sides.
Reimplemented from TimeShower.
|
virtual |
Top-level routine to do a full time-like shower in resonance decay.
Add new system, automatically with two empty beam slots.
Loop over allowed range to find all final-state particles.
Let prepare routine do the setup.
Begin evolution down in pT from hard pT scale.
Do a final-state emission (if allowed).
Keep on evolving until nothing is left to be done.
Return number of emissions that were performed.
Reimplemented from TimeShower.
|
virtual |
Top-level routine for QED radiation in hadronic decay to two leptons.
Top-level routine for QED radiation in hadronic decay to two leptons. Intentionally only does photon radiation, i.e. no photon branchings.
Add new system, automatically with two empty beam slots.
Reset scales to allow setting up dipoles.
Prepare all dipoles for evolution.
Begin evolution down in pT from hard pT scale.
Do a final-state emission (if allowed).
Done with branching
Keep on evolving until nothing is left to be done.
Undo scale setting.
Return number of emissions that were performed.
Reimplemented from TimeShower.
|
virtual |
Update dipole list after each ISR emission.
Update dipole list after each ISR emission (so not used for resonances).
Store dipoles in other systems.
Reset dipole-ends.
No dipoles for 2 -> 1 processes.
Loop through final state of system to find possible dipole ends.
Find dipole end formed by colour index.
Find dipole end formed by anticolour index.
Now find non-QCD dipoles and/or update the existing dipoles.
End loop over system final state. Have now found the dipole ends.
Setup decay dipoles.
Now update masses and allowed emissions.
Reimplemented from TimeShower.
void updateAfterFF | ( | int | iSysSelNow, |
int | iSysSelRec, | ||
Event & | event, | ||
int | iRadBef, | ||
int | iRecBef, | ||
int | iRad, | ||
int | iEmt, | ||
int | iRec, | ||
int | flavour, | ||
int | colType, | ||
double | pTsel | ||
) |
Update dipole list after final-final splitting.
Gluon emission: update both dipole ends and add two new ones.
Strive to match colour to anticolour inside closed system.
When recoiler was uncoloured particle, in resonance decays, assign recoil to coloured particle.
Set dipole mass properties.
Gluon branching to q qbar: update current dipole and other of gluon.
Update dipoles for second step in 1->3 splitting.
Nothing to be done if dipole end has already been updated.
Strive to match colour to anticolour inside closed system.
Potentially also update recoiler!
Nothing to be done if dipole end has already been updated.
Always update the production pT.
Just update old radiator/recoiler to current outgoing particles.
Update radiator-recoiler end.
Update recoiler-radiator end.
Now update other dipoles that also involved the radiator or recoiler.
Nothing to be done if dipole end has already been updated.
Now update or construct new dipoles if the radiator or emission allow for new types of emissions.
Now check if a new dipole end a-b should be added: First check if the dipole end is already existing.
If the dipole end exists, attempt to update the dipole end (a) for the current a-b dipole.
If no dipole exists and idEmtAfter != 0, create new dipole end (a).
Copy or set lifetime for new final state.
Finally update the list of all partons in all systems.
Loop through final state of system to find possible new dipole ends.
Find dipole end formed by colour index.
Find dipole end formed by anticolour index.
Now find non-QCD dipoles and/or update the existing dipoles.
Setup decay dipoles.
Now remove inactive dipoles.
Now update all dipoles.
Bookkeep shower-induced resonances.
else direInfoPtr->addResPos(iRad);
Done.
void updateAfterFI | ( | int | iSysSelNow, |
int | iSysSelRec, | ||
Event & | event, | ||
int | iRadBef, | ||
int | iRecBef, | ||
int | iRad, | ||
int | iEmt, | ||
int | iRec, | ||
int | flavour, | ||
int | colType, | ||
double | pTsel, | ||
double | xNew | ||
) |
Update dipole list after final-final splitting.
For initial-state recoiler also update beam and sHat info.
Strive to match colour to anticolour inside closed system.
When recoiler was uncoloured particle, in resonance decays, assign recoil to coloured particle.
Set dipole mass properties.
Set dipole mass properties.
Gluon branching to q qbar: update current dipole and other of gluon.
Update dipoles for second step in 1->3 splitting.
Nothing to be done if dipole end has already been updated.
Strive to match colour to anticolour inside closed system.
Nothing to be done if dipole end has already been updated.
Always update the production pT.
Just update old radiator/recoiler to current outgoing particles.
Update radiator-recoiler end.
Update recoiler-radiator end.
Now update other dipoles that also involved the radiator or recoiler. Note: For 1->3 splittings, this step has already been done earlier!
Nothing to be done if dipole end has already been updated.
Now update or construct new dipoles if the radiator or emission allow for new types of emissions. Careful: do not include initial recoiler as potential radiator.
Now check if a new dipole end a-b should be added: First check if the dipole end is already existing.
If the dipole end exists, attempt to update the dipole end (a) for the current a-b dipole.
If no dipole exists and idEmtAfter != 0, create new dipole end (a).
Copy or set lifetime for new final state.
Finally update the list of all partons in all systems.
Loop through final state of system to find possible new dipole ends.
Find dipole end formed by colour index.
Find dipole end formed by anticolour index.
Now find non-QCD dipoles and/or update the existing dipoles.
Setup decay dipoles.
Now update all dipoles.
Bookkeep shower-induced resonances.
Done.
bool validEvent | ( | const Event & | state, |
bool | isProcess = false , |
||
int | iSysCheck = -1 |
||
) |
Check colour/flavour correctness of state.
Check for NaNs or INFs.
Done if the state is already broken.
Check if event is coloured
Check colour of quarks
No corresponding anticolour in final state
No corresponding colour in initial state
Check anticolour of antiquarks
No corresponding colour in final state
No corresponding anticolour in initial state
No uncontracted colour (anticolour) charge of gluons
No corresponding anticolour in final state
No corresponding colour in initial state
No corresponding colour in final state
No corresponding anticolour in initial state
Check charge sum in initial and final state
Check if particles are on mass shell
Check that overall pT is vanishing.
if ( i ==3 || i == 4 ) pSum -= event[i].p();
Check for negative energies.
Done with loop over systems.
bool validMomentum | ( | const Vec4 & | p, |
int | id, | ||
int | status | ||
) |
Check that particle has sensible momentum.
Check colour/flavour correctness of state.
Check for NaNs or INFs.
Check if particles is on mass shell
Do not check on-shell condition for massive intermediate (!) resonances. Assuming all non-SM particles are heavy here!
Check for negative energies.
Done
bool validMotherDaughter | ( | const Event & | state | ) |
Check that mother-daughter-relations are correctly set.
Loop through the event and check that there are beam particles.
Check that mother and daughter lists not empty where not expected to.
Check that the particle appears in the daughters list of each mother.
Check that the particle appears in the mothers list of each daughter.
Mother-daughter relations not correct if any lists do not match.
Done.