9 #include "Pythia8/ParticleDecays.h"
23 const int ParticleDecays::NTRYDECAY = 10;
26 const int ParticleDecays::NTRYPICK = 100;
29 const int ParticleDecays::NTRYMEWT = 1000;
32 const int ParticleDecays::NTRYDALITZ = 1000;
35 const double ParticleDecays::MSAFEDALITZ = 1.000001;
39 const double ParticleDecays::WTCORRECTION[11] = { 1., 1., 1.,
40 2., 5., 15., 60., 250., 1250., 7000., 50000. };
46 void ParticleDecays::init(Info* infoPtrIn, Settings& settings,
47 ParticleData* particleDataPtrIn, Rndm* rndmPtrIn,
48 Couplings* couplingsPtrIn, TimeShower* timesDecPtrIn,
49 StringFlav* flavSelPtrIn, DecayHandler* decayHandlePtrIn,
50 vector<int> handledParticles) {
54 particleDataPtr = particleDataPtrIn;
56 couplingsPtr = couplingsPtrIn;
57 flavSelPtr = flavSelPtrIn;
60 timesDecPtr = timesDecPtrIn;
63 decayHandlePtr = decayHandlePtrIn;
66 if (decayHandlePtr != 0)
67 for (
int i = 0; i < int(handledParticles.size()); ++i)
68 particleDataPtr->doExternalDecay(handledParticles[i],
true);
71 mSafety = settings.parm(
"ParticleDecays:mSafety");
74 limitTau0 = settings.flag(
"ParticleDecays:limitTau0");
75 tau0Max = settings.parm(
"ParticleDecays:tau0Max");
76 limitTau = settings.flag(
"ParticleDecays:limitTau");
77 tauMax = settings.parm(
"ParticleDecays:tauMax");
78 limitRadius = settings.flag(
"ParticleDecays:limitRadius");
79 rMax = settings.parm(
"ParticleDecays:rMax");
80 limitCylinder = settings.flag(
"ParticleDecays:limitCylinder");
81 xyMax = settings.parm(
"ParticleDecays:xyMax");
82 zMax = settings.parm(
"ParticleDecays:zMax");
83 limitDecay = limitTau0 || limitTau || limitRadius || limitCylinder;
86 mixB = settings.flag(
"ParticleDecays:mixB");
87 xBdMix = settings.parm(
"ParticleDecays:xBdMix");
88 xBsMix = settings.parm(
"ParticleDecays:xBsMix");
91 sigmaSoft = settings.parm(
"ParticleDecays:sigmaSoft");
94 multIncrease = settings.parm(
"ParticleDecays:multIncrease");
95 multIncreaseWeak = settings.parm(
"ParticleDecays:multIncreaseWeak");
96 multRefMass = settings.parm(
"ParticleDecays:multRefMass");
97 multGoffset = settings.parm(
"ParticleDecays:multGoffset");
98 colRearrange = settings.parm(
"ParticleDecays:colRearrange");
101 stopMass = settings.parm(
"StringFragmentation:stopMass");
104 sRhoDal = pow2(particleDataPtr->m0(113));
105 wRhoDal = pow2(particleDataPtr->mWidth(113));
108 doFSRinDecays = settings.flag(
"ParticleDecays:FSRinDecays");
109 doGammaRad = settings.flag(
"ParticleDecays:allowPhotonRadiation");
112 tauMode = settings.mode(
"TauDecays:mode");
115 if (tauMode) tauDecayer.init(infoPtr, &settings,
116 particleDataPtr, rndmPtr, couplingsPtr);
124 bool ParticleDecays::decay(
int iDec,
Event& event) {
127 Particle& decayer =
event[iDec];
130 if (limitDecay && !checkVertex(decayer))
return true;
133 if (decayer.isResonance()) {
134 infoPtr->errorMsg(
"Warning in ParticleDecays::decay: "
135 "resonance left undecayed");
140 idDec = decayer.id();
144 iProd.push_back( iDec );
145 idProd.push_back( idDec );
146 mProd.push_back( decayer.m() );
149 bool hasOscillated = (abs(idDec) == 511 || abs(idDec) == 531)
150 ? oscillateB(decayer) :
false;
151 if (hasOscillated) {idDec = - idDec; idProd[0] = idDec;}
154 decDataPtr = &decayer.particleDataEntry();
157 bool doneExternally =
false;
158 if (decDataPtr->doExternalDecay()) {
159 motherProd.resize(0);
160 motherProd.push_back( 0 );
162 pProd.push_back(decayer.p());
163 doneExternally = decayHandlePtr->chainDecay( idProd, motherProd,
164 mProd, pProd, iDec, event);
165 if (!doneExternally) doneExternally = decayHandlePtr->decay( idProd,
166 mProd, pProd, iDec, event);
169 if (doneExternally) {
170 int oldSize =
event.size();
171 mult = idProd.size() - 1;
172 while (
int(motherProd.size()) <= mult) motherProd.push_back( 0 );
173 int status = (hasOscillated) ? 94 : 93;
174 for (
int i = 1; i <= mult; ++i) {
175 int iMo = (motherProd[i] == 0) ? iDec : oldSize + motherProd[i] - 1;
176 int iPos =
event.append( idProd[i], status, iMo, 0, 0, 0,
177 0, 0, pProd[i], mProd[i]);
178 iProd.push_back( iPos);
181 event[iMo].statusNeg();
182 if ( event[iMo].daughter1() == 0)
event[iMo].daughter1( iPos);
183 else event[iMo].daughter2( iPos);
189 if (decayer.idAbs() == 15 && !doneExternally && tauMode) {
190 doneExternally = tauDecayer.decay(iDec, event);
191 if (doneExternally)
return true;
195 if (!doneExternally) {
198 if (!decDataPtr->preparePick(idDec, decayer.m()))
return false;
199 bool foundChannel =
false;
200 bool hasStored =
false;
201 for (
int iTryChannel = 0; iTryChannel < NTRYDECAY; ++iTryChannel) {
204 if (hasStored)
event.popBack(mult);
208 DecayChannel& channel = decDataPtr->pickChannel();
209 meMode = channel.meMode();
210 keepPartons = (meMode > 90 && meMode <= 100);
211 mult = channel.multiplicity();
214 bool foundMode =
false;
216 for (
int iTryMode = 0; iTryMode < NTRYDECAY; ++iTryMode) {
223 for (
int i = 0; i < mult; ++i) {
224 int idNow = channel.product(i);
225 int idAbs = abs(idNow);
226 if ( idAbs < 10 || idAbs == 21 || idAbs == 81 || idAbs == 82
227 || idAbs == 83 || (idAbs > 1000 && idAbs < 10000
228 && (idAbs/10)%10 == 0) ) hasPartons =
true;
229 if (idDec < 0 && particleDataPtr->hasAnti(idNow)) idNow = -idNow;
230 double mNow = particleDataPtr->mSel(idNow);
231 idProd.push_back( idNow);
232 mProd.push_back( mNow);
236 if (hasPartons && !keepPartons && !pickHadrons())
continue;
241 for (
int i = 0; i <= mult; ++i) {
245 if (hasPartons && keepPartons && !setColours(event))
continue;
249 double mDiff = mProd[0];
250 for (
int i = 1; i <= mult; ++i) mDiff -= mProd[i];
251 if (mDiff < mSafety)
continue;
258 if (!foundMode)
continue;
261 int status = (hasOscillated) ? 92 : 91;
262 for (
int i = 1; i <= mult; ++i) {
263 int iPos =
event.append( idProd[i], status, iDec, 0, 0, 0,
264 cols[i], acols[i], Vec4(0., 0., 0., 0.), mProd[i], scale);
265 iProd.push_back( iPos);
270 if ( (meMode == 11 || meMode == 12 || meMode == 13)
271 && !dalitzMass() )
continue;
274 bool decayed =
false;
275 if (mult == 1) decayed = oneBody(event);
276 else if (mult == 2) decayed = twoBody(event);
277 else if (mult == 3) decayed = threeBody(event);
278 else decayed = mGenerator(event);
279 if (!decayed)
continue;
282 if (meMode == 11 || meMode == 12 || meMode == 13)
283 dalitzKinematics(event);
292 event[iDec].statusNeg();
293 event[iDec].daughters( iProd[1], iProd[mult]);
297 if (hasStored)
event.popBack(mult);
298 infoPtr->errorMsg(
"Error in ParticleDecays::decay: "
299 "failed to find workable decay channel");
307 if (event[iDec].hasVertex() || event[iDec].tau() > 0.) {
308 Vec4 vDec =
event[iDec].vDec();
309 for (
int i = event[iDec].daughter1(); i <=
event[iDec].daughter2(); ++i)
310 event[i].vProd( vDec );
314 for (
int i = iProd[1]; i <= iProd[mult]; ++i) {
315 event[i].tau( event[i].tau0() * rndmPtr->exp() );
316 if (!event[i].isFinal() && (
event[i].hasVertex() ||
event[i].tau() > 0.)) {
317 Vec4 vDecR =
event[i].vDec();
318 for (
int iR = event[i].daughter1(); iR <=
event[i].daughter2(); ++iR)
319 event[iR].vProd( vDecR );
325 if (hasPartons && keepPartons && doFSRinDecays)
326 timesDecPtr->shower( iProd[1], iProd.back(), event, mProd[0]);
329 else if (doGammaRad && mult == 2 && event[iProd[1]].isLepton()
330 &&
event[iProd[2]].isLepton())
331 timesDecPtr->showerQED( iProd[1], iProd[2], event, mProd[0]);
334 else if (event[iDec].idAbs() > 4900000 && event[iDec].idAbs() < 5000000
335 && doFSRinDecays && mult == 2 && event[iProd[1]].isLepton()) {
336 event[iProd[1]].scale(mProd[0]);
337 event[iProd[2]].scale(mProd[0]);
338 timesDecPtr->shower( iProd[1], iProd.back(), event, mProd[0]);
350 bool ParticleDecays::checkVertex(Particle& decayer) {
353 if (limitTau0 && decayer.tau0() > tau0Max)
return false;
354 if (limitTau && decayer.tau() > tauMax)
return false;
355 if (limitRadius && pow2(decayer.xDec()) + pow2(decayer.yDec())
356 + pow2(decayer.zDec()) > pow2(rMax))
return false;
357 if (limitCylinder && (pow2(decayer.xDec()) + pow2(decayer.yDec())
358 > pow2(xyMax) || abs(decayer.zDec()) > zMax) )
return false;
369 bool ParticleDecays::oscillateB(Particle& decayer) {
372 if (!mixB)
return false;
373 double xBmix = (abs(decayer.id()) == 511) ? xBdMix : xBsMix;
374 double tau = decayer.tau();
375 double tau0 = decayer.tau0();
376 double probosc = pow2(sin(0.5 * xBmix * tau / tau0));
377 return (probosc > rndmPtr->flat());
385 bool ParticleDecays::oneBody(
Event& event) {
388 Particle& decayer =
event[iProd[0]];
389 Particle& prod =
event[iProd[1]];
392 prod.p( decayer.p() );
393 prod.m( decayer.m() );
394 prod.mother2( iProd[0] );
405 bool ParticleDecays::twoBody(
Event& event) {
408 Particle& decayer =
event[iProd[0]];
409 Particle& prod1 =
event[iProd[1]];
410 Particle& prod2 =
event[iProd[2]];
413 double m0 = mProd[0];
414 double m1 = mProd[1];
415 double m2 = mProd[2];
418 if (m1 + m2 + mSafety > m0)
return false;
419 double e1 = 0.5 * (m0*m0 + m1*m1 - m2*m2) / m0;
420 double e2 = 0.5 * (m0*m0 + m2*m2 - m1*m1) / m0;
421 double pAbs = 0.5 * sqrtpos( (m0 - m1 - m2) * (m0 + m1 + m2)
422 * (m0 + m1 - m2) * (m0 - m1 + m2) ) / m0;
426 int iMother =
event[iProd[0]].mother1();
429 if (iMother <= 0 || iMother >= iProd[0]) meMode = 0;
431 int iDaughter1 =
event[iMother].daughter1();
432 int iDaughter2 =
event[iMother].daughter2();
433 if (iDaughter2 != iDaughter1 + 1) meMode = 0;
435 int idMother = abs( event[iMother].
id() );
436 if (idMother <= 100 || idMother%10 !=1
437 || (idMother/1000)%10 != 0) meMode = 0;
439 int iSister = (iProd[0] == iDaughter1) ? iDaughter2 : iDaughter1;
440 idSister = abs( event[iSister].
id() );
441 if ( (idSister <= 100 || idSister%10 !=1
442 || (idSister/1000)%10 != 0) && idSister != 22) meMode = 0;
449 double wtME, wtMEmax;
457 double cosTheta = 2. * rndmPtr->flat() - 1.;
458 double sinTheta = sqrt(1. - cosTheta*cosTheta);
459 double phi = 2. * M_PI * rndmPtr->flat();
460 double pX = pAbs * sinTheta * cos(phi);
461 double pY = pAbs * sinTheta * sin(phi);
462 double pZ = pAbs * cosTheta;
465 prod1.p( pX, pY, pZ, e1);
466 prod2.p( -pX, -pY, -pZ, e2);
467 prod1.bst( decayer.p(), decayer.m() );
468 prod2.bst( decayer.p(), decayer.m() );
474 double p10 = decayer.p() *
event[iMother].p();
475 double p12 = decayer.p() * prod1.p();
476 double p02 =
event[iMother].p() * prod1.p();
477 double s0 = pow2(event[iMother].m());
478 double s1 = pow2(decayer.m());
479 double s2 = pow2(prod1.m());
480 if (idSister != 22) wtME = pow2(p10 * p12 - s1 * p02);
481 else wtME = s1 * (2. * p10 * p12 * p02 - s1 * p02*p02
482 - s0 * p12*p12 - s2 * p10*p10 + s1 * s0 * s2);
483 wtME = max( wtME, 1e-6 * s1*s1 * s0 * s2);
484 wtMEmax = (p10*p10 - s1 * s0) * (p12*p12 - s1 * s2);
488 if(loop > NTRYMEWT) {
489 infoPtr->errorMsg(
"ParticleDecays::twoBody: "
490 "caught in infinite ME weight loop");
495 }
while ( wtME < rndmPtr->flat() * wtMEmax );
506 bool ParticleDecays::threeBody(
Event& event) {
509 Particle& decayer =
event[iProd[0]];
510 Particle& prod1 =
event[iProd[1]];
511 Particle& prod2 =
event[iProd[2]];
512 Particle& prod3 =
event[iProd[3]];
515 double m0 = mProd[0];
516 double m1 = mProd[1];
517 double m2 = mProd[2];
518 double m3 = mProd[3];
519 double mSum = m1 + m2 + m3;
520 double mDiff = m0 - mSum;
521 if (mDiff < mSafety)
return false;
524 double m23Min = m2 + m3;
525 double m23Max = m0 - m1;
526 double p1Max = 0.5 * sqrtpos( (m0 - m1 - m23Min) * (m0 + m1 + m23Min)
527 * (m0 + m1 - m23Min) * (m0 - m1 + m23Min) ) / m0;
528 double p23Max = 0.5 * sqrtpos( (m23Max - m2 - m3) * (m23Max + m2 + m3)
529 * (m23Max + m2 - m3) * (m23Max - m2 + m3) ) / m23Max;
530 double wtPSmax = 0.5 * p1Max * p23Max;
533 double wtME, wtMEmax, wtPS, m23, p1Abs, p23Abs;
540 m23 = m23Min + rndmPtr->flat() * mDiff;
543 p1Abs = 0.5 * sqrtpos( (m0 - m1 - m23) * (m0 + m1 + m23)
544 * (m0 + m1 - m23) * (m0 - m1 + m23) ) / m0;
545 p23Abs = 0.5 * sqrtpos( (m23 - m2 - m3) * (m23 + m2 + m3)
546 * (m23 + m2 - m3) * (m23 - m2 + m3) ) / m23;
547 wtPS = p1Abs * p23Abs;
550 }
while ( wtPS < rndmPtr->flat() * wtPSmax );
553 double cosTheta = 2. * rndmPtr->flat() - 1.;
554 double sinTheta = sqrt(1. - cosTheta*cosTheta);
555 double phi = 2. * M_PI * rndmPtr->flat();
556 double pX = p23Abs * sinTheta * cos(phi);
557 double pY = p23Abs * sinTheta * sin(phi);
558 double pZ = p23Abs * cosTheta;
559 double e2 = sqrt( m2*m2 + p23Abs*p23Abs);
560 double e3 = sqrt( m3*m3 + p23Abs*p23Abs);
561 prod2.p( pX, pY, pZ, e2);
562 prod3.p( -pX, -pY, -pZ, e3);
565 cosTheta = 2. * rndmPtr->flat() - 1.;
566 sinTheta = sqrt(1. - cosTheta*cosTheta);
567 phi = 2. * M_PI * rndmPtr->flat();
568 pX = p1Abs * sinTheta * cos(phi);
569 pY = p1Abs * sinTheta * sin(phi);
570 pZ = p1Abs * cosTheta;
571 double e1 = sqrt( m1*m1 + p1Abs*p1Abs);
572 double e23 = sqrt( m23*m23 + p1Abs*p1Abs);
573 prod1.p( pX, pY, pZ, e1);
576 Vec4 p23( -pX, -pY, -pZ, e23);
577 prod2.bst( p23, m23 );
578 prod3.bst( p23, m23 );
582 double p1p2 = prod1.p() * prod2.p();
583 double p1p3 = prod1.p() * prod3.p();
584 double p2p3 = prod2.p() * prod3.p();
585 wtME = pow2(m1 * m2 * m3) - pow2(m1 * p2p3) - pow2(m2 * p1p3)
586 - pow2(m3 * p1p2) + 2. * p1p2 * p1p3 * p2p3;
587 wtMEmax = pow3(m0 * m0) / 150.;
590 }
else if (meMode == 21) {
591 double x1 = 2. * prod1.e() / m0;
592 wtME = x1 * (3. - 2. * x1);
593 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
594 wtMEmax = xMax * (3. - 2. * xMax);
597 }
else if ((meMode == 22 || meMode == 23) && prod1.isLepton()) {
598 wtME = m0 * prod1.e() * (prod2.p() * prod3.p());
599 wtMEmax = min( pow4(m0) / 16., m0 * (m0 - m1 - m2) * (m0 - m1 - m3)
603 }
else if (meMode == 22 || meMode == 23) {
604 double x1 = 2. * prod1.pAbs() / m0;
605 wtME = x1 * (3. - 2. * x1);
606 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
607 wtMEmax = xMax * (3. - 2. * xMax);
610 }
else if (meMode == 31) {
611 double x1 = 2. * prod1.e() / m0;
613 double x1Max = 1. - pow2(mSum / m0);
614 wtMEmax = pow3(x1Max);
617 }
else if (meMode == 92) {
618 double x1 = 2. * prod1.e() / m0;
619 double x2 = 2. * prod2.e() / m0;
620 double x3 = 2. * prod3.e() / m0;
621 wtME = pow2( (1. - x1) / (x2 * x3) ) + pow2( (1. - x2) / (x1 * x3) )
622 + pow2( (1. - x3) / (x1 * x2) );
625 if (prod1.id() == 22 && sqrt(1. - x1) * m0 < 2. * stopMass) wtME = 0.;
626 if (prod2.id() == 22 && sqrt(1. - x2) * m0 < 2. * stopMass) wtME = 0.;
627 if (prod3.id() == 22 && sqrt(1. - x3) * m0 < 2. * stopMass) wtME = 0.;
631 }
while ( wtME < rndmPtr->flat() * wtMEmax );
634 prod1.bst( decayer.p(), decayer.m() );
635 prod2.bst( decayer.p(), decayer.m() );
636 prod3.bst( decayer.p(), decayer.m() );
647 bool ParticleDecays::mGenerator(
Event& event) {
650 double m0 = mProd[0];
651 double mSum = mProd[1];
652 for (
int i = 2; i <= mult; ++i) mSum += mProd[i];
653 double mDiff = m0 - mSum;
654 if (mDiff < mSafety)
return false;
658 for (
int i = 0; i <= mult; ++i) mInv.push_back( mProd[i]);
661 double wtPS, wtME, wtMEmax;
662 double wtPSmax = 1. / WTCORRECTION[mult];
663 double mMax = mDiff + mProd[mult];
665 for (
int i = mult - 1; i > 0; --i) {
668 double mNow = mProd[i];
669 wtPSmax *= 0.5 * sqrtpos( (mMax - mMin - mNow) * (mMax + mMin + mNow)
670 * (mMax + mMin - mNow) * (mMax - mMin + mNow) ) / mMax;
684 rndmOrd.push_back(1.);
685 for (
int i = 1; i < mult - 1; ++i) {
686 double rndm = rndmPtr->flat();
687 rndmOrd.push_back(rndm);
688 for (
int j = i - 1; j > 0; --j) {
689 if (rndm > rndmOrd[j]) swap( rndmOrd[j], rndmOrd[j+1] );
693 rndmOrd.push_back(0.);
696 for (
int i = mult - 1; i > 0; --i) {
697 mInv[i] = mInv[i+1] + mProd[i] + (rndmOrd[i-1] - rndmOrd[i]) * mDiff;
698 wtPS *= 0.5 * sqrtpos( (mInv[i] - mInv[i+1] - mProd[i])
699 * (mInv[i] + mInv[i+1] + mProd[i]) * (mInv[i] + mInv[i+1] - mProd[i])
700 * (mInv[i] - mInv[i+1] + mProd[i]) ) / mInv[i];
704 }
while ( wtPS < rndmPtr->flat() * wtPSmax );
707 pInv.resize(mult + 1);
708 for (
int i = 1; i < mult; ++i) {
709 double pAbs = 0.5 * sqrtpos( (mInv[i] - mInv[i+1] - mProd[i])
710 * (mInv[i] + mInv[i+1] + mProd[i]) * (mInv[i] + mInv[i+1] - mProd[i])
711 * (mInv[i] - mInv[i+1] + mProd[i]) ) / mInv[i];
714 double cosTheta = 2. * rndmPtr->flat() - 1.;
715 double sinTheta = sqrt(1. - cosTheta*cosTheta);
716 double phi = 2. * M_PI * rndmPtr->flat();
717 double pX = pAbs * sinTheta * cos(phi);
718 double pY = pAbs * sinTheta * sin(phi);
719 double pZ = pAbs * cosTheta;
722 double eHad = sqrt( mProd[i]*mProd[i] + pAbs*pAbs);
723 double eInv = sqrt( mInv[i+1]*mInv[i+1] + pAbs*pAbs);
724 event[iProd[i]].p( pX, pY, pZ, eHad);
725 pInv[i+1].p( -pX, -pY, -pZ, eInv);
729 event[iProd[mult]].p( pInv[mult] );
730 for (
int iFrame = mult - 1; iFrame > 1; --iFrame)
731 for (
int i = iFrame; i <= mult; ++i)
732 event[iProd[i]].bst( pInv[iFrame], mInv[iFrame]);
735 if (meMode == 21 && event[iProd[1]].isLepton()) {
736 double x1 = 2. *
event[iProd[1]].e() / m0;
737 wtME = x1 * (3. - 2. * x1);
738 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
739 wtMEmax = xMax * (3. - 2. * xMax);
743 }
else if ((meMode == 22 || meMode == 23) &&
event[iProd[1]].isLepton()) {
744 Vec4 pRest =
event[iProd[3]].p();
745 for (
int i = 4; i <= mult; ++i) pRest += event[iProd[i]].p();
746 wtME = m0 *
event[iProd[1]].e() * (
event[iProd[2]].p() * pRest);
747 for (
int i = 4; i <= mult; ++i) wtME
748 *= exp(- event[iProd[i]].pAbs2() / pow2(sigmaSoft) );
749 wtMEmax = pow4(m0) / 16.;
752 }
else if (meMode == 22 || meMode == 23) {
753 double x1 = 2. *
event[iProd[1]].pAbs() / m0;
754 wtME = x1 * (3. - 2. * x1);
755 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
756 wtMEmax = xMax * (3. - 2. * xMax);
759 }
else if (meMode == 31) {
760 double x1 = 2. *
event[iProd[1]].e() / m0;
762 double x1Max = 1. - pow2(mSum / m0);
763 wtMEmax = pow3(x1Max);
767 }
while ( wtME < rndmPtr->flat() * wtMEmax );
770 pInv[1].p( event[iProd[0]].p() );
771 for (
int i = 1; i <= mult; ++i) event[iProd[i]].bst( pInv[1], mInv[1] );
782 bool ParticleDecays::dalitzMass() {
786 for (
int i = 1; i <= mult - 2; ++i) mSum1 += mProd[i];
787 if (meMode == 13) mSum1 *= MSAFEDALITZ;
788 double mSum2 = MSAFEDALITZ * (mProd[mult -1] + mProd[mult]);
789 double mDiff = mProd[0] - mSum1 - mSum2;
792 if (mDiff < mSafety)
return false;
793 if (idProd[mult - 1] + idProd[mult] != 0
794 || mProd[mult - 1] != mProd[mult]) {
795 infoPtr->errorMsg(
"Error in ParticleDecays::dalitzMass:"
796 " inconsistent flavour/mass assignments");
799 if ( meMode == 13 && (idProd[1] + idProd[2] != 0
800 || mProd[1] != mProd[2]) ) {
801 infoPtr->errorMsg(
"Error in ParticleDecays::dalitzMass:"
802 " inconsistent flavour/mass assignments");
807 if (meMode == 11 || meMode == 12) {
810 double sGamMin = pow2(mSum2);
811 double sGamMax = pow2(mProd[0] - mSum1);
816 if (++loop > NTRYDALITZ)
return false;
817 sGam = sGamMin * pow( sGamMax / sGamMin, rndmPtr->flat() );
818 wtGam = (1. + 0.5 * sGamMin / sGam) * sqrt(1. - sGamMin / sGam)
819 * pow3(1. - sGam / sGamMax) * sRhoDal * (sRhoDal + wRhoDal)
820 / ( pow2(sGam - sRhoDal) + sRhoDal * wRhoDal );
821 }
while ( wtGam < rndmPtr->flat() );
825 mProd[mult] = sqrt(sGam);
831 double s0 = pow2(mProd[0]);
832 double s12Min = pow2(mSum1);
833 double s12Max = pow2(mProd[0] - mSum2);
834 double s34Min = pow2(mSum2);
835 double s34Max = pow2(mProd[0] - mSum1);
838 double s12, s34, wt12, wt34, wtPAbs, wtAll;
841 if (++loop > NTRYDALITZ)
return false;
842 s12 = s12Min * pow( s12Max / s12Min, rndmPtr->flat() );
843 wt12 = (1. + 0.5 * s12Min / s12) * sqrt(1. - s12Min / s12)
844 * sRhoDal * (sRhoDal + wRhoDal)
845 / ( pow2(s12 - sRhoDal) + sRhoDal * wRhoDal );
846 s34 = s34Min * pow( s34Max / s34Min, rndmPtr->flat() );
847 wt34 = (1. + 0.5 * s34Min / s34) * sqrt(1. - s34Min / s34)
848 * sRhoDal * (sRhoDal + wRhoDal)
849 / ( pow2(s34 - sRhoDal) + sRhoDal * wRhoDal );
850 wtPAbs = sqrtpos( pow2(1. - (s12 + s34)/ s0)
851 - 4. * s12 * s34 / (s0 * s0) );
852 wtAll = wt12 * wt34 * pow3(wtPAbs);
853 if (wtAll > 1.) infoPtr->errorMsg(
854 "Error in ParticleDecays::dalitzMass: weight > 1");
855 }
while (wtAll < rndmPtr->flat());
859 mProd[1] = sqrt(s12);
860 mProd[2] = sqrt(s34);
872 bool ParticleDecays::dalitzKinematics(
Event& event) {
875 int nDal = (meMode < 13) ? 1 : 2;
879 for (
int iDal = 0; iDal < nDal; ++iDal) {
882 Particle& decayer =
event[iProd[0]];
883 Particle& prodA = (iDal == 0) ? event[iProd[mult - 1]]
885 Particle& prodB = (iDal == 0) ? event[iProd[mult]]
889 Vec4 pDec = decayer.p();
890 int iGam = (meMode < 13) ? mult - 1 : 2 - iDal;
891 Vec4 pGam =
event[iProd[iGam]].p();
892 pGam.bstback( pDec, decayer.m() );
893 double phiGam = pGam.phi();
894 pGam.rot( 0., -phiGam);
895 double thetaGam = pGam.theta();
896 pGam.rot( -thetaGam, 0.);
899 double mGam = (meMode < 13) ? mProd[mult - 1] : mProd[2 - iDal];
900 double mA = prodA.m();
901 double mB = prodB.m();
902 double mGamMin = MSAFEDALITZ * (mA + mB);
903 double mGamRat = pow2(mGamMin / mGam);
904 double pGamAbs = 0.5 * sqrtpos( (mGam - mA - mB) * (mGam + mA + mB) );
907 double cosTheta, cos2Theta;
909 cosTheta = 2. * rndmPtr->flat() - 1.;
910 cos2Theta = cosTheta * cosTheta;
911 }
while ( 1. + cos2Theta + mGamRat * (1. - cos2Theta)
912 < 2. * rndmPtr->flat() );
913 double sinTheta = sqrt(1. - cosTheta*cosTheta);
914 double phi = 2. * M_PI * rndmPtr->flat();
915 double pX = pGamAbs * sinTheta * cos(phi);
916 double pY = pGamAbs * sinTheta * sin(phi);
917 double pZ = pGamAbs * cosTheta;
918 double eA = sqrt( mA*mA + pGamAbs*pGamAbs);
919 double eB = sqrt( mB*mB + pGamAbs*pGamAbs);
920 prodA.p( pX, pY, pZ, eA);
921 prodB.p( -pX, -pY, -pZ, eB);
924 prodA.bst( pGam, mGam);
925 prodB.bst( pGam, mGam);
926 prodA.rot( thetaGam, phiGam);
927 prodB.rot( thetaGam, phiGam);
928 prodA.bst( pDec, decayer.m() );
929 prodB.bst( pDec, decayer.m() );
941 bool ParticleDecays::pickHadrons() {
948 bool closedGLoop =
false;
949 for (
int i = 1; i <= mult; ++i) {
950 int idAbs = abs(idProd[i]);
951 if ( idAbs < 9 || (idAbs > 1000 && idAbs < 10000 && (idAbs/10)%10 == 0)
952 || idAbs == 81 || idAbs == 82 || idAbs == 83) {
954 idPartons.push_back(idProd[i]);
955 if (idAbs == 83) closedGLoop =
true;
959 idProd[nKnown] = idProd[i];
960 mProd[nKnown] = mProd[i];
966 for (
int i = 0; i < nPartons; ++i) {
967 int idPart = idPartons[i];
970 int idAbs = abs(idDec);
971 if ( (idAbs/1000)%10 == 0 ) {
972 idNew = -(idAbs/10)%10;
973 if ((idAbs/100)%2 == 1) idNew = -idNew;
974 }
else if ( (idAbs/100)%10 >= (idAbs/10)%10 )
975 idNew = 100 * ((idAbs/10)%100) + 3;
976 else idNew = 1000 * ((idAbs/10)%10) + 100 * ((idAbs/100)%10) + 1;
977 if (idDec < 0) idNew = -idNew;
980 }
else if (idPart == 82 || idPart == 83) {
983 int idDummy = -flavSelPtr->pickLightQ();
984 FlavContainer flavDummy(idDummy, idPart - 82);
985 do idNew = flavSelPtr->pick(flavDummy).id;
987 mFlav = particleDataPtr->constituentMass(idNew);
988 }
while (2. * mFlav + stopMass > mProd[0]);
989 }
else if (idPart == -82 || idPart == -83) {
990 idNew = -idPartons[i-1];
992 idPartons[i] = idNew;
996 int nMin = max( 2, nKnown + nPartons / 2);
997 if (meMode == 23) nMin = 3;
998 if (meMode > 41 && meMode <= 50) nMin = meMode - 40;
999 if (meMode > 51 && meMode <= 60) nMin = meMode - 50;
1001 if (meMode == 0) nFix = nMin;
1002 if (meMode == 11) nFix = 3;
1003 if (meMode == 12) nFix = 4;
1004 if (meMode > 61 && meMode <= 70) nFix = meMode - 60;
1005 if (meMode > 71 && meMode <= 80) nFix = meMode - 70;
1006 if (nFix > 0 && nKnown + nPartons/2 > nFix)
return false;
1009 int nFilled, nTotal, nNew, nSpec, nLeft;
1012 bool diquarkClash =
false;
1013 bool usedChannel =
false;
1018 if (nTry > NTRYPICK)
return false;
1021 nFilled = nKnown + 1;
1022 idProd.resize(nFilled);
1023 mProd.resize(nFilled);
1028 for (
int i = 1; i < nFilled; ++i) mDiff -= mProd[i];
1029 diquarkClash =
false;
1030 usedChannel =
false;
1033 if ( (meMode == 22 || meMode == 23) && nLeft > 1) {
1034 FlavContainer flav1( idPartons[nPartons - 2] );
1035 FlavContainer flav2( idPartons[nPartons - 1] );
1037 do idHad = flavSelPtr->combine( flav1, flav2);
1039 double mHad = particleDataPtr->mSel(idHad);
1041 idProd.push_back( idHad);
1042 mProd.push_back( mHad);
1053 double geom = rndmPtr->flat();
1058 }
while (geom < 1. && nTotal < 10);
1061 }
else if (nFix == 0) {
1062 double multIncreaseNow = (meMode == 23)
1063 ? multIncreaseWeak : multIncrease;
1064 double mDiffPS = mDiff;
1065 for (
int i = 0; i < nLeft; ++i)
1066 mDiffPS -= particleDataPtr->constituentMass( idPartons[i] );
1067 double average = 0.5 * (nKnown + nSpec) + 0.25 * nPartons
1068 + multIncreaseNow * log( max( 1.1, mDiffPS / multRefMass ) );
1069 if (closedGLoop) average += multGoffset;
1074 for (
int nNow = nMin + 1; nNow <= 10; ++nNow) {
1075 value *= average / nNow;
1080 sum *= rndmPtr->flat();
1084 value *= average / nTotal;
1086 }
while (sum > 0. && nTotal < 10);
1092 nNew = nTotal - nKnown - nSpec;
1096 for (
int i = 0; i < nLeft; ++i) {
1097 flavEnds.push_back( FlavContainer(idPartons[i]) );
1098 if (abs(idPartons[i]) > 100) flavSelPtr->assignPopQ( flavEnds[i] );
1102 if (nNew > nLeft/2) {
1103 FlavContainer flavNew;
1105 for (
int i = 0; i < nNew - nLeft/2; ++i) {
1107 int iEnd = int( (nLeft - 1.) * rndmPtr->flat() );
1110 flavNew = flavSelPtr->pick( flavEnds[iEnd] );
1111 idHad = flavSelPtr->combine( flavEnds[iEnd], flavNew);
1112 }
while (idHad == 0);
1114 idProd.push_back( idHad);
1115 flavEnds[iEnd].anti(flavNew);
1122 if ( abs(flavEnds[0].
id) > 8 && abs(flavEnds[1].
id) > 8)
1123 diquarkClash =
true;
1125 do idHad = flavSelPtr->combine( flavEnds[0], flavEnds[1]);
1127 idProd.push_back( idHad);
1136 if ( rndmPtr->flat() < colRearrange) iEnd2 = 3;
1138 ( (flavEnds[iEnd1].id > 0 && flavEnds[iEnd1].id < 9)
1139 || flavEnds[iEnd1].
id < -10 ) ? 1 : -1;
1140 if ( (flavEnds[iEnd2].
id < 0 && flavEnds[iEnd2].
id > -9)
1141 || flavEnds[iEnd2].id > 10 ) relColSign *= -1;
1142 if (relColSign == 1) iEnd2 = 2;
1143 if (iEnd2 == 2) iEnd3 = 1;
1144 if (iEnd2 == 3) iEnd4 = 1;
1148 if ( abs(flavEnds[iEnd1].
id) > 8 && abs(flavEnds[iEnd2].
id) > 8)
1149 diquarkClash =
true;
1151 do idHad = flavSelPtr->combine( flavEnds[iEnd1], flavEnds[iEnd2]);
1153 idProd.push_back( idHad);
1155 if ( abs(flavEnds[iEnd3].
id) > 8 && abs(flavEnds[iEnd4].
id) > 8)
1156 diquarkClash =
true;
1158 do idHad = flavSelPtr->combine( flavEnds[iEnd3], flavEnds[iEnd4]);
1160 idProd.push_back( idHad);
1165 for (
int i = nFilled; i < int(idProd.size()) ; ++i) {
1166 double mHad = particleDataPtr->mSel(idProd[i]);
1167 mProd.push_back( mHad);
1173 if ( ( (meMode > 51 && meMode <= 60) || (meMode > 71 && meMode <= 80) )
1174 && mDiff > mSafety && !diquarkClash ) {
1175 int idMatch[10], idPNow;
1176 usedChannel =
false;
1177 bool matched =
false;
1179 for (
int i = 0; i < decDataPtr->sizeChannels(); ++i) {
1180 DecayChannel& channel = decDataPtr->channel(i);
1181 if (channel.multiplicity() != nTotal)
continue;
1182 for (
int k = 0; k < nTotal; ++k) idMatch[k] = channel.product(k);
1185 for (
int j = 0; j < nTotal; ++j) {
1187 idPNow = idProd[j + 1];
1188 if (idPNow == -311 || idPNow == 130 || idPNow == 310) idPNow = 311;
1189 for (
int k = 0; k < nTotal - j; ++k)
1190 if (idMatch[k] == idPNow || (idMatch[k] == -311 && idPNow == 311)) {
1192 idMatch[k] = idMatch[nTotal - 1 - j];
1196 if (!matched)
break;
1199 if (matched) {usedChannel =
true;
break;}
1204 }
while (mDiff < mSafety || diquarkClash || usedChannel);
1207 mult = idProd.size() - 1;
1210 if (meMode == 11 || meMode == 12) {
1213 for (
int i = 1; i <= mult; ++i) {
1214 if (idProd[i] == 11 || idProd[i] == 13 || idProd[i] == 15) iL1 = i;
1215 if (idProd[i] == -11 || idProd[i] == -13 || idProd[i] == -15) iL2 = i;
1217 if (iL1 > 0 && iL2 > 0) {
1218 int idL1 = idProd[iL1];
1219 int idL2 = idProd[iL2];
1220 double mL1 = mProd[iL1];
1221 double mL2 = mProd[iL2];
1223 for (
int i = 1; i <= mult; ++i)
if (i != iL1 && i != iL2) {
1225 idProd[iMove] = idProd[i];
1226 mProd[iMove] = mProd[i];
1228 idProd[mult - 1] = idL1;
1229 idProd[mult] = idL2;
1230 mProd[mult - 1] = mL1;
1244 bool ParticleDecays::setColours(
Event& event) {
1247 if (meMode == 91 && idProd[1] > 0 && idProd[1] < 9) {
1248 int newCol =
event.nextColTag();
1251 }
else if (meMode == 91 && idProd[1] < 0 && idProd[1] > -9) {
1252 int newCol =
event.nextColTag();
1257 }
else if (meMode == 91 && idProd[1] == 21) {
1258 int newCol1 =
event.nextColTag();
1259 int newCol2 =
event.nextColTag();
1266 }
else if (meMode == 92 && idProd[1] == 21 && idProd[2] == 21
1267 && idProd[3] == 21) {
1268 int newCol1 =
event.nextColTag();
1269 int newCol2 =
event.nextColTag();
1270 int newCol3 =
event.nextColTag();
1279 }
else if (meMode == 92) {
1280 int iGlu1 = (idProd[1] == 21) ? 1 : 3;
1281 int iGlu2 = (idProd[2] == 21) ? 2 : 3;
1282 int newCol1 =
event.nextColTag();
1283 int newCol2 =
event.nextColTag();
1284 cols[iGlu1] = newCol1;
1285 acols[iGlu1] = newCol2;
1286 cols[iGlu2] = newCol2;
1287 acols[iGlu2] = newCol1;
1290 }
else return false;