8 #ifndef Pythia8_LowEnergyProcess_H 9 #define Pythia8_LowEnergyProcess_H 11 #include "Pythia8/Basics.h" 12 #include "Pythia8/Event.h" 13 #include "Pythia8/FragmentationSystems.h" 14 #include "Pythia8/SigmaLowEnergy.h" 15 #include "Pythia8/MiniStringFragmentation.h" 16 #include "Pythia8/HadronWidths.h" 17 #include "Pythia8/NucleonExcitations.h" 18 #include "Pythia8/PythiaStdlib.h" 19 #include "Pythia8/StringFragmentation.h" 49 double bSlope(
int id1In,
int id2In,
double eCMIn,
double mAIn,
double mBIn,
50 int typeIn = 2) { id1 = id1In; id2 = id2In; eCM = eCMIn, sCM = eCM * eCM;
51 mA = mAIn; mB = mBIn; type = typeIn;
return bSlope();}
59 double probStoUD, fracEtass, fracEtaPss, xPowMes, xPowBar, xDiqEnhance,
60 sigmaQ, mStringMin, sProton, probDoubleAnn;
64 bool isBaryon1, isBaryon2;
65 int type, sizeOld, id1, id2, idc1, idac1, idc2, idac2, nHadron,
66 id1sv = {}, id2sv = {};
67 double m1, m2, eCM, sCM, mThr1, mThr2, z1, z2, mT1, mT2, mA, mB,
68 mc1, mac1, px1, py1, pTs1, mTsc1, mTsac1, mTc1, mTac1,
69 mc2, mac2, px2, py2, pTs2, mTsc2, mTsac2, mTc2, mTac2, bA, bB;
105 bool simpleHadronization();
114 bool splitA(
double mMax,
double redMpT,
bool splitFlavour =
true);
115 bool splitB(
double mMax,
double redMpT,
bool splitFlavour =
true);
118 pair< int, int> splitFlav(
int id);
121 double splitZ(
int iq1,
int iq2,
double mRat1,
double mRat2);
124 double mThreshold(
int iq1,
int iq2);
127 double mDiffThr(
int idNow,
double mNow);
Definition: LowEnergyProcess.h:28
bool collide(int i1, int i2, int typeIn, Event &event, Vec4 vtx=Vec4(), Vec4 vtx1=Vec4(), Vec4 vtx2=Vec4())
Produce outgoing primary hadrons from collision of incoming pair.
Definition: LowEnergyProcess.cc:105
double bSlope(int id1In, int id2In, double eCMIn, double mAIn, double mBIn, int typeIn=2)
Give access to b slope in elastic and diffractive interactions.
Definition: LowEnergyProcess.h:49
Definition: NucleonExcitations.h:23
Definition: PhysicsBase.h:27
The Event class holds all info on the generated event.
Definition: Event.h:453
Gets cross sections for hadron-hadron collisions at low energies.
Definition: SigmaLowEnergy.h:22
Definition: StringFragmentation.h:119
LowEnergyProcess()=default
Constructor.
The ColConfig class describes the colour configuration of the whole event.
Definition: FragmentationSystems.h:60
Event leEvent
Event record to handle hadronization.
Definition: LowEnergyProcess.h:46
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The StringFlav class is used to select quark and hadron flavours.
Definition: FragmentationFlavZpT.h:84
Definition: MiniStringFragmentation.h:30
void init(StringFlav *flavSelPtrIn, StringFragmentation *stringFragPtrIn, MiniStringFragmentation *ministringFragPtrIn, SigmaLowEnergy *sigmaLowEnergyPtrIn, NucleonExcitations *nucleonExcitationsPtrIn)
Initialize the class.
Definition: LowEnergyProcess.cc:50