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( Info* infoPtrIn, Settings& settings,
51 ParticleData* particleDataPtrIn, Rndm* rndmPtrIn) {
55 particleDataPtr = particleDataPtrIn;
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");
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;
81 m0Sb = particleDataPtr->m0(idRSb);
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);
88 m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
89 particleDataPtr->m0( idR, m0RHad);
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);
104 m0St = particleDataPtr->m0(idRSt);
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);
111 m0RHad += particleDataPtr->constituentMass( (idR%1000)/100 );
112 particleDataPtr->m0( idR, m0RHad);
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);
127 m0Go = particleDataPtr->m0(idRGo);
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);
137 m0RHad += particleDataPtr->constituentMass( (idR%10000)/1000 );
138 particleDataPtr->m0( idR, m0RHad);
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);
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;
171 bool RHadrons::produce( ColConfig& colConfig,
Event& event) {
174 if (!allowSomeR)
return true;
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);
193 nRHad = iRHadron.size();
196 if (nRHad == 0)
return true;
200 infoPtr->errorMsg(
"Error in RHadrons::produce: "
201 "cannot handle more than two R-hadrons");
204 if (nRHad > 1 && rndmPtr->flat() > 0.5) swap( iBefRHad[0], iBefRHad[1]);
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");
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");
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");
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");
248 int iSys1 = colConfig.findSinglet( iBefRHad[0]);
249 int iSys2 = colConfig.findSinglet( iBefRHad[1]);
250 if (iSys2 == 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");
262 for (iRHad = 0; iRHad < nRHad; ++iRHad) {
263 iBef = iBefRHad[iRHad];
264 iSys = colConfig.findSinglet( iBef);
266 infoPtr->errorMsg(
"Error in RHadrons::produce: "
267 "sparticle not in any colour singlet");
270 systemPtr = &colConfig[iSys];
273 if (systemPtr->hasJunction) {
274 infoPtr->errorMsg(
"Error in RHadrons::produce: "
275 "cannot handle system with junction");
278 if (systemPtr->isClosed) {
279 infoPtr->errorMsg(
"Error in RHadrons::produce: "
280 "cannot handle closed colour loop");
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;
301 bool RHadrons::decay(
Event& event) {
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();
314 pair<int,int> idPair = (isTriplet[iRHad])
315 ? fromIdWithSquark( idRHad) : fromIdWithGluino( idRHad);
316 int id1 = idPair.first;
317 int id2 = idPair.second;
321 double fracR = mRBef / mRHad;
323 infoPtr->errorMsg(
"Error in RHadrons::decay: "
324 "too low R-hadron mass for decay");
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;
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.);
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);
348 int col1 =
event.nextColTag();
349 int col2 =
event.nextColTag();
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.);
361 event[iRNow].statusNeg();
362 event[iRNow].daughters( iR0, iR2);
363 iAftRHad[iRHad] = iR0;
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);
383 bool RHadrons::splitOffJunction( ColConfig& colConfig,
Event& event) {
386 vector<int> leg1, leg2, leg3;
391 for (
int i = 0; i < systemPtr->size(); ++i) {
393 int iP = systemPtr->iParton[i];
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);
405 if (iLegSP == 0)
return false;
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;
416 Vec4 pSP =
event[iBef].p();
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;
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,
434 leg1.insert( leg1.begin(), iNewG);
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);
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);
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);
457 int iNewSP =
event.copy( iBef, 101);
458 event[iNewSP].mother2(0);
459 event[iBef].daughter1(iNewG);
460 event[iNewSP].rotbst( MSP);
461 leg1[iPosSP] = iNewSP;
462 if (iBefRHad[0] == iBef) iBefRHad[0] = iNewSP;
463 else if (nRHad > 1 && iBefRHad[1] == iBef) iBefRHad[1] = iNewSP;
465 for (
int i = 0; i < int(leg2.size()); ++i) {
466 int iCopy =
event.copy( leg2[i], 101);
467 event[iCopy].rotbst( MRec);
468 if (iBefRHad[0] == leg2[i]) iBefRHad[0] = iCopy;
469 else if (nRHad > 1 && iBefRHad[1] == leg2[i]) iBefRHad[1] = iCopy;
472 for (
int i = 0; i < int(leg3.size()); ++i) {
473 int iCopy =
event.copy( leg3[i], 101);
474 event[iCopy].rotbst( MRec);
475 if (iBefRHad[0] == leg3[i]) iBefRHad[0] = iCopy;
476 else if (nRHad > 1 && iBefRHad[1] == leg3[i]) iBefRHad[1] = iCopy;
485 double eGspl =
event[leg1[0]].p() *
event[iBef].p();
486 for (
int i = 1; i < iPosSP; ++i) {
487 double eGnow =
event[leg1[i]].p() *
event[iBef].p();
493 int iG = leg1[iGspl];
496 int idNewQ = flavSelPtr->pickLightQ();
497 int iNewQ =
event.append( idNewQ, 101, iG, 0, 0, 0, event[iG].col(), 0,
498 0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
499 int iNewQb =
event.append( -idNewQ, 101, iG, 0, 0, 0, 0, event[iG].acol(),
500 0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
501 event[iG].statusNeg();
502 event[iG].daughters( iNewQ, iNewQb);
503 if (event[leg1.back()].col() == 0) swap( iNewQ, iNewQb);
506 vector<int> iNewSys1, iNewSys2;
507 iNewSys1.push_back( iNewQb);
508 for (
int i = iGspl + 1; i < int(leg1.size()); ++i)
509 iNewSys1.push_back( leg1[i]);
510 iNewSys2.push_back( -10);
511 for (
int i = 0; i < iGspl; ++i) iNewSys2.push_back( leg1[i]);
512 iNewSys2.push_back( iNewQ);
513 iNewSys2.push_back( -11);
514 for (
int i = 0; i < int(leg2.size()); ++i) iNewSys2.push_back( leg2[i]);
515 iNewSys2.push_back( -12);
516 for (
int i = 0; i < int(leg3.size()); ++i) iNewSys2.push_back( leg3[i]);
519 colConfig.erase(iSys);
520 colConfig.insert( iNewSys1, event);
521 colConfig.insert( iNewSys2, event);
532 bool RHadrons::openClosedLoop( ColConfig& colConfig,
Event& event) {
537 for (
int i = 0; i < systemPtr->size(); ++i) {
538 int iGNow = systemPtr->iParton[i];
539 if (event[iGNow].
id() == 21) {
540 double eGnow =
event[iGNow].p() *
event[iBef].p();
547 if (iGspl == -1)
return false;
550 int iG = systemPtr->iParton[iGspl];
551 int idNewQ = flavSelPtr->pickLightQ();
552 int iNewQ =
event.append( idNewQ, 101, iG, 0, 0, 0, event[iG].col(), 0,
553 0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
554 int iNewQb =
event.append( -idNewQ, 101, iG, 0, 0, 0, 0, event[iG].acol(),
555 0.5 * event[iG].p(), 0.5 * event[iG].m(), 0.);
556 event[iG].statusNeg();
557 event[iG].daughters( iNewQ, iNewQb);
560 int iNext = iGspl + 1;
561 if (iNext == systemPtr->size()) iNext = 0;
562 if (event[ systemPtr->iParton[iNext]].acol() !=
event[iNewQ].col())
563 swap( iNewQ, iNewQb);
565 iNewSys.push_back( iNewQ);
566 for (
int i = iGspl + 1; i < systemPtr->size(); ++i)
567 iNewSys.push_back( systemPtr->iParton[i]);
568 for (
int i = 0; i < iGspl; ++i)
569 iNewSys.push_back( systemPtr->iParton[i]);
570 iNewSys.push_back( iNewQb);
573 colConfig.erase(iSys);
574 colConfig.insert( iNewSys, event);
586 bool RHadrons::splitSystem( ColConfig& colConfig,
Event& event) {
591 for (
int i = 0; i < int(systemPtr->size()); ++i) {
592 int iTmp = systemPtr->iParton[i];
593 if ( givesRHadron( event[iTmp].
id()) ) {
594 if (iFirst == -1) iFirst = i;
598 int nLeg = iSecond - iFirst;
601 int idNewQ = flavSelPtr->pickLightQ();
602 double mNewQ = particleDataPtr->constituentMass( idNewQ);
603 vector<int> iNewSys1, iNewSys2;
609 int i1Old = systemPtr->iParton[iFirst];
610 int i2Old = systemPtr->iParton[iSecond];
613 double mHat = (
event[i1Old].p() +
event[i2Old].p()).mCalc();
614 double mMax = mHat -
event[i1Old].m() -
event[i2Old].m();
615 if (mMax < 2. * (mNewQ + MSAFETY))
return false;
616 double mEff = min( 2. * (mNewQ + mOffsetCloudRH), mMax - 2. * MSAFETY);
617 double frac = mEff / mHat;
618 Vec4 pEff = frac * (
event[i1Old].p() +
event[i2Old].p());
623 if ( !newKin( event[i1Old].p(), event[i2Old].p(),
624 event[i1Old].m() / (1. - frac), event[i2Old].m() / (1. - frac),
625 p1New, p2New) )
return false;
630 int i1New, i2New, i3New, i4New;
631 int newCol =
event.nextColTag();
632 i1New =
event.copy( i1Old, 101);
633 if (event[i2Old].acol() == event[i1Old].col()) {
634 i3New =
event.append( -idNewQ, 101, i1Old, 0, 0, 0,
635 0, event[i2Old].acol(), 0.5 * pEff, 0.5 * mEff, 0.);
636 i2New =
event.copy( i2Old, 101);
637 event[i2New].acol( newCol);
638 i4New =
event.append( idNewQ, 101, i2Old, 0, 0, 0,
639 newCol, 0, 0.5 * pEff, 0.5 * mEff, 0.);
641 i3New =
event.append( idNewQ, 101, i1Old, 0, 0, 0,
642 event[i2Old].col(), 0, 0.5 * pEff, 0.5 * mEff, 0.);
643 i2New =
event.copy( i2Old, 101);
644 event[i2New].col( newCol);
645 i4New =
event.append( -idNewQ, 101, i2Old, 0, 0, 0,
646 0, newCol, 0.5 * pEff, 0.5 * mEff, 0.);
651 event[i1New].p( p1New);
652 event[i2New].p( p2New);
653 event[i1Old].daughters( i1New, i3New);
654 event[i1New].mother2( 0);
655 event[i2Old].daughters( i2New, i4New);
656 event[i2New].mother2( 0);
661 for (
int i = 0; i < iFirst; ++i)
662 iNewSys1.push_back( systemPtr->iParton[i]);
663 iNewSys1.push_back( i1New);
664 iNewSys1.push_back( i3New);
665 iNewSys2.push_back( i4New);
666 iNewSys2.push_back( i2New);
667 for (
int i = iSecond + 1; i < int(systemPtr->size()); ++i)
668 iNewSys2.push_back( systemPtr->iParton[i]);
672 }
else if (nLeg == 2) {
675 int iOld = systemPtr->iParton[iFirst + 1];
676 int i1New =
event.append( idNewQ, 101, iOld, 0, 0, 0,
677 event[iOld].col(), 0, 0.5 * event[iOld].p(),
678 0.5 * event[iOld].m(), 0.);
679 int i2New =
event.append( -idNewQ, 101, iOld, 0, 0, 0,
680 0, event[iOld].acol(), 0.5 * event[iOld].p(),
681 0.5 * event[iOld].m(), 0.);
682 event[iOld].statusNeg();
683 event[iOld].daughters( i1New, i2New);
686 if (event[systemPtr->iParton[iFirst]].col() ==
event[i2New].acol())
688 for (
int i = 0; i <= iFirst; ++i)
689 iNewSys1.push_back( systemPtr->iParton[i]);
690 iNewSys1.push_back( i1New);
691 iNewSys2.push_back( i2New);
692 for (
int i = iSecond; i < int(systemPtr->size()); ++i)
693 iNewSys2.push_back( systemPtr->iParton[i]);
704 for (
int i = iFirst + 1; i < iSecond - 1; ++i) {
705 int i1Tmp = systemPtr->iParton[i];
706 int i2Tmp = systemPtr->iParton[i + 1];
707 double mTmp = (
event[i1Tmp].p() +
event[i2Tmp].p()).mCalc();
715 double mEff = min( mNewQ + mOffsetCloudRH, 0.4 * mMin);
719 if ( !newKin( event[i1Old].p(), event[i2Old].p(),
720 mEff, mEff, p1New, p2New,
false) )
return false;
724 if (event[systemPtr->iParton[0]].acol() == 0) {
725 i1New =
event.append( -idNewQ, 101, i1Old, 0, 0, 0,
726 0, event[i1Old].acol(), p1New, mEff, 0.);
727 i2New =
event.append( idNewQ, 101, i2Old, 0, 0, 0,
728 event[i2Old].col(), 0, p2New, mEff, 0.);
730 i1New =
event.append( idNewQ, 101, i1Old, 0, 0, 0,
731 event[i1Old].col(), 0, p1New, mEff, 0.);
732 i2New =
event.append( -idNewQ, 101, i2Old, 0, 0, 0,
733 0, event[i2Old].acol(), p2New, mEff, 0.);
735 event[i1Old].statusNeg();
736 event[i2Old].statusNeg();
737 event[i1Old].daughters( i1New, 0);
738 event[i2Old].daughters( i2New, 0);
741 for (
int i = 0; i < iMin; ++i)
742 iNewSys1.push_back( systemPtr->iParton[i]);
743 iNewSys1.push_back( i1New);
744 iNewSys2.push_back( i2New);
745 for (
int i = iMin + 2; i < int(systemPtr->size()); ++i)
746 iNewSys2.push_back( systemPtr->iParton[i]);
750 colConfig.erase(iSys);
751 colConfig.insert( iNewSys1, event);
752 colConfig.insert( iNewSys2, event);
763 bool RHadrons::produceSquark( ColConfig& colConfig,
Event& event) {
772 int idAbsTop =
event[ systemPtr->iParton[0] ].idAbs();
773 bool sqAtTop = (allowRSb && idAbsTop == idRSb)
774 || (allowRSt && idAbsTop == idRSt);
777 int iBeg =
event.size();
778 iCreRHad[iRHad] = iBeg;
779 if (sqAtTop)
for (
int i = 0; i < systemPtr->size(); ++i)
780 event.copy( systemPtr->iParton[i], 102);
781 else for (
int i = systemPtr->size() - 1; i >= 0; --i)
782 event.copy( systemPtr->iParton[i], 102);
783 int iEnd =
event.size() - 1;
786 int idOldH =
event[iBeg].id();
787 int idOldL =
event[iEnd].id();
790 FlavContainer flavOld( idOldH%10);
791 int idNewQ = flavSelPtr->pick(flavOld).id;
792 int idRHad = toIdWithSquark( idOldH, idNewQ);
794 infoPtr->errorMsg(
"Error in RHadrons::produceSquark: "
795 "cannot form R-hadron code");
800 double mRHad = particleDataPtr->m0(idRHad) +
event[iBeg].m()
801 - ( (abs(idOldH) == idRSb) ? m0Sb : m0St );
802 double z = zSelPtr->zFrag( idOldH, idNewQ, mRHad*mRHad);
805 Vec4 pOldH =
event[iBeg].p();
806 int iOldL = iBeg + 1;
807 Vec4 pOldL =
event[iOldL].p();
808 double mOldL =
event[iOldL].m();
809 double mNewH = mRHad / z;
810 double sSys = (pOldH + pOldL).m2Calc();
811 double sRem = (1. - z) * (sSys - mNewH*mNewH);
812 double sMin = pow2(mOldL + mCollapseRH);
815 while ( ( sRem < sMin || sSys < pow2(mNewH + mOldL + MSAFETY) )
818 pOldL +=
event[iOldL].p();
819 mOldL =
event[iOldL].m();
820 sSys = (pOldH + pOldL).m2Calc();
821 sRem = (1. - z) * (sSys - mNewH*mNewH);
822 sMin = pow2(mOldL + mCollapseRH);
826 if ( sRem > sMin && sSys > pow2(mNewH + mOldL + MSAFETY) ) {
828 if ( !newKin( pOldH, pOldL, mNewH, mOldL, pNewH, pNewL) ) {
829 infoPtr->errorMsg(
"Error in RHadrons::produceSquark: "
830 "failed to construct kinematics with reduced system");
835 iRNow =
event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
836 z * pNewH, mRHad, 0.);
840 bool hasCol = (idNewQ > 0 && idNewQ < 10) || idNewQ < -10;
841 int col = (hasCol) ? event[iOldL].acol() : 0;
842 int acol = (hasCol) ? 0 : event[iOldL].col();
843 iNewQ =
event.append( idNewQ, 105, iBeg, iOldL, 0, 0, col, acol,
844 (1. - z) * pNewH, (1. - z) * mNewH, 0.);
845 iNewL =
event.copy( iOldL, 105);
846 event[iNewL].mothers( iBeg, iOldL);
847 event[iNewL].p( pNewL);
855 FlavContainer flav1( idOldL);
856 FlavContainer flav2( -idNewQ);
858 int idNewL = flavSelPtr->combine( flav1, flav2);
859 while (++iTry < NTRYMAX && idNewL == 0)
860 idNewL = flavSelPtr->combine( flav1, flav2);
862 infoPtr->errorMsg(
"Error in RHadrons::produceSquark: "
863 "cannot form light hadron code");
866 double mNewL = particleDataPtr->mSel( idNewL);
869 if ( sSys > pow2(mRHad + mNewL + MSAFETY) ) {
871 if ( !newKin( pOldH, pOldL, mRHad, mNewL, pRHad, pNewL) ) {
872 infoPtr->errorMsg(
"Error in RHadrons::produceSquark: "
873 "failed to construct kinematics for two-hadron decay");
878 iRNow =
event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
880 event.append( idNewL, 105, iBeg, iOldL, 0, 0, 0, 0,
890 idRHad = toIdWithSquark( idOldH, idOldL);
892 infoPtr->errorMsg(
"Error in RHadrons::produceSquark: "
893 "cannot form R-hadron code");
898 iRNow =
event.append( idRHad, 104, iBeg, iOldL, 0, 0, 0, 0,
899 systemPtr->pSum, systemPtr->mass, 0.);
906 iRHadron[iRHad] = iRNow;
907 int iLast =
event.size() - 1;
908 for (
int i = iBeg; i <= iOldL; ++i) {
909 event[i].statusNeg();
910 event[i].daughters( iRNow, iLast);
914 colConfig.erase(iSys);
917 iNewSys.push_back( iNewQ);
918 iNewSys.push_back( iNewL);
919 for (
int i = iOldL + 1; i <= iEnd; ++i) iNewSys.push_back( i);
920 colConfig.insert( iNewSys, event);
924 event[iRNow].tau( event[iBef].tau() );
925 if (event[iBef].hasVertex())
event[iRNow].vProd( event[iBef].vProd() );
936 bool RHadrons::produceGluino( ColConfig& colConfig,
Event& event) {
943 vector<int> iSide1, iSide2, iSideTmp, iNewSys1, iNewSys2;
944 Vec4 pGlui, pSide1, pSide2;
947 bool isGBall = (rndmPtr->flat() < probGluinoballRH);
950 for (
int i = 0; i < systemPtr->size(); ++i) {
951 int iTmp = systemPtr->iParton[i];
952 if (iGlui == 0 && event[ iTmp ].
id() == idRGo) {
954 pGlui =
event[ iTmp ].p();
955 }
else if (iGlui == 0) {
956 iSideTmp.push_back( iTmp);
957 pSide1 +=
event[ iTmp ].p();
959 iSide2.push_back( iTmp);
960 pSide2 +=
event[ iTmp ].p();
965 for (
int i =
int(iSideTmp.size()) - 1; i >= 0; --i)
966 iSide1.push_back( iSideTmp[i]);
967 double m1H = (pGlui + pSide1).mCalc() -
event[ iSide1.back() ].m();
968 double m2H = (pGlui + pSide2).mCalc() -
event[ iSide2.back() ].m();
970 swap( iSide1, iSide2);
971 swap( pSide1, pSide2);
975 for (
int iSide = 1; iSide < 3; ++iSide) {
987 int iBeg =
event.size();
989 iCreRHad[iRHad] = iBeg;
990 event.copy( iGlui, 102);
991 for (
int i = 0; i < int(iSide1.size()); ++i)
992 event.copy( iSide1[i], 102);
994 event.copy( iGlui, 102);
995 for (
int i = 0; i < int(iSide2.size()); ++i)
996 event.copy( iSide2[i], 102);
998 int iEnd =
event.size() - 1;
1001 int idOldL =
event[iEnd].id();
1005 FlavContainer flavOld( idOldL);
1006 idNewQ = -flavSelPtr->pick(flavOld).id;
1007 }
while (iSide == 2 && isDiq1 && abs(idNewQ) > 10);
1008 if (iSide == 1) isDiq1 = (abs(idNewQ) > 10);
1010 bool hasCol = (
event[iEnd].col() > 0);
1017 idRHad = (hasCol) ? 1009002 : -1009002;
1018 if (hasCol) colR =
event[iBeg].col();
1019 else acolR =
event[iBeg].acol();
1021 double mNewQ = particleDataPtr->constituentMass( idNewQ);
1022 if (isGBall) mNewQ *= 0.5;
1023 mRHad =
event[iBeg].m() + mOffsetCloudRH + mNewQ;
1025 idRHad = toIdWithGluino( idSave, idNewQ);
1027 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1028 "cannot form R-hadron code");
1032 mRHad = particleDataPtr->m0(idRHad) +
event[iBeg].m() - m0Go;
1036 double z = zSelPtr->zFrag( idRGo, idNewQ, mRHad*mRHad);
1039 Vec4 pOldH =
event[iBeg].p();
1040 int iOldL = iBeg + 1;
1041 Vec4 pOldL =
event[iOldL].p();
1042 double mOldL =
event[iOldL].m();
1043 double mNewH = mRHad / z;
1044 double sSys = (pOldH + pOldL).m2Calc();
1045 double sRem = (1. - z) * (sSys - mNewH*mNewH);
1046 double sMin = pow2(mOldL + mCollapseRH);
1049 while ( ( sRem < sMin || sSys < pow2(mNewH + mOldL + MSAFETY) )
1052 pOldL +=
event[iOldL].p();
1053 mOldL =
event[iOldL].m();
1054 sSys = (pOldH + pOldL).m2Calc();
1055 sRem = (1. - z) * (sSys - mNewH*mNewH);
1056 sMin = pow2(mOldL + mCollapseRH);
1060 if ( sRem > sMin && sSys > pow2(mNewH + mOldL + MSAFETY) ) {
1062 if ( !newKin( pOldH, pOldL, mNewH, mOldL, pNewH, pNewL) ) {
1063 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1064 "failed to construct kinematics with reduced system");
1069 iRNow =
event.append( idRHad, statusRHad, iBeg, iOldL,
1070 0, 0, colR, acolR, z * pNewH, mRHad, 0.);
1073 if (!isGBall) idNewQ = -idNewQ;
1074 int colN = (hasCol) ? 0 : event[iOldL].acol();
1075 int acolN = (hasCol) ? event[iOldL].col() : 0;
1076 iNewQ =
event.append( idNewQ, 105, iBeg, iOldL, 0, 0,
1077 colN, acolN, (1. - z) * pNewH, (1. - z) * mNewH, 0.);
1078 iNewL =
event.copy( iOldL, 105);
1079 event[iNewL].mothers( iBeg, iOldL);
1080 event[iNewL].p( pNewL);
1084 iNewSys1.push_back( iNewQ);
1085 iNewSys1.push_back( iNewL);
1086 for (
int i = iOldL + 1; i <= iEnd; ++i) iNewSys1.push_back( i);
1088 iNewSys2.push_back( iNewQ);
1089 iNewSys2.push_back( iNewL);
1090 for (
int i = iOldL + 1; i <= iEnd; ++i) iNewSys2.push_back( i);
1098 if (nBody == 0 && isGBall && iSide == 1) {
1099 idQLeap =
event[iOldL].id();
1100 Vec4 pNewH =
event[iBeg].p() + pOldL;
1103 iRNow =
event.append( idRHad, statusRHad, iBeg, iEnd,
1104 0, 0, colR, acolR, pNewH, pNewH.mCalc(), 0.);
1110 if (nBody == 0 && (!isGBall || (iSide == 2 && idQLeap != 0) )) {
1111 if (isGBall) idNewQ = -idQLeap;
1112 FlavContainer flav1( idOldL);
1113 FlavContainer flav2( -idNewQ);
1115 int idNewL = flavSelPtr->combine( flav1, flav2);
1116 while (++iTry < NTRYMAX && idNewL == 0)
1117 idNewL = flavSelPtr->combine( flav1, flav2);
1119 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1120 "cannot form light hadron code");
1123 double mNewL = particleDataPtr->mSel( idNewL);
1126 if ( sSys > pow2(mRHad + mNewL + MSAFETY) ) {
1128 if ( !newKin( pOldH, pOldL, mRHad, mNewL, pRHad, pNewL) ) {
1129 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1130 "failed to construct kinematics for two-hadron decay");
1135 iRNow =
event.append( idRHad, statusRHad, iBeg, iOldL, 0, 0,
1136 colR, acolR, pRHad, mRHad, 0.);
1137 event.append( idNewL, 105, iBeg, iOldL, 0, 0, 0, 0,
1148 if (nBody == 0 && (!isGBall || (iSide == 2 && idQLeap != 0) )) {
1149 if (isGBall) idSave = idQLeap;
1150 if (iSide == 1) idSave = idOldL;
1151 else idRHad = toIdWithGluino( idSave, idOldL);
1153 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1154 "cannot form R-hadron code");
1159 iRNow =
event.append( idRHad, statusRHad, iBeg, iOldL, 0, 0,
1160 colR, acolR, pOldH + pOldL, (pOldH + pOldL).mCalc(), 0.);
1168 int iLast =
event.size() - 1;
1169 for (
int i = iBeg; i <= iOldL; ++i) {
1170 event[i].statusNeg();
1171 event[i].daughters( iRNow, iLast);
1180 infoPtr->errorMsg(
"Error in RHadrons::produceGluino: "
1181 "could not handle gluinoball kinematics");
1184 iRHadron[iRHad] = iGlui;
1185 colConfig.erase(iSys);
1189 if (iNewSys1.size() > 0) colConfig.insert( iNewSys1, event);
1190 if (iNewSys2.size() > 0) colConfig.insert( iNewSys2, event);
1194 }
else if (idQLeap == 0) {
1195 int iG1 = iNewSys1[0];
1196 int iG2 = iNewSys2[0];
1197 int colG =
event[iG1].col() +
event[iG2].col();
1198 int acolG =
event[iG1].acol() +
event[iG2].acol();
1199 Vec4 pG =
event[iG1].p() +
event[iG2].p();
1200 int iG12 =
event.append( 21, 107, iG1, iG2, 0, 0, colG, acolG,
1201 pG, pG.mCalc(), 0.);
1206 event[iG1].statusNeg();
1207 event[iG2].statusNeg();
1208 event[iG1].daughter1( iG12);
1209 event[iG2].daughter1( iG12);
1210 int colBridge =
event.nextColTag();
1211 if (event[iG1].col() == 0) {
1212 event[iG1].col( colBridge);
1213 event[iG2].acol( colBridge);
1215 event[iG1].acol( colBridge);
1216 event[iG2].col( colBridge);
1220 vector<int> iNewSys12;
1221 for (
int i =
int(iNewSys1.size()) - 1; i > 0; --i)
1222 iNewSys12.push_back( iNewSys1[i]);
1223 iNewSys12.push_back( iG12);
1224 for (
int i = 1; i < int(iNewSys2.size()); ++i)
1225 iNewSys12.push_back( iNewSys2[i]);
1226 colConfig.insert( iNewSys12, event);
1231 int iG2 = iNewSys2[0];
1232 event[iG2].id( idQLeap);
1233 colConfig.insert( iNewSys2, event);
1237 event[iGlui].tau( event[iBef].tau() );
1238 if (event[iBef].hasVertex())
event[iGlui].vProd( event[iBef].vProd() );
1250 int RHadrons::toIdWithSquark(
int id1,
int id2) {
1253 int id1Abs = abs(id1);
1254 int id2Abs = abs(id2);
1255 if (id2Abs < 10 && id1 > 0 && id2 > 0)
return 0;
1256 if (id2Abs < 10 && id1 < 0 && id2 < 0)
return 0;
1257 if (id2Abs > 10 && id1 > 0 && id2 < 0)
return 0;
1258 if (id2Abs > 10 && id1 < 0 && id2 > 0)
return 0;
1261 bool isSt = (id1Abs == idRSt);
1262 int idRHad = 1000000;
1263 if (id2Abs < 10) idRHad += ((isSt) ? 600 : 500) + 10 * id2Abs + 2;
1264 else idRHad += ((isSt) ? 6000 : 5000) + 10 * (id2Abs/100) + id2Abs%10;
1265 if (id1 < 0) idRHad = -idRHad;
1277 pair<int,int> RHadrons::fromIdWithSquark(
int idRHad) {
1280 int idLight = (abs(idRHad) - 1000000) / 10;
1281 int idSq = (idLight < 100) ? idLight/10 : idLight/100;
1282 int id1 = (idSq == 6) ? idRSt : idRSb;
1283 if (idRHad < 0) id1 = -id1;
1286 int id2 = (idLight < 100) ? idLight%10 : idLight%100;
1287 if (id2 > 10) id2 = 100 * id2 + abs(idRHad)%10;
1288 if ((id2 < 10 && idRHad > 0) || (id2 > 10 && idRHad < 0)) id2 = -id2;
1291 return make_pair( id1, id2);
1299 int RHadrons::toIdWithGluino(
int id1,
int id2) {
1302 int id1Abs = abs(id1);
1303 int id2Abs = abs(id2);
1304 if (id1Abs == 21 && id2Abs == 21)
return 1000993;
1305 int idMax = max( id1Abs, id2Abs);
1306 int idMin = min( id1Abs, id2Abs);
1307 if (idMin > 10)
return 0;
1308 if (idMax > 10 && id1 > 0 && id2 < 0)
return 0;
1309 if (idMax > 10 && id1 < 0 && id2 > 0)
return 0;
1310 if (idMax < 10 && id1 > 0 && id2 > 0)
return 0;
1311 if (idMax < 10 && id1 < 0 && id2 < 0)
return 0;
1316 idRHad = 1009003 + 100 * idMax + 10 * idMin;
1317 if (idMin != idMax && idMax%2 == 1) {
1318 if (id1Abs == idMax && id1 > 0) idRHad = -idRHad;
1319 if (id2Abs == idMax && id2 > 0) idRHad = -idRHad;
1321 if (idMin != idMax && idMax%2 == 0) {
1322 if (id1Abs == idMax && id1 < 0) idRHad = -idRHad;
1323 if (id2Abs == idMax && id2 < 0) idRHad = -idRHad;
1328 int idA = idMax/1000;
1329 int idB = (idMax/100)%10;
1331 if (idC > idB) swap( idB, idC);
1332 if (idB > idA) swap( idA, idB);
1333 if (idC > idB) swap( idB, idC);
1334 idRHad = 1090004 + 1000 * idA + 100 * idB + 10 * idC;
1335 if (id1 < 0) idRHad = -idRHad;
1348 pair<int,int> RHadrons::fromIdWithGluino(
int idRHad) {
1351 int idLight = (abs(idRHad) - 1000000) / 10;
1352 int id1, id2, idTmp, idA, idB, idC;
1355 if (idLight < 100) {
1356 id1 = (rndmPtr->flat() < 0.5) ? 1 : 2;
1360 }
else if (idLight < 1000) {
1361 id1 = (idLight / 10) % 10;
1362 id2 = -(idLight % 10);
1373 idA = (idLight / 100) % 10;
1374 idB = (idLight / 10) % 10;
1376 double rndmQ = 3. * rndmPtr->flat();
1377 if (idA > 3) rndmQ = 0.5;
1380 id2 = 1000 * idB + 100 * idC + 3;
1381 if (idB != idC && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2;
1382 }
else if (rndmQ < 2.) {
1384 id2 = 1000 * idA + 100 * idC + 3;
1385 if (idA != idC && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2;
1388 id2 = 1000 * idA + 100 * idB +3;
1389 if (idA != idB && rndmPtr->flat() > diquarkSpin1RH) id2 -= 2;
1401 return make_pair( id1, id2);
1411 bool RHadrons::newKin( Vec4 pOld1, Vec4 pOld2,
double mNew1,
double mNew2,
1412 Vec4& pNew1, Vec4& pNew2,
bool checkMargin) {
1415 double sSum = (pOld1 + pOld2).m2Calc();
1416 double sOld1 = pOld1.m2Calc();
1417 double sOld2 = pOld2.m2Calc();
1418 double sNew1 = mNew1 * mNew1;
1419 double sNew2 = mNew2 * mNew2;
1422 if (checkMargin && pow2(mNew1 + mNew2 + MSAFETY) > sSum)
return false;
1425 double lamOld = sqrt( pow2(sSum - sOld1 - sOld2) - 4. * sOld1 * sOld2 );
1426 double lamNew = sqrt( pow2(sSum - sNew1 - sNew2) - 4. * sNew1 * sNew2 );
1427 double move1 = (lamNew * (sSum - sOld1 + sOld2)
1428 - lamOld * (sSum - sNew1 + sNew2)) / (2. * sSum * lamOld);
1429 double move2 = (lamNew * (sSum + sOld1 - sOld2)
1430 - lamOld * (sSum + sNew1 - sNew2)) / (2. * sSum * lamOld);
1433 pNew1 = (1. + move1) * pOld1 - move2 * pOld2;
1434 pNew2 = (1. + move2) * pOld2 - move1 * pOld1;