14 using namespace Pythia8;
43 void MyBeamShape::pick() {
46 deltaPxA = deltaPyA = deltaPzA = deltaPxB = deltaPyB = deltaPzB
47 = vertexX = vertexY = vertexZ = vertexT = 0.;
50 if (allowMomentumSpread) {
51 double totalDev, gauss;
55 gauss = rndmPtr->gauss();
56 deltaPxA = sigmaPxA * gauss;
57 totalDev += gauss * gauss;
60 gauss = rndmPtr->gauss();
61 deltaPyA = sigmaPyA * gauss;
62 totalDev += gauss * gauss;
64 }
while (totalDev > maxDevA * maxDevA);
69 deltaPzA = sigmaPzA * ( 1. - sqrt(rndmPtr->flat()) );
70 if (rndmPtr->flat() < 0.5) deltaPzA = -deltaPzA;
77 gauss = rndmPtr->gauss();
78 deltaPxB = sigmaPxB * gauss;
79 totalDev += gauss * gauss;
82 gauss = rndmPtr->gauss();
83 deltaPyB = sigmaPyB * gauss;
84 totalDev += gauss * gauss;
86 }
while (totalDev > maxDevB * maxDevB);
91 deltaPzB = sigmaPzB * ( 1. - sqrt(rndmPtr->flat()) );
92 if (rndmPtr->flat() < 0.5) deltaPzB = -deltaPzB;
97 if (allowVertexSpread) {
98 double totalDev, gauss;
101 if (sigmaVertexX > 0.) {
102 gauss = rndmPtr->gauss();
103 vertexX = sigmaVertexX * gauss;
104 totalDev += gauss * gauss;
106 if (sigmaVertexY > 0.) {
107 gauss = rndmPtr->gauss();
108 vertexY = sigmaVertexY * gauss;
109 totalDev += gauss * gauss;
111 }
while (totalDev > maxDevVertex * maxDevVertex);
116 if (sigmaVertexZ > 0.) {
117 vertexZ = sigmaVertexZ * ( 1. - sqrt(rndmPtr->flat()) );
118 if (rndmPtr->flat() < 0.5) vertexZ = -vertexZ;
122 vertexT = (2. * rndmPtr->flat() - 1.) * (sigmaVertexZ - abs(vertexZ));
163 void stupidRndm::init() {
175 double stupidRndm::flat() {
181 value -= double(
int(value));
182 if (value < 0.) value += 1.;
183 }
while (value <= 0. || value >= 1.);
199 Scaling(
int idBeamIn = 2212) :
PDF(idBeamIn) {}
204 void xfUpdate(
int id,
double x,
double Q2);
212 void Scaling::xfUpdate(
int id,
double x,
double ) {
215 double dv = 4. * x * pow3(1. - x);
219 double gl = 2. * pow5(1. - x);
220 double sea = 0.4 * pow5(1. - x);
224 xu = uv + 0.18 * sea;
225 xd = dv + 0.18 * sea;
256 pythia.readString(
"HardQCD:all = on");
257 pythia.readString(
"PhaseSpace:pTHatMin = 20.");
261 pythia.readString(
"Beams:frameType = 3");
262 pythia.readString(
"Beams:pxA = 1.");
263 pythia.readString(
"Beams:pxB = 1.");
270 pythia.setBeamShapePtr( myBeamShape);
273 pythia.readString(
"Beams:allowMomentumSpread = on");
274 pythia.readString(
"Beams:sigmapxA = 0.1");
275 pythia.readString(
"Beams:sigmapyA = 0.1");
276 pythia.readString(
"Beams:sigmapzA = 5.");
277 pythia.readString(
"Beams:sigmapxB = 0.1");
278 pythia.readString(
"Beams:sigmapyB = 0.1");
279 pythia.readString(
"Beams:sigmapzB = 5.");
282 pythia.readString(
"Beams:allowVertexSpread = on");
283 pythia.readString(
"Beams:sigmaVertexX = 0.3");
284 pythia.readString(
"Beams:sigmaVertexY = 0.3");
285 pythia.readString(
"Beams:sigmaVertexZ = 50.");
290 pythia.readString(
"PartonLevel:all = off");
294 pythia.setRndmEnginePtr( badRndm);
299 pythia.setPDFPtr( pdfAPtr, pdfBPtr);
305 double eCMnom = pythia.info.eCM();
308 Hist eCM(
"center-of-mass energy deviation", 100, -20., 20.);
309 Hist pXsum(
"net pX offset", 100, -1.0, 1.0);
310 Hist pYsum(
"net pY offset", 100, -1.0, 1.0);
311 Hist pZsum(
"net pZ offset", 100, -10., 10.);
312 Hist pZind(
"individual abs(pZ) offset", 100, -10., 10.);
313 Hist vtxX(
"vertex x position", 100, -1.0, 1.0);
314 Hist vtxY(
"vertex y position", 100, -1.0, 1.0);
315 Hist vtxZ(
"vertex z position", 100, -100., 100.);
316 Hist vtxT(
"vertex time", 100, -100., 100.);
317 Hist vtxZT(
"vertex |x| + |t|", 100, 0., 100.);
321 for (
int iEvent = 0; iEvent < nEvent; ++iEvent) {
322 if (!pythia.next()) {
326 pythia.process.list();
328 if (++iAbort < nAbort)
continue;
329 cout <<
" Event generation aborted prematurely, owing to error!\n";
334 double eCMnow = pythia.info.eCM();
335 eCM.fill( eCMnow - eCMnom);
336 pXsum.fill( pythia.process[0].px() - 2. );
337 pYsum.fill( pythia.process[0].py() );
338 pZsum.fill( pythia.process[0].pz() );
339 pZind.fill( pythia.process[1].pz() - 7000. );
340 pZind.fill( -pythia.process[2].pz() - 7000. );
341 vtxX.fill( pythia.process[0].xProd() );
342 vtxY.fill( pythia.process[0].yProd() );
343 vtxZ.fill( pythia.process[0].zProd() );
344 vtxT.fill( pythia.process[0].tProd() );
345 double absSum = abs(pythia.process[0].zProd())
346 + abs(pythia.process[0].tProd());
347 vtxZT.fill( absSum );
352 cout << eCM << pXsum << pYsum << pZsum << pZind
353 << vtxX << vtxY << vtxZ << vtxT << vtxZT;
356 Hist rndmDist(
"standard random number distribution", 100, 0., 1.);
357 Hist rndmCorr(
"standard random number correlation", 100, 0., 1.);
359 double rndmOld = pythia.rndm.flat();
360 for (
int i = 0; i < 100000; ++i) {
361 rndmNow = pythia.rndm.flat();
362 rndmDist.fill(rndmNow);
363 rndmCorr.fill( abs(rndmNow - rndmOld) );
366 cout << rndmDist << rndmCorr;
369 Hist rndmDist2(
"stupid random number distribution", 100, 0., 1.);
370 Hist rndmCorr2(
"stupid random number correlation", 100, 0., 1.);
371 rndmOld = pythia.rndm.flat();
372 for (
int i = 0; i < 100000; ++i) {
373 rndmNow = pythia.rndm.flat();
374 rndmDist2.fill(rndmNow);
375 rndmCorr2.fill( abs(rndmNow - rndmOld) );
378 cout << rndmDist2 << rndmCorr2;