9#include <framework/logging/Logger.h>
10#include <framework/particledb/EvtGenDatabasePDG.h>
13#include "generators/evtgen/EvtGenModelRegister.h"
14#include "generators/evtgen/models/EvtbTosllNPR.h"
16#include "EvtGenBase/EvtGenKine.hh"
17#include "EvtGenBase/EvtPDL.hh"
18#include "EvtGenBase/EvtParticle.hh"
19#include "EvtGenBase/EvtScalarParticle.hh"
20#include "EvtGenBase/EvtRandom.hh"
22#include "EvtGenBase/EvtAmp.hh"
23#include "EvtGenBase/EvtDiracSpinor.hh"
24#include "EvtGenBase/EvtId.hh"
25#include "EvtGenBase/EvtPDL.hh"
26#include "EvtGenBase/EvtParticle.hh"
27#include "EvtGenBase/EvtTensor4C.hh"
28#include "EvtGenBase/EvtVector4C.hh"
30#include "EvtGenModels/EvtbTosllAmp.hh"
31#include "EvtGenModels/EvtbTosllFF.hh"
40using namespace std::complex_literals;
62static void EvtLeptonVandACurrents(EvtVector4C&, EvtVector4C&,
const EvtDiracSpinor&,
const EvtDiracSpinor&);
82 void init(
const std::vector<EvtPointf>& _v);
88 std::pair<double, double>
operator()(
double)
const;
91 std::vector<EvtPointf>
m_v;
94 std::vector<double>
m_I;
114static double PhaseSpacePole2(
double M,
double m1,
double m2,
double m3, EvtVector4R* p4,
const EvtLinSample* ls)
117 double m13sqmin = (m1 + m3) * (m1 + m3), m13sqmax = (M - m2) * (M - m2);
118 double d13 = m13sqmax - m13sqmin;
119 double M2 = M * M, m32 = m3 * m3, m12 = m1 * m1, m22 = m2 * m2;
120 double c0 = M2 - m32, c1 = m12 - m22, c2 = m32 + m12;
121 double m12sq, m13sq, weight;
123 double z0 = EvtRandom::random(), z1 = EvtRandom::random();
124 m13sq = m13sqmin + z0 * d13;
128 double E3s = c0 - m12sq, E1s = m12sq + c1;
129 double w = 2 * m12sq, e = 4 * m12sq;
130 double A = (m13sq - c2) * w - E3s * E1s;
131 double B = (E3s * E3s - e * m32) * (E1s * E1s - e * m12);
136 double E2 = (M2 + m22 - m13sq) * iM;
137 double E3 = (M2 + m32 - m12sq) * iM;
138 double E1 = M - E2 - E3;
139 double p1 =
sqrt(E1 * E1 - m12);
140 double p3 =
sqrt(E3 * E3 - m32);
141 double cost13 = (2.0 * E1 * E3 + (m12 + m32 - m13sq)) / (2.0 * p1 * p3);
143 double px = p1 * cost13;
145 double v1x = -px - p3;
146 double py = p1 *
sqrt(1.0 - cost13 * cost13);
149 double phi = (EvtRandom::random() - 0.5) * (2 * M_PI);
150 double x = cos(phi), y = sin(phi);
151 phi = (EvtRandom::random() - 0.5) * (2 * M_PI);
152 double c = cos(phi), s = sin(phi);
153 double u = EvtRandom::random() * 2, z = u - 1, t =
sqrt(1 - z * z);
154 double Rx = s * y - c * x, Ry = c * y + s * x, ux = u * x, uy = u * y;
155 double R00 = ux * Rx + c, R01 = s - ux * Ry, R10 = uy * Rx - s,
156 R11 = c - uy * Ry, R20 = -t * Rx, R21 = t * Ry;
157 double pyx = R01 * py, pyy = R11 * py, pyz = R21 * py;
159 p4[0].set(E1, R00 * v0x + pyx, R10 * v0x + pyy, R20 * v0x + pyz);
160 p4[1].set(E2, R00 * v1x - pyx, R10 * v1x - pyy, R20 * v1x - pyz);
161 p4[2].set(E3, R00 * v2x, R10 * v2x, R20 * v2x);
168 static EvtVector4R p4[3];
169 static double mass[3];
170 double m_b = p->mass();
171 for (
int i = 0; i < 3; i++)
172 mass[i] = p->getDaug(i)->mass();
173 EvtId* daughters = getDaugs();
174 double weight = PhaseSpacePole2(m_b, mass[2], mass[1], mass[0], p4,
m_ls);
175 p->getDaug(0)->init(daughters[0], p4[2]);
176 p->getDaug(1)->init(daughters[1], p4[1]);
177 p->getDaug(2)->init(daughters[2], p4[0]);
191 for (
unsigned i = 0; i <
m_v.size() - 1; i++)
192 m_I.push_back((
m_v[i + 1].y +
m_v[i + 0].y) * (
m_v[i + 1].x -
m_v[i + 0].x) +
m_I.back());
193 double norm = 1 /
m_I.back();
194 for (
unsigned i = 0; i <
m_v.size(); i++)
204 int j = upper_bound(
m_I.begin(),
m_I.end(), r) -
m_I.begin();
205 double dI =
m_I[j] -
m_I[j - 1];
206 r = (r -
m_I[j - 1]) / dI;
207 double f0 =
m_v[j - 1].y, f1 =
m_v[j].y, x0 =
m_v[j - 1].x, x1 =
m_v[j].x,
208 df = f1 - f0, dx = x1 - x0, z;
209 if (fabs(df) > f0 * 1e-3) {
210 z = (f1 * x0 - x1 * f0 + dx * sqrt(df * (f0 + f1) * r + f0 * f0)) / df;
213 z = x0 + (r * dx) * (f1 - r * df) / f0;
215 z = x1 - (r * dx) * (f0 + r * df) / f1;
218 return std::make_pair(z, f0 + (df / dx) * (z - x0));
241 BW_t(
double m,
double gtot,
double ghad,
double bll,
double alpha = 1 / 137.035999084)
246 A = (9 / (alpha * alpha)) * bll * gtot * ghad;
252 double R(
double s)
const
254 double z = s -
Mv2, z2 = z * z;
255 return A * s / (z2 +
Mv2 *
Gv2);
261 double Rp(
double s)
const
263 double z = s -
Mv2, z2 = z * z, t = 1 / (z2 +
Mv2 *
Gv2);
264 t *= 1 - 2 * s * z * t;
271 double Rpp(
double s)
const
273 double z = s -
Mv2, z2 = z * z, t = 1 / (z2 +
Mv2 *
Gv2);
274 t = t * t * (8 * s * z2 * t - 4 * z - 2 * s);
284 double zp = sp -
Mv2, zp2 = zp * zp;
285 double z = s -
Mv2, z2 = z * z;
286 double GM2 =
Mv2 *
Gv2, GM = sqrt(GM2);
288 return A * (0.5 * log(ds * ds / (zp2 + GM2)) - z * atan(zp / GM) / GM) / (z2 + GM2);
299 double Mv = sqrt(
Mv2), Gv = sqrt(
Gv2);
301 double smin = (Mv - w * Gv) * (Mv - w * Gv);
302 double smax = (Mv + w * Gv) * (Mv + w * Gv);
303 en.push_back(Mv * Mv);
304 for (
int i = 0, n = 30; i < n; i++) {
305 double s = smin + (smax - smin) * (i + 0.5) / n;
309 double smin1 = (Mv - w * Gv) * (Mv - w * Gv);
310 double smax1 = (Mv + w * Gv) * (Mv + w * Gv);
311 for (
int i = 0, n = 30; i < n; i++) {
312 double s = smin1 + (smax1 - smin1) * (i + 0.5) / n;
313 if (s >= smin && s <= smax)
318 double smin2 = (Mv - w * Gv) * (Mv - w * Gv);
319 double smax2 = (Mv + w * Gv) * (Mv + w * Gv);
320 for (
int i = 0, n = 30; i < n; i++) {
321 double s = smin2 + (smax2 - smin2) * (i + 0.5) / n;
322 if (s >= smin1 && s <= smax1)
327 double smin3 = (Mv - w * Gv) * (Mv - w * Gv);
328 double smax3 = (Mv + w * Gv) * (Mv + w * Gv);
329 for (
int i = 0, n = 30; i < n; i++) {
330 double s = smin3 + (smax3 - smin3) * (i + 0.5) / n;
331 if (s >= smin2 && s <= smax2)
336 double smin4 = (Mv - w * Gv) * (Mv - w * Gv);
337 double smax4 = (Mv + w * Gv) * (Mv + w * Gv);
338 for (
int i = 0, n = 20; i < n; i++) {
339 double s = smin4 + (smax4 - smin4) * (i + 0.5) / n;
340 if (s >= smin3 && s <= smax3)
351static double Rcont(
double s)
354 const double mb = 4.8 , mb2 = mb * mb;
355 if (s < 0.6 * mb2)
return 0;
357 return (1.02 / ((0.69 - 0.6) * mb2)) * (s - 0.6 * mb2);
366static double uR(
double s,
double wjpsi,
double wpsi2s,
double wrest,
const std::vector<BW_t*>& res)
369 for (
unsigned i = 2; i < res.size(); i++) sum += res[i]->
R(s);
370 return Rcont(s) + wjpsi * res[0]->R(s) + wpsi2s * res[1]->R(s) + wrest * sum;
377static double linint(
double a,
double b,
double s,
double sp)
379 double x = std::abs(sp - s);
380 return ((a * s - b) * log(x) + b * log(sp)) / s;
386static double cnstint(
double c,
double s,
double sp)
388 double x = std::abs(sp - s);
389 return c * (log(x) - log(sp)) / s;
395static double DRIcont(
double s)
399 const double mb = 4.8, mb2 = mb * mb;
400 double s0 = 0.6 * mb2, s1 = 0.69 * mb2;
401 double c = 1.02, a = c / (s1 - s0), b = c * s0 / (s1 - s0);
402 return (linint(a, b, s, s1) - linint(a, b, s, s0)) - cnstint(c, s, s1);
410static double DRInt(
double s,
double wjpsi,
double wpsi2s,
double wrest,
const std::vector<BW_t*>& res)
413 const double ct = 4 * Mpi * Mpi;
414 auto Integ = [ct, s](
const BW_t * r) {
return r->DRIntegral(s, 1e8) - r->DRIntegral(s, ct);};
416 for (
unsigned i = 2; i < res.size(); i++) sum += Integ(res[i]);
417 return DRIcont(s) + wjpsi * Integ(res[0]) + wpsi2s * Integ(res[1]) + wrest * sum;
425static std::vector<double> getgrid(
const std::vector<BW_t*>& reso)
427 const double m_c = 1.3, m_s = 0.2;
433 std::vector<double> en;
434 if (reso.size() > 0) {
435 for (
const auto t : reso)
437 const BW_t& jpsi = *(reso[0]);
438 en.push_back(jpsi.
Mv2 + 0.1);
439 en.push_back(jpsi.
Mv2 - 0.1);
440 en.push_back(jpsi.
Mv2 + 0.08);
441 en.push_back(jpsi.
Mv2 - 0.08);
442 en.push_back(jpsi.
Mv2 + 0.065);
443 en.push_back(jpsi.
Mv2 - 0.065);
447 for (
unsigned i = 0, n = 200; i < n; i++) {
448 double s = smax / n * (i + 0.5);
452 double t0, q0 = 4 * m_e * m_e;
453 for (t0 = q0; t0 < 1.004 * q0; t0 *= 1.0008) en.push_back(t0);
454 for (t0 = q0; t0 < 1.020 * q0; t0 *= 1.004) en.push_back(t0);
455 for (t0 = q0; t0 < 1.100 * q0; t0 *= 1.02) en.push_back(t0);
456 for (t0 = q0; t0 < 1.500 * q0; t0 *= 1.10) en.push_back(t0);
457 for (t0 = q0; t0 < 0.150; t0 *= 1.5) en.push_back(t0);
459 double q1 = 4 * m_mu * m_mu;
460 for (t0 = q1; t0 < 1.004 * q1; t0 *= 1.0008) en.push_back(t0);
461 for (t0 = q1; t0 < 1.020 * q1; t0 *= 1.004) en.push_back(t0);
462 for (t0 = q1; t0 < 1.100 * q1; t0 *= 1.02) en.push_back(t0);
463 for (t0 = q1; t0 < 1.500 * q1; t0 *= 1.1) en.push_back(t0);
466 for (t0 = 4 * m_e * m_e; t0 < 0.5; t0 *= 1.5) {
467 en.push_back(4 * MD0 * MD0 + t0);
468 en.push_back(4 * MD0 * MD0 - t0);
471 en.push_back(4 * m_c * m_c);
472 for (t0 = 0.00125; t0 < 0.05; t0 *= 1.5) {
473 en.push_back(4 * m_c * m_c + t0);
474 en.push_back(4 * m_c * m_c - t0);
477 en.push_back(4 * m_s * m_s);
478 for (t0 = 0.00125; t0 < 0.1; t0 *= 1.5) {
479 en.push_back(4 * m_s * m_s + t0);
480 en.push_back(4 * m_s * m_s - t0);
483 std::sort(en.begin(), en.end());
485 for (std::vector<double>::iterator it = en.begin(); it != en.end(); ++it) {
487 en.erase(it, en.end());
492 en.erase(std::unique(en.begin(), en.end()), en.end());
498 EvtId pId = getParentId(), mId = getDaug(0), l1Id = getDaug(1), l2Id = getDaug(2);
500 std::vector<double> v;
501 std::vector<EvtPointf> pp;
503 EvtScalarParticle* scalar_part;
504 EvtParticle* root_part;
506 scalar_part =
new EvtScalarParticle;
509 scalar_part->noLifeTime();
513 p_init.set(EvtPDL::getMass(pId), 0.0, 0.0, 0.0);
514 scalar_part->init(pId, p_init);
515 root_part = (EvtParticle*)scalar_part;
516 root_part->setDiagonalSpinDensity();
518 EvtParticle* daughter, *lep1, *lep2;
527 amp.init(pId, 3, listdaug);
529 root_part->makeDaughters(3, listdaug);
530 daughter = root_part->getDaug(0);
531 lep1 = root_part->getDaug(1);
532 lep2 = root_part->getDaug(2);
535 daughter->noLifeTime();
541 rho.setDiag(root_part->getSpinStates());
545 double m = root_part->mass();
547 EvtVector4R p4meson, p4lepton1, p4lepton2, p4w;
549 double maxprobfound = 0.0;
551 mass[1] = EvtPDL::getMeanMass(l1Id);
552 mass[2] = EvtPDL::getMeanMass(l2Id);
553 std::vector<double> mH;
554 mH.push_back(EvtPDL::getMeanMass(mId));
556 double g = EvtPDL::getWidth(mId);
558 mH.push_back(EvtPDL::getMinMass(mId));
560 std::min(EvtPDL::getMaxMass(mId), m - mass[1] - mass[2]));
561 mH.push_back(mH.front() - g);
562 mH.push_back(mH.front() + g);
565 double q2min = (mass[1] + mass[2]) * (mass[1] + mass[2]);
567 double m0 = EvtPDL::getMinMass(mId);
568 double q2max = (m - m0) * (m - m0);
572 for (
double q2 : v) {
574 if (!(q2min <= q2 && q2 < q2max))
576 double Erho = (m * m + m0 * m0 - q2) / (2.0 * m);
578 m0 = EvtPDL::getMinMass(mId);
579 Erho = (m * m + m0 * m0 - q2) / (2.0 * m);
581 double Prho = sqrt((Erho - m0) * (Erho + m0));
583 p4meson.set(Erho, 0.0, 0.0, -Prho);
584 p4w.set(m - Erho, 0.0, 0.0, Prho);
587 double Elep = sqrt(q2) * 0.5,
588 Plep = sqrt((Elep - mass[1]) * (Elep + mass[1]));
590 const int nj = 3 + 2 + 4 + 8 + 32;
591 double cmin = -1, cmax = 1, dc = (cmax - cmin) / (nj - 1);
593 for (
int j = 0; j < nj; j++) {
594 double c = cmin + dc * j;
597 double Py = Plep * sqrt(1.0 - c * c), Pz = Plep * c;
598 p4lepton1.set(Elep, 0.0, Py, Pz);
599 p4lepton2.set(Elep, 0.0, -Py, -Pz);
601 p4lepton1 = boostTo(p4lepton1, p4w);
602 p4lepton2 = boostTo(p4lepton2, p4w);
605 daughter->init(mId, p4meson);
606 lep1->init(l1Id, p4lepton1);
607 lep2->init(l2Id, p4lepton2);
609 double prob = rho.normalizedProb(amp.getSpinDensity());
610 maxprob = std::max(maxprob, prob);
612 pp.push_back({ (float)q2, (
float)maxprob });
613 maxprobfound = std::max(maxprobfound, maxprob);
615 root_part->deleteTree();
619 setProbMax(maxprobfound);
622static EvtId IdB0, IdaB0, IdBp, IdBm, IdBs, IdaBs, IdRhop, IdRhom, IdRho0,
623 IdOmega, IdKst0, IdaKst0, IdKstp, IdKstm, IdK0, IdaK0, IdKp, IdKm;
624static bool cafirst =
true;
630 IdB0 = EvtPDL::getId(
"B0");
631 IdaB0 = EvtPDL::getId(
"anti-B0");
632 IdBp = EvtPDL::getId(
"B+");
633 IdBm = EvtPDL::getId(
"B-");
634 IdBs = EvtPDL::getId(
"B_s0");
635 IdaBs = EvtPDL::getId(
"anti-B_s0");
636 IdRhop = EvtPDL::getId(
"rho+");
637 IdRhom = EvtPDL::getId(
"rho-");
638 IdRho0 = EvtPDL::getId(
"rho0");
639 IdOmega = EvtPDL::getId(
"omega");
640 IdKst0 = EvtPDL::getId(
"K*0");
641 IdaKst0 = EvtPDL::getId(
"anti-K*0");
642 IdKstp = EvtPDL::getId(
"K*+");
643 IdKstm = EvtPDL::getId(
"K*-");
644 IdK0 = EvtPDL::getId(
"K0");
645 IdaK0 = EvtPDL::getId(
"anti-K0");
646 IdKp = EvtPDL::getId(
"K+");
647 IdKm = EvtPDL::getId(
"K-");
658 checkSpinParent(EvtSpinType::SCALAR);
659 checkSpinDaughter(1, EvtSpinType::DIRAC);
660 checkSpinDaughter(2, EvtSpinType::DIRAC);
662 EvtId mId = getDaug(0);
663 if (mId != IdKst0 && mId != IdaKst0 && mId != IdKstp && mId != IdKstm &&
664 mId != IdK0 && mId != IdaK0 && mId != IdKp && mId != IdKm) {
665 EvtGenReport(EVTGEN_ERROR,
"EvtGen")
666 <<
"EvtbTosllNPR generator expected a K(*) 1st daughter, found: "
667 << EvtPDL::name(getDaug(0)) << endl;
668 EvtGenReport(EVTGEN_ERROR,
"EvtGen")
669 <<
"Will terminate execution!" << endl;
673 auto getInteger = [
this](
int i) ->
int {
674 double a = getArg(i);
675 if (a -
static_cast<int>(a) != 0)
677 EvtGenReport(EVTGEN_ERROR,
"EvtGen")
678 <<
"Parameters is not integer in the BTOSLLNPR decay model: " << i
680 EvtGenReport(EVTGEN_ERROR,
"EvtGen")
681 <<
"Will terminate execution!" << endl;
684 return static_cast<int>(a);
687 double phi = 0.0, n1 = 0.0, d1 = 1.0, eta = 0.0, theta_Jpsi = 0.0, theta_psi2S = 0.0, Nr = 1.0;
689 int i = 0, coordsyst = getInteger(i++);
690 auto getcomplex = [
this, &i, coordsyst]() -> EvtComplex {
691 double a0 = getArg(i++);
692 double a1 = getArg(i++);
693 return (coordsyst) ? EvtComplex(a0 * cos(a1), a0 * sin(a1))
694 : EvtComplex(a0, a1);
696 auto getreal = [
this, &i]() ->
double {
return getArg(i++); };
698 int id = getInteger(i++);
700 m_dc7 = getcomplex();
702 m_dc9 = getcomplex();
706 m_c7p = getcomplex();
708 m_c9p = getcomplex();
726 theta_Jpsi = getreal();
728 theta_psi2S = getreal();
747 BW_t* jpsi =
new BW_t(p_jpsi->Mass(), p_jpsi->Width(), p_jpsi->Width() * 0.877, 0.05971);
748 m_rs.push_back(jpsi);
751 BW_t* psi2s =
new BW_t(p_psi2s->Mass(), p_psi2s->Width(), p_psi2s->Width() * 0.9785, 0.00793);
752 m_rs.push_back(psi2s);
755 BW_t* psi3770 =
new BW_t(p_psi3770->Mass(), p_psi3770->Width(), p_psi3770->Width(), 1e-5);
756 m_rs.push_back(psi3770);
759 BW_t* psi4040 =
new BW_t(p_psi4040->Mass(), p_psi4040->Width(), p_psi4040->Width(), 1.07e-5);
760 m_rs.push_back(psi4040);
763 BW_t* psi4160 =
new BW_t(p_psi4160->Mass(), p_psi4160->Width(), p_psi4160->Width(), 6.9e-6);
764 m_rs.push_back(psi4160);
766 BW_t* psi4230 =
new BW_t(4220 / 1000., 60 / 1000., 60 / 1000., 1e-5);
767 m_rs.push_back(psi4230);
769 std::vector<double> v = getgrid(
m_rs);
770 std::complex<double> eiphi = exp(1i * phi), pJpsi = Nr * exp(1i * theta_Jpsi), ppsi2S = Nr * exp(1i * theta_psi2S);
772 const double C1 = -0.257, C2 = 1.009, C3 = -0.005, C5 = 0;
773 double sc = 1.0 / (4.0 / 3.0 * C1 + C2 + 6.0 * C3 + 60.0 * C5);
774 const double m_b = 4.8, m_c = 1.3;
777 double re = -8.0 / 9.0 * log(m_c / m_b) - 4.0 / 9.0 +
778 t / 3.0 * DRInt(t, real(pJpsi), real(ppsi2S), Nr,
m_rs) - M_PI / 3 * uR(t, imag(pJpsi), imag(ppsi2S), 0,
m_rs);
780 t / 3.0 * DRInt(t, imag(pJpsi), imag(ppsi2S), 0,
m_rs) + M_PI / 3 * uR(t, real(pJpsi), real(ppsi2S), Nr,
m_rs);
781 std::complex<double> r(re, im);
783 double z = t - jpsi->
Mv2;
784 r += eiphi * (n1 * sc / (z * z + d1));
786 m_reso.push_back(std::make_pair(t, r));
796static std::complex<double> hfun(
double q2,
double mq)
798 const double mu = 4.8;
800 return std::complex<double>(8. / 27 + 4. / 9 * log(mu * mu / q2), 4. / 9 * M_PI);
801 double z = 4 * mq * mq / q2;
802 std::complex<double> t;
804 t = atan(1 / sqrt(z - 1));
806 t = std::complex<double>(log((1 +
sqrt(1 - z)) /
sqrt(z)), -M_PI / 2);
808 return -4 / 9. * (2 * log(mq / mu) - 2 / 3. - z) - (4 / 9. * (2 + z) *
sqrt(fabs(z - 1))) * t;
811static EvtComplex C9(
double q2, std::complex<double> hReso = std::complex<double>(0, 0))
814 double C1 = -0.257, C2 = 1.009, C3 = -0.005, C4 = -0.078, C5 = 0, C6 = 0.001;
815 std::complex<double> Y = (hReso) *
816 (4 / 3. * C1 + C2 + 6 * C3 + 60 * C5);
817 Y -= 0.5 * hfun(q2, m_b) * (7 * C3 + 4 / 3. * C4 + 76 * C5 + 64 / 3. * C6);
818 Y -= 0.5 * hfun(q2, 0.0) * (C3 + 4 / 3. * C4 + 16 * C5 + 64 / 3. * C6);
819 Y += 4 / 3. * C3 + 64 / 9. * C5 + 64 / 27. * C6;
820 return EvtComplex(4.211 + real(Y), imag(Y));
823static std::complex<double> interpol(
const hvec_t& v,
double s)
829 int j = std::lower_bound(v.begin(), v.end(), s, [](
const helem_t& a,
double b) { return a.first < b; }) - v.begin();
831 double dx = v[j].first - v[j - 1].first;
832 auto dy = v[j].second - v[j - 1].second;
833 return v[j - 1].second + (s - v[j - 1].first) * (dy / dx);
838 EvtId dId = parent->getDaug(0)->getId();
839 if (dId == IdKst0 || dId == IdaKst0 || dId == IdKstp || dId == IdKstm) {
842 if (dId == IdK0 || dId == IdaK0 || dId == IdKp || dId == IdKm) {
847static inline double poly(
double x,
int n,
const double* c)
850 while (n--) t = c[n] + x * t;
854static void getVectorFF(
double q2,
double mB,
double mV,
double& a1,
double& a2,
double& a0,
double& v,
double& t1,
double& t2,
857 static const double alfa[7][4] = {
860 { 5.415000, 0.376313, -1.165970, 2.424430 },
861 { 5.366000, 0.369196, -1.365840, 0.128191 },
862 { 5.829000, 0.297250, 0.392378, 1.189160 },
863 { 5.829000, 0.265375, 0.533638, 0.483166 },
864 { 5.415000, 0.312055, -1.008930, 1.527200 },
865 { 5.829000, 0.312055, 0.496846, 1.614310 },
866 { 5.829000, 0.667412, 1.318120, 3.823340 }
869 double mBaV = mB + mV, mBsV = mB - mV;
870 double tp = mBaV * mBaV;
871 double s0 =
sqrt(2 * mBaV *
sqrt(mB * mV));
872 double z0 = (mBaV - s0) / (mBaV + s0);
874 double s =
sqrt(tp - q2), z = (s - s0) / (s + s0), dz = z - z0, ff[7];
875 for (
int j = 0; j < 7; j++) {
876 double mR = alfa[j][0], mR2 = mR * mR;
877 ff[j] = (mR2 / (mR2 - q2)) * poly(dz, 3, alfa[j] + 1);
882 double lambda = (mBaV * mBaV - q2) * (mBsV * mBsV - q2);
890 a2 = mBaV * ((mBaV * (mBaV * mBsV - q2)) * a1 - (16 * mB * mV * mV) * a12) / lambda;
896 t3 = mBsV * (mBaV * (mB * mB + 3 * mV * mV - q2) * t2 - (8 * mB * mV * mV) * t23) / lambda;
899static EvtTensor4C asymProd(
const EvtVector4R& a,
const EvtVector4R& b)
905 double t01 = a.get(2) * b.get(3) - a.get(3) * b.get(2);
908 double t02 = a.get(3) * b.get(1) - a.get(1) * b.get(3);
911 double t03 = a.get(1) * b.get(2) - a.get(2) * b.get(1);
914 double t12 = a.get(0) * b.get(3) - a.get(3) * b.get(0);
917 double t13 = a.get(2) * b.get(0) - a.get(0) * b.get(2);
920 double t23 = a.get(0) * b.get(1) - a.get(1) * b.get(0);
930static void addScaled(EvtTensor4C& out, EvtComplex scale,
const EvtTensor4C& in)
932 for (
int i = 0; i < 4; i++)
933 for (
int j = 0; j < 4; j++) {
934 out.set(i, j, out.get(i, j) + scale * in.get(i, j));
941static EvtComplex cprod(
const EvtComplex& c1,
const EvtComplex& c2)
943 return EvtComplex(real(c1) * real(c2) + imag(c1) * imag(c2),
944 real(c1) * imag(c2) - imag(c1) * real(c2));
950static void EvtLeptonVandACurrents(EvtVector4C& v, EvtVector4C& a,
const EvtDiracSpinor& x,
const EvtDiracSpinor& y)
952 EvtComplex w0 = cprod(x.get_spinor(0), y.get_spinor(0)), w1 = cprod(x.get_spinor(1), y.get_spinor(1)), w2 = cprod(x.get_spinor(2),
953 y.get_spinor(2)), w3 = cprod(x.get_spinor(3), y.get_spinor(3));
954 EvtComplex w02 = w0 + w2, w13 = w1 + w3;
955 EvtComplex W1 = w02 + w13, W2 = w02 - w13;
956 EvtComplex q0 = cprod(x.get_spinor(0), y.get_spinor(2)), q1 = cprod(x.get_spinor(1), y.get_spinor(3)), q2 = cprod(x.get_spinor(2),
957 y.get_spinor(0)), q3 = cprod(x.get_spinor(3), y.get_spinor(1));
958 EvtComplex q02 = q0 + q2, q13 = q1 + q3;
959 EvtComplex Q1 = q02 + q13, Q2 = q02 - q13;
960 EvtComplex e0 = cprod(x.get_spinor(0), y.get_spinor(3)), e1 = cprod(x.get_spinor(1), y.get_spinor(2)), e2 = cprod(x.get_spinor(2),
961 y.get_spinor(1)), e3 = cprod(x.get_spinor(3), y.get_spinor(0));
962 EvtComplex e20 = e0 + e2, e13 = e1 + e3;
963 EvtComplex E1 = e13 + e20, E2 = e13 - e20;
964 v.set(W1, E1, EvtComplex(-imag(E2), real(E2)), Q2);
965 EvtComplex t0 = cprod(x.get_spinor(0), y.get_spinor(1)), t1 = cprod(x.get_spinor(1), y.get_spinor(0)), t2 = cprod(x.get_spinor(2),
966 y.get_spinor(3)), t3 = cprod(x.get_spinor(3), y.get_spinor(2));
967 EvtComplex t20 = t0 + t2, t13 = t1 + t3;
968 EvtComplex T1 = t13 + t20, T2 = t13 - t20;
969 a.set(Q1, T1, EvtComplex(-imag(T2), real(T2)), W2);
975 EvtId pId = parent->getId();
977 EvtVector4R q = parent->getDaug(1)->getP4() + parent->getDaug(2)->getP4();
978 double q2 = q.mass2();
979 double mesonmass = parent->getDaug(0)->mass();
980 double parentmass = parent->mass();
982 double a1, a2, a0, v, t1, t2, t3;
983 getVectorFF(q2, parentmass, mesonmass, a1, a2, a0, v, t1, t2, t3);
985 const EvtVector4R& p4meson = parent->getDaug(0)->getP4();
986 double pmass = parent->mass(), ipmass = 1 / pmass;
987 const EvtVector4R p4b(pmass, 0.0, 0.0, 0.0);
988 EvtVector4R pbhat = p4b * ipmass;
989 EvtVector4R qhat = q * ipmass;
990 EvtVector4R pkstarhat = p4meson * ipmass;
991 EvtVector4R phat = pbhat + pkstarhat;
993 EvtComplex c7 = -0.304;
994 EvtComplex c9 = C9(q2, interpol(
m_reso, q2));
995 EvtComplex c10 = -4.103;
1000 EvtComplex I(0.0, 1.0);
1002 double mb = 4.8 * ipmass, ms = 0.093 * ipmass;
1003 double mH = mesonmass * ipmass, oamH = 1 + mH, osmH = 1 - mH,
1004 osmH2 = oamH * osmH, iosmH2 = 1 / osmH2;
1005 double is = pmass * pmass / q2;
1009 double cs0 = (a1 - a2 - a0) * is;
1013 EvtComplex a = (c9 +
m_c9p) * v + (c7 +
m_c7p) * (4 * mb * is * t1);
1014 EvtComplex b = (c9 -
m_c9p) * a1 +
1015 (c7 -
m_c7p) * (2 * mb * is * osmH2 * t2);
1016 EvtComplex c = (c9 -
m_c9p) * a2 +
1017 (c7 -
m_c7p) * (2 * mb * (t3 * iosmH2 + t2 * is));
1018 EvtComplex d = (c9 -
m_c9p) * cs0 - (c7 -
m_c7p) * (2 * mb * is * t3);
1019 EvtComplex e = (c10 +
m_c10p) * v;
1020 EvtComplex f = (c10 -
m_c10p) * a1;
1021 EvtComplex g = (c10 -
m_c10p) * a2;
1022 EvtComplex h = (c10 -
m_c10p) * cs0;
1023 double sscale = a0 / (mb + ms);
1024 EvtComplex hs =
m_cS * sscale, hp =
m_cP * sscale;
1026 int charge1 = EvtPDL::chg3(parent->getDaug(1)->getId());
1028 EvtParticle* lepPos = parent->getDaug(2 - (charge1 > 0));
1029 EvtParticle* lepNeg = parent->getDaug(1 + (charge1 > 0));
1031 EvtDiracSpinor lp0(lepPos->spParent(0)), lp1(lepPos->spParent(1));
1032 EvtDiracSpinor lm0(lepNeg->spParent(0)), lm1(lepNeg->spParent(1));
1034 EvtVector4C l11, l12, l21, l22, a11, a12, a21, a22;
1035 EvtComplex s11, s12, s21, s22, p11, p12, p21, p22;
1036 EvtTensor4C tt0 = asymProd(pbhat, pkstarhat);
1038 EvtTensor4C T1(tt0), T2(tt0);
1039 const EvtTensor4C& gmn = EvtTensor4C::g();
1040 EvtTensor4C tt1(EvtGenFunctions::directProd(pbhat, phat));
1041 EvtTensor4C tt2(EvtGenFunctions::directProd(pbhat, qhat));
1049 if (pId == IdBm || pId == IdaB0 || pId == IdaBs) {
1051 addScaled(T1, -b, gmn);
1052 addScaled(T1, c, tt1);
1053 addScaled(T1, d, tt2);
1055 addScaled(T2, -f, gmn);
1056 addScaled(T2, g, tt1);
1057 addScaled(T2, h, tt2);
1059 EvtLeptonVandACurrents(l11, a11, lp0, lm0);
1060 EvtLeptonVandACurrents(l21, a21, lp1, lm0);
1061 EvtLeptonVandACurrents(l12, a12, lp0, lm1);
1062 EvtLeptonVandACurrents(l22, a22, lp1, lm1);
1064 s11 = EvtLeptonSCurrent(lp0, lm0);
1065 p11 = EvtLeptonPCurrent(lp0, lm0);
1066 s21 = EvtLeptonSCurrent(lp1, lm0);
1067 p21 = EvtLeptonPCurrent(lp1, lm0);
1068 s12 = EvtLeptonSCurrent(lp0, lm1);
1069 p12 = EvtLeptonPCurrent(lp0, lm1);
1070 s22 = EvtLeptonSCurrent(lp1, lm1);
1071 p22 = EvtLeptonPCurrent(lp1, lm1);
1072 }
else if (pId == IdBp || pId == IdB0 || pId == IdBs) {
1074 addScaled(T1, -b, gmn);
1075 addScaled(T1, c, tt1);
1076 addScaled(T1, d, tt2);
1078 addScaled(T2, -f, gmn);
1079 addScaled(T2, g, tt1);
1080 addScaled(T2, h, tt2);
1082 EvtLeptonVandACurrents(l11, a11, lp1, lm1);
1083 EvtLeptonVandACurrents(l21, a21, lp0, lm1);
1084 EvtLeptonVandACurrents(l12, a12, lp1, lm0);
1085 EvtLeptonVandACurrents(l22, a22, lp0, lm0);
1087 s11 = EvtLeptonSCurrent(lp1, lm1);
1088 p11 = EvtLeptonPCurrent(lp1, lm1);
1089 s21 = EvtLeptonSCurrent(lp0, lm1);
1090 p21 = EvtLeptonPCurrent(lp0, lm1);
1091 s12 = EvtLeptonSCurrent(lp1, lm0);
1092 p12 = EvtLeptonPCurrent(lp1, lm0);
1093 s22 = EvtLeptonSCurrent(lp0, lm0);
1094 p22 = EvtLeptonPCurrent(lp0, lm0);
1096 EvtGenReport(EVTGEN_ERROR,
"EvtGen") <<
"Wrong lepton number\n";
1100 for (
int i = 0; i < 3; i++) {
1101 EvtVector4C eps = parent->getDaug(0)->epsParent(i).conj();
1103 EvtVector4C E1 = T1.cont1(eps);
1104 EvtVector4C E2 = T2.cont1(eps);
1106 EvtComplex epsq = I * (eps * q), epsqs = epsq * hs, epsqp = epsq * hp;
1108 amp.vertex(i, 0, 0, l11 * E1 + a11 * E2 - s11 * epsqs - p11 * epsqp);
1109 amp.vertex(i, 0, 1, l12 * E1 + a12 * E2 - s12 * epsqs - p12 * epsqp);
1110 amp.vertex(i, 1, 0, l21 * E1 + a21 * E2 - s21 * epsqs - p21 * epsqp);
1111 amp.vertex(i, 1, 1, l22 * E1 + a22 * E2 - s22 * epsqs - p22 * epsqp);
1115void getScalarFF(
double q2,
double& fp,
double& f0,
double& fT)
1121 static const double MK = 0.495644, MB = 5.279495108051249,
1122 MBs0 = 5.729495108051249, MBsstar = 5.415766151925566,
1123 logsB = 1.3036556717286227;
1124 static const int N = 3;
1125 static const double a0B[] = { 0.2546162729845652, 0.211016713841977,
1126 0.02690776720598433, 0.2546162729845652,
1127 -0.7084710010870424, 0.3096901516968882,
1128 0.2549078414069112, -0.6549384905373116,
1131 *apB = a0B + N, *aTB = apB + N;
1133 double pole0 = 1 / (1 - q2 / (MBs0 * MBs0));
1134 double polep = 1 / (1 - q2 / (MBsstar * MBsstar));
1135 double B = MB + MK, A =
sqrt(B * B - q2), z = (A - B) / (A + B),
1136 zN = z * z * z / N, zn = z;
1140 for (
int i = 1; i < N; i++) {
1144 double izN = i * zN;
1154 f0 *= pole0 * logsB;
1155 fp *= polep * logsB;
1156 fT *= polep * logsB;
1162 EvtId pId = parent->getId();
1164 EvtVector4R q = parent->getDaug(1)->getP4() + parent->getDaug(2)->getP4();
1165 double q2 = q.mass2();
1166 double dmass = parent->getDaug(0)->mass();
1169 getScalarFF(q2, fp, f0, ft);
1170 const EvtVector4R& k = parent->getDaug(0)->getP4();
1171 double pmass = parent->mass();
1172 const EvtVector4R p(pmass, 0.0, 0.0, 0.0);
1173 EvtVector4R pk = p + k;
1175 EvtComplex c7 = -0.304;
1176 EvtComplex c9 = C9(q2, interpol(
m_reso, q2));
1177 EvtComplex c10 = -4.103;
1182 double mb = 4.8 , ms = 0.093 ;
1184 int charge1 = EvtPDL::chg3(parent->getDaug(1)->getId());
1186 EvtParticle* lepPos = (charge1 > 0) ? parent->getDaug(1)
1187 : parent->getDaug(2);
1188 EvtParticle* lepNeg = (charge1 < 0) ? parent->getDaug(1)
1189 : parent->getDaug(2);
1191 EvtDiracSpinor lp0(lepPos->spParent(0)), lp1(lepPos->spParent(1));
1192 EvtDiracSpinor lm0(lepNeg->spParent(0)), lm1(lepNeg->spParent(1));
1194 EvtVector4C l11, l12, l21, l22, a11, a12, a21, a22;
1195 EvtComplex s11, s12, s21, s22, p11, p12, p21, p22;
1197 if (pId == IdBm || pId == IdaB0 || pId == IdaBs) {
1198 EvtLeptonVandACurrents(l11, a11, lp0, lm0);
1199 EvtLeptonVandACurrents(l21, a21, lp1, lm0);
1200 EvtLeptonVandACurrents(l12, a12, lp0, lm1);
1201 EvtLeptonVandACurrents(l22, a22, lp1, lm1);
1203 s11 = EvtLeptonSCurrent(lp0, lm0);
1204 p11 = EvtLeptonPCurrent(lp0, lm0);
1205 s21 = EvtLeptonSCurrent(lp1, lm0);
1206 p21 = EvtLeptonPCurrent(lp1, lm0);
1207 s12 = EvtLeptonSCurrent(lp0, lm1);
1208 p12 = EvtLeptonPCurrent(lp0, lm1);
1209 s22 = EvtLeptonSCurrent(lp1, lm1);
1210 p22 = EvtLeptonPCurrent(lp1, lm1);
1211 }
else if (pId == IdBp || pId == IdB0 || pId == IdBs) {
1212 EvtLeptonVandACurrents(l11, a11, lp1, lm1);
1213 EvtLeptonVandACurrents(l21, a21, lp0, lm1);
1214 EvtLeptonVandACurrents(l12, a12, lp1, lm0);
1215 EvtLeptonVandACurrents(l22, a22, lp0, lm0);
1217 s11 = EvtLeptonSCurrent(lp1, lm1);
1218 p11 = EvtLeptonPCurrent(lp1, lm1);
1219 s21 = EvtLeptonSCurrent(lp0, lm1);
1220 p21 = EvtLeptonPCurrent(lp0, lm1);
1221 s12 = EvtLeptonSCurrent(lp1, lm0);
1222 p12 = EvtLeptonPCurrent(lp1, lm0);
1223 s22 = EvtLeptonSCurrent(lp0, lm0);
1224 p22 = EvtLeptonPCurrent(lp0, lm0);
1226 EvtGenReport(EVTGEN_ERROR,
"EvtGen") <<
"Wrong lepton number\n";
1228 double dm2 = pmass * pmass - dmass * dmass, t0 = dm2 / q2,
1229 t1 = 2 * mb * ft / (pmass + dmass);
1233 EvtVector4C E1 = (c9 * fp + c7 * t1) * pk +
1234 (t0 * (c9 * (f0 - fp) - c7 * t1)) * q;
1235 EvtVector4C E2 = (c10 * fp) * pk + (t0 * (f0 - fp)) * q;
1236 double s = dm2 / (mb - ms) * f0;
1237 amp.vertex(0, 0, l11 * E1 + a11 * E2 + (
m_cS * s11 +
m_cP * p11) * s);
1238 amp.vertex(0, 1, l12 * E1 + a12 * E2 + (
m_cS * s12 +
m_cP * p12) * s);
1239 amp.vertex(1, 0, l21 * E1 + a21 * E2 + (
m_cS * s21 +
m_cP * p21) * s);
1240 amp.vertex(1, 1, l22 * E1 + a22 * E2 + (
m_cS * s22 +
m_cP * p22) * s);
static EvtGenDatabasePDG * Instance()
Instance method that loads the EvtGen table.
Description: Implementation of the B->K(*)l^+l^- decays with New Physics contributions and c\bar{c}...
std::vector< BW_t * > m_rs
vector with c\bar{c} resonance lineshapes
void CalcSAmp(EvtParticle *, EvtAmp &)
The method to evaluate the scalar kaon decay amplitude.
EvtDecayBase * clone() override
The function which makes a copy of the model.
void decay(EvtParticle *p) override
The method to calculate a decay amplitude.
EvtLinSample * m_ls
piece-wise interpolation of maximum of the matrix element depend on q^2
EvtComplex m_dc10
delta C_10eff – addition to NNLO SM value
hvec_t m_reso
tabulated resonance contribution
std::string getName() override
The function which returns the name of the model.
void CalcAmp(EvtParticle *, EvtAmp &)
The method to evaluate the decay amplitude.
EvtComplex m_cS
(C_S - C'_S) – scalar right and left polarizations
void CalcVAmp(EvtParticle *, EvtAmp &)
The method to evaluate the vector kaon decay amplitude.
EvtComplex m_cP
(C_P - C'_P) – pseudo-scalar right and left polarizations
virtual ~EvtbTosllNPR()
The destructor to clean up objects created during the initialization.
EvtComplex m_c9p
C'_9eff – right hand polarizations.
EvtComplex m_dc7
delta C_7eff – addition to NNLO SM value
void init() override
Initialization method.
int m_flag
flag is set nonzero to include resonances
EvtComplex m_dc9
delta C_9eff – addition to NNLO SM value
EvtComplex m_c7p
C'_7eff – right hand polarizations.
EvtComplex m_c10p
c'_10eff – right hand polarizations
void initProbMax() override
The method to evaluate the maximum amplitude.
#define B2_EVTGEN_REGISTER_MODEL(classname)
Class to register B2_EVTGEN_REGISTER_MODEL.
double sqrt(double a)
sqrt for double
Description: The class to treat resonances, their lineshapes, and dispersion relation integral.
BW_t(double m, double gtot, double ghad, double bll, double alpha=1/137.035999084)
m is the mass of a resonance gtot is the total width of a resonance ghad is the hadronic width of a r...
double A
resonance total amplitude
void addknots(std::vector< double > &en) const
Based on the resonance parameters place knots for linear piece-wise approximation with enough precisi...
double Mv2
resonance mass squared
double R(double s) const
resonance shape
double Rp(double s) const
the first derivative of a resonance shape
double DRIntegral(double s, double sp) const
Dispersion relation indefinite integral: s is the parmeter, sp is the integration variable.
double Rpp(double s) const
the second derivative of a resonance shape
double Gv2
resonance width squared
Description: The simple helper class to generate arbitrary distribution based on a linear piece-wise ...
std::pair< double, double > operator()(double) const
The operator to return <x,weight> depend on argument in [0,1] range which represent fraction of the t...
EvtLinSample()
The default constructor.
void init(const std::vector< EvtPointf > &_v)
The method to initialize data.
std::vector< double > m_I
cumulative integral data
EvtLinSample(const std::vector< EvtPointf > &_v)
The primary constructor with an actual shape to generate.
std::vector< EvtPointf > m_v
shape data
Description: struct for 2-dimensional point in the float format instead of double to save memory.