StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
RHadrons.cc
1 // RHadrons.cc is a part of the PYTHIA event generator.
2 // Copyright (C) 2012 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL version 2, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // Function definitions (not found in the header) for the RHadrons class.
7 
8 #include "RHadrons.h"
9 
10 namespace Pythia8 {
11 
12 //==========================================================================
13 
14 // The RHadrons class.
15 
16 //--------------------------------------------------------------------------
17 
18 // Constants: could be changed here if desired, but normally should not.
19 // These are of technical nature, as described for each.
20 
21 const int RHadrons::IDRHADSB[14] = { 1000512, 1000522, 1000532,
22  1000542, 1000552, 1005113, 1005211, 1005213, 1005223, 1005311,
23  1005313, 1005321, 1005323, 1005333 };
24 
25 const int RHadrons::IDRHADST[14] = { 1000612, 1000622, 1000632,
26  1000642, 1000652, 1006113, 1006211, 1006213, 1006223, 1006311,
27  1006313, 1006321, 1006323, 1006333 };
28 
29 // Gluino code and list of gluino R-hadron codes.
30 const int RHadrons::IDRHADGO[38] = { 1000993, 1009113, 1009213,
31  1009223, 1009313, 1009323, 1009333, 1009413, 1009423, 1009433,
32  1009443, 1009513, 1009523, 1009533, 1009543, 1009553, 1091114,
33  1092114, 1092214, 1092224, 1093114, 1093214, 1093224, 1093314,
34  1093324, 1093334, 1094114, 1094214, 1094224, 1094314, 1094324,
35  1094334, 1095114, 1095214, 1095224, 1095314, 1095324, 1095334 };
36 
37 // Allow a few tries for flavour selection since failure is allowed.
38 const int RHadrons::NTRYMAX = 10;
39 
40 // Safety margin (in GeV) when constructing kinematics of system.
41 const double RHadrons::MSAFETY = 0.1;
42 
43 // Maximal energy to borrow for gluon to insert on leg in to junction.
44 const double RHadrons::EGBORROWMAX = 4.;
45 
46 //--------------------------------------------------------------------------
47 
48 // Main routine to initialize R-hadron handling.
49 
50 bool RHadrons::init( Info* infoPtrIn, Settings& settings,
51  ParticleData* particleDataPtrIn, Rndm* rndmPtrIn) {
52 
53  // Store input pointers for future use.
54  infoPtr = infoPtrIn;
55  particleDataPtr = particleDataPtrIn;
56  rndmPtr = rndmPtrIn;
57 
58  // Flags and parameters related to R-hadron formation and decay.
59  allowRH = settings.flag("RHadrons:allow");
60  maxWidthRH = settings.parm("RHadrons:maxWidth");
61  idRSb = settings.mode("RHadrons:idSbottom");
62  idRSt = settings.mode("RHadrons:idStop");
63  idRGo = settings.mode("RHadrons:idGluino");
64  setMassesRH = settings.flag("RHadrons:setMasses");
65  probGluinoballRH = settings.parm("RHadrons:probGluinoball");
66  mOffsetCloudRH = settings.parm("RHadrons:mOffsetCloud");
67  mCollapseRH = settings.parm("RHadrons:mCollapse");
68  diquarkSpin1RH = settings.parm("RHadrons:diquarkSpin1");
69 
70  // Check whether sbottom, stop or gluino should form R-hadrons.
71  allowRSb = allowRH && idRSb > 0
72  && (particleDataPtr->mWidth(idRSb) < maxWidthRH);
73  allowRSt = allowRH && idRSt > 0
74  && (particleDataPtr->mWidth(idRSt) < maxWidthRH);
75  allowRGo = allowRH && idRGo > 0
76  && (particleDataPtr->mWidth(idRGo) < maxWidthRH);
77  allowSomeR = allowRSb || allowRSt || allowRGo;
78 
79  // Set nominal masses of sbottom R-mesons and R-baryons.
80  if (allowRSb) {
81  m0Sb = particleDataPtr->m0(idRSb);
82  if (setMassesRH) {
83  for (int i = 0; i < 14; ++i) {
84  int idR = IDRHADSB[i];
85  double m0RHad = m0Sb + mOffsetCloudRH;
86  m0RHad += particleDataPtr->constituentMass( (idR%100)/10);
87  if (i > 4)
88  m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
89  particleDataPtr->m0( idR, m0RHad);
90  }
91  }
92 
93  // Set widths and lifetimes of sbottom R-hadrons.
94  double mWidthRHad = particleDataPtr->mWidth(idRSb);
95  double tau0RHad = particleDataPtr->tau0( idRSb);
96  for (int i = 0; i < 14; ++i) {
97  particleDataPtr->mWidth( IDRHADSB[i], mWidthRHad);
98  particleDataPtr->tau0( IDRHADSB[i], tau0RHad);
99  }
100  }
101 
102  // Set nominal masses of stop R-mesons and R-baryons.
103  if (allowRSt) {
104  m0St = particleDataPtr->m0(idRSt);
105  if (setMassesRH) {
106  for (int i = 0; i < 14; ++i) {
107  int idR = IDRHADST[i];
108  double m0RHad = m0St + mOffsetCloudRH;
109  m0RHad += particleDataPtr->constituentMass( (idR%100)/10);
110  if (i > 4)
111  m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
112  particleDataPtr->m0( idR, m0RHad);
113  }
114  }
115 
116  // Set widths and lifetimes of stop R-hadrons.
117  double mWidthRHad = particleDataPtr->mWidth(idRSt);
118  double tau0RHad = particleDataPtr->tau0( idRSt);
119  for (int i = 0; i < 14; ++i) {
120  particleDataPtr->mWidth( IDRHADST[i], mWidthRHad);
121  particleDataPtr->tau0( IDRHADST[i], tau0RHad);
122  }
123  }
124 
125  // Set nominal masses of gluino R-glueballs, R-mesons and R-baryons.
126  if (allowRGo) {
127  m0Go = particleDataPtr->m0(idRGo);
128  if (setMassesRH) {
129  particleDataPtr->m0( IDRHADGO[0], m0Go + 2. * mOffsetCloudRH
130  + particleDataPtr->constituentMass(21) );
131  for (int i = 1; i < 38; ++i) {
132  int idR = IDRHADGO[i];
133  double m0RHad = m0Go + 2. * mOffsetCloudRH;
134  m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
135  m0RHad += particleDataPtr->constituentMass( (idR%100)/10);
136  if (i > 15)
137  m0RHad += particleDataPtr->constituentMass( (idR%10000)/1000 );
138  particleDataPtr->m0( idR, m0RHad);
139  }
140  }
141 
142  // Set widths and lifetimes of gluino R-hadrons.
143  double mWidthRHad = particleDataPtr->mWidth(idRGo);
144  double tau0RHad = particleDataPtr->tau0( idRGo);
145  for (int i = 0; i < 38; ++i) {
146  particleDataPtr->mWidth( IDRHADGO[i], mWidthRHad);
147  particleDataPtr->tau0( IDRHADGO[i], tau0RHad);
148  }
149  }
150 
151  // Done.
152  return true;
153 
154 }
155 //--------------------------------------------------------------------------
156 
157 // Check if a given particle can produce R-hadrons.
158 
159 bool RHadrons::givesRHadron( int id) {
160  if (allowRSb && abs(id) == idRSb) return true;
161  if (allowRSt && abs(id) == idRSt) return true;
162  if (allowRGo && id == idRGo) return true;
163  return false;
164 
165 }
166 
167 //--------------------------------------------------------------------------
168 
169 // Produce R-hadrons by fragmenting them off from existing strings.
170 
171 bool RHadrons::produce( ColConfig& colConfig, Event& event) {
172 
173  // Check whether some sparticles are allowed to hadronize.
174  if (!allowSomeR) return true;
175 
176  // Reset arrays and values.
177  iBefRHad.resize(0);
178  iCreRHad.resize(0);
179  iRHadron.resize(0);
180  iAftRHad.resize(0);
181  isTriplet.resize(0);
182  nRHad = 0;
183 
184  // Loop over event and identify hadronizing sparticles.
185  for (int i = 0; i < event.size(); ++i)
186  if (event[i].isFinal() && givesRHadron(event[i].id())) {
187  iBefRHad.push_back(i);
188  iCreRHad.push_back(i);
189  iRHadron.push_back(0);
190  iAftRHad.push_back(0);
191  isTriplet.push_back(true);
192  }
193  nRHad = iRHadron.size();
194 
195  // Done if no hadronizing sparticles.
196  if (nRHad == 0) return true;
197 
198  // Max two R-hadrons. Randomize order of processing.
199  if (nRHad > 2) {
200  infoPtr->errorMsg("Error in RHadrons::produce: "
201  "cannot handle more than two R-hadrons");
202  return false;
203  }
204  if (nRHad > 1 && rndmPtr->flat() > 0.5) swap( iBefRHad[0], iBefRHad[1]);
205 
206  // Split a system with both a sparticle and a junction.
207  iBef = iBefRHad[0];
208  iSys = colConfig.findSinglet( iBef);
209  systemPtr = &colConfig[iSys];
210  if (systemPtr->hasJunction && !splitOffJunction( colConfig, event)) {
211  infoPtr->errorMsg("Error in RHadrons::produce: "
212  "cannot handle system with junction");
213  return false;
214  }
215  if (nRHad == 2) {
216  iBef = iBefRHad[1];
217  iSys = colConfig.findSinglet( iBefRHad[1]);
218  systemPtr = &colConfig[iSys];
219  if (systemPtr->hasJunction && !splitOffJunction( colConfig, event)) {
220  infoPtr->errorMsg("Error in RHadrons::produce: "
221  "cannot handle system with junction");
222  return false;
223  }
224  }
225 
226  // Open up a closed gluon/gluino loop.
227  iBef = iBefRHad[0];
228  iSys = colConfig.findSinglet( iBef);
229  systemPtr = &colConfig[iSys];
230  if (systemPtr->isClosed && !openClosedLoop( colConfig, event)) {
231  infoPtr->errorMsg("Error in RHadrons::produce: "
232  "cannot open up closed gluon/gluino loop");
233  return false;
234  }
235  if (nRHad == 2) {
236  iBef = iBefRHad[1];
237  iSys = colConfig.findSinglet( iBefRHad[1]);
238  systemPtr = &colConfig[iSys];
239  if (systemPtr->isClosed && !openClosedLoop( colConfig, event)) {
240  infoPtr->errorMsg("Error in RHadrons::produce: "
241  "cannot open up closed gluon/gluino loop");
242  return false;
243  }
244  }
245 
246  // Split up a colour singlet system that should give two R-hadrons.
247  if (nRHad == 2) {
248  int iSys1 = colConfig.findSinglet( iBefRHad[0]);
249  int iSys2 = colConfig.findSinglet( iBefRHad[1]);
250  if (iSys2 == iSys1) {
251  iSys = iSys1;
252  systemPtr = &colConfig[iSys];
253  if ( !splitSystem( colConfig, event) ) {
254  infoPtr->errorMsg("Error in RHadrons::produce: "
255  "failed to handle two sparticles in same system");
256  return false;
257  }
258  }
259  }
260 
261  // Loop over R-hadrons to be formed. Find its colour singlet system.
262  for (iRHad = 0; iRHad < nRHad; ++iRHad) {
263  iBef = iBefRHad[iRHad];
264  iSys = colConfig.findSinglet( iBef);
265  if (iSys < 0) {
266  infoPtr->errorMsg("Error in RHadrons::produce: "
267  "sparticle not in any colour singlet");
268  return false;
269  }
270  systemPtr = &colConfig[iSys];
271 
272  // For now don't handle systems involving junctions or loops.
273  if (systemPtr->hasJunction) {
274  infoPtr->errorMsg("Error in RHadrons::produce: "
275  "cannot handle system with junction");
276  return false;
277  }
278  if (systemPtr->isClosed) {
279  infoPtr->errorMsg("Error in RHadrons::produce: "
280  "cannot handle closed colour loop");
281  return false;
282  }
283 
284  // Handle formation of R-hadron separately for gluino and squark.
285  if (event[iBef].id() == idRGo) isTriplet[iRHad] = false;
286  bool formed = (isTriplet[iRHad]) ? produceSquark( colConfig, event)
287  : produceGluino( colConfig, event);
288  if (!formed) return false;
289 
290  // End of loop over R-hadrons. Done.
291  }
292  return true;
293 
294 }
295 
296 //--------------------------------------------------------------------------
297 
298 // Decay R-hadrons by resolving them into string systems and letting the
299 // heavy unstable particle decay as normal.
300 
301 bool RHadrons::decay( Event& event) {
302 
303  // Loop over R-hadrons to decay.
304  for (iRHad = 0; iRHad < nRHad; ++iRHad) {
305  int iRNow = iRHadron[iRHad];
306  int iRBef = iBefRHad[iRHad];
307  int idRHad = event[iRNow].id();
308  double mRHad = event[iRNow].m();
309  double mRBef = event[iRBef].m();
310  int iR0 = 0;
311  int iR2 = 0;
312 
313  // Find flavour content of squark or gluino R-hadron.
314  pair<int,int> idPair = (isTriplet[iRHad])
315  ? fromIdWithSquark( idRHad) : fromIdWithGluino( idRHad);
316  int id1 = idPair.first;
317  int id2 = idPair.second;
318 
319  // Sharing of momentum: the squark/gluino should be restored
320  // to original mass, but error if negative-mass spectators.
321  double fracR = mRBef / mRHad;
322  if (fracR >= 1.) {
323  infoPtr->errorMsg("Error in RHadrons::decay: "
324  "too low R-hadron mass for decay");
325  return false;
326  }
327 
328  // Squark: new colour needed in the breakup.
329  if (isTriplet[iRHad]) {
330  int colNew = event.nextColTag();
331  int col = (event[iRBef].col() != 0) ? colNew : 0;
332  int acol = (col == 0) ? colNew : 0;
333 
334  // Store the constituents of a squark R-hadron.
335  iR0 = event.append( id1, 106, iRNow, 0, 0, 0, col, acol,
336  fracR * event[iRNow].p(), fracR * mRHad, 0.);
337  iR2 = event.append( id2, 106, iRNow, 0, 0, 0, acol, col,
338  (1. - fracR) * event[iRNow].p(), (1. - fracR) * mRHad, 0.);
339 
340  // Gluino: set mass sharing between two spectators.
341  } else {
342  double m1Eff = particleDataPtr->constituentMass(id1) + mOffsetCloudRH;
343  double m2Eff = particleDataPtr->constituentMass(id2) + mOffsetCloudRH;
344  double frac1 = (1. - fracR) * m1Eff / ( m1Eff + m2Eff);
345  double frac2 = (1. - fracR) * m2Eff / ( m1Eff + m2Eff);
346 
347  // Two new colours needed in the breakups.
348  int col1 = event.nextColTag();
349  int col2 = event.nextColTag();
350 
351  // Store the constituents of a gluino R-hadron.
352  iR0 = event.append( idRGo, 106, iRNow, 0, 0, 0, col2, col1,
353  fracR * event[iRNow].p(), fracR * mRHad, 0.);
354  event.append( id1, 106, iRNow, 0, 0, 0, col1, 0,
355  frac1 * event[iRNow].p(), frac1 * mRHad, 0.);
356  iR2 = event.append( id2, 106, iRNow, 0, 0, 0, 0, col2,
357  frac2 * event[iRNow].p(), frac2 * mRHad, 0.);
358  }
359 
360  // Mark R-hadron as decayed and update history.
361  event[iRNow].statusNeg();
362  event[iRNow].daughters( iR0, iR2);
363  iAftRHad[iRHad] = iR0;
364 
365  // Set secondary vertex for decay products, but no lifetime.
366  Vec4 vDec = event[iRNow].vProd() + event[iRNow].tau()
367  * event[iR0].p() / event[iR0].m();
368  for (int iRd = iR0; iRd <= iR2; ++iRd) event[iRd].vProd( vDec);
369 
370  // End loop over R-hadron decays, based on velocity of squark.
371 
372  }
373 
374  // Done.
375  return true;
376 
377 }
378 
379 //--------------------------------------------------------------------------
380 
381 // Split a system that contains both a sparticle and a junction.
382 
383 bool RHadrons::splitOffJunction( ColConfig& colConfig, Event& event) {
384 
385  // Classify system into three legs, and find sparticle location.
386  vector<int> leg1, leg2, leg3;
387  int iLegSP = 0;
388  int iPosSP = 0;
389  int iLeg = 0;
390  int iPos = 0;
391  for (int i = 0; i < systemPtr->size(); ++i) {
392  ++iPos;
393  int iP = systemPtr->iParton[i];
394  if (iP < 0) {
395  ++iLeg;
396  iPos = -1;
397  } else if (iLeg == 1) leg1.push_back( iP);
398  else if (iLeg == 2) leg2.push_back( iP);
399  else if (iLeg == 3) leg3.push_back( iP);
400  if (iP == iBef) {
401  iLegSP = iLeg;
402  iPosSP = iPos;
403  }
404  }
405  if (iLegSP == 0) return false;
406 
407  // Swap so leg 1 contains sparticle. If not innermost sparticle then
408  // skip for now, and wait for this other (gluino!) to be split off.
409  if (iLegSP == 2) swap(leg2, leg1);
410  else if (iLegSP == 3) swap(leg3, leg1);
411  for (int i = 0; i < iPosSP; ++i)
412  if (event[leg1[i]].id() != 21) return true;
413 
414  // No gluon between sparticle and junction: find kinetic energy of system.
415  if (iPosSP == 0) {
416  Vec4 pSP = event[iBef].p();
417  Vec4 pRec = 0.;
418  for (int i = 0; i < int(leg2.size()); ++i) pRec += event[leg2[i]].p();
419  for (int i = 0; i < int(leg3.size()); ++i) pRec += event[leg3[i]].p();
420  double mSys = (pSP + pRec).mCalc();
421  double mSP = pSP.mCalc();
422  double mRec = pRec.mCalc();
423  double eKin = mSys - mSP - mRec;
424 
425  // Insert new gluon that borrows part of kinetic energy.
426  double mNewG = EGBORROWMAX * eKin / (EGBORROWMAX + eKin) ;
427  Vec4 pNewG = (mNewG / mSys) * (pSP + pRec);
428  int newCol = event.nextColTag();
429  bool isCol = (event[leg1.back()].col() > 0);
430  int colNG = (isCol)? newCol : event[iBef].acol();
431  int acolNG = (isCol) ? event[iBef].col() : newCol;
432  int iNewG = event.append( 21, 101, iBef, 0, 0, 0, colNG, acolNG,
433  pNewG, mNewG, 0.);
434  leg1.insert( leg1.begin(), iNewG);
435  ++iPosSP;
436 
437  // Boosts for remainder systems that gave up energy.
438  double mNewSys = mSys - mNewG;
439  double pAbsOld = 0.5 * sqrtpos( pow2(mSys*mSys - mSP*mSP - mRec*mRec)
440  - pow2(2. * mSP * mRec) ) / mSys;
441  double pAbsNew = 0.5 * sqrtpos( pow2(mNewSys*mNewSys - mSP*mSP - mRec*mRec)
442  - pow2(2. * mSP * mRec) ) / mNewSys;
443  RotBstMatrix MtoCM, MfromCM, MSP, MRec;
444  MtoCM.toCMframe( pSP, pRec);
445  MfromCM = MtoCM;
446  MfromCM.invert();
447  MSP = MtoCM;
448  MSP.bst( 0., 0., -pAbsOld / sqrt(mSP * mSP + pAbsOld * pAbsOld));
449  MSP.bst( 0., 0., pAbsNew / sqrt(mSP * mSP + pAbsNew * pAbsNew));
450  MSP.rotbst( MfromCM);
451  MRec = MtoCM;
452  MRec.bst( 0., 0., pAbsOld / sqrt(mRec * mRec + pAbsOld * pAbsOld));
453  MRec.bst( 0., 0., -pAbsNew / sqrt(mRec * mRec + pAbsNew * pAbsNew));
454  MRec.rotbst( MfromCM);
455 
456  // Copy down recoling partons and boost their momenta.
457  int iNewSP = event.copy( iBef, 101);
458  event[iNewSP].rotbst( MSP);
459  leg1[iPosSP] = iNewSP;
460  if (iBefRHad[0] == iBef) iBefRHad[0] = iNewSP;
461  else if (nRHad > 1 && iBefRHad[1] == iBef) iBefRHad[1] = iNewSP;
462  iBef = iNewSP;
463  for (int i = 0; i < int(leg2.size()); ++i) {
464  int iCopy = event.copy( leg2[i], 101);
465  event[iCopy].rotbst( MRec);
466  if (iBefRHad[0] == leg2[i]) iBefRHad[0] = iCopy;
467  else if (nRHad > 1 && iBefRHad[1] == leg2[i]) iBefRHad[1] = iCopy;
468  leg2[i] = iCopy;
469  }
470  for (int i = 0; i < int(leg3.size()); ++i) {
471  int iCopy = event.copy( leg3[i], 101);
472  event[iCopy].rotbst( MRec);
473  if (iBefRHad[0] == leg3[i]) iBefRHad[0] = iCopy;
474  else if (nRHad > 1 && iBefRHad[1] == leg3[i]) iBefRHad[1] = iCopy;
475  leg3[i] = iCopy;
476  }
477 
478  // Now always at least one gluon between sparticle and junction.
479  }
480 
481  // Find gluon with largest energy in sparticle rest frame.
482  int iGspl = 0;
483  double eGspl = event[leg1[0]].p() * event[iBef].p();
484  for (int i = 1; i < iPosSP; ++i) {
485  double eGnow = event[leg1[i]].p() * event[iBef].p();
486  if (eGnow > eGspl) {
487  iGspl = i;
488  eGspl = eGnow;
489  }
490  }
491  int iG = leg1[iGspl];
492 
493  // Split this gluon into a collinear quark.antiquark pair.
494  int idNewQ = flavSelPtr->pickLightQ();
495  int iNewQ = event.append( idNewQ, 101, iG, 0, 0, 0, event[iG].col(), 0,
496  0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
497  int iNewQb = event.append( -idNewQ, 101, iG, 0, 0, 0, 0, event[iG].acol(),
498  0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
499  event[iG].statusNeg();
500  event[iG].daughters( iNewQ, iNewQb);
501  if (event[leg1.back()].col() == 0) swap( iNewQ, iNewQb);
502 
503  // Collect two new systems after split.
504  vector<int> iNewSys1, iNewSys2;
505  iNewSys1.push_back( iNewQb);
506  for (int i = iGspl + 1; i < int(leg1.size()); ++i)
507  iNewSys1.push_back( leg1[i]);
508  iNewSys2.push_back( -10);
509  for (int i = 0; i < iGspl; ++i) iNewSys2.push_back( leg1[i]);
510  iNewSys2.push_back( iNewQ);
511  iNewSys2.push_back( -11);
512  for (int i = 0; i < int(leg2.size()); ++i) iNewSys2.push_back( leg2[i]);
513  iNewSys2.push_back( -12);
514  for (int i = 0; i < int(leg3.size()); ++i) iNewSys2.push_back( leg3[i]);
515 
516  // Remove old system and insert two new ones.
517  colConfig.erase(iSys);
518  colConfig.insert( iNewSys1, event);
519  colConfig.insert( iNewSys2, event);
520 
521  // Done.
522  return true;
523 
524 }
525 
526 //--------------------------------------------------------------------------
527 
528 // Open up a closed gluon/gluino loop.
529 
530 bool RHadrons::openClosedLoop( ColConfig& colConfig, Event& event) {
531 
532  // Find gluon with largest energy in gluino rest frame.
533  int iGspl = -1;
534  double eGspl = 0.;
535  for (int i = 0; i < systemPtr->size(); ++i) {
536  int iGNow = systemPtr->iParton[i];
537  if (event[iGNow].id() == 21) {
538  double eGnow = event[iGNow].p() * event[iBef].p();
539  if (eGnow > eGspl) {
540  iGspl = i;
541  eGspl = eGnow;
542  }
543  }
544  }
545  if (iGspl == -1) return false;
546 
547  // Split this gluon into a collinear quark.antiquark pair.
548  int iG = systemPtr->iParton[iGspl];
549  int idNewQ = flavSelPtr->pickLightQ();
550  int iNewQ = event.append( idNewQ, 101, iG, 0, 0, 0, event[iG].col(), 0,
551  0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
552  int iNewQb = event.append( -idNewQ, 101, iG, 0, 0, 0, 0, event[iG].acol(),
553  0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
554  event[iG].statusNeg();
555  event[iG].daughters( iNewQ, iNewQb);
556 
557  // Order partons in new system. Note order of colour flow.
558  int iNext = iGspl + 1;
559  if (iNext == systemPtr->size()) iNext = 0;
560  if (event[ systemPtr->iParton[iNext]].acol() != event[iNewQ].col())
561  swap( iNewQ, iNewQb);
562  vector<int> iNewSys;
563  iNewSys.push_back( iNewQ);
564  for (int i = iGspl + 1; i < systemPtr->size(); ++i)
565  iNewSys.push_back( systemPtr->iParton[i]);
566  for (int i = 0; i < iGspl; ++i)
567  iNewSys.push_back( systemPtr->iParton[i]);
568  iNewSys.push_back( iNewQb);
569 
570  // Erase the old system and insert the new one instead.
571  colConfig.erase(iSys);
572  colConfig.insert( iNewSys, event);
573 
574  // Done.
575  return true;
576 
577 }
578 
579 //--------------------------------------------------------------------------
580 
581 // Split a single colour singlet that contains two sparticles.
582 // To fix : if nLeg >= 3 && mMin large handle as in nLeg == 1??
583 
584 bool RHadrons::splitSystem( ColConfig& colConfig, Event& event) {
585 
586  // First and second R-hadron mother. Number of legs in between.
587  int iFirst = -1;
588  int iSecond = -1;
589  for (int i = 0; i < int(systemPtr->size()); ++i) {
590  int iTmp = systemPtr->iParton[i];
591  if ( givesRHadron( event[iTmp].id()) ) {
592  if (iFirst == -1) iFirst = i;
593  else iSecond = i;
594  }
595  }
596  int nLeg = iSecond - iFirst;
597 
598  // New flavour pair for breaking the string, and its mass.
599  int idNewQ = flavSelPtr->pickLightQ();
600  double mNewQ = particleDataPtr->constituentMass( idNewQ);
601  vector<int> iNewSys1, iNewSys2;
602 
603  // If sparticles are neighbours then need new q-qbar pair in between.
604  if (nLeg == 1) {
605 
606  // The location of the two sparticles.
607  int i1Old = systemPtr->iParton[iFirst];
608  int i2Old = systemPtr->iParton[iSecond];
609 
610  // Take a fraction of momentum to give breakup quark pair.
611  double mHat = (event[i1Old].p() + event[i2Old].p()).mCalc();
612  double mMax = mHat - event[i1Old].m() - event[i2Old].m();
613  if (mMax < 2. * (mNewQ + MSAFETY)) return false;
614  double mEff = min( 2. * (mNewQ + mOffsetCloudRH), mMax - 2. * MSAFETY);
615  double frac = mEff / mHat;
616  Vec4 pEff = frac * (event[i1Old].p() + event[i2Old].p());
617 
618  // New kinematics by (1) go to same mHat with bigger masses, and
619  // (2) rescale back down to original masses and reduced mHat.
620  Vec4 p1New, p2New;
621  if ( !newKin( event[i1Old].p(), event[i2Old].p(),
622  event[i1Old].m() / (1. - frac), event[i2Old].m() / (1. - frac),
623  p1New, p2New) ) return false;
624  p1New *= 1. - frac;
625  p2New *= 1. - frac;
626 
627  // Fill in new partons after branching, with correct colour/flavour sign.
628  int i1New, i2New, i3New, i4New;
629  int newCol = event.nextColTag();
630  i1New = event.copy( i1Old, 101);
631  if (event[i2Old].acol() == event[i1Old].col()) {
632  i3New = event.append( -idNewQ, 101, i1Old, 0, 0, 0,
633  0, event[i2Old].acol(), 0.5 * pEff, 0.5 * mEff, 0.);
634  i2New = event.copy( i2Old, 101);
635  event[i2New].acol( newCol);
636  i4New = event.append( idNewQ, 101, i2Old, 0, 0, 0,
637  newCol, 0, 0.5 * pEff, 0.5 * mEff, 0.);
638  } else {
639  i3New = event.append( idNewQ, 101, i1Old, 0, 0, 0,
640  event[i2Old].col(), 0, 0.5 * pEff, 0.5 * mEff, 0.);
641  i2New = event.copy( i2Old, 101);
642  event[i2New].col( newCol);
643  i4New = event.append( -idNewQ, 101, i2Old, 0, 0, 0,
644  0, newCol, 0.5 * pEff, 0.5 * mEff, 0.);
645  }
646 
647  // Modify momenta and history. For iBotCopyId tracing asymmetric
648  // bookkeeping where one sparticle is radiator and the other recoiler.
649  event[i1New].p( p1New);
650  event[i2New].p( p2New);
651  event[i1Old].daughters( i1New, i3New);
652  event[i1New].mother2( 0);
653  event[i2Old].daughters( i2New, i4New);
654  event[i2New].mother2( 0);
655  iBefRHad[0] = i1New;
656  iBefRHad[1] = i2New;
657 
658  // Partons in the two new systems.
659  for (int i = 0; i < iFirst; ++i)
660  iNewSys1.push_back( systemPtr->iParton[i]);
661  iNewSys1.push_back( i1New);
662  iNewSys1.push_back( i3New);
663  iNewSys2.push_back( i4New);
664  iNewSys2.push_back( i2New);
665  for (int i = iSecond + 1; i < int(systemPtr->size()); ++i)
666  iNewSys2.push_back( systemPtr->iParton[i]);
667 
668  // If one gluon between sparticles then split it into a
669  // collinear q-qbar pair.
670  } else if (nLeg == 2) {
671 
672  // Gluon to be split and its two daughters.
673  int iOld = systemPtr->iParton[iFirst + 1];
674  int i1New = event.append( idNewQ, 101, iOld, 0, 0, 0,
675  event[iOld].col(), 0, 0.5 * event[iOld].p(),
676  0.5 * event[iOld].m(), 0.);
677  int i2New = event.append( -idNewQ, 101, iOld, 0, 0, 0,
678  0, event[iOld].acol(), 0.5 * event[iOld].p(),
679  0.5 * event[iOld].m(), 0.);
680  event[iOld].statusNeg();
681  event[iOld].daughters( i1New, i2New);
682 
683  // Partons in the two new systems.
684  if (event[systemPtr->iParton[iFirst]].col() == event[i2New].acol())
685  swap( i1New, i2New);
686  for (int i = 0; i <= iFirst; ++i)
687  iNewSys1.push_back( systemPtr->iParton[i]);
688  iNewSys1.push_back( i1New);
689  iNewSys2.push_back( i2New);
690  for (int i = iSecond; i < int(systemPtr->size()); ++i)
691  iNewSys2.push_back( systemPtr->iParton[i]);
692 
693  // If several gluons between sparticles then find lowest-mass gluon pair
694  // and replace it by a q-qbar pair.
695  } else {
696 
697  // Find lowest-mass gluon pair and adjust effective quark mass.
698  int iMin = 0;
699  int i1Old = 0;
700  int i2Old = 0;
701  double mMin = 1e20;
702  for (int i = iFirst + 1; i < iSecond - 1; ++i) {
703  int i1Tmp = systemPtr->iParton[i];
704  int i2Tmp = systemPtr->iParton[i + 1];
705  double mTmp = (event[i1Tmp].p() + event[i2Tmp].p()).mCalc();
706  if (mTmp < mMin) {
707  iMin = i;
708  i1Old = i1Tmp;
709  i2Old = i2Tmp;
710  mMin = mTmp;
711  }
712  }
713  double mEff = min( mNewQ + mOffsetCloudRH, 0.4 * mMin);
714 
715  // New kinematics by sharing mHat appropriately.
716  Vec4 p1New, p2New;
717  if ( !newKin( event[i1Old].p(), event[i2Old].p(),
718  mEff, mEff, p1New, p2New, false) ) return false;
719 
720  // Insert new partons and update history.
721  int i1New, i2New;
722  if (event[systemPtr->iParton[0]].acol() == 0) {
723  i1New = event.append( -idNewQ, 101, i1Old, 0, 0, 0,
724  0, event[i1Old].acol(), p1New, mEff, 0.);
725  i2New = event.append( idNewQ, 101, i2Old, 0, 0, 0,
726  event[i2Old].col(), 0, p2New, mEff, 0.);
727  } else {
728  i1New = event.append( idNewQ, 101, i1Old, 0, 0, 0,
729  event[i1Old].col(), 0, p1New, mEff, 0.);
730  i2New = event.append( -idNewQ, 101, i2Old, 0, 0, 0,
731  0, event[i2Old].acol(), p2New, mEff, 0.);
732  }
733  event[i1Old].statusNeg();
734  event[i2Old].statusNeg();
735  event[i1Old].daughters( i1New, 0);
736  event[i2Old].daughters( i2New, 0);
737 
738  // Partons in the two new systems.
739  for (int i = 0; i < iMin; ++i)
740  iNewSys1.push_back( systemPtr->iParton[i]);
741  iNewSys1.push_back( i1New);
742  iNewSys2.push_back( i2New);
743  for (int i = iMin + 2; i < int(systemPtr->size()); ++i)
744  iNewSys2.push_back( systemPtr->iParton[i]);
745  }
746 
747  // Erase the old system and insert the two new ones instead.
748  colConfig.erase(iSys);
749  colConfig.insert( iNewSys1, event);
750  colConfig.insert( iNewSys2, event);
751 
752  // Done.
753  return true;
754 
755 }
756 
757 //--------------------------------------------------------------------------
758 
759 // Produce a R-hadron from a squark and another string end.
760 
761 bool RHadrons::produceSquark( ColConfig& colConfig, Event& event) {
762 
763  // Initial values.
764  int nBody = 0;
765  int iRNow = 0;
766  int iNewQ = 0;
767  int iNewL = 0;
768 
769  // Check at which end of the string the squark is located.
770  int idAbsTop = event[ systemPtr->iParton[0] ].idAbs();
771  bool sqAtTop = (allowRSb && idAbsTop == idRSb)
772  || (allowRSt && idAbsTop == idRSt);
773 
774  // Copy down system consecutively from squark end.
775  int iBeg = event.size();
776  iCreRHad[iRHad] = iBeg;
777  if (sqAtTop) for (int i = 0; i < systemPtr->size(); ++i)
778  event.copy( systemPtr->iParton[i], 102);
779  else for (int i = systemPtr->size() - 1; i >= 0; --i)
780  event.copy( systemPtr->iParton[i], 102);
781  int iEnd = event.size() - 1;
782 
783  // Input flavours. From now on H = heavy and L = light string ends.
784  int idOldH = event[iBeg].id();
785  int idOldL = event[iEnd].id();
786 
787  // Pick new flavour to form R-hadron.
788  FlavContainer flavOld( idOldH%10);
789  int idNewQ = flavSelPtr->pick(flavOld).id;
790  int idRHad = toIdWithSquark( idOldH, idNewQ);
791  if (idRHad == 0) {
792  infoPtr->errorMsg("Error in RHadrons::produceSquark: "
793  "cannot form R-hadron code");
794  return false;
795  }
796 
797  // Target mass of R-hadron and z value of fragmentation function.
798  double mRHad = particleDataPtr->m0(idRHad) + event[iBeg].m()
799  - ( (abs(idOldH) == idRSb) ? m0Sb : m0St );
800  double z = zSelPtr->zFrag( idOldH, idNewQ, mRHad*mRHad);
801 
802  // Basic kinematics of string piece where break is to occur.
803  Vec4 pOldH = event[iBeg].p();
804  int iOldL = iBeg + 1;
805  Vec4 pOldL = event[iOldL].p();
806  double mOldL = event[iOldL].m();
807  double mNewH = mRHad / z;
808  double sSys = (pOldH + pOldL).m2Calc();
809  double sRem = (1. - z) * (sSys - mNewH*mNewH);
810  double sMin = pow2(mOldL + mCollapseRH);
811 
812  // If too low remaining mass in system then add one more parton to it.
813  while ( ( sRem < sMin || sSys < pow2(mNewH + mOldL + MSAFETY) )
814  && iOldL < iEnd ) {
815  ++iOldL;
816  pOldL += event[iOldL].p();
817  mOldL = event[iOldL].m();
818  sSys = (pOldH + pOldL).m2Calc();
819  sRem = (1. - z) * (sSys - mNewH*mNewH);
820  sMin = pow2(mOldL + mCollapseRH);
821  }
822 
823  // If enough mass then split off R-hadron and reduced system.
824  if ( sRem > sMin && sSys > pow2(mNewH + mOldL + MSAFETY) ) {
825  Vec4 pNewH, pNewL;
826  if ( !newKin( pOldH, pOldL, mNewH, mOldL, pNewH, pNewL) ) {
827  infoPtr->errorMsg("Error in RHadrons::produceSquark: "
828  "failed to construct kinematics with reduced system");
829  return false;
830  }
831 
832  // Insert R-hadron with new momentum.
833  iRNow = event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
834  z * pNewH, mRHad, 0.);
835 
836  // Reduced system with new string endpoint and modified recoiler.
837  idNewQ = -idNewQ;
838  bool hasCol = (idNewQ > 0 && idNewQ < 10) || idNewQ < -10;
839  int col = (hasCol) ? event[iOldL].acol() : 0;
840  int acol = (hasCol) ? 0 : event[iOldL].col();
841  iNewQ = event.append( idNewQ, 105, iBeg, iOldL, 0, 0, col, acol,
842  (1. - z) * pNewH, (1. - z) * mNewH, 0.);
843  iNewL = event.copy( iOldL, 105);
844  event[iNewL].mothers( iBeg, iOldL);
845  event[iNewL].p( pNewL);
846 
847  // Done with processing of split to R-hadron and reduced system.
848  nBody = 3;
849  }
850 
851  // Two-body final state: form light hadron from endpoint and new flavour.
852  if (nBody == 0) {
853  FlavContainer flav1( idOldL);
854  FlavContainer flav2( -idNewQ);
855  int iTry = 0;
856  int idNewL = flavSelPtr->combine( flav1, flav2);
857  while (++iTry < NTRYMAX && idNewL == 0)
858  idNewL = flavSelPtr->combine( flav1, flav2);
859  if (idNewL == 0) {
860  infoPtr->errorMsg("Error in RHadrons::produceSquark: "
861  "cannot form light hadron code");
862  return false;
863  }
864  double mNewL = particleDataPtr->mass( idNewL);
865 
866  // Check that energy enough for light hadron and R-hadron.
867  if ( sSys > pow2(mRHad + mNewL + MSAFETY) ) {
868  Vec4 pRHad, pNewL;
869  if ( !newKin( pOldH, pOldL, mRHad, mNewL, pRHad, pNewL) ) {
870  infoPtr->errorMsg("Error in RHadrons::produceSquark: "
871  "failed to construct kinematics for two-hadron decay");
872  return false;
873  }
874 
875  // Insert R-hadron and light hadron.
876  iRNow = event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
877  pRHad, mRHad, 0.);
878  event.append( idNewL, 105, iBeg, iOldL, 0, 0, 0, 0,
879  pNewL, mNewL, 0.);
880 
881  // Done for two-body case.
882  nBody = 2;
883  }
884  }
885 
886  // Final case: for very low mass collapse to one single R-hadron.
887  if (nBody == 0) {
888  idRHad = toIdWithSquark( idOldH, idOldL);
889  if (idRHad == 0) {
890  infoPtr->errorMsg("Error in RHadrons::produceSquark: "
891  "cannot form R-hadron code");
892  return false;
893  }
894 
895  // Insert R-hadron with new momentum.
896  iRNow = event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
897  systemPtr->pSum, systemPtr->mass, 0.);
898 
899  // Done with one-body case.
900  nBody = 1;
901  }
902 
903  // History bookkeeping: the R-hadron and processed partons.
904  iRHadron[iRHad] = iRNow;
905  int iLast = event.size() - 1;
906  for (int i = iBeg; i <= iOldL; ++i) {
907  event[i].statusNeg();
908  event[i].daughters( iRNow, iLast);
909  }
910 
911  // Remove old string system and, if needed, insert a new one.
912  colConfig.erase(iSys);
913  if (nBody == 3) {
914  vector<int> iNewSys;
915  iNewSys.push_back( iNewQ);
916  iNewSys.push_back( iNewL);
917  for ( int i = iOldL + 1; i <= iEnd; ++i) iNewSys.push_back( i);
918  colConfig.insert( iNewSys, event);
919  }
920 
921  // Copy lifetime and vertex from sparticle to R-hadron.
922  event[iRNow].tau( event[iBef].tau() );
923  if (event[iBef].hasVertex()) event[iRNow].vProd( event[iBef].vProd() );
924 
925  // Done with production of a R-hadron from a squark.
926  return true;
927 
928 }
929 
930 //--------------------------------------------------------------------------
931 
932 // Produce a R-hadron from a gluino and two string ends (or a gluon).
933 
934 bool RHadrons::produceGluino( ColConfig& colConfig, Event& event) {
935 
936  // Initial values.
937  int iGlui = 0;
938  int idSave = 0;
939  int idQLeap = 0;
940  bool isDiq1 = false;
941  vector<int> iSide1, iSide2, iSideTmp, iNewSys1, iNewSys2;
942  Vec4 pGlui, pSide1, pSide2;
943 
944  // Decide whether to produce a gluinoball or not.
945  bool isGBall = (rndmPtr->flat() < probGluinoballRH);
946 
947  // Extract one string system on either side of the gluino.
948  for (int i = 0; i < systemPtr->size(); ++i) {
949  int iTmp = systemPtr->iParton[i];
950  if (iGlui == 0 && event[ iTmp ].id() == idRGo) {
951  iGlui = iTmp;
952  pGlui = event[ iTmp ].p();
953  } else if (iGlui == 0) {
954  iSideTmp.push_back( iTmp);
955  pSide1 += event[ iTmp ].p();
956  } else {
957  iSide2.push_back( iTmp);
958  pSide2 += event[ iTmp ].p();
959  }
960  }
961 
962  // Order sides from gluino outwards and with lowest relative mass first.
963  for (int i = int(iSideTmp.size()) - 1; i >= 0; --i)
964  iSide1.push_back( iSideTmp[i]);
965  double m1H = (pGlui + pSide1).mCalc() - event[ iSide1.back() ].m();
966  double m2H = (pGlui + pSide2).mCalc() - event[ iSide2.back() ].m();
967  if (m2H < m1H) {
968  swap( iSide1, iSide2);
969  swap( pSide1, pSide2);
970  }
971 
972  // Begin loop over two sides of gluino, with lowest relative mass first.
973  for (int iSide = 1; iSide < 3; ++iSide) {
974 
975  // Begin processing of lower-mass string side.
976  int idRHad = 0;
977  double mRHad = 0.;
978  int nBody = 0;
979  int iRNow = 0;
980  int iNewQ = 0;
981  int iNewL = 0;
982  int statusRHad = 0;
983 
984  // Copy down system consecutively from gluino end.
985  int iBeg = event.size();
986  if (iSide == 1) {
987  iCreRHad[iRHad] = iBeg;
988  event.copy( iGlui, 102);
989  for (int i = 0; i < int(iSide1.size()); ++i)
990  event.copy( iSide1[i], 102);
991  } else {
992  event.copy( iGlui, 102);
993  for (int i = 0; i < int(iSide2.size()); ++i)
994  event.copy( iSide2[i], 102);
995  }
996  int iEnd = event.size() - 1;
997 
998  // Pick new flavour to help form R-hadron. Initial colour values.
999  int idOldL = event[iEnd].id();
1000  int idNewQ = 21;
1001  if (!isGBall) {
1002  do {
1003  FlavContainer flavOld( idOldL);
1004  idNewQ = -flavSelPtr->pick(flavOld).id;
1005  } while (iSide == 2 && isDiq1 && abs(idNewQ) > 10);
1006  if (iSide == 1) isDiq1 = (abs(idNewQ) > 10);
1007  }
1008  bool hasCol = (event[iEnd].col() > 0);
1009  int colR = 0;
1010  int acolR = 0;
1011 
1012  // Target intermediary mass or R-hadron mass.
1013  if (iSide == 1) {
1014  idSave = idNewQ;
1015  idRHad = (hasCol) ? 1009002 : -1009002;
1016  if (hasCol) colR = event[iBeg].col();
1017  else acolR = event[iBeg].acol();
1018  statusRHad = 103;
1019  double mNewQ = particleDataPtr->constituentMass( idNewQ);
1020  if (isGBall) mNewQ *= 0.5;
1021  mRHad = event[iBeg].m() + mOffsetCloudRH + mNewQ;
1022  } else {
1023  idRHad = toIdWithGluino( idSave, idNewQ);
1024  if (idRHad == 0) {
1025  infoPtr->errorMsg("Error in RHadrons::produceGluino: "
1026  "cannot form R-hadron code");
1027  return false;
1028  }
1029  statusRHad = 104;
1030  mRHad = particleDataPtr->m0(idRHad) + event[iBeg].m() - m0Go;
1031  }
1032 
1033  // z value of fragmentation function.
1034  double z = zSelPtr->zFrag( idRGo, idNewQ, mRHad*mRHad);
1035 
1036  // Basic kinematics of string piece where break is to occur.
1037  Vec4 pOldH = event[iBeg].p();
1038  int iOldL = iBeg + 1;
1039  Vec4 pOldL = event[iOldL].p();
1040  double mOldL = event[iOldL].m();
1041  double mNewH = mRHad / z;
1042  double sSys = (pOldH + pOldL).m2Calc();
1043  double sRem = (1. - z) * (sSys - mNewH*mNewH);
1044  double sMin = pow2(mOldL + mCollapseRH);
1045 
1046  // If too low remaining mass in system then add one more parton to it.
1047  while ( ( sRem < sMin || sSys < pow2(mNewH + mOldL + MSAFETY) )
1048  && iOldL < iEnd ) {
1049  ++iOldL;
1050  pOldL += event[iOldL].p();
1051  mOldL = event[iOldL].m();
1052  sSys = (pOldH + pOldL).m2Calc();
1053  sRem = (1. - z) * (sSys - mNewH*mNewH);
1054  sMin = pow2(mOldL + mCollapseRH);
1055  }
1056 
1057  // If enough mass then split off R-hadron and reduced system.
1058  if ( sRem > sMin && sSys > pow2(mNewH + mOldL + MSAFETY) ) {
1059  Vec4 pNewH, pNewL;
1060  if ( !newKin( pOldH, pOldL, mNewH, mOldL, pNewH, pNewL) ) {
1061  infoPtr->errorMsg("Error in RHadrons::produceGluino: "
1062  "failed to construct kinematics with reduced system");
1063  return false;
1064  }
1065 
1066  // Insert intermediary or R-hadron with new momentum, less colour.
1067  iRNow = event.append( idRHad, statusRHad, iBeg, iOldL,
1068  0, 0, colR, acolR, z * pNewH, mRHad, 0.);
1069 
1070  // Reduced system with new string endpoint and modified recoiler.
1071  if (!isGBall) idNewQ = -idNewQ;
1072  int colN = (hasCol) ? 0 : event[iOldL].acol();
1073  int acolN = (hasCol) ? event[iOldL].col() : 0;
1074  iNewQ = event.append( idNewQ, 105, iBeg, iOldL, 0, 0,
1075  colN, acolN, (1. - z) * pNewH, (1. - z) * mNewH, 0.);
1076  iNewL = event.copy( iOldL, 105);
1077  event[iNewL].mothers( iBeg, iOldL);
1078  event[iNewL].p( pNewL);
1079 
1080  // Collect partons of new string system.
1081  if (iSide == 1) {
1082  iNewSys1.push_back( iNewQ);
1083  iNewSys1.push_back( iNewL);
1084  for ( int i = iOldL + 1; i <= iEnd; ++i) iNewSys1.push_back( i);
1085  } else {
1086  iNewSys2.push_back( iNewQ);
1087  iNewSys2.push_back( iNewL);
1088  for ( int i = iOldL + 1; i <= iEnd; ++i) iNewSys2.push_back( i);
1089  }
1090 
1091  // Done with processing of split to R-hadron and reduced system.
1092  nBody = 3;
1093  }
1094 
1095  // For side-1 low-mass glueball system reabsorb full momentum.
1096  if (nBody == 0 && isGBall && iSide == 1) {
1097  idQLeap = event[iOldL].id();
1098  Vec4 pNewH = event[iBeg].p() + pOldL;
1099 
1100  // Insert intermediary R-hadron with new momentum, less colour.
1101  iRNow = event.append( idRHad, statusRHad, iBeg, iEnd,
1102  0, 0, colR, acolR, pNewH, pNewH.mCalc(), 0.);
1103  nBody = 1;
1104  }
1105 
1106  // Two-body final state: form light hadron from endpoint and new flavour.
1107  // Also for side 2 if gluinoball formation gives quark from side 1.
1108  if (nBody == 0 && (!isGBall || (iSide == 2 && idQLeap != 0) )) {
1109  if (isGBall) idNewQ = -idQLeap;
1110  FlavContainer flav1( idOldL);
1111  FlavContainer flav2( -idNewQ);
1112  int iTry = 0;
1113  int idNewL = flavSelPtr->combine( flav1, flav2);
1114  while (++iTry < NTRYMAX && idNewL == 0)
1115  idNewL = flavSelPtr->combine( flav1, flav2);
1116  if (idNewL == 0) {
1117  infoPtr->errorMsg("Error in RHadrons::produceGluino: "
1118  "cannot form light hadron code");
1119  return false;
1120  }
1121  double mNewL = particleDataPtr->mass( idNewL);
1122 
1123  // Check that energy enough for light hadron and R-hadron.
1124  if ( sSys > pow2(mRHad + mNewL + MSAFETY) ) {
1125  Vec4 pRHad, pNewL;
1126  if ( !newKin( pOldH, pOldL, mRHad, mNewL, pRHad, pNewL) ) {
1127  infoPtr->errorMsg("Error in RHadrons::produceGluino: "
1128  "failed to construct kinematics for two-hadron decay");
1129  return false;
1130  }
1131 
1132  // Insert intermediary or R-hadron and light hadron.
1133  iRNow = event.append( idRHad, statusRHad, iBeg, iOldL, 0, 0,
1134  colR, acolR, pRHad, mRHad, 0.);
1135  event.append( idNewL, 105, iBeg, iOldL, 0, 0, 0, 0,
1136  pNewL, mNewL, 0.);
1137 
1138  // Done for two-body case.
1139  nBody = 2;
1140  // For this case gluinoball should be handled as normal flavour.
1141  isGBall = false;
1142  }
1143  }
1144 
1145  // Final case: for very low mass collapse to one single R-hadron.
1146  if (nBody == 0 && (!isGBall || (iSide == 2 && idQLeap != 0) )) {
1147  if (isGBall) idSave = idQLeap;
1148  if (iSide == 1) idSave = idOldL;
1149  else idRHad = toIdWithGluino( idSave, idOldL);
1150  if (idRHad == 0) {
1151  infoPtr->errorMsg("Error in RHadrons::produceGluino: "
1152  "cannot form R-hadron code");
1153  return false;
1154  }
1155 
1156  // Insert R-hadron with new momentum.
1157  iRNow = event.append( idRHad, statusRHad, iBeg, iOldL, 0, 0,
1158  colR, acolR, pOldH + pOldL, (pOldH + pOldL).mCalc(), 0.);
1159 
1160  // Done with one-body case.
1161  // Even if hoped-for, it was not possible to create a gluinoball.
1162  isGBall = false;
1163  }
1164 
1165  // History bookkeeping: the processed partons.
1166  int iLast = event.size() - 1;
1167  for (int i = iBeg; i <= iOldL; ++i) {
1168  event[i].statusNeg();
1169  event[i].daughters( iRNow, iLast);
1170  }
1171 
1172  // End of loop over two sides of the gluino.
1173  iGlui = iRNow;
1174  }
1175 
1176  // History bookkeeping: insert R-hadron; delete old string system.
1177  if (iGlui == 0) {
1178  infoPtr->errorMsg("Error in RHadrons::produceGluino: "
1179  "could not handle gluinoball kinematics");
1180  return false;
1181  }
1182  iRHadron[iRHad] = iGlui;
1183  colConfig.erase(iSys);
1184 
1185  // Non-gluinoball: insert (up to) two new string systems.
1186  if (!isGBall) {
1187  if (iNewSys1.size() > 0) colConfig.insert( iNewSys1, event);
1188  if (iNewSys2.size() > 0) colConfig.insert( iNewSys2, event);
1189 
1190  // Gluinoball with enough energy in first string:
1191  // join two temporary gluons into one.
1192  } else if (idQLeap == 0) {
1193  int iG1 = iNewSys1[0];
1194  int iG2 = iNewSys2[0];
1195  int colG = event[iG1].col() + event[iG2].col();
1196  int acolG = event[iG1].acol() + event[iG2].acol();
1197  Vec4 pG = event[iG1].p() + event[iG2].p();
1198  int iG12 = event.append( 21, 105, iG1, iG2, 0, 0, colG, acolG,
1199  pG, pG.mCalc(), 0.);
1200 
1201  // Temporary gluons no longer needed, but new colour to avoid warnings.
1202  event[iG1].id( 21);
1203  event[iG2].id( 21);
1204  event[iG1].statusNeg();
1205  event[iG2].statusNeg();
1206  int colBridge = event.nextColTag();
1207  if (event[iG1].col() == 0) {
1208  event[iG1].col( colBridge);
1209  event[iG2].acol( colBridge);
1210  } else {
1211  event[iG1].acol( colBridge);
1212  event[iG2].col( colBridge);
1213  }
1214 
1215  // Form the remnant system from which the R-hadron has been split off.
1216  vector<int> iNewSys12;
1217  for (int i = int(iNewSys1.size()) - 1; i > 0; --i)
1218  iNewSys12.push_back( iNewSys1[i]);
1219  iNewSys12.push_back( iG12);
1220  for (int i = 1; i < int(iNewSys2.size()); ++i)
1221  iNewSys12.push_back( iNewSys2[i]);
1222  colConfig.insert( iNewSys12, event);
1223 
1224  // Gluinoball where side 1 was fully eaten, and its flavour content
1225  // should leap over to the other side, to a gluon there.
1226  } else {
1227  int iG2 = iNewSys2[0];
1228  event[iG2].id( idQLeap);
1229  colConfig.insert( iNewSys2, event);
1230  }
1231 
1232  // Copy lifetime and vertex from sparticle to R-hadron.
1233  event[iGlui].tau( event[iBef].tau() );
1234  if (event[iBef].hasVertex()) event[iGlui].vProd( event[iBef].vProd() );
1235 
1236  // Done with production of a R-hadron from a gluino.
1237  return true;
1238 
1239 }
1240 
1241 //--------------------------------------------------------------------------
1242 
1243 // Form a R-hadron code from a squark and a (di)quark code.
1244 // First argument is the (anti)squark, second the (anti)(di)quark.
1245 
1246 int RHadrons::toIdWithSquark( int id1, int id2) {
1247 
1248  // Check that physical combination; return 0 if not.
1249  int id1Abs = abs(id1);
1250  int id2Abs = abs(id2);
1251  if (id2Abs < 10 && id1 > 0 && id2 > 0) return 0;
1252  if (id2Abs < 10 && id1 < 0 && id2 < 0) return 0;
1253  if (id2Abs > 10 && id1 > 0 && id2 < 0) return 0;
1254  if (id2Abs > 10 && id1 < 0 && id2 > 0) return 0;
1255 
1256  // Form R-hadron code. Flip sign for antisquark.
1257  bool isSt = (id1Abs == idRSt);
1258  int idRHad = 1000000;
1259  if (id2Abs < 10) idRHad += ((isSt) ? 600 : 500) + 10 * id2Abs + 2;
1260  else idRHad += ((isSt) ? 6000 : 5000) + 10 * (id2Abs/100) + id2Abs%10;
1261  if (id1 < 0) idRHad = -idRHad;
1262 
1263  // Done.
1264  return idRHad;
1265 
1266 }
1267 
1268 //--------------------------------------------------------------------------
1269 
1270 // Resolve a R-hadron code into a squark and a (di)quark.
1271 // Return a pair, where first is the squark and the second the (di)quark.
1272 
1273 pair<int,int> RHadrons::fromIdWithSquark( int idRHad) {
1274 
1275  // Find squark flavour content.
1276  int idLight = (abs(idRHad) - 1000000) / 10;
1277  int idSq = (idLight < 100) ? idLight/10 : idLight/100;
1278  int id1 = (idSq == 6) ? idRSt : idRSb;
1279  if (idRHad < 0) id1 = -id1;
1280 
1281  // Find light (di)quark flavour content.
1282  int id2 = (idLight < 100) ? idLight%10 : idLight%100;
1283  if (id2 > 10) id2 = 100 * id2 + abs(idRHad)%10;
1284  if ((id2 < 10 && idRHad > 0) || (id2 > 10 && idRHad < 0)) id2 = -id2;
1285 
1286  // Done.
1287  return make_pair( id1, id2);
1288 
1289 }
1290 
1291 //--------------------------------------------------------------------------
1292 
1293 // Form a R-hadron code from two quark/diquark endpoints and a gluino.
1294 
1295 int RHadrons::toIdWithGluino( int id1, int id2) {
1296 
1297  // Check that physical combination; return 0 if not. Handle gluinoball.
1298  int id1Abs = abs(id1);
1299  int id2Abs = abs(id2);
1300  if (id1Abs == 21 && id2Abs == 21) return 1000993;
1301  int idMax = max( id1Abs, id2Abs);
1302  int idMin = min( id1Abs, id2Abs);
1303  if (idMin > 10) return 0;
1304  if (idMax > 10 && id1 > 0 && id2 < 0) return 0;
1305  if (idMax > 10 && id1 < 0 && id2 > 0) return 0;
1306  if (idMax < 10 && id1 > 0 && id2 > 0) return 0;
1307  if (idMax < 10 && id1 < 0 && id2 < 0) return 0;
1308 
1309  // Form R-meson code. Flip sign for antiparticle.
1310  int idRHad = 0;
1311  if (idMax < 10) {
1312  idRHad = 1009003 + 100 * idMax + 10 * idMin;
1313  if (idMin != idMax && idMax%2 == 1) {
1314  if (id1Abs == idMax && id1 > 0) idRHad = -idRHad;
1315  if (id2Abs == idMax && id2 > 0) idRHad = -idRHad;
1316  }
1317  if (idMin != idMax && idMax%2 == 0) {
1318  if (id1Abs == idMax && id1 < 0) idRHad = -idRHad;
1319  if (id2Abs == idMax && id2 < 0) idRHad = -idRHad;
1320  }
1321 
1322  // Form R-baryon code. Flip sign for antiparticle.
1323  } else {
1324  int idA = idMax/1000;
1325  int idB = (idMax/100)%10;
1326  int idC = idMin;
1327  if (idC > idB) swap( idB, idC);
1328  if (idB > idA) swap( idA, idB);
1329  if (idC > idB) swap( idB, idC);
1330  idRHad = 1090004 + 1000 * idA + 100 * idB + 10 * idC;
1331  if (id1 < 0) idRHad = -idRHad;
1332  }
1333 
1334  // Done.
1335  return idRHad;
1336 
1337 }
1338 
1339 //--------------------------------------------------------------------------
1340 
1341 // Resolve a R-hadron code into two quark/diquark endpoints (and a gluino).
1342 // Return a pair, where first carries colour and second anticolour.
1343 
1344 pair<int,int> RHadrons::fromIdWithGluino( int idRHad) {
1345 
1346  // Find light flavour content of R-hadron.
1347  int idLight = (abs(idRHad) - 1000000) / 10;
1348  int id1, id2, idTmp, idA, idB, idC;
1349 
1350  // Gluinoballs: split g into d dbar or u ubar.
1351  if (idLight < 100) {
1352  id1 = (rndmPtr->flat() < 0.5) ? 1 : 2;
1353  id2 = -id1;
1354 
1355  // Gluino-meson: split into q + qbar.
1356  } else if (idLight < 1000) {
1357  id1 = (idLight / 10) % 10;
1358  id2 = -(idLight % 10);
1359  // Flip signs when first quark of down-type.
1360  if (id1%2 == 1) {
1361  idTmp = id1;
1362  id1 = -id2;
1363  id2 = -idTmp;
1364  }
1365 
1366  // Gluino-baryon: split to q + qq (diquark).
1367  // Pick diquark at random, except if c or b involved.
1368  } else {
1369  idA = (idLight / 100) % 10;
1370  idB = (idLight / 10) % 10;
1371  idC = idLight % 10;
1372  double rndmQ = 3. * rndmPtr->flat();
1373  if (idA > 3) rndmQ = 0.5;
1374  if (rndmQ < 1.) {
1375  id1 = idA;
1376  id2 = 1000 * idB + 100 * idC + 3;
1377  if (idB != idC && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2;
1378  } else if (rndmQ < 2.) {
1379  id1 = idB;
1380  id2 = 1000 * idA + 100 * idC + 3;
1381  if (idA != idC && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2;
1382  } else {
1383  id1 = idC;
1384  id2 = 1000 * idA + 100 * idB +3;
1385  if (idA != idB && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2;
1386  }
1387  }
1388 
1389  // Flip signs for anti-R-hadron.
1390  if (idRHad < 0) {
1391  idTmp = id1;
1392  id1 = -id2;
1393  id2 = -idTmp;
1394  }
1395 
1396  // Done.
1397  return make_pair( id1, id2);
1398 
1399 }
1400 
1401 //--------------------------------------------------------------------------
1402 
1403 // Construct modified four-vectors to match modified masses:
1404 // minimal reshuffling of momentum along common axis.
1405 // Note that last two arguments are used to provide output!
1406 
1407 bool RHadrons::newKin( Vec4 pOld1, Vec4 pOld2, double mNew1, double mNew2,
1408  Vec4& pNew1, Vec4& pNew2, bool checkMargin) {
1409 
1410  // Squared masses in initial and final kinematics.
1411  double sSum = (pOld1 + pOld2).m2Calc();
1412  double sOld1 = pOld1.m2Calc();
1413  double sOld2 = pOld2.m2Calc();
1414  double sNew1 = mNew1 * mNew1;
1415  double sNew2 = mNew2 * mNew2;
1416 
1417  // Check that kinematically possible.
1418  if (checkMargin && pow2(mNew1 + mNew2 + MSAFETY) > sSum) return false;
1419 
1420  // Transfer coefficients to give four-vectors with the new masses.
1421  double lamOld = sqrt( pow2(sSum - sOld1 - sOld2) - 4. * sOld1 * sOld2 );
1422  double lamNew = sqrt( pow2(sSum - sNew1 - sNew2) - 4. * sNew1 * sNew2 );
1423  double move1 = (lamNew * (sSum - sOld1 + sOld2)
1424  - lamOld * (sSum - sNew1 + sNew2)) / (2. * sSum * lamOld);
1425  double move2 = (lamNew * (sSum + sOld1 - sOld2)
1426  - lamOld * (sSum + sNew1 - sNew2)) / (2. * sSum * lamOld);
1427 
1428  // Construct final vectors. Done.
1429  pNew1 = (1. + move1) * pOld1 - move2 * pOld2;
1430  pNew2 = (1. + move2) * pOld2 - move1 * pOld1;
1431  return true;
1432 
1433 }
1434 
1435 //==========================================================================
1436 
1437 } // end namespace Pythia8
Definition: AgUStep.h:26