StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StringFragmentation.cc
1 // StringFragmentation.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, 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 StringEnd and
7 // StringFragmentation classes.
8 
9 #include "Pythia8/StringFragmentation.h"
10 
11 namespace Pythia8 {
12 
13 //==========================================================================
14 
15 // The StringEnd class.
16 
17 //--------------------------------------------------------------------------
18 
19 // Constants: could be changed here if desired, but normally should not.
20 
21 // Avoid unphysical solutions to equation system.
22 const double StringEnd::TINY = 1e-6;
23 
24 // Assume two (eX, eY) regions are related if pT2 differs by less.
25 const double StringEnd::PT2SAME = 0.01;
26 
27 // Fictitious typical mass and pT used to look ahead for approximately where
28 // the next hadron will be produced, to quantify string close-packing there.
29 const double StringEnd::MEANMMIN = 0.2;
30 const double StringEnd::MEANM = 0.8;
31 const double StringEnd::MEANPT = 0.4;
32 
33 //--------------------------------------------------------------------------
34 
35 // Set up initial endpoint values from input.
36 
37 void StringEnd::setUp(bool fromPosIn, int iEndIn, int idOldIn, int iMaxIn,
38  double pxIn, double pyIn, double GammaIn, double xPosIn, double xNegIn,
39  int colIn) {
40 
41  // Simple transcription from input.
42  fromPos = fromPosIn;
43  iEnd = iEndIn;
44  iMax = iMaxIn;
45  flavOld = FlavContainer(idOldIn);
46  pxOld = pxIn;
47  pyOld = pyIn;
48  GammaOld = GammaIn;
49  iPosOld = (fromPos) ? 0 : iMax;
50  iNegOld = (fromPos) ? iMax : 0;
51  xPosOld = xPosIn;
52  xNegOld = xNegIn;
53  colOld = colIn;
54 
55 }
56 
57 //--------------------------------------------------------------------------
58 
59 // Fragment off one hadron from the string system, in flavour and pT.
60 
61 void StringEnd::newHadron(double nNSP) {
62 
63  // In case we are using the thermal model or Gaussian with
64  // mT2 suppression we have to pick the pT first.
65  if (thermalModel || mT2suppression) {
66 
67  // Pick its transverse momentum.
68  pair<double, double> pxy = pTSelPtr->pxy(flavNew.id, nNSP);
69  pxNew = pxy.first;
70  pyNew = pxy.second;
71  pxHad = pxOld + pxNew;
72  pyHad = pyOld + pyNew;
73  double pT2Had = pow2(pxHad) + pow2(pyHad);
74 
75  // Pick new flavour and form a new hadron.
76  do {
77  flavNew = flavSelPtr->pick( flavOld, sqrt(pT2Had), nNSP);
78  idHad = flavSelPtr->getHadronID( flavOld, flavNew);
79  } while (idHad == 0);
80 
81  // Get its mass and thereby define its transverse mass.
82  mHad = flavSelPtr->getHadronMassWin(idHad);
83  mT2Had = pow2(mHad) + pow2(pxHad) + pow2(pyHad);
84  }
85 
86  // In case of the Gaussian without mT2 suppression we pick
87  // the new flavour first to make the width flavour dependent.
88  else {
89 
90  // Pick new flavour and form a new hadron.
91  do {
92  flavNew = flavSelPtr->pick( flavOld);
93  idHad = flavSelPtr->combine( flavOld, flavNew);
94  } while (idHad == 0);
95 
96  // Pick its transverse momentum.
97  pair<double, double> pxy = pTSelPtr->pxy(flavNew.id, nNSP);
98  pxNew = pxy.first;
99  pyNew = pxy.second;
100  pxHad = pxOld + pxNew;
101  pyHad = pyOld + pyNew;
102 
103  // Pick its mass and thereby define its transverse mass.
104  mHad = particleDataPtr->mSel(idHad);
105  mT2Had = pow2(mHad) + pow2(pxHad) + pow2(pyHad);
106  }
107 
108 }
109 
110 //--------------------------------------------------------------------------
111 
112 // Fragment off one hadron from the string system, in momentum space,
113 // by taking steps from positive end.
114 
115 Vec4 StringEnd::kinematicsHadron( StringSystem& system,
116  vector<StringVertex>& stringVertices, bool useInputZ, double zHadIn) {
117 
118  // Pick fragmentation step z and calculate new Gamma.
119  if (useInputZ) zHad = zHadIn;
120  else zHad = zSelPtr->zFrag( flavOld.id, flavNew.id, mT2Had);
121  GammaNew = (1. - zHad) * (GammaOld + mT2Had / zHad);
122 
123  // Set up references that are direction-neutral;
124  // ...Dir for direction of iteration and ...Inv for its inverse.
125  int& iDirOld = (fromPos) ? iPosOld : iNegOld;
126  int& iInvOld = (fromPos) ? iNegOld : iPosOld;
127  int& iDirNew = (fromPos) ? iPosNew : iNegNew;
128  int& iInvNew = (fromPos) ? iNegNew : iPosNew;
129  double& xDirOld = (fromPos) ? xPosOld : xNegOld;
130  double& xInvOld = (fromPos) ? xNegOld : xPosOld;
131  double& xDirNew = (fromPos) ? xPosNew : xNegNew;
132  double& xInvNew = (fromPos) ? xNegNew : xPosNew;
133  double& xDirHad = (fromPos) ? xPosHad : xNegHad;
134  double& xInvHad = (fromPos) ? xNegHad : xPosHad;
135 
136  // Start search for new breakup in the old region.
137  iDirNew = iDirOld;
138  iInvNew = iInvOld;
139  Vec4 pTNew;
140 
141  // Each step corresponds to trying a new string region.
142  for (int iStep = 0; ; ++iStep) {
143 
144  // Referance to current string region.
145  StringRegion& region = system.region( iPosNew, iNegNew);
146  colNew = fromPos ? region.colPos : region.colNeg;
147 
148  // Now begin special section for rapid processing of low region.
149  if (iStep == 0 && iPosOld + iNegOld == iMax) {
150 
151  // A first step within a low region is easy.
152  if (mT2Had < zHad * xDirOld * (1. - xInvOld) * region.w2) {
153 
154  // Translate into x coordinates.
155  xDirHad = zHad * xDirOld;
156  xInvHad = mT2Had / (xDirHad * region.w2);
157  xDirNew = xDirOld - xDirHad;
158  xInvNew = xInvOld + xInvHad;
159 
160  // Store breakup vertex information from the fragmentation process.
161  stringVertices.push_back( StringVertex( fromPos, iPosNew, iNegNew,
162  xPosNew, xNegNew) );
163 
164  // Find and return four-momentum of the produced particle.
165  return region.pHad( xPosHad, xNegHad, pxHad, pyHad);
166 
167  // A first step out of a low region also OK, if there are more regions.
168  // Negative energy signals failure, i.e. in last region.
169  } else {
170  --iInvNew;
171  if (iInvNew < 0) return Vec4(0., 0., 0., -1.);
172 
173  // Momentum taken by stepping out of region. Continue to next region.
174  xInvHad = 1. - xInvOld;
175  xDirHad = 0.;
176  pSoFar = region.pHad( xPosHad, xNegHad, pxOld, pyOld);
177  continue;
178  }
179 
180  // Else, for first step, take into account starting pT.
181  } else if (iStep == 0) {
182  pSoFar = region.pHad( 0., 0., pxOld, pyOld);
183  pTNew = region.pHad( 0., 0., pxNew, pyNew);
184  }
185 
186 
187  // Now begin normal treatment of nontrivial regions.
188  // Set up four-vectors in a region not visited before.
189  if (!region.isSetUp) region.setUp(
190  system.regionLowPos(iPosNew).pPos,
191  system.regionLowNeg(iNegNew).pNeg,
192  system.regionLowPos(iPosNew).colPos,
193  system.regionLowNeg(iNegNew).colNeg, true);
194 
195  // If new region is vanishingly small, continue immediately to next.
196  // Negative energy signals failure to do this, i.e. moved too low.
197  if (region.isEmpty) {
198  xDirHad = (iDirNew == iDirOld) ? xDirOld : 1.;
199  xInvHad = 0.;
200  pSoFar += region.pHad( xPosHad, xNegHad, 0., 0.);
201  ++iDirNew;
202  if (iDirNew + iInvNew > iMax) return Vec4(0., 0., 0., -1.);
203  continue;
204  }
205 
206  // Reexpress pTNew w.r.t. base vectors in new region, if possible.
207  // Recall minus sign from normalization e_x^2 = e_y^2 = -1.
208  double pxNewTemp = -pTNew * region.eX;
209  double pyNewTemp = -pTNew * region.eY;
210  if (abs( pxNewTemp * pxNewTemp + pyNewTemp * pyNewTemp
211  - pxNew * pxNew - pyNew * pyNew) < PT2SAME) {
212  pxNew = pxNewTemp;
213  pyNew = pyNewTemp;
214  }
215 
216  // Four-momentum taken so far, including new pT.
217  Vec4 pTemp = pSoFar + region.pHad( 0., 0., pxNew, pyNew);
218 
219  // Derive coefficients for m2 expression.
220  // cM2 * x+ + cM3 * x- + cM4 * x+ * x- = m^2 - cM1;
221  double cM1 = pTemp.m2Calc();
222  double cM2 = 2. * (pTemp * region.pPos);
223  double cM3 = 2. * (pTemp * region.pNeg);
224  double cM4 = region.w2;
225  if (!fromPos) swap( cM2, cM3);
226 
227  // Derive coefficients for Gamma expression.
228  // cGam2 * x+ + cGam3 * x- + cGam4 * x+ * x- = Gamma_new - cGam1;
229  double cGam1 = 0.;
230  double cGam2 = 0.;
231  double cGam3 = 0.;
232  double cGam4 = 0.;
233  for (int iInv = iInvNew; iInv <= iMax - iDirNew; ++iInv) {
234  double xInv = 1.;
235  if (iInv == iInvNew) xInv = (iInvNew == iInvOld) ? xInvOld : 0.;
236  for (int iDir = iDirNew; iDir <= iMax - iInv; ++iDir) {
237  double xDir = (iDir == iDirOld) ? xDirOld : 1.;
238  int iPos = (fromPos) ? iDir : iInv;
239  int iNeg = (fromPos) ? iInv : iDir;
240  StringRegion& regionGam = system.region( iPos, iNeg);
241  if (!regionGam.isSetUp) regionGam.setUp(
242  system.regionLowPos(iPos).pPos,
243  system.regionLowNeg(iNeg).pNeg,
244  system.regionLowPos(iPos).colPos,
245  system.regionLowNeg(iNeg).colNeg, true);
246  double w2 = regionGam.w2;
247  cGam1 += xDir * xInv * w2;
248  if (iDir == iDirNew) cGam2 -= xInv * w2;
249  if (iInv == iInvNew) cGam3 += xDir * w2;
250  if (iDir == iDirNew && iInv == iInvNew) cGam4 -= w2;
251  }
252  }
253 
254  // Solve (m2, Gamma) equation system => r2 * x-^2 + r1 * x- + r0 = 0.
255  double cM0 = pow2(mHad) - cM1;
256  double cGam0 = GammaNew - cGam1;
257  double r2 = cM3 * cGam4 - cM4 * cGam3;
258  double r1 = cM4 * cGam0 - cM0 * cGam4 + cM3 * cGam2 - cM2 * cGam3;
259  double r0 = cM2 * cGam0 - cM0 * cGam2;
260  double root = sqrtpos( r1*r1 - 4. * r2 * r0 );
261  if (abs(r2) < TINY || root < TINY) return Vec4(0., 0., 0., -1.);
262  xInvHad = 0.5 * (root / abs(r2) - r1 / r2);
263  if (abs(cM2 + cM4 * xInvHad) < TINY) return Vec4(0., 0., 0., -1.);
264  xDirHad = (cM0 - cM3 * xInvHad) / (cM2 + cM4 * xInvHad);
265 
266  // Define position of new trial vertex.
267  xDirNew = (iDirNew == iDirOld) ? xDirOld - xDirHad : 1. - xDirHad;
268  xInvNew = (iInvNew == iInvOld) ? xInvOld + xInvHad : xInvHad;
269 
270  // Step up to new region if new x- > 1.
271  if (xInvNew > 1.) {
272  xInvHad = (iInvNew == iInvOld) ? 1. - xInvOld : 1.;
273  xDirHad = 0.;
274  pSoFar += region.pHad( xPosHad, xNegHad, 0., 0.);
275  --iInvNew;
276  if (iInvNew < 0) return Vec4(0., 0., 0., -1.);
277  continue;
278 
279  // Step down to new region if new x+ < 0.
280  } else if (xDirNew < 0.) {
281  xDirHad = (iDirNew == iDirOld) ? xDirOld : 1.;
282  xInvHad = 0.;
283  pSoFar += region.pHad( xPosHad, xNegHad, 0., 0.);
284  ++iDirNew;
285  if (iDirNew + iInvNew > iMax) return Vec4(0., 0., 0., -1.);
286  continue;
287  }
288 
289  // Store breakup vertex information from the fragmentation process.
290  stringVertices.push_back( StringVertex( fromPos, iPosNew, iNegNew,
291  xPosNew, xNegNew) );
292 
293  // Else we have found the correct region, and can set the new
294  // colour index and return four-momentum.
295  colNew = fromPos ? region.colPos : region.colNeg;
296  if (useInputZ) return Vec4( 0., 0., 0., 0.);
297  else return pSoFar + region.pHad( xPosHad, xNegHad, pxNew, pyNew);
298 
299  // End of "infinite" loop of stepping to new region.
300  }
301 
302 }
303 
304 //--------------------------------------------------------------------------
305 
306 // Generate momentum for some possible next hadron, based on mean values
307 // to get an estimate for rapidity and pT.
308 
309 Vec4 StringEnd::kinematicsHadronTmp( StringSystem system, Vec4 pRem,
310  double phi, double mult) {
311 
312  // Now estimate the energy the next hadron will take.
313  double mRem = pRem.mCalc();
314  double meanM = (mRem > 0.0) ? max( MEANMMIN, min( MEANM, mRem) ) : MEANM;
315  double meanMT2 = pow2(meanM) + pow2(MEANPT);
316  double GammaNow = (1.0 + aLund) / bLund;
317  // Modify Gamma value in case of earlier fails.
318  if (mult > 0.0) GammaNow *= mult;
319  double tmp = ( GammaNow + meanMT2 - GammaOld ) / GammaOld;
320  double zPlus = (-0.5 * tmp + sqrt(0.25 * pow2(tmp) + meanMT2 / GammaOld));
321  double zMinus = (-0.5 * tmp - sqrt(0.25 * pow2(tmp) + meanMT2 / GammaOld));
322  // Special case of first hadron.
323  if (GammaOld < 1e-10) {
324  zPlus = pow(1.0 + meanMT2 / GammaNow, -1.0);
325  zMinus = -1.0;
326  }
327  bool zPlusOk = (zPlus < 1.0) && (zPlus > 0.0);
328  bool zMinusOk = (zMinus < 1.0) && (zMinus > 0.0);
329  // Negative energy signals failure.
330  if ( (!zPlusOk) && (!zMinusOk) ) return Vec4(0., 0., 0., -1.);
331  double zHadTmp = (zPlusOk ? zPlus : zMinus);
332 
333  double pxHadTmp = cos(phi) * MEANPT;
334  double pyHadTmp = sin(phi) * MEANPT;
335 
336  // First make a copy of all variables to not overwrite anything.
337  int iPosOldTmp = iPosOld, iNegOldTmp = iNegOld;
338  int iPosNewTmp = iPosNew, iNegNewTmp = iNegNew;
339  double xPosOldTmp = xPosOld, xNegOldTmp = xNegOld;
340  double xPosNewTmp = xPosNew, xNegNewTmp = xNegNew;
341  double xPosHadTmp = xPosHad, xNegHadTmp = xNegHad;
342  double pxNewTmp = pxNew, pxOldTmp = pxOld;
343  double pyNewTmp = pyNew, pyOldTmp = pyOld;
344 
345  Vec4 pSoFarTmp = pSoFar;
346 
347  // Set up references that are direction-neutral;
348  // ...Dir for direction of iteration and ...Inv for its inverse.
349  int& iDirOld = (fromPos) ? iPosOldTmp : iNegOldTmp;
350  int& iInvOld = (fromPos) ? iNegOldTmp : iPosOldTmp;
351  int& iDirNew = (fromPos) ? iPosNewTmp : iNegNewTmp;
352  int& iInvNew = (fromPos) ? iNegNewTmp : iPosNewTmp;
353  double& xDirOld = (fromPos) ? xPosOldTmp : xNegOldTmp;
354  double& xInvOld = (fromPos) ? xNegOldTmp : xPosOldTmp;
355  double& xDirNew = (fromPos) ? xPosNewTmp : xNegNewTmp;
356  double& xInvNew = (fromPos) ? xNegNewTmp : xPosNewTmp;
357  double& xDirHad = (fromPos) ? xPosHadTmp : xNegHadTmp;
358  double& xInvHad = (fromPos) ? xNegHadTmp : xPosHadTmp;
359 
360  // Start search for new breakup in the old region.
361  iDirNew = iDirOld;
362  iInvNew = iInvOld;
363  Vec4 pTNew;
364 
365  // Each step corresponds to trying a new string region.
366  for (int iStep = 0; ; ++iStep) {
367 
368  // Referance to current string region.
369  StringRegion region = system.region( iPosNewTmp, iNegNewTmp);
370 
371  // Now begin special section for rapid processing of low region.
372  if (iStep == 0 && iPosOldTmp + iNegOldTmp == iMax) {
373 
374  // A first step within a low region is easy. Make sure we use this
375  // region in case it's the last one.
376  if ( (meanMT2 < zHadTmp * xDirOld * (1. - xInvOld) * region.w2)
377  || (iInvNew - 1 < 0) ) {
378 
379  if (iInvNew - 1 < 0)
380  zHadTmp = meanMT2 / xDirOld / (1. - xInvOld) / region.w2;
381 
382  // Translate into x coordinates.
383  xDirHad = zHad * xDirOld;
384  xInvHad = meanMT2 / (xDirHad * region.w2);
385  xDirNew = xDirOld - xDirHad;
386  xInvNew = xInvOld + xInvHad;
387 
388  // Find and return four-momentum of the produced particle.
389  return region.pHad( xPosHadTmp, xNegHadTmp, pxHadTmp, pyHadTmp);
390 
391  // A first step out of a low region also OK, if there are more regions.
392  // Negative energy signals failure, i.e. in last region.
393  } else {
394  --iInvNew;
395  // Should be covered by the above check.
396  if (iInvNew < 0) return Vec4(0., 0., 0., -1.);
397 
398  // Momentum taken by stepping out of region. Continue to next region.
399  xInvHad = 1. - xInvOld;
400  xDirHad = 0.;
401  pSoFarTmp = region.pHad( xPosHadTmp, xNegHadTmp, pxOldTmp, pyOldTmp);
402  continue;
403  }
404 
405  // Else, for first step, take into account starting pT.
406  } else if (iStep == 0) {
407  pSoFarTmp = region.pHad( 0., 0., pxOldTmp, pyOldTmp);
408  pTNew = region.pHad( 0., 0., pxNewTmp, pyNewTmp);
409  }
410 
411  // Now begin normal treatment of nontrivial regions.
412  // Set up four-vectors in a region not visited before.
413  if (!region.isSetUp) region.setUp(
414  system.regionLowPos(iPosNewTmp).pPos,
415  system.regionLowNeg(iNegNewTmp).pNeg,
416  system.regionLowPos(iPosNewTmp).colPos,
417  system.regionLowNeg(iNegNewTmp).colNeg, true);
418 
419  // If new region is vanishingly small, continue immediately to next.
420  // Negative energy signals failure to do this, i.e. moved too low.
421  if (region.isEmpty) {
422  xDirHad = (iDirNew == iDirOld) ? xDirOld : 1.;
423  xInvHad = 0.;
424  pSoFarTmp += region.pHad( xPosHadTmp, xNegHadTmp, 0., 0.);
425  ++iDirNew;
426  if (iDirNew + iInvNew > iMax) return Vec4(0., 0., 0., -1.);
427  continue;
428  }
429 
430  // Reexpress pTNew w.r.t. base vectors in new region, if possible.
431  // Recall minus sign from normalization e_x^2 = e_y^2 = -1.
432  double pxNewRegNow = -pTNew * region.eX;
433  double pyNewRegNow = -pTNew * region.eY;
434  if (abs( pxNewRegNow * pxNewRegNow + pyNewRegNow * pyNewRegNow
435  - pxNewTmp * pxNewTmp - pyNewTmp * pyNewTmp) < PT2SAME) {
436  pxNewTmp = pxNewRegNow;
437  pyNewTmp = pyNewRegNow;
438  }
439 
440  // Four-momentum taken so far, including new pT.
441  Vec4 pTemp = pSoFarTmp + region.pHad( 0., 0., pxNewTmp, pyNewTmp);
442 
443  // Derive coefficients for m2 expression.
444  // cM2 * x+ + cM3 * x- + cM4 * x+ * x- = m^2 - cM1;
445  double cM1 = pTemp.m2Calc();
446  double cM2 = 2. * (pTemp * region.pPos);
447  double cM3 = 2. * (pTemp * region.pNeg);
448  double cM4 = region.w2;
449  if (!fromPos) swap( cM2, cM3);
450 
451  // Derive coefficients for Gamma expression.
452  // cGam2 * x+ + cGam3 * x- + cGam4 * x+ * x- = Gamma_new - cGam1;
453  double cGam1 = 0.;
454  double cGam2 = 0.;
455  double cGam3 = 0.;
456  double cGam4 = 0.;
457  for (int iInv = iInvNew; iInv <= iMax - iDirNew; ++iInv) {
458  double xInv = 1.;
459  if (iInv == iInvNew) xInv = (iInvNew == iInvOld) ? xInvOld : 0.;
460  for (int iDir = iDirNew; iDir <= iMax - iInv; ++iDir) {
461  double xDir = (iDir == iDirOld) ? xDirOld : 1.;
462  int iPos = (fromPos) ? iDir : iInv;
463  int iNeg = (fromPos) ? iInv : iDir;
464  StringRegion regionGam = system.region( iPos, iNeg);
465  if (!regionGam.isSetUp) regionGam.setUp(
466  system.regionLowPos(iPos).pPos,
467  system.regionLowNeg(iNeg).pNeg,
468  system.regionLowPos(iPos).colPos,
469  system.regionLowNeg(iNeg).colNeg, true);
470  double w2 = regionGam.w2;
471  cGam1 += xDir * xInv * w2;
472  if (iDir == iDirNew) cGam2 -= xInv * w2;
473  if (iInv == iInvNew) cGam3 += xDir * w2;
474  if (iDir == iDirNew && iInv == iInvNew) cGam4 -= w2;
475  }
476  }
477 
478  // Solve (m2, Gamma) equation system => r2 * x-^2 + r1 * x- + r0 = 0.
479  double cM0 = pow2(meanM) - cM1;
480  double cGam0 = GammaNow - cGam1;
481  double r2 = cM3 * cGam4 - cM4 * cGam3;
482  double r1 = cM4 * cGam0 - cM0 * cGam4 + cM3 * cGam2 - cM2 * cGam3;
483  double r0 = cM2 * cGam0 - cM0 * cGam2;
484  double root = sqrtpos( r1*r1 - 4. * r2 * r0 );
485  if (abs(r2) < TINY || root < TINY) return Vec4(0., 0., 0., -1.);
486  xInvHad = 0.5 * (root / abs(r2) - r1 / r2);
487  xDirHad = (cM0 - cM3 * xInvHad) / (cM2 + cM4 * xInvHad);
488 
489  // Define position of new trial vertex.
490  xDirNew = (iDirNew == iDirOld) ? xDirOld - xDirHad : 1. - xDirHad;
491  xInvNew = (iInvNew == iInvOld) ? xInvOld + xInvHad : xInvHad;
492 
493  // Step up to new region if new x- > 1.
494  if (xInvNew > 1.) {
495  xInvHad = (iInvNew == iInvOld) ? 1. - xInvOld : 1.;
496  xDirHad = 0.;
497  pSoFarTmp += region.pHad( xPosHadTmp, xNegHadTmp, 0., 0.);
498  --iInvNew;
499  if (iInvNew < 0) return Vec4(0., 0., 0., -1.);
500  continue;
501 
502  // Step down to new region if new x+ < 0.
503  } else if (xDirNew < 0.) {
504  xDirHad = (iDirNew == iDirOld) ? xDirOld : 1.;
505  xInvHad = 0.;
506  pSoFarTmp += region.pHad( xPosHadTmp, xNegHadTmp, 0., 0.);
507  ++iDirNew;
508  if (iDirNew + iInvNew > iMax) return Vec4(0., 0., 0., -1.);
509  continue;
510  }
511 
512  // Else we have found the correct region, and can return four-momentum.
513  return pSoFarTmp + region.pHad( xPosHadTmp, xNegHadTmp, pxNewTmp,
514  pyNewTmp);
515 
516  // End of "infinite" loop of stepping to new region.
517  }
518 
519 }
520 
521 //--------------------------------------------------------------------------
522 
523 // Update string end information after a hadron has been removed.
524 
525 void StringEnd::update() {
526 
527  flavOld.anti(flavNew);
528  iPosOld = iPosNew;
529  iNegOld = iNegNew;
530  pxOld = -pxNew;
531  pyOld = -pyNew;
532  GammaOld = GammaNew;
533  xPosOld = xPosNew;
534  xNegOld = xNegNew;
535  colOld = colNew;
536 
537 }
538 
539 //==========================================================================
540 
541 // The StringFragmentation class.
542 
543 //--------------------------------------------------------------------------
544 
545 // Constants: could be changed here if desired, but normally should not.
546 // These are of technical nature, as described for each.
547 
548 // Maximum number of tries to (flavour-, energy) join the two string ends.
549 const int StringFragmentation::NTRYFLAV = 10;
550 const int StringFragmentation::NTRYJOIN = 30;
551 
552 // The last few times gradually increase the stop mass to make it easier.
553 const int StringFragmentation::NSTOPMASS = 15;
554 const double StringFragmentation::FACSTOPMASS = 1.05;
555 
556 // For closed string, pick a Gamma by taking a step with fictitious mass.
557 const double StringFragmentation::CLOSEDM2MAX = 25.;
558 const double StringFragmentation::CLOSEDM2FRAC = 0.1;
559 
560 // Do not allow too large argument to exp function.
561 const double StringFragmentation::EXPMAX = 50.;
562 
563 // Matching criterion that p+ and p- not the same (can happen in gg loop).
564 const double StringFragmentation::MATCHPOSNEG = 1e-4;
565 
566 // For pull on junction, do not trace too far down each leg.
567 const double StringFragmentation::EJNWEIGHTMAX = 10.;
568 
569 // Iterate junction rest frame boost until convergence or too many tries.
570 const double StringFragmentation::CONVJNREST = 1e-5;
571 const int StringFragmentation::NTRYJNREST = 20;
572 
573 // Fail and try again when two legs combined to diquark (3 loops).
574 const int StringFragmentation::NTRYJNMATCH = 20;
575 const double StringFragmentation::EEXTRAJNMATCH = 0.5;
576 const double StringFragmentation::MDIQUARKMIN = -2.0;
577 
578 // Consider junction-leg parton as massless if m2 tiny.
579 const double StringFragmentation::M2MAXJRF = 1e-4;
580 
581 // Protect against numerical precision giving zero or negative m2.
582 const double StringFragmentation::M2MINJRF = 1e-4;
583 
584 // Iterate junction rest frame equation until convergence or too many tries.
585 const double StringFragmentation::CONVJRFEQ = 1e-12;
586 const int StringFragmentation::NTRYJRFEQ = 40;
587 
588 // Retry smearing of breakup vertex if too big shifts.
589 const int StringFragmentation::NTRYSMEAR = 100;
590 
591 // Check that breakup vertex does not have negative tau^2 or t within roundoff.
592 const double StringFragmentation::CHECKPOS = 1e-10;
593 
594 //--------------------------------------------------------------------------
595 
596 // Initialize and save pointers.
597 
598 void StringFragmentation::init(StringFlav* flavSelPtrIn,
599  StringPT* pTSelPtrIn, StringZ* zSelPtrIn, FragModPtr fragModPtrIn) {
600 
601  // Save pointers.
602  flavSelPtr = flavSelPtrIn;
603  pTSelPtr = pTSelPtrIn;
604  zSelPtr = zSelPtrIn;
605  flavRopePtr = fragModPtrIn;
606 
607  // Initialize the StringFragmentation class.
608  stopMass = zSelPtr->stopMass();
609  stopNewFlav = zSelPtr->stopNewFlav();
610  stopSmear = zSelPtr->stopSmear();
611  eNormJunction = parm("StringFragmentation:eNormJunction");
612  eBothLeftJunction = parm("StringFragmentation:eBothLeftJunction");
613  eMaxLeftJunction = parm("StringFragmentation:eMaxLeftJunction");
614  eMinLeftJunction = parm("StringFragmentation:eMinLeftJunction");
615 
616  // Calculation and definition of hadron space-time production vertices.
617  hadronVertex = mode("HadronVertex:mode");
618  setVertices = flag("Fragmentation:setVertices")
619  || flag("HadronLevel:Rescatter");
620  kappaVtx = parm("HadronVertex:kappa");
621  smearOn = flag("HadronVertex:smearOn");
622  xySmear = parm("HadronVertex:xySmear");
623  maxSmear = parm("HadronVertex:maxSmear");
624  constantTau = flag("HadronVertex:constantTau");
625  maxTau = parm("HadronVertex:maxTau");
626 
627  // Tracing of colours for primary hadrons.
628  traceColours = flag("StringFragmentation:TraceColours");
629 
630  // Joining of nearby partons along the string.
631  mJoin = parm("FragmentationSystems:mJoin");
632 
633  // Initialize the b parameter of the z spectrum, used when joining jets.
634  bLund = zSelPtr->bAreaLund();
635 
636  // Charm and bottom quark masses used for space-time offset.
637  mc = particleDataPtr->m0(4);
638  mb = particleDataPtr->m0(5);
639 
640  // MPI pT0, used for calculating effective number of strings.
641  pT20 = pow2(parm("MultipartonInteractions:pT0Ref"));
642 
643  // Initialize the hadrons instance of an event record.
644  hadrons.init( "(string fragmentation)", particleDataPtr);
645 
646  // Send on pointers to the two StringEnd instances.
647  posEnd.init( particleDataPtr, flavSelPtr, pTSelPtr, zSelPtr, *settingsPtr);
648  negEnd.init( particleDataPtr, flavSelPtr, pTSelPtr, zSelPtr, *settingsPtr);
649 
650  // Check for number of nearby string pieces (nNSP) or not.
651  closePacking = flag("StringPT:closePacking");
652 
653 }
654 
655 //--------------------------------------------------------------------------
656 
657 // Perform the fragmentation.
658 
659 bool StringFragmentation::fragment( int iSub, ColConfig& colConfig,
660  Event& event) {
661 
662  // Find partons and their total four-momentum.
663  iParton = colConfig[iSub].iParton;
664  iPos = iParton[0];
665  if (iPos < 0) iPos = iParton[1];
666  int idPos = event[iPos].id();
667  iNeg = iParton.back();
668  int idNeg = event[iNeg].id();
669  pSum = colConfig[iSub].pSum;
670 
671  // Rapidity pairs [yMin, yMax] of string piece ends.
672  vector< vector< pair<double,double> > > rapPairs = colConfig.rapPairs;
673 
674  // Reset the local event record and vertex arrays.
675  hadrons.clear();
676  stringVertices.clear();
677  legMinVertices.clear();
678  legMidVertices.clear();
679  posEnd.hadSoFar = 0;
680  negEnd.hadSoFar = 0;
681 
682  // For closed gluon string: pick first breakup region.
683  isClosed = colConfig[iSub].isClosed;
684  if (isClosed) iParton = findFirstRegion(iSub, colConfig, event);
685 
686  // For junction topology: fragment off two of the string legs.
687  // Then iParton overwritten to remaining leg + leftover diquark.
688  pJunctionHadrons = 0.;
689  hasJunction = colConfig[iSub].hasJunction;
690  if (hasJunction && !fragmentToJunction(event)) return false;
691  int junctionHadrons = hadrons.size();
692  if (hasJunction) {
693  idPos = event[ iParton[0] ].id();
694  idNeg = event.back().id();
695  pSum -= pJunctionHadrons;
696  }
697 
698  // Set up kinematics of string evolution ( = motion).
699  system.setUp(iParton, event);
700  stopMassNow = stopMass;
701  int nExtraJoin = 0;
702 
703  // Fallback loop, when joining in the middle fails. Bailout if stuck.
704  for ( int iTry = 0; ; ++iTry) {
705  if (iTry > NTRYJOIN) {
706  infoPtr->errorMsg("Error in StringFragmentation::fragment: "
707  "stuck in joining");
708  if (hasJunction) ++nExtraJoin;
709  if (nExtraJoin > 0) event.popBack(nExtraJoin);
710  return false;
711  }
712 
713  // After several failed tries join some (extra) nearby partons.
714  if (iTry == NTRYJOIN / 3) nExtraJoin = extraJoin( 2., event);
715  if (iTry == 2 * NTRYJOIN / 3) nExtraJoin += extraJoin( 4., event);
716 
717  // After several failed tries gradually allow larger stop mass.
718  if (iTry > NTRYJOIN - NSTOPMASS) stopMassNow *= FACSTOPMASS;
719 
720  // Set up flavours of two string ends, and reset other info.
721  setStartEnds(idPos, idNeg, system);
722  pRem = pSum;
723 
724  // Begin fragmentation loop, interleaved from the two ends.
725  bool fromPos;
726 
727  // Variables used to help identifying baryons from junction splittings.
728  bool usedPosJun = false, usedNegJun = false;
729 
730  // Keep track of the momentum of hadrons taken from left and right.
731  Vec4 hadMomPos, hadMomNeg;
732 
733  // Inform the UserHooks about the string to he hadronised.
734  if ( userHooksPtr && userHooksPtr->canChangeFragPar() )
735  userHooksPtr->setStringEnds(&posEnd, &negEnd, iParton);
736 
737  for ( ; ; ) {
738 
739  // Take a step either from the positive or the negative end.
740  fromPos = (rndmPtr->flat() < 0.5);
741  StringEnd& nowEnd = (fromPos) ? posEnd : negEnd;
742 
743  // Check how many nearby string pieces there are for the next hadron.
744  double nNSP = (closePacking) ? nearStringPieces(nowEnd, rapPairs) : 0.;
745 
746  // The FlavourRope treatment changes the fragmentation parameters.
747  if (flavRopePtr) {
748  if (!flavRopePtr->doChangeFragPar(flavSelPtr, zSelPtr, pTSelPtr,
749  (fromPos ? hadMomPos.m2Calc() : hadMomNeg.m2Calc()), iParton,
750  (fromPos ? idPos : idNeg)) )
751  infoPtr->errorMsg("Error in StringFragmentation::fragment: "
752  "FlavourRope failed to change fragmentation parameters.");
753  }
754 
755  // Possibility for a user to change the fragmentation parameters.
756  if ( (userHooksPtr != 0) && userHooksPtr->canChangeFragPar() ) {
757  if ( !userHooksPtr->doChangeFragPar( flavSelPtr, zSelPtr, pTSelPtr,
758  (fromPos ? idPos : idNeg),
759  (fromPos ? hadMomPos.m2Calc() : hadMomNeg.m2Calc()),
760  iParton, &nowEnd) )
761  infoPtr->errorMsg("Error in StringFragmentation::fragment: "
762  "failed to change hadronisation parameters.");
763  }
764 
765  // Construct trial hadron and check that energy remains.
766  nowEnd.newHadron(nNSP);
767 
768  if ( energyUsedUp(fromPos) ) break;
769 
770  // Construct kinematics of the new hadron and store it.
771  Vec4 pHad = nowEnd.kinematicsHadron(system, stringVertices);
772  int statusHad = (fromPos) ? 83 : 84;
773  nowEnd.hadSoFar += 1;
774 
775  // Change status code if hadron from junction.
776  if (abs(nowEnd.idHad) > 1000 && abs(nowEnd.idHad) < 10000) {
777  if (fromPos && event[iPos].statusAbs() == 74 && !usedPosJun) {
778  statusHad = 87;
779  usedPosJun = true;
780  }
781  if (!fromPos && event[iNeg].statusAbs() == 74 && !usedNegJun) {
782  statusHad = 88;
783  usedNegJun = true;
784  }
785  if (!fromPos && hasJunction && !usedNegJun) {
786  statusHad = 88;
787  usedNegJun = true;
788  }
789  }
790 
791  // Possibility for a user to veto the hadron production.
792  if ( (userHooksPtr != 0) && userHooksPtr->canChangeFragPar() ) {
793  // Provide full particle info for veto decision.
794  if ( userHooksPtr->doVetoFragmentation( Particle( nowEnd.idHad,
795  statusHad, iPos, iNeg, 0, 0, 0, 0, pHad, nowEnd.mHad), &nowEnd ) )
796  continue;
797  }
798 
799  // Bookkeeping of momentum taken away.
800  if (fromPos) hadMomPos += pHad;
801  else hadMomNeg += pHad;
802 
803  // Append produced hadron.
804  int colHadOld = nowEnd.colOld;
805  int colHadNew = nowEnd.colNew;
806  if ( !nowEnd.fromPos ) swap(colHadOld, colHadNew);
807  hadrons.append( nowEnd.idHad, statusHad, iPos, iNeg,
808  0, 0, colHadOld, colHadNew, pHad, nowEnd.mHad);
809  if (pHad.e() < 0.) break;
810 
811  // Update string end and remaining momentum.
812  nowEnd.update();
813  pRem -= pHad;
814 
815  // End of fragmentation loop.
816  }
817 
818  // Check how many nearby string pieces there are for the last hadron.
819  double nNSP = (closePacking) ? nearStringPieces(
820  ((rndmPtr->flat() < 0.5) ? posEnd : negEnd), rapPairs) : 0.;
821 
822  // When done, join in the middle. If this works, then really done.
823  if ( finalTwo(fromPos, event, usedPosJun, usedNegJun, nNSP) ) break;
824 
825  // Else remove produced particles (except from first two junction legs)
826  // and start all over.
827  int newHadrons = hadrons.size() - junctionHadrons;
828  hadrons.popBack(newHadrons);
829  stringVertices.clear();
830  posEnd.hadSoFar = 0;
831  negEnd.hadSoFar = 0;
832  }
833 
834  // Junctions & extra joins: remove fictitious end, restore original partons.
835  if (hasJunction) ++nExtraJoin;
836  if (nExtraJoin > 0) {
837  event.popBack(nExtraJoin);
838  iParton = colConfig[iSub].iParton;
839  }
840 
841  // Store the hadrons in the normal event record, ordered from one end.
842  store(event);
843 
844  // Store hadron production space-time vertices.
845  bool saneVertices = (setVertices) ? setHadronVertices( event) : true;
846 
847  // Done.
848  return saneVertices;
849 
850 }
851 
852 //--------------------------------------------------------------------------
853 
854 // Find region where to put first string break for closed gluon loop.
855 
856 vector<int> StringFragmentation::findFirstRegion(int iSub,
857  ColConfig& colConfig, Event& event) {
858 
859  // Partons and their total four-momentum.
860  vector<int> iPartonIn = colConfig[iSub].iParton;
861 
862  // Evaluate mass-squared for all adjacent gluon pairs.
863  vector<double> m2Pair;
864  double m2Sum = 0.;
865  int size = iPartonIn.size();
866  for (int i = 0; i < size; ++i) {
867  double m2Now = 0.5 * event[ iPartonIn[i] ].p()
868  * event[ iPartonIn[(i + 1)%size] ].p();
869  m2Pair.push_back(m2Now);
870  m2Sum += m2Now;
871  }
872 
873  // Pick breakup region with probability proportional to mass-squared.
874  double m2Reg = m2Sum * rndmPtr->flat();
875  int iReg = -1;
876  do m2Reg -= m2Pair[++iReg];
877  while (m2Reg > 0. && iReg < size - 1);
878 
879  // Create reordered parton list, with breakup string region duplicated.
880  vector<int> iPartonOut;
881  for (int i = 0; i < size + 2; ++i)
882  iPartonOut.push_back( iPartonIn[(i + iReg)%size] );
883 
884  // Done.
885  return iPartonOut;
886 
887 }
888 
889 //--------------------------------------------------------------------------
890 
891 // Set flavours and momentum position for initial string endpoints.
892 
893 void StringFragmentation::setStartEnds(int idPosIn, int idNegIn,
894 StringSystem systemNow, int legNow) {
895 
896  // Variables characterizing string endpoints: defaults for open string.
897  int idPos = idPosIn;
898  int idNeg = idNegIn;
899  double px = 0.;
900  double py = 0.;
901  double Gamma = 0.;
902  double xPosFromPos = 1.;
903  double xNegFromPos = 0.;
904  double xPosFromNeg = 0.;
905  double xNegFromNeg = 1.;
906 
907  // For closed gluon loop need to pick an initial flavour.
908  if (isClosed) {
909  do {
910  int idTry = flavSelPtr->pickLightQ();
911  FlavContainer flavTry(idTry, 1);
912  flavTry = flavSelPtr->pick( flavTry);
913  flavTry = flavSelPtr->pick( flavTry);
914  idPos = flavTry.id;
915  idNeg = -idPos;
916  } while (idPos == 0);
917 
918  // Also need pT and breakup vertex position in region.
919  pair<double, double> pxy = pTSelPtr->pxy(idPos);
920  px = pxy.first;
921  py = pxy.second;
922  double m2Region = systemNow.regionLowPos(0).w2;
923  double m2Temp = min( CLOSEDM2MAX, CLOSEDM2FRAC * m2Region);
924  do {
925  double zTemp = zSelPtr->zFrag( idPos, idNeg, m2Temp);
926  xPosFromPos = 1. - zTemp;
927  xNegFromPos = m2Temp / (zTemp * m2Region);
928  } while (xNegFromPos > 1.);
929  Gamma = xPosFromPos * xNegFromPos * m2Region;
930  xPosFromNeg = xPosFromPos;
931  xNegFromNeg = xNegFromPos;
932  }
933 
934  // Initialize two string endpoints.
935  posEnd.setUp( true, iPos, idPos, systemNow.iMax, px, py,
936  Gamma, xPosFromPos, xNegFromPos,
937  systemNow.regionLowPos(0).colPos);
938  negEnd.setUp( false, iNeg, idNeg, systemNow.iMax, -px, -py,
939  Gamma, xPosFromNeg, xNegFromNeg,
940  systemNow.regionLowNeg(0).colPos);
941  // Store breakup vertex information from the first and last points.
942  if (setVertices) {
943  if (legNow == legMin) legMinVertices.push_back(
944  StringVertex( true, 0, systemNow.iMax, xPosFromPos, xNegFromPos) );
945  else if (legNow == legMid) legMidVertices.push_back(
946  StringVertex( true, 0, systemNow.iMax, xPosFromPos, xNegFromPos) );
947  else {
948  stringVertices.push_back(
949  StringVertex (true, 0, systemNow.iMax, xPosFromPos, xNegFromPos) );
950  stringVertices.push_back(
951  StringVertex( false, systemNow.iMax, 0, xPosFromNeg, xNegFromNeg) );
952  }
953  }
954  // For closed gluon loop can allow popcorn on one side but not both.
955  if (isClosed) {
956  flavSelPtr->assignPopQ(posEnd.flavOld);
957  flavSelPtr->assignPopQ(negEnd.flavOld);
958  if (rndmPtr->flat() < 0.5) posEnd.flavOld.nPop = 0;
959  else negEnd.flavOld.nPop = 0;
960  posEnd.flavOld.rank = 1;
961  negEnd.flavOld.rank = 1;
962  }
963 
964  // Done.
965 
966 }
967 
968 //--------------------------------------------------------------------------
969 
970 // Check remaining energy-momentum whether it is OK to continue.
971 
972 bool StringFragmentation::energyUsedUp(bool fromPos) {
973 
974  // If remaining negative energy then abort right away.
975  if (pRem.e() < 0.) return true;
976 
977  // Calculate W2_minimum and done if remaining W2 is below it.
978  double wMin = stopMassNow
979  + particleDataPtr->constituentMass(posEnd.flavOld.id)
980  + particleDataPtr->constituentMass(negEnd.flavOld.id);
981  if (fromPos) wMin += stopNewFlav
982  * particleDataPtr->constituentMass(posEnd.flavNew.id);
983  else wMin += stopNewFlav
984  * particleDataPtr->constituentMass(negEnd.flavNew.id);
985  wMin *= 1. + (2. * rndmPtr->flat() - 1.) * stopSmear;
986  w2Rem = pRem.m2Calc();
987  if (w2Rem < pow2(wMin)) return true;
988 
989  // Else still enough energy left to continue iteration.
990  return false;
991 
992 }
993 
994 //--------------------------------------------------------------------------
995 
996 // Store hadron production vertices in the event record.
997 
998 bool StringFragmentation::setHadronVertices( Event& event) {
999 
1000  // Order breakup points from one end to the other.
1001  bool saneVertices = true;
1002  int vertexSize = stringVertices.size();
1003  vector<StringVertex> orderedVertices;
1004  for (int i = 0; i < vertexSize; ++i) if (stringVertices[i].fromPos)
1005  orderedVertices.push_back( stringVertices[i] );
1006  for (int i = vertexSize - 1; i >= 0; --i) if (!stringVertices[i].fromPos)
1007  orderedVertices.push_back( stringVertices[i] );
1008 
1009  // Obtain space-time picture for breakup points.
1010  vector<Vec4> longitudinal;
1011  int finalSpacePos = 0;
1012  int finalVertexPos = 0;
1013  vector<int> iPartonIn = (hasJunction) ? iPartonMax : iParton;
1014  int id1 = event[ iPartonIn.front() ].idAbs();
1015  int id2 = (hasJunction) ? idDiquark : event[ iPartonIn.back() ].idAbs();
1016  int iHadJunc = (hasJunction)
1017  ? legMinVertices.size() + legMidVertices.size() - 2 : 0;
1018 
1019  // Loop over ordered vertices to obtain longitudinal space-time locations.
1020  for (int i = 0; i < vertexSize; ++i) {
1021  int iPosIn = orderedVertices[i].iRegPos;
1022  int iNegIn = orderedVertices[i].iRegNeg;
1023  double xPosIn = orderedVertices[i].xRegPos;
1024  double xNegIn = orderedVertices[i].xRegNeg;
1025 
1026  // Special case for complicated solutions in finalTwo: calculated after.
1027  if (iPosIn == -1 && iNegIn == -1) {
1028  longitudinal.push_back( Vec4( 0., 0., 0., 0.) );
1029  finalSpacePos = i + 1;
1030  finalVertexPos = i;
1031 
1032  // The normal cases.
1033  } else {
1034  StringRegion currentRegion = system.region( iPosIn, iNegIn);
1035  Vec4 gluonOffset = currentRegion.gluonOffset( iPartonIn, event, iPosIn,
1036  iNegIn) / kappaVtx;
1037  Vec4 noOffset = (xPosIn * currentRegion.pPos +
1038  xNegIn * currentRegion.pNeg) / kappaVtx;
1039  Vec4 pRegion = (currentRegion.pPos + currentRegion.pNeg) / kappaVtx;
1040  Vec4 fromBreaks = noOffset + gluonOffset;
1041 
1042  // Correction added to the space-time location of breakup points
1043  // if negative squared invariant time.
1044  if (fromBreaks.m2Calc() < 0.) {
1045  double cPlus = (-pRegion * fromBreaks + sqrt( pow2(pRegion
1046  * fromBreaks) - pRegion.m2Calc() * fromBreaks.m2Calc()))
1047  / pRegion.m2Calc();
1048  Vec4 pCorrection = noOffset + cPlus * pRegion;
1049  Vec4 fracCorrection;
1050  bool betterFrac = false;
1051  if (xPosIn < 0. || xPosIn > 1. || xNegIn < 0. || xNegIn > 1.) {
1052  xPosIn = max( 0., min( 1., xPosIn) );
1053  xNegIn = max( 0., min( 1., xNegIn) );
1054  fracCorrection = (xPosIn * currentRegion.pPos
1055  + xNegIn * currentRegion.pNeg) / kappaVtx;
1056  betterFrac = abs(noOffset.e() - pCorrection.e())
1057  > abs(noOffset.e() - fracCorrection.e());
1058  }
1059  noOffset = (betterFrac) ? fracCorrection : pCorrection;
1060  fromBreaks = noOffset + gluonOffset;
1061  }
1062 
1063  // Store vertex and check positivity.
1064  longitudinal.push_back(fromBreaks);
1065  if (fromBreaks.m2Calc() < -CHECKPOS * max(1., pow2(fromBreaks.e())))
1066  infoPtr->errorMsg("Warning in StringFragmentation::setVertices: "
1067  "negative tau^2 from breaks");
1068  }
1069  }
1070 
1071  // Breakup longitudinal space--time location for the special finalTwo case.
1072  if (finalSpacePos != 0) {
1073  double xPosIn = orderedVertices[finalVertexPos].xRegPos;
1074  double xNegIn = orderedVertices[finalVertexPos].xRegNeg;
1075  Vec4 v1 = longitudinal[finalSpacePos - 1];
1076  Vec4 v2 = longitudinal[finalSpacePos + 1];
1077  Vec4 vl = v1 - v2;
1078  Vec4 vk = pPosFinalReg + pNegFinalReg;
1079  double r = 0.;
1080  if (vl.m2Calc() > 0.) {
1081  Vec4 va = ( 0.5 * (v1 + v2) + xPosIn * pPosFinalReg +
1082  xNegIn * pNegFinalReg ) / kappaVtx;
1083  r = (va * vk - sqrt(pow2(va * vk) - vk.m2Calc() * va.m2Calc()) )
1084  / vk.m2Calc();
1085  } else r = 0.5 * sqrt(-vl.m2Calc() / vk.m2Calc());
1086  longitudinal[finalSpacePos] = ( 0.5 * (v1 + v2) + (xPosIn - r)
1087  * pPosFinalReg + (xNegIn - r) * pNegFinalReg ) / kappaVtx;
1088  }
1089 
1090  // Longitudinal offset of breakup points for massive quarks.
1091  for (int i = 0; i < vertexSize; ++i) {
1092  int iPosIn = orderedVertices[i].iRegPos;
1093  int iNegIn = orderedVertices[i].iRegNeg;
1094  if (iPosIn != -1) {
1095  StringRegion currentRegion = system.region( iPosIn, iNegIn);
1096  if ( currentRegion.massiveOffset( iPosIn, iNegIn, system.iMax,
1097  id1, id2, mc, mb) ) {
1098  Vec4 v1, v2;
1099 
1100  // Initial endpoint correction.
1101  if (i == 0 && (id1 == 4 || id1 == 5)) {
1102  int iPosIn2 = orderedVertices[i + 1].iRegPos;
1103  int iNegIn2 = orderedVertices[i + 1].iRegNeg;
1104  v2 = longitudinal[i + 1];
1105  double mHad = event[event.size() + iHadJunc - hadrons.size()].m();
1106  double pPosMass = particleDataPtr->m0(id1);
1107  if (iPosIn == iPosIn2 && iNegIn == iNegIn2) {
1108  v1 = longitudinal[i];
1109  longitudinal[i] = v1 + (pPosMass / mHad) * (v2 - v1);
1110  if (longitudinal[i].m2Calc()
1111  < -CHECKPOS * max(1., pow2(longitudinal[i].e())))
1112  infoPtr->errorMsg("Warning in StringFragmentation::set"
1113  "Vertices: negative tau^2 for endpoint massive correction");
1114  } else {
1115  StringRegion region2 = system.region( iPosIn2, iNegIn2);
1116  Vec4 gluonOffset = currentRegion.gluonOffset( iPartonIn, event,
1117  iPosIn, iNegIn);
1118  v1 = (region2.pPos + gluonOffset) / kappaVtx;
1119  longitudinal[i] = v1 + (pPosMass / mHad) * (v2 - v1);
1120  if (longitudinal[i].m2Calc()
1121  < -CHECKPOS * max(1., pow2(longitudinal[i].e())))
1122  infoPtr->errorMsg("Warning in StringFragmentation::set"
1123  "Vertices: negative tau^2 for endpoint massive correction");
1124  continue;
1125  }
1126  }
1127 
1128  // Final endpoint correction.
1129  if (i == vertexSize - 1 && (id2 == 4 || id2 == 5) && !hasJunction) {
1130  int iPosIn2 = orderedVertices[i - 1].iRegPos;
1131  int iNegIn2 = orderedVertices[i - 1].iRegNeg;
1132  double mHad = event[i - 1 + event.size() + iHadJunc
1133  - hadrons.size()].m();
1134  double pNegMass = particleDataPtr->m0(id2);
1135  if (iPosIn == iPosIn2 && iNegIn == iNegIn2) {
1136  v1 = longitudinal[i];
1137  v2 = longitudinal[i - 1] + currentRegion.massOffset / kappaVtx;
1138  longitudinal[i] = v1 + (pNegMass / mHad) * (v2 - v1);
1139  if (longitudinal[i].m2Calc()
1140  < -CHECKPOS * max(1., pow2(longitudinal[i].e())))
1141  infoPtr->errorMsg("Warning in StringFragmentation::set"
1142  "Vertices: negative tau^2 for endpoint massive correction");
1143 
1144  } else {
1145  StringRegion region2 = system.region( iPosIn2, iNegIn2);
1146  Vec4 gluonOffset = currentRegion.gluonOffset( iPartonIn, event,
1147  iPosIn, iNegIn);
1148  v1 = (region2.pNeg + gluonOffset) / kappaVtx;
1149  v2 = longitudinal[i - 1];
1150  longitudinal[i] = v1 + (pNegMass / mHad) * (v2 - v1);
1151  if (longitudinal[i].m2Calc()
1152  < -CHECKPOS * max(1., pow2(longitudinal[i].e())))
1153  infoPtr->errorMsg("Warning in StringFragmentation::set"
1154  "Vertices: negative tau^2 for endpoint massive correction");
1155  continue;
1156  }
1157  }
1158 
1159  // Add mass offset for all breakup points.
1160  Vec4 massOffset = currentRegion.massOffset / kappaVtx;
1161  Vec4 position = longitudinal[i] - massOffset;
1162  if (position.m2Calc() < 0.) {
1163  double cMinus = 0.;
1164  if (position.m2Calc() > -CHECKPOS * max(1., pow2(position.e())))
1165  position.e( position.pAbs() );
1166  else {
1167  if (massOffset.m2Calc() > 1e-6)
1168  cMinus = (longitudinal[i] * massOffset
1169  - sqrt(pow2(longitudinal[i] * massOffset)
1170  - longitudinal[i].m2Calc() * massOffset.m2Calc()))
1171  / massOffset.m2Calc();
1172  else cMinus = 0.5 * longitudinal[i].m2Calc()
1173  / (longitudinal[i] * massOffset);
1174  position = longitudinal[i] - cMinus * massOffset;
1175  }
1176  }
1177  longitudinal[i] = position;
1178  }
1179  }
1180  }
1181 
1182  // Begin smearing in transverse space. Endpoint vertices unchanged.
1183  vector<Vec4> spaceTime;
1184  for (int i = 0; i < vertexSize; ++i) {
1185  Vec4& longi = longitudinal[i];
1186  Vec4 positionTot = longi;
1187 
1188  // Find two spacelike transverse four-vector directions.
1189  if (smearOn && (isClosed || (i > 0 && i < vertexSize - 1))) {
1190  Vec4 eX, eY;
1191  int iPosIn = orderedVertices[i].iRegPos;
1192  int iNegIn = orderedVertices[i].iRegNeg;
1193  if (iPosIn == -1 && iNegIn == -1) {
1194  eX = eXFinalReg;
1195  eY = eYFinalReg;
1196  } else{
1197  StringRegion currentRegion = system.region(iPosIn, iNegIn);
1198  eX = currentRegion.eX;
1199  eY = currentRegion.eY;
1200  }
1201 
1202  // Loop over tries; give up if struck.
1203  double longiLen = sqrt(longi.pAbs2() + pow2(longi.e()) + pow2(xySmear));
1204  for (int iTry = 0; ; ++iTry) {
1205  if (iTry == NTRYSMEAR) {
1206  infoPtr->errorMsg("Warning in StringFragmentation::set"
1207  "Vertices: failed to smear vertex (normal string)");
1208  positionTot = longi;
1209  break;
1210  }
1211 
1212  // Smearing calculated randomly following a Gaussian.
1213  double transX = rndmPtr -> gauss();
1214  double transY = rndmPtr -> gauss();
1215  Vec4 transPos = xySmear * (transX * eX + transY * eY) / sqrt(2.);
1216  positionTot = transPos + longi;
1217 
1218  // Keep proper or actual time constant when including the smearing.
1219  if (constantTau)
1220  positionTot.e( sqrt(longi.m2Calc() + positionTot.pAbs2()) );
1221  if ( sqrt(transPos.pAbs2() + pow2(positionTot.e() - longi.e()))
1222  < maxSmear * longiLen) break;
1223  }
1224  }
1225  spaceTime.push_back(positionTot);
1226  }
1227 
1228  // Space-time location of the breakup points in two initial junction legs.
1229  vector<Vec4> legMinSpaceTime, legMidSpaceTime;
1230  if (hasJunction) {
1231  int hadSoFar = 0;
1232  // Loop over the two lowest-energy legs.
1233  for (int legLoop = 0; legLoop < 2; ++legLoop) {
1234  vector<StringVertex> legVertices = (legLoop == 0) ? legMinVertices
1235  : legMidVertices;
1236  StringSystem& systemNow = (legLoop == 0) ? systemMin : systemMid;
1237  vector<int> iPartonNow = (legLoop == 0) ? iPartonMinLeg : iPartonMidLeg;
1238  vector<Vec4> longitudinalPos;
1239  if (int(legVertices.size()) < 2) continue;
1240 
1241  // Loop over all breakup points of the specific leg.
1242  for (int i = 0; i < int(legVertices.size()); i++) {
1243 
1244  // Obtain longitudinal space-time location of breakup vertices.
1245  int iPosIn = legVertices[i].iRegPos;
1246  int iNegIn = legVertices[i].iRegNeg;
1247  double xPosIn = legVertices[i].xRegPos;
1248  double xNegIn = legVertices[i].xRegNeg;
1249  StringRegion currentRegion = systemNow.region( iPosIn, iNegIn);
1250  Vec4 gluonOffset = currentRegion.gluonOffsetJRF( iPartonNow, event,
1251  iPosIn, iNegIn, MtoJRF) / kappaVtx;
1252  Vec4 noOffset = (xPosIn * currentRegion.pPos
1253  + xNegIn * currentRegion.pNeg) / kappaVtx;
1254  Vec4 pRegion = (currentRegion.pPos + currentRegion.pNeg) / kappaVtx;
1255  Vec4 fromBreaks = noOffset + gluonOffset;
1256 
1257  // Correction added to the space-time location of breakup points
1258  // if negative squared invariant time.
1259  if (fromBreaks.m2Calc() < 0.) {
1260  double cPlus = (-pRegion * fromBreaks + sqrt( pow2(pRegion
1261  * fromBreaks) - pRegion.m2Calc() * fromBreaks.m2Calc()))
1262  / pRegion.m2Calc();
1263  Vec4 pCorrection = noOffset + cPlus * pRegion;
1264  Vec4 fracCorrection;
1265  bool betterFrac = false;
1266  if (xPosIn < 0. || xPosIn > 1. || xNegIn < 0. || xNegIn > 1.) {
1267  xPosIn = max( 0., min( 1., xPosIn) );
1268  xNegIn = max( 0., min( 1., xNegIn) );
1269  fracCorrection = (xPosIn * currentRegion.pPos
1270  + xNegIn * currentRegion.pNeg) / kappaVtx;
1271  betterFrac = abs(noOffset.e() - pCorrection.e())
1272  > abs(noOffset.e() - fracCorrection.e());
1273  }
1274  noOffset = (betterFrac) ? fracCorrection : pCorrection;
1275  fromBreaks = noOffset + gluonOffset;
1276  }
1277 
1278  // Store vertex and check positivity.
1279  longitudinalPos.push_back(fromBreaks);
1280  if (fromBreaks.m2Calc() < -CHECKPOS * max(1., pow2(fromBreaks.e())))
1281  infoPtr->errorMsg("Warning in StringFragmentation::setVertices: "
1282  "negative tau^2 from breaks");
1283  }
1284 
1285  // Longitudinal offset of breakup points for massive quarks.
1286  for (int i = 0; i < int(legVertices.size()); ++i) {
1287  int iPosIn = legVertices[i].iRegPos;
1288  int iNegIn = legVertices[i].iRegNeg;
1289  StringRegion currentRegion = systemNow.region( iPosIn, iNegIn);
1290  int id = event[ iPartonNow.front() ].idAbs();
1291 
1292  if ( currentRegion.massiveOffset( iPosIn, iNegIn, systemNow.iMax,
1293  id, 2, mc, mb) ) {
1294  Vec4 v1, v2;
1295  // Initial endpoint correction.
1296  if (i == 0 && (id == 4 || id == 5)) {
1297  int iPosIn2 = legVertices[i + 1].iRegPos;
1298  int iNegIn2 = legVertices[i + 1].iRegNeg;
1299  v2 = longitudinalPos[i + 1];
1300  double mHad = event[hadSoFar + event.size() - hadrons.size()].m();
1301  double pPosMass = particleDataPtr->m0(id);
1302  if (iPosIn == iPosIn2 && iNegIn == iNegIn2) {
1303  v1 = longitudinalPos[i];
1304  longitudinalPos[i] = v1 + (pPosMass / mHad) * (v2 - v1);
1305  if (longitudinalPos[i].m2Calc()
1306  < -CHECKPOS * max(1., pow2(longitudinalPos[i].e())))
1307  infoPtr->errorMsg("Warning in StringFragmentation::set"
1308  "Vertices: negative tau^2 for endpoint massive correction");
1309  } else {
1310  StringRegion region2 = systemNow.region( iPosIn2, iNegIn2);
1311  Vec4 gluonOffset = currentRegion.gluonOffsetJRF( iPartonNow,
1312  event, iPosIn2, iNegIn2, MtoJRF);
1313  v1 = (region2.pPos + gluonOffset) / kappaVtx;
1314  longitudinalPos[i] = v1 + (pPosMass / mHad) * (v2 - v1);
1315  if (longitudinalPos[i].m2Calc()
1316  < -CHECKPOS * max(1., pow2(longitudinalPos[i].e())))
1317  infoPtr->errorMsg("Warning in StringFragmentation::set"
1318  "Vertices: negative tau^2 for endpoint massive correction");
1319  continue;
1320  }
1321  }
1322 
1323  // Add mass offset for all breakup points.
1324  Vec4 massOffset = currentRegion.massOffset / kappaVtx;
1325  Vec4 position = longitudinalPos[i] - massOffset;
1326 
1327  // Correction for non-physical situations.
1328  if (position.m2Calc() < 0.) {
1329  double cMinus = 0.;
1330  if (position.m2Calc() > -CHECKPOS * max(1., pow2(position.e())))
1331  position.e( position.pAbs() );
1332  else {
1333  if (massOffset.m2Calc() > 1e-6)
1334  cMinus = (longitudinalPos[i] * massOffset
1335  - sqrt(pow2(longitudinalPos[i] * massOffset)
1336  - longitudinalPos[i].m2Calc() * massOffset.m2Calc()))
1337  / massOffset.m2Calc();
1338  else cMinus = 0.5 * longitudinalPos[i].m2Calc()
1339  / (longitudinalPos[i] * massOffset);
1340  position = longitudinalPos[i] - cMinus * massOffset;
1341  }
1342  }
1343  longitudinalPos[i] = position;
1344  }
1345  }
1346 
1347  // Begin smearing in transverse space.
1348  for (int i = 0; i < int(legVertices.size()); ++i) {
1349  Vec4& longi = longitudinalPos[i];
1350  Vec4 positionTot = longi;
1351 
1352  // Find two spacelike transverse four-vector directions.
1353  if (smearOn && i > 0) {
1354  int iPosIn = legVertices[i].iRegPos;
1355  int iNegIn = legVertices[i].iRegNeg;
1356  StringRegion currentRegion = systemNow.region( iPosIn, iNegIn);
1357  Vec4 eX = currentRegion.eX;
1358  Vec4 eY = currentRegion.eY;
1359 
1360  // Loop over tries; give up if struck.
1361  double longiLen = sqrt(longi.pAbs2() + pow2(longi.e())
1362  + pow2(xySmear));
1363  for (int iTry = 0; ; ++iTry) {
1364  if (iTry == NTRYSMEAR) {
1365  infoPtr->errorMsg("Warning in StringFragmentation::set"
1366  "Vertices: failed to smear vertex (junction string)");
1367  positionTot = longi;
1368  break;
1369  }
1370 
1371  // Smearing calculated randomly following a Gaussian.
1372  double transX = rndmPtr->gauss();
1373  double transY = rndmPtr->gauss();
1374  Vec4 transPos = xySmear * (transX * eX + transY * eY) / sqrt(2.);
1375  positionTot = transPos + longi;
1376 
1377  // Keep proper or actual time constant when including the smearing.
1378  if (constantTau)
1379  positionTot.e( sqrt(longi.m2Calc() + positionTot.pAbs2()) );
1380  if ( sqrt(transPos.pAbs2() + pow2(positionTot.e() - longi.e()))
1381  < maxSmear * longiLen) break;
1382  }
1383  }
1384 
1385  // Boost from the rest frame of the junction to the original frame.
1386  positionTot.rotbst(MfromJRF);
1387 
1388  // Recalculate time to compensate for numerical precision loss
1389  // in iterative calculation of MfromJRF. Store final result.
1390  if (positionTot.m2Calc() < 0.)
1391  positionTot.e( positionTot.pAbs() );
1392  if (legLoop == 0) legMinSpaceTime.push_back(positionTot);
1393  else legMidSpaceTime.push_back(positionTot);
1394 
1395  // End of loops over all breakup point of two lowest-energy legs.
1396  }
1397  hadSoFar = hadSoFar + legVertices.size() - 1;
1398  }
1399  }
1400 
1401  // Calculate hadron production points of the two initial junction legs.
1402  if (hasJunction) {
1403  int hadSoFar = 0;
1404 
1405  // Loop over the two lowest-energy legs.
1406  for (int legLoop = 0; legLoop < 2; legLoop++) {
1407  vector<Vec4>& finalLocation = (legLoop == 0) ? legMinSpaceTime
1408  : legMidSpaceTime;
1409  vector<int> iPartonNow = (legLoop == 0) ? iPartonMinLeg : iPartonMidLeg;
1410  int id = event[iPartonNow.front()].idAbs();
1411 
1412  // Calculate hadron production points from breakup vertices
1413  // using one of the three definitions.
1414  for (int i = 0; i < int(finalLocation.size()) - 1; ++i) {
1415  Vec4 middlePoint = 0.5 * (finalLocation[i] + finalLocation[i + 1]);
1416  if (abs(middlePoint.mCalc()) > maxTau) saneVertices = false;
1417  int iHad = i + hadSoFar + event.size() - hadrons.size();
1418  Vec4 pHad = event[iHad].p();
1419  Vec4 prodPoints = Vec4( 0., 0., 0., 0.);
1420 
1421  // Reduced oscillation period if system contains massive quarks.
1422  double mHad = event[iHad].m();
1423  double redOsc = (i == 0 && (id == 4 || id == 5))
1424  ? 1. - pow2(particleDataPtr->m0(id) / mHad) : 1.;
1425 
1426  // Hadron production points calculation depending on definition.
1427  if (hadronVertex == 0) prodPoints = middlePoint;
1428  else if (hadronVertex == 1)
1429  prodPoints = middlePoint + 0.5 * redOsc * pHad / kappaVtx;
1430  else {
1431  prodPoints = middlePoint - 0.5 * redOsc * pHad / kappaVtx;
1432  if (prodPoints.m2Calc() < 0. || prodPoints.e() < 0.) {
1433  double midpProd = redOsc * middlePoint * pHad;
1434  double tau0fac = 2. * (midpProd - sqrt( pow2(midpProd)
1435  - middlePoint.m2Calc() * pow2(redOsc * mHad)))
1436  / pow2(redOsc * mHad);
1437  prodPoints = middlePoint - 0.5 * tau0fac * redOsc * pHad
1438  / kappaVtx;
1439  }
1440  }
1441  event[iHad].vProd( event[iHad].vProd() + prodPoints * FM2MM );
1442  }
1443 
1444  // End of the two legs loop. Number of hadrons with stored vertices.
1445  if (finalLocation.size() > 0) hadSoFar += finalLocation.size() - 1;
1446  }
1447  }
1448 
1449  // Normal string system or last leg: calculate hadron production points
1450  // from breakup vertices using one of the three definitions.
1451  for (int i = 0; i < int(spaceTime.size()) - 1; ++i) {
1452  Vec4 middlePoint = 0.5 * (spaceTime[i] + spaceTime[i + 1]);
1453  if (abs(middlePoint.mCalc()) > maxTau) saneVertices = false;
1454  int iHad = i + iHadJunc + event.size() - hadrons.size();
1455  Vec4 pHad = event[iHad].p();
1456  Vec4 prodPoints = Vec4( 0., 0., 0., 0.);
1457 
1458  // Reduced oscillation period if system contains massive quarks.
1459  double redOsc = 1.;
1460  double mHad = event[iHad].m();
1461  if (i == 0 && (id1 == 4 || id1 == 5))
1462  redOsc = 1. - pow2( particleDataPtr->m0(id1) / mHad);
1463  if (i == int(spaceTime.size()) - 1 && (id2 == 4 || id2 == 5))
1464  redOsc = 1. - pow2( particleDataPtr->m0(id2) / mHad);
1465 
1466  // Hadron production points calculation depending on definition.
1467  if (hadronVertex == 0) prodPoints = middlePoint;
1468  else if (hadronVertex == 1)
1469  prodPoints = middlePoint + 0.5 * redOsc * pHad / kappaVtx;
1470  else {
1471  prodPoints = middlePoint - 0.5 * redOsc * pHad / kappaVtx;
1472  if (prodPoints.m2Calc() < 0. || prodPoints.e() < 0.) {
1473  double midpProd = redOsc * middlePoint * pHad;
1474  double tau0fac = 2. * ( midpProd - sqrt( pow2(midpProd)
1475  - middlePoint.m2Calc() * pow2(redOsc * mHad))) / pow2(redOsc * mHad);
1476  prodPoints = middlePoint - 0.5 * tau0fac * redOsc * pHad / kappaVtx;
1477  }
1478  }
1479  event[iHad].vProd( event[iHad].vProd() + prodPoints * FM2MM );
1480  }
1481 
1482  // Done.
1483  if (!saneVertices) infoPtr->errorMsg("Error in StringFragmentation::set"
1484  "Vertices: too large |tau| so make new try");
1485  return saneVertices;
1486 
1487 }
1488 
1489 //--------------------------------------------------------------------------
1490 
1491 // Produce the final two partons to complete the system.
1492 
1493 bool StringFragmentation::finalTwo(bool fromPos, Event& event,
1494  bool usedPosJun, bool usedNegJun, double nNSP) {
1495 
1496  // Check whether we went too far in p+-.
1497  if (pRem.e() < 0. || w2Rem < 0. || (hadrons.size() > 0
1498  && hadrons.back().e() < 0.) ) return false;
1499  if ( posEnd.iPosOld > negEnd.iPosOld || negEnd.iNegOld > posEnd.iNegOld)
1500  return false;
1501  if ( posEnd.iPosOld == negEnd.iPosOld && posEnd.xPosOld < negEnd.xPosOld)
1502  return false;
1503  if ( posEnd.iNegOld == negEnd.iNegOld && posEnd.xNegOld > negEnd.xNegOld)
1504  return false;
1505 
1506  // Impossible to join two diquarks.
1507  FlavContainer flav1 = (fromPos) ? posEnd.flavNew.anti() : posEnd.flavOld;
1508  FlavContainer flav2 = (fromPos) ? negEnd.flavOld : negEnd.flavNew.anti();
1509  if (flav1.isDiquark() && flav2.isDiquark()) return false;
1510 
1511  // Construct preliminary hadron pT as if no region changes.
1512  double pHadPrelim[2] = { 0.0, 0.0 };
1513  if (fromPos) {
1514  pHadPrelim[0] = negEnd.pxOld-posEnd.pxNew;
1515  pHadPrelim[1] = negEnd.pyOld-posEnd.pyNew;
1516  } else {
1517  pHadPrelim[0] = posEnd.pxOld-negEnd.pxNew;
1518  pHadPrelim[1] = posEnd.pyOld-negEnd.pyNew;
1519  }
1520  double pThadPrelim = sqrt( pow2(pHadPrelim[0]) + pow2(pHadPrelim[1]) );
1521 
1522  // Construct the final hadron from the leftover flavours. Break if stuck.
1523  int idHad;
1524  for (int iTry = 0; iTry < NTRYFLAV; ++iTry) {
1525  idHad = flavSelPtr->getHadronID( flav1, flav2, pThadPrelim, nNSP, true);
1526  if (idHad != 0) break;
1527  }
1528  if (idHad == 0) return false;
1529 
1530  // Store the final particle and its new pT, and construct its mass.
1531  if (fromPos) {
1532  negEnd.idHad = idHad;
1533  negEnd.pxNew = -posEnd.pxNew;
1534  negEnd.pyNew = -posEnd.pyNew;
1535  negEnd.mHad = flavSelPtr->getHadronMassWin(idHad);
1536  } else {
1537  posEnd.idHad = idHad;
1538  posEnd.pxNew = -negEnd.pxNew;
1539  posEnd.pyNew = -negEnd.pyNew;
1540  posEnd.mHad = flavSelPtr->getHadronMassWin(idHad);
1541  }
1542 
1543  // String region in which to do the joining.
1544  StringRegion region = finalRegion();
1545  if (region.isEmpty) return false;
1546 
1547  // Project remaining momentum along longitudinal and transverse directions.
1548  region.project( pRem);
1549  double pxRem = region.px() - posEnd.pxOld - negEnd.pxOld;
1550  double pyRem = region.py() - posEnd.pyOld - negEnd.pyOld;
1551  double xPosRem = region.xPos();
1552  double xNegRem = region.xNeg();
1553 
1554  // Share extra pT kick evenly between final two hadrons.
1555  posEnd.pxOld += 0.5 * pxRem;
1556  posEnd.pyOld += 0.5 * pyRem;
1557  negEnd.pxOld += 0.5 * pxRem;
1558  negEnd.pyOld += 0.5 * pyRem;
1559 
1560  // Construct new pT and mT of the final two particles.
1561  posEnd.pxHad = posEnd.pxOld + posEnd.pxNew;
1562  posEnd.pyHad = posEnd.pyOld + posEnd.pyNew;
1563  posEnd.mT2Had = pow2(posEnd.mHad) + pow2(posEnd.pxHad)
1564  + pow2(posEnd.pyHad);
1565  negEnd.pxHad = negEnd.pxOld + negEnd.pxNew;
1566  negEnd.pyHad = negEnd.pyOld + negEnd.pyNew;
1567  negEnd.mT2Had = pow2(negEnd.mHad) + pow2(negEnd.pxHad)
1568  + pow2(negEnd.pyHad);
1569 
1570  // Construct remaining system transverse mass.
1571  double wT2Rem = w2Rem + pow2( posEnd.pxHad + negEnd.pxHad)
1572  + pow2( posEnd.pyHad + negEnd.pyHad);
1573 
1574  // Check that kinematics possible.
1575  if ( sqrt(wT2Rem) < sqrt(posEnd.mT2Had) + sqrt(negEnd.mT2Had) )
1576  return false;
1577  double lambda2 = pow2( wT2Rem - posEnd.mT2Had - negEnd.mT2Had)
1578  - 4. * posEnd.mT2Had * negEnd.mT2Had;
1579  if (lambda2 <= 0.) return false;
1580 
1581  // Construct kinematics, as viewed in the transverse rest frame.
1582  double lambda = sqrt(lambda2);
1583  double probReverse = 1. / (1. + exp( min( EXPMAX, bLund * lambda) ) );
1584  double xpzPos = 0.5 * lambda/ wT2Rem;
1585  if (probReverse > rndmPtr->flat()) xpzPos = -xpzPos;
1586  double xmDiff = (posEnd.mT2Had - negEnd.mT2Had) / wT2Rem;
1587  double xePos = 0.5 * (1. + xmDiff);
1588  double xeNeg = 0.5 * (1. - xmDiff );
1589 
1590  // Translate this into kinematics in the string frame.
1591  Vec4 pHadPos = region.pHad( (xePos + xpzPos) * xPosRem,
1592  (xePos - xpzPos) * xNegRem, posEnd.pxHad, posEnd.pyHad);
1593  Vec4 pHadNeg = region.pHad( (xeNeg - xpzPos) * xPosRem,
1594  (xeNeg + xpzPos) * xNegRem, negEnd.pxHad, negEnd.pyHad);
1595 
1596  // Project pHadPos momentum fraction on the positive region
1597  // to obtain breakup vertices with respect to that region.
1598  if (setVertices) {
1599  StringRegion posRegion = system.region( posEnd.iPosOld, posEnd.iNegOld);
1600  posRegion.project(pHadPos);
1601  double xFromPosPos = posEnd.xPosOld - posRegion.xPos();
1602  double xFromPosNeg = posEnd.xNegOld + posRegion.xNeg();
1603 
1604  // Same, but projecting pHadNeg fractions on the negative region.
1605  StringRegion negRegion = system.region( negEnd.iPosOld, negEnd.iNegOld);
1606  negRegion.project(pHadNeg);
1607  double xFromNegPos = negEnd.xPosOld + negRegion.xPos();
1608  double xFromNegNeg = negEnd.xNegOld - negRegion.xNeg();
1609 
1610  // Store energy-momentum coordinates for the final breakup vertex.
1611  // If projections give valid results, store them as breakup fractions.
1612  if (xFromPosPos > 0. && xFromPosPos < 1. && xFromPosNeg > 0.
1613  && xFromPosNeg < 1.) stringVertices.push_back( StringVertex(
1614  fromPos, posEnd.iPosOld, posEnd.iNegOld, xFromPosPos, xFromPosNeg) );
1615  else if (xFromNegPos > 0. && xFromNegPos < 1. && xFromNegNeg > 0.
1616  && xFromNegNeg < 1.) stringVertices.push_back( StringVertex(
1617  fromPos, negEnd.iPosOld, negEnd.iNegOld, xFromNegPos, xFromNegNeg) );
1618 
1619  // If above procedures do not work, calculate a new zHad and use
1620  // the kinematicsHadron method, first from the positive end.
1621  else {
1622  double gammaPosOld = posEnd.GammaOld;
1623  double gammaNegOld = negEnd.GammaOld;
1624  double zNewReg = 0.;
1625  if (posEnd.hadSoFar == 0) zNewReg = wT2Rem / (wT2Rem + gammaNegOld);
1626  else zNewReg = 0.5 * ( sqrt( pow2(wT2Rem + gammaNegOld - gammaPosOld)
1627  + 4. * wT2Rem * gammaPosOld) - (wT2Rem + gammaNegOld - gammaPosOld) )
1628  / gammaPosOld;
1629  double zHad = (xePos + xpzPos) * zNewReg;
1630  Vec4 proof = posEnd.kinematicsHadron(system, stringVertices, true, zHad);
1631 
1632  // Try negative-end kinematicsHadron method if positive-end one failed.
1633  if (proof.e() < -1e-8) {
1634  if (negEnd.hadSoFar == 0) zNewReg = wT2Rem / (wT2Rem + gammaPosOld);
1635  else zNewReg = 0.5 * ( sqrt( pow2(wT2Rem + gammaPosOld - gammaNegOld)
1636  + 4. * wT2Rem * gammaNegOld) - (wT2Rem + gammaPosOld - gammaNegOld) )
1637  / gammaNegOld;
1638  zHad = (xeNeg + xpzPos) * zNewReg;
1639  proof = negEnd.kinematicsHadron( system, stringVertices, true, zHad);
1640 
1641  // As last resort use the final region created by the finalTwo method.
1642  if (proof.e() < -1.) {
1643  pPosFinalReg = region.pPos;
1644  pNegFinalReg = region.pNeg;
1645  eXFinalReg = region.eX;
1646  eYFinalReg = region.eY;
1647  stringVertices.push_back( StringVertex( true, -1, -1,
1648  1. - (xePos + xpzPos) * xPosRem, (xePos - xpzPos) * xNegRem) );
1649  }
1650  }
1651  }
1652  }
1653 
1654  // Update status codes for junction baryons.
1655  int statusHadPos = 83;
1656  int statusHadNeg = 84;
1657  if (fromPos) {
1658  if (abs(posEnd.idHad) > 1000 && abs(posEnd.idHad) < 10000) {
1659  if (event[iPos].statusAbs() == 74 && !usedPosJun) {
1660  statusHadPos = 87;
1661  usedPosJun = true;
1662  }
1663  }
1664  if (abs(idHad) > 1000 && abs(idHad) < 10000) {
1665  if ((!usedNegJun && (event[iNeg].statusAbs() == 74 || hasJunction))
1666  || (!usedPosJun && event[iPos].statusAbs() == 74)) {
1667  statusHadNeg = 88;
1668  }
1669  }
1670  } else {
1671  if (abs(negEnd.idHad) > 1000 && abs(negEnd.idHad) < 10000) {
1672  if (!usedNegJun && (event[iNeg].statusAbs() == 74 || hasJunction)) {
1673  statusHadNeg = 88;
1674  usedNegJun = true;
1675  }
1676  }
1677  if (abs(idHad) > 1000 && abs(idHad) < 10000) {
1678  if ((!usedNegJun && (event[iNeg].statusAbs() == 74 || hasJunction))
1679  || (!usedPosJun && event[iPos].statusAbs() == 74)) {
1680  statusHadPos = 87;
1681  }
1682  }
1683  }
1684 
1685  int colMid = (fromPos? negEnd.colOld: posEnd.colOld);
1686  // Possibility for a user to veto the hadron production.
1687  if ( (userHooksPtr != 0) && userHooksPtr->canChangeFragPar() ) {
1688  // Provide full particle info for veto decision.
1689  if ( userHooksPtr->doVetoFragmentation(
1690  Particle( posEnd.idHad, statusHadPos, posEnd.iEnd, negEnd.iEnd,
1691  0, 0, posEnd.colOld, colMid, pHadPos, posEnd.mHad),
1692  Particle( negEnd.idHad, statusHadNeg, posEnd.iEnd, negEnd.iEnd,
1693  0, 0, colMid, negEnd.colOld, pHadNeg, negEnd.mHad),
1694  &posEnd, &negEnd ) ) return false;
1695  }
1696 
1697  // Add produced particles to the event record.
1698  hadrons.append( posEnd.idHad, statusHadPos, posEnd.iEnd, negEnd.iEnd,
1699  0, 0, posEnd.colOld, colMid, pHadPos, posEnd.mHad);
1700  hadrons.append( negEnd.idHad, statusHadNeg, posEnd.iEnd, negEnd.iEnd,
1701  0, 0, colMid, negEnd.colOld, pHadNeg, negEnd.mHad);
1702 
1703  // It worked.
1704  return true;
1705 
1706 }
1707 
1708 //--------------------------------------------------------------------------
1709 
1710 // Construct a special joining region for the final two hadrons.
1711 
1712 StringRegion StringFragmentation::finalRegion() {
1713 
1714  // Simple case when both string ends are in the same region.
1715  if (posEnd.iPosOld == negEnd.iPosOld && posEnd.iNegOld == negEnd.iNegOld)
1716  return system.region( posEnd.iPosOld, posEnd.iNegOld);
1717 
1718  // Start out with empty region. (Empty used for error returns.)
1719  StringRegion region;
1720 
1721  // Add up all remaining p+.
1722  Vec4 pPosJoin;
1723  int colPos = system.regionLowPos(posEnd.iPosOld).colPos;
1724  int colNeg = system.regionLowNeg(negEnd.iNegOld).colNeg;
1725  if ( posEnd.iPosOld == negEnd.iPosOld) {
1726  double xPosJoin = posEnd.xPosOld - negEnd.xPosOld;
1727  if (xPosJoin < 0.) return region;
1728  pPosJoin = system.regionLowPos(posEnd.iPosOld).pHad( xPosJoin, 0., 0., 0.);
1729  } else {
1730  for (int iPosNow = posEnd.iPosOld; iPosNow <= negEnd.iPosOld; ++iPosNow) {
1731  if (iPosNow == posEnd.iPosOld) pPosJoin
1732  += system.regionLowPos(iPosNow).pHad( posEnd.xPosOld, 0., 0., 0.);
1733  else if (iPosNow == negEnd.iPosOld) pPosJoin
1734  += system.regionLowPos(iPosNow).pHad( 1. - negEnd.xPosOld, 0., 0., 0.);
1735  else pPosJoin += system.regionLowPos(iPosNow).pHad( 1., 0., 0., 0.);
1736  }
1737  }
1738 
1739  // Add up all remaining p-.
1740  Vec4 pNegJoin;
1741  if ( negEnd.iNegOld == posEnd.iNegOld) {
1742  double xNegJoin = negEnd.xNegOld - posEnd.xNegOld;
1743  if (xNegJoin < 0.) return region;
1744  pNegJoin = system.regionLowNeg(negEnd.iNegOld).pHad( 0., xNegJoin, 0., 0.);
1745  } else {
1746  for (int iNegNow = negEnd.iNegOld; iNegNow <= posEnd.iNegOld; ++iNegNow) {
1747  if (iNegNow == negEnd.iNegOld) pNegJoin
1748  += system.regionLowNeg(iNegNow).pHad( 0., negEnd.xNegOld, 0., 0.);
1749  else if (iNegNow == posEnd.iNegOld) pNegJoin
1750  += system.regionLowNeg(iNegNow).pHad( 0., 1. - posEnd.xNegOld, 0., 0.);
1751  else pNegJoin += system.regionLowNeg(iNegNow).pHad( 0., 1., 0., 0.);
1752  }
1753  }
1754 
1755  // For a closed gluon loop pPosJoin == pNegJoin and the above does not work.
1756  // So reshuffle; "perfect" for g g systems, OK in general.
1757  Vec4 pTest = pPosJoin - pNegJoin;
1758  if ( abs(pTest.px()) + abs(pTest.py()) + abs(pTest.pz()) + abs(pTest.e())
1759  < MATCHPOSNEG * (pPosJoin.e() + pNegJoin.e()) ) {
1760  Vec4 delta
1761  = system.regionLowPos(posEnd.iPosOld + 1).pHad( 1., 0., 0., 0.)
1762  - system.regionLowNeg(negEnd.iNegOld + 1).pHad( 0., 1., 0., 0.);
1763  // If reshuffle did not help then pick random axis to break tie.
1764  // (Needed for low-mass q-g-qbar with q-qbar perfectly parallel.)
1765  if ( abs(delta.px()) + abs(delta.py()) + abs(delta.pz()) + abs(delta.e())
1766  < MATCHPOSNEG * (pPosJoin.e() + pNegJoin.e()) ) {
1767  double cthe = 2. * rndmPtr->flat() - 1.;
1768  double sthe = sqrtpos(1. - cthe * cthe);
1769  double phi = 2. * M_PI * rndmPtr->flat();
1770  delta = 0.5 * min( pPosJoin.e(), pNegJoin.e())
1771  * Vec4( sthe * sin(phi), sthe * cos(phi), cthe, 0.);
1772  infoPtr->errorMsg("Warning in StringFragmentation::finalRegion: "
1773  "random axis needed to break tie");
1774  }
1775  pPosJoin -= delta;
1776  pNegJoin += delta;
1777  }
1778 
1779  // Construct a new region from remaining p+ and p-.
1780  region.setUp( pPosJoin, pNegJoin, colPos, colNeg);
1781  if (region.isEmpty) return region;
1782 
1783  // Project the existing pTold vectors onto the new directions.
1784  Vec4 pTposOld = system.region( posEnd.iPosOld, posEnd.iNegOld).pHad(
1785  0., 0., posEnd.pxOld, posEnd.pyOld);
1786  region.project( pTposOld);
1787  posEnd.pxOld = region.px();
1788  posEnd.pyOld = region.py();
1789  Vec4 pTnegOld = system.region( negEnd.iPosOld, negEnd.iNegOld).pHad(
1790  0., 0., negEnd.pxOld, negEnd.pyOld);
1791  region.project( pTnegOld);
1792  negEnd.pxOld = region.px();
1793  negEnd.pyOld = region.py();
1794 
1795  // Done.
1796  return region;
1797 
1798 }
1799 
1800 //--------------------------------------------------------------------------
1801 
1802 // Store the hadrons in the normal event record, ordered from one end.
1803 
1804 void StringFragmentation::store(Event& event) {
1805 
1806  // Starting position.
1807  int iFirst = event.size();
1808 
1809  // Remove colour indices to avoid confusion by default.
1810  if ( !traceColours )
1811  for (int i = 0; i < hadrons.size(); ++i) {
1812  hadrons[i].col(0);
1813  hadrons[i].acol(0);
1814  }
1815 
1816  // Copy straight over from first two junction legs.
1817  if (hasJunction) {
1818  for (int i = 0; i < hadrons.size(); ++i)
1819  if (hadrons[i].status() == 85 || hadrons[i].status() == 86)
1820  event.append( hadrons[i] );
1821  }
1822 
1823  // Loop downwards, copying all from the positive end.
1824  for (int i = 0; i < hadrons.size(); ++i)
1825  if (hadrons[i].status() == 83 || hadrons[i].status() == 87)
1826  event.append( hadrons[i] );
1827 
1828  // Loop upwards, copying all from the negative end.
1829  for (int i = hadrons.size() - 1; i >= 0 ; --i)
1830  if (hadrons[i].status() == 84 || hadrons[i].status() == 88)
1831  event.append( hadrons[i] );
1832 
1833  int iLast = event.size() - 1;
1834 
1835  // Set decay vertex when this is displaced.
1836  if (event[posEnd.iEnd].hasVertex()) {
1837  Vec4 vDec = event[posEnd.iEnd].vDec();
1838  for (int i = iFirst; i <= iLast; ++i) event[i].vProd( vDec );
1839  }
1840 
1841  // Set lifetime of hadrons.
1842  for (int i = iFirst; i <= iLast; ++i)
1843  event[i].tau( event[i].tau0() * rndmPtr->exp() );
1844 
1845  // Mark original partons as hadronized and set their daughter range.
1846  for (int i = 0; i < int(iParton.size()); ++i)
1847  if (iParton[i] >= 0) {
1848  event[ iParton[i] ].statusNeg();
1849  event[ iParton[i] ].daughters(iFirst, iLast);
1850  }
1851 
1852 }
1853 
1854 //--------------------------------------------------------------------------
1855 
1856 // Fragment off two of the string legs in to a junction.
1857 
1858 bool StringFragmentation::fragmentToJunction(Event& event) {
1859 
1860  // Identify range of partons on the three legs.
1861  // (Each leg begins with an iParton[i] = -(10 + 10*junctionNumber + leg),
1862  // and partons then appear ordered from the junction outwards.)
1863  int legBeg[3] = { 0, 0, 0};
1864  int legEnd[3] = { 0, 0, 0};
1865  int leg = -1;
1866  // PS (4/10/2011) Protect against invalid systems
1867  if (iParton[0] > 0) {
1868  infoPtr->errorMsg("Error in StringFragmentation::fragment"
1869  "ToJunction: iParton[0] not a valid junctionNumber");
1870  return false;
1871  }
1872  for (int i = 0; i < int(iParton.size()); ++i) {
1873  if (iParton[i] < 0) {
1874  if (leg == 2) {
1875  infoPtr->errorMsg("Error in StringFragmentation::fragment"
1876  "ToJunction: unprocessed multi-junction system");
1877  return false;
1878  }
1879  legBeg[++leg] = i + 1;
1880  }
1881  else legEnd[leg] = i;
1882  }
1883 
1884  // Iterate from system rest frame towards the junction rest frame (JRF).
1885  RotBstMatrix Mstep;
1886  MtoJRF.reset();
1887  MtoJRF.bstback(pSum);
1888  Vec4 pWTinJRF[3];
1889  int iter = 0;
1890  double errInCM = 0.;
1891 
1892  do {
1893  ++iter;
1894  // Find weighted sum of momenta on the three sides of the junction.
1895  for (leg = 0; leg < 3; ++ leg) {
1896  pWTinJRF[leg] = 0.;
1897  double eWeight = 0.;
1898  for (int i = legBeg[leg]; i <= legEnd[leg]; ++i) {
1899  Vec4 pTemp = event[ iParton[i] ].p();
1900  pTemp.rotbst(MtoJRF);
1901  pWTinJRF[leg] += pTemp * exp(-eWeight);
1902  eWeight += pTemp.e() / eNormJunction;
1903  if (eWeight > EJNWEIGHTMAX) break;
1904  }
1905  }
1906 
1907  // Store original deviation from 120 degree topology.
1908  if (iter == 1) errInCM = pow2(costheta(pWTinJRF[0], pWTinJRF[1]) + 0.5)
1909  + pow2(costheta(pWTinJRF[0], pWTinJRF[2]) + 0.5)
1910  + pow2(costheta(pWTinJRF[1], pWTinJRF[2]) + 0.5);
1911 
1912  // Check numerical instabilities in boost, use rest frame if it fails.
1913  if ( (pWTinJRF[0] + pWTinJRF[1]).m2Calc() < M2MINJRF
1914  || (pWTinJRF[0] + pWTinJRF[2]).m2Calc() < M2MINJRF
1915  || (pWTinJRF[1] + pWTinJRF[2]).m2Calc() < M2MINJRF ) {
1916  infoPtr->errorMsg("Warning in StringFragmentation::fragmentTo"
1917  "Junction: Negative invariant masses in junction rest frame");
1918  MtoJRF.reset();
1919  MtoJRF.bstback(pSum);
1920  break;
1921  }
1922  // Find new JRF from the set of weighted momenta.
1923  Mstep = junctionRestFrame( pWTinJRF[0], pWTinJRF[1], pWTinJRF[2]);
1924  // Fortran code will not take full step after the first few
1925  // iterations. How implement this in terms of an M matrix??
1926  MtoJRF.rotbst( Mstep );
1927  } while (iter < 3 || (Mstep.deviation() > CONVJNREST && iter < NTRYJNREST) );
1928 
1929  // If final deviation from 120 degrees is bigger than in CM then revert.
1930  double errInJRF = pow2(costheta(pWTinJRF[0], pWTinJRF[1]) + 0.5)
1931  + pow2(costheta(pWTinJRF[0], pWTinJRF[2]) + 0.5)
1932  + pow2(costheta(pWTinJRF[1], pWTinJRF[2]) + 0.5);
1933  if (errInJRF > errInCM + CONVJNREST) {
1934  infoPtr->errorMsg("Warning in StringFragmentation::fragmentTo"
1935  "Junction: bad convergence junction rest frame");
1936  MtoJRF.reset();
1937  MtoJRF.bstback(pSum);
1938  }
1939 
1940  // Opposite operation: boost from JRF to original system.
1941  MfromJRF = MtoJRF;
1942  MfromJRF.invert();
1943 
1944  // Sum leg four-momenta in original frame and in JRF.
1945  Vec4 pInLeg[3], pInJRF[3];
1946  for (leg = 0; leg < 3; ++leg) {
1947  pInLeg[leg] = 0.;
1948  for (int i = legBeg[leg]; i <= legEnd[leg]; ++i)
1949  pInLeg[leg] += event[ iParton[i] ].p();
1950  pInJRF[leg] = pInLeg[leg];
1951  pInJRF[leg].rotbst(MtoJRF);
1952  }
1953 
1954  // Pick the two legs with lowest energy in JRF.
1955  legMin = (pInJRF[0].e() < pInJRF[1].e()) ? 0 : 1;
1956  int legMax = 1 - legMin;
1957  if (pInJRF[2].e() < min(pInJRF[0].e(), pInJRF[1].e()) ) legMin = 2;
1958  else if (pInJRF[2].e() > max(pInJRF[0].e(), pInJRF[1].e()) ) legMax = 2;
1959  legMid = 3 - legMin - legMax;
1960 
1961  // Save info on which status codes belong with the three legs.
1962  int iJunction = (-iParton[0]) / 10 - 1;
1963  event.statusJunction( iJunction, legMin, 85);
1964  event.statusJunction( iJunction, legMid, 86);
1965  event.statusJunction( iJunction, legMax, 83);
1966 
1967  // Temporarily copy the partons on the low-energy legs, into the JRF,
1968  // in reverse order, so (anti)quark leg end first.
1969  vector<int> iPartonMin;
1970  iPartonMinLeg.clear();
1971  for (int i = legEnd[legMin]; i >= legBeg[legMin]; --i) {
1972  if (setVertices) iPartonMinLeg.push_back( iParton[i] );
1973  int iNew = event.append( event[ iParton[i] ] );
1974  event[iNew].rotbst(MtoJRF);
1975  iPartonMin.push_back( iNew );
1976  }
1977  vector<int> iPartonMid;
1978  iPartonMidLeg.clear();
1979  for (int i = legEnd[legMid]; i >= legBeg[legMid]; --i) {
1980  if (setVertices) iPartonMidLeg.push_back( iParton[i] );
1981  int iNew = event.append( event[ iParton[i] ] );
1982  event[iNew].rotbst(MtoJRF);
1983  iPartonMid.push_back( iNew );
1984  }
1985 
1986  // Find final weighted sum of momenta on each of the two legs.
1987  double eWeight = 0.;
1988  pWTinJRF[legMin] = 0.;
1989  for (int i = iPartonMin.size() - 1; i >= 0; --i) {
1990  pWTinJRF[legMin] += event[ iPartonMin[i] ].p() * exp(-eWeight);
1991  eWeight += event[ iPartonMin[i] ].e() / eNormJunction;
1992  if (eWeight > EJNWEIGHTMAX) break;
1993  }
1994  eWeight = 0.;
1995  pWTinJRF[legMid] = 0.;
1996  for (int i = iPartonMid.size() - 1; i >= 0; --i) {
1997  pWTinJRF[legMid] += event[ iPartonMid[i] ].p() * exp(-eWeight);
1998  eWeight += event[ iPartonMid[i] ].e() / eNormJunction;
1999  if (eWeight > EJNWEIGHTMAX) break;
2000  }
2001 
2002  // Define fictitious opposing partons in JRF and store as string ends.
2003  Vec4 pOppose = pWTinJRF[legMin];
2004  pOppose.flip3();
2005  int idOppose = (rndmPtr->flat() > 0.5) ? 2 : 1;
2006  int colOppose = event[iPartonMin.back()].acol();
2007  int acolOppose = 0;
2008  if (event[ iPartonMin[0] ].col() > 0) {
2009  idOppose = -idOppose;
2010  colOppose = 0;
2011  acolOppose = event[iPartonMin.back()].col();
2012  }
2013  int iOppose = event.append( idOppose, 77, 0, 0, 0, 0, colOppose, acolOppose,
2014  pOppose, 0.);
2015  iPartonMin.push_back( iOppose);
2016  pOppose = pWTinJRF[legMid];
2017  pOppose.flip3();
2018  idOppose = (rndmPtr->flat() > 0.5) ? 2 : 1;
2019  colOppose = event[iPartonMid.back()].acol();
2020  acolOppose = 0;
2021  if (event[ iPartonMid[0] ].col() > 0) {
2022  idOppose = -idOppose;
2023  colOppose = 0;
2024  acolOppose = event[iPartonMid.back()].col();
2025  }
2026  iOppose = event.append( idOppose, 77, 0, 0, 0, 0, colOppose, acolOppose,
2027  pOppose, 0.);
2028  iPartonMid.push_back( iOppose);
2029 
2030  // Set up kinematics of string evolution in low-energy temporary systems.
2031  systemMin.setUp(iPartonMin, event);
2032  systemMid.setUp(iPartonMid, event);
2033 
2034  // Outer fallback loop, when too little energy left for third leg.
2035  int idMin = 0;
2036  int idMid = 0;
2037  Vec4 pDiquark;
2038  for ( int iTryOuter = 0; ; ++iTryOuter) {
2039 
2040  // Middle fallback loop, when much unused energy in leg remnants.
2041  double eLeftMin = 0.;
2042  double eLeftMid = 0.;
2043  for ( int iTryMiddle = 0; ; ++iTryMiddle) {
2044 
2045  // Loop over the two lowest-energy legs.
2046  for (int legLoop = 0; legLoop < 2; ++ legLoop) {
2047  int legNow = (legLoop == 0) ? legMin : legMid;
2048 
2049  // Read in properties specific to this leg.
2050  StringSystem& systemNow = (legLoop == 0) ? systemMin : systemMid;
2051  vector<int>& iPartonNow = (legLoop == 0) ? iPartonMin : iPartonMid;
2052  int idPos = event[ iPartonNow[0] ].id();
2053  idOppose = event[ iPartonNow.back() ].id();
2054  double eInJRF = pInJRF[legNow].e();
2055  int statusHad = (legLoop == 0) ? 85 : 86;
2056 
2057  // Inner fallback loop, when a diquark comes in to junction.
2058  double eUsed = 0.;
2059  vector<StringVertex> junctionVertices;
2060  for ( int iTryInner = 0; ; ++iTryInner) {
2061 
2062  if (iTryInner > 2 * NTRYJNMATCH) {
2063  infoPtr->errorMsg("Error in StringFragmentation::fragment"
2064  "ToJunction: caught in junction flavour loop");
2065  event.popBack( iPartonMin.size() + iPartonMid.size() );
2066  return false;
2067  }
2068 
2069  bool needBaryon = (abs(idPos) > 10 && iTryInner > NTRYJNMATCH);
2070  double eExtra = (iTryInner > NTRYJNMATCH) ? EEXTRAJNMATCH : 0.;
2071 
2072  // Set up two string ends, and begin fragmentation loop.
2073  setStartEnds(idPos, idOppose, systemNow, legNow);
2074  eUsed = 0.;
2075  int nHadrons = 0;
2076  bool noNegE = true;
2077 
2078  // Keep track of hadron momentum.
2079  Vec4 hadMom;
2080 
2081  // Inform the UserHooks about the string to he hadronised.
2082  if ( userHooksPtr && userHooksPtr->canChangeFragPar() )
2083  userHooksPtr->setStringEnds(&posEnd, 0, iPartonNow);
2084 
2085  for ( ; ; ++nHadrons) {
2086 
2087  // The FlavourRope treatment changes the fragmentation parameters.
2088  if (flavRopePtr) {
2089  if (!flavRopePtr->doChangeFragPar(flavSelPtr, zSelPtr, pTSelPtr,
2090  hadMom.m2Calc(), (legLoop == 0 ? iPartonMin : iPartonMid ),
2091  idPos )) infoPtr->errorMsg("Error in StringFragmentation::"
2092  "fragmentToJunction: FlavourRope failed to change "
2093  "fragmentation parameters.");
2094  }
2095 
2096  // Possibility for a user to change the fragmentation parameters.
2097  if ( (userHooksPtr != 0) && userHooksPtr->canChangeFragPar() ) {
2098  if ( !userHooksPtr->doChangeFragPar( flavSelPtr, zSelPtr,
2099  pTSelPtr, idPos, hadMom.m2Calc(), iPartonNow, &posEnd) )
2100  infoPtr->errorMsg("Error in StringFragmentation::fragment"
2101  "ToJunction: failed to change hadronisation parameters.");
2102  }
2103 
2104  // Construct trial hadron from positive end.
2105  posEnd.newHadron();
2106  Vec4 pHad = posEnd.kinematicsHadron(systemNow, junctionVertices);
2107 
2108  // Possibility for a user to veto the hadron production.
2109  if ( (userHooksPtr != 0) && userHooksPtr->canChangeFragPar() ) {
2110  // Provide full particle info for veto decision.
2111  if ( userHooksPtr->doVetoFragmentation( Particle( posEnd.idHad,
2112  statusHad, iPos, iNeg, 0, 0, 0, 0,
2113  pHad, posEnd.mHad), &posEnd ) )
2114  continue;
2115  }
2116 
2117  // Negative energy signals failure in construction.
2118  if (pHad.e() < 0. ) { noNegE = false; break; }
2119 
2120  // Break if passed system midpoint ( = junction) in energy.
2121  // Exceptions: small systems, and/or with diquark end.
2122  bool delayedBreak = false;
2123  if (eUsed + pHad.e() + eExtra > eInJRF) {
2124  if (nHadrons > 0 || !needBaryon) {
2125  junctionVertices.pop_back();
2126  break;
2127  }
2128  delayedBreak = true;
2129  }
2130 
2131  // Else construct kinematics of the new hadron and store it.
2132  hadrons.append( posEnd.idHad, statusHad, iPos, iNeg,
2133  0, 0, posEnd.colOld, posEnd.colNew, pHad, posEnd.mHad);
2134 
2135  // Update hadron, string end and remaining momentum.
2136  hadMom += pHad;
2137  posEnd.update();
2138  eUsed += pHad.e();
2139 
2140  // Delayed break in small systems, and/or with diquark end.
2141  if (delayedBreak) {
2142  ++nHadrons;
2143  break;
2144  }
2145  }
2146 
2147  // Possible to produce zero hadrons if the endpoint is not a diquark.
2148  if (iTryInner > NTRYJNMATCH && !noNegE && nHadrons == 0 &&
2149  abs(idPos) < 10) break;
2150 
2151  // End of fragmentation loop. Inner loopback if ends on a diquark.
2152  if ( noNegE && abs(posEnd.flavOld.id) < 10 ) break;
2153  hadrons.popBack(nHadrons);
2154  junctionVertices.clear();
2155  if (legNow == legMin) legMinVertices.clear();
2156  else legMidVertices.clear();
2157  }
2158 
2159  // End of one-leg fragmentation. Store end quark and remnant energy.
2160  if (legNow == legMin) {
2161  idMin = posEnd.flavOld.id;
2162  eLeftMin = eInJRF - eUsed;
2163  } else {
2164  idMid = posEnd.flavOld.id;
2165  eLeftMid = eInJRF - eUsed;
2166  }
2167 
2168  // Store breakup vertices in two different vectors depending on leg.
2169  if (setVertices) {
2170  for (int i = 0; i < int(junctionVertices.size()); ++i) {
2171  if (legNow == legMin)
2172  legMinVertices.push_back( junctionVertices[i]);
2173  else legMidVertices.push_back( junctionVertices[i]);
2174  }
2175  }
2176  }
2177 
2178  // End of both-leg fragmentation.
2179  // Middle loopback if too much energy left.
2180  double eTrial = eBothLeftJunction + rndmPtr->flat() * eMaxLeftJunction;
2181  if (iTryMiddle > NTRYJNMATCH
2182  || ( min( eLeftMin, eLeftMid) < eBothLeftJunction
2183  && max( eLeftMin, eLeftMid) < eTrial ) ) break;
2184  hadrons.clear();
2185  legMinVertices.clear();
2186  legMidVertices.clear();
2187  }
2188 
2189  // Boost hadrons away from the JRF to the original frame.
2190  for (int i = 0; i < hadrons.size(); ++i) {
2191  hadrons[i].rotbst(MfromJRF);
2192 
2193  // Recalculate energy to compensate for numerical precision loss
2194  // in iterative calculation of MfromJRF.
2195  hadrons[i].e( hadrons[i].eCalc() );
2196  pJunctionHadrons += hadrons[i].p();
2197  }
2198 
2199  // Outer loopback if combined diquark mass too negative
2200  // or too little energy left in third leg.
2201  pDiquark = pInLeg[legMin] + pInLeg[legMid] - pJunctionHadrons;
2202  double m2Left = m2( pInLeg[legMax], pDiquark);
2203  if (iTryOuter > NTRYJNMATCH || (pDiquark.mCalc() > MDIQUARKMIN
2204  && m2Left > eMinLeftJunction * pInLeg[legMax].e()) ) break;
2205  hadrons.clear();
2206  legMinVertices.clear();
2207  legMidVertices.clear();
2208  pJunctionHadrons = 0.;
2209  }
2210 
2211  // Now found solution; no more loopback. Remove temporary parton copies.
2212  event.popBack( iPartonMin.size() + iPartonMid.size() );
2213 
2214  // Construct and store an effective diquark string end from the
2215  // two remnant quark ends, for temporary usage.
2216  idDiquark = flavSelPtr->makeDiquark( idMin, idMid);
2217  int iDiquark = event.append( idDiquark, 78, 0, 0, 0, 0, 0, 0,
2218  pDiquark, pDiquark.mCalc());
2219 
2220  // Find the partons on the last leg, again in reverse order.
2221  iPartonMax.clear();
2222  for (int i = legEnd[legMax]; i >= legBeg[legMax]; --i)
2223  iPartonMax.push_back( iParton[i] );
2224  iPartonMax.push_back( iDiquark );
2225 
2226  // Recluster gluons nearby to diquark end when taken too much energy.
2227  int iPsize = iPartonMax.size();
2228  double m0Diquark = event[iDiquark].m0();
2229  while (iPsize > 2) {
2230  Vec4 pGluNear = event[ iPartonMax[iPsize - 2] ].p();
2231  if ( (pDiquark + 0.5 * pGluNear).mCalc() > m0Diquark + mJoin ) break;
2232  pDiquark += pGluNear;
2233  event[iDiquark].p( pDiquark );
2234  event[iDiquark].m( pDiquark.mCalc() );
2235  iPartonMax.pop_back();
2236  --iPsize;
2237  iPartonMax[iPsize - 1] = iDiquark;
2238  }
2239  if ( idDiquark > 0 )
2240  event[iDiquark].acol(event[iPartonMax[iPsize - 2]].col());
2241  else
2242  event[iDiquark].col(event[iPartonMax[iPsize - 2]].acol());
2243 
2244  // Modify parton list to remaining leg + remnant of the first two.
2245  iParton = iPartonMax;
2246 
2247  // Done.
2248  return true;
2249 }
2250 
2251 //--------------------------------------------------------------------------
2252 
2253 // Find the boost matrix to the rest frame of a junction,
2254 // given the three respective endpoint four-momenta.
2255 
2256 RotBstMatrix StringFragmentation::junctionRestFrame(Vec4& p0, Vec4& p1,
2257  Vec4& p2) {
2258 
2259  // Calculate masses and other invariants.
2260  Vec4 pSumJun = p0 + p1 + p2;
2261  double sHat = pSumJun.m2Calc();
2262  double pp[3][3];
2263  pp[0][0] = p0.m2Calc();
2264  pp[1][1] = p1.m2Calc();
2265  pp[2][2] = p2.m2Calc();
2266  pp[0][1] = pp[1][0] = p0 * p1;
2267  pp[0][2] = pp[2][0] = p0 * p2;
2268  pp[1][2] = pp[2][1] = p1 * p2;
2269 
2270  // Requirement (eiMax)_j = pi*pj/mj < (eiMax)_k = pi*pk/mk, used below,
2271  // here rewritten as pi*pj * mk < pi*pk * mj and squared.
2272  double eMax01 = pow2(pp[0][1]) * pp[2][2];
2273  double eMax02 = pow2(pp[0][2]) * pp[1][1];
2274  double eMax12 = pow2(pp[1][2]) * pp[0][0];
2275 
2276  // Initially pick i to be the most massive parton. but allow other tries.
2277  int i = (pp[1][1] > pp[0][0]) ? 1 : 0;
2278  if (pp[2][2] > max(pp[0][0], pp[1][1])) i = 2;
2279  int j, k;
2280  double ei = 0.;
2281  double ej = 0.;
2282  double ek = 0.;
2283  for (int iTry = 0; iTry < 3; ++iTry) {
2284 
2285  // Pick j to give minimal eiMax, and k the third vector.
2286  if (i == 0) j = (eMax02 < eMax01) ? 2 : 1;
2287  else if (i == 1) j = (eMax12 < eMax01) ? 2 : 0;
2288  else j = (eMax12 < eMax02) ? 1 : 0;
2289  k = 3 - i - j;
2290 
2291  // Alternative names according to i, j, k conventions.
2292  double m2i = pp[i][i];
2293  double m2j = pp[j][j];
2294  double m2k = pp[k][k];
2295  double pipj = pp[i][j];
2296  double pipk = pp[i][k];
2297  double pjpk = pp[j][k];
2298 
2299  // Trivial to find new parton energies if all three partons are massless.
2300  if (m2i < M2MAXJRF) {
2301  ei = sqrt( 2. * pipk * pipj / (3. * pjpk) );
2302  ej = sqrt( 2. * pjpk * pipj / (3. * pipk) );
2303  ek = sqrt( 2. * pipk * pjpk / (3. * pipj) );
2304 
2305  // Else find three-momentum range for parton i and values at extremes.
2306  } else {
2307  // Minimum when i is at rest.
2308  double piMin = 0.;
2309  double eiMin = sqrt(m2i);
2310  double ejMin = pipj / eiMin;
2311  double ekMin = pipk / eiMin;
2312  double pjMin = sqrtpos( ejMin*ejMin - m2j );
2313  double pkMin = sqrtpos( ekMin*ekMin - m2k );
2314  double fMin = ejMin * ekMin + 0.5 * pjMin * pkMin - pjpk;
2315  // Maximum estimated when j + k is at rest, alternatively j at rest.
2316  double eiMax = (pipj + pipk)
2317  / sqrt( max( M2MINJRF, m2j + m2k + 2. * pjpk) );
2318  if (m2j > M2MAXJRF) eiMax = min( eiMax, pipj / sqrt(m2j) );
2319  double piMax = sqrtpos( eiMax*eiMax - m2i );
2320  double temp = eiMax*eiMax - 0.25 *piMax*piMax;
2321  double pjMax = (eiMax * sqrtpos( pipj*pipj - m2j * temp )
2322  - 0.5 * piMax * pipj) / temp;
2323  double pkMax = (eiMax * sqrtpos( pipk*pipk - m2k * temp )
2324  - 0.5 * piMax * pipk) / temp;
2325  double ejMax = sqrt(pjMax*pjMax + m2j);
2326  double ekMax = sqrt(pkMax*pkMax + m2k);
2327  double fMax = ejMax * ekMax + 0.5 * pjMax * pkMax - pjpk;
2328 
2329  // If unexpected values at upper endpoint then pick another parton.
2330  if (fMax > 0.) {
2331  int iPrel = (i + 1)%3;
2332  if (pp[iPrel][iPrel] > M2MAXJRF) {i = iPrel; continue;}
2333  ++iTry;
2334  iPrel = (i + 2)%3;
2335  if (iTry < 3 && pp[iPrel][iPrel] > M2MAXJRF) {i = iPrel; continue;}
2336  }
2337 
2338  // Start binary + linear search to find solution inside range.
2339  int iterMin = 0;
2340  int iterMax = 0;
2341  double pi = 0.5 * (piMin + piMax);
2342  for (int iter = 0; iter < NTRYJRFEQ; ++iter) {
2343 
2344  // Derive momentum of other two partons and distance to root.
2345  ei = sqrt(pi*pi + m2i);
2346  temp = ei*ei - 0.25 * pi*pi;
2347  double pj = (ei * sqrtpos( pipj*pipj - m2j * temp )
2348  - 0.5 * pi * pipj) / temp;
2349  double pk = (ei * sqrtpos( pipk*pipk - m2k * temp )
2350  - 0.5 * pi * pipk) / temp;
2351  ej = sqrt(pj*pj + m2j);
2352  ek = sqrt(pk*pk + m2k);
2353  double fNow = ej * ek + 0.5 * pj * pk - pjpk;
2354 
2355  // Replace lower or upper bound by new value.
2356  if (fNow > 0.) { ++iterMin; piMin = pi; fMin = fNow;}
2357  else {++iterMax; piMax = pi; fMax = fNow;}
2358 
2359  // Pick next i momentum to explore, hopefully closer to root.
2360  if (2 * iter < NTRYJRFEQ
2361  && (iterMin < 2 || iterMax < 2 || 4 * iter < NTRYJRFEQ))
2362  { pi = 0.5 * (piMin + piMax); continue;}
2363  if (fMin < 0. || fMax > 0. || abs(fNow) < CONVJRFEQ * sHat) break;
2364  pi = piMin + (piMax - piMin) * fMin / (fMin - fMax);
2365  }
2366 
2367  // If arrived here then either succeeded or exhausted possibilities.
2368  } break;
2369  }
2370 
2371  // Now we know the energies in the junction rest frame.
2372  double eNew[3] = { 0., 0., 0.};
2373  eNew[i] = ei;
2374  eNew[j] = ej;
2375  eNew[k] = ek;
2376 
2377  // Boost (copy of) partons to their rest frame.
2378  RotBstMatrix Mmove;
2379  Vec4 p0cm = p0;
2380  Vec4 p1cm = p1;
2381  Vec4 p2cm = p2;
2382  Mmove.bstback(pSumJun);
2383  p0cm.rotbst(Mmove);
2384  p1cm.rotbst(Mmove);
2385  p2cm.rotbst(Mmove);
2386 
2387  // Construct difference vectors and the boost to junction rest frame.
2388  Vec4 pDir01 = p0cm / p0cm.e() - p1cm / p1cm.e();
2389  Vec4 pDir02 = p0cm / p0cm.e() - p2cm / p2cm.e();
2390  double pDiff01 = pDir01.pAbs2();
2391  double pDiff02 = pDir02.pAbs2();
2392  double pDiff0102 = dot3(pDir01, pDir02);
2393  double eDiff01 = eNew[0] / p0cm.e() - eNew[1] / p1cm.e();
2394  double eDiff02 = eNew[0] / p0cm.e() - eNew[2] / p2cm.e();
2395  double denom = pDiff01 * pDiff02 - pDiff0102*pDiff0102;
2396  double coef01 = (eDiff01 * pDiff02 - eDiff02 * pDiff0102) / denom;
2397  double coef02 = (eDiff02 * pDiff01 - eDiff01 * pDiff0102) / denom;
2398  Vec4 vJunction = coef01 * pDir01 + coef02 * pDir02;
2399  vJunction.e( sqrt(1. + vJunction.pAbs2()) );
2400 
2401  // Add two boosts, giving final result.
2402  Mmove.bst(vJunction);
2403  return Mmove;
2404 
2405 }
2406 
2407 //--------------------------------------------------------------------------
2408 
2409 // When string fragmentation has failed several times,
2410 // try to join some more nearby partons.
2411 
2412 int StringFragmentation::extraJoin(double facExtra, Event& event) {
2413 
2414  // Keep on looping while pairs found below joining threshold.
2415  int nJoin = 0;
2416  int iPsize = iParton.size();
2417  while (iPsize > 2) {
2418 
2419  // Look for the pair of neighbour partons (along string) with
2420  // the smallest invariant mass (subtracting quark masses).
2421  int iJoinMin = -1;
2422  double mJoinMin = 2. * facExtra * mJoin;
2423  for (int i = 0; i < iPsize - 1; ++i) {
2424  Particle& parton1 = event[ iParton[i] ];
2425  Particle& parton2 = event[ iParton[i + 1] ];
2426  Vec4 pSumNow;
2427  pSumNow += (parton2.isGluon()) ? 0.5 * parton1.p() : parton1.p();
2428  pSumNow += (parton2.isGluon()) ? 0.5 * parton2.p() : parton2.p();
2429  double mJoinNow = pSumNow.mCalc();
2430  if (!parton1.isGluon()) mJoinNow -= parton1.m0();
2431  if (!parton2.isGluon()) mJoinNow -= parton2.m0();
2432  if (mJoinNow < mJoinMin) { iJoinMin = i; mJoinMin = mJoinNow; }
2433  }
2434 
2435  // Decide whether to join, if not finished.
2436  if (iJoinMin < 0 || mJoinMin > facExtra * mJoin) return nJoin;
2437  ++nJoin;
2438 
2439  // Create new joined parton.
2440  int iJoin1 = iParton[iJoinMin];
2441  int iJoin2 = iParton[iJoinMin + 1];
2442  int idNew = (event[iJoin1].isGluon()) ? event[iJoin2].id()
2443  : event[iJoin1].id();
2444  int colNew = event[iJoin1].col();
2445  int acolNew = event[iJoin2].acol();
2446  if (colNew == acolNew) {
2447  colNew = event[iJoin2].col();
2448  acolNew = event[iJoin1].acol();
2449  }
2450  Vec4 pNew = event[iJoin1].p() + event[iJoin2].p();
2451 
2452  // Append joined parton to event record and reduce parton list.
2453  int iNew = event.append( idNew, 73, min(iJoin1, iJoin2),
2454  max(iJoin1, iJoin2), 0, 0, colNew, acolNew, pNew, pNew.mCalc() );
2455  iParton[iJoinMin] = iNew;
2456  for (int i = iJoinMin + 1; i < iPsize - 1; ++i)
2457  iParton[i] = iParton[i + 1];
2458  iParton.pop_back();
2459  --iPsize;
2460 
2461  // Done.
2462  }
2463  return nJoin;
2464 }
2465 
2466 //--------------------------------------------------------------------------
2467 
2468 // When the next hadron is generated the temperature can be dependent on the
2469 // number of nearby string pieces with respect to the string piece the hadron
2470 // will be produced on.
2471 
2472 double StringFragmentation::nearStringPieces(StringEnd end,
2473  vector< vector< pair<double,double> > >& rapPairs) {
2474 
2475  // No modification for junctions.
2476  if (hasJunction) return 1;
2477 
2478  // Get temporary hadron momentum.
2479  double phi = 2.0 * M_PI * rndmPtr->flat();
2480  double mult = -1.0;
2481  int nTryMax = 100;
2482  double multStep = 5.0 / ((double)nTryMax/2);
2483  double multNow = 1.0 + multStep;
2484  Vec4 pHad = Vec4(0., 0., 0., -1.);
2485  for (int i = 1; i <= nTryMax; i++) {
2486  pHad = end.kinematicsHadronTmp(system, pRem, phi, mult);
2487  // If valid momentum found, done.
2488  if (pHad.e() > 0.0) break;
2489  // Set mult as multiplicative factor. Alternate between adding and
2490  // subtracting multStep.
2491  mult = 1.0;
2492  if (i%2 == 0) {
2493  mult *= multNow;
2494  multNow += multStep;
2495  } else mult /= multNow;
2496  }
2497 
2498  // In case of failure, use remnant momentum.
2499  if (pHad.e() < 0.0) pHad = pRem;
2500 
2501  // Now loop through the list of rapidity pairs and count strings
2502  // sitting at the hadron rapidity.
2503  Particle hadron = Particle();
2504  hadron.p(pHad); hadron.m(pHad.mCalc());
2505  double yHad = hadron.y();
2506  int nString = -1;
2507  for (int iSub = 0; iSub < int(rapPairs.size()); iSub++) {
2508  vector< pair<double,double> > pairNow = rapPairs[iSub];
2509  for (int iPair = 0; iPair < int(pairNow.size()); iPair++) {
2510  double y1 = pairNow[iPair].first;
2511  double y2 = pairNow[iPair].second;
2512  if ( (y1 < yHad) && (yHad < y2) ) nString++;
2513  }
2514  }
2515  // Effective number of strings takes pT into account.
2516  double pT2Had = pHad.pT2();
2517  double nStringEff = double(nString) / (1.0 + pT2Had / pT20);
2518 
2519  return nStringEff + 1.0;
2520 }
2521 
2522 //==========================================================================
2523 
2524 } // end namespace Pythia8
Definition: AgUStep.h:26