PYTHIA  8.312
Public Member Functions | List of all members
PythiaCascade Class Reference

Intended flow: More...

#include <PythiaCascade.h>

Public Member Functions

 PythiaCascade ()=default
 Default constructor, all setup is done in init().
 
bool init (double eMaxIn=1e9, bool listFinalIn=false, bool rapidDecaysIn=false, double smallTau0In=1e-10, bool reuseMPI=true, string initFile="pythiaCascade.mpi")
 Initialize PythiaCascade for a given maximal incoming energy. More...
 
double nCollAvg (int A)
 
bool sigmaSetuphN (int idNowIn, Vec4 pNowIn, double mNowIn)
 
double sigmahA (int A)
 
EventnextColl (int Znow, int Anow, Vec4 vNow=Vec4())
 
EventnextDecay (int idNowIn, Vec4 pNowIn, double mNowIn, Vec4 vNow=Vec4())
 
void compress ()
 
void stat ()
 Summary of aborts, errors and warnings.
 
ParticleDataparticleData ()
 
Rndmrndm ()
 

Detailed Description

Intended flow:

Wrapper class for Pythia usage in cosmic ray or detector cascades. Here it is assumed that the user wants to control the full evolution. The complete event record is then kept somewhere separate from PYTHIA, and specific input on production and decay generates the next subevent.- init sets up all generation, given a maximum energy. This may be time-consuming, less so if MPI initialization data can be reused. - sigmaSetuphN prepares a possible collision for a given hadron by calculating the hadron-nucleon cross section (if possible). Moderately time-consuming. - sigmahA uses the hN cross section calculated by sigaSetuphN and returns hadron-ion cross section for a collision with a specified nucleus. It can be called several times to cover a mix of target nuclei, with minimal time usage. - nextColl performs the hadron-nucleus collision, as a sequences of hadron-nucleon ones. Can be quite time-consuming. - nextDecay can be used anytime to decay a particle. Each individual decay is rather fast, but there may be many of them. - stat can be used at end of run to give a summary of error messages. Negligible time usage. - references to particleData() and rndm() can be used in the main

Member Function Documentation

void compress ( )
inline

Compress the event record by removing initial and intermediate particles. Keep line 0, since the += operator for event records only adds from 1 on.

Reference to the main event record. Original and new size.

Loop through all particles and move up the final ones to the top. Remove history information. Update event record size.

Shrink event record to new size.

bool init ( double  eMaxIn = 1e9,
bool  listFinalIn = false,
bool  rapidDecaysIn = false,
double  smallTau0In = 1e-10,
bool  reuseMPI = true,
string  initFile = "pythiaCascade.mpi" 
)
inline

Initialize PythiaCascade for a given maximal incoming energy.

Set eMax to the maximal incoming energy (in GeV) on fixed target. Keep listFinal = false if you want to get back the full event record, else only the "final" particles of the collision or decay are returned. Keep rapidDecays = false if you want to do every decay yourself, else all particle types with a tau0 below smallTau0 will be decayed immediately. The tau0 units are mm/c, so default gives c

  • tau0 = 1e-10 mm = 100 fm. Note that time dilation effects are not included, and they can be large. Set reuseMPI = false in case you cannot reuse the pythiaCascade.mpi file, or any other relevant stored file. The saved initialization must have been done with the same or a larger eMax than the one now being used.

Store input for future usage.

Proton mass.

Main Pythia object for managing the cascade evolution in a nucleus. Can also do decays, but no hard processes.

Redure statistics printout to relevant ones.

Initialize. Return if failure.

Secondary Pythia object for performing individual collisions, or decays. Variable incoming beam type and energy.

Set up for fixed-target collisions.

Must use the soft and low-energy QCD processes.

Primary (single) decay to be done by pythiaColl, to circumvent limitTau0.

Secondary decays to be done by pythiaMain, respecting limitTau0.

Reduce printout and relax energy-momentum conservation. (Unusually large errors unfortunate consequence of large boosts.)

Redure statistics printout to relevant ones.

Reuse MPI initialization file if it exists; else create a new one.

Initialize.

double nCollAvg ( int  A)
inline

Calculate the average number of collisions. Average number of hadron-nucleon collisions in a hadron-nucleus one. If not in table then interpolate, knowing that <n> - 1 propto A^{2/3}.

Event& nextColl ( int  Znow,
int  Anow,
Vec4  vNow = Vec4() 
)
inline

Generate a collision, and return the event record. Input (Z, A) of nucleus, and optionally collision vertex.

References to the two event records. Clear main event record.

Restrict to allowed range 1 <= A <= 208.

Insert incoming particle in cleared main event record.

Set up for collisions on a nucleus.

Drop rate of geometric series. (Deuterium has corrected nCollAvg.)

Loop over varying number of hit nucleons in target nucleus.

Pick incoming projectile: trivial for first subcollision, else ...

... find highest-pLongitudinal particle from latest subcollision.

No further subcollision if no particle with enough energy.

Choose process; only SD or ND at perturbative energies.

Pick one p or n from target.

Do a projectile-nucleon subcollision. Return empty event if failure.

Insert target nucleon. Mothers are (0,iProj) to mark who it interacted with. Always use proton mass for simplicity.

Update full energy of the event with the proton mass.

Insert secondary produced particles (but skip intermediate partons) into main event record and shift to correct production vertex.

Update daughters of colliding hadrons and other history.

End of loop over interactions in a nucleus.

Optionally do decays of short-lived particles.

Optionally compress event record.

Return generated collision.

Event& nextDecay ( int  idNowIn,
Vec4  pNowIn,
double  mNowIn,
Vec4  vNow = Vec4() 
)
inline

Generate a particle decay, and return the event record. You can allow sequential decays, if they occur rapidly enough.

Save incoming quantities. (Not needed, but by analogy with collisions.)

References to the event records. Clear them.

Insert incoming particle in cleared collision event record.

Decay incoming particle. Return empty event if fail. Copy event record.

Optionally do secondary decays of short-lived particles.

Optionally compress event record.

Return generated collision.

ParticleData& particleData ( )
inline

Possibility to access particle data and random numbers from pythiaMain.

double sigmahA ( int  A)
inline

Calculate the hadron-nucleus cross section for a given nucleon number A, using hN cross section from sigmaSetuphN. Interpolate where not (yet) available.

Restrict to allowed range 1 <= A <= 208.

Correction factor for number of h-nucleon collisions per h-nucleus one.

Done.

bool sigmaSetuphN ( int  idNowIn,
Vec4  pNowIn,
double  mNowIn 
)
inline

Calculate the hadron-proton collision cross section. Return false if not possible to find.

Cannot handle low-energy hadrons.

Cannot handle hadrons above maximum energy set at initialization.

Save incoming quantities for reuse in later methods.

Calculate hadron-nucleon cross section. Check if cross section vanishes.

Done.


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