13#include "EvtGenBase/EvtCPUtil.hh"
14#include "EvtGenBase/EvtPDL.hh"
15#include "EvtGenBase/EvtTensor4C.hh"
17#include <generators/evtgen/EvtGenModelRegister.h>
18#include "generators/evtgen/models/EvtB0toKsKK.h"
29 EvtB0toKsKK::~EvtB0toKsKK() {}
50 if ((getParentId() != EvtPDL::getId(
"B0") &&
51 getParentId() != EvtPDL::getId(
"anti-B0")) ||
52 (getDaug(0) != EvtPDL::getId(
"K_S0")) ||
53 (getDaug(1) != EvtPDL::getId(
"K+") &&
54 getDaug(1) != EvtPDL::getId(
"K-")) ||
55 (getDaug(2) != EvtPDL::getId(
"K+") &&
56 getDaug(2) != EvtPDL::getId(
"K-"))) {
57 std::cout <<
"ERROR: Invalid decay" << std::endl;
58 std::cout <<
"USAGE: K_S0 K+ K-" << std::endl;
63 (1.0 + getArg(3)) * EvtComplex(getArg(1) * cos(getArg(2) * M_PI / 180.0),
64 getArg(1) * sin(getArg(2) * M_PI / 180.0));
66 (1.0 + getArg(7)) * EvtComplex(getArg(5) * cos(getArg(6) * M_PI / 180.0),
67 getArg(5) * sin(getArg(6) * M_PI / 180.0));
69 (1.0 + getArg(11)) * EvtComplex(getArg(9) * cos(getArg(10) * M_PI / 180.0),
70 getArg(9) * sin(getArg(10) * M_PI / 180.0));
72 (1.0 + getArg(15)) * EvtComplex(getArg(13) * cos(getArg(14) * M_PI / 180.0),
73 getArg(13) * sin(getArg(14) * M_PI / 180.0));
75 (1.0 + getArg(19)) * EvtComplex(getArg(17) * cos(getArg(18) * M_PI / 180.0),
76 getArg(17) * sin(getArg(18) * M_PI / 180.0));
78 (1.0 + getArg(24)) * EvtComplex(getArg(22) * cos(getArg(23) * M_PI / 180.0),
79 getArg(22) * sin(getArg(23) * M_PI / 180.0));
81 (1.0 + getArg(29)) * EvtComplex(getArg(27) * cos(getArg(28) * M_PI / 180.0),
82 getArg(27) * sin(getArg(28) * M_PI / 180.0));
85 (1.0 - getArg(3)) * EvtComplex(getArg(1) * cos(getArg(2) * M_PI / 180.0),
86 getArg(1) * sin(getArg(2) * M_PI / 180.0));
88 (1.0 - getArg(7)) * EvtComplex(getArg(5) * cos(getArg(6) * M_PI / 180.0),
89 getArg(5) * sin(getArg(6) * M_PI / 180.0));
91 (1.0 - getArg(11)) * EvtComplex(getArg(9) * cos(getArg(10) * M_PI / 180.0),
92 getArg(9) * sin(getArg(10) * M_PI / 180.0));
94 (1.0 - getArg(15)) * EvtComplex(getArg(13) * cos(getArg(14) * M_PI / 180.0),
95 getArg(13) * sin(getArg(14) * M_PI / 180.0));
97 (1.0 - getArg(19)) * EvtComplex(getArg(17) * cos(getArg(18) * M_PI / 180.0),
98 getArg(17) * sin(getArg(18) * M_PI / 180.0));
100 (1.0 - getArg(24)) * EvtComplex(getArg(22) * cos(getArg(23) * M_PI / 180.0),
101 getArg(22) * sin(getArg(23) * M_PI / 180.0));
103 (1.0 - getArg(29)) * EvtComplex(getArg(27) * cos(getArg(28) * M_PI / 180.0),
104 getArg(27) * sin(getArg(28) * M_PI / 180.0));
110 std::cout << setiosflags(std::ios::left) << std::endl
111 <<
"B0 Channel " << std::setw(20) <<
"Relative amplitude" << std::setw(20) <<
"Relative phase" << std::endl
112 <<
"f0ks " << std::setw(20) << real(
a_f0ks_) << std::setw(20) << imag(
a_f0ks_) << std::endl
113 <<
"phiks " << std::setw(20) << real(
a_phiks_) << std::setw(20) << imag(
a_phiks_) << std::endl
114 <<
"fxks " << std::setw(20) << real(
a_fxks_) << std::setw(20) << imag(
a_fxks_) << std::endl
115 <<
"chic0ks " << std::setw(20) << real(
a_chic0ks_) << std::setw(20) << imag(
a_chic0ks_) << std::endl
116 <<
"kpkmnr " << std::setw(20) << real(
a_kpkmnr_) << std::setw(20) << imag(
a_kpkmnr_) << std::endl
117 <<
"kskpnr " << std::setw(20) << real(
a_kskpnr_) << std::setw(20) << imag(
a_kskpnr_) << std::endl
118 <<
"kskmnr " << std::setw(20) << real(
a_kskmnr_) << std::setw(20) << imag(
a_kskmnr_) << std::endl
121 std::cout << setiosflags(std::ios::left) << std::endl
122 <<
"B0B Channel " << std::setw(20) <<
"Relative amplitude" << std::setw(20) <<
"Relative phase" << std::endl
123 <<
"f0ks " << std::setw(20) << real(
abar_f0ks_) << std::setw(20) << imag(
abar_f0ks_) << std::endl
125 <<
"fxks " << std::setw(20) << real(
abar_fxks_) << std::setw(20) << imag(
abar_fxks_) << std::endl
143 setProbMax(1100000.0);
149 static EvtId B0 = EvtPDL::getId(
"B0");
150 static EvtId B0B = EvtPDL::getId(
"anti-B0");
157 EvtCPUtil::getInstance()->OtherB(
p, t, other_b, 0.5);
162 p->initializePhaseSpace(getNDaug(), getDaugs());
163 EvtVector4R p4ks, p4kp, p4km;
164 if (
p->getId() == B0) {
165 p4ks =
p->getDaug(0)->getP4();
166 p4kp =
p->getDaug(1)->getP4();
167 p4km =
p->getDaug(2)->getP4();
174 p4ks =
p->getDaug(0)->getP4();
175 p4kp =
p->getDaug(2)->getP4();
176 p4km =
p->getDaug(1)->getP4();
217 (1.0 + getArg(3)) * EvtComplex(getArg(1) * cos(getArg(2) * M_PI / 180.0),
218 getArg(1) * sin(getArg(2) * M_PI / 180.0));
220 (1.0 + getArg(7)) * EvtComplex(getArg(5) * cos(getArg(6) * M_PI / 180.0),
221 getArg(5) * sin(getArg(6) * M_PI / 180.0));
223 (1.0 + getArg(11)) * EvtComplex(getArg(9) * cos(getArg(10) * M_PI / 180.0),
224 getArg(9) * sin(getArg(10) * M_PI / 180.0));
226 (1.0 + getArg(15)) * EvtComplex(getArg(13) * cos(getArg(14) * M_PI / 180.0),
227 getArg(13) * sin(getArg(14) * M_PI / 180.0));
229 (1.0 + getArg(19)) * EvtComplex(getArg(17) * cos(getArg(18) * M_PI / 180.0),
230 getArg(17) * sin(getArg(18) * M_PI / 180.0));
232 (1.0 + getArg(24)) * EvtComplex(getArg(22) * cos(getArg(23) * M_PI / 180.0),
233 getArg(22) * sin(getArg(23) * M_PI / 180.0));
235 (1.0 + getArg(29)) * EvtComplex(getArg(27) * cos(getArg(28) * M_PI / 180.0),
236 getArg(27) * sin(getArg(28) * M_PI / 180.0));
239 (1.0 - getArg(3)) * EvtComplex(getArg(1) * cos(getArg(2) * M_PI / 180.0),
240 getArg(1) * sin(getArg(2) * M_PI / 180.0));
242 (1.0 - getArg(7)) * EvtComplex(getArg(5) * cos(getArg(6) * M_PI / 180.0),
243 getArg(5) * sin(getArg(6) * M_PI / 180.0));
245 (1.0 - getArg(11)) * EvtComplex(getArg(9) * cos(getArg(10) * M_PI / 180.0),
246 getArg(9) * sin(getArg(10) * M_PI / 180.0));
248 (1.0 - getArg(15)) * EvtComplex(getArg(13) * cos(getArg(14) * M_PI / 180.0),
249 getArg(13) * sin(getArg(14) * M_PI / 180.0));
251 (1.0 - getArg(19)) * EvtComplex(getArg(17) * cos(getArg(18) * M_PI / 180.0),
252 getArg(17) * sin(getArg(18) * M_PI / 180.0));
254 (1.0 - getArg(24)) * EvtComplex(getArg(22) * cos(getArg(23) * M_PI / 180.0),
255 getArg(22) * sin(getArg(23) * M_PI / 180.0));
257 (1.0 - getArg(29)) * EvtComplex(getArg(27) * cos(getArg(28) * M_PI / 180.0),
258 getArg(27) * sin(getArg(28) * M_PI / 180.0));
261 const double pCP_f0ks = getArg(4) * M_PI / 180.0;
262 const double pCP_phiks = getArg(8) * M_PI / 180.0;
263 const double pCP_fxks = getArg(12) * M_PI / 180.0;
264 const double pCP_chic0ks = getArg(16) * M_PI / 180.0;
265 const double pCP_kpkmnr = getArg(20) * M_PI / 180.0;
266 const double pCP_kskpnr = getArg(25) * M_PI / 180.0;
267 const double pCP_kskmnr = getArg(30) * M_PI / 180.0;
271 EvtComplex(cos(+2.0 * pCP_f0ks), sin(+2.0 * pCP_f0ks));
273 EvtComplex(cos(+2.0 * pCP_phiks), sin(+2.0 * pCP_phiks));
275 EvtComplex(cos(+2.0 * pCP_fxks), sin(+2.0 * pCP_fxks));
277 EvtComplex(cos(+2.0 * pCP_chic0ks), sin(+2.0 * pCP_chic0ks));
279 EvtComplex(cos(+2.0 * pCP_kpkmnr), sin(+2.0 * pCP_kpkmnr));
281 EvtComplex(cos(+2.0 * pCP_kskpnr), sin(+2.0 * pCP_kskpnr));
283 EvtComplex(cos(+2.0 * pCP_kskmnr), sin(+2.0 * pCP_kskmnr));
285 if (other_b == B0B) {
287 EvtComplex(cos(-2.0 * pCP_f0ks), sin(-2.0 * pCP_f0ks));
289 EvtComplex(cos(-2.0 * pCP_phiks), sin(-2.0 * pCP_phiks));
291 EvtComplex(cos(-2.0 * pCP_fxks), sin(-2.0 * pCP_fxks));
293 EvtComplex(cos(-2.0 * pCP_chic0ks), sin(-2.0 * pCP_chic0ks));
295 EvtComplex(cos(-2.0 * pCP_kpkmnr), sin(-2.0 * pCP_kpkmnr));
297 EvtComplex(cos(-2.0 * pCP_kskpnr), sin(-2.0 * pCP_kskpnr));
299 EvtComplex(cos(-2.0 * pCP_kskmnr), sin(-2.0 * pCP_kskmnr));
303 EvtComplex Amp_f0ks =
A_f0ks(p4ks, p4kp, p4km);
304 EvtComplex Amp_phiks =
A_phiks(p4ks, p4kp, p4km);
305 EvtComplex Amp_fxks =
A_fxks(p4ks, p4kp, p4km);
306 EvtComplex Amp_chic0ks =
A_chic0ks(p4ks, p4kp, p4km);
311 EvtComplex Ampbar_f0ks =
A_f0ks(p4ks, p4km, p4kp);
312 EvtComplex Ampbar_phiks =
A_phiks(p4ks, p4km, p4kp);
313 EvtComplex Ampbar_fxks =
A_fxks(p4ks, p4km, p4kp);
314 EvtComplex Ampbar_chic0ks =
A_chic0ks(p4ks, p4km, p4kp);
319 const EvtComplex A_B0toKsKK =
328 const EvtComplex Abar_B0toKsKK =
338 const double dm = getArg(0);
341 if (other_b == B0B) {
343 A_B0toKsKK * cos(dm * t / (2.0 * EvtConst::c)) +
344 EvtComplex(0.0, 1.0) * Abar_B0toKsKK * sin(dm * t / (2.0 * EvtConst::c));
349 EvtComplex(0.0, 1.0) * sin(dm * t / (2.0 * EvtConst::c)) +
350 Abar_B0toKsKK * cos(dm * t / (2.0 * EvtConst::c));
353 if (abs2(amp) > 50000.0)
361 const EvtVector4R& p4c)
363 const double s = (p4a + p4b + p4c) * (p4a + p4b + p4c);
365 const EvtVector4R
umu = (p4a + p4b + p4c) /
sqrt(s);
373 const double ma = p4a.mass();
374 const double mb = p4b.mass();
375 const double mres = (p4a + p4b).mass();
376 const double ma2 = ma * ma;
377 const double mb2 = mb * mb;
378 const double mres2 = mres * mres;
382 (
sqrt(mres2 - ((ma + mb) * (ma + mb))) *
sqrt(mres2 - ((ma - mb) * (ma - mb))));
384 const EvtVector4R
Smu =
385 N1 * ((p4a - p4b) - ((p4a + p4b) * (ma2 - mb2) / mres2));
391 const EvtVector4R& p4c)
393 const double mc = p4c.mass();
394 const double mres = (p4a + p4b).mass();
395 const double mc2 = mc * mc;
396 const double mres2 = mres * mres;
397 const double s = (p4a + p4b + p4c) * (p4a + p4b + p4c);
401 (
sqrt(s - ((mres + mc) * (mres + mc))) *
sqrt(s - ((mres - mc) * (mres - mc))));
403 const EvtVector4R
Lmu =
404 N2 * ((p4c - (p4a + p4b)) - ((p4c + (p4a + p4b)) * (mc2 - mres2) / s));
410 const EvtVector4R& p4b,
411 const EvtVector4R& p4c)
413 const double s = (p4a + p4b + p4c) * (p4a + p4b + p4c);
415 const EvtVector4R
umu = (p4a + p4b + p4c) /
sqrt(s);
418 EvtTensor4C::g() - EvtGenFunctions::directProd(
umu,
umu);
424 const EvtVector4R& p4c)
426 const double mres = (p4a + p4b).mass();
427 const EvtVector4R wmu = (p4a + p4b) / mres;
429 const EvtTensor4C
Tmunu =
431 (EvtGenFunctions::directProd(
Smu(p4a, p4b, p4c),
Smu(p4a, p4b, p4c)) +
432 (1.0 / 3.0 * (EvtTensor4C::g() - EvtGenFunctions::directProd(wmu, wmu))));
438 const EvtTensor4C& t2)
441 for (
unsigned int i = 0; i < 4; i++)
442 for (
unsigned int j = 0; j < 4; j++) {
443 const EvtComplex element =
444 t1.get(i, 0) * t2.get(0, j) +
445 t1.get(i, 1) * t2.get(1, j) +
446 t1.get(i, 2) * t2.get(2, j) +
447 t1.get(i, 3) * t2.get(3, j);
448 t.set(i, j, element);
455 const EvtTensor4C tmp =
Multiply(t, EvtTensor4C::g());
456 const EvtTensor4C t_raised =
Multiply(EvtTensor4C::g(), tmp);
463 vector.set(1, -vector.get(1));
464 vector.set(2, -vector.get(2));
465 vector.set(3, -vector.get(3));
469 const EvtVector4R& p4c)
471 const EvtTensor4C
Mmunu =
473 (EvtGenFunctions::directProd(
Lmu(p4a, p4b, p4c),
Lmu(p4a, p4b, p4c)) +
482 const double d = 1.0;
484 const double z =
q *
q * d * d;
489 bwbf =
sqrt(2.0 * z / (1.0 + z));
491 bwbf =
sqrt(13.0 * z * z / (((z - 3.0) * (z - 3.0)) + (9.0 * z)));
493 std::cout <<
"ERROR: BWBF not implemented for L>2" << std::endl;
500 const unsigned int& L)
503 const double d = 1.0;
505 const double z =
q *
q * d * d;
506 const double z0 = q0 * q0 * d * d;
511 bwbf =
sqrt((1.0 + z0) / (1.0 + z));
513 bwbf =
sqrt((((z0 - 3.0) * (z0 - 3.0)) + (9.0 * z0)) / (((z - 3.0) * (z - 3.0)) + (9.0 * z)));
515 std::cout <<
"ERROR: BWBF not implemented for L>2" << std::endl;
522 const double& Gamma0,
523 const double& q,
const double& q0,
524 const unsigned int& L)
526 const double s = m * m;
527 const double s0 = m0 * m0;
530 Gamma0 * m0 / m * pow(
q / q0, (2 * L) + 1) *
533 const EvtComplex
BreitWigner = 1.0 / EvtComplex(s0 - s, -m0 * Gamma);
539 const EvtVector4R& boost)
541 const double bx = boost.get(1) / boost.get(0);
542 const double by = boost.get(2) / boost.get(0);
543 const double bz = boost.get(3) / boost.get(0);
545 const double bx2 = bx * bx;
546 const double by2 = by * by;
547 const double bz2 = bz * bz;
549 const double b2 = bx2 + by2 + bz2;
554 const double gamma = 1.0 /
sqrt(1 - b2);
556 const double gb2 = (gamma - 1.0) / b2;
558 const double gb2xy = gb2 * bx * by;
559 const double gb2xz = gb2 * bx * bz;
560 const double gb2yz = gb2 * by * bz;
562 const double gbx = gamma * bx;
563 const double gby = gamma * by;
564 const double gbz = gamma * bz;
567 (gamma * p4.get(0)) - (gbx * p4.get(1)) - (gby * p4.get(2)) - (gbz * p4.get(3));
569 (-gbx * p4.get(0)) + ((1.0 + (gb2 * bx2)) * p4.get(1)) +
570 (gb2xy * p4.get(2)) + (gb2xz * p4.get(3));
572 (-gby * p4.get(0)) + (gb2xy * p4.get(1)) +
573 ((1.0 + (gb2 * by2)) * p4.get(2)) + (gb2yz * p4.get(3));
575 (-gbz * p4.get(0)) + (gb2xz * p4.get(1)) +
576 (gb2yz * p4.get(2)) + ((1.0 + (gb2 * bz2)) * p4.get(3));
578 const EvtVector4R p4_b(e_b, x_b, y_b, z_b);
585 const double sab = mab * mab;
586 const double M_p_mc2 = (M + mc) * (M + mc);
587 const double M_m_mc2 = (M - mc) * (M - mc);
589 const double p =
sqrt((sab - M_p_mc2) * (sab - M_m_mc2)) / (2.0 * mab);
596 const double mab2 = mab * mab;
597 const double ma_p_mb2 = (ma + mb) * (ma + mb);
598 const double ma_m_mb2 = (ma - mb) * (ma - mb);
600 const double q =
sqrt((mab2 - ma_p_mb2) * (mab2 - ma_m_mb2)) / (2.0 * mab);
607 const double k2 = 1.0 - (4.0 * m_h * m_h / s);
610 k = EvtComplex(0.0,
sqrt(fabs(k2)));
619 const double g_pipi = 0.165;
620 const double g_kk = 4.21 * g_pipi;
622 static EvtId pip = EvtPDL::getId(
"pi+");
623 static EvtId kp = EvtPDL::getId(
"K+");
625 const double s = m * m;
626 const double s0 = m0 * m0;
628 const EvtComplex rho_pipi = 2.0 *
Flatte_k(s, EvtPDL::getMeanMass(pip)) / m;
629 const EvtComplex rho_kk = 2.0 *
Flatte_k(s, EvtPDL::getMeanMass(kp)) / m;
631 const EvtComplex Gamma =
632 EvtComplex(0.0, 1.0) * ((g_pipi * rho_pipi) + (g_kk * rho_kk));
634 const EvtComplex
Flatte = 1.0 / (s0 - s - Gamma);
640 const EvtVector4R& p4kp,
641 const EvtVector4R& p4km)
643 static EvtId f0 = EvtPDL::getId(
"f_0");
645 const double f0_m = EvtPDL::getMeanMass(f0);
647 const EvtVector4R p4kpkm = p4kp + p4km;
648 const double mkpkm = p4kpkm.mass();
651 const EvtComplex H_f0ks = 1.0;
654 const EvtComplex D_f0ks = 1.0;
657 const EvtComplex L_f0ks =
Flatte(mkpkm, f0_m);
660 const EvtComplex
A_f0ks = D_f0ks * H_f0ks * L_f0ks;
666 double s13_min(
const double& s12,
const double& M,
667 const double& m1,
const double& m2,
const double& m3)
669 const double E2star = (s12 - (m2 * m2) + (m1 * m1)) / 2.0 /
sqrt(s12);
670 const double E3star = (M * M - s12 - (m3 * m3)) / 2.0 /
sqrt(s12);
672 const double s23_min =
673 ((E2star + E3star) * (E2star + E3star)) -
674 ((
sqrt((E2star * E2star) - (m1 * m1)) +
sqrt((E3star * E3star) - (m3 * m3))) *
675 (
sqrt((E2star * E2star) - (m1 * m1)) +
sqrt((E3star * E3star) - (m3 * m3))));
681 double s13_max(
const double& s12,
const double& M,
682 const double& m1,
const double& m2,
const double& m3)
684 const double E2star = (s12 - (m2 * m2) + (m1 * m1)) / 2.0 /
sqrt(s12);
685 const double E3star = (M * M - s12 - (m3 * m3)) / 2.0 /
sqrt(s12);
687 const double s23_max =
688 ((E2star + E3star) * (E2star + E3star)) -
689 ((
sqrt((E2star * E2star) - (m1 * m1)) -
sqrt((E3star * E3star) - (m3 * m3))) *
690 (
sqrt((E2star * E2star) - (m1 * m1)) -
sqrt((E3star * E3star) - (m3 * m3))));
696 const EvtVector4R& p4kp,
697 const EvtVector4R& p4km)
699 static EvtId phi = EvtPDL::getId(
"phi");
701 const double phi_m = EvtPDL::getMeanMass(phi);
702 const double phi_w = EvtPDL::getWidth(phi);
704 const EvtVector4R p4kpkm = p4kp + p4km;
705 const double mkpkm = p4kpkm.mass();
708 const EvtVector4R S_mu =
Smu(p4kp, p4km, p4ks);
709 const EvtVector4R L_mu =
Lmu(p4kp, p4km, p4ks);
711 const EvtTensor4C g_munu_tilde =
gmunu_tilde(p4kp, p4km, p4ks);
712 const EvtTensor4C g__munu_tilde =
RaiseIndices(g_munu_tilde);
714 EvtComplex H_phiks = 0.0;
715 for (
unsigned int i = 0; i < 4; i++)
716 for (
unsigned int j = 0; j < 4; j++)
717 H_phiks += g__munu_tilde.get(i, j) * S_mu.get(i) * L_mu.get(j);
720 const EvtVector4R p4kp_kpkm =
Boost(p4kp, p4kpkm);
722 const EvtComplex D_phiks =
BWBF(p4kp_kpkm.d3mag(), 1);
725 const double q0_phi =
q(phi_m, p4kp.mass(), p4km.mass());
726 const EvtComplex L_phiks =
727 BreitWigner(mkpkm, phi_m, phi_w, p4kp_kpkm.d3mag(), q0_phi, 1);
730 const EvtComplex
A_phiks = D_phiks * H_phiks * L_phiks;
771 const EvtVector4R& p4kp,
772 const EvtVector4R& p4km)
774 static EvtId fx = EvtPDL::getId(
"f_0(1500)");
776 const double fx_m = EvtPDL::getMeanMass(fx);
777 const double fx_w = EvtPDL::getWidth(fx);
779 const EvtVector4R p4kpkm = p4kp + p4km;
780 const double mkpkm = p4kpkm.mass();
783 const EvtComplex H_fxks = 1.0;
786 const EvtComplex D_fxks = 1.0;
789 const EvtVector4R p4kp_kpkm =
Boost(p4kp, p4kpkm);
790 const double q0_fx =
q(fx_m, p4kp.mass(), p4km.mass());
791 const EvtComplex L_fxks =
792 BreitWigner(mkpkm, fx_m, fx_w, p4kp_kpkm.d3mag(), q0_fx, 0);
795 const EvtComplex
A_fxks = D_fxks * H_fxks * L_fxks;
801 const EvtVector4R& p4kp,
802 const EvtVector4R& p4km)
804 static EvtId chic0 = EvtPDL::getId(
"chi_c0");
806 const double chic0_m = EvtPDL::getMeanMass(chic0);
807 const double chic0_w = EvtPDL::getWidth(chic0);
809 const EvtVector4R p4kpkm = p4kp + p4km;
810 const double mkpkm = p4kpkm.mass();
813 const EvtComplex H_chic0ks = 1.0;
816 const EvtComplex D_chic0ks = 1.0;
819 const EvtVector4R p4kp_kpkm =
Boost(p4kp, p4kpkm);
820 const double q0_chic0 =
q(chic0_m, p4kp.mass(), p4km.mass());
821 const EvtComplex L_chic0ks =
822 BreitWigner(mkpkm, chic0_m, chic0_w, p4kp_kpkm.d3mag(), q0_chic0, 0);
825 const EvtComplex
A_chic0ks = D_chic0ks * H_chic0ks * L_chic0ks;
831 const EvtVector4R& p4k2,
832 const double& alpha_kk)
834 const EvtVector4R p4kk = p4k1 + p4k2;
835 const double m2kk = p4kk.mass2();
837 const EvtComplex
A_kknr = exp(-alpha_kk * m2kk);
Register Decay model EvtB0toKsKK.
EvtComplex abar_kskpnr_
Variable member abar_kskpnr_
double alpha_kskmnr
Variable member alpha_kskmnr.
EvtComplex a_chic0ks_
Variable member a_chic0ks_.
EvtComplex a_fxks_
Variable member a_fxks_
std::ofstream debugfile_
debuging stream
EvtComplex a_f0ks_
<Variable names for form factors
double alpha_kpkmnr
Variable member alpha_kpkmnr.
EvtComplex abar_fxks_
Variable member abar_fxks_
EvtComplex abar_f0ks_
Variable member abar_f0ks_
double alpha_kskpnr
Variable member alpha_kskpnr.
EvtComplex abar_chic0ks_
Variable member abar_chic0ks_.
EvtComplex a_kskmnr_
Variable member a_kskmnr_.
EvtComplex abar_kskmnr_
Variable member abar_kskmnr_
EvtComplex a_kpkmnr_
Variable member a_kpkmnr_.
EvtComplex abar_phiks_
Variable member abar_phiks_.
EvtComplex a_kskpnr_
Variable member a_kskpnr_.
EvtComplex a_phiks_
Variable member a_phiks_
EvtComplex abar_kpkmnr_
Variable member abar_kpkmnr_
void init()
Initialize standard stream objects
EvtVector4R Smu(const EvtVector4R &p4a, const EvtVector4R &p4b, const EvtVector4R &p4c)
Function 4Vector Smu.
EvtComplex A_f0ks(const EvtVector4R &p4ks, const EvtVector4R &p4kp, const EvtVector4R &p4km)
A_f0ks is amplitude of f0.
EvtComplex A_phiks(const EvtVector4R &p4ks, const EvtVector4R &p4kp, const EvtVector4R &p4km)
A_phiks is amplitude of phi.
EvtVector4R umu(const EvtVector4R &p4a, const EvtVector4R &p4b, const EvtVector4R &p4c)
Function 4Vector umu.
EvtComplex Flatte(const double &m, const double &m0)
Constant Flatte.
double s13_max(const double &s12, const double &M, const double &m1, const double &m2, const double &m3)
maximum s13
double q(const double &mab, const double &ma, const double &mb)
Constants q.
EvtTensor4C gmunu_tilde(const EvtVector4R &p4a, const EvtVector4R &p4b, const EvtVector4R &p4c)
Function Tensor gmunu
EvtComplex BreitWigner(const double &m, const double &m0, const double &Gamma0, const double &q, const double &q0, const unsigned int &L)
BreitWigner Shape.
EvtTensor4C Mmunu(const EvtVector4R &p4a, const EvtVector4R &p4b, const EvtVector4R &p4c)
Function Tensor Mmunu.
EvtTensor4C RaiseIndices(const EvtTensor4C &t)
Function RaiseIndices
void RaiseIndex(EvtVector4R &vector)
Member function RaiseIndices.
EvtDecayBase * clone()
Clone the decay of B0toKsKK.
EvtVector4R Lmu(const EvtVector4R &p4a, const EvtVector4R &p4b, const EvtVector4R &p4c)
Function 4Vector Lmu.
EvtComplex A_fxks(const EvtVector4R &p4ks, const EvtVector4R &p4kp, const EvtVector4R &p4km)
A_fxks is amplitude of fxks.
double p(const double &mab, const double &M, const double &mc)
Constants p
EvtComplex Flatte_k(const double &s, const double &m_h)
Constant Flatte_k.
EvtComplex A_kknr(const EvtVector4R &p4k1, const EvtVector4R &p4k2, const double &alpha_kk)
A_kknr is amplitude of kknr.
void initProbMax()
Initialize standard stream objects for probability function
double BWBF(const double &q, const unsigned int &L)
Meson radius
EvtTensor4C Tmunu(const EvtVector4R &p4a, const EvtVector4R &p4b, const EvtVector4R &p4c)
Function Tensor Tmunu
EvtVector4R Boost(const EvtVector4R &p4, const EvtVector4R &boost)
Parameter for boost frame
std::string getName()
Get function Name
#define B2_EVTGEN_REGISTER_MODEL(classname)
Class to register B2_EVTGEN_REGISTER_MODEL.
double sqrt(double a)
sqrt for double
void decay(EvtParticle *p)
Member of particle in EvtGen.
EvtTensor4C Multiply(const EvtTensor4C &t1, const EvtTensor4C &t2)
Function Tensor Multiply
double s13_min(const double &s12, const double &M, const double &m1, const double &m2, const double &m3)
minimum s13
EvtComplex A_chic0ks(const EvtVector4R &p4ks, const EvtVector4R &p4kp, const EvtVector4R &p4km)
A_chic0ks is amplitude of chic0ks.
Abstract base class for different kinds of events.