7 #ifndef Pythia8_PythiaCascade_H 8 #define Pythia8_PythiaCascade_H 10 #include "Pythia8/Pythia.h" 89 bool init(
double eMaxIn = 1e9,
bool listFinalIn =
false,
90 bool rapidDecaysIn =
false,
double smallTau0In = 1e-10,
91 bool reuseMPI =
true,
string initFile =
"pythiaCascade.mpi") {
95 listFinal = listFinalIn;
96 rapidDecays = rapidDecaysIn;
97 smallTau0 = smallTau0In;
104 pythiaMain.
readString(
"ProcessLevel:all = off");
109 pythiaMain.
settings.
flag(
"ParticleDecays:limitTau0", rapidDecays);
110 pythiaMain.
settings.parm(
"ParticleDecays:tau0Max", smallTau0);
113 pythiaMain.
readString(
"Stat:showProcessLevel = off");
114 pythiaMain.
readString(
"Stat:showPartonLevel = off");
117 if (!pythiaMain.
init())
return false;
121 pythiaColl.
readString(
"Beams:allowVariableEnergy = on");
122 pythiaColl.
readString(
"Beams:allowIDAswitch = on");
126 pythiaColl.
settings.parm(
"Beams:pzA", eMax);
131 pythiaColl.
readString(
"LowEnergyQCD:all = on");
141 pythiaColl.
readString(
"HadronLevel:Decay = off");
146 pythiaColl.
readString(
"Check:epTolErr = 0.01");
147 pythiaColl.
readString(
"Check:epTolWarn = 0.0001");
148 pythiaColl.
readString(
"Check:mTolErr = 0.01");
151 pythiaColl.
readString(
"Stat:showProcessLevel = off");
152 pythiaColl.
readString(
"Stat:showPartonLevel = off");
156 pythiaColl.
readString(
"MultipartonInteractions:reuseInit = 3");
157 else pythiaColl.
readString(
"MultipartonInteractions:reuseInit = 1");
158 pythiaColl.
settings.word(
"MultipartonInteractions:initFile", initFile);
161 if (!pythiaColl.
init())
return false;
173 for (
int i = 0; i < nA; ++i) {
175 return (sigmaNow < tabBorder) ? 1. + tabSlopeLo[i] * sigmaNow
176 : 1. + tabOffset[i] + tabSlope[i] * sigmaNow;
177 }
else if (A < tabA[i]) {
178 double nColl1 = (sigmaNow < tabBorder) ? tabSlopeLo[i - 1] * sigmaNow
179 : tabOffset[i - 1] + tabSlope[i - 1] * sigmaNow;
180 double nColl2 = (sigmaNow < tabBorder) ? tabSlopeLo[i] * sigmaNow
181 : tabOffset[i] + tabSlope[i] * sigmaNow;
182 double wt1 = (tabA[i] - A) / (tabA[i] - tabA[i - 1]);
183 return 1. + wt1 * pow( A / tabA[i - 1], 2./3.) * nColl1
184 + (1. - wt1) * pow( A / tabA[i], 2./3.) * nColl2;
188 return numeric_limits<double>::quiet_NaN();
199 if (pNowIn.e() - mNowIn < eKinMin)
return false;
202 if (pNowIn.e() > eMax) {
203 logger.ERROR_MSG(
"too high energy");
214 eCMNow = (pNow +
Vec4(0, 0, 0, mp)).mCalc();
215 sigmaNow = pythiaColl.
getSigmaTotal(idNow, 2212, eCMNow, mNow, mp);
216 if (sigmaNow <= 0.) {
217 if (eCMNow - mNow - mp > eKinMin)
218 logger.ERROR_MSG(
"vanishing cross section");
236 if (A < 1 || A > 208) {
237 logger.ERROR_MSG(
"A is outside of valid range (1 <= A <= 208)");
263 if (Anow < 1 || Anow > 208) {
264 logger.ERROR_MSG(
"A is outside of valid range (1 <= A <= 208)");
269 eventMain.
append(90, -11, 0, 0, 1, 1, 0, 0, pNow, mNow);
270 int iHad = eventMain.
append(idNow, 12, 0, 0, 0, 0, 0, 0, pNow, mNow);
271 eventMain[iHad].vProd(vNow);
275 int nn = Anow - Znow;
278 Vec4 dirNow = pNow / pNow.pAbs();
282 double probMore = 1. - 1. /
nCollAvg(Anow);
285 for (
int iColl = 1; iColl <= Anow; ++iColl) {
286 if (iColl > 1 && rndm.
flat() > probMore)
break;
296 for (
int i = sizeOld; i < sizeNew; ++i)
297 if ( eventMain[i].isFinal() && eventMain[i].isHadron()) {
298 double pp =
dot3(dirNow, eventMain[i].p());
306 if ( iProj == 0 || eventMain[iProj].e() - eventMain[iProj].
m()
310 double eCMSub = (eventMain[iProj].p() +
Vec4(0, 0, 0, mp)).mCalc();
311 if (eCMSub > 10.) procType = (rndm.
flat() < probSD) ? 4 : 1;
315 int idProj = eventMain[iProj].id();
316 bool doProton = rndm.
flat() < (np / double(np + nn));
317 if (doProton) np -= 1;
319 int idNuc = (doProton) ? 2212 : 2112;
324 if (!pythiaColl.
next(procType)) {
331 int statusNuc = (iColl == 1) ? -181 : -182;
332 int iNuc = eventMain.
append( idNuc, statusNuc, 0, iProj, 0, 0, 0, 0,
334 eventMain[iNuc].vProdAdd(vNow);
337 eventMain[0].e( eventMain[0].e() + mp);
338 eventMain[0].m( eventMain[0].p().mCalc() );
342 sizeOld = eventMain.
size();
343 for (
int iSub = 3; iSub < eventColl.
size(); ++iSub) {
344 if (!eventColl[iSub].isFinal())
continue;
345 int iNew = eventMain.
append(eventColl[iSub]);
346 eventMain[iNew].mothers(iNuc, iProj);
347 eventMain[iNew].vProdAdd(vNow);
349 sizeNew = eventMain.
size();
352 eventMain[iProj].daughters(sizeOld, sizeNew - 1);
353 eventMain[iNuc].daughters(sizeOld, sizeNew - 1);
354 eventMain[iProj].statusNeg();
355 eventMain[iProj].tau(0.);
392 eventColl.
append(90, -11, 0, 0, 1, 1, 0, 0, pNow, mNow);
393 int iHad = eventColl.
append(idNow, 12, 0, 0, 0, 0, 0, 0, pNow, mNow);
394 eventColl[iHad].vProd(vNow);
397 if (!pythiaColl.
moreDecays(iHad))
return eventMain;
398 eventMain = eventColl;
421 int sizeOld = eventMain.
size();
426 for (
int i = 1; i < sizeOld; ++i)
if (eventMain[i].isFinal()) {
427 eventMain[sizeNew] = eventMain[i];
428 eventMain[sizeNew].mothers( 0, 0);
429 eventMain[sizeNew].daughters( 0, 0);
434 eventMain.
popBack( sizeOld - sizeNew);
455 Rndm& rndm() {
return pythiaMain.
rndm;}
464 Pythia pythiaMain, pythiaColl;
470 bool listFinal, rapidDecays;
472 double eMax, smallTau0, mp, mNow, eCMNow, sigmaNow;
476 static constexpr
double probSD = 0.3;
480 static constexpr
double eKinMin = 0.2;
485 static constexpr
int nA = 16;
486 static constexpr
double tabBorder = 20.;
487 static const double tabA[nA];
488 static const double tabOffset[nA];
489 static const double tabSlope[nA];
490 static const double tabSlopeLo[nA];
495 const double PythiaCascade::tabA[] = {1, 2, 4, 9, 12, 14, 16, 27, 40, 56,
496 63, 84, 107, 129, 197, 208};
497 const double PythiaCascade::tabOffset[] = {0., 0.03, 0.08, 0.15, 0.20, 0.20,
498 0.20, 0.26, 0.30, 0.34, 0.40, 0.40, 0.40, 0.50, 0.50, 0.60};
499 const double PythiaCascade::tabSlope[] = {0., 0.0016, 0.0033, 0.0075, 0.0092,
500 0.0105, 0.012, 0.017, 0.022, 0.027, 0.028, 0.034, 0.040, 0.044, 0.055,
502 const double PythiaCascade::tabSlopeLo[] = {0., 0.0031, 0.0073, 0.015, 0.0192,
503 0.0205, 0.022, 0.03, 0.037, 0.044, 0.048, 0.054, 0.06, 0.069, 0.08, 0.085};
double nCollAvg(int A)
Definition: PythiaCascade.h:172
bool flag(string keyIn)
Give back current value, with check that key exists.
Definition: Settings.cc:1478
PythiaCascade()=default
Default constructor, all setup is done in init().
Rndm rndm
Random number generator.
Definition: Pythia.h:338
ParticleData particleData
ParticleData: the particle data table/database.
Definition: Pythia.h:335
The Event class holds all info on the generated event.
Definition: Event.h:453
Settings settings
Settings: databases of flags/modes/parms/words to control run.
Definition: Pythia.h:332
void clear()
Clear event record.
Definition: Event.h:473
bool setKinematics(double eCMIn)
Switch beam kinematics.
Definition: Pythia.cc:1420
Event event
The event record for the complete event history.
Definition: Pythia.h:323
bool init()
Initialize.
Definition: Pythia.cc:491
double dot3(const Vec4 &v1, const Vec4 &v2)
Scalar and cross product of 3-vector parts.
Definition: Basics.cc:625
void stat()
Summary of aborts, errors and warnings.
Definition: PythiaCascade.h:442
double sigmahA(int A)
Definition: PythiaCascade.h:233
bool setBeamIDs(int idAin, int idBin=0)
Switch to new beam particle identities; for similar hadrons only.
Definition: Pythia.cc:1398
bool moreDecays()
Special routine to allow more decays if on/off switches changed.
Definition: Pythia.h:253
int append(Particle entryIn)
Put a new particle at the end of the event record; return index.
Definition: Event.h:507
Event & nextDecay(int idNowIn, Vec4 pNowIn, double mNowIn, Vec4 vNow=Vec4())
Definition: PythiaCascade.h:376
bool sigmaSetuphN(int idNowIn, Vec4 pNowIn, double mNowIn)
Definition: PythiaCascade.h:196
Intended flow:
Definition: PythiaCascade.h:61
void stat()
Main routine to provide final statistics on generation.
Definition: Pythia.cc:1748
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.
Definition: PythiaCascade.h:89
void popBack(int nRemove=1)
Remove last n entries.
Definition: Event.h:561
int size() const
Event record size.
Definition: Event.h:504
double getSigmaTotal()
Get total cross section for two hadrons in the event record or standalone.
Definition: Pythia.h:266
Header for classes to set beam momentum and interaction vertex spread.
Definition: Analysis.h:20
The Pythia class contains the top-level routines to generate an event.
Definition: Pythia.h:71
double m(const Vec4 &v1)
Invariant mass and its square.
Definition: Basics.cc:595
void compress()
Definition: PythiaCascade.h:417
ParticleData & particleData()
Definition: PythiaCascade.h:453
bool next()
Generate the next event.
Definition: Pythia.h:229
bool readString(string, bool warn=true, int subrun=SUBRUNDEFAULT)
Read in one update for a setting or particle data from a single line.
Definition: Pythia.cc:365
double flat()
Generate next random number uniformly between 0 and 1.
Definition: Basics.cc:189
void errorStatistics() const
Print error statistics.
Definition: Logger.h:113
This class holds a map of all ParticleDataEntries.
Definition: ParticleData.h:422
Event & nextColl(int Znow, int Anow, Vec4 vNow=Vec4())
Definition: PythiaCascade.h:255