9 #include "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");
111 sophisticatedTau = settings.mode(
"ParticleDecays:sophisticatedTau");
114 if (sophisticatedTau) tauDecayer.init(infoPtr, &settings,
115 particleDataPtr, rndmPtr, couplingsPtr);
123 bool ParticleDecays::decay(
int iDec,
Event& event) {
126 Particle& decayer =
event[iDec];
129 if (limitDecay && !checkVertex(decayer))
return true;
132 idDec = decayer.id();
136 iProd.push_back( iDec );
137 idProd.push_back( idDec );
138 mProd.push_back( decayer.m() );
141 bool hasOscillated = (abs(idDec) == 511 || abs(idDec) == 531)
142 ? oscillateB(decayer) :
false;
143 if (hasOscillated) {idDec = - idDec; idProd[0] = idDec;}
146 decDataPtr = &decayer.particleDataEntry();
149 bool doneExternally =
false;
150 if (decDataPtr->doExternalDecay()) {
152 pProd.push_back(decayer.p());
153 doneExternally = decayHandlePtr->decay(idProd, mProd, pProd,
157 if (doneExternally) {
158 mult = idProd.size() - 1;
159 int status = (hasOscillated) ? 94 : 93;
160 for (
int i = 1; i <= mult; ++i) {
161 int iPos =
event.append( idProd[i], status, iDec, 0, 0, 0,
162 0, 0, pProd[i], mProd[i]);
163 iProd.push_back( iPos);
167 event[iDec].statusNeg();
168 event[iDec].daughters( iProd[1], iProd[mult]);
173 if (decayer.idAbs() == 15 && !doneExternally && sophisticatedTau) {
174 doneExternally = tauDecayer.decay(iDec, event);
175 if (doneExternally)
return true;
179 if (!doneExternally) {
182 if (!decDataPtr->preparePick(idDec))
return false;
183 bool foundChannel =
false;
184 bool hasStored =
false;
185 for (
int iTryChannel = 0; iTryChannel < NTRYDECAY; ++iTryChannel) {
188 if (hasStored)
event.popBack(mult);
192 DecayChannel& channel = decDataPtr->pickChannel();
193 meMode = channel.meMode();
194 keepPartons = (meMode > 90 && meMode <= 100);
195 mult = channel.multiplicity();
198 bool foundMode =
false;
199 for (
int iTryMode = 0; iTryMode < NTRYDECAY; ++iTryMode) {
206 for (
int i = 0; i < mult; ++i) {
207 int idNow = channel.product(i);
208 int idAbs = abs(idNow);
209 if ( idAbs < 10 || idAbs == 21 || idAbs == 81 || idAbs == 82
210 || idAbs == 83 || (idAbs > 1000 && idAbs < 10000
211 && (idAbs/10)%10 == 0) ) hasPartons =
true;
212 if (idDec < 0 && particleDataPtr->hasAnti(idNow)) idNow = -idNow;
213 double mNow = particleDataPtr->mass(idNow);
214 idProd.push_back( idNow);
215 mProd.push_back( mNow);
219 if (hasPartons && !keepPartons && !pickHadrons())
continue;
224 for (
int i = 0; i <= mult; ++i) {
228 if (hasPartons && keepPartons && !setColours(event))
continue;
232 double mDiff = mProd[0];
233 for (
int i = 1; i <= mult; ++i) mDiff -= mProd[i];
234 if (mDiff < mSafety)
continue;
241 if (!foundMode)
continue;
244 int status = (hasOscillated) ? 92 : 91;
245 for (
int i = 1; i <= mult; ++i) {
246 int iPos =
event.append( idProd[i], status, iDec, 0, 0, 0,
247 cols[i], acols[i], Vec4(0., 0., 0., 0.), mProd[i], scale);
248 iProd.push_back( iPos);
253 if ( (meMode == 11 || meMode == 12 || meMode == 13)
254 && !dalitzMass() )
continue;
257 bool decayed =
false;
258 if (mult == 1) decayed = oneBody(event);
259 else if (mult == 2) decayed = twoBody(event);
260 else if (mult == 3) decayed = threeBody(event);
261 else decayed = mGenerator(event);
262 if (!decayed)
continue;
265 if (meMode == 11 || meMode == 12 || meMode == 13)
266 dalitzKinematics(event);
275 event[iDec].statusNeg();
276 event[iDec].daughters( iProd[1], iProd[mult]);
280 if (hasStored)
event.popBack(mult);
281 infoPtr->errorMsg(
"Error in ParticleDecays::decay: "
282 "failed to find workable decay channel");
290 if (event[iDec].hasVertex() || event[iDec].tau() > 0.) {
291 Vec4 vDec =
event[iDec].vDec();
292 for (
int i = 1; i <= mult; ++i) event[iProd[i]].vProd( vDec );
296 for (
int i = 1; i <= mult; ++i)
297 event[iProd[i]].tau( event[iProd[i]].tau0() * rndmPtr->exp() );
301 if (hasPartons && keepPartons && doFSRinDecays)
302 timesDecPtr->shower( iProd[1], iProd.back(), event, mProd[0]);
305 else if (event[iDec].idAbs() > 4900000 &&
event[iDec].idAbs() < 5000000
306 && doFSRinDecays && mult == 2 &&
event[iProd[1]].isLepton()) {
307 event[iProd[1]].scale(mProd[0]);
308 event[iProd[2]].scale(mProd[0]);
309 timesDecPtr->shower( iProd[1], iProd.back(), event, mProd[0]);
321 bool ParticleDecays::checkVertex(Particle& decayer) {
324 if (limitTau0 && decayer.tau0() > tau0Max)
return false;
325 if (limitTau && decayer.tau() > tauMax)
return false;
326 if (limitRadius && pow2(decayer.xDec()) + pow2(decayer.yDec())
327 + pow2(decayer.zDec()) > pow2(rMax))
return false;
328 if (limitCylinder && (pow2(decayer.xDec()) + pow2(decayer.yDec())
329 > pow2(xyMax) || abs(decayer.zDec()) > zMax) )
return false;
340 bool ParticleDecays::oscillateB(Particle& decayer) {
343 double xBmix = (abs(decayer.id()) == 511) ? xBdMix : xBsMix;
344 double tau = decayer.tau();
345 double tau0 = decayer.tau0();
346 double probosc = pow2(sin(0.5 * xBmix * tau / tau0));
347 return (probosc > rndmPtr->flat());
355 bool ParticleDecays::oneBody(
Event& event) {
358 Particle& decayer =
event[iProd[0]];
359 Particle& prod =
event[iProd[1]];
362 prod.p( decayer.p() );
363 prod.m( decayer.m() );
364 prod.mother2( iProd[0] );
375 bool ParticleDecays::twoBody(
Event& event) {
378 Particle& decayer =
event[iProd[0]];
379 Particle& prod1 =
event[iProd[1]];
380 Particle& prod2 =
event[iProd[2]];
383 double m0 = mProd[0];
384 double m1 = mProd[1];
385 double m2 = mProd[2];
388 if (m1 + m2 + mSafety > m0)
return false;
389 double e1 = 0.5 * (m0*m0 + m1*m1 - m2*m2) / m0;
390 double e2 = 0.5 * (m0*m0 + m2*m2 - m1*m1) / m0;
391 double pAbs = 0.5 * sqrtpos( (m0 - m1 - m2) * (m0 + m1 + m2)
392 * (m0 + m1 - m2) * (m0 - m1 + m2) ) / m0;
396 int iMother =
event[iProd[0]].mother1();
399 if (iMother <= 0 || iMother >= iProd[0]) meMode = 0;
401 int iDaughter1 =
event[iMother].daughter1();
402 int iDaughter2 =
event[iMother].daughter2();
403 if (iDaughter2 != iDaughter1 + 1) meMode = 0;
405 int idMother = abs( event[iMother].
id() );
406 if (idMother <= 100 || idMother%10 !=1
407 || (idMother/1000)%10 != 0) meMode = 0;
409 int iSister = (iProd[0] == iDaughter1) ? iDaughter2 : iDaughter1;
410 idSister = abs( event[iSister].
id() );
411 if ( (idSister <= 100 || idSister%10 !=1
412 || (idSister/1000)%10 != 0) && idSister != 22) meMode = 0;
419 double wtME, wtMEmax;
427 double cosTheta = 2. * rndmPtr->flat() - 1.;
428 double sinTheta = sqrt(1. - cosTheta*cosTheta);
429 double phi = 2. * M_PI * rndmPtr->flat();
430 double pX = pAbs * sinTheta * cos(phi);
431 double pY = pAbs * sinTheta * sin(phi);
432 double pZ = pAbs * cosTheta;
435 prod1.p( pX, pY, pZ, e1);
436 prod2.p( -pX, -pY, -pZ, e2);
437 prod1.bst( decayer.p(), decayer.m() );
438 prod2.bst( decayer.p(), decayer.m() );
444 double p10 = decayer.p() *
event[iMother].p();
445 double p12 = decayer.p() * prod1.p();
446 double p02 =
event[iMother].p() * prod1.p();
447 double s0 = pow2(event[iMother].m());
448 double s1 = pow2(decayer.m());
449 double s2 = pow2(prod1.m());
450 if (idSister != 22) wtME = pow2(p10 * p12 - s1 * p02);
451 else wtME = s1 * (2. * p10 * p12 * p02 - s1 * p02*p02
452 - s0 * p12*p12 - s2 * p10*p10 + s1 * s0 * s2);
453 wtME = max( wtME, 1e-6 * s1*s1 * s0 * s2);
454 wtMEmax = (p10*p10 - s1 * s0) * (p12*p12 - s1 * s2);
458 if(loop > NTRYMEWT) {
459 infoPtr->errorMsg(
"ParticleDecays::twoBody: "
460 "caught in infinite ME weight loop");
465 }
while ( wtME < rndmPtr->flat() * wtMEmax );
476 bool ParticleDecays::threeBody(
Event& event) {
479 Particle& decayer =
event[iProd[0]];
480 Particle& prod1 =
event[iProd[1]];
481 Particle& prod2 =
event[iProd[2]];
482 Particle& prod3 =
event[iProd[3]];
485 double m0 = mProd[0];
486 double m1 = mProd[1];
487 double m2 = mProd[2];
488 double m3 = mProd[3];
489 double mSum = m1 + m2 + m3;
490 double mDiff = m0 - mSum;
491 if (mDiff < mSafety)
return false;
494 double m23Min = m2 + m3;
495 double m23Max = m0 - m1;
496 double p1Max = 0.5 * sqrtpos( (m0 - m1 - m23Min) * (m0 + m1 + m23Min)
497 * (m0 + m1 - m23Min) * (m0 - m1 + m23Min) ) / m0;
498 double p23Max = 0.5 * sqrtpos( (m23Max - m2 - m3) * (m23Max + m2 + m3)
499 * (m23Max + m2 - m3) * (m23Max - m2 + m3) ) / m23Max;
500 double wtPSmax = 0.5 * p1Max * p23Max;
503 double wtME, wtMEmax, wtPS, m23, p1Abs, p23Abs;
510 m23 = m23Min + rndmPtr->flat() * mDiff;
513 p1Abs = 0.5 * sqrtpos( (m0 - m1 - m23) * (m0 + m1 + m23)
514 * (m0 + m1 - m23) * (m0 - m1 + m23) ) / m0;
515 p23Abs = 0.5 * sqrtpos( (m23 - m2 - m3) * (m23 + m2 + m3)
516 * (m23 + m2 - m3) * (m23 - m2 + m3) ) / m23;
517 wtPS = p1Abs * p23Abs;
520 }
while ( wtPS < rndmPtr->flat() * wtPSmax );
523 double cosTheta = 2. * rndmPtr->flat() - 1.;
524 double sinTheta = sqrt(1. - cosTheta*cosTheta);
525 double phi = 2. * M_PI * rndmPtr->flat();
526 double pX = p23Abs * sinTheta * cos(phi);
527 double pY = p23Abs * sinTheta * sin(phi);
528 double pZ = p23Abs * cosTheta;
529 double e2 = sqrt( m2*m2 + p23Abs*p23Abs);
530 double e3 = sqrt( m3*m3 + p23Abs*p23Abs);
531 prod2.p( pX, pY, pZ, e2);
532 prod3.p( -pX, -pY, -pZ, e3);
535 cosTheta = 2. * rndmPtr->flat() - 1.;
536 sinTheta = sqrt(1. - cosTheta*cosTheta);
537 phi = 2. * M_PI * rndmPtr->flat();
538 pX = p1Abs * sinTheta * cos(phi);
539 pY = p1Abs * sinTheta * sin(phi);
540 pZ = p1Abs * cosTheta;
541 double e1 = sqrt( m1*m1 + p1Abs*p1Abs);
542 double e23 = sqrt( m23*m23 + p1Abs*p1Abs);
543 prod1.p( pX, pY, pZ, e1);
546 Vec4 p23( -pX, -pY, -pZ, e23);
547 prod2.bst( p23, m23 );
548 prod3.bst( p23, m23 );
552 double p1p2 = prod1.p() * prod2.p();
553 double p1p3 = prod1.p() * prod3.p();
554 double p2p3 = prod2.p() * prod3.p();
555 wtME = pow2(m1 * m2 * m3) - pow2(m1 * p2p3) - pow2(m2 * p1p3)
556 - pow2(m3 * p1p2) + 2. * p1p2 * p1p3 * p2p3;
557 wtMEmax = pow3(m0 * m0) / 150.;
560 }
else if (meMode == 21) {
561 double x1 = 2. * prod1.e() / m0;
562 wtME = x1 * (3. - 2. * x1);
563 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
564 wtMEmax = xMax * (3. - 2. * xMax);
567 }
else if ((meMode == 22 || meMode == 23) && prod1.isLepton()) {
568 wtME = m0 * prod1.e() * (prod2.p() * prod3.p());
569 wtMEmax = min( pow4(m0) / 16., m0 * (m0 - m1 - m2) * (m0 - m1 - m3)
573 }
else if (meMode == 22 || meMode == 23) {
574 double x1 = 2. * prod1.pAbs() / m0;
575 wtME = x1 * (3. - 2. * x1);
576 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
577 wtMEmax = xMax * (3. - 2. * xMax);
580 }
else if (meMode == 31) {
581 double x1 = 2. * prod1.e() / m0;
583 double x1Max = 1. - pow2(mSum / m0);
584 wtMEmax = pow3(x1Max);
587 }
else if (meMode == 92) {
588 double x1 = 2. * prod1.e() / m0;
589 double x2 = 2. * prod2.e() / m0;
590 double x3 = 2. * prod3.e() / m0;
591 wtME = pow2( (1. - x1) / (x2 * x3) ) + pow2( (1. - x2) / (x1 * x3) )
592 + pow2( (1. - x3) / (x1 * x2) );
595 if (prod1.id() == 22 && sqrt(1. - x1) * m0 < 2. * stopMass) wtME = 0.;
596 if (prod2.id() == 22 && sqrt(1. - x2) * m0 < 2. * stopMass) wtME = 0.;
597 if (prod3.id() == 22 && sqrt(1. - x3) * m0 < 2. * stopMass) wtME = 0.;
601 }
while ( wtME < rndmPtr->flat() * wtMEmax );
604 prod1.bst( decayer.p(), decayer.m() );
605 prod2.bst( decayer.p(), decayer.m() );
606 prod3.bst( decayer.p(), decayer.m() );
617 bool ParticleDecays::mGenerator(
Event& event) {
620 double m0 = mProd[0];
621 double mSum = mProd[1];
622 for (
int i = 2; i <= mult; ++i) mSum += mProd[i];
623 double mDiff = m0 - mSum;
624 if (mDiff < mSafety)
return false;
628 for (
int i = 0; i <= mult; ++i) mInv.push_back( mProd[i]);
631 double wtPS, wtME, wtMEmax;
632 double wtPSmax = 1. / WTCORRECTION[mult];
633 double mMax = mDiff + mProd[mult];
635 for (
int i = mult - 1; i > 0; --i) {
638 double mNow = mProd[i];
639 wtPSmax *= 0.5 * sqrtpos( (mMax - mMin - mNow) * (mMax + mMin + mNow)
640 * (mMax + mMin - mNow) * (mMax - mMin + mNow) ) / mMax;
654 rndmOrd.push_back(1.);
655 for (
int i = 1; i < mult - 1; ++i) {
656 double rndm = rndmPtr->flat();
657 rndmOrd.push_back(rndm);
658 for (
int j = i - 1; j > 0; --j) {
659 if (rndm > rndmOrd[j]) swap( rndmOrd[j], rndmOrd[j+1] );
663 rndmOrd.push_back(0.);
666 for (
int i = mult - 1; i > 0; --i) {
667 mInv[i] = mInv[i+1] + mProd[i] + (rndmOrd[i-1] - rndmOrd[i]) * mDiff;
668 wtPS *= 0.5 * sqrtpos( (mInv[i] - mInv[i+1] - mProd[i])
669 * (mInv[i] + mInv[i+1] + mProd[i]) * (mInv[i] + mInv[i+1] - mProd[i])
670 * (mInv[i] - mInv[i+1] + mProd[i]) ) / mInv[i];
674 }
while ( wtPS < rndmPtr->flat() * wtPSmax );
677 pInv.resize(mult + 1);
678 for (
int i = 1; i < mult; ++i) {
679 double pAbs = 0.5 * sqrtpos( (mInv[i] - mInv[i+1] - mProd[i])
680 * (mInv[i] + mInv[i+1] + mProd[i]) * (mInv[i] + mInv[i+1] - mProd[i])
681 * (mInv[i] - mInv[i+1] + mProd[i]) ) / mInv[i];
684 double cosTheta = 2. * rndmPtr->flat() - 1.;
685 double sinTheta = sqrt(1. - cosTheta*cosTheta);
686 double phi = 2. * M_PI * rndmPtr->flat();
687 double pX = pAbs * sinTheta * cos(phi);
688 double pY = pAbs * sinTheta * sin(phi);
689 double pZ = pAbs * cosTheta;
692 double eHad = sqrt( mProd[i]*mProd[i] + pAbs*pAbs);
693 double eInv = sqrt( mInv[i+1]*mInv[i+1] + pAbs*pAbs);
694 event[iProd[i]].p( pX, pY, pZ, eHad);
695 pInv[i+1].p( -pX, -pY, -pZ, eInv);
699 event[iProd[mult]].p( pInv[mult] );
700 for (
int iFrame = mult - 1; iFrame > 1; --iFrame)
701 for (
int i = iFrame; i <= mult; ++i)
702 event[iProd[i]].bst( pInv[iFrame], mInv[iFrame]);
705 if (meMode == 21 && event[iProd[1]].isLepton()) {
706 double x1 = 2. *
event[iProd[1]].e() / m0;
707 wtME = x1 * (3. - 2. * x1);
708 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
709 wtMEmax = xMax * (3. - 2. * xMax);
713 }
else if ((meMode == 22 || meMode == 23) &&
event[iProd[1]].isLepton()) {
714 Vec4 pRest =
event[iProd[3]].p();
715 for (
int i = 4; i <= mult; ++i) pRest += event[iProd[i]].p();
716 wtME = m0 *
event[iProd[1]].e() * (
event[iProd[2]].p() * pRest);
717 for (
int i = 4; i <= mult; ++i) wtME
718 *= exp(- event[iProd[i]].pAbs2() / pow2(sigmaSoft) );
719 wtMEmax = pow4(m0) / 16.;
722 }
else if (meMode == 22 || meMode == 23) {
723 double x1 = 2. *
event[iProd[1]].pAbs() / m0;
724 wtME = x1 * (3. - 2. * x1);
725 double xMax = min( 0.75, 2. * (1. - mSum / m0) );
726 wtMEmax = xMax * (3. - 2. * xMax);
729 }
else if (meMode == 31) {
730 double x1 = 2. *
event[iProd[1]].e() / m0;
732 double x1Max = 1. - pow2(mSum / m0);
733 wtMEmax = pow3(x1Max);
737 }
while ( wtME < rndmPtr->flat() * wtMEmax );
740 pInv[1].p( event[iProd[0]].p() );
741 for (
int i = 1; i <= mult; ++i) event[iProd[i]].bst( pInv[1], mInv[1] );
752 bool ParticleDecays::dalitzMass() {
756 for (
int i = 1; i <= mult - 2; ++i) mSum1 += mProd[i];
757 if (meMode == 13) mSum1 *= MSAFEDALITZ;
758 double mSum2 = MSAFEDALITZ * (mProd[mult -1] + mProd[mult]);
759 double mDiff = mProd[0] - mSum1 - mSum2;
762 if (mDiff < mSafety)
return false;
763 if (idProd[mult - 1] + idProd[mult] != 0
764 || mProd[mult - 1] != mProd[mult]) {
765 infoPtr->errorMsg(
"Error in ParticleDecays::dalitzMass:"
766 " inconsistent flavour/mass assignments");
769 if ( meMode == 13 && (idProd[1] + idProd[2] != 0
770 || mProd[1] != mProd[2]) ) {
771 infoPtr->errorMsg(
"Error in ParticleDecays::dalitzMass:"
772 " inconsistent flavour/mass assignments");
777 if (meMode == 11 || meMode == 12) {
780 double sGamMin = pow2(mSum2);
781 double sGamMax = pow2(mProd[0] - mSum1);
786 if (++loop > NTRYDALITZ)
return false;
787 sGam = sGamMin * pow( sGamMax / sGamMin, rndmPtr->flat() );
788 wtGam = (1. + 0.5 * sGamMin / sGam) * sqrt(1. - sGamMin / sGam)
789 * pow3(1. - sGam / sGamMax) * sRhoDal * (sRhoDal + wRhoDal)
790 / ( pow2(sGam - sRhoDal) + sRhoDal * wRhoDal );
791 }
while ( wtGam < rndmPtr->flat() );
795 mProd[mult] = sqrt(sGam);
801 double s0 = pow2(mProd[0]);
802 double s12Min = pow2(mSum1);
803 double s12Max = pow2(mProd[0] - mSum2);
804 double s34Min = pow2(mSum2);
805 double s34Max = pow2(mProd[0] - mSum1);
808 double s12, s34, wt12, wt34, wtPAbs, wtAll;
811 if (++loop > NTRYDALITZ)
return false;
812 s12 = s12Min * pow( s12Max / s12Min, rndmPtr->flat() );
813 wt12 = (1. + 0.5 * s12Min / s12) * sqrt(1. - s12Min / s12)
814 * sRhoDal * (sRhoDal + wRhoDal)
815 / ( pow2(s12 - sRhoDal) + sRhoDal * wRhoDal );
816 s34 = s34Min * pow( s34Max / s34Min, rndmPtr->flat() );
817 wt34 = (1. + 0.5 * s34Min / s34) * sqrt(1. - s34Min / s34)
818 * sRhoDal * (sRhoDal + wRhoDal)
819 / ( pow2(s34 - sRhoDal) + sRhoDal * wRhoDal );
820 wtPAbs = sqrtpos( pow2(1. - (s12 + s34)/ s0)
821 - 4. * s12 * s34 / (s0 * s0) );
822 wtAll = wt12 * wt34 * pow3(wtPAbs);
823 if (wtAll > 1.) infoPtr->errorMsg(
824 "Error in ParticleDecays::dalitzMass: weight > 1");
825 }
while (wtAll < rndmPtr->flat());
829 mProd[1] = sqrt(s12);
830 mProd[2] = sqrt(s34);
842 bool ParticleDecays::dalitzKinematics(
Event& event) {
845 int nDal = (meMode < 13) ? 1 : 2;
849 for (
int iDal = 0; iDal < nDal; ++iDal) {
852 Particle& decayer =
event[iProd[0]];
853 Particle& prodA = (iDal == 0) ? event[iProd[mult - 1]]
855 Particle& prodB = (iDal == 0) ? event[iProd[mult]]
859 Vec4 pDec = decayer.p();
860 int iGam = (meMode < 13) ? mult - 1 : 2 - iDal;
861 Vec4 pGam =
event[iProd[iGam]].p();
862 pGam.bstback( pDec, decayer.m() );
863 double phiGam = pGam.phi();
864 pGam.rot( 0., -phiGam);
865 double thetaGam = pGam.theta();
866 pGam.rot( -thetaGam, 0.);
869 double mGam = (meMode < 13) ? mProd[mult - 1] : mProd[2 - iDal];
870 double mA = prodA.m();
871 double mB = prodB.m();
872 double mGamMin = MSAFEDALITZ * (mA + mB);
873 double mGamRat = pow2(mGamMin / mGam);
874 double pGamAbs = 0.5 * sqrtpos( (mGam - mA - mB) * (mGam + mA + mB) );
877 double cosTheta, cos2Theta;
879 cosTheta = 2. * rndmPtr->flat() - 1.;
880 cos2Theta = cosTheta * cosTheta;
881 }
while ( 1. + cos2Theta + mGamRat * (1. - cos2Theta)
882 < 2. * rndmPtr->flat() );
883 double sinTheta = sqrt(1. - cosTheta*cosTheta);
884 double phi = 2. * M_PI * rndmPtr->flat();
885 double pX = pGamAbs * sinTheta * cos(phi);
886 double pY = pGamAbs * sinTheta * sin(phi);
887 double pZ = pGamAbs * cosTheta;
888 double eA = sqrt( mA*mA + pGamAbs*pGamAbs);
889 double eB = sqrt( mB*mB + pGamAbs*pGamAbs);
890 prodA.p( pX, pY, pZ, eA);
891 prodB.p( -pX, -pY, -pZ, eB);
894 prodA.bst( pGam, mGam);
895 prodB.bst( pGam, mGam);
896 prodA.rot( thetaGam, phiGam);
897 prodB.rot( thetaGam, phiGam);
898 prodA.bst( pDec, decayer.m() );
899 prodB.bst( pDec, decayer.m() );
911 bool ParticleDecays::pickHadrons() {
918 bool closedGLoop =
false;
919 for (
int i = 1; i <= mult; ++i) {
920 int idAbs = abs(idProd[i]);
921 if ( idAbs < 9 || (idAbs > 1000 && idAbs < 10000 && (idAbs/10)%10 == 0)
922 || idAbs == 81 || idAbs == 82 || idAbs == 83) {
924 idPartons.push_back(idProd[i]);
925 if (idAbs == 83) closedGLoop =
true;
929 idProd[nKnown] = idProd[i];
930 mProd[nKnown] = mProd[i];
936 for (
int i = 0; i < nPartons; ++i) {
937 int idPart = idPartons[i];
940 int idAbs = abs(idDec);
941 if ( (idAbs/1000)%10 == 0 ) {
942 idNew = -(idAbs/10)%10;
943 if ((idAbs/100)%2 == 1) idNew = -idNew;
944 }
else if ( (idAbs/100)%10 >= (idAbs/10)%10 )
945 idNew = 100 * ((idAbs/10)%100) + 3;
946 else idNew = 1000 * ((idAbs/10)%10) + 100 * ((idAbs/100)%10) + 1;
947 if (idDec < 0) idNew = -idNew;
950 }
else if (idPart == 82 || idPart == 83) {
953 int idDummy = -flavSelPtr->pickLightQ();
954 FlavContainer flavDummy(idDummy, idPart - 82);
955 do idNew = flavSelPtr->pick(flavDummy).id;
957 mFlav = particleDataPtr->constituentMass(idNew);
958 }
while (2. * mFlav + stopMass > mProd[0]);
959 }
else if (idPart == -82 || idPart == -83) {
960 idNew = -idPartons[i-1];
962 idPartons[i] = idNew;
966 int nMin = max( 2, nKnown + nPartons / 2);
967 if (meMode == 23) nMin = 3;
968 if (meMode > 41 && meMode <= 50) nMin = meMode - 40;
969 if (meMode > 51 && meMode <= 60) nMin = meMode - 50;
971 if (meMode == 0) nFix = nMin;
972 if (meMode == 11) nFix = 3;
973 if (meMode == 12) nFix = 4;
974 if (meMode > 61 && meMode <= 70) nFix = meMode - 60;
975 if (meMode > 71 && meMode <= 80) nFix = meMode - 70;
976 if (nFix > 0 && nKnown + nPartons/2 > nFix)
return false;
979 int nFilled, nTotal, nNew, nSpec, nLeft;
982 bool diquarkClash =
false;
983 bool usedChannel =
false;
988 if (nTry > NTRYPICK)
return false;
991 nFilled = nKnown + 1;
992 idProd.resize(nFilled);
993 mProd.resize(nFilled);
998 for (
int i = 1; i < nFilled; ++i) mDiff -= mProd[i];
999 diquarkClash =
false;
1000 usedChannel =
false;
1003 if ( (meMode == 22 || meMode == 23) && nLeft > 1) {
1004 FlavContainer flav1( idPartons[nPartons - 2] );
1005 FlavContainer flav2( idPartons[nPartons - 1] );
1007 do idHad = flavSelPtr->combine( flav1, flav2);
1009 double mHad = particleDataPtr->mass(idHad);
1011 idProd.push_back( idHad);
1012 mProd.push_back( mHad);
1023 double geom = rndmPtr->flat();
1028 }
while (geom < 1. && nTotal < 10);
1031 }
else if (nFix == 0) {
1032 double multIncreaseNow = (meMode == 23)
1033 ? multIncreaseWeak : multIncrease;
1034 double mDiffPS = mDiff;
1035 for (
int i = 0; i < nLeft; ++i)
1036 mDiffPS -= particleDataPtr->constituentMass( idPartons[i] );
1037 double average = 0.5 * (nKnown + nSpec) + 0.25 * nPartons
1038 + multIncreaseNow * log( max( 1.1, mDiffPS / multRefMass ) );
1039 if (closedGLoop) average += multGoffset;
1044 for (
int nNow = nMin + 1; nNow <= 10; ++nNow) {
1045 value *= average / nNow;
1050 sum *= rndmPtr->flat();
1054 value *= average / nTotal;
1056 }
while (sum > 0. && nTotal < 10);
1062 nNew = nTotal - nKnown - nSpec;
1066 for (
int i = 0; i < nLeft; ++i) {
1067 flavEnds.push_back( FlavContainer(idPartons[i]) );
1068 if (abs(idPartons[i]) > 100) flavSelPtr->assignPopQ( flavEnds[i] );
1072 if (nNew > nLeft/2) {
1073 FlavContainer flavNew;
1075 for (
int i = 0; i < nNew - nLeft/2; ++i) {
1077 int iEnd = int( (nLeft - 1.) * rndmPtr->flat() );
1080 flavNew = flavSelPtr->pick( flavEnds[iEnd] );
1081 idHad = flavSelPtr->combine( flavEnds[iEnd], flavNew);
1082 }
while (idHad == 0);
1084 idProd.push_back( idHad);
1085 flavEnds[iEnd].anti(flavNew);
1092 if ( abs(flavEnds[0].
id) > 8 && abs(flavEnds[1].
id) > 8)
1093 diquarkClash =
true;
1095 do idHad = flavSelPtr->combine( flavEnds[0], flavEnds[1]);
1097 idProd.push_back( idHad);
1106 if ( rndmPtr->flat() < colRearrange) iEnd2 = 3;
1108 ( (flavEnds[iEnd1].id > 0 && flavEnds[iEnd1].id < 9)
1109 || flavEnds[iEnd1].
id < -10 ) ? 1 : -1;
1110 if ( (flavEnds[iEnd2].
id < 0 && flavEnds[iEnd2].
id > -9)
1111 || flavEnds[iEnd2].id > 10 ) relColSign *= -1;
1112 if (relColSign == 1) iEnd2 = 2;
1113 if (iEnd2 == 2) iEnd3 = 1;
1114 if (iEnd2 == 3) iEnd4 = 1;
1118 if ( abs(flavEnds[iEnd1].
id) > 8 && abs(flavEnds[iEnd2].
id) > 8)
1119 diquarkClash =
true;
1121 do idHad = flavSelPtr->combine( flavEnds[iEnd1], flavEnds[iEnd2]);
1123 idProd.push_back( idHad);
1125 if ( abs(flavEnds[iEnd3].
id) > 8 && abs(flavEnds[iEnd4].
id) > 8)
1126 diquarkClash =
true;
1128 do idHad = flavSelPtr->combine( flavEnds[iEnd3], flavEnds[iEnd4]);
1130 idProd.push_back( idHad);
1135 for (
int i = nFilled; i < int(idProd.size()) ; ++i) {
1136 double mHad = particleDataPtr->mass(idProd[i]);
1137 mProd.push_back( mHad);
1143 if ( (meMode > 61 && meMode <= 80) && mDiff > mSafety && !diquarkClash ) {
1144 int idMatch[10], idPNow;
1145 usedChannel =
false;
1146 bool matched =
false;
1148 for (
int i = 0; i < decDataPtr->sizeChannels(); ++i) {
1149 DecayChannel& channel = decDataPtr->channel(i);
1150 if (channel.multiplicity() != nTotal)
continue;
1151 for (
int k = 0; k < nTotal; ++k) idMatch[k] = channel.product(k);
1154 for (
int j = 0; j < nTotal; ++j) {
1156 idPNow = idProd[j + 1];
1157 if (idPNow == -311 || idPNow == 130 || idPNow == 310) idPNow = 311;
1158 for (
int k = 0; k < nTotal - j; ++k)
1159 if (idMatch[k] == idPNow || (idMatch[k] == -311 && idPNow == 311)) {
1161 idMatch[k] = idMatch[nTotal - 1 - j];
1165 if (!matched)
break;
1168 if (matched) {usedChannel =
true;
break;}
1173 }
while (mDiff < mSafety || diquarkClash || usedChannel);
1176 mult = idProd.size() - 1;
1179 if (meMode == 11 || meMode == 12) {
1182 for (
int i = 1; i <= mult; ++i) {
1183 if (idProd[i] == 11 || idProd[i] == 13 || idProd[i] == 15) iL1 = i;
1184 if (idProd[i] == -11 || idProd[i] == -13 || idProd[i] == -15) iL2 = i;
1186 if (iL1 > 0 && iL2 > 0) {
1187 int idL1 = idProd[iL1];
1188 int idL2 = idProd[iL2];
1189 double mL1 = mProd[iL1];
1190 double mL2 = mProd[iL2];
1192 for (
int i = 1; i <= mult; ++i)
if (i != iL1 && i != iL2) {
1194 idProd[iMove] = idProd[i];
1195 mProd[iMove] = mProd[i];
1197 idProd[mult - 1] = idL1;
1198 idProd[mult] = idL2;
1199 mProd[mult - 1] = mL1;
1213 bool ParticleDecays::setColours(
Event& event) {
1216 if (meMode == 91 && idProd[1] > 0 && idProd[1] < 9) {
1217 int newCol =
event.nextColTag();
1220 }
else if (meMode == 91 && idProd[1] < 0 && idProd[1] > -9) {
1221 int newCol =
event.nextColTag();
1226 }
else if (meMode == 91 && idProd[1] == 21) {
1227 int newCol1 =
event.nextColTag();
1228 int newCol2 =
event.nextColTag();
1235 }
else if (meMode == 92 && idProd[1] == 21 && idProd[2] == 21
1236 && idProd[3] == 21) {
1237 int newCol1 =
event.nextColTag();
1238 int newCol2 =
event.nextColTag();
1239 int newCol3 =
event.nextColTag();
1248 }
else if (meMode == 92) {
1249 int iGlu1 = (idProd[1] == 21) ? 1 : 3;
1250 int iGlu2 = (idProd[2] == 21) ? 2 : 3;
1251 int newCol1 =
event.nextColTag();
1252 int newCol2 =
event.nextColTag();
1253 cols[iGlu1] = newCol1;
1254 acols[iGlu1] = newCol2;
1255 cols[iGlu2] = newCol2;
1256 acols[iGlu2] = newCol1;
1259 }
else return false;