8 #include "Pythia8/TauDecays.h"
22 const int TauDecays::NTRYCHANNEL = 10;
25 const int TauDecays::NTRYDECAY = 10000;
29 const double TauDecays::WTCORRECTION[11] = { 1., 1., 1.,
30 2., 5., 15., 60., 250., 1250., 7000., 50000. };
39 void TauDecays::init() {
42 hmeTwoFermions2W2TwoFermions
43 .initPointers(particleDataPtr, coupSMPtr, settingsPtr);
44 hmeTwoFermions2GammaZ2TwoFermions
45 .initPointers(particleDataPtr, coupSMPtr, settingsPtr);
47 .initPointers(particleDataPtr, coupSMPtr, settingsPtr);
49 .initPointers(particleDataPtr, coupSMPtr, settingsPtr);
51 .initPointers(particleDataPtr, coupSMPtr);
53 .initPointers(particleDataPtr, coupSMPtr, settingsPtr);
56 hmeTau2Meson .initPointers(particleDataPtr, coupSMPtr);
57 hmeTau2TwoLeptons .initPointers(particleDataPtr, coupSMPtr);
58 hmeTau2TwoMesonsViaVector .initPointers(particleDataPtr, coupSMPtr);
59 hmeTau2TwoMesonsViaVectorScalar.initPointers(particleDataPtr, coupSMPtr);
60 hmeTau2ThreePions .initPointers(particleDataPtr, coupSMPtr);
61 hmeTau2ThreeMesonsWithKaons .initPointers(particleDataPtr, coupSMPtr);
62 hmeTau2ThreeMesonsGeneric .initPointers(particleDataPtr, coupSMPtr);
63 hmeTau2TwoPionsGamma .initPointers(particleDataPtr, coupSMPtr);
64 hmeTau2FourPions .initPointers(particleDataPtr, coupSMPtr);
65 hmeTau2FivePions .initPointers(particleDataPtr, coupSMPtr);
66 hmeTau2PhaseSpace .initPointers(particleDataPtr, coupSMPtr);
69 tauExt = mode(
"TauDecays:externalMode");
70 tauMode = mode(
"TauDecays:mode");
71 tauMother = mode(
"TauDecays:tauMother");
72 tauPol = parm(
"TauDecays:tauPolarization");
75 limitTau0 = flag(
"ParticleDecays:limitTau0");
76 tau0Max = parm(
"ParticleDecays:tau0Max");
77 limitTau = flag(
"ParticleDecays:limitTau");
78 tauMax = parm(
"ParticleDecays:tauMax");
79 limitRadius = flag(
"ParticleDecays:limitRadius");
80 rMax = parm(
"ParticleDecays:rMax");
81 limitCylinder = flag(
"ParticleDecays:limitCylinder");
82 xyMax = parm(
"ParticleDecays:xyMax");
83 zMax = parm(
"ParticleDecays:zMax");
84 limitDecay = limitTau0 || limitTau || limitRadius || limitCylinder;
95 bool TauDecays::decay(
int idxOut1,
Event& event) {
98 out1 = HelicityParticle(event[idxOut1]);
99 int idxOut1Top = out1.iTopCopyId();
100 vector<int> sistersOut1 =
event[idxOut1Top].sisterList();
101 int idxOut2Top = idxOut1Top;
102 if (sistersOut1.size() == 1) idxOut2Top = sistersOut1[0];
105 int tau(-1), tnu(-1), lep(-1), lnu(-1);
106 for (
int i = 0; i < int(sistersOut1.size()); ++i) {
107 int sn = out1.id() == 15 ? -1 : 1;
108 int id =
event[sistersOut1[i]].id();
109 if (
id == sn * 15 && tau == -1) tau = sistersOut1[i];
110 else if (
id == sn * 16 && tnu == -1) tnu = sistersOut1[i];
111 else if ((
id == sn * 11 || (
id == sn * 13)) && lep == -1)
112 lep = sistersOut1[i];
113 else if ((
id == sn * 12 || (
id == sn * 14)) && lnu == -1)
114 lnu = sistersOut1[i];
116 if (tau > 0) idxOut2Top = tau;
117 else if (tnu > 0) idxOut2Top = tnu;
118 else if (lep > 0) idxOut2Top = lep;
119 else if (lnu > 0) idxOut2Top = lnu;
121 int idxOut2 =
event[idxOut2Top].iBotCopyId();
122 out2 = HelicityParticle(event[idxOut2]);
125 int idxMediator =
event[idxOut1Top].mother1();
126 mediator = HelicityParticle(event[idxMediator]);
127 mediator.direction = -1;
128 if (mediator.m() < out1.m() + out2.m()) {
129 Vec4 p = out1.p() + out2.p();
131 mediator.m(p.mCalc());
135 int idxMediatorTop = mediator.iTopCopyId();
136 int idxIn1 =
event[idxMediatorTop].mother1();
137 int idxIn2 =
event[idxMediatorTop].mother2();
138 in1 = HelicityParticle(event[idxIn1]);
140 in2 = HelicityParticle(event[idxIn2]);
145 particles.push_back(in1);
146 particles.push_back(in2);
147 particles.push_back(out1);
148 particles.push_back(out2);
152 if (idxOut1 != idxOut2 &&
153 (abs(out2.id()) == 11 || abs(out2.id()) == 13 || abs(out2.id()) == 15)) {
156 if (limitTau0 && out2.tau0() > tau0Max) correlated =
false;
157 else if (limitTau && out2.tau() > tauMax) correlated =
false;
158 else if (limitRadius && pow2(out2.xDec()) + pow2(out2.yDec())
159 + pow2(out2.zDec()) > pow2(rMax)) correlated =
false;
160 else if (limitCylinder && (pow2(out2.xDec()) + pow2(out2.yDec())
161 > pow2(xyMax) || abs(out2.zDec()) > zMax))
164 else if (!out2.canDecay()) correlated =
false;
165 else if (!out2.mayDecay()) correlated =
false;
171 if (tauMode == 4) known = internalMechanism(event);
172 else if (tauMode == 5) known = externalMechanism(event);
174 if ((tauMode == 2 && abs(mediator.id()) == tauMother) || tauMode == 3) {
177 double sign = out1.id() == -15 ? -1 : 1;
178 particles[2].rho[0][0] = (1 - sign * tauPol) / 2;
179 particles[2].rho[1][1] = (1 + sign * tauPol) / 2;
181 if (!externalMechanism(event)) known = internalMechanism(event);
188 particles[1] = mediator;
189 if (abs(mediator.id()) == 22)
190 hardME = hmeGamma2TwoFermions.initChannel(particles);
191 else if (abs(mediator.id()) == 23 || abs(mediator.id()) == 32)
192 hardME = hmeZ2TwoFermions.initChannel(particles);
193 else if (abs(mediator.id()) == 24 || abs(mediator.id()) == 34)
194 hardME = hmeW2TwoFermions.initChannel(particles);
195 else if (correlated) {
196 Vec4 p = out1.p() + out2.p();
197 particles[1] = HelicityParticle(22, -22, idxIn1, idxIn2, idxOut1,
198 idxOut2, 0, 0, p, p.mCalc(), 0, particleDataPtr);
199 hardME = hmeGamma2TwoFermions.initChannel(particles);
200 infoPtr->errorMsg(
"Warning in TauDecays::decay: unknown correlated "
201 "tau production, assuming from unpolarized photon");
203 infoPtr->errorMsg(
"Warning in TauDecays::decay: unknown uncorrelated "
204 "tau production, assuming unpolarized");
209 if (correlated && !out2.isFinal()) event[out2.index()].undoDecay();
212 HelicityParticle* tau;
214 if (correlated) idx = (rndmPtr->flat() < 0.5) ? 2 : 3;
215 tau = &particles[idx];
218 if (hardME) hardME->calculateRho(idx, particles);
219 vector<HelicityParticle> children = createChildren(*tau);
220 if (children.size() == 0)
return false;
223 bool accepted =
false;
226 isotropicDecay(children);
227 double decayWeight = decayME->decayWeight(children);
228 double decayWeightMax = decayME->decayWeightMax(children);
229 accepted = (rndmPtr->flat() < decayWeight / decayWeightMax);
230 if (decayWeight > decayWeightMax)
231 infoPtr->errorMsg(
"Warning in TauDecays::decay: maximum "
232 "decay weight exceeded in tau decay");
233 if (tries > NTRYDECAY) {
234 infoPtr->errorMsg(
"Warning in TauDecays::decay: maximum "
235 "number of decay attempts exceeded");
240 writeDecay(event,children);
244 idx = (idx == 2) ? 3 : 2;
246 decayME->calculateD(children);
248 tau->D = children[0].D;
250 tau = &particles[idx];
252 hardME->calculateRho(idx, particles);
256 children = createChildren(*tau);
257 if (children.size() == 0)
return false;
261 isotropicDecay(children);
262 double decayWeight = decayME->decayWeight(children);
263 double decayWeightMax = decayME->decayWeightMax(children);
264 accepted = (rndmPtr->flat() < decayWeight / decayWeightMax);
265 if (decayWeight > decayWeightMax)
266 infoPtr->errorMsg(
"Warning in TauDecays::decay: maximum "
267 "decay weight exceeded in correlated tau decay");
268 if (tries > NTRYDECAY) {
269 infoPtr->errorMsg(
"Warning in TauDecays::decay: maximum "
270 "number of decay attempts exceeded");
275 writeDecay(event,children);
288 bool TauDecays::internalMechanism(
Event&) {
294 if (abs(mediator.id()) == 22 || abs(mediator.id()) == 23 ||
295 abs(mediator.id()) == 32) {
297 if (abs(in1.id()) <= 18 && abs(in2.id()) <= 18 && in1.daughter2() == 0 &&
298 in2.daughter2() == 0 && in1.daughter1() == in2.daughter1()) {
299 particles.push_back(mediator);
300 hardME = hmeTwoFermions2GammaZ2TwoFermions.initChannel(particles);
302 }
else known =
false;
305 }
else if (abs(mediator.id()) == 24 || abs(mediator.id()) == 34) {
307 if (abs(in1.id()) <= 18 && abs(in2.id()) <= 18 && in1.daughter2() == 0 &&
308 in2.daughter2() == 0 && in1.daughter1() == in2.daughter1()) {
309 particles.push_back(mediator);
310 hardME = hmeTwoFermions2W2TwoFermions.initChannel(particles);
312 }
else known =
false;
315 }
else if (abs(mediator.id()) == 25 || abs(mediator.id()) == 35 ||
316 abs(mediator.id()) == 36 || abs(mediator.id()) == 37) {
317 particles[1] = mediator;
318 hardME = hmeHiggs2TwoFermions.initChannel(particles);
321 }
else if ((abs(mediator.id()) == 411 || abs(mediator.id()) == 431
322 || abs(mediator.id()) == 511 || abs(mediator.id()) == 521
323 || abs(mediator.id()) == 531 || abs(mediator.id()) == 541
324 || (abs(mediator.id()) > 5100 && abs(mediator.id()) < 5600) )
325 && abs(out2.id()) == 16) {
326 int idBmother = (mediator.id() > 0) ? -5 : 5;
327 if (abs(mediator.id()) > 5100) idBmother = -idBmother;
328 particles[0] = HelicityParticle( idBmother, 0, 0, 0, 0, 0, 0, 0,
329 0., 0., 0., 0., 0., 0., particleDataPtr);
330 particles[1] = HelicityParticle( -idBmother, 0, 0, 0, 0, 0, 0, 0,
331 0., 0., 0., 0., 0., 0., particleDataPtr);
332 particles[0].index(-1);
333 particles[1].index(-1);
336 if (mediator.daughter1() + 2 == mediator.daughter2()) {
337 particles[0].p(mediator.p());
338 particles[1].direction = 1;
339 particles[1].id(-particles[1].
id());
340 particles[1].p(particles[0].p() - particles[2].p() - particles[3].p());
345 particles[0].p(mediator.p()/2);
346 particles[1].p(mediator.p()/2);
348 hardME = hmeTwoFermions2W2TwoFermions.initChannel(particles);
351 }
else known =
false;
361 bool TauDecays::externalMechanism(
Event &event) {
367 if (tauExt == 0) correlated =
false;
369 double spinup = particles[2].pol();
370 if (abs(spinup) > 1.001) spinup =
event[particles[2].iTopCopyId()].pol();
371 if (abs(spinup) > 1.001) known =
false;
373 particles[2].rho[0][0] = (1 - spinup) / 2;
374 particles[2].rho[1][1] = (1 + spinup) / 2;
378 }
else if (tauExt == 1) {
379 double spinup = mediator.pol();
380 if (abs(spinup) > 1.001) spinup =
event[mediator.iTopCopyId()].pol();
381 if (abs(spinup) > 1.001) spinup = 0;
382 if (mediator.rho.size() > 1) {
383 mediator.rho[0][0] = (1 - spinup) / mediator.spinStates();
384 mediator.rho[1][1] = (1 + spinup) / mediator.spinStates();
386 particles[1] = mediator;
387 if (abs(mediator.id()) == 22)
388 hardME = hmeGamma2TwoFermions.initChannel(particles);
389 else if (abs(mediator.id()) == 23 || abs(mediator.id()) == 32)
390 hardME = hmeZ2TwoFermions.initChannel(particles);
391 else if (abs(mediator.id()) == 24 || abs(mediator.id()) == 34)
392 hardME = hmeZ2TwoFermions.initChannel(particles);
393 else if (abs(mediator.id()) == 25 || abs(mediator.id()) == 35 ||
394 abs(mediator.id()) == 36 || abs(mediator.id()) == 37)
395 hardME = hmeHiggs2TwoFermions.initChannel(particles);
399 }
else known =
false;
410 vector<HelicityParticle> TauDecays::createChildren(HelicityParticle parent) {
414 vector<HelicityParticle> children;
417 parent.direction = -1;
420 ParticleDataEntry decayData = parent.particleDataEntry();
423 if (!decayData.preparePick(parent.id()))
return children;
426 bool decayed =
false;
428 while (!decayed && decayTries < NTRYCHANNEL) {
431 DecayChannel decayChannel = decayData.pickChannel();
432 meMode = decayChannel.meMode();
433 int decayMult = decayChannel.multiplicity();
436 bool allowed =
false;
437 int channelTries = 0;
438 while (!allowed && channelTries < NTRYCHANNEL) {
440 children.push_back(parent);
441 for (
int i = 0; i < decayMult; ++i) {
443 int childId = decayChannel.product(i);
445 if (parent.id() < 0 && particleDataPtr->hasAnti(childId))
448 double childMass = particleDataPtr->mSel(childId);
450 children.push_back( HelicityParticle(childId, 91, parent.index(), 0, 0,
451 0, 0, 0, 0., 0., 0., 0.,childMass, 0., particleDataPtr) );
456 double massDiff = parent.m();
457 for (
int i = 0; i < decayMult; ++i) massDiff -= children[i].m();
472 if (parent.idAbs() == 13 && children.size() == 4 && meMode == 22)
473 swap(children[1], children[3]);
477 if (children.size() == 3) {
479 decayME = hmeTau2Meson.initChannel(children);
480 else decayME = hmeTau2PhaseSpace.initChannel(children);
484 else if (children.size() == 4) {
486 if (meMode == 1531 || (parent.idAbs() == 13 && meMode == 22))
487 decayME = hmeTau2TwoLeptons.initChannel(children);
489 else if (meMode == 1532)
490 decayME = hmeTau2TwoMesonsViaVector.initChannel(children);
492 else if (meMode == 1533)
493 decayME = hmeTau2TwoMesonsViaVectorScalar.initChannel(children);
495 else decayME = hmeTau2PhaseSpace.initChannel(children);
499 else if (children.size() == 5) {
502 decayME = hmeTau2ThreePions.initChannel(children);
504 else if (meMode == 1542)
505 decayME = hmeTau2ThreeMesonsWithKaons.initChannel(children);
507 else if (meMode == 1543)
508 decayME = hmeTau2ThreeMesonsGeneric.initChannel(children);
510 else if (meMode == 1544)
511 decayME = hmeTau2TwoPionsGamma.initChannel(children);
513 else decayME = hmeTau2PhaseSpace.initChannel(children);
517 else if (children.size() == 6) {
520 decayME = hmeTau2FourPions.initChannel(children);
522 else decayME = hmeTau2PhaseSpace.initChannel(children);
526 else if (children.size() == 7) {
529 decayME = hmeTau2FivePions.initChannel(children);
531 else decayME = hmeTau2PhaseSpace.initChannel(children);
535 else decayME = hmeTau2PhaseSpace.initChannel(children);
549 void TauDecays::isotropicDecay(vector<HelicityParticle>& children) {
552 int decayMult = children.size() - 1;
553 double m0 = children[0].m();
554 double mSum = children[1].m();
555 for (
int i = 2; i <= decayMult; ++i) mSum += children[i].m();
556 double mDiff = m0 - mSum;
560 for (
int i = 0; i <= decayMult; ++i) mInv.push_back( children[i].m());
564 double wtPSmax = 1. / WTCORRECTION[decayMult];
565 double mMax = mDiff + children[decayMult].m();
567 for (
int i = decayMult - 1; i > 0; --i) {
568 mMax += children[i].m();
569 mMin += children[i+1].m();
570 double mNow = children[i].m();
571 wtPSmax *= 0.5 * sqrtpos( (mMax - mMin - mNow) * (mMax + mMin + mNow)
572 * (mMax + mMin - mNow) * (mMax - mMin + mNow) ) / mMax;
576 vector<double> rndmOrd;
582 rndmOrd.push_back(1.);
583 for (
int i = 1; i < decayMult - 1; ++i) {
584 double random = rndmPtr->flat();
585 rndmOrd.push_back(random);
586 for (
int j = i - 1; j > 0; --j) {
587 if (random > rndmOrd[j]) swap( rndmOrd[j], rndmOrd[j+1] );
591 rndmOrd.push_back(0.);
594 for (
int i = decayMult - 1; i > 0; --i) {
595 mInv[i] = mInv[i+1] + children[i].m()
596 + (rndmOrd[i-1] - rndmOrd[i]) * mDiff;
597 wtPS *= 0.5 * sqrtpos( (mInv[i] - mInv[i+1] - children[i].m())
598 * (mInv[i] + mInv[i+1] + children[i].m())
599 * (mInv[i] + mInv[i+1] - children[i].m())
600 * (mInv[i] - mInv[i+1] + children[i].m()) ) / mInv[i];
604 }
while ( wtPS < rndmPtr->flat() * wtPSmax );
607 vector<Vec4> pInv(decayMult + 1);
608 for (
int i = 1; i < decayMult; ++i) {
610 pair<Vec4, Vec4> ps =
611 rndmPtr->phaseSpace2(mInv[i], mInv[i+1], children[i].m());
612 pInv[i+1].p(ps.first);
613 children[i].p(ps.second);
617 children[decayMult].p( pInv[decayMult] );
618 for (
int iFrame = decayMult - 1; iFrame > 1; --iFrame)
619 for (
int i = iFrame; i <= decayMult; ++i)
620 children[i].bst( pInv[iFrame], mInv[iFrame]);
623 pInv[1].p( children[0].p() );
624 for (
int i = 1; i <= decayMult; ++i) children[i].bst( pInv[1], mInv[1] );
636 void TauDecays::writeDecay(
Event& event, vector<HelicityParticle>& children) {
639 int decayMult = children.size() - 1;
640 Vec4 decayVertex = children[0].vDec();
641 for (
int i = 1; i <= decayMult; i++) {
643 children[i].tau(children[i].tau0() * rndmPtr->exp());
645 children[i].vProd(decayVertex);
647 children[i].index(event.append(children[i]));
651 event[children[0].index()].statusNeg();
652 event[children[0].index()].daughters(children[1].index(),
653 children[decayMult].index());