1 | // main23.cc is a part of the PYTHIA event generator. |
---|
2 | // Copyright (C) 2012 Torbjorn Sjostrand. |
---|
3 | // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details. |
---|
4 | // Please respect the MCnet Guidelines, see GUIDELINES for details. |
---|
5 | |
---|
6 | // Example how to write a derived class for beam momentum and vertex spread, |
---|
7 | // with an instance handed to Pythia for internal generation. |
---|
8 | // Also how to write a derived class for external random numbers, |
---|
9 | // and how to write a derived class for external parton distributions. |
---|
10 | // Warning: the parameters are not realistic. |
---|
11 | |
---|
12 | #include "Pythia.h" |
---|
13 | |
---|
14 | using namespace Pythia8; |
---|
15 | |
---|
16 | //========================================================================== |
---|
17 | |
---|
18 | // A derived class to set beam momentum and interaction vertex spread. |
---|
19 | |
---|
20 | class MyBeamShape : public BeamShape { |
---|
21 | |
---|
22 | public: |
---|
23 | |
---|
24 | // Constructor. |
---|
25 | MyBeamShape() {} |
---|
26 | |
---|
27 | // Initialize beam parameters. |
---|
28 | // In this particular example we will reuse the existing settings names |
---|
29 | // but with modified meaning, so init() in the base class can be kept. |
---|
30 | //virtual void init( Settings& settings, Rndm* rndmPtrIn); |
---|
31 | |
---|
32 | // Set the two beam momentum deviations and the beam vertex. |
---|
33 | virtual void pick(); |
---|
34 | |
---|
35 | }; |
---|
36 | |
---|
37 | //-------------------------------------------------------------------------- |
---|
38 | |
---|
39 | // Set the two beam momentum deviations and the beam vertex. |
---|
40 | // Note that momenta are in units of GeV and vertices in mm, |
---|
41 | // always with c = 1, so that e.g. time is in mm/c. |
---|
42 | |
---|
43 | void MyBeamShape::pick() { |
---|
44 | |
---|
45 | // Reset all values. |
---|
46 | deltaPxA = deltaPyA = deltaPzA = deltaPxB = deltaPyB = deltaPzB |
---|
47 | = vertexX = vertexY = vertexZ = vertexT = 0.; |
---|
48 | |
---|
49 | // Set beam A transverse momentum deviation by a two-dimensional Gaussian. |
---|
50 | if (allowMomentumSpread) { |
---|
51 | double totalDev, gauss; |
---|
52 | do { |
---|
53 | totalDev = 0.; |
---|
54 | if (sigmaPxA > 0.) { |
---|
55 | gauss = rndmPtr->gauss(); |
---|
56 | deltaPxA = sigmaPxA * gauss; |
---|
57 | totalDev += gauss * gauss; |
---|
58 | } |
---|
59 | if (sigmaPyA > 0.) { |
---|
60 | gauss = rndmPtr->gauss(); |
---|
61 | deltaPyA = sigmaPyA * gauss; |
---|
62 | totalDev += gauss * gauss; |
---|
63 | } |
---|
64 | } while (totalDev > maxDevA * maxDevA); |
---|
65 | |
---|
66 | // Set beam A longitudinal momentum as a triangular shape. |
---|
67 | // Reuse sigmaPzA to represent maximum deviation in this case. |
---|
68 | if (sigmaPzA > 0.) { |
---|
69 | deltaPzA = sigmaPzA * ( 1. - sqrt(rndmPtr->flat()) ); |
---|
70 | if (rndmPtr->flat() < 0.5) deltaPzA = -deltaPzA; |
---|
71 | } |
---|
72 | |
---|
73 | // Set beam B transverse momentum deviation by a two-dimensional Gaussian. |
---|
74 | do { |
---|
75 | totalDev = 0.; |
---|
76 | if (sigmaPxB > 0.) { |
---|
77 | gauss = rndmPtr->gauss(); |
---|
78 | deltaPxB = sigmaPxB * gauss; |
---|
79 | totalDev += gauss * gauss; |
---|
80 | } |
---|
81 | if (sigmaPyB > 0.) { |
---|
82 | gauss = rndmPtr->gauss(); |
---|
83 | deltaPyB = sigmaPyB * gauss; |
---|
84 | totalDev += gauss * gauss; |
---|
85 | } |
---|
86 | } while (totalDev > maxDevB * maxDevB); |
---|
87 | |
---|
88 | // Set beam B longitudinal momentum as a triangular shape. |
---|
89 | // Reuse sigmaPzB to represent maximum deviation in this case. |
---|
90 | if (sigmaPzB > 0.) { |
---|
91 | deltaPzB = sigmaPzB * ( 1. - sqrt(rndmPtr->flat()) ); |
---|
92 | if (rndmPtr->flat() < 0.5) deltaPzB = -deltaPzB; |
---|
93 | } |
---|
94 | } |
---|
95 | |
---|
96 | // Set beam vertex location by a two-dimensional Gaussian. |
---|
97 | if (allowVertexSpread) { |
---|
98 | double totalDev, gauss; |
---|
99 | do { |
---|
100 | totalDev = 0.; |
---|
101 | if (sigmaVertexX > 0.) { |
---|
102 | gauss = rndmPtr->gauss(); |
---|
103 | vertexX = sigmaVertexX * gauss; |
---|
104 | totalDev += gauss * gauss; |
---|
105 | } |
---|
106 | if (sigmaVertexY > 0.) { |
---|
107 | gauss = rndmPtr->gauss(); |
---|
108 | vertexY = sigmaVertexY * gauss; |
---|
109 | totalDev += gauss * gauss; |
---|
110 | } |
---|
111 | } while (totalDev > maxDevVertex * maxDevVertex); |
---|
112 | |
---|
113 | // Set beam B longitudinal momentum as a triangular shape. |
---|
114 | // This corresponds to two step-function beams colliding. |
---|
115 | // Reuse sigmaVertexZ to represent maximum deviation in this case. |
---|
116 | if (sigmaVertexZ > 0.) { |
---|
117 | vertexZ = sigmaVertexZ * ( 1. - sqrt(rndmPtr->flat()) ); |
---|
118 | if (rndmPtr->flat() < 0.5) vertexZ = -vertexZ; |
---|
119 | |
---|
120 | // Set beam collision time flat between +-(sigmaVertexZ - |vertexZ|). |
---|
121 | // This corresponds to two step-function beams colliding (with v = c). |
---|
122 | vertexT = (2. * rndmPtr->flat() - 1.) * (sigmaVertexZ - abs(vertexZ)); |
---|
123 | } |
---|
124 | |
---|
125 | // Add offset to beam vertex. |
---|
126 | vertexX += offsetX; |
---|
127 | vertexY += offsetY; |
---|
128 | vertexZ += offsetZ; |
---|
129 | vertexT += offsetT; |
---|
130 | } |
---|
131 | |
---|
132 | } |
---|
133 | |
---|
134 | //========================================================================== |
---|
135 | |
---|
136 | // A derived class to generate random numbers. |
---|
137 | // A guranteed VERY STUPID generator, just to show principles. |
---|
138 | |
---|
139 | class stupidRndm : public RndmEngine { |
---|
140 | |
---|
141 | public: |
---|
142 | |
---|
143 | // Constructor. |
---|
144 | stupidRndm() { init();} |
---|
145 | |
---|
146 | // Routine for generating a random number. |
---|
147 | double flat(); |
---|
148 | |
---|
149 | private: |
---|
150 | |
---|
151 | // Initialization. |
---|
152 | void init(); |
---|
153 | |
---|
154 | // State of the generator. |
---|
155 | double value, exp10; |
---|
156 | |
---|
157 | }; |
---|
158 | |
---|
159 | //-------------------------------------------------------------------------- |
---|
160 | |
---|
161 | // Initialization method for the random numbers. |
---|
162 | |
---|
163 | void stupidRndm::init() { |
---|
164 | |
---|
165 | // Initial values. |
---|
166 | value = 0.5; |
---|
167 | exp10 = exp(10.); |
---|
168 | |
---|
169 | } |
---|
170 | |
---|
171 | //-------------------------------------------------------------------------- |
---|
172 | |
---|
173 | // Generation method for the random numbers. |
---|
174 | |
---|
175 | double stupidRndm::flat() { |
---|
176 | |
---|
177 | // Update counter. Add to current value. |
---|
178 | do { |
---|
179 | value *= exp10; |
---|
180 | value += M_PI; |
---|
181 | value -= double(int(value)); |
---|
182 | if (value < 0.) value += 1.; |
---|
183 | } while (value <= 0. || value >= 1.); |
---|
184 | |
---|
185 | // Return new value. |
---|
186 | return value; |
---|
187 | |
---|
188 | } |
---|
189 | |
---|
190 | //========================================================================== |
---|
191 | |
---|
192 | // A simple scaling PDF. Not realistic; only to check that it works. |
---|
193 | |
---|
194 | class Scaling : public PDF { |
---|
195 | |
---|
196 | public: |
---|
197 | |
---|
198 | // Constructor. |
---|
199 | Scaling(int idBeamIn = 2212) : PDF(idBeamIn) {} |
---|
200 | |
---|
201 | private: |
---|
202 | |
---|
203 | // Update PDF values. |
---|
204 | void xfUpdate(int id, double x, double Q2); |
---|
205 | |
---|
206 | }; |
---|
207 | |
---|
208 | //-------------------------------------------------------------------------- |
---|
209 | |
---|
210 | // No dependence on Q2, so leave out name for last argument. |
---|
211 | |
---|
212 | void Scaling::xfUpdate(int id, double x, double ) { |
---|
213 | |
---|
214 | // Valence quarks, carrying 60% of the momentum. |
---|
215 | double dv = 4. * x * pow3(1. - x); |
---|
216 | double uv = 2. * dv; |
---|
217 | |
---|
218 | // Gluons and sea quarks carrying the rest. |
---|
219 | double gl = 2. * pow5(1. - x); |
---|
220 | double sea = 0.4 * pow5(1. - x); |
---|
221 | |
---|
222 | // Update values |
---|
223 | xg = gl; |
---|
224 | xu = uv + 0.18 * sea; |
---|
225 | xd = dv + 0.18 * sea; |
---|
226 | xubar = 0.18 * sea; |
---|
227 | xdbar = 0.18 * sea; |
---|
228 | xs = 0.08 * sea; |
---|
229 | xc = 0.04 * sea; |
---|
230 | xb = 0.02 * sea; |
---|
231 | |
---|
232 | // Subdivision of valence and sea. |
---|
233 | xuVal = uv; |
---|
234 | xuSea = xubar; |
---|
235 | xdVal = dv; |
---|
236 | xdSea = xdbar; |
---|
237 | |
---|
238 | // idSav = 9 to indicate that all flavours reset. id change dummy. |
---|
239 | idSav = 9; |
---|
240 | id = 0; |
---|
241 | |
---|
242 | } |
---|
243 | |
---|
244 | //========================================================================== |
---|
245 | |
---|
246 | int main() { |
---|
247 | |
---|
248 | // Number of events to generate. Max number of errors. |
---|
249 | int nEvent = 10000; |
---|
250 | int nAbort = 5; |
---|
251 | |
---|
252 | // Pythia generator. |
---|
253 | Pythia pythia; |
---|
254 | |
---|
255 | // Process selection. |
---|
256 | pythia.readString("HardQCD:all = on"); |
---|
257 | pythia.readString("PhaseSpace:pTHatMin = 20."); |
---|
258 | |
---|
259 | // LHC with acollinear beams in the x plane. |
---|
260 | // Use that default is pp with pz = +-7000 GeV, so this need not be set. |
---|
261 | pythia.readString("Beams:frameType = 3"); |
---|
262 | pythia.readString("Beams:pxA = 1."); |
---|
263 | pythia.readString("Beams:pxB = 1."); |
---|
264 | |
---|
265 | // A class to generate beam parameters according to own parametrization. |
---|
266 | BeamShape* myBeamShape = new MyBeamShape(); |
---|
267 | |
---|
268 | // Hand pointer to Pythia. |
---|
269 | // If you comment this out you get internal Gaussian-style implementation. |
---|
270 | pythia.setBeamShapePtr( myBeamShape); |
---|
271 | |
---|
272 | // Set up beam spread parameters - reused by 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."); |
---|
280 | |
---|
281 | // Set up beam vertex parameters - reused by MyBeamShape. |
---|
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."); |
---|
286 | // In MyBeamShape the time width is not an independent parameter. |
---|
287 | //pythia.readString("Beams:sigmaTime = 50."); |
---|
288 | |
---|
289 | // Optionally simplify generation. |
---|
290 | pythia.readString("PartonLevel:all = off"); |
---|
291 | |
---|
292 | // A class to do random numbers externally. Hand pointer to Pythia. |
---|
293 | RndmEngine* badRndm = new stupidRndm(); |
---|
294 | pythia.setRndmEnginePtr( badRndm); |
---|
295 | |
---|
296 | // Two classes to do the two PDFs externally. Hand pointers to Pythia. |
---|
297 | PDF* pdfAPtr = new Scaling(2212); |
---|
298 | PDF* pdfBPtr = new Scaling(2212); |
---|
299 | pythia.setPDFPtr( pdfAPtr, pdfBPtr); |
---|
300 | |
---|
301 | // Initialization. |
---|
302 | pythia.init(); |
---|
303 | |
---|
304 | // Read out nominal energy. |
---|
305 | double eCMnom = pythia.info.eCM(); |
---|
306 | |
---|
307 | // Histograms. |
---|
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.); |
---|
318 | |
---|
319 | // Begin event loop. Generate event. |
---|
320 | int iAbort = 0; |
---|
321 | for (int iEvent = 0; iEvent < nEvent; ++iEvent) { |
---|
322 | if (!pythia.next()) { |
---|
323 | |
---|
324 | // List faulty events and quit if many failures. |
---|
325 | pythia.info.list(); |
---|
326 | pythia.process.list(); |
---|
327 | //pythia.event.list(); |
---|
328 | if (++iAbort < nAbort) continue; |
---|
329 | cout << " Event generation aborted prematurely, owing to error!\n"; |
---|
330 | break; |
---|
331 | } |
---|
332 | |
---|
333 | // Fill histograms. |
---|
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 ); |
---|
348 | |
---|
349 | // End of event loop. Statistics. Histograms. |
---|
350 | } |
---|
351 | pythia.stat(); |
---|
352 | cout << eCM << pXsum << pYsum << pZsum << pZind |
---|
353 | << vtxX << vtxY << vtxZ << vtxT << vtxZT; |
---|
354 | |
---|
355 | // Study standard Pythia random number generator. |
---|
356 | Hist rndmDist("standard random number distribution", 100, 0., 1.); |
---|
357 | Hist rndmCorr("standard random number correlation", 100, 0., 1.); |
---|
358 | double rndmNow; |
---|
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) ); |
---|
364 | rndmOld = rndmNow; |
---|
365 | } |
---|
366 | cout << rndmDist << rndmCorr; |
---|
367 | |
---|
368 | // Study bad "new" random number generator. |
---|
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) ); |
---|
376 | rndmOld = rndmNow; |
---|
377 | } |
---|
378 | cout << rndmDist2 << rndmCorr2; |
---|
379 | |
---|
380 | // Done. |
---|
381 | return 0; |
---|
382 | } |
---|