11 using namespace Pythia8;
19 void fillParticle(
int id,
double ee,
double thetaIn,
double phiIn,
26 double mm = pdt.mass(
id);
27 double pp = sqrtpos(ee*ee - mm*mm);
30 double cThe, sThe, phi;
36 cThe = 2. * rndm.flat() - 1.;
37 sThe = sqrtpos(1. - cThe * cThe);
38 phi = 2. * M_PI * rndm.flat();
42 event.append(
id, 1, 0, 0, pp * sThe * cos(phi), pp * sThe * sin(phi),
59 double mm = pdt.m0(
id);
60 double pp = sqrtpos(ee*ee - mm*mm);
61 event.append(
id, 23, 101, 0, 0., 0., pp, ee, mm);
62 event.append( -
id, 23, 0, 101, 0., 0., -pp, ee, mm);
65 }
else if (type == 2) {
66 event.append( 21, 23, 101, 102, 0., 0., ee, ee);
67 event.append( 21, 23, 102, 101, 0., 0., -ee, ee);
70 }
else if (type == 3) {
71 event.append( 21, 23, 101, 102, 0., 0., ee, ee);
72 event.append( 21, 23, 102, 103, 0.8 * ee, 0., -0.6 * ee, ee);
73 event.append( 21, 23, 103, 101, -0.8 * ee, 0., -0.6 * ee, ee);
76 }
else if (type == 4 || type == 5) {
79 event.append( 1000022, -21, 0, 0, 2, 4, 0, 0,
80 0., 0., 1.01 * ee, 1.01 * ee);
83 double rt75 = sqrt(0.75);
84 event.append( 2, 23, 1, 0, 0, 0, 101, 0,
85 0., 0., 1.01 * ee, 1.01 * ee);
86 event.append( 2, 23, 1, 0, 0, 0, 102, 0,
87 rt75 * ee, 0., -0.5 * ee, ee );
88 event.append( 1, 23, 1, 0, 0, 0, 103, 0,
89 -rt75 * ee, 0., -0.5 * ee, ee );
93 int colNow[4] = {0, 101, 102, 103};
95 pQ[1] =
Vec4(0., 0., 1., 0.);
96 pQ[2] =
Vec4( rt75, 0., -0.5, 0.);
97 pQ[3] =
Vec4(-rt75, 0., -0.5, 0.);
100 double cosThetaMin =0.;
103 for (
int nglu = 0; nglu < 5; ++nglu) {
104 int iq = 1 + int( 2.99999 * rndm.flat() );
105 double px, py, pz, e, prod;
107 e = ee * rndm.flat();
108 double cThe = 2. * rndm.flat() - 1.;
109 double phi = 2. * M_PI * rndm.flat();
110 px = e * sqrt(1. - cThe*cThe) * cos(phi);
111 py = e * sqrt(1. - cThe*cThe) * sin(phi);
113 prod = ( px * pQ[iq].px() + py * pQ[iq].py() + pz * pQ[iq].pz() )
115 }
while (prod < cosThetaMin);
116 int colNew = 104 + nglu;
117 event.append( 21, 23, 1, 0, 0, 0, colNew, colNow[iq],
122 event[1].daughters(1, event.size() - 1);
127 }
else if (type >= 6) {
130 event.append( 2212, -12, 0, 0, 3, 5, 0, 0, 0., 0., ee, ee, 0.);
131 event.append(-2212, -12, 0, 0, 6, 8, 0, 0, 0., 0., ee, ee, 0.);
135 double cThe = cos(theta);
136 double sThe = sin(theta);
139 int acol = (type == 6) ? 103 : 106;
143 event.append( 2, 23, 1, 0, 0, 0, 101, 0,
144 ee * sThe, 0., ee * cThe, ee, 0.);
145 event.append( 1, 23, 1, 0, 0, 0, 102, 0,
146 -ee * sThe, 0., ee * cThe, ee, 0.);
147 event.append( 2, -21, 1, 0, 0, 0, 103, 0,
148 0., 0., ee , ee, 0.);
149 event.append( -2, 23, 2, 0, 0, 0, 0, 104,
150 ee * sThe, 0., -ee * cThe, ee, 0.);
151 event.append( -1, 23, 2, 0, 0, 0, 0, 105,
152 -ee * sThe, 0., -ee * cThe, ee, 0.);
153 event.append( -2, -21, 2, 0, 0, 0, 0, acol,
154 0., 0., -ee , ee, 0.);
158 event.append( 21, 23, 5, 8, 0, 0, 103, 106, 0., ee, 0., ee, 0.);
159 }
else if (type == 8) {
160 event.append( 21, 23, 5, 8, 0, 0, 103, 108, 0., ee, 0., ee, 0.);
161 event.append( 21, 23, 5, 8, 0, 0, 108, 106, 0.,-ee, 0., ee, 0.);
162 }
else if (type == 9) {
163 event.append( 21, 23, 5, 8, 0, 0, 103, 107, 0., ee, 0., ee, 0.);
164 event.append( 21, 23, 5, 8, 0, 0, 107, 108, ee, 0., 0., ee, 0.);
165 event.append( 21, 23, 5, 8, 0, 0, 108, 106, 0.,-ee, 0., ee, 0.);
166 }
else if (type == 10) {
167 event.append( 21, 23, 5, 8, 0, 0, 103, 107, 0., ee, 0., ee, 0.);
168 event.append( 21, 23, 5, 8, 0, 0, 107, 108, ee, 0., 0., ee, 0.);
169 event.append( 21, 23, 5, 8, 0, 0, 108, 109, 0.,-ee, 0., ee, 0.);
170 event.append( 21, 23, 5, 8, 0, 0, 109, 106,-ee, 0., 0., ee, 0.);
205 Event&
event = pythia.event;
209 pythia.readString(
"ProcessLevel:all = off");
215 pythia.readString(
"Next:numberShowInfo = 0");
216 pythia.readString(
"Next:numberShowProcess = 0");
217 pythia.readString(
"Next:numberShowEvent = 0");
223 Hist epCons(
"deviation from energy-momentum conservation",100,0.,1e-4);
224 Hist chgCons(
"deviation from charge conservation",57,-9.5,9.5);
225 Hist nFinal(
"final particle multiplicity",100,-0.5,99.5);
226 Hist dnparticledp(
"dn/dp for particles",100,0.,ee);
227 Hist status85(
"multiplicity status code 85",50,-0.5,49.5);
228 Hist status86(
"multiplicity status code 86",50,-0.5,49.5);
229 Hist status83(
"multiplicity status code 83",50,-0.5,49.5);
230 Hist status84(
"multiplicity status code 84",50,-0.5,49.5);
231 Hist dndtheta(
"particle flow in event plane",100,-M_PI,M_PI);
232 Hist dedtheta(
"energy flow in event plane",100,-M_PI,M_PI);
233 Hist dpartondtheta(
"parton flow in event plane",100,-M_PI,M_PI);
234 Hist dndyAnti(
"dn/dy primaries antijunction",100, -10., 10.);
235 Hist dndyJun(
"dn/dy primaries junction",100, -10., 10.);
236 Hist dndySum(
"dn/dy all primaries",100, -10., 10.);
239 for (
int iEvent = 0; iEvent < nEvent; ++iEvent) {
242 if (type == 0) fillParticle( idGun, eeGun, -1., 0., event, pdt,
246 else fillPartons( type, ee, event, pdt, pythia.rndm);
249 if (!pythia.next()) {
250 cout <<
" Event generation aborted prematurely, owing to error!\n";
255 if (iEvent < nList) {
258 event.listJunctions();
262 Vec4 pSum = -
event[0].p();
263 double chargeSum = 0.;
264 if (type == 0) chargeSum = -
event[1].charge();
265 if (type == 4 || type == 5) chargeSum = -1;
273 for (
int i = 0; i <
event.size(); ++i) {
274 int status =
event[i].statusAbs();
277 int id =
event[i].id();
278 if (
id == 0 || !pdt.isParticle(
id))
279 cout <<
" Error! Unknown code id = " <<
id <<
"\n";
282 double eCalc =
event[i].eCalc();
283 if (abs(eCalc/event[i].e() - 1.) > 1e-6) cout <<
" e mismatch, i = "
284 << i <<
" e_nominal = " <<
event[i].e() <<
" e-from-p = "
285 << eCalc <<
" m-from-e " <<
event[i].mCalc() <<
"\n";
288 if (status == 71 || status == 72) {
289 double thetaXZ =
event[i].thetaXZ();
290 dpartondtheta.fill(thetaXZ);
294 if (status == 85) ++n85;
295 if (status == 86) ++n86;
296 if (status == 83) ++n83;
297 if (status == 84) ++n84;
300 if (status > 80 && status < 90) {
301 double eAbs =
event[i].e();
302 if (eAbs < 0.) {cout <<
" e < 0 line " << i;
event.list();}
303 double thetaXZ =
event[i].thetaXZ();
304 dndtheta.fill(thetaXZ);
305 dedtheta.fill(thetaXZ, eAbs);
308 double y =
event[i].y();
311 int motherId =
event[
event[i].mother1()].id();
312 if (motherId > 0 ) dndyJun.fill(event[i].y());
313 else dndyAnti.fill(event[i].y());
318 if (event[i].isFinal()) {
319 pSum +=
event[i].p();
320 chargeSum +=
event[i].charge();
322 double pAbs =
event[i].pAbs();
323 dnparticledp.fill(pAbs);
328 double epDev = abs(pSum.e()) + abs(pSum.px()) + abs(pSum.py())
331 chgCons.fill(chargeSum);
337 if (epDev > 1e-3 || abs(chargeSum) > 0.1)
event.list();
344 cout << epCons << chgCons << nFinal << dnparticledp
345 << dndtheta << dedtheta << dpartondtheta << dndySum;
346 if (type >= 4) cout << status85 << status86 << status83
348 if (type >= 6) cout << dndyJun << dndyAnti;