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 sophisticatedTau = settings.mode(
"ParticleDecays:sophisticatedTau");
115 if (sophisticatedTau) 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()) {
160 pProd.push_back(decayer.p());
161 doneExternally = decayHandlePtr->decay(idProd, mProd, pProd,
165 if (doneExternally) {
166 mult = idProd.size() - 1;
167 int status = (hasOscillated) ? 94 : 93;
168 for (
int i = 1; i <= mult; ++i) {
169 int iPos =
event.append( idProd[i], status, iDec, 0, 0, 0,
170 0, 0, pProd[i], mProd[i]);
171 iProd.push_back( iPos);
175 event[iDec].statusNeg();
176 event[iDec].daughters( iProd[1], iProd[mult]);
181 if (decayer.idAbs() == 15 && !doneExternally && sophisticatedTau) {
182 doneExternally = tauDecayer.decay(iDec, event);
183 if (doneExternally)
return true;
187 if (!doneExternally) {
190 if (!decDataPtr->preparePick(idDec, decayer.m()))
return false;
191 bool foundChannel =
false;
192 bool hasStored =
false;
193 for (
int iTryChannel = 0; iTryChannel < NTRYDECAY; ++iTryChannel) {
196 if (hasStored)
event.popBack(mult);
200 DecayChannel& channel = decDataPtr->pickChannel();
201 meMode = channel.meMode();
202 keepPartons = (meMode > 90 && meMode <= 100);
203 mult = channel.multiplicity();
206 bool foundMode =
false;
207 for (
int iTryMode = 0; iTryMode < NTRYDECAY; ++iTryMode) {
214 for (
int i = 0; i < mult; ++i) {
215 int idNow = channel.product(i);
216 int idAbs = abs(idNow);
217 if ( idAbs < 10 || idAbs == 21 || idAbs == 81 || idAbs == 82
218 || idAbs == 83 || (idAbs > 1000 && idAbs < 10000
219 && (idAbs/10)%10 == 0) ) hasPartons =
true;
220 if (idDec < 0 && particleDataPtr->hasAnti(idNow)) idNow = -idNow;
221 double mNow = particleDataPtr->mSel(idNow);
222 idProd.push_back( idNow);
223 mProd.push_back( mNow);
227 if (hasPartons && !keepPartons && !pickHadrons())
continue;
232 for (
int i = 0; i <= mult; ++i) {
236 if (hasPartons && keepPartons && !setColours(event))
continue;
240 double mDiff = mProd[0];
241 for (
int i = 1; i <= mult; ++i) mDiff -= mProd[i];
242 if (mDiff < mSafety)
continue;
249 if (!foundMode)
continue;
252 int status = (hasOscillated) ? 92 : 91;
253 for (
int i = 1; i <= mult; ++i) {
254 int iPos =
event.append( idProd[i], status, iDec, 0, 0, 0,
255 cols[i], acols[i], Vec4(0., 0., 0., 0.), mProd[i], scale);
256 iProd.push_back( iPos);
261 if ( (meMode == 11 || meMode == 12 || meMode == 13)
262 && !dalitzMass() )
continue;
265 bool decayed =
false;
266 if (mult == 1) decayed = oneBody(event);
267 else if (mult == 2) decayed = twoBody(event);
268 else if (mult == 3) decayed = threeBody(event);
269 else decayed = mGenerator(event);
270 if (!decayed)
continue;
273 if (meMode == 11 || meMode == 12 || meMode == 13)
274 dalitzKinematics(event);
283 event[iDec].statusNeg();
284 event[iDec].daughters( iProd[1], iProd[mult]);
288 if (hasStored)
event.popBack(mult);
289 infoPtr->errorMsg(
"Error in ParticleDecays::decay: "
290 "failed to find workable decay channel");
298 if (event[iDec].hasVertex() || event[iDec].tau() > 0.) {
299 Vec4 vDec =
event[iDec].vDec();
300 for (
int i = 1; i <= mult; ++i) event[iProd[i]].vProd( vDec );
304 for (
int i = 1; i <= mult; ++i)
305 event[iProd[i]].tau( event[iProd[i]].tau0() * rndmPtr->exp() );
309 if (hasPartons && keepPartons && doFSRinDecays)
310 timesDecPtr->shower( iProd[1], iProd.back(), event, mProd[0]);
313 else if (doGammaRad && mult == 2 && event[iProd[1]].isLepton()
314 &&
event[iProd[2]].isLepton())
315 timesDecPtr->showerQED( iProd[1], iProd[2], event, mProd[0]);
318 else if (event[iDec].idAbs() > 4900000 && event[iDec].idAbs() < 5000000
319 && doFSRinDecays && mult == 2 && event[iProd[1]].isLepton()) {
320 event[iProd[1]].scale(mProd[0]);
321 event[iProd[2]].scale(mProd[0]);
322 timesDecPtr->shower( iProd[1], iProd.back(), event, mProd[0]);
334 bool ParticleDecays::checkVertex(Particle& decayer) {
337 if (limitTau0 && decayer.tau0() > tau0Max)
return false;
338 if (limitTau && decayer.tau() > tauMax)
return false;
339 if (limitRadius && pow2(decayer.xDec()) + pow2(decayer.yDec())
340 + pow2(decayer.zDec()) > pow2(rMax))
return false;
341 if (limitCylinder && (pow2(decayer.xDec()) + pow2(decayer.yDec())
342 > pow2(xyMax) || abs(decayer.zDec()) > zMax) )
return false;
353 bool ParticleDecays::oscillateB(Particle& decayer) {
356 if (!mixB)
return false;
357 double xBmix = (abs(decayer.id()) == 511) ? xBdMix : xBsMix;
358 double tau = decayer.tau();
359 double tau0 = decayer.tau0();
360 double probosc = pow2(sin(0.5 * xBmix * tau / tau0));
361 return (probosc > rndmPtr->flat());
369 bool ParticleDecays::oneBody(
Event& event) {
372 Particle& decayer =
event[iProd[0]];
373 Particle& prod =
event[iProd[1]];
376 prod.p( decayer.p() );
377 prod.m( decayer.m() );
378 prod.mother2( iProd[0] );
389 bool ParticleDecays::twoBody(
Event& event) {
392 Particle& decayer =
event[iProd[0]];
393 Particle& prod1 =
event[iProd[1]];
394 Particle& prod2 =
event[iProd[2]];
397 double m0 = mProd[0];
398 double m1 = mProd[1];
399 double m2 = mProd[2];
402 if (m1 + m2 + mSafety > m0)
return false;
403 double e1 = 0.5 * (m0*m0 + m1*m1 - m2*m2) / m0;
404 double e2 = 0.5 * (m0*m0 + m2*m2 - m1*m1) / m0;
405 double pAbs = 0.5 * sqrtpos( (m0 - m1 - m2) * (m0 + m1 + m2)
406 * (m0 + m1 - m2) * (m0 - m1 + m2) ) / m0;
410 int iMother =
event[iProd[0]].mother1();
413 if (iMother <= 0 || iMother >= iProd[0]) meMode = 0;
415 int iDaughter1 =
event[iMother].daughter1();
416 int iDaughter2 =
event[iMother].daughter2();
417 if (iDaughter2 != iDaughter1 + 1) meMode = 0;
419 int idMother = abs( event[iMother].
id() );
420 if (idMother <= 100 || idMother%10 !=1
421 || (idMother/1000)%10 != 0) meMode = 0;
423 int iSister = (iProd[0] == iDaughter1) ? iDaughter2 : iDaughter1;
424 idSister = abs( event[iSister].
id() );
425 if ( (idSister <= 100 || idSister%10 !=1
426 || (idSister/1000)%10 != 0) && idSister != 22) meMode = 0;
433 double wtME, wtMEmax;
441 double cosTheta = 2. * rndmPtr->flat() - 1.;
442 double sinTheta = sqrt(1. - cosTheta*cosTheta);
443 double phi = 2. * M_PI * rndmPtr->flat();
444 double pX = pAbs * sinTheta * cos(phi);
445 double pY = pAbs * sinTheta * sin(phi);
446 double pZ = pAbs * cosTheta;
449 prod1.p( pX, pY, pZ, e1);
450 prod2.p( -pX, -pY, -pZ, e2);
451 prod1.bst( decayer.p(), decayer.m() );
452 prod2.bst( decayer.p(), decayer.m() );
458 double p10 = decayer.p() *
event[iMother].p();
459 double p12 = decayer.p() * prod1.p();
460 double p02 =
event[iMother].p() * prod1.p();
461 double s0 = pow2(event[iMother].m());
462 double s1 = pow2(decayer.m());
463 double s2 = pow2(prod1.m());
464 if (idSister != 22) wtME = pow2(p10 * p12 - s1 * p02);
465 else wtME = s1 * (2. * p10 * p12 * p02 - s1 * p02*p02
466 - s0 * p12*p12 - s2 * p10*p10 + s1 * s0 * s2);
467 wtME = max( wtME, 1e-6 * s1*s1 * s0 * s2);
468 wtMEmax = (p10*p10 - s1 * s0) * (p12*p12 - s1 * s2);
472 if(loop > NTRYMEWT) {
473 infoPtr->errorMsg(
"ParticleDecays::twoBody: "
474 "caught in infinite ME weight loop");
479 }
while ( wtME < rndmPtr->flat() * wtMEmax );
490 bool ParticleDecays::threeBody(
Event& event) {
493 Particle& decayer =
event[iProd[0]];
494 Particle& prod1 =
event[iProd[1]];
495 Particle& prod2 =
event[iProd[2]];
496 Particle& prod3 =
event[iProd[3]];
499 double m0 = mProd[0];
500 double m1 = mProd[1];
501 double m2 = mProd[2];
502 double m3 = mProd[3];
503 double mSum = m1 + m2 + m3;
504 double mDiff = m0 - mSum;
505 if (mDiff < mSafety)
return false;
508 double m23Min = m2 + m3;
509 double m23Max = m0 - m1;
510 double p1Max = 0.5 * sqrtpos( (m0 - m1 - m23Min) * (m0 + m1 + m23Min)
511 * (m0 + m1 - m23Min) * (m0 - m1 + m23Min) ) / m0;
512 double p23Max = 0.5 * sqrtpos( (m23Max - m2 - m3) * (m23Max + m2 + m3)
513 * (m23Max + m2 - m3) * (m23Max - m2 + m3) ) / m23Max;
514 double wtPSmax = 0.5 * p1Max * p23Max;
517 double wtME, wtMEmax, wtPS, m23, p1Abs, p23Abs;
524 m23 = m23Min + rndmPtr->flat() * mDiff;
527 p1Abs = 0.5 * sqrtpos( (m0 - m1 - m23) * (m0 + m1 + m23)
528 * (m0 + m1 - m23) * (m0 - m1 + m23) ) / m0;
529 p23Abs = 0.5 * sqrtpos( (m23 - m2 - m3) * (m23 + m2 + m3)
530 * (m23 + m2 - m3) * (m23 - m2 + m3) ) / m23;
531 wtPS = p1Abs * p23Abs;
534 }
while ( wtPS < rndmPtr->flat() * wtPSmax );
537 double cosTheta = 2. * rndmPtr->flat() - 1.;
538 double sinTheta = sqrt(1. - cosTheta*cosTheta);
539 double phi = 2. * M_PI * rndmPtr->flat();
540 double pX = p23Abs * sinTheta * cos(phi);
541 double pY = p23Abs * sinTheta * sin(phi);
542 double pZ = p23Abs * cosTheta;
543 double e2 = sqrt( m2*m2 + p23Abs*p23Abs);
544 double e3 = sqrt( m3*m3 + p23Abs*p23Abs);
545 prod2.p( pX, pY, pZ, e2);
546 prod3.p( -pX, -pY, -pZ, e3);
549 cosTheta = 2. * rndmPtr->flat() - 1.;
550 sinTheta = sqrt(1. - cosTheta*cosTheta);
551 phi = 2. * M_PI * rndmPtr->flat();
552 pX = p1Abs * sinTheta * cos(phi);
553 pY = p1Abs * sinTheta * sin(phi);
554 pZ = p1Abs * cosTheta;
555 double e1 = sqrt( m1*m1 + p1Abs*p1Abs);
556 double e23 = sqrt( m23*m23 + p1Abs*p1Abs);
557 prod1.p( pX, pY, pZ, e1);
560 Vec4 p23( -pX, -pY, -pZ, e23);
561 prod2.bst( p23, m23 );
562 prod3.bst( p23, m23 );
566 double p1p2 = prod1.p() * prod2.p();
567 double p1p3 = prod1.p() * prod3.p();
568 double p2p3 = prod2.p() * prod3.p();
569 wtME = pow2(m1 * m2 * m3) - pow2(m1 * p2p3) - pow2(m2 * p1p3)
570 - pow2(m3 * p1p2) + 2. * p1p2 * p1p3 * p2p3;
571 wtMEmax = pow3(m0 * m0) / 150.;
574 }
else if (meMode == 21) {
575 double x1 = 2. * prod1.e() / m0;
576 wtME = x1 * (3. - 2. * x1);
577 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
578 wtMEmax = xMax * (3. - 2. * xMax);
581 }
else if ((meMode == 22 || meMode == 23) && prod1.isLepton()) {
582 wtME = m0 * prod1.e() * (prod2.p() * prod3.p());
583 wtMEmax = min( pow4(m0) / 16., m0 * (m0 - m1 - m2) * (m0 - m1 - m3)
587 }
else if (meMode == 22 || meMode == 23) {
588 double x1 = 2. * prod1.pAbs() / m0;
589 wtME = x1 * (3. - 2. * x1);
590 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
591 wtMEmax = xMax * (3. - 2. * xMax);
594 }
else if (meMode == 31) {
595 double x1 = 2. * prod1.e() / m0;
597 double x1Max = 1. - pow2(mSum / m0);
598 wtMEmax = pow3(x1Max);
601 }
else if (meMode == 92) {
602 double x1 = 2. * prod1.e() / m0;
603 double x2 = 2. * prod2.e() / m0;
604 double x3 = 2. * prod3.e() / m0;
605 wtME = pow2( (1. - x1) / (x2 * x3) ) + pow2( (1. - x2) / (x1 * x3) )
606 + pow2( (1. - x3) / (x1 * x2) );
609 if (prod1.id() == 22 && sqrt(1. - x1) * m0 < 2. * stopMass) wtME = 0.;
610 if (prod2.id() == 22 && sqrt(1. - x2) * m0 < 2. * stopMass) wtME = 0.;
611 if (prod3.id() == 22 && sqrt(1. - x3) * m0 < 2. * stopMass) wtME = 0.;
615 }
while ( wtME < rndmPtr->flat() * wtMEmax );
618 prod1.bst( decayer.p(), decayer.m() );
619 prod2.bst( decayer.p(), decayer.m() );
620 prod3.bst( decayer.p(), decayer.m() );
631 bool ParticleDecays::mGenerator(
Event& event) {
634 double m0 = mProd[0];
635 double mSum = mProd[1];
636 for (
int i = 2; i <= mult; ++i) mSum += mProd[i];
637 double mDiff = m0 - mSum;
638 if (mDiff < mSafety)
return false;
642 for (
int i = 0; i <= mult; ++i) mInv.push_back( mProd[i]);
645 double wtPS, wtME, wtMEmax;
646 double wtPSmax = 1. / WTCORRECTION[mult];
647 double mMax = mDiff + mProd[mult];
649 for (
int i = mult - 1; i > 0; --i) {
652 double mNow = mProd[i];
653 wtPSmax *= 0.5 * sqrtpos( (mMax - mMin - mNow) * (mMax + mMin + mNow)
654 * (mMax + mMin - mNow) * (mMax - mMin + mNow) ) / mMax;
668 rndmOrd.push_back(1.);
669 for (
int i = 1; i < mult - 1; ++i) {
670 double rndm = rndmPtr->flat();
671 rndmOrd.push_back(rndm);
672 for (
int j = i - 1; j > 0; --j) {
673 if (rndm > rndmOrd[j]) swap( rndmOrd[j], rndmOrd[j+1] );
677 rndmOrd.push_back(0.);
680 for (
int i = mult - 1; i > 0; --i) {
681 mInv[i] = mInv[i+1] + mProd[i] + (rndmOrd[i-1] - rndmOrd[i]) * mDiff;
682 wtPS *= 0.5 * sqrtpos( (mInv[i] - mInv[i+1] - mProd[i])
683 * (mInv[i] + mInv[i+1] + mProd[i]) * (mInv[i] + mInv[i+1] - mProd[i])
684 * (mInv[i] - mInv[i+1] + mProd[i]) ) / mInv[i];
688 }
while ( wtPS < rndmPtr->flat() * wtPSmax );
691 pInv.resize(mult + 1);
692 for (
int i = 1; i < mult; ++i) {
693 double pAbs = 0.5 * sqrtpos( (mInv[i] - mInv[i+1] - mProd[i])
694 * (mInv[i] + mInv[i+1] + mProd[i]) * (mInv[i] + mInv[i+1] - mProd[i])
695 * (mInv[i] - mInv[i+1] + mProd[i]) ) / mInv[i];
698 double cosTheta = 2. * rndmPtr->flat() - 1.;
699 double sinTheta = sqrt(1. - cosTheta*cosTheta);
700 double phi = 2. * M_PI * rndmPtr->flat();
701 double pX = pAbs * sinTheta * cos(phi);
702 double pY = pAbs * sinTheta * sin(phi);
703 double pZ = pAbs * cosTheta;
706 double eHad = sqrt( mProd[i]*mProd[i] + pAbs*pAbs);
707 double eInv = sqrt( mInv[i+1]*mInv[i+1] + pAbs*pAbs);
708 event[iProd[i]].p( pX, pY, pZ, eHad);
709 pInv[i+1].p( -pX, -pY, -pZ, eInv);
713 event[iProd[mult]].p( pInv[mult] );
714 for (
int iFrame = mult - 1; iFrame > 1; --iFrame)
715 for (
int i = iFrame; i <= mult; ++i)
716 event[iProd[i]].bst( pInv[iFrame], mInv[iFrame]);
719 if (meMode == 21 && event[iProd[1]].isLepton()) {
720 double x1 = 2. *
event[iProd[1]].e() / m0;
721 wtME = x1 * (3. - 2. * x1);
722 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
723 wtMEmax = xMax * (3. - 2. * xMax);
727 }
else if ((meMode == 22 || meMode == 23) &&
event[iProd[1]].isLepton()) {
728 Vec4 pRest =
event[iProd[3]].p();
729 for (
int i = 4; i <= mult; ++i) pRest += event[iProd[i]].p();
730 wtME = m0 *
event[iProd[1]].e() * (
event[iProd[2]].p() * pRest);
731 for (
int i = 4; i <= mult; ++i) wtME
732 *= exp(- event[iProd[i]].pAbs2() / pow2(sigmaSoft) );
733 wtMEmax = pow4(m0) / 16.;
736 }
else if (meMode == 22 || meMode == 23) {
737 double x1 = 2. *
event[iProd[1]].pAbs() / m0;
738 wtME = x1 * (3. - 2. * x1);
739 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
740 wtMEmax = xMax * (3. - 2. * xMax);
743 }
else if (meMode == 31) {
744 double x1 = 2. *
event[iProd[1]].e() / m0;
746 double x1Max = 1. - pow2(mSum / m0);
747 wtMEmax = pow3(x1Max);
751 }
while ( wtME < rndmPtr->flat() * wtMEmax );
754 pInv[1].p( event[iProd[0]].p() );
755 for (
int i = 1; i <= mult; ++i) event[iProd[i]].bst( pInv[1], mInv[1] );
766 bool ParticleDecays::dalitzMass() {
770 for (
int i = 1; i <= mult - 2; ++i) mSum1 += mProd[i];
771 if (meMode == 13) mSum1 *= MSAFEDALITZ;
772 double mSum2 = MSAFEDALITZ * (mProd[mult -1] + mProd[mult]);
773 double mDiff = mProd[0] - mSum1 - mSum2;
776 if (mDiff < mSafety)
return false;
777 if (idProd[mult - 1] + idProd[mult] != 0
778 || mProd[mult - 1] != mProd[mult]) {
779 infoPtr->errorMsg(
"Error in ParticleDecays::dalitzMass:"
780 " inconsistent flavour/mass assignments");
783 if ( meMode == 13 && (idProd[1] + idProd[2] != 0
784 || mProd[1] != mProd[2]) ) {
785 infoPtr->errorMsg(
"Error in ParticleDecays::dalitzMass:"
786 " inconsistent flavour/mass assignments");
791 if (meMode == 11 || meMode == 12) {
794 double sGamMin = pow2(mSum2);
795 double sGamMax = pow2(mProd[0] - mSum1);
800 if (++loop > NTRYDALITZ)
return false;
801 sGam = sGamMin * pow( sGamMax / sGamMin, rndmPtr->flat() );
802 wtGam = (1. + 0.5 * sGamMin / sGam) * sqrt(1. - sGamMin / sGam)
803 * pow3(1. - sGam / sGamMax) * sRhoDal * (sRhoDal + wRhoDal)
804 / ( pow2(sGam - sRhoDal) + sRhoDal * wRhoDal );
805 }
while ( wtGam < rndmPtr->flat() );
809 mProd[mult] = sqrt(sGam);
815 double s0 = pow2(mProd[0]);
816 double s12Min = pow2(mSum1);
817 double s12Max = pow2(mProd[0] - mSum2);
818 double s34Min = pow2(mSum2);
819 double s34Max = pow2(mProd[0] - mSum1);
822 double s12, s34, wt12, wt34, wtPAbs, wtAll;
825 if (++loop > NTRYDALITZ)
return false;
826 s12 = s12Min * pow( s12Max / s12Min, rndmPtr->flat() );
827 wt12 = (1. + 0.5 * s12Min / s12) * sqrt(1. - s12Min / s12)
828 * sRhoDal * (sRhoDal + wRhoDal)
829 / ( pow2(s12 - sRhoDal) + sRhoDal * wRhoDal );
830 s34 = s34Min * pow( s34Max / s34Min, rndmPtr->flat() );
831 wt34 = (1. + 0.5 * s34Min / s34) * sqrt(1. - s34Min / s34)
832 * sRhoDal * (sRhoDal + wRhoDal)
833 / ( pow2(s34 - sRhoDal) + sRhoDal * wRhoDal );
834 wtPAbs = sqrtpos( pow2(1. - (s12 + s34)/ s0)
835 - 4. * s12 * s34 / (s0 * s0) );
836 wtAll = wt12 * wt34 * pow3(wtPAbs);
837 if (wtAll > 1.) infoPtr->errorMsg(
838 "Error in ParticleDecays::dalitzMass: weight > 1");
839 }
while (wtAll < rndmPtr->flat());
843 mProd[1] = sqrt(s12);
844 mProd[2] = sqrt(s34);
856 bool ParticleDecays::dalitzKinematics(
Event& event) {
859 int nDal = (meMode < 13) ? 1 : 2;
863 for (
int iDal = 0; iDal < nDal; ++iDal) {
866 Particle& decayer =
event[iProd[0]];
867 Particle& prodA = (iDal == 0) ? event[iProd[mult - 1]]
869 Particle& prodB = (iDal == 0) ? event[iProd[mult]]
873 Vec4 pDec = decayer.p();
874 int iGam = (meMode < 13) ? mult - 1 : 2 - iDal;
875 Vec4 pGam =
event[iProd[iGam]].p();
876 pGam.bstback( pDec, decayer.m() );
877 double phiGam = pGam.phi();
878 pGam.rot( 0., -phiGam);
879 double thetaGam = pGam.theta();
880 pGam.rot( -thetaGam, 0.);
883 double mGam = (meMode < 13) ? mProd[mult - 1] : mProd[2 - iDal];
884 double mA = prodA.m();
885 double mB = prodB.m();
886 double mGamMin = MSAFEDALITZ * (mA + mB);
887 double mGamRat = pow2(mGamMin / mGam);
888 double pGamAbs = 0.5 * sqrtpos( (mGam - mA - mB) * (mGam + mA + mB) );
891 double cosTheta, cos2Theta;
893 cosTheta = 2. * rndmPtr->flat() - 1.;
894 cos2Theta = cosTheta * cosTheta;
895 }
while ( 1. + cos2Theta + mGamRat * (1. - cos2Theta)
896 < 2. * rndmPtr->flat() );
897 double sinTheta = sqrt(1. - cosTheta*cosTheta);
898 double phi = 2. * M_PI * rndmPtr->flat();
899 double pX = pGamAbs * sinTheta * cos(phi);
900 double pY = pGamAbs * sinTheta * sin(phi);
901 double pZ = pGamAbs * cosTheta;
902 double eA = sqrt( mA*mA + pGamAbs*pGamAbs);
903 double eB = sqrt( mB*mB + pGamAbs*pGamAbs);
904 prodA.p( pX, pY, pZ, eA);
905 prodB.p( -pX, -pY, -pZ, eB);
908 prodA.bst( pGam, mGam);
909 prodB.bst( pGam, mGam);
910 prodA.rot( thetaGam, phiGam);
911 prodB.rot( thetaGam, phiGam);
912 prodA.bst( pDec, decayer.m() );
913 prodB.bst( pDec, decayer.m() );
925 bool ParticleDecays::pickHadrons() {
932 bool closedGLoop =
false;
933 for (
int i = 1; i <= mult; ++i) {
934 int idAbs = abs(idProd[i]);
935 if ( idAbs < 9 || (idAbs > 1000 && idAbs < 10000 && (idAbs/10)%10 == 0)
936 || idAbs == 81 || idAbs == 82 || idAbs == 83) {
938 idPartons.push_back(idProd[i]);
939 if (idAbs == 83) closedGLoop =
true;
943 idProd[nKnown] = idProd[i];
944 mProd[nKnown] = mProd[i];
950 for (
int i = 0; i < nPartons; ++i) {
951 int idPart = idPartons[i];
954 int idAbs = abs(idDec);
955 if ( (idAbs/1000)%10 == 0 ) {
956 idNew = -(idAbs/10)%10;
957 if ((idAbs/100)%2 == 1) idNew = -idNew;
958 }
else if ( (idAbs/100)%10 >= (idAbs/10)%10 )
959 idNew = 100 * ((idAbs/10)%100) + 3;
960 else idNew = 1000 * ((idAbs/10)%10) + 100 * ((idAbs/100)%10) + 1;
961 if (idDec < 0) idNew = -idNew;
964 }
else if (idPart == 82 || idPart == 83) {
967 int idDummy = -flavSelPtr->pickLightQ();
968 FlavContainer flavDummy(idDummy, idPart - 82);
969 do idNew = flavSelPtr->pick(flavDummy).id;
971 mFlav = particleDataPtr->constituentMass(idNew);
972 }
while (2. * mFlav + stopMass > mProd[0]);
973 }
else if (idPart == -82 || idPart == -83) {
974 idNew = -idPartons[i-1];
976 idPartons[i] = idNew;
980 int nMin = max( 2, nKnown + nPartons / 2);
981 if (meMode == 23) nMin = 3;
982 if (meMode > 41 && meMode <= 50) nMin = meMode - 40;
983 if (meMode > 51 && meMode <= 60) nMin = meMode - 50;
985 if (meMode == 0) nFix = nMin;
986 if (meMode == 11) nFix = 3;
987 if (meMode == 12) nFix = 4;
988 if (meMode > 61 && meMode <= 70) nFix = meMode - 60;
989 if (meMode > 71 && meMode <= 80) nFix = meMode - 70;
990 if (nFix > 0 && nKnown + nPartons/2 > nFix)
return false;
993 int nFilled, nTotal, nNew, nSpec, nLeft;
996 bool diquarkClash =
false;
997 bool usedChannel =
false;
1002 if (nTry > NTRYPICK)
return false;
1005 nFilled = nKnown + 1;
1006 idProd.resize(nFilled);
1007 mProd.resize(nFilled);
1012 for (
int i = 1; i < nFilled; ++i) mDiff -= mProd[i];
1013 diquarkClash =
false;
1014 usedChannel =
false;
1017 if ( (meMode == 22 || meMode == 23) && nLeft > 1) {
1018 FlavContainer flav1( idPartons[nPartons - 2] );
1019 FlavContainer flav2( idPartons[nPartons - 1] );
1021 do idHad = flavSelPtr->combine( flav1, flav2);
1023 double mHad = particleDataPtr->mSel(idHad);
1025 idProd.push_back( idHad);
1026 mProd.push_back( mHad);
1037 double geom = rndmPtr->flat();
1042 }
while (geom < 1. && nTotal < 10);
1045 }
else if (nFix == 0) {
1046 double multIncreaseNow = (meMode == 23)
1047 ? multIncreaseWeak : multIncrease;
1048 double mDiffPS = mDiff;
1049 for (
int i = 0; i < nLeft; ++i)
1050 mDiffPS -= particleDataPtr->constituentMass( idPartons[i] );
1051 double average = 0.5 * (nKnown + nSpec) + 0.25 * nPartons
1052 + multIncreaseNow * log( max( 1.1, mDiffPS / multRefMass ) );
1053 if (closedGLoop) average += multGoffset;
1058 for (
int nNow = nMin + 1; nNow <= 10; ++nNow) {
1059 value *= average / nNow;
1064 sum *= rndmPtr->flat();
1068 value *= average / nTotal;
1070 }
while (sum > 0. && nTotal < 10);
1076 nNew = nTotal - nKnown - nSpec;
1080 for (
int i = 0; i < nLeft; ++i) {
1081 flavEnds.push_back( FlavContainer(idPartons[i]) );
1082 if (abs(idPartons[i]) > 100) flavSelPtr->assignPopQ( flavEnds[i] );
1086 if (nNew > nLeft/2) {
1087 FlavContainer flavNew;
1089 for (
int i = 0; i < nNew - nLeft/2; ++i) {
1091 int iEnd = int( (nLeft - 1.) * rndmPtr->flat() );
1094 flavNew = flavSelPtr->pick( flavEnds[iEnd] );
1095 idHad = flavSelPtr->combine( flavEnds[iEnd], flavNew);
1096 }
while (idHad == 0);
1098 idProd.push_back( idHad);
1099 flavEnds[iEnd].anti(flavNew);
1106 if ( abs(flavEnds[0].
id) > 8 && abs(flavEnds[1].
id) > 8)
1107 diquarkClash =
true;
1109 do idHad = flavSelPtr->combine( flavEnds[0], flavEnds[1]);
1111 idProd.push_back( idHad);
1120 if ( rndmPtr->flat() < colRearrange) iEnd2 = 3;
1122 ( (flavEnds[iEnd1].id > 0 && flavEnds[iEnd1].id < 9)
1123 || flavEnds[iEnd1].
id < -10 ) ? 1 : -1;
1124 if ( (flavEnds[iEnd2].
id < 0 && flavEnds[iEnd2].
id > -9)
1125 || flavEnds[iEnd2].id > 10 ) relColSign *= -1;
1126 if (relColSign == 1) iEnd2 = 2;
1127 if (iEnd2 == 2) iEnd3 = 1;
1128 if (iEnd2 == 3) iEnd4 = 1;
1132 if ( abs(flavEnds[iEnd1].
id) > 8 && abs(flavEnds[iEnd2].
id) > 8)
1133 diquarkClash =
true;
1135 do idHad = flavSelPtr->combine( flavEnds[iEnd1], flavEnds[iEnd2]);
1137 idProd.push_back( idHad);
1139 if ( abs(flavEnds[iEnd3].
id) > 8 && abs(flavEnds[iEnd4].
id) > 8)
1140 diquarkClash =
true;
1142 do idHad = flavSelPtr->combine( flavEnds[iEnd3], flavEnds[iEnd4]);
1144 idProd.push_back( idHad);
1149 for (
int i = nFilled; i < int(idProd.size()) ; ++i) {
1150 double mHad = particleDataPtr->mSel(idProd[i]);
1151 mProd.push_back( mHad);
1157 if ( (meMode > 61 && meMode <= 80) && mDiff > mSafety && !diquarkClash ) {
1158 int idMatch[10], idPNow;
1159 usedChannel =
false;
1160 bool matched =
false;
1162 for (
int i = 0; i < decDataPtr->sizeChannels(); ++i) {
1163 DecayChannel& channel = decDataPtr->channel(i);
1164 if (channel.multiplicity() != nTotal)
continue;
1165 for (
int k = 0; k < nTotal; ++k) idMatch[k] = channel.product(k);
1168 for (
int j = 0; j < nTotal; ++j) {
1170 idPNow = idProd[j + 1];
1171 if (idPNow == -311 || idPNow == 130 || idPNow == 310) idPNow = 311;
1172 for (
int k = 0; k < nTotal - j; ++k)
1173 if (idMatch[k] == idPNow || (idMatch[k] == -311 && idPNow == 311)) {
1175 idMatch[k] = idMatch[nTotal - 1 - j];
1179 if (!matched)
break;
1182 if (matched) {usedChannel =
true;
break;}
1187 }
while (mDiff < mSafety || diquarkClash || usedChannel);
1190 mult = idProd.size() - 1;
1193 if (meMode == 11 || meMode == 12) {
1196 for (
int i = 1; i <= mult; ++i) {
1197 if (idProd[i] == 11 || idProd[i] == 13 || idProd[i] == 15) iL1 = i;
1198 if (idProd[i] == -11 || idProd[i] == -13 || idProd[i] == -15) iL2 = i;
1200 if (iL1 > 0 && iL2 > 0) {
1201 int idL1 = idProd[iL1];
1202 int idL2 = idProd[iL2];
1203 double mL1 = mProd[iL1];
1204 double mL2 = mProd[iL2];
1206 for (
int i = 1; i <= mult; ++i)
if (i != iL1 && i != iL2) {
1208 idProd[iMove] = idProd[i];
1209 mProd[iMove] = mProd[i];
1211 idProd[mult - 1] = idL1;
1212 idProd[mult] = idL2;
1213 mProd[mult - 1] = mL1;
1227 bool ParticleDecays::setColours(
Event& event) {
1230 if (meMode == 91 && idProd[1] > 0 && idProd[1] < 9) {
1231 int newCol =
event.nextColTag();
1234 }
else if (meMode == 91 && idProd[1] < 0 && idProd[1] > -9) {
1235 int newCol =
event.nextColTag();
1240 }
else if (meMode == 91 && idProd[1] == 21) {
1241 int newCol1 =
event.nextColTag();
1242 int newCol2 =
event.nextColTag();
1249 }
else if (meMode == 92 && idProd[1] == 21 && idProd[2] == 21
1250 && idProd[3] == 21) {
1251 int newCol1 =
event.nextColTag();
1252 int newCol2 =
event.nextColTag();
1253 int newCol3 =
event.nextColTag();
1262 }
else if (meMode == 92) {
1263 int iGlu1 = (idProd[1] == 21) ? 1 : 3;
1264 int iGlu2 = (idProd[2] == 21) ? 2 : 3;
1265 int newCol1 =
event.nextColTag();
1266 int newCol2 =
event.nextColTag();
1267 cols[iGlu1] = newCol1;
1268 acols[iGlu1] = newCol2;
1269 cols[iGlu2] = newCol2;
1270 acols[iGlu2] = newCol1;
1273 }
else return false;