a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
MPll Class Reference

A class for the \(M \to P l^+ l^-\) decay. More...

#include <MPll.h>

Detailed Description

A class for the \(M \to P l^+ l^-\) decay.

Author
HEPfit Collaboration

This class is used to compute all the functions needed in order to build the observables relative to the \(M \to P l^+ l^-\) decays, where \(M\) is a generic meson and \(P\) is a pseudoscalar meson.

MPll parameters

The mandatory parameters of MPll are summarized below:

Label LaTeX symbol Description
r_1_fplus, r_2_fplus, m_fit2_fplus \(r_1^{f_+}, r_2^{f_+}, m_{fit}^{2,f_+}\) The fit parameters for the LCSR form factor \(f_+\) of the \(B\to K\).
r_1_fT, r_2_fT, m_fit2_fT \(r_1^{f_T}, r_2^{f_T}, m_{fit}^{2,f_T}\) The fit parameters for the LCSR form factor \(f_T\) of the \(B\to K\).
r_2_f0, m_fit2_f0 \(r_2^{f_0}, m_{fit}^{2,f_0}\) The fit parameters for the LCSR form factor \(f_0\) of the \(B\to K\).
absh_0_MP, absh_0_1_MP \(\mathrm{Abs}h_0^{(0)}, \mathrm{Abs}h_0^{(1)}\) The constant and linear terms of the absolute value of the hadronic parameter \(h_0\) of the \(B\to K\).
argh_0_MP, argh_0_1_MP \(\mathrm{Arg}h_0^{(0)}, \mathrm{Arg}h_0^{(1)}\) The constant and linear terms of the argument of the hadronic parameter \(h_0\) of the \(B\to K\).

This kind of decays can be described by means of the \(\Delta B = 1 \) weak effective Hamiltonian

\[ \mathcal{H}_\mathrm{eff}^{\Delta B = 1} = \mathcal{H}_\mathrm{eff}^\mathrm{had} + \mathcal{H}_\mathrm{eff}^\mathrm{sl+\gamma}, \]

where the first term is the hadronic contribution

\[ \mathcal{H}_\mathrm{eff}^\mathrm{had} = \frac{4G_F}{\sqrt{2}}\Bigg[\sum_{p=u,c}\lambda_p\bigg(C_1 Q^{p}_1 + C_2 Q^{p}_2\bigg) -\lambda_t \bigg(\sum_{i=3}^{6} C_i P_i + C_{8}Q_{8g} \bigg)\Bigg] \,, \]

involving current-current, chromodynamic penguin and chromomagnetic dipole operators, while the second one, given by

\[ \mathcal{H}_\mathrm{eff}^\mathrm{sl+\gamma} = - \frac{4G_F}{\sqrt{2}}\lambda_t \bigg( C_7Q_{7\gamma} + C_9Q_{9V} + C_{10}Q_{10A} \bigg) \,, \]

includes the electromagnetic penguin plus the semileptonic operators.

Considering the matrix element of \(\mathcal{H}_\mathrm{eff}^{\Delta B = 1}\) between the initial state \(M\) and the final state \(P l^+ l^-\), only the contribution of \(\mathcal{H}_\mathrm{eff}^\mathrm{sl+\gamma}\) clearly factorizes into the product of hadronic form factors and leptonic tensors at all orders in strong interactions. Following [144], we implemented the amplitude in the helicity basis; hence we made use of the helicity form factors \( \tilde{V}_0(q^2), \tilde{T}_0(q^2)\) and \(\tilde{S}(q^2) \), which are related to the ones in the transverse basis through the following relations :

\[ \tilde{V}_0(q^2) = i \frac{\sqrt{\lambda(q^2)}}{2m_M\sqrt{q^2}}f_+(q^2)\,,\\ \tilde{T}_0(q^2) = i \frac{\sqrt{\lambda(q^2)q^2}}{2m_M^2(m_M+m_P)}f_T(q^2)\,,\\ \tilde{S}(q^2) = -\frac{m_M^2-m_P^2}{2m_M(m_b+m_s)}\frac{1+m_s/m_b}{1-m_s/m_b}f_0(q^2)\,, \]

where \(\lambda(q^2) = 4m_M^2|\vec{k}|^2\), with \(\vec{k}\) as the 3-momentum of the meson \(P\) in the \(M\) rest frame.

