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