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

Gets cross sections for hadron-hadron collisions at low energies. More...

#include <SigmaLowEnergy.h>

Inheritance diagram for SigmaLowEnergy:
PhysicsBase

Public Member Functions

void init (NucleonExcitations *nucleonExcitationsPtrIn)
 Initialize. More...
 
double sigmaTotal (int idA, int idB, double eCM, double mA, double mB)
 Get the total cross section for the specified collision. More...
 
double sigmaTotal (int idAIn, int idBIn, double eCMIn)
 
double sigmaPartial (int idA, int idB, double eCM, double mA, double mB, int proc)
 Gets the partial cross section for the specified process. More...
 
double sigmaPartial (int idAIn, int idBIn, double eCMIn, int proc)
 
bool sigmaPartial (int idA, int idB, double eCM, double mA, double mB, vector< int > &procsOut, vector< double > &sigmasOut)
 
int pickProcess (int idA, int idB, double eCM, double mA, double mB)
 Picks a process randomly according to their partial cross sections.
 
int pickResonance (int idA, int idB, double eCM)
 Picks a resonance according to their partial cross sections. More...
 
bool hasExcitation (int idAIn, int idBIn) const
 
void updateResonances ()
 Update the list of internal resonances. More...
 
- 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
 

Static Public Attributes

static constexpr double TINYSIGMA = 1.e-9
 Cross sections below threshold are assumed numerically equal to zero.
 

Friends

class LowEnergyProcess
 LowEnergyProcess should have access to nqEffAQM for slope in t.
 

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

Gets cross sections for hadron-hadron collisions at low energies.

Member Function Documentation

bool hasExcitation ( int  idAIn,
int  idBIn 
) const
inline

For NN / N Nbar / Nbar Nbar collisions explicit excitation states. replace a generic smeared-out low-mass diffraction component.

void init ( NucleonExcitations nucleonExcitationsPtrIn)

Initialize.

The SigmaLowEnergy class.

Initialize.

Flag to allow or suppress inelastic processes.

Mode for calculating total cross sections for pi pi and pi K.

Suppression factors in Additive Quark Model (AQM).

Mixing for eta and eta'; specifically fraction of s sbar.

Some standard masses.

Store pointer.

Fill explicit resonances.

int pickResonance ( int  idA,
int  idB,
double  eCM 
)

Picks a resonance according to their partial cross sections.

Set canonical ordering.

Fail if no resonances exist.

Calculate cross section for each resonance.

Return zero if no resonances are available.

Pick resonance at random.

double sigmaPartial ( int  idA,
int  idB,
double  eCM,
double  mA,
double  mB,
int  proc 
)

Gets the partial cross section for the specified process.

Get the partial cross section for the specified collision and process. proc | 0: total; | 1: nondiff; | 2 : el; | 3: SD (XB); | 4: SD (AX); | 5: DD; | 6: CD (AXB, not implemented) | 7: excitation | 8: annihilation | 9: resonant | >100: resonant through the specified resonance particle

Energy cannot be less than the hadron masses.

For K0S/K0L, take average of K0 and K0bar.

Total cross section.

Get all partial cross sections.

Return partial cross section.

bool sigmaPartial ( int  idAIn,
int  idBIn,
double  eCMIn,
double  mAIn,
double  mBIn,
vector< int > &  procsOut,
vector< double > &  sigmasOut 
)

Gets all partial cross sections for the specified collision. This is used when all cross sections are needed to determine which process to execute. Returns false if no processes are available.

Gets all partial cross sections for the specified collision. Returns whether any processes have positive cross sections.

No cross sections below threshold.

If K_S/K_L + X then take the average K0 and K0bar.

If X + K_S/K_L then take the average K0 and K0bar.

Store current configuration.

Calculate total, annihilation and resonant cross sections (if needed).

Abort early if total cross section vanishes.

If inelastic is off, return only elastic cros ssection.

Calculate diffractive cross sections.

Calculate elastic cross sections.

Calculate excitation cross sections for NN.

Calculate nondiffractive as difference between total and the others.

Give warning if sigmaND is very negative.

Rescale pi pi and pi K cross sections to match PelĂ ez total cross section.

Rescale for pi pi below 1.42 GeV, and for pi K below 1.8 GeV

Get cross section from parameterization.

Write results to output arrays.

double sigmaTotal ( int  idA,
int  idB,
double  eCM,
double  mA,
double  mB 
)

Get the total cross section for the specified collision.

Energy cannot be less than the hadron masses.

For K0S/K0L, take average of K0 and K0bar.

Fix particle ordering.

Get custom cross section if applicable.

Special handling for pi pi and pi K cross sections

Calculate total.

void updateResonances ( )

Update the list of internal resonances.

Use setConfig to ensure canonical ordering.


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