The effect of the operators of \(\mathcal{H}_\mathrm{eff}^\mathrm{had}\) due to exchange of soft gluon can be reabsorbed in the following parameterization,

\[ h_0(q^2) = \frac{\epsilon^*_\mu(\lambda)}{m_M^2} \int d^4x e^{iqx} \langle \bar P \vert T\{j^{\mu}_\mathrm{em} (x) \mathcal{H}_\mathrm{eff}^\mathrm{had} (0)\} \vert \bar M \rangle = h_0^{(0)} + \frac{q^2}{1\,\mathrm{GeV}^2} h_0^{(1)}\,. \]

The amplitude can be therefore parametrized in terms of the following helicity amplitudes:

\[ H_V = -i\, N \Big\{C_{9} \tilde{V}_{L,0} +C_{9}' \tilde{V}_{R,0} + \frac{m_M^2}{q^2} \Big[\frac{2\, m_b}{m_M} (C_{7} \tilde{T}_{L,0} + C_{7}' \tilde{T}_{R,0}) - 16 \pi^2 h_0 \Big] \Big\} \,, \\ H_A = -i\, N (C_{10} \tilde{V}_{L,0} + C_{10}'\tilde{V}_{R,0}) \,, \\ H_S = i\, N \frac{ m_b}{m_W} (C_S \tilde{S}_L + C_S' \tilde{S}_R)\,, \\ H_P = i\, N \Big\{ \frac{ m_b}{m_W} (C_P \tilde{S}_L + C_P' \tilde{S}_R) + \frac{2\,m_\ell m_b}{q^2} \left[C_{10} \Big(\tilde{S}_L - \frac{m_s}{m_b} \tilde{S}_R \Big) + C_{10}' \Big(\tilde{S}_R - \frac{m_s}{m_b} \tilde{S}_L\Big) \right] \Big\} \,, \]

where \( N = - \frac{4 G_F m_M}{\sqrt{2}}\frac{e^2}{16\pi^2}\lambda_t\) and we have defined

\[ \tilde{V}_{L,0}(q^2) = -\tilde{V}_{R,0}(q^2)=\tilde{V}_0(q^2)\,,\\ \tilde{T}_{L,0}(q^2) = -\tilde{T}_{R,0}(q^2)=\tilde{V}_0(q^2)\,,\\ \tilde{S}_L(q^2) = -\tilde{S}_R(q^2)=\tilde{S}(q^2)\,. \]

Squaring the amplitude and summing over the spins it is possible to obtain the fully differential decay rate, which is

\[ \frac{d^{(4)} \Gamma}{dq^2\,d(\cos\theta_l)} = \frac{9}{32\,\pi} \Big( I^c_1 +I^c_2\cos2\theta_l + I_6^c \cos\theta_l \Big) \]

The angular coefficients involved in the differential decay rate are related to the helicity amplitudes according to the following relations:

\[ I_1^c = F \left\{ \frac{1}{2}\left(|H_V^0|^2+|H_A^0|^2\right)+ |H_P|^2+\frac{2m_\ell^2}{q^2}\left(|H_V^0|^2-|H_A^0|^2\right) + \beta^2 |H_S|^2 \right\}\,,\\ I_2^c = -F\, \frac{\beta^2}{2}\left(|H_V^0|^2+|H_A^0|^2\right)\,,\\ I_6^c = 2 F \frac{\beta\, m_\ell}{\sqrt{q^2}} {\rm Re} \left[ H_S^* H_V^0 \right]\,,\\ \]

where

\[ F=\frac{ \lambda^{1/2}\beta\, q^2}{3 \times 2^{5} \,\pi^3\, m_M^3}\,, \qquad \beta = \sqrt{1 - \frac{4 m_\ell^2}{q^2} }\,. \]

The final observables are hence build employing CP-averages \(\Sigma_i\) or CP-asymmetries \(\Delta_i\) of such angular coefficients; however, since on the experimental side the observables are averaged over \( q^2 \) bins, an integration of the coeffiecients over such bins has to be performed before they are combined in order to build the observables.

