]>
Commit | Line | Data |
---|---|---|
b584e2f5 | 1 | // HadronLevel.cc is a part of the PYTHIA event generator. |
2 | // Copyright (C) 2010 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 | // Function definitions (not found in the header) for the HadronLevel class. | |
7 | ||
8 | #include "HadronLevel.h" | |
9 | ||
10 | namespace Pythia8 { | |
11 | ||
12 | //========================================================================== | |
13 | ||
14 | // The HadronLevel class. | |
15 | ||
16 | //-------------------------------------------------------------------------- | |
17 | ||
18 | // Constants: could be changed here if desired, but normally should not. | |
19 | // These are of technical nature, as described for each. | |
20 | ||
21 | // For breaking J-J string, pick a Gamma by taking a step with fictitious mass. | |
22 | const double HadronLevel::JJSTRINGM2MAX = 25.; | |
23 | const double HadronLevel::JJSTRINGM2FRAC = 0.1; | |
24 | ||
25 | // Iterate junction rest frame boost until convergence or too many tries. | |
26 | const double HadronLevel::CONVJNREST = 1e-5; | |
27 | const int HadronLevel::NTRYJNREST = 20; | |
28 | ||
29 | // Typical average transvere primary hadron mass <mThad>. | |
30 | const double HadronLevel::MTHAD = 0.9; | |
31 | ||
32 | //-------------------------------------------------------------------------- | |
33 | ||
34 | // Find settings. Initialize HadronLevel classes as required. | |
35 | ||
36 | bool HadronLevel::init(Info* infoPtrIn, Settings& settings, | |
37 | ParticleData& particleData, Rndm* rndmPtrIn, TimeShower* timesDecPtr, | |
38 | DecayHandler* decayHandlePtr, vector<int> handledParticles) { | |
39 | ||
40 | // Save pointers. | |
41 | infoPtr = infoPtrIn; | |
42 | rndmPtr = rndmPtrIn; | |
43 | ||
44 | // Main flags. | |
45 | doHadronize = settings.flag("HadronLevel:Hadronize"); | |
46 | doDecay = settings.flag("HadronLevel:Decay"); | |
47 | doBoseEinstein = settings.flag("HadronLevel:BoseEinstein"); | |
48 | ||
49 | // Boundary mass between string and ministring handling. | |
50 | mStringMin = settings.parm("HadronLevel:mStringMin"); | |
51 | ||
52 | // For junction processing. | |
53 | eNormJunction = settings.parm("StringFragmentation:eNormJunction"); | |
54 | ||
55 | // Particles that should decay or not before Bose-Einstein stage. | |
56 | widthSepBE = settings.parm("BoseEinstein:widthSep"); | |
57 | ||
58 | // Initialize string and ministring fragmentation. | |
59 | stringFrag.init(infoPtr, settings, &particleData, rndmPtr, &flavSel, &pTSel, | |
60 | &zSel); | |
61 | ministringFrag.init(infoPtr, settings, &particleData, rndmPtr, &flavSel); | |
62 | ||
63 | // Initialize particle decays. | |
64 | decays.init(infoPtr, settings, &particleData, rndmPtr, timesDecPtr, &flavSel, | |
65 | decayHandlePtr, handledParticles); | |
66 | ||
67 | // Initialize BoseEinstein. | |
68 | boseEinstein.init(infoPtr, settings, particleData); | |
69 | ||
70 | // Initialize auxiliary administrative classes. | |
71 | colConfig.init(infoPtr, settings, &flavSel); | |
72 | ||
73 | // Initialize auxiliary fragmentation classes. | |
74 | flavSel.init(settings, rndmPtr); | |
75 | pTSel.init(settings, rndmPtr); | |
76 | zSel.init(settings, particleData, rndmPtr); | |
77 | ||
78 | // Done. | |
79 | return true; | |
80 | ||
81 | } | |
82 | ||
83 | //-------------------------------------------------------------------------- | |
84 | ||
85 | // Hadronize and decay the next parton-level. | |
86 | ||
87 | bool HadronLevel::next( Event& event) { | |
88 | ||
89 | // Colour-octet onia states must be decayed to singlet + gluon. | |
90 | if (!decayOctetOnia(event)) return false; | |
91 | ||
92 | // Possibility of hadronization inside decay, but then no BE second time. | |
93 | bool moreToDo; | |
94 | bool doBoseEinsteinNow = doBoseEinstein; | |
95 | do { | |
96 | moreToDo = false; | |
97 | ||
98 | // First part: string fragmentation. | |
99 | if (doHadronize) { | |
100 | ||
101 | // Find the complete colour singlet configuration of the event. | |
102 | if (!findSinglets( event)) return false; | |
103 | ||
104 | // Process all colour singlet (sub)system | |
105 | for (int iSub = 0; iSub < colConfig.size(); ++iSub) { | |
106 | ||
107 | // Collect sequentially all partons in a colour singlet subsystem. | |
108 | colConfig.collect(iSub, event); | |
109 | ||
110 | // String fragmentation of each colour singlet (sub)system. | |
111 | if ( colConfig[iSub].massExcess > mStringMin ) { | |
112 | if (!stringFrag.fragment( iSub, colConfig, event)) return false; | |
113 | ||
114 | // Low-mass string treated separately. Tell if diffractive system. | |
115 | } else { | |
116 | bool isDiff = infoPtr->isDiffractiveA() | |
117 | || infoPtr->isDiffractiveB(); | |
118 | if (!ministringFrag.fragment( iSub, colConfig, event, isDiff)) | |
119 | return false; | |
120 | } | |
121 | } | |
122 | } | |
123 | ||
124 | // Second part: sequential decays of short-lived particles (incl. K0). | |
125 | if (doDecay) { | |
126 | ||
127 | // Loop through all entries to find those that should decay. | |
128 | int iDec = 0; | |
129 | do { | |
130 | Particle& decayer = event[iDec]; | |
131 | if ( decayer.isFinal() && decayer.canDecay() && decayer.mayDecay() | |
132 | && (decayer.mWidth() > widthSepBE || decayer.idAbs() == 311) ) { | |
133 | decays.decay( iDec, event); | |
134 | if (decays.moreToDo()) moreToDo = true; | |
135 | } | |
136 | } while (++iDec < event.size()); | |
137 | } | |
138 | ||
139 | // Third part: include Bose-Einstein effects among current particles. | |
140 | if (doBoseEinsteinNow) { | |
141 | if (!boseEinstein.shiftEvent(event)) return false; | |
142 | doBoseEinsteinNow = false; | |
143 | } | |
144 | ||
145 | // Fourth part: sequential decays also of long-lived particles. | |
146 | if (doDecay) { | |
147 | ||
148 | // Loop through all entries to find those that should decay. | |
149 | int iDec = 0; | |
150 | do { | |
151 | Particle& decayer = event[iDec]; | |
152 | if ( decayer.isFinal() && decayer.canDecay() && decayer.mayDecay() ) { | |
153 | decays.decay( iDec, event); | |
154 | if (decays.moreToDo()) moreToDo = true; | |
155 | } | |
156 | } while (++iDec < event.size()); | |
157 | } | |
158 | ||
159 | // Normally done first time around, but sometimes not (e.g. Upsilon). | |
160 | } while (moreToDo); | |
161 | ||
162 | // Done. | |
163 | return true; | |
164 | ||
165 | } | |
166 | ||
167 | //-------------------------------------------------------------------------- | |
168 | ||
169 | // Allow more decays if on/off switches changed. | |
170 | // Note: does not do sequential hadronization, e.g. for Upsilon. | |
171 | ||
172 | bool HadronLevel::moreDecays( Event& event) { | |
173 | ||
174 | // Colour-octet onia states must be decayed to singlet + gluon. | |
175 | if (!decayOctetOnia(event)) return false; | |
176 | ||
177 | // Loop through all entries to find those that should decay. | |
178 | int iDec = 0; | |
179 | do { | |
180 | if ( event[iDec].isFinal() && event[iDec].canDecay() | |
181 | && event[iDec].mayDecay() ) decays.decay( iDec, event); | |
182 | } while (++iDec < event.size()); | |
183 | ||
184 | // Done. | |
185 | return true; | |
186 | ||
187 | } | |
188 | ||
189 | //-------------------------------------------------------------------------- | |
190 | ||
191 | // Decay colour-octet onium states. | |
192 | ||
193 | bool HadronLevel::decayOctetOnia(Event& event) { | |
194 | ||
195 | // Onium states to be decayed. | |
196 | int idOnium[6] = { 9900443, 9900441, 9910441, | |
197 | 9900553, 9900551, 9910551 }; | |
198 | ||
199 | // Loop over particles and identify onia. | |
200 | for (int iDec = 0; iDec < event.size(); ++iDec) | |
201 | if (event[iDec].isFinal()) { | |
202 | int id = event[iDec].id(); | |
203 | bool isOnium = false; | |
204 | for (int j = 0; j < 6; ++j) if (id == idOnium[j]) isOnium = true; | |
205 | ||
206 | // Decay any onia encountered. | |
207 | if (isOnium) { | |
208 | if (!decays.decay( iDec, event)) return false; | |
209 | ||
210 | // Set colour flow by hand: gluon inherits octet-onium state. | |
211 | int iGlu = event.size() - 1; | |
212 | event[iGlu].cols( event[iDec].col(), event[iDec].acol() ); | |
213 | } | |
214 | } | |
215 | ||
216 | // Done. | |
217 | return true; | |
218 | ||
219 | } | |
220 | ||
221 | //-------------------------------------------------------------------------- | |
222 | ||
223 | // Trace colour flow in the event to form colour singlet subsystems. | |
224 | ||
225 | bool HadronLevel::findSinglets(Event& event) { | |
226 | ||
227 | // Find a list of final partons and of all colour ends and gluons. | |
228 | iColEnd.resize(0); | |
229 | iAcolEnd.resize(0); | |
230 | iColAndAcol.resize(0); | |
231 | for (int i = 0; i < event.size(); ++i) if (event[i].isFinal()) { | |
232 | if (event[i].col() > 0 && event[i].acol() > 0) iColAndAcol.push_back(i); | |
233 | else if (event[i].col() > 0) iColEnd.push_back(i); | |
234 | else if (event[i].acol() > 0) iAcolEnd.push_back(i); | |
235 | } | |
236 | ||
237 | // Begin arrange the partons into separate colour singlets. | |
238 | colConfig.clear(); | |
239 | iPartonJun.resize(0); | |
240 | iPartonAntiJun.resize(0); | |
241 | ||
242 | // Junctions: loop over them, and identify kind. | |
243 | for (int iJun = 0; iJun < event.sizeJunction(); ++iJun) | |
244 | if (event.remainsJunction(iJun)) { | |
245 | event.remainsJunction(iJun, false); | |
246 | int kindJun = event.kindJunction(iJun); | |
247 | iParton.resize(0); | |
248 | ||
249 | // First kind: three (anti)colours out from junction. | |
250 | if (kindJun == 1 || kindJun == 2) { | |
251 | for (int iCol = 0; iCol < 3; ++iCol) { | |
252 | int indxCol = event.colJunction(iJun, iCol); | |
253 | iParton.push_back( -(10 + 10 * iJun + iCol) ); | |
254 | if (kindJun == 1 && !traceFromAcol(indxCol, event, iJun, iCol)) | |
255 | return false; | |
256 | if (kindJun == 2 && !traceFromCol(indxCol, event, iJun, iCol)) | |
257 | return false; | |
258 | } | |
259 | } | |
260 | ||
261 | // Keep in memory a junction hooked up with an antijunction, | |
262 | // else store found single-junction system. | |
263 | int nNeg = 0; | |
264 | for (int i = 0; i < int(iParton.size()); ++i) if (iParton[i] < 0) | |
265 | ++nNeg; | |
266 | if (nNeg > 3 && kindJun == 1) { | |
267 | for (int i = 0; i < int(iParton.size()); ++i) | |
268 | iPartonJun.push_back(iParton[i]); | |
269 | } else if (nNeg > 3 && kindJun == 2) { | |
270 | for (int i = 0; i < int(iParton.size()); ++i) | |
271 | iPartonAntiJun.push_back(iParton[i]); | |
272 | } else { | |
273 | // A junction may be eliminated by insert if two quarks are nearby. | |
274 | int nJunOld = event.sizeJunction(); | |
275 | if (!colConfig.insert(iParton, event)) return false; | |
276 | if (event.sizeJunction() < nJunOld) --iJun; | |
277 | } | |
278 | } | |
279 | ||
280 | // Split junction-antijunction system into two, and store those. | |
281 | // (Only one system in extreme cases, and then second empty.) | |
282 | if (iPartonJun.size() > 0 && iPartonAntiJun.size() > 0) { | |
283 | if (!splitJunctionPair(event)) return false; | |
284 | if (!colConfig.insert(iPartonJun, event)) return false; | |
285 | if (iPartonAntiJun.size() > 0) | |
286 | if (!colConfig.insert(iPartonAntiJun, event)) return false; | |
287 | // Error if only one of junction and antijuction left here. | |
288 | } else if (iPartonJun.size() > 0 || iPartonAntiJun.size() > 0) { | |
289 | infoPtr->errorMsg("Error in HadronLevel::findSinglets: " | |
290 | "unmatched (anti)junction"); | |
291 | return false; | |
292 | } | |
293 | ||
294 | // Open strings: pick up each colour end and trace to its anticolor end. | |
295 | for (int iEnd = 0; iEnd < int(iColEnd.size()); ++iEnd) { | |
296 | iParton.resize(0); | |
297 | iParton.push_back( iColEnd[iEnd] ); | |
298 | int indxCol = event[ iColEnd[iEnd] ].col(); | |
299 | if (!traceFromCol(indxCol, event)) return false; | |
300 | ||
301 | // Store found open string system. Analyze its properties. | |
302 | if (!colConfig.insert(iParton, event)) return false; | |
303 | } | |
304 | ||
305 | // Closed strings : begin at any gluon and trace until back at it. | |
306 | while (iColAndAcol.size() > 0) { | |
307 | iParton.resize(0); | |
308 | iParton.push_back( iColAndAcol[0] ); | |
309 | int indxCol = event[ iColAndAcol[0] ].col(); | |
310 | int indxAcol = event[ iColAndAcol[0] ].acol(); | |
311 | iColAndAcol[0] = iColAndAcol.back(); | |
312 | iColAndAcol.pop_back(); | |
313 | if (!traceInLoop(indxCol, indxAcol, event)) return false; | |
314 | ||
315 | // Store found closed string system. Analyze its properties. | |
316 | if (!colConfig.insert(iParton, event)) return false; | |
317 | } | |
318 | ||
319 | // Done. | |
320 | return true; | |
321 | ||
322 | } | |
323 | ||
324 | //-------------------------------------------------------------------------- | |
325 | ||
326 | // Trace a colour line, from a colour to an anticolour. | |
327 | ||
328 | bool HadronLevel::traceFromCol(int indxCol, Event& event, int iJun, | |
329 | int iCol) { | |
330 | ||
331 | // Junction kind, if any. | |
332 | int kindJun = (iJun >= 0) ? event.kindJunction(iJun) : 0; | |
333 | ||
334 | // Begin to look for a matching anticolour. | |
335 | int loop = 0; | |
336 | int loopMax = iColAndAcol.size() + 2; | |
337 | bool hasFound = false; | |
338 | do { | |
339 | ++loop; | |
340 | hasFound= false; | |
341 | ||
342 | // First check list of matching anticolour ends. | |
343 | for (int i = 0; i < int(iAcolEnd.size()); ++i) | |
344 | if (event[ iAcolEnd[i] ].acol() == indxCol) { | |
345 | iParton.push_back( iAcolEnd[i] ); | |
346 | indxCol = 0; | |
347 | iAcolEnd[i] = iAcolEnd.back(); | |
348 | iAcolEnd.pop_back(); | |
349 | hasFound = true; | |
350 | break; | |
351 | } | |
352 | ||
353 | // Then check list of intermediate gluons. | |
354 | if (!hasFound) | |
355 | for (int i = 0; i < int(iColAndAcol.size()); ++i) | |
356 | if (event[ iColAndAcol[i] ].acol() == indxCol) { | |
357 | iParton.push_back( iColAndAcol[i] ); | |
358 | ||
359 | // Update to new colour. Remove gluon. | |
360 | indxCol = event[ iColAndAcol[i] ].col(); | |
361 | if (kindJun > 0) event.endColJunction(iJun, iCol, indxCol); | |
362 | iColAndAcol[i] = iColAndAcol.back(); | |
363 | iColAndAcol.pop_back(); | |
364 | hasFound = true; | |
365 | break; | |
366 | } | |
367 | ||
368 | // In a pinch, check list of end colours on other (anti)junction. | |
369 | if (!hasFound && kindJun == 2 && event.sizeJunction() > 1) | |
370 | for (int iAntiJun = 0; iAntiJun < event.sizeJunction(); ++iAntiJun) | |
371 | if (iAntiJun != iJun && event.kindJunction(iAntiJun) == 1) | |
372 | for (int iColAnti = 0; iColAnti < 3; ++iColAnti) | |
373 | if (event.endColJunction(iAntiJun, iColAnti) == indxCol) { | |
374 | iParton.push_back( -(10 + 10 * iAntiJun + iColAnti) ); | |
375 | indxCol = 0; | |
376 | hasFound = true; | |
377 | break; | |
378 | } | |
379 | ||
380 | // Keep on tracing via gluons until reached end of leg. | |
381 | } while (hasFound && indxCol > 0 && loop < loopMax); | |
382 | ||
383 | // Something went wrong in colour tracing. | |
384 | if (!hasFound || loop == loopMax) { | |
385 | infoPtr->errorMsg("Error in HadronLevel::traceFromCol: " | |
386 | "colour tracing failed"); | |
387 | return false; | |
388 | } | |
389 | ||
390 | // Done. | |
391 | return true; | |
392 | ||
393 | } | |
394 | ||
395 | //-------------------------------------------------------------------------- | |
396 | ||
397 | // Trace a colour line, from an anticolour to a colour. | |
398 | ||
399 | bool HadronLevel::traceFromAcol(int indxCol, Event& event, int iJun, | |
400 | int iCol) { | |
401 | ||
402 | // Junction kind, if any. | |
403 | int kindJun = (iJun >= 0) ? event.kindJunction(iJun) : 0; | |
404 | ||
405 | // Begin to look for a matching colour. | |
406 | int loop = 0; | |
407 | int loopMax = iColAndAcol.size() + 2; | |
408 | bool hasFound = false; | |
409 | do { | |
410 | ++loop; | |
411 | hasFound= false; | |
412 | ||
413 | // First check list of matching colour ends. | |
414 | for (int i = 0; i < int(iColEnd.size()); ++i) | |
415 | if (event[ iColEnd[i] ].col() == indxCol) { | |
416 | iParton.push_back( iColEnd[i] ); | |
417 | indxCol = 0; | |
418 | iColEnd[i] = iColEnd.back(); | |
419 | iColEnd.pop_back(); | |
420 | hasFound = true; | |
421 | break; | |
422 | } | |
423 | ||
424 | // Then check list of intermediate gluons. | |
425 | if (!hasFound) | |
426 | for (int i = 0; i < int(iColAndAcol.size()); ++i) | |
427 | if (event[ iColAndAcol[i] ].col() == indxCol) { | |
428 | iParton.push_back( iColAndAcol[i] ); | |
429 | ||
430 | // Update to new colour. Remove gluon. | |
431 | indxCol = event[ iColAndAcol[i] ].acol(); | |
432 | if (kindJun > 0) event.endColJunction(iJun, iCol, indxCol); | |
433 | iColAndAcol[i] = iColAndAcol.back(); | |
434 | iColAndAcol.pop_back(); | |
435 | hasFound = true; | |
436 | break; | |
437 | } | |
438 | ||
439 | // In a pinch, check list of colours on other (anti)junction. | |
440 | if (!hasFound && kindJun == 1 && event.sizeJunction() > 1) | |
441 | for (int iAntiJun = 0; iAntiJun < event.sizeJunction(); ++iAntiJun) | |
442 | if (iAntiJun != iJun && event.kindJunction(iAntiJun) == 2) | |
443 | for (int iColAnti = 0; iColAnti < 3; ++iColAnti) | |
444 | if (event.endColJunction(iAntiJun, iColAnti) == indxCol) { | |
445 | iParton.push_back( -(10 + 10 * iAntiJun + iColAnti) ); | |
446 | indxCol = 0; | |
447 | hasFound = true; | |
448 | break; | |
449 | } | |
450 | ||
451 | // Keep on tracing via gluons until reached end of leg. | |
452 | } while (hasFound && indxCol > 0 && loop < loopMax); | |
453 | ||
454 | // Something went wrong in colour tracing. | |
455 | if (!hasFound || loop == loopMax) { | |
456 | infoPtr->errorMsg("Error in HadronLevel::traceFromAcol: " | |
457 | "colour tracing failed"); | |
458 | return false; | |
459 | } | |
460 | ||
461 | // Done. | |
462 | return true; | |
463 | ||
464 | } | |
465 | ||
466 | //-------------------------------------------------------------------------- | |
467 | ||
468 | // Trace a colour loop, from a colour back to the anticolour of the same. | |
469 | ||
470 | bool HadronLevel::traceInLoop(int indxCol, int indxAcol, Event& event) { | |
471 | ||
472 | // Move around until back where begun. | |
473 | int loop = 0; | |
474 | int loopMax = iColAndAcol.size() + 2; | |
475 | bool hasFound = false; | |
476 | do { | |
477 | ++loop; | |
478 | hasFound= false; | |
479 | ||
480 | // Check list of gluons. | |
481 | for (int i = 0; i < int(iColAndAcol.size()); ++i) | |
482 | if (event[ iColAndAcol[i] ].acol() == indxCol) { | |
483 | iParton.push_back( iColAndAcol[i] ); | |
484 | indxCol = event[ iColAndAcol[i] ].col(); | |
485 | iColAndAcol[i] = iColAndAcol.back(); | |
486 | iColAndAcol.pop_back(); | |
487 | hasFound = true; | |
488 | break; | |
489 | } | |
490 | } while (hasFound && indxCol != indxAcol && loop < loopMax); | |
491 | ||
492 | // Something went wrong in colour tracing. | |
493 | if (!hasFound || loop == loopMax) { | |
494 | infoPtr->errorMsg("Error in HadronLevel::traceInLoop: " | |
495 | "colour tracing failed"); | |
496 | return false; | |
497 | } | |
498 | ||
499 | // Done. | |
500 | return true; | |
501 | ||
502 | } | |
503 | ||
504 | //-------------------------------------------------------------------------- | |
505 | ||
506 | // Split junction-antijunction system into two, or simplify other way. | |
507 | ||
508 | bool HadronLevel::splitJunctionPair(Event& event) { | |
509 | ||
510 | // Construct separate index arrays for the three junction legs. | |
511 | int identJun = (-iPartonJun[0])/10; | |
512 | iJunLegA.resize(0); | |
513 | iJunLegB.resize(0); | |
514 | iJunLegC.resize(0); | |
515 | int leg = -1; | |
516 | for (int i = 0; i < int(iPartonJun.size()); ++ i) { | |
517 | if ( (-iPartonJun[i])/10 == identJun) ++leg; | |
518 | if (leg == 0) iJunLegA.push_back( iPartonJun[i] ); | |
519 | else if (leg == 1) iJunLegB.push_back( iPartonJun[i] ); | |
520 | else iJunLegC.push_back( iPartonJun[i] ); | |
521 | } | |
522 | ||
523 | // Construct separate index arrays for the three antijunction legs. | |
524 | int identAnti = (-iPartonAntiJun[0])/10; | |
525 | iAntiLegA.resize(0); | |
526 | iAntiLegB.resize(0); | |
527 | iAntiLegC.resize(0); | |
528 | leg = -1; | |
529 | for (int i = 0; i < int(iPartonAntiJun.size()); ++ i) { | |
530 | if ( (-iPartonAntiJun[i])/10 == identAnti) ++leg; | |
531 | if (leg == 0) iAntiLegA.push_back( iPartonAntiJun[i] ); | |
532 | else if (leg == 1) iAntiLegB.push_back( iPartonAntiJun[i] ); | |
533 | else iAntiLegC.push_back( iPartonAntiJun[i] ); | |
534 | } | |
535 | ||
536 | // Find interjunction legs, i.e. between junction and antijunction. | |
537 | int nMatch = 0; | |
538 | int legJun[3], legAnti[3], nGluLeg[3]; | |
539 | if (iJunLegA.back() < 0) { legJun[nMatch] = 0; | |
540 | legAnti[nMatch] = (-iJunLegA.back())%10; ++nMatch;} | |
541 | if (iJunLegB.back() < 0) { legJun[nMatch] = 1; | |
542 | legAnti[nMatch] = (-iJunLegB.back())%10; ++nMatch;} | |
543 | if (iJunLegC.back() < 0) { legJun[nMatch] = 2; | |
544 | legAnti[nMatch] = (-iJunLegC.back())%10; ++nMatch;} | |
545 | ||
546 | // Loop over interjunction legs. | |
547 | for (int iMatch = 0; iMatch < nMatch; ++iMatch) { | |
548 | vector<int>& iJunLeg = (legJun[iMatch] == 0) ? iJunLegA | |
549 | : ( (legJun[iMatch] == 1) ? iJunLegB : iJunLegC ); | |
550 | vector<int>& iAntiLeg = (legAnti[iMatch] == 0) ? iAntiLegA | |
551 | : ( (legAnti[iMatch] == 1) ? iAntiLegB : iAntiLegC ); | |
552 | ||
553 | // Find number of gluons on each. Do nothing for now if none. | |
554 | nGluLeg[iMatch] = iJunLeg.size() + iAntiLeg.size() - 4; | |
555 | if (nGluLeg[iMatch] == 0) continue; | |
556 | ||
557 | // Else pick up the gluons on the interjunction leg in order. | |
558 | iGluLeg.resize(0); | |
559 | for (int i = 1; i < int(iJunLeg.size()) - 1; ++i) | |
560 | iGluLeg.push_back( iJunLeg[i] ); | |
561 | for (int i = int(iAntiLeg.size()) - 2; i > 0; --i) | |
562 | iGluLeg.push_back( iAntiLeg[i] ); | |
563 | ||
564 | // Remove those gluons from the junction/antijunction leg lists. | |
565 | iJunLeg.resize(1); | |
566 | iAntiLeg.resize(1); | |
567 | ||
568 | // Pick a new quark at random; for simplicity no diquarks. | |
569 | int idQ = flavSel.pickLightQ(); | |
570 | int colQ, acolQ; | |
571 | ||
572 | // If one gluon on leg, split it into a collinear q-qbar pair. | |
573 | if (iGluLeg.size() == 1) { | |
574 | ||
575 | // Store the new q qbar pair, sharing gluon colour and momentum. | |
576 | colQ = event[ iGluLeg[0] ].col(); | |
577 | acolQ = event[ iGluLeg[0] ].acol(); | |
578 | Vec4 pQ = 0.5 * event[ iGluLeg[0] ].p(); | |
579 | double mQ = 0.5 * event[ iGluLeg[0] ].m(); | |
580 | int iQ = event.append( idQ, 75, iGluLeg[0], 0, 0, 0, colQ, 0, pQ, mQ ); | |
581 | int iQbar = event.append( -idQ, 75, iGluLeg[0], 0, 0, 0, 0, acolQ, | |
582 | pQ, mQ ); | |
583 | ||
584 | // Mark split gluon and update junction and antijunction legs. | |
585 | event[ iGluLeg[0] ].statusNeg(); | |
586 | event[ iGluLeg[0] ].daughters( iQ, iQbar); | |
587 | iJunLeg.push_back(iQ); | |
588 | iAntiLeg.push_back(iQbar); | |
589 | ||
590 | // If several gluons on the string, decide which g-g region to split up. | |
591 | } else { | |
592 | ||
593 | // Evaluate mass-squared for all adjacent gluon pairs. | |
594 | m2Pair.resize(0); | |
595 | double m2Sum = 0.; | |
596 | for (int i = 0; i < int(iGluLeg.size()) - 1; ++i) { | |
597 | double m2Now = 0.5 * event[ iGluLeg[i] ].p() | |
598 | * event[ iGluLeg[i + 1] ].p(); | |
599 | m2Pair.push_back(m2Now); | |
600 | m2Sum += m2Now; | |
601 | } | |
602 | ||
603 | // Pick breakup region with probability proportional to mass-squared. | |
604 | double m2Reg = m2Sum * rndmPtr->flat(); | |
605 | int iReg = -1; | |
606 | do m2Reg -= m2Pair[++iReg]; | |
607 | while (m2Reg > 0. && iReg < int(iGluLeg.size()) - 1); | |
608 | m2Reg = m2Pair[iReg]; | |
609 | ||
610 | // Pick breaking point of string in chosen region (symmetrically). | |
611 | double m2Temp = min( JJSTRINGM2MAX, JJSTRINGM2FRAC * m2Reg); | |
612 | double xPos = 0.5; | |
613 | double xNeg = 0.5; | |
614 | do { | |
615 | double zTemp = zSel.zFrag( idQ, 0, m2Temp); | |
616 | xPos = 1. - zTemp; | |
617 | xNeg = m2Temp / (zTemp * m2Reg); | |
618 | } while (xNeg > 1.); | |
619 | if (rndmPtr->flat() > 0.5) swap(xPos, xNeg); | |
620 | ||
621 | // Pick up two "mother" gluons of breakup. Mark them decayed. | |
622 | Particle& gJun = event[ iGluLeg[iReg] ]; | |
623 | Particle& gAnti = event[ iGluLeg[iReg + 1] ]; | |
624 | gJun.statusNeg(); | |
625 | gAnti.statusNeg(); | |
626 | int dau1 = event.size(); | |
627 | gJun.daughters(dau1, dau1 + 3); | |
628 | gAnti.daughters(dau1, dau1 + 3); | |
629 | int mother1 = min( iGluLeg[iReg], iGluLeg[iReg + 1]); | |
630 | int mother2 = max( iGluLeg[iReg], iGluLeg[iReg + 1]); | |
631 | ||
632 | // Can keep one of old colours but need one new so unambiguous. | |
633 | colQ = gJun.acol(); | |
634 | acolQ = event.nextColTag(); | |
635 | ||
636 | // Store copied gluons with reduced momenta. | |
637 | int iGjun = event.append( 21, 75, mother1, mother2, 0, 0, | |
638 | gJun.col(), gJun.acol(), (1. - 0.5 * xPos) * gJun.p(), | |
639 | (1. - 0.5 * xPos) * gJun.m()); | |
640 | int iGanti = event.append( 21, 75, mother1, mother2, 0, 0, | |
641 | acolQ, gAnti.acol(), (1. - 0.5 * xNeg) * gAnti.p(), | |
642 | (1. - 0.5 * xNeg) * gAnti.m()); | |
643 | ||
644 | // Store the new q qbar pair with remaining momenta. | |
645 | int iQ = event.append( idQ, 75, mother1, mother2, 0, 0, | |
646 | colQ, 0, 0.5 * xNeg * gAnti.p(), 0.5 * xNeg * gAnti.m() ); | |
647 | int iQbar = event.append( -idQ, 75, mother1, mother2, 0, 0, | |
648 | 0, acolQ, 0.5 * xPos * gJun.p(), 0.5 * xPos * gJun.m() ); | |
649 | ||
650 | // Update junction and antijunction legs with gluons and quarks. | |
651 | for (int i = 0; i < iReg; ++i) | |
652 | iJunLeg.push_back( iGluLeg[i] ); | |
653 | iJunLeg.push_back(iGjun); | |
654 | iJunLeg.push_back(iQ); | |
655 | for (int i = int(iGluLeg.size()) - 1; i > iReg + 1; --i) | |
656 | iAntiLeg.push_back( iGluLeg[i] ); | |
657 | iAntiLeg.push_back(iGanti); | |
658 | iAntiLeg.push_back(iQbar); | |
659 | } | |
660 | ||
661 | // Update end colours for both g -> q qbar and g g -> g g q qbar. | |
662 | event.endColJunction(identJun - 1, legJun[iMatch], colQ); | |
663 | event.endColJunction(identAnti - 1, legAnti[iMatch], acolQ); | |
664 | } | |
665 | ||
666 | // Update list of interjunction legs after splittings above. | |
667 | int iMatchUp = 0; | |
668 | while (iMatchUp < nMatch) { | |
669 | if (nGluLeg[iMatchUp] > 0) { | |
670 | for (int i = iMatchUp; i < nMatch - 1; ++i) { | |
671 | legJun[i] = legJun[i + 1]; | |
672 | legAnti[i] = legAnti[i + 1]; | |
673 | nGluLeg[i] = nGluLeg[i + 1]; | |
674 | } --nMatch; | |
675 | } else ++iMatchUp; | |
676 | } | |
677 | ||
678 | // Should not ever have three empty interjunction legs. | |
679 | if (nMatch == 3) { | |
680 | infoPtr->errorMsg("Error in HadronLevel::splitJunctionPair: " | |
681 | "three empty junction-junction legs"); | |
682 | return false; | |
683 | } | |
684 | ||
685 | // If two legs are empty, then collapse system to a single string. | |
686 | if (nMatch == 2) { | |
687 | int legJunLeft = 3 - legJun[0] - legJun[1]; | |
688 | int legAntiLeft = 3 - legAnti[0] - legAnti[1]; | |
689 | vector<int>& iJunLeg = (legJunLeft == 0) ? iJunLegA | |
690 | : ( (legJunLeft == 1) ? iJunLegB : iJunLegC ); | |
691 | vector<int>& iAntiLeg = (legAntiLeft == 0) ? iAntiLegA | |
692 | : ( (legAntiLeft == 1) ? iAntiLegB : iAntiLegC ); | |
693 | iPartonJun.resize(0); | |
694 | for (int i = int(iJunLeg.size()) - 1; i > 0; --i) | |
695 | iPartonJun.push_back( iJunLeg[i] ); | |
696 | for (int i = 1; i < int(iAntiLeg.size()); ++i) | |
697 | iPartonJun.push_back( iAntiLeg[i] ); | |
698 | ||
699 | // Match up the colours where the strings are joined. | |
700 | int iColJoin = iJunLeg[1]; | |
701 | int iAcolJoin = iAntiLeg[1]; | |
702 | event[iAcolJoin].acol( event[iColJoin].col() ); | |
703 | ||
704 | // Other string system empty. Remove junctions from their list. Done. | |
705 | iPartonAntiJun.resize(0); | |
706 | event.eraseJunction( max(identJun, identAnti) - 1); | |
707 | event.eraseJunction( min(identJun, identAnti) - 1); | |
708 | return true; | |
709 | } | |
710 | ||
711 | // If one leg is empty then, depending on string length, either | |
712 | // (a) annihilate junction and antijunction into two simple strings, or | |
713 | // (b) split the empty leg by borrowing energy from nearby legs. | |
714 | if (nMatch == 1) { | |
715 | ||
716 | // Identify the two external legs of either junction. | |
717 | vector<int>& iJunLeg0 = (legJun[0] == 0) ? iJunLegB : iJunLegA; | |
718 | vector<int>& iJunLeg1 = (legJun[0] == 2) ? iJunLegB : iJunLegC; | |
719 | vector<int>& iAntiLeg0 = (legAnti[0] == 0) ? iAntiLegB : iAntiLegA; | |
720 | vector<int>& iAntiLeg1 = (legAnti[0] == 2) ? iAntiLegB : iAntiLegC; | |
721 | ||
722 | // Simplified procedure: mainly study first parton on each leg. | |
723 | Vec4 pJunLeg0 = event[ iJunLeg0[1] ].p(); | |
724 | Vec4 pJunLeg1 = event[ iJunLeg1[1] ].p(); | |
725 | Vec4 pAntiLeg0 = event[ iAntiLeg0[1] ].p(); | |
726 | Vec4 pAntiLeg1 = event[ iAntiLeg1[1] ].p(); | |
727 | ||
728 | // Starting frame hopefully intermediate to two junction directions. | |
729 | Vec4 pStart = pJunLeg0 / pJunLeg0.e() + pJunLeg1 / pJunLeg1.e() | |
730 | + pAntiLeg0 / pAntiLeg0.e() + pAntiLeg1 / pAntiLeg1.e(); | |
731 | ||
732 | // Loop over iteration to junction/antijunction rest frames (JRF/ARF). | |
733 | RotBstMatrix MtoJRF, MtoARF; | |
734 | Vec4 pInJRF[3], pInARF[3]; | |
735 | for (int iJun = 0; iJun < 2; ++iJun) { | |
736 | int offset = (iJun == 0) ? 0 : 2; | |
737 | ||
738 | // Iterate from system rest frame towards the junction rest frame. | |
739 | RotBstMatrix MtoRF, Mstep; | |
740 | MtoRF.bstback(pStart); | |
741 | Vec4 pInRF[4]; | |
742 | int iter = 0; | |
743 | do { | |
744 | ++iter; | |
745 | ||
746 | // Find rest-frame momenta on the three sides of the junction. | |
747 | // Only consider first parton on each leg, for simplicity. | |
748 | pInRF[0 + offset] = pJunLeg0; | |
749 | pInRF[1 + offset] = pJunLeg1; | |
750 | pInRF[2 - offset] = pAntiLeg0; | |
751 | pInRF[3 - offset] = pAntiLeg1; | |
752 | for (int i = 0; i < 4; ++i) pInRF[i].rotbst(MtoRF); | |
753 | ||
754 | // For third side add both legs beyond other junction, weighted. | |
755 | double wt2 = 1. - exp( -pInRF[2].e() / eNormJunction); | |
756 | double wt3 = 1. - exp( -pInRF[3].e() / eNormJunction); | |
757 | pInRF[2] = wt2 * pInRF[2] + wt3 * pInRF[3]; | |
758 | ||
759 | // Find new junction rest frame from the set of momenta. | |
760 | Mstep = stringFrag.junctionRestFrame( pInRF[0], pInRF[1], pInRF[2]); | |
761 | MtoRF.rotbst( Mstep ); | |
762 | } while (iter < 3 || (Mstep.deviation() > CONVJNREST | |
763 | && iter < NTRYJNREST) ); | |
764 | ||
765 | // Store final boost and rest-frame (weighted) momenta. | |
766 | if (iJun == 0) { | |
767 | MtoJRF = MtoRF; | |
768 | for (int i = 0; i < 3; ++i) pInJRF[i] = pInRF[i]; | |
769 | } else { | |
770 | MtoARF = MtoRF; | |
771 | for (int i = 0; i < 3; ++i) pInARF[i] = pInRF[i]; | |
772 | } | |
773 | } | |
774 | ||
775 | // Opposite operations: boost from JRF/ARF to original system. | |
776 | RotBstMatrix MfromJRF = MtoJRF; | |
777 | MfromJRF.invert(); | |
778 | RotBstMatrix MfromARF = MtoARF; | |
779 | MfromARF.invert(); | |
780 | ||
781 | // Velocity vectors of junctions and momentum of legs in lab frame. | |
782 | Vec4 vJun(0., 0., 0., 1.); | |
783 | vJun.rotbst(MfromJRF); | |
784 | Vec4 vAnti(0., 0., 0., 1.); | |
785 | vAnti.rotbst(MfromARF); | |
786 | Vec4 pLabJ[3], pLabA[3]; | |
787 | for (int i = 0; i < 3; ++i) { | |
788 | pLabJ[i] = pInJRF[i]; | |
789 | pLabJ[i].rotbst(MfromJRF); | |
790 | pLabA[i] = pInARF[i]; | |
791 | pLabA[i].rotbst(MfromARF); | |
792 | } | |
793 | ||
794 | // Calculate Lambda-measure length of three possible topologies. | |
795 | double vJvA = vJun * vAnti; | |
796 | double vJvAe2y = vJvA + sqrt(vJvA*vJvA - 1.); | |
797 | double LambdaJA = (2. * pInJRF[0].e()) * (2. * pInJRF[1].e()) | |
798 | * (2. * pInARF[0].e()) * (2. * pInARF[1].e()) * vJvAe2y; | |
799 | double Lambda00 = (2. * pLabJ[0] * pLabA[0]) | |
800 | * (2. * pLabJ[1] * pLabA[1]); | |
801 | double Lambda01 = (2. * pLabJ[0] * pLabA[1]) | |
802 | * (2. * pLabJ[1] * pLabA[0]); | |
803 | ||
804 | // Case when either topology without junctions is the shorter one. | |
805 | if (LambdaJA > min( Lambda00, Lambda01)) { | |
806 | vector<int>& iAntiMatch0 = (Lambda00 < Lambda01) | |
807 | ? iAntiLeg0 : iAntiLeg1; | |
808 | vector<int>& iAntiMatch1 = (Lambda00 < Lambda01) | |
809 | ? iAntiLeg1 : iAntiLeg0; | |
810 | ||
811 | // Define two quark-antiquark strings. | |
812 | iPartonJun.resize(0); | |
813 | for (int i = int(iJunLeg0.size()) - 1; i > 0; --i) | |
814 | iPartonJun.push_back( iJunLeg0[i] ); | |
815 | for (int i = 1; i < int(iAntiMatch0.size()); ++i) | |
816 | iPartonJun.push_back( iAntiMatch0[i] ); | |
817 | iPartonAntiJun.resize(0); | |
818 | for (int i = int(iJunLeg1.size()) - 1; i > 0; --i) | |
819 | iPartonAntiJun.push_back( iJunLeg1[i] ); | |
820 | for (int i = 1; i < int(iAntiMatch1.size()); ++i) | |
821 | iPartonAntiJun.push_back( iAntiMatch1[i] ); | |
822 | ||
823 | // Match up the colours where the strings are joined. | |
824 | int iColJoin = iJunLeg0[1]; | |
825 | int iAcolJoin = iAntiMatch0[1]; | |
826 | event[iAcolJoin].acol( event[iColJoin].col() ); | |
827 | iColJoin = iJunLeg1[1]; | |
828 | iAcolJoin = iAntiMatch1[1]; | |
829 | event[iAcolJoin].acol( event[iColJoin].col() ); | |
830 | ||
831 | // Remove junctions from their list. Done. | |
832 | event.eraseJunction( max(identJun, identAnti) - 1); | |
833 | event.eraseJunction( min(identJun, identAnti) - 1); | |
834 | return true; | |
835 | } | |
836 | ||
837 | // Case where junction and antijunction to be separated. | |
838 | // Shuffle (p+/p-) momentum of order <mThad> between systems, | |
839 | // times 2/3 for 120 degree in JRF, times 1/2 for two legs, | |
840 | // but not more than half of what nearest parton carries. | |
841 | double eShift = MTHAD / (3. * sqrt(vJvAe2y)); | |
842 | double fracJ0 = min(0.5, eShift / pInJRF[0].e()); | |
843 | double fracJ1 = min(0.5, eShift / pInJRF[0].e()); | |
844 | Vec4 pFromJun = fracJ0 * pJunLeg0 + fracJ1 * pJunLeg1; | |
845 | double fracA0 = min(0.5, eShift / pInARF[0].e()); | |
846 | double fracA1 = min(0.5, eShift / pInARF[0].e()); | |
847 | Vec4 pFromAnti = fracA0 * pAntiLeg0 + fracA1 * pAntiLeg1; | |
848 | ||
849 | // Copy partons with scaled-down momenta and update legs. | |
850 | int iNew = event.copy(iJunLeg0[1], 76); | |
851 | event[iNew].rescale5(1. - fracJ0); | |
852 | iJunLeg0[1] = iNew; | |
853 | iNew = event.copy(iJunLeg1[1], 76); | |
854 | event[iNew].rescale5(1. - fracJ1); | |
855 | iJunLeg1[1] = iNew; | |
856 | iNew = event.copy(iAntiLeg0[1], 76); | |
857 | event[iNew].rescale5(1. - fracA0); | |
858 | iAntiLeg0[1] = iNew; | |
859 | iNew = event.copy(iAntiLeg1[1], 76); | |
860 | event[iNew].rescale5(1. - fracA1); | |
861 | iAntiLeg1[1] = iNew; | |
862 | ||
863 | // Pick a new quark at random; for simplicity no diquarks. | |
864 | int idQ = flavSel.pickLightQ(); | |
865 | ||
866 | // Update junction colours for new quark and antiquark. | |
867 | int colQ = event.nextColTag(); | |
868 | int acolQ = event.nextColTag(); | |
869 | event.endColJunction(identJun - 1, legJun[0], colQ); | |
870 | event.endColJunction(identAnti - 1, legAnti[0], acolQ); | |
871 | ||
872 | // Store new quark and antiquark with momentum from other junction. | |
873 | int mother1 = min(iJunLeg0[1], iJunLeg1[1]); | |
874 | int mother2 = max(iJunLeg0[1], iJunLeg1[1]); | |
875 | int iNewJ = event.append( idQ, 76, mother1, mother2, 0, 0, | |
876 | colQ, 0, pFromAnti, pFromAnti.mCalc() ); | |
877 | mother1 = min(iAntiLeg0[1], iAntiLeg1[1]); | |
878 | mother2 = max(iAntiLeg0[1], iAntiLeg1[1]); | |
879 | int iNewA = event.append( -idQ, 76, mother1, mother2, 0, 0, | |
880 | 0, acolQ, pFromJun, pFromJun.mCalc() ); | |
881 | ||
882 | // Bookkeep new quark and antiquark on third legs. | |
883 | if (legJun[0] == 0) iJunLegA[1] = iNewJ; | |
884 | else if (legJun[0] == 1) iJunLegB[1] = iNewJ; | |
885 | else iJunLegC[1] = iNewJ; | |
886 | if (legAnti[0] == 0) iAntiLegA[1] = iNewA; | |
887 | else if (legAnti[0] == 1) iAntiLegB[1] = iNewA; | |
888 | else iAntiLegC[1] = iNewA; | |
889 | ||
890 | // Done with splitting junction from antijunction. | |
891 | } | |
892 | ||
893 | // Put together new junction parton list. | |
894 | iPartonJun.resize(0); | |
895 | for (int i = 0; i < int(iJunLegA.size()); ++i) | |
896 | iPartonJun.push_back( iJunLegA[i] ); | |
897 | for (int i = 0; i < int(iJunLegB.size()); ++i) | |
898 | iPartonJun.push_back( iJunLegB[i] ); | |
899 | for (int i = 0; i < int(iJunLegC.size()); ++i) | |
900 | iPartonJun.push_back( iJunLegC[i] ); | |
901 | ||
902 | // Put together new antijunction parton list. | |
903 | iPartonAntiJun.resize(0); | |
904 | for (int i = 0; i < int(iAntiLegA.size()); ++i) | |
905 | iPartonAntiJun.push_back( iAntiLegA[i] ); | |
906 | for (int i = 0; i < int(iAntiLegB.size()); ++i) | |
907 | iPartonAntiJun.push_back( iAntiLegB[i] ); | |
908 | for (int i = 0; i < int(iAntiLegC.size()); ++i) | |
909 | iPartonAntiJun.push_back( iAntiLegC[i] ); | |
910 | ||
911 | // Now the two junction systems are separated and can be stored. | |
912 | return true; | |
913 | ||
914 | } | |
915 | ||
916 | //========================================================================== | |
917 | ||
918 | } // end namespace Pythia8 | |
919 |