8 #include "Pythia8/RHadrons.h"
21 const int RHadrons::IDRHADSB[14] = { 1000512, 1000522, 1000532,
22 1000542, 1000552, 1005113, 1005211, 1005213, 1005223, 1005311,
23 1005313, 1005321, 1005323, 1005333 };
25 const int RHadrons::IDRHADST[14] = { 1000612, 1000622, 1000632,
26 1000642, 1000652, 1006113, 1006211, 1006213, 1006223, 1006311,
27 1006313, 1006321, 1006323, 1006333 };
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 };
38 const int RHadrons::NTRYMAX = 10;
41 const double RHadrons::MSAFETY = 0.1;
44 const double RHadrons::EGBORROWMAX = 4.;
50 bool RHadrons::init() {
53 allowRH = flag(
"RHadrons:allow");
54 maxWidthRH = parm(
"RHadrons:maxWidth");
55 idRSb = mode(
"RHadrons:idSbottom");
56 idRSt = mode(
"RHadrons:idStop");
57 idRGo = mode(
"RHadrons:idGluino");
58 setMassesRH = flag(
"RHadrons:setMasses");
59 probGluinoballRH = parm(
"RHadrons:probGluinoball");
60 mOffsetCloudRH = parm(
"RHadrons:mOffsetCloud");
61 mCollapseRH = parm(
"RHadrons:mCollapse");
62 diquarkSpin1RH = parm(
"RHadrons:diquarkSpin1");
65 allowRSb = allowRH && idRSb > 0
66 && (particleDataPtr->mWidth(idRSb) < maxWidthRH);
67 allowRSt = allowRH && idRSt > 0
68 && (particleDataPtr->mWidth(idRSt) < maxWidthRH);
69 allowRGo = allowRH && idRGo > 0
70 && (particleDataPtr->mWidth(idRGo) < maxWidthRH);
71 allowSomeR = allowRSb || allowRSt || allowRGo;
75 m0Sb = particleDataPtr->m0(idRSb);
77 for (
int i = 0; i < 14; ++i) {
78 int idR = IDRHADSB[i];
79 double m0RHad = m0Sb + mOffsetCloudRH;
80 m0RHad += particleDataPtr->constituentMass( (idR%100)/10);
82 m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
83 particleDataPtr->m0( idR, m0RHad);
88 double mWidthRHad = particleDataPtr->mWidth(idRSb);
89 double tau0RHad = particleDataPtr->tau0( idRSb);
90 for (
int i = 0; i < 14; ++i) {
91 particleDataPtr->mWidth( IDRHADSB[i], mWidthRHad);
92 particleDataPtr->tau0( IDRHADSB[i], tau0RHad);
98 m0St = particleDataPtr->m0(idRSt);
100 for (
int i = 0; i < 14; ++i) {
101 int idR = IDRHADST[i];
102 double m0RHad = m0St + mOffsetCloudRH;
103 m0RHad += particleDataPtr->constituentMass( (idR%100)/10);
105 m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
106 particleDataPtr->m0( idR, m0RHad);
111 double mWidthRHad = particleDataPtr->mWidth(idRSt);
112 double tau0RHad = particleDataPtr->tau0( idRSt);
113 for (
int i = 0; i < 14; ++i) {
114 particleDataPtr->mWidth( IDRHADST[i], mWidthRHad);
115 particleDataPtr->tau0( IDRHADST[i], tau0RHad);
121 m0Go = particleDataPtr->m0(idRGo);
123 particleDataPtr->m0( IDRHADGO[0], m0Go + 2. * mOffsetCloudRH
124 + particleDataPtr->constituentMass(21) );
125 for (
int i = 1; i < 38; ++i) {
126 int idR = IDRHADGO[i];
127 double m0RHad = m0Go + 2. * mOffsetCloudRH;
128 m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
129 m0RHad += particleDataPtr->constituentMass( (idR%100)/10);
131 m0RHad += particleDataPtr->constituentMass( (idR%10000)/1000 );
132 particleDataPtr->m0( idR, m0RHad);
137 double mWidthRHad = particleDataPtr->mWidth(idRGo);
138 double tau0RHad = particleDataPtr->tau0( idRGo);
139 for (
int i = 0; i < 38; ++i) {
140 particleDataPtr->mWidth( IDRHADGO[i], mWidthRHad);
141 particleDataPtr->tau0( IDRHADGO[i], tau0RHad);
153 bool RHadrons::givesRHadron(
int id) {
154 if (allowRSb && abs(
id) == idRSb)
return true;
155 if (allowRSt && abs(
id) == idRSt)
return true;
156 if (allowRGo &&
id == idRGo)
return true;
165 bool RHadrons::produce( ColConfig& colConfig,
Event& event) {
168 if (!allowSomeR)
return true;
179 for (
int i = 0; i <
event.size(); ++i)
180 if (event[i].isFinal() && givesRHadron(event[i].
id())) {
181 iBefRHad.push_back(i);
182 iCreRHad.push_back(i);
183 iRHadron.push_back(0);
184 iAftRHad.push_back(0);
185 isTriplet.push_back(
true);
187 nRHad = iRHadron.size();
190 if (nRHad == 0)
return true;
194 infoPtr->errorMsg(
"Error in RHadrons::produce: "
195 "cannot handle more than two R-hadrons");
198 if (nRHad > 1 && rndmPtr->flat() > 0.5) swap( iBefRHad[0], iBefRHad[1]);
202 iSys = colConfig.findSinglet( iBef);
203 systemPtr = &colConfig[iSys];
204 if (systemPtr->hasJunction && !splitOffJunction( colConfig, event)) {
205 infoPtr->errorMsg(
"Error in RHadrons::produce: "
206 "cannot handle system with junction");
211 iSys = colConfig.findSinglet( iBef);
212 systemPtr = &colConfig[iSys];
213 if (systemPtr->hasJunction && !splitOffJunction( colConfig, event)) {
214 infoPtr->errorMsg(
"Error in RHadrons::produce: "
215 "cannot handle system with junction");
222 iSys = colConfig.findSinglet( iBef);
223 systemPtr = &colConfig[iSys];
224 if (systemPtr->isClosed && !openClosedLoop( colConfig, event)) {
225 infoPtr->errorMsg(
"Error in RHadrons::produce: "
226 "cannot open up closed gluon/gluino loop");
231 iSys = colConfig.findSinglet( iBef);
232 systemPtr = &colConfig[iSys];
233 if (systemPtr->isClosed && !openClosedLoop( colConfig, event)) {
234 infoPtr->errorMsg(
"Error in RHadrons::produce: "
235 "cannot open up closed gluon/gluino loop");
242 int iSys1 = colConfig.findSinglet( iBefRHad[0]);
243 int iSys2 = colConfig.findSinglet( iBefRHad[1]);
244 if (iSys2 == iSys1) {
246 systemPtr = &colConfig[iSys];
247 if ( !splitSystem( colConfig, event) ) {
248 infoPtr->errorMsg(
"Error in RHadrons::produce: "
249 "failed to handle two sparticles in same system");
256 for (iRHad = 0; iRHad < nRHad; ++iRHad) {
257 iBef = iBefRHad[iRHad];
258 iSys = colConfig.findSinglet( iBef);
260 infoPtr->errorMsg(
"Error in RHadrons::produce: "
261 "sparticle not in any colour singlet");
264 systemPtr = &colConfig[iSys];
267 if (systemPtr->hasJunction) {
268 infoPtr->errorMsg(
"Error in RHadrons::produce: "
269 "cannot handle system with junction");
272 if (systemPtr->isClosed) {
273 infoPtr->errorMsg(
"Error in RHadrons::produce: "
274 "cannot handle closed colour loop");
279 if (event[iBef].
id() == idRGo) isTriplet[iRHad] =
false;
280 bool formed = (isTriplet[iRHad]) ? produceSquark( colConfig, event)
281 : produceGluino( colConfig, event);
282 if (!formed)
return false;
295 bool RHadrons::decay(
Event& event) {
298 for (iRHad = 0; iRHad < nRHad; ++iRHad) {
299 int iRNow = iRHadron[iRHad];
300 int iRBef = iBefRHad[iRHad];
301 int idRHad =
event[iRNow].id();
302 double mRHad =
event[iRNow].m();
303 double mRBef =
event[iRBef].m();
308 pair<int,int> idPair = (isTriplet[iRHad])
309 ? fromIdWithSquark( idRHad) : fromIdWithGluino( idRHad);
310 int id1 = idPair.first;
311 int id2 = idPair.second;
315 double fracR = mRBef / mRHad;
317 infoPtr->errorMsg(
"Error in RHadrons::decay: "
318 "too low R-hadron mass for decay");
323 if (isTriplet[iRHad]) {
324 int colNew =
event.nextColTag();
325 int col = (
event[iRBef].col() != 0) ? colNew : 0;
326 int acol = (col == 0) ? colNew : 0;
329 iR0 =
event.append( id1, 106, iRNow, 0, 0, 0, col, acol,
330 fracR * event[iRNow].p(), fracR * mRHad, 0.);
331 iR2 =
event.append( id2, 106, iRNow, 0, 0, 0, acol, col,
332 (1. - fracR) * event[iRNow].p(), (1. - fracR) * mRHad, 0.);
336 double m1Eff = particleDataPtr->constituentMass(id1) + mOffsetCloudRH;
337 double m2Eff = particleDataPtr->constituentMass(id2) + mOffsetCloudRH;
338 double frac1 = (1. - fracR) * m1Eff / ( m1Eff + m2Eff);
339 double frac2 = (1. - fracR) * m2Eff / ( m1Eff + m2Eff);
342 int col1 =
event.nextColTag();
343 int col2 =
event.nextColTag();
346 iR0 =
event.append( idRGo, 106, iRNow, 0, 0, 0, col2, col1,
347 fracR * event[iRNow].p(), fracR * mRHad, 0.);
348 event.append( id1, 106, iRNow, 0, 0, 0, col1, 0,
349 frac1 * event[iRNow].p(), frac1 * mRHad, 0.);
350 iR2 =
event.append( id2, 106, iRNow, 0, 0, 0, 0, col2,
351 frac2 * event[iRNow].p(), frac2 * mRHad, 0.);
355 event[iRNow].statusNeg();
356 event[iRNow].daughters( iR0, iR2);
357 iAftRHad[iRHad] = iR0;
360 Vec4 vDec =
event[iRNow].vProd() +
event[iRNow].tau()
361 *
event[iR0].p() /
event[iR0].m();
362 for (
int iRd = iR0; iRd <= iR2; ++iRd) event[iRd].vProd( vDec);
377 bool RHadrons::splitOffJunction( ColConfig& colConfig,
Event& event) {
380 vector<int> leg1, leg2, leg3;
385 for (
int i = 0; i < systemPtr->size(); ++i) {
387 int iP = systemPtr->iParton[i];
391 }
else if (iLeg == 1) leg1.push_back( iP);
392 else if (iLeg == 2) leg2.push_back( iP);
393 else if (iLeg == 3) leg3.push_back( iP);
399 if (iLegSP == 0)
return false;
403 if (iLegSP == 2) swap(leg2, leg1);
404 else if (iLegSP == 3) swap(leg3, leg1);
405 for (
int i = 0; i < iPosSP; ++i)
406 if (event[leg1[i]].
id() != 21)
return true;
410 Vec4 pSP =
event[iBef].p();
412 for (
int i = 0; i < int(leg2.size()); ++i) pRec += event[leg2[i]].p();
413 for (
int i = 0; i < int(leg3.size()); ++i) pRec += event[leg3[i]].p();
414 double mSys = (pSP + pRec).mCalc();
415 double mSP = pSP.mCalc();
416 double mRec = pRec.mCalc();
417 double eKin = mSys - mSP - mRec;
420 double mNewG = EGBORROWMAX * eKin / (EGBORROWMAX + eKin) ;
421 Vec4 pNewG = (mNewG / mSys) * (pSP + pRec);
422 int newCol =
event.nextColTag();
423 bool isCol = (
event[leg1.back()].col() > 0);
424 int colNG = (isCol)? newCol : event[iBef].acol();
425 int acolNG = (isCol) ? event[iBef].col() : newCol;
426 int iNewG =
event.append( 21, 101, iBef, 0, 0, 0, colNG, acolNG,
428 leg1.insert( leg1.begin(), iNewG);
432 double mNewSys = mSys - mNewG;
433 double pAbsOld = 0.5 * sqrtpos( pow2(mSys*mSys - mSP*mSP - mRec*mRec)
434 - pow2(2. * mSP * mRec) ) / mSys;
435 double pAbsNew = 0.5 * sqrtpos( pow2(mNewSys*mNewSys - mSP*mSP - mRec*mRec)
436 - pow2(2. * mSP * mRec) ) / mNewSys;
437 RotBstMatrix MtoCM, MfromCM, MSP, MRec;
438 MtoCM.toCMframe( pSP, pRec);
442 MSP.bst( 0., 0., -pAbsOld / sqrt(mSP * mSP + pAbsOld * pAbsOld));
443 MSP.bst( 0., 0., pAbsNew / sqrt(mSP * mSP + pAbsNew * pAbsNew));
444 MSP.rotbst( MfromCM);
446 MRec.bst( 0., 0., pAbsOld / sqrt(mRec * mRec + pAbsOld * pAbsOld));
447 MRec.bst( 0., 0., -pAbsNew / sqrt(mRec * mRec + pAbsNew * pAbsNew));
448 MRec.rotbst( MfromCM);
451 int iNewSP =
event.copy( iBef, 101);
452 event[iNewSP].mother2(0);
453 event[iBef].daughter1(iNewG);
454 event[iNewSP].rotbst( MSP);
455 leg1[iPosSP] = iNewSP;
456 if (iBefRHad[0] == iBef) iBefRHad[0] = iNewSP;
457 else if (nRHad > 1 && iBefRHad[1] == iBef) iBefRHad[1] = iNewSP;
459 for (
int i = 0; i < int(leg2.size()); ++i) {
460 int iCopy =
event.copy( leg2[i], 101);
461 event[iCopy].rotbst( MRec);
462 if (iBefRHad[0] == leg2[i]) iBefRHad[0] = iCopy;
463 else if (nRHad > 1 && iBefRHad[1] == leg2[i]) iBefRHad[1] = iCopy;
466 for (
int i = 0; i < int(leg3.size()); ++i) {
467 int iCopy =
event.copy( leg3[i], 101);
468 event[iCopy].rotbst( MRec);
469 if (iBefRHad[0] == leg3[i]) iBefRHad[0] = iCopy;
470 else if (nRHad > 1 && iBefRHad[1] == leg3[i]) iBefRHad[1] = iCopy;
479 double eGspl =
event[leg1[0]].p() *
event[iBef].p();
480 for (
int i = 1; i < iPosSP; ++i) {
481 double eGnow =
event[leg1[i]].p() *
event[iBef].p();
487 int iG = leg1[iGspl];
490 int idNewQ = flavSelPtr->pickLightQ();
491 int iNewQ =
event.append( idNewQ, 101, iG, 0, 0, 0, event[iG].col(), 0,
492 0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
493 int iNewQb =
event.append( -idNewQ, 101, iG, 0, 0, 0, 0, event[iG].acol(),
494 0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
495 event[iG].statusNeg();
496 event[iG].daughters( iNewQ, iNewQb);
497 if (event[leg1.back()].col() == 0) swap( iNewQ, iNewQb);
500 vector<int> iNewSys1, iNewSys2;
501 iNewSys1.push_back( iNewQb);
502 for (
int i = iGspl + 1; i < int(leg1.size()); ++i)
503 iNewSys1.push_back( leg1[i]);
504 iNewSys2.push_back( -10);
505 for (
int i = 0; i < iGspl; ++i) iNewSys2.push_back( leg1[i]);
506 iNewSys2.push_back( iNewQ);
507 iNewSys2.push_back( -11);
508 for (
int i = 0; i < int(leg2.size()); ++i) iNewSys2.push_back( leg2[i]);
509 iNewSys2.push_back( -12);
510 for (
int i = 0; i < int(leg3.size()); ++i) iNewSys2.push_back( leg3[i]);
513 colConfig.erase(iSys);
514 colConfig.insert( iNewSys1, event);
515 colConfig.insert( iNewSys2, event);
526 bool RHadrons::openClosedLoop( ColConfig& colConfig,
Event& event) {
531 for (
int i = 0; i < systemPtr->size(); ++i) {
532 int iGNow = systemPtr->iParton[i];
533 if (event[iGNow].
id() == 21) {
534 double eGnow =
event[iGNow].p() *
event[iBef].p();
541 if (iGspl == -1)
return false;
544 int iG = systemPtr->iParton[iGspl];
545 int idNewQ = flavSelPtr->pickLightQ();
546 int iNewQ =
event.append( idNewQ, 101, iG, 0, 0, 0, event[iG].col(), 0,
547 0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
548 int iNewQb =
event.append( -idNewQ, 101, iG, 0, 0, 0, 0, event[iG].acol(),
549 0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
550 event[iG].statusNeg();
551 event[iG].daughters( iNewQ, iNewQb);
554 int iNext = iGspl + 1;
555 if (iNext == systemPtr->size()) iNext = 0;
556 if (event[ systemPtr->iParton[iNext]].acol() !=
event[iNewQ].col())
557 swap( iNewQ, iNewQb);
559 iNewSys.push_back( iNewQ);
560 for (
int i = iGspl + 1; i < systemPtr->size(); ++i)
561 iNewSys.push_back( systemPtr->iParton[i]);
562 for (
int i = 0; i < iGspl; ++i)
563 iNewSys.push_back( systemPtr->iParton[i]);
564 iNewSys.push_back( iNewQb);
567 colConfig.erase(iSys);
568 colConfig.insert( iNewSys, event);
580 bool RHadrons::splitSystem( ColConfig& colConfig,
Event& event) {
585 for (
int i = 0; i < int(systemPtr->size()); ++i) {
586 int iTmp = systemPtr->iParton[i];
587 if ( givesRHadron( event[iTmp].
id()) ) {
588 if (iFirst == -1) iFirst = i;
592 int nLeg = iSecond - iFirst;
595 int idNewQ = flavSelPtr->pickLightQ();
596 double mNewQ = particleDataPtr->constituentMass( idNewQ);
597 vector<int> iNewSys1, iNewSys2;
603 int i1Old = systemPtr->iParton[iFirst];
604 int i2Old = systemPtr->iParton[iSecond];
607 double mHat = (
event[i1Old].p() +
event[i2Old].p()).mCalc();
608 double mMax = mHat -
event[i1Old].m() -
event[i2Old].m();
609 if (mMax < 2. * (mNewQ + MSAFETY))
return false;
610 double mEff = min( 2. * (mNewQ + mOffsetCloudRH), mMax - 2. * MSAFETY);
611 double frac = mEff / mHat;
612 Vec4 pEff = frac * (
event[i1Old].p() +
event[i2Old].p());
617 if ( !newKin( event[i1Old].p(), event[i2Old].p(),
618 event[i1Old].m() / (1. - frac), event[i2Old].m() / (1. - frac),
619 p1New, p2New) )
return false;
624 int i1New, i2New, i3New, i4New;
625 int newCol =
event.nextColTag();
626 i1New =
event.copy( i1Old, 101);
627 if (event[i2Old].acol() == event[i1Old].col()) {
628 i3New =
event.append( -idNewQ, 101, i1Old, 0, 0, 0,
629 0, event[i2Old].acol(), 0.5 * pEff, 0.5 * mEff, 0.);
630 i2New =
event.copy( i2Old, 101);
631 event[i2New].acol( newCol);
632 i4New =
event.append( idNewQ, 101, i2Old, 0, 0, 0,
633 newCol, 0, 0.5 * pEff, 0.5 * mEff, 0.);
635 i3New =
event.append( idNewQ, 101, i1Old, 0, 0, 0,
636 event[i2Old].col(), 0, 0.5 * pEff, 0.5 * mEff, 0.);
637 i2New =
event.copy( i2Old, 101);
638 event[i2New].col( newCol);
639 i4New =
event.append( -idNewQ, 101, i2Old, 0, 0, 0,
640 0, newCol, 0.5 * pEff, 0.5 * mEff, 0.);
645 event[i1New].p( p1New);
646 event[i2New].p( p2New);
647 event[i1Old].daughters( i1New, i3New);
648 event[i1New].mother2( 0);
649 event[i2Old].daughters( i2New, i4New);
650 event[i2New].mother2( 0);
655 for (
int i = 0; i < iFirst; ++i)
656 iNewSys1.push_back( systemPtr->iParton[i]);
657 iNewSys1.push_back( i1New);
658 iNewSys1.push_back( i3New);
659 iNewSys2.push_back( i4New);
660 iNewSys2.push_back( i2New);
661 for (
int i = iSecond + 1; i < int(systemPtr->size()); ++i)
662 iNewSys2.push_back( systemPtr->iParton[i]);
666 }
else if (nLeg == 2) {
669 int iOld = systemPtr->iParton[iFirst + 1];
670 int i1New =
event.append( idNewQ, 101, iOld, 0, 0, 0,
671 event[iOld].col(), 0, 0.5 * event[iOld].p(),
672 0.5 * event[iOld].m(), 0.);
673 int i2New =
event.append( -idNewQ, 101, iOld, 0, 0, 0,
674 0, event[iOld].acol(), 0.5 * event[iOld].p(),
675 0.5 * event[iOld].m(), 0.);
676 event[iOld].statusNeg();
677 event[iOld].daughters( i1New, i2New);
680 if (event[systemPtr->iParton[iFirst]].col() ==
event[i2New].acol())
682 for (
int i = 0; i <= iFirst; ++i)
683 iNewSys1.push_back( systemPtr->iParton[i]);
684 iNewSys1.push_back( i1New);
685 iNewSys2.push_back( i2New);
686 for (
int i = iSecond; i < int(systemPtr->size()); ++i)
687 iNewSys2.push_back( systemPtr->iParton[i]);
698 for (
int i = iFirst + 1; i < iSecond - 1; ++i) {
699 int i1Tmp = systemPtr->iParton[i];
700 int i2Tmp = systemPtr->iParton[i + 1];
701 double mTmp = (
event[i1Tmp].p() +
event[i2Tmp].p()).mCalc();
709 double mEff = min( mNewQ + mOffsetCloudRH, 0.4 * mMin);
713 if ( !newKin( event[i1Old].p(), event[i2Old].p(),
714 mEff, mEff, p1New, p2New,
false) )
return false;
718 if (event[systemPtr->iParton[0]].acol() == 0) {
719 i1New =
event.append( -idNewQ, 101, i1Old, 0, 0, 0,
720 0, event[i1Old].acol(), p1New, mEff, 0.);
721 i2New =
event.append( idNewQ, 101, i2Old, 0, 0, 0,
722 event[i2Old].col(), 0, p2New, mEff, 0.);
724 i1New =
event.append( idNewQ, 101, i1Old, 0, 0, 0,
725 event[i1Old].col(), 0, p1New, mEff, 0.);
726 i2New =
event.append( -idNewQ, 101, i2Old, 0, 0, 0,
727 0, event[i2Old].acol(), p2New, mEff, 0.);
729 event[i1Old].statusNeg();
730 event[i2Old].statusNeg();
731 event[i1Old].daughters( i1New, 0);
732 event[i2Old].daughters( i2New, 0);
735 for (
int i = 0; i < iMin; ++i)
736 iNewSys1.push_back( systemPtr->iParton[i]);
737 iNewSys1.push_back( i1New);
738 iNewSys2.push_back( i2New);
739 for (
int i = iMin + 2; i < int(systemPtr->size()); ++i)
740 iNewSys2.push_back( systemPtr->iParton[i]);
744 colConfig.erase(iSys);
745 colConfig.insert( iNewSys1, event);
746 colConfig.insert( iNewSys2, event);
757 bool RHadrons::produceSquark( ColConfig& colConfig,
Event& event) {
766 int idAbsTop =
event[ systemPtr->iParton[0] ].idAbs();
767 bool sqAtTop = (allowRSb && idAbsTop == idRSb)
768 || (allowRSt && idAbsTop == idRSt);
771 int iBeg =
event.size();
772 iCreRHad[iRHad] = iBeg;
773 if (sqAtTop)
for (
int i = 0; i < systemPtr->size(); ++i)
774 event.copy( systemPtr->iParton[i], 102);
775 else for (
int i = systemPtr->size() - 1; i >= 0; --i)
776 event.copy( systemPtr->iParton[i], 102);
777 int iEnd =
event.size() - 1;
780 int idOldH =
event[iBeg].id();
781 int idOldL =
event[iEnd].id();
784 FlavContainer flavOld( idOldH%10);
785 int idNewQ = flavSelPtr->pick(flavOld).id;
786 int idRHad = toIdWithSquark( idOldH, idNewQ);
788 infoPtr->errorMsg(
"Error in RHadrons::produceSquark: "
789 "cannot form R-hadron code");
794 double mRHad = particleDataPtr->m0(idRHad) +
event[iBeg].m()
795 - ( (abs(idOldH) == idRSb) ? m0Sb : m0St );
796 double z = zSelPtr->zFrag( idOldH, idNewQ, mRHad*mRHad);
799 Vec4 pOldH =
event[iBeg].p();
800 int iOldL = iBeg + 1;
801 Vec4 pOldL =
event[iOldL].p();
802 double mOldL =
event[iOldL].m();
803 double mNewH = mRHad / z;
804 double sSys = (pOldH + pOldL).m2Calc();
805 double sRem = (1. - z) * (sSys - mNewH*mNewH);
806 double sMin = pow2(mOldL + mCollapseRH);
809 while ( ( sRem < sMin || sSys < pow2(mNewH + mOldL + MSAFETY) )
812 pOldL +=
event[iOldL].p();
813 mOldL =
event[iOldL].m();
814 sSys = (pOldH + pOldL).m2Calc();
815 sRem = (1. - z) * (sSys - mNewH*mNewH);
816 sMin = pow2(mOldL + mCollapseRH);
820 if ( sRem > sMin && sSys > pow2(mNewH + mOldL + MSAFETY) ) {
822 if ( !newKin( pOldH, pOldL, mNewH, mOldL, pNewH, pNewL) ) {
823 infoPtr->errorMsg(
"Error in RHadrons::produceSquark: "
824 "failed to construct kinematics with reduced system");
829 iRNow =
event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
830 z * pNewH, mRHad, 0.);
834 bool hasCol = (idNewQ > 0 && idNewQ < 10) || idNewQ < -10;
835 int col = (hasCol) ? event[iOldL].acol() : 0;
836 int acol = (hasCol) ? 0 : event[iOldL].col();
837 iNewQ =
event.append( idNewQ, 105, iBeg, iOldL, 0, 0, col, acol,
838 (1. - z) * pNewH, (1. - z) * mNewH, 0.);
839 iNewL =
event.copy( iOldL, 105);
840 event[iNewL].mothers( iBeg, iOldL);
841 event[iNewL].p( pNewL);
849 FlavContainer flav1( idOldL);
850 FlavContainer flav2( -idNewQ);
852 int idNewL = flavSelPtr->combine( flav1, flav2);
853 while (++iTry < NTRYMAX && idNewL == 0)
854 idNewL = flavSelPtr->combine( flav1, flav2);
856 infoPtr->errorMsg(
"Error in RHadrons::produceSquark: "
857 "cannot form light hadron code");
860 double mNewL = particleDataPtr->mSel( idNewL);
863 if ( sSys > pow2(mRHad + mNewL + MSAFETY) ) {
865 if ( !newKin( pOldH, pOldL, mRHad, mNewL, pRHad, pNewL) ) {
866 infoPtr->errorMsg(
"Error in RHadrons::produceSquark: "
867 "failed to construct kinematics for two-hadron decay");
872 iRNow =
event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
874 event.append( idNewL, 105, iBeg, iOldL, 0, 0, 0, 0,
884 idRHad = toIdWithSquark( idOldH, idOldL);
886 infoPtr->errorMsg(
"Error in RHadrons::produceSquark: "
887 "cannot form R-hadron code");
892 iRNow =
event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
893 systemPtr->pSum, systemPtr->mass, 0.);
900 iRHadron[iRHad] = iRNow;
901 int iLast =
event.size() - 1;
902 for (
int i = iBeg; i <= iOldL; ++i) {
903 event[i].statusNeg();
904 event[i].daughters( iRNow, iLast);
908 colConfig.erase(iSys);
911 iNewSys.push_back( iNewQ);
912 iNewSys.push_back( iNewL);
913 for (
int i = iOldL + 1; i <= iEnd; ++i) iNewSys.push_back( i);
914 colConfig.insert( iNewSys, event);
918 event[iRNow].tau( event[iBef].tau() );
919 if (event[iBef].hasVertex())
event[iRNow].vProd( event[iBef].vProd() );
930 bool RHadrons::produceGluino( ColConfig& colConfig,
Event& event) {
937 vector<int> iSide1, iSide2, iSideTmp, iNewSys1, iNewSys2;
938 Vec4 pGlui, pSide1, pSide2;
941 bool isGBall = (rndmPtr->flat() < probGluinoballRH);
944 for (
int i = 0; i < systemPtr->size(); ++i) {
945 int iTmp = systemPtr->iParton[i];
946 if (iGlui == 0 && event[ iTmp ].
id() == idRGo) {
948 pGlui =
event[ iTmp ].p();
949 }
else if (iGlui == 0) {
950 iSideTmp.push_back( iTmp);
951 pSide1 +=
event[ iTmp ].p();
953 iSide2.push_back( iTmp);
954 pSide2 +=
event[ iTmp ].p();
957 int iGluiFix = iGlui;
960 for (
int i =
int(iSideTmp.size()) - 1; i >= 0; --i)
961 iSide1.push_back( iSideTmp[i]);
962 double m1H = (pGlui + pSide1).mCalc() -
event[ iSide1.back() ].m();
963 double m2H = (pGlui + pSide2).mCalc() -
event[ iSide2.back() ].m();
965 swap( iSide1, iSide2);
966 swap( pSide1, pSide2);
970 for (
int iSide = 1; iSide < 3; ++iSide) {
982 int iBeg =
event.size();
984 iCreRHad[iRHad] = iBeg;
985 event.copy( iGlui, 102);
986 for (
int i = 0; i < int(iSide1.size()); ++i)
987 event.copy( iSide1[i], 102);
989 event.copy( iGlui, 102);
990 for (
int i = 0; i < int(iSide2.size()); ++i)
991 event.copy( iSide2[i], 102);
993 int iEnd =
event.size() - 1;
996 int idOldL =
event[iEnd].id();
1000 FlavContainer flavOld( idOldL);
1001 idNewQ = -flavSelPtr->pick(flavOld).id;
1002 }
while (iSide == 2 && isDiq1 && abs(idNewQ) > 10);
1003 if (iSide == 1) isDiq1 = (abs(idNewQ) > 10);
1005 bool hasCol = (
event[iEnd].col() > 0);
1012 idRHad = (hasCol) ? 1009002 : -1009002;
1013 if (hasCol) colR =
event[iBeg].col();
1014 else acolR =
event[iBeg].acol();
1016 double mNewQ = particleDataPtr->constituentMass( idNewQ);
1017 if (isGBall) mNewQ *= 0.5;
1018 mRHad =
event[iBeg].m() + mOffsetCloudRH + mNewQ;
1020 idRHad = toIdWithGluino( idSave, idNewQ);
1022 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1023 "cannot form R-hadron code");
1027 mRHad = particleDataPtr->m0(idRHad) +
event[iGluiFix].m() - m0Go;
1031 double z = zSelPtr->zFrag( idRGo, idNewQ, mRHad*mRHad);
1034 Vec4 pOldH =
event[iBeg].p();
1035 int iOldL = iBeg + 1;
1036 Vec4 pOldL =
event[iOldL].p();
1037 double mOldL =
event[iOldL].m();
1038 double mNewH = mRHad / z;
1039 double sSys = (pOldH + pOldL).m2Calc();
1040 double sRem = (1. - z) * (sSys - mNewH*mNewH);
1041 double sMin = pow2(mOldL + mCollapseRH);
1044 while ( ( sRem < sMin || sSys < pow2(mNewH + mOldL + MSAFETY) )
1047 pOldL +=
event[iOldL].p();
1048 mOldL =
event[iOldL].m();
1049 sSys = (pOldH + pOldL).m2Calc();
1050 sRem = (1. - z) * (sSys - mNewH*mNewH);
1051 sMin = pow2(mOldL + mCollapseRH);
1055 if ( sRem > sMin && sSys > pow2(mNewH + mOldL + MSAFETY) ) {
1057 if ( !newKin( pOldH, pOldL, mNewH, mOldL, pNewH, pNewL) ) {
1058 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1059 "failed to construct kinematics with reduced system");
1064 iRNow =
event.append( idRHad, statusRHad, iBeg, iOldL,
1065 0, 0, colR, acolR, z * pNewH, mRHad, 0.);
1068 if (!isGBall) idNewQ = -idNewQ;
1069 int colN = (hasCol) ? 0 : event[iOldL].acol();
1070 int acolN = (hasCol) ? event[iOldL].col() : 0;
1071 iNewQ =
event.append( idNewQ, 105, iBeg, iOldL, 0, 0,
1072 colN, acolN, (1. - z) * pNewH, (1. - z) * mNewH, 0.);
1073 iNewL =
event.copy( iOldL, 105);
1074 event[iNewL].mothers( iBeg, iOldL);
1075 event[iNewL].p( pNewL);
1079 iNewSys1.push_back( iNewQ);
1080 iNewSys1.push_back( iNewL);
1081 for (
int i = iOldL + 1; i <= iEnd; ++i) iNewSys1.push_back( i);
1083 iNewSys2.push_back( iNewQ);
1084 iNewSys2.push_back( iNewL);
1085 for (
int i = iOldL + 1; i <= iEnd; ++i) iNewSys2.push_back( i);
1093 if (nBody == 0 && isGBall && iSide == 1) {
1094 idQLeap =
event[iOldL].id();
1095 Vec4 pNewH =
event[iBeg].p() + pOldL;
1098 iRNow =
event.append( idRHad, statusRHad, iBeg, iEnd,
1099 0, 0, colR, acolR, pNewH, pNewH.mCalc(), 0.);
1105 if (nBody == 0 && (!isGBall || (iSide == 2 && idQLeap != 0) )) {
1106 if (isGBall) idNewQ = -idQLeap;
1107 FlavContainer flav1( idOldL);
1108 FlavContainer flav2( -idNewQ);
1110 int idNewL = flavSelPtr->combine( flav1, flav2);
1111 while (++iTry < NTRYMAX && idNewL == 0)
1112 idNewL = flavSelPtr->combine( flav1, flav2);
1114 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1115 "cannot form light hadron code");
1118 double mNewL = particleDataPtr->mSel( idNewL);
1121 if ( sSys > pow2(mRHad + mNewL + MSAFETY) ) {
1123 if ( !newKin( pOldH, pOldL, mRHad, mNewL, pRHad, pNewL) ) {
1124 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1125 "failed to construct kinematics for two-hadron decay");
1130 iRNow =
event.append( idRHad, statusRHad, iBeg, iOldL, 0, 0,
1131 colR, acolR, pRHad, mRHad, 0.);
1132 event.append( idNewL, 105, iBeg, iOldL, 0, 0, 0, 0,
1143 if (nBody == 0 && (!isGBall || (iSide == 2 && idQLeap != 0) )) {
1144 if (isGBall) idSave = idQLeap;
1145 if (iSide == 1) idSave = idOldL;
1146 else idRHad = toIdWithGluino( idSave, idOldL);
1148 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1149 "cannot form R-hadron code");
1154 iRNow =
event.append( idRHad, statusRHad, iBeg, iOldL, 0, 0,
1155 colR, acolR, pOldH + pOldL, (pOldH + pOldL).mCalc(), 0.);
1163 int iLast =
event.size() - 1;
1164 for (
int i = iBeg; i <= iOldL; ++i) {
1165 event[i].statusNeg();
1166 event[i].daughters( iRNow, iLast);
1175 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1176 "could not handle gluinoball kinematics");
1179 iRHadron[iRHad] = iGlui;
1180 colConfig.erase(iSys);
1184 if (iNewSys1.size() > 0) colConfig.insert( iNewSys1, event);
1185 if (iNewSys2.size() > 0) colConfig.insert( iNewSys2, event);
1189 }
else if (idQLeap == 0) {
1190 int iG1 = iNewSys1[0];
1191 int iG2 = iNewSys2[0];
1192 int colG =
event[iG1].col() +
event[iG2].col();
1193 int acolG =
event[iG1].acol() +
event[iG2].acol();
1194 Vec4 pG =
event[iG1].p() +
event[iG2].p();
1195 int iG12 =
event.append( 21, 107, iG1, iG2, 0, 0, colG, acolG,
1196 pG, pG.mCalc(), 0.);
1201 event[iG1].statusNeg();
1202 event[iG2].statusNeg();
1203 event[iG1].daughter1( iG12);
1204 event[iG2].daughter1( iG12);
1205 int colBridge =
event.nextColTag();
1206 if (event[iG1].col() == 0) {
1207 event[iG1].col( colBridge);
1208 event[iG2].acol( colBridge);
1210 event[iG1].acol( colBridge);
1211 event[iG2].col( colBridge);
1215 vector<int> iNewSys12;
1216 for (
int i =
int(iNewSys1.size()) - 1; i > 0; --i)
1217 iNewSys12.push_back( iNewSys1[i]);
1218 iNewSys12.push_back( iG12);
1219 for (
int i = 1; i < int(iNewSys2.size()); ++i)
1220 iNewSys12.push_back( iNewSys2[i]);
1221 colConfig.insert( iNewSys12, event);
1226 int iG2 = iNewSys2[0];
1227 event[iG2].id( idQLeap);
1228 colConfig.insert( iNewSys2, event);
1232 event[iGlui].tau( event[iBef].tau() );
1233 if (event[iBef].hasVertex())
event[iGlui].vProd( event[iBef].vProd() );
1245 int RHadrons::toIdWithSquark(
int id1,
int id2) {
1248 int id1Abs = abs(id1);
1249 int id2Abs = abs(id2);
1250 if (id2Abs < 10 && id1 > 0 && id2 > 0)
return 0;
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;
1256 bool isSt = (id1Abs == idRSt);
1257 int idRHad = 1000000;
1258 if (id2Abs < 10) idRHad += ((isSt) ? 600 : 500) + 10 * id2Abs + 2;
1259 else idRHad += ((isSt) ? 6000 : 5000) + 10 * (id2Abs/100) + id2Abs%10;
1260 if (id1 < 0) idRHad = -idRHad;
1272 pair<int,int> RHadrons::fromIdWithSquark(
int idRHad) {
1275 int idLight = (abs(idRHad) - 1000000) / 10;
1276 int idSq = (idLight < 100) ? idLight/10 : idLight/100;
1277 int id1 = (idSq == 6) ? idRSt : idRSb;
1278 if (idRHad < 0) id1 = -id1;
1281 int id2 = (idLight < 100) ? idLight%10 : idLight%100;
1282 if (id2 > 10) id2 = 100 * id2 + abs(idRHad)%10;
1283 if ((id2 < 10 && idRHad > 0) || (id2 > 10 && idRHad < 0)) id2 = -id2;
1286 return make_pair( id1, id2);
1294 int RHadrons::toIdWithGluino(
int id1,
int id2) {
1297 int id1Abs = abs(id1);
1298 int id2Abs = abs(id2);
1299 if (id1Abs == 21 && id2Abs == 21)
return 1000993;
1300 int idMax = max( id1Abs, id2Abs);
1301 int idMin = min( id1Abs, id2Abs);
1302 if (idMin > 10)
return 0;
1303 if (idMax > 10 && id1 > 0 && id2 < 0)
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;
1311 idRHad = 1009003 + 100 * idMax + 10 * idMin;
1312 if (idMin != idMax && idMax%2 == 1) {
1313 if (id1Abs == idMax && id1 > 0) idRHad = -idRHad;
1314 if (id2Abs == idMax && id2 > 0) idRHad = -idRHad;
1316 if (idMin != idMax && idMax%2 == 0) {
1317 if (id1Abs == idMax && id1 < 0) idRHad = -idRHad;
1318 if (id2Abs == idMax && id2 < 0) idRHad = -idRHad;
1323 int idA = idMax/1000;
1324 int idB = (idMax/100)%10;
1326 if (idC > idB) swap( idB, idC);
1327 if (idB > idA) swap( idA, idB);
1328 if (idC > idB) swap( idB, idC);
1329 idRHad = 1090004 + 1000 * idA + 100 * idB + 10 * idC;
1330 if (id1 < 0) idRHad = -idRHad;
1343 pair<int,int> RHadrons::fromIdWithGluino(
int idRHad) {
1346 int idLight = (abs(idRHad) - 1000000) / 10;
1347 int id1, id2, idTmp, idA, idB, idC;
1350 if (idLight < 100) {
1351 id1 = (rndmPtr->flat() < 0.5) ? 1 : 2;
1355 }
else if (idLight < 1000) {
1356 id1 = (idLight / 10) % 10;
1357 id2 = -(idLight % 10);
1368 idA = (idLight / 100) % 10;
1369 idB = (idLight / 10) % 10;
1371 double rndmQ = 3. * rndmPtr->flat();
1372 if (idA > 3) rndmQ = 0.5;
1375 id2 = 1000 * idB + 100 * idC + 3;
1376 if (idB != idC && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2;
1377 }
else if (rndmQ < 2.) {
1379 id2 = 1000 * idA + 100 * idC + 3;
1380 if (idA != idC && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2;
1383 id2 = 1000 * idA + 100 * idB +3;
1384 if (idA != idB && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2;
1396 return make_pair( id1, id2);
1406 bool RHadrons::newKin( Vec4 pOld1, Vec4 pOld2,
double mNew1,
double mNew2,
1407 Vec4& pNew1, Vec4& pNew2,
bool checkMargin) {
1410 double sSum = (pOld1 + pOld2).m2Calc();
1411 double sOld1 = pOld1.m2Calc();
1412 double sOld2 = pOld2.m2Calc();
1413 double sNew1 = mNew1 * mNew1;
1414 double sNew2 = mNew2 * mNew2;
1417 if (checkMargin && pow2(mNew1 + mNew2 + MSAFETY) > sSum)
return false;
1420 double lamOld = sqrt( pow2(sSum - sOld1 - sOld2) - 4. * sOld1 * sOld2 );
1421 double lamNew = sqrt( pow2(sSum - sNew1 - sNew2) - 4. * sNew1 * sNew2 );
1422 double move1 = (lamNew * (sSum - sOld1 + sOld2)
1423 - lamOld * (sSum - sNew1 + sNew2)) / (2. * sSum * lamOld);
1424 double move2 = (lamNew * (sSum + sOld1 - sOld2)
1425 - lamOld * (sSum + sNew1 - sNew2)) / (2. * sSum * lamOld);
1428 pNew1 = (1. + move1) * pOld1 - move2 * pOld2;
1429 pNew2 = (1. + move2) * pOld2 - move1 * pOld1;