The class is organized as follows: after the parameters are updated in updateParameters() and the cache is checked in checkCache(), the form factor are build in the transverse basis in the functions f_plus(), f_0() and f_T() [31]. They are consequentely translated in the helicity basis through the functions V_L(), V_R(), T_L(), T_R(), S_L() and S_R(). Form factors and parameters are combined together in the functions H_V(), H_A(), H_S() and H_P() in order to build the helicity amplitudes, which are consequentely combined to create the angular coefficients in the function I(). Those coefficients are used to create the CP averaged coefficients in the function Sigma() ad the CP asymmetric coefficients in the function Delta(). Form factors, CP averaged and asymmetric coefficients and hadronic contributions are integrated in the functions integrateSigma() and integrateDelta() in order to be further used to build the observables.

Definition at line 172 of file MPll.h.

Public Member Functions

gslpp::complex DeltaC9_KD (double q2)
 
gslpp::complex funct_g (double q2)
 
double getSigma (int i, double q_2)
 The value of \( \Sigma_{i} \) from \(q_{min}\) to \(q_{max}\). More...
 
double getwidth ()
 The width of the meson M. More...
 
gslpp::complex H_A (double q2)
 The helicity amplitude \( H_A^{\lambda} \) . More...
 
gslpp::complex h_lambda (double q2)
 The non-pertubative ccbar contributions to the helicity amplitudes. More...
 
gslpp::complex H_P (double q2)
 The helicity amplitude \( H_P^{\lambda} \) . More...
 
gslpp::complex H_S (double q2)
 The helicity amplitude \( H_S^{\lambda} \) . More...
 
gslpp::complex H_V (double q2)
 The helicity amplitude \( H_V^{\lambda} \) . More...
 
std::vector< std::string > initializeMPllParameters ()
 A method for initializing the parameters necessary for MPll. More...
 
double integrateDelta (int i, double q_min, double q_max)
 The integral of \( \Delta_{i} \) from \(q_{min}\) to \(q_{max}\). More...
 
double integrateSigma (int i, double q_min, double q_max)
 The integral of \( \Sigma_{i} \) from \(q_{min}\) to \(q_{max}\). More...
 
 MPll (const StandardModel &SM_i, QCD::meson meson_i, QCD::meson pseudoscalar_i, QCD::lepton lep_i)
 Constructor. More...
 
virtual ~MPll ()
 Destructor. More...
 

Private Attributes

double ale
 
bool dispersion
 
bool FixedWCbtos
 
double GF
 
QCD::lepton lep
 
double Mb
 
double mb_pole
 
double Mc
 
double mc_pole
 
QCD::meson meson
 
double mJ2
 
double Mlep
 
double MM
 
double MP
 
std::vector< std::string > mpllParameters
 
double Ms
 
double mu_b
 
double mu_h
 
std::unique_ptr< F_1myF_1
 
std::unique_ptr< F_2myF_2
 
const StandardModelmySM
 
QCD::meson pseudoscalar
 
double spectator_charge
 
double width
 

Constructor & Destructor Documentation

◆ MPll()

MPll::MPll ( const StandardModel SM_i,
QCD::meson  meson_i,
QCD::meson  pseudoscalar_i,
QCD::lepton  lep_i 
)

Constructor.

Parameters
[in]SM_ia reference to an object of type StandardModel
[in]meson_iinitial meson of the decay
[in]pseudoscalar_ifinal pseudoscalar meson of the decay
[in]lep_ifinal leptons of the decay

Definition at line 19 of file MPll.cpp.

20 : mySM(SM_i), myF_1(new F_1()), myF_2(new F_2()),
21 #if LATTICE
22 fplus_lat_cache(3, 0.),
23 fT_lat_cache(3, 0.),
24 f0_lat_cache(3, 0.),
25 #else
26 fplus_cache(2, 0.),
27 fT_cache(2, 0.),
28 f0_cache(2, 0.),
29 #endif
30 k2_cache(2, 0.),
31 SL_cache(2, 0.),
32 N_cache(3, 0.),
33 Ycache(2, 0.),
34 H_V0cache(2, 0.),
35 H_Scache(2, 0.),
36 H_P_cache(4, 0.),
37 T_cache(5, 0.)
38 {
39  lep = lep_i;
40  meson = meson_i;
41  pseudoscalar = pseudoscalar_i;
42  dispersion = true;
43  FixedWCbtos = false;
44  mJ2 = 3.096 * 3.096;
45 
46  I0_updated = 0;
47  I2_updated = 0;
48  I8_updated = 0;
49 
50  VL_updated = 0;
51  TL_updated = 0;
52  SL_updated = 0;
53 
54  deltaTparpupdated = 0;
55  deltaTparmupdated = 0;
56 
57  w_sigma = gsl_integration_cquad_workspace_alloc(100);
58  w_delta = gsl_integration_cquad_workspace_alloc(100);
59  w_DTPPR = gsl_integration_cquad_workspace_alloc(100);
60 
61  acc_Re_deltaC7_QCDF = gsl_interp_accel_alloc();
62  acc_Im_deltaC7_QCDF = gsl_interp_accel_alloc();
63  acc_Re_deltaC9_QCDF = gsl_interp_accel_alloc();
64  acc_Im_deltaC9_QCDF = gsl_interp_accel_alloc();
65 
66  spline_Re_deltaC7_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
67  spline_Im_deltaC7_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
68  spline_Re_deltaC9_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
69  spline_Im_deltaC9_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
70 }

◆ ~MPll()

MPll::~MPll ( )
virtual

Destructor.

Definition at line 72 of file MPll.cpp.

73 {
74 }

Member Function Documentation

◆ DeltaC9_KD()

gslpp::complex MPll::DeltaC9_KD ( double  q2)

Definition at line 1090 of file MPll.cpp.

1094 {

◆ funct_g()

gslpp::complex MPll::funct_g ( double  q2)

Definition at line 1082 of file MPll.cpp.

1086 {
1087  if (q2 < 4. * Mc * Mc)
1088  return -8. / 9. * log(Mc / Mb) + 8. / 27. + 16. / 9. * Mc * Mc / q2 - 4. / 9. * (2. + 4. * Mc * Mc / q2) * (sqrt(4. * Mc * Mc / q2 - 1.) * atan(1. / sqrt(4. * Mc * Mc / q2 - 1.)));

◆ getSigma()

double MPll::getSigma ( int  i,
double  q_2 
)

The value of \( \Sigma_{i} \) from \(q_{min}\) to \(q_{max}\).

Parameters
[in]iindex of the angular coefficient \( I_{i} \)
[in]

Definition at line 1211 of file MPll.cpp.

1216 {
1217  updateParameters();
1218 
1219  switch (i) {
1220  case 0:
1221  return getSigma1c(q_2);
1222  break;
1223  case 2:
1224  return getSigma2c(q_2);
1225  break;
1226  case 8:
1227  return getSigma6c(q_2);
1228  break;
1229  default:
1230  std::stringstream out;

◆ getwidth()

double MPll::getwidth ( )
inline

The width of the meson M.

Returns
\( \Gamma_M \)

Definition at line 262 of file MPll.h.

262  {
263  updateParameters();
264  return width;
265  }

◆ H_A()

gslpp::complex MPll::H_A ( double  q2)

The helicity amplitude \( H_A^{\lambda} \) .

Parameters
[in]q2\(q^2\) of the decay
Returns
\( H_A^{\lambda} \)

Definition at line 1110 of file MPll.cpp.

1114 {

◆ h_lambda()

gslpp::complex MPll::h_lambda ( double  q2)

The non-pertubative ccbar contributions to the helicity amplitudes.

Parameters
helhelicity
q2\(q^2\)
Returns
\(h_{hel}(q^2)\)

Definition at line 1096 of file MPll.cpp.

1100 {
1101  if (!dispersion) {
1102  return (twoMboMM * h_0 * T_L(q2) + h_1 * q2 / MM2 * V_L(q2)) / sixteenM_PI2 + h_2 * q2 * q2;
1103  } else {

◆ H_P()

gslpp::complex MPll::H_P ( double  q2)

The helicity amplitude \( H_P^{\lambda} \) .

Parameters
[in]q2\(q^2\) of the decay
Returns
\( H_P^{\lambda} \)

Definition at line 1120 of file MPll.cpp.

1124 {

◆ H_S()

gslpp::complex MPll::H_S ( double  q2)

The helicity amplitude \( H_S^{\lambda} \) .

Parameters
[in]q2\(q^2\) of the decay
Returns
\( H_S^{\lambda} \)

Definition at line 1115 of file MPll.cpp.

1119 {

◆ H_V()

gslpp::complex MPll::H_V ( double  q2)

The helicity amplitude \( H_V^{\lambda} \) .

Parameters
[in]q2\(q^2\) of the decay
Returns
\( H_V^{\lambda} \)

Definition at line 1105 of file MPll.cpp.

1109 {

◆ initializeMPllParameters()

std::vector< std::string > MPll::initializeMPllParameters ( )

A method for initializing the parameters necessary for MPll.

Returns
the vector of MPll specific parameters

Definition at line 76 of file MPll.cpp.

77 {
80 
81 #if NFPOLARBASIS_MPLL
83 #if LATTICE
84  << "b_0_fplus" << "b_1_fplus" << "b_2_fplus" << "m_fit_fplus_lat"
85  << "b_0_fT" << "b_1_fT" << "b_2_fT" << "m_fit_fT_lat"
86  << "b_0_f0" << "b_1_f0" << "b_2_f0" << "m_fit_f0_lat"
87 #else
88  << "r_1_fplus" << "r_2_fplus" << "m_fit2_fplus" << "r_1_fT" << "r_2_fT" << "m_fit2_fT" << "r_2_f0" << "m_fit2_f0"
89 #endif
90  << "absh_0_MP" << "argh_0_MP" << "absh_1_MP" << "argh_1_MP" << "absh_2_MP" << "argh_2_MP";
91 #else
93 #if LATTICE
94  << "b_0_fplus" << "b_1_fplus" << "b_2_fplus" << "m_fit_fplus_lat"
95  << "b_0_fT" << "b_1_fT" << "b_2_fT" << "m_fit_fT_lat"
96  << "b_0_f0" << "b_1_f0" << "b_2_f0" << "m_fit_f0_lat"
97 #else
98  << "r_1_fplus" << "r_2_fplus" << "m_fit2_fplus" << "r_1_fT" << "r_2_fT" << "m_fit2_fT" << "r_2_f0" << "m_fit2_f0"
99 #endif
100  << "reh_0_MP" << "imh_0_MP" << "reh_1_MP" << "imh_1_MP" << "reh_2_MP" << "imh_2_MP";
101 #endif
102  else {
103  std::stringstream out;
104  out << pseudoscalar;
105  throw std::runtime_error("MPll: pseudoscalar " + out.str() + " not implemented");
106  }
107 
108  if (dispersion) {
109  mpllParameters.clear();
111 #if LATTICE
112  << "b_0_fplus" << "b_1_fplus" << "b_2_fplus" << "m_fit_fplus_lat"
113  << "b_0_fT" << "b_1_fT" << "b_2_fT" << "m_fit_fT_lat"
114  << "b_0_f0" << "b_1_f0" << "b_2_f0" << "m_fit_f0_lat"
115 #else
116  << "r_1_fplus" << "r_2_fplus" << "m_fit2_fplus" << "r_1_fT" << "r_2_fT" << "m_fit2_fT" << "r_2_f0" << "m_fit2_f0"
117 #endif
118  << "r1_BK" << "r2_BK" << "deltaC9_BK" << "phDC9_BK";
119  }
120 
121  if (FixedWCbtos) mpllParameters.insert(mpllParameters.end(), { "C7_SM", "C9_SM", "C10_SM" });
124  return mpllParameters;
125 }

◆ integrateDelta()

double MPll::integrateDelta ( int  i,
double  q_min,
double  q_max 
)

The integral of \( \Delta_{i} \) from \(q_{min}\) to \(q_{max}\).

Parameters
[in]iindex of the angular coefficient \( I_{i} \)
[in]q_minminimum q^2 of the integral
[in]q_maxmaximum q^2 of the integral
Returns
\( <\Delta_{i}> \)

Definition at line 1232 of file MPll.cpp.

1232  : index " + out.str() + " not implemented");
1233  }
1234 }
1235 
1236 double MPll::integrateDelta(int i, double q_min, double q_max)
1237 {
1238  updateParameters();
1239 
1240  std::pair<double, double > qbin = std::make_pair(q_min, q_max);
1241 
1242  old_handler = gsl_set_error_handler_off();
1243 
1244  switch (i) {
1245  case 0:
1246  if (delta0Cached[qbin] == 0) {
1247  FD = convertToGslFunction(boost::bind(&MPll::getDelta1c, &(*this), _1));
1248  if (gsl_integration_cquad(&FD, q_min, q_max, 1.e-2, 1.e-1, w_delta, &avaDelta, &errDelta, NULL) != 0) return std::numeric_limits<double>::quiet_NaN();
1249  cacheDelta0[qbin] = NN*avaDelta;
1250  delta0Cached[qbin] = 1;
1251  }
1252  return cacheDelta0[qbin];
1253  break;
1254  case 2:
1255  if (delta2Cached[qbin] == 0) {
1256  FD = convertToGslFunction(boost::bind(&MPll::getDelta2c, &(*this), _1));
1257  if (gsl_integration_cquad(&FD, q_min, q_max, 1.e-2, 1.e-1, w_delta, &avaDelta, &errDelta, NULL) != 0) return std::numeric_limits<double>::quiet_NaN();
1258  cacheDelta2[qbin] = NN*avaDelta;
1259  delta2Cached[qbin] = 1;
1260  }
1261  return cacheDelta2[qbin];
1262  break;
1263  default:
1264  std::stringstream out;
1265  out << i;
1266  throw std::runtime_error("MPll::integrateDelta: index " + out.str() + " not implemented");
1267  }

◆ integrateSigma()

double MPll::integrateSigma ( int  i,
double  q_min,
double  q_max 
)

The integral of \( \Sigma_{i} \) from \(q_{min}\) to \(q_{max}\).

Parameters
[in]iindex of the angular coefficient \( I_{i} \)
[in]q_minminimum q^2 of the integral
[in]q_maxmaximum q^2 of the integral
Returns
\( <\Sigma_{i}> \)

Definition at line 1174 of file MPll.cpp.

1179 {
1180  updateParameters();
1181 
1182  std::pair<double, double > qbin = std::make_pair(q_min, q_max);
1183 
1184  old_handler = gsl_set_error_handler_off();
1185 
1186  switch (i) {
1187  case 0:
1188  if (sigma0Cached[qbin] == 0) {
1189  FS = convertToGslFunction(boost::bind(&MPll::getSigma1c, &(*this), _1));
1190  if (gsl_integration_cquad(&FS, q_min, q_max, 1.e-2, 1.e-1, w_sigma, &avaSigma, &errSigma, NULL) != 0) return std::numeric_limits<double>::quiet_NaN();
1191  cacheSigma0[qbin] = NN*avaSigma;
1192  sigma0Cached[qbin] = 1;
1193  }
1194  return cacheSigma0[qbin];
1195  break;
1196  case 2:
1197  if (sigma2Cached[qbin] == 0) {
1198  FS = convertToGslFunction(boost::bind(&MPll::getSigma2c, &(*this), _1));
1199  if (gsl_integration_cquad(&FS, q_min, q_max, 1.e-2, 1.e-1, w_sigma, &avaSigma, &errSigma, NULL) != 0) return std::numeric_limits<double>::quiet_NaN();
1200  cacheSigma2[qbin] = NN*avaSigma;
1201  sigma2Cached[qbin] = 1;
1202  }
1203  return cacheSigma2[qbin];
1204  break;
1205  default:
1206  std::stringstream out;
1207  out << i;
1208  throw std::runtime_error("MPll::integrateSigma: index " + out.str() + " not implemented");
1209  }

Member Data Documentation

◆ ale

double MPll::ale
private

alpha electromagnetic

Definition at line 287 of file MPll.h.

◆ dispersion

bool MPll::dispersion
private

Definition at line 282 of file MPll.h.

◆ FixedWCbtos

bool MPll::FixedWCbtos
private

Definition at line 283 of file MPll.h.

◆ GF

double MPll::GF
private

Fermi constant

Definition at line 286 of file MPll.h.

◆ lep

QCD::lepton MPll::lep
private

Final leptons type

Definition at line 276 of file MPll.h.

◆ Mb

double MPll::Mb
private

b quark mass

Definition at line 291 of file MPll.h.

◆ mb_pole

double MPll::mb_pole
private

b quark pole mass

Definition at line 295 of file MPll.h.

◆ Mc

double MPll::Mc
private

c quark mass

Definition at line 294 of file MPll.h.

◆ mc_pole

double MPll::mc_pole
private

c quark pole mass

Definition at line 296 of file MPll.h.

◆ meson

QCD::meson MPll::meson
private

Initial meson type

Definition at line 277 of file MPll.h.

◆ mJ2

double MPll::mJ2
private

Definition at line 284 of file MPll.h.

◆ Mlep

double MPll::Mlep
private

muon mass

Definition at line 288 of file MPll.h.

◆ MM

double MPll::MM
private

initial meson mass

Definition at line 289 of file MPll.h.

◆ MP

double MPll::MP
private

final pseudoscalar meson mass

Definition at line 290 of file MPll.h.

◆ mpllParameters

std::vector<std::string> MPll::mpllParameters
private

The string of mandatory MPll parameters

Definition at line 279 of file MPll.h.

◆ Ms

double MPll::Ms
private

s quark mass

Definition at line 297 of file MPll.h.

◆ mu_b

double MPll::mu_b
private

b mass scale

Definition at line 292 of file MPll.h.

◆ mu_h

double MPll::mu_h
private

\(\sqrt{\mu_b*\lambda_{QCD}}\)

Definition at line 293 of file MPll.h.

◆ myF_1

std::unique_ptr<F_1> MPll::myF_1
private

Definition at line 280 of file MPll.h.

◆ myF_2

std::unique_ptr<F_2> MPll::myF_2
private

Definition at line 281 of file MPll.h.

◆ mySM

const StandardModel& MPll::mySM
private

Model type

Definition at line 275 of file MPll.h.

◆ pseudoscalar

QCD::meson MPll::pseudoscalar
private

Final pseudoscalar meson type

Definition at line 278 of file MPll.h.

◆ spectator_charge

double MPll::spectator_charge
private

spectator quark charge

Definition at line 298 of file MPll.h.

◆ width

double MPll::width
private

initial meson width

Definition at line 299 of file MPll.h.


The documentation for this class was generated from the following files:
MPll::dispersion
bool dispersion
Definition: MPll.h:282
MPll::pseudoscalar
QCD::meson pseudoscalar
Definition: MPll.h:278
MPll::meson
QCD::meson meson
Definition: MPll.h:277
QCD::initializeMeson
void initializeMeson(QCD::meson meson_i) const
A method to initialize a meson.
Definition: QCD.cpp:236
make_vector
Definition: std_make_vector.h:15
gslpp::log
complex log(const complex &z)
Definition: gslpp_complex.cpp:342
h_0
A class for the correction in .
Definition: MVllObservables.h:1587
MPll::myF_1
std::unique_ptr< F_1 > myF_1
Definition: MPll.h:280
MPll::mpllParameters
std::vector< std::string > mpllParameters
Definition: MPll.h:279
MPll::integrateDelta
double integrateDelta(int i, double q_min, double q_max)
The integral of from to .
Definition: MPll.cpp:1232
F_2
Definition: F_2.h:15
MPll::lep
QCD::lepton lep
Definition: MPll.h:276
gslpp::sqrt
complex sqrt(const complex &z)
Definition: gslpp_complex.cpp:385
MPll::Mc
double Mc
Definition: MPll.h:294
MPll::FixedWCbtos
bool FixedWCbtos
Definition: MPll.h:283
MPll::Mb
double Mb
Definition: MPll.h:291
QCD::K_P
Definition: QCD.h:340
StandardModel::getFlavour
const Flavour & getFlavour() const
Definition: StandardModel.h:1020
F_1
Definition: F_1.h:15
MPll::mJ2
double mJ2
Definition: MPll.h:284
QCD::K_0
Definition: QCD.h:339
MPll::mySM
const StandardModel & mySM
Definition: MPll.h:275
MPll::width
double width
Definition: MPll.h:299
convertToGslFunction
gsl_function convertToGslFunction(const F &f)
Definition: gslpp_function_adapter.h:24
MPll::myF_2
std::unique_ptr< F_2 > myF_2
Definition: MPll.h:281
Flavour::getFlagUseDispersionRelation
bool getFlagUseDispersionRelation() const
Definition: Flavour.h:250
Flavour::getFlagFixedWCbtos
bool getFlagFixedWCbtos() const
Definition: Flavour.h:260