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 [Jager:2012uw], 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() [Ball:2004ye]. 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 173 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_nunu (double q2, QCD::lepton lep)
 The helicity amplitude \( H^{\lambda} \) for the invisible decay. 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...
 
double integrateSigmaTree (double q_min, double q_max)
 The integral of \( \Sigma_{tree} \) from \(q_{min}\) to \(q_{max}\) (arxiv/2301.06990) More...
 
 MPll (const StandardModel &SM_i, QCD::meson meson_i, QCD::meson pseudoscalar_i, QCD::lepton lep_i)
 Constructor. More...
 
double S_L (double q2)
 The helicity form factor \( S_L \). More...
 
gslpp::complex T_L (double q2)
 The helicity form factor \( T_L^0 \). More...
 
gslpp::complex V_L (double q2)
 The helicity form factor \( V_L^0 \). More...
 
virtual ~MPll ()
 Destructor. More...
 

Private Member Functions

gslpp::complex Cpar (double q2)
 The correction \( C_{\parallel} \) from [Beneke:2001at] . More...
 
gslpp::complex deltaC7_QCDF (double q2, bool spline=false)
 QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et. al (arXiv:0810.4077).. More...
 
gslpp::complex DeltaC9 (double q2)
 The total QCDF correction \( \Delta C_9 \) computed integrating over \( u \). More...
 
gslpp::complex deltaC9_QCDF (double q2, bool spline=false)
 QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et. al (arXiv:0810.4077).. More...
 
gslpp::complex deltaTpar (double q2)
 The total correction \( \Delta \mathcal{T}^{\parallel} \) from [Beneke:2001at] . More...
 
gslpp::complex F19 (double q2)
 The correction \( F_{19} \) from [Asatrian:2001de]. More...
 
gslpp::complex F27 (double q2)
 The correction \( F_{27} \) from [Asatrian:2001de]. More...
 
gslpp::complex F29 (double q2)
 The correction \( F_{29} \) from [Asatrian:2001de]. More...
 
gslpp::complex F87 (double q2)
 The correction \( F_{87} \) from [Asatrian:2001de]. More...
 
double F89 (double q2)
 The correction \( F_{89} \) from [Asatrian:2001de]. More...
 
gslpp::complex fDeltaC9 (double q2)
 The total QCDF correction \( \Delta C_9 \) computed fitting over \( u \). More...
 
void fit_DeltaC9_mumu ()
 The fitting routine for the QCDF correction \( \Delta C_9^0 \) in the muon channel. More...
 
double getDelta1c (double q2)
 The CP asymmetry \( \Delta_{1s} \) . More...
 
double getDelta2c (double q2)
 The CP asymmetry \( \Delta_{2s} \) . More...
 
double getintegratedSigmaTree ()
 The integral of \( \Sigma_{tree} \) from 0 to \(q_{cut}\). More...
 
double getSigma1c (double q2)
 The CP average \( \Sigma_{1s} \) . More...
 
double getSigma2c (double q2)
 The CP average \( \Sigma_{2s} \) . More...
 
double getSigma6c (double q2)
 The CP average \( \Sigma_{6s} \) . More...
 
gslpp::complex I1 (double u, double q2)
 The \( I_1 \) function from [Beneke:2001at] . More...
 
double imDC9fit (double *x, double *p)
 The fit function for the imaginary part of the QCDF correction \( \Delta C_9^{\lambda} \). More...
 
double Integrand_ImTpar_pm (double up)
 The sum of Integrand_ImTparplus() and Integrand_ImTparminus(). More...
 
double Integrand_ImTparminus (double up)
 The imaginary part of the integral involving \( T^{\parallel}_- \) at fixed \( q^2 \), according to [Beneke:2001at] . More...
 
double Integrand_ImTparplus (double up)
 The imaginary part of the integral involving \( T^{\parallel}_+ \) at fixed \( q^2 \), according to [Beneke:2001at] . More...
 
double Integrand_ReTpar_pm (double up)
 The sum of Integrand_ReTparplus() and Integrand_ReTparminus(). More...
 
double Integrand_ReTparminus (double up)
 The real part of the integral involving \( T^{\parallel}_- \) at fixed \( q^2 \), according to [Beneke:2001at] . More...
 
double Integrand_ReTparplus (double up)
 The real part of the integral involving \( T^{\parallel}_+ \) at fixed \( q^2 \), according to [Beneke:2001at] . More...
 
double reDC9fit (double *x, double *p)
 The fit function for the real part of the QCDF correction \( \Delta C_9^{\lambda} \). More...
 
double SigmaTree (double q2)
 
void spline_QCDF_func ()
 
gslpp::complex Tparminus (double u, double q2)
 The \( T^{\parallel}_- \) function from [Khodjamirian:2012rm] . More...
 
gslpp::complex Tparplus (double u, double q2)
 The \( T^{\parallel}_+ \) function from [Khodjamirian:2012rm] . More...
 

Private Attributes

double ale
 
double alpha_s_mub
 
bool dispersion = false
 
bool FixedWCbtos = false
 
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
 
bool MPll_DM_flag
 
bool MPll_GRvDV_flag
 
bool MPll_Lattice_flag
 
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
 
bool NeutrinoTree_flag
 
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 20 of file MPll.cpp.

21: mySM(SM_i), myF_1(new F_1()), myF_2(new F_2()),
22fplus_lat_cache(3, 0.),
23fT_lat_cache(3, 0.),
24f0_lat_cache(3, 0.),
25fplus_cache(2, 0.),
26fT_cache(2, 0.),
27k2_cache(2, 0.),
28SL_cache(2, 0.),
29N_cache(3, 0.),
30Ycache(2, 0.),
31H_V0cache(2, 0.),
32H_Scache(2, 0.),
33H_P_cache(4, 0.),
34Itree_cache(3, 0.),
35T_cache(5, 0.)
36{
37 lep = lep_i;
38 meson = meson_i;
39 pseudoscalar = pseudoscalar_i;
40 dispersion = true;
41 FixedWCbtos = false;
42 MPll_Lattice_flag = false;
43 MPll_GRvDV_flag = false;
44 MPll_DM_flag = false;
45 NeutrinoTree_flag = false;
46 mJ2 = 3.096 * 3.096;
47
48 I0_updated = 0;
49 I2_updated = 0;
50 I8_updated = 0;
51 Itree_updated = 0;
52
53 VL_updated = 0;
54 TL_updated = 0;
55 SL_updated = 0;
56
57 deltaTparpupdated = 0;
58 deltaTparmupdated = 0;
59
60 w_sigma = gsl_integration_cquad_workspace_alloc(100);
61 w_delta = gsl_integration_cquad_workspace_alloc(100);
62 w_sigmaTree = gsl_integration_cquad_workspace_alloc(100);
63 w_DTPPR = gsl_integration_cquad_workspace_alloc(100);
64
65 acc_Re_deltaC7_QCDF = gsl_interp_accel_alloc();
66 acc_Im_deltaC7_QCDF = gsl_interp_accel_alloc();
67 acc_Re_deltaC9_QCDF = gsl_interp_accel_alloc();
68 acc_Im_deltaC9_QCDF = gsl_interp_accel_alloc();
69
70 spline_Re_deltaC7_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
71 spline_Im_deltaC7_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
72 spline_Re_deltaC9_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
73 spline_Im_deltaC9_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
74}
Definition: F_1.h:15
Definition: F_2.h:15
bool FixedWCbtos
Definition: MPll.h:323
bool MPll_GRvDV_flag
Definition: MPll.h:325
QCD::meson pseudoscalar
Definition: MPll.h:318
std::unique_ptr< F_2 > myF_2
Definition: MPll.h:321
bool MPll_Lattice_flag
Definition: MPll.h:324
bool MPll_DM_flag
Definition: MPll.h:326
QCD::meson meson
Definition: MPll.h:317
const StandardModel & mySM
Definition: MPll.h:315
QCD::lepton lep
Definition: MPll.h:316
double mJ2
Definition: MPll.h:328
std::unique_ptr< F_1 > myF_1
Definition: MPll.h:320
bool NeutrinoTree_flag
Definition: MPll.h:327
bool dispersion
Definition: MPll.h:322

◆ ~MPll()

MPll::~MPll ( )
virtual

Destructor.

Definition at line 76 of file MPll.cpp.

77{
78}

Member Function Documentation

◆ Cpar()

gslpp::complex MPll::Cpar ( double  q2)
private

The correction \( C_{\parallel} \) from [Beneke:2001at] .

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

Definition at line 1109 of file MPll.cpp.

1110{
1111 return -(C_2L_bar * F27(q2) + C_8L * F87(q2) + MM / 2. / Mb *
1112 (C_2L_bar * F29(q2) + 2. * C_1L_bar * (F19(q2) + F29(q2) / 6.) + C_8L * F89(q2)));
1113}
double Mb
Definition: MPll.h:335
double F89(double q2)
The correction from .
Definition: MPll.cpp:1102
gslpp::complex F87(double q2)
The correction from .
Definition: MPll.cpp:1095
gslpp::complex F19(double q2)
The correction from .
Definition: MPll.cpp:1063
double MM
Definition: MPll.h:333
gslpp::complex F27(double q2)
The correction from .
Definition: MPll.cpp:1077
gslpp::complex F29(double q2)
The correction from .
Definition: MPll.cpp:1086

◆ deltaC7_QCDF()

gslpp::complex MPll::deltaC7_QCDF ( double  q2,
bool  spline = false 
)
private

QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et. al (arXiv:0810.4077)..

Parameters
q2\(q^2\) of the decay
Returns
\( \Delta C_{7}^{QCDF} \)

Definition at line 1198 of file MPll.cpp.

1199{
1200#if SPLINE
1201 if (spline) return gsl_spline_eval(spline_Re_deltaC7_QCDF, q2, acc_Re_deltaC7_QCDF);
1202#endif
1203
1204 double muh = mu_b / mb_pole;
1205 double z = mc_pole * mc_pole / mb_pole / mb_pole;
1206 double sh = q2 / mb_pole / mb_pole;
1207 double sh2 = sh*sh;
1208
1209 //gslpp::complex A_Sdl = A_Seidel(q2, mb_pole*mb_pole); /* hep-ph/0403185v2.*/
1210 //gslpp::complex Fu_17 = -A_Sdl; /* sign different from hep-ph/0403185v2 but consistent with hep-ph/0412400 */
1211 //gslpp::complex Fu_27 = 6. * A_Sdl; /* sign different from hep-ph/0403185v2 but consistent with hep-ph/0412400 */
1212 gslpp::complex F_17 = myF_1->F_17re(muh, z, sh, 20) + gslpp::complex::i() * myF_1->F_17im(muh, z, sh, 20); /*q^2 = 0 gives nan. Independent of how small q^2 is. arXiv:0810.4077*/
1213 gslpp::complex F_27 = myF_2->F_27re(muh, z, sh, 20) + gslpp::complex::i() * myF_2->F_27im(muh, z, sh, 20); /*q^2 = 0 gives nan. Independent of how small q^2 is. arXiv:0810.4077*/
1214 gslpp::complex F_87 = F87_0 + F87_1 * sh + F87_2 * sh2 + F87_3 * sh * sh2 - 8. / 9. * log(sh) * (sh + sh2 + sh * sh2);
1215
1216 gslpp::complex delta = C_1 * F_17 + C_2 * F_27;
1217 gslpp::complex delta_t = C_8 * F_87 + delta;
1218 //gslpp::complex delta_u = delta + C_1 * Fu_17 + C_2 * Fu_27;
1219
1220 return -alpha_s_mub / (4. * M_PI) * (delta_t /*- lambda_u / lambda_t * delta_u */);
1221}
double mc_pole
Definition: MPll.h:340
double alpha_s_mub
Definition: MPll.h:344
double mb_pole
Definition: MPll.h:339
double mu_b
Definition: MPll.h:336

◆ DeltaC9()

gslpp::complex MPll::DeltaC9 ( double  q2)
private

The total QCDF correction \( \Delta C_9 \) computed integrating over \( u \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \Delta C_9 \)

Definition at line 1193 of file MPll.cpp.

1194{
1195 return deltaTpar(q2);
1196}
gslpp::complex deltaTpar(double q2)
The total correction from .
Definition: MPll.cpp:1115

◆ DeltaC9_KD()

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

Definition at line 1315 of file MPll.cpp.

1316{
1317 return ((Delta_C9 + r_1 * q2 / mJ2) / (1. - r_2 * q2 / mJ2) - (3. * (-0.267) + 1.117) * funct_g(q2))*exp_Phase;
1318 /* C_1 = -0.267 and C_2 = 1.117 in KMPW */
1319}
gslpp::complex funct_g(double q2)
Definition: MPll.cpp:1307

◆ deltaC9_QCDF()

gslpp::complex MPll::deltaC9_QCDF ( double  q2,
bool  spline = false 
)
private

QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et. al (arXiv:0810.4077)..

Parameters
q2\(q^2\) of the decay
Returns
\( \Delta C_{9}^{QCDF} \)

Definition at line 1223 of file MPll.cpp.

1224{
1225#if SPLINE
1226 if (spline) return gsl_spline_eval(spline_Re_deltaC9_QCDF, q2, acc_Re_deltaC9_QCDF);
1227#endif
1228 double muh = mu_b / mb_pole;
1229 double z = mc_pole * mc_pole / mb_pole / mb_pole;
1230 double sh = q2 / mb_pole / mb_pole;
1231 double sh2 = sh*sh;
1232
1233 //gslpp::complex B_Sdl = B_Seidel(q2, mb_pole*mb_pole); /* hep-ph/0403185v2.*/
1234 //gslpp::complex C_Sdl = C_Seidel(q2); /* hep-ph/0403185v2.*/
1235 //gslpp::complex Fu_19 = -(B_Sdl + 4. * C_Sdl); /* sign different from hep-ph/0403185v2 but consistent with hep-ph/0412400 */
1236 //gslpp::complex Fu_29 = -(-6. * B_Sdl + 3. * C_Sdl); /* sign different from hep-ph/0403185v2 but consistent with hep-ph/0412400 */
1237 gslpp::complex F_19 = myF_1->F_19re(muh, z, sh, 20) + gslpp::complex::i() * myF_1->F_19im(muh, z, sh, 20); /*q^2 = 0 gives nan. Independent of how small q^2 is. arXiv:0810.4077*/
1238 gslpp::complex F_29 = myF_2->F_29re(muh, z, sh, 20) + gslpp::complex::i() * myF_2->F_29im(muh, z, sh, 20); /*q^2 = 0 gives nan. Independent of how small q^2 is. arXiv:0810.4077*/
1239 gslpp::complex F_89 = (F89_0 + F89_1 * sh + F89_2 * sh2 + F89_3 * sh * sh2 + 16. / 9. * log(sh) * (1. + sh + sh2 + sh * sh2));
1240
1241 gslpp::complex delta = C_1 * F_19 + C_2 * F_29;
1242 gslpp::complex delta_t = C_8 * F_89 + delta;
1243 //gslpp::complex delta_u = delta + C_1 * Fu_19 + C_2 * Fu_29;
1244
1245 return -alpha_s_mub / (4. * M_PI) * (delta_t /*- lambda_u / lambda_t * delta_u*/);
1246}

◆ deltaTpar()

gslpp::complex MPll::deltaTpar ( double  q2)
private

The total correction \( \Delta \mathcal{T}^{\parallel} \) from [Beneke:2001at] .

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \Delta \mathcal{T}^{\parallel} \)

Definition at line 1115 of file MPll.cpp.

1116{
1117 tmpq2 = q2;
1118
1119 //old_handler = gsl_set_error_handler_off();
1120
1121 if (deltaTparpCached[q2] == 0) {
1122
1123 DTPPR = convertToGslFunction(bind(&MPll::Integrand_ReTpar_pm, &(*this), _1));
1124 if (gsl_integration_cquad(&DTPPR, 0., 1., 1.e-2, 1.e-1, w_DTPPR, &avaDTPPR, &errDTPPR, NULL) != 0) return std::numeric_limits<double>::quiet_NaN();
1125 double ReTppint = avaDTPPR;
1126
1127 DTPPR = convertToGslFunction(bind(&MPll::Integrand_ImTpar_pm, &(*this), _1));
1128 if (gsl_integration_cquad(&DTPPR, 0., 1., 1.e-2, 1.e-1, w_DTPPR, &avaDTPPR, &errDTPPR, NULL) != 0) return std::numeric_limits<double>::quiet_NaN();
1129 double ImTppint = avaDTPPR;
1130
1131 cacheDeltaTparp[q2] = (ReTppint + gslpp::complex::i() * ImTppint);
1132 deltaTparpCached[q2] = 1;
1133 }
1134
1135 //gsl_set_error_handler(old_handler);
1136
1137 return deltaT_0 * Cpar(q2) + deltaT_1par * cacheDeltaTparp[q2] / f_plus(q2);
1138}
double Integrand_ReTpar_pm(double up)
The sum of Integrand_ReTparplus() and Integrand_ReTparminus().
Definition: MPll.cpp:1058
gslpp::complex Cpar(double q2)
The correction from .
Definition: MPll.cpp:1109
double Integrand_ImTpar_pm(double up)
The sum of Integrand_ImTparplus() and Integrand_ImTparminus().
Definition: MPll.cpp:1053

◆ F19()

gslpp::complex MPll::F19 ( double  q2)
private

The correction \( F_{19} \) from [Asatrian:2001de].

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

Definition at line 1063 of file MPll.cpp.

1064{
1065 double s = q2 / Mb2;
1066 double s2 = s*s;
1067 double Ls = log(s);
1068 double Lc = log(Mc / Mb);
1069 double Lm = log(mu_b / Mb);
1070 gslpp::complex i = gslpp::complex::i();
1071 return (-1424. / 729. + 16. / 243. * i * M_PI + 64. / 27. * Lc)*Lm - 16. / 243. * Lm * Ls + (16. / 1215. - 32. / 135. / Mc2 * Mb2) * Lm * s
1072 + (4. / 2835. - 8. / 315. / Mc2 * Mb2 / Mc2 * Mb2) * Lm * s2 + (16. / 76545. - 32. / 8505 / Mc2 * Mb2 / Mc2 * Mb2 / Mc2 * Mb2) * Lm * s * s2 - 256. / 243. * Lm * Lm
1073 + (-11.65 + 0.18223 * i + (-24.709 - 0.13474 * i) * s + (-43.588 - 0.4738 * i) * s2 + (-86.22 - 1.3542 * i) * s * s2
1074 + (-0.080959 - 0.051864 * i + (-0.036585 + 0.024753 * i) * s + (-0.021692 + 0.036925 * i) * s2 + (0.013282 + 0.052023 * i) * s * s2) * Ls);
1075}
double Mc
Definition: MPll.h:338
Test Observable.

◆ F27()

gslpp::complex MPll::F27 ( double  q2)
private

The correction \( F_{27} \) from [Asatrian:2001de].

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

Definition at line 1077 of file MPll.cpp.

1078{
1079 double s = q2 / Mb2;
1080 double s2 = s*s;
1081 double Ls = log(s);
1082 gslpp::complex i = gslpp::complex::i();
1083 return F27_0 + F27_1 * s + F27_2 * s2 + F27_3 * s * s2 + F27_L1_1 * Ls * s + F27_L1_2 * Ls * s2 + F27_L1_3 * Ls * s * s2;
1084}

◆ F29()

gslpp::complex MPll::F29 ( double  q2)
private

The correction \( F_{29} \) from [Asatrian:2001de].

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

Definition at line 1086 of file MPll.cpp.

1087{
1088 double s = q2 / Mb2;
1089 double s2 = s*s;
1090 double Ls = log(s);
1091 gslpp::complex i = gslpp::complex::i();
1092 return F29_0 + F29_L1 * Ls + F29_1 * s + F29_2 * s2 + F29_3 * s * s2 + F29_L1_1 * Ls * s + F29_L1_2 * Ls * s2 + F29_L1_3 * Ls * s2 *s;
1093}

◆ F87()

gslpp::complex MPll::F87 ( double  q2)
private

The correction \( F_{87} \) from [Asatrian:2001de].

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

Definition at line 1095 of file MPll.cpp.

1096{
1097 double s = q2 / Mb2;
1098 double s2 = s*s;
1099 return F87_0 + F87_1 * s + F87_2 * s2 + F87_3 * s * s2 - 0.888889 * log(s)*(1. + s + s2 + s * s2);
1100}

◆ F89()

double MPll::F89 ( double  q2)
private

The correction \( F_{89} \) from [Asatrian:2001de].

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

Definition at line 1102 of file MPll.cpp.

1103{
1104 double s = q2 / Mb2;
1105 double s2 = s*s;
1106 return F89_0 + F89_1 * s + F89_2 * s2 + F89_3 * s * s2 + 1.77778 * log(s)*(1. + s + s2 + s * s2);
1107}

◆ fDeltaC9()

gslpp::complex MPll::fDeltaC9 ( double  q2)
private

The total QCDF correction \( \Delta C_9 \) computed fitting over \( u \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \Delta C_9 \)

Definition at line 1184 of file MPll.cpp.

1185{
1186 if (q2 < MPllSWITCH) return (reDC9fit(&q2, const_cast<double *> (refres->GetParams()))
1187 + gslpp::complex::i() * imDC9fit(&q2, const_cast<double *> (imfres->GetParams())));
1188 else return (reDC9fit(&q2, const_cast<double *> (refres->GetParams()))
1189 + gslpp::complex::i() * imDC9fit(&q2, const_cast<double *> (imfres->GetParams()))) / q2;
1190
1191}
double imDC9fit(double *x, double *p)
The fit function for the imaginary part of the QCDF correction .
Definition: MPll.cpp:1145
double reDC9fit(double *x, double *p)
The fit function for the real part of the QCDF correction .
Definition: MPll.cpp:1140

◆ fit_DeltaC9_mumu()

void MPll::fit_DeltaC9_mumu ( )
private

The fitting routine for the QCDF correction \( \Delta C_9^0 \) in the muon channel.

Definition at line 1152 of file MPll.cpp.

1153{
1154 int dim = 0;
1155 for (double i = 0.1; i < MPllSWITCH; i += 0.4) {
1156 double q2tmp = i;
1157 myq2.push_back(q2tmp);
1158 ReDeltaC9.push_back((deltaTpar(q2tmp)).real());
1159 ImDeltaC9.push_back((deltaTpar(q2tmp)).imag());
1160 dim++;
1161 }
1162 for (double i = MPllSWITCH; i < 8.2; i += 0.4) {
1163 double q2tmp = i;
1164 myq2.push_back(q2tmp);
1165 ReDeltaC9.push_back(q2tmp * (deltaTpar(q2tmp)).real());
1166 ImDeltaC9.push_back(q2tmp * (deltaTpar(q2tmp)).imag());
1167 dim++;
1168 }
1169
1170 gr1 = TGraph(dim, myq2.data(), ReDeltaC9.data());
1171 gr2 = TGraph(dim, myq2.data(), ImDeltaC9.data());
1172
1173 reffit = TF1("reffit", this, &MPll::reDC9fit, 0, 8.1, 7);
1174 imffit = TF1("imffit", this, &MPll::imDC9fit, 0, 8.1, 8);
1175
1176 refres = gr1.Fit(&reffit, "SQN0+rob=0.99");
1177 imfres = gr2.Fit(&imffit, "SQN0+rob=0.99");
1178
1179 ReDeltaC9.clear();
1180 ImDeltaC9.clear();
1181 myq2.clear();
1182}

◆ funct_g()

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

Definition at line 1307 of file MPll.cpp.

1308{
1309 if (q2 < 4. * Mc * Mc)
1310 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.)));
1311 else
1312 return -8. / 9. * log(Mc / Mb) + 8. / 27. + 16. / 9. * Mc * Mc / q2 - 4. / 9. * (2. + 4. * Mc * Mc / q2) * (sqrt(1. - 4. * Mc * Mc / q2) * (log(1. + sqrt(1. - 4. * Mc * Mc / q2) / sqrt(4. * Mc * Mc / q2)) - gslpp::complex::i() * M_PI_2));
1313}

◆ getDelta1c()

double MPll::getDelta1c ( double  q2)
inlineprivate

The CP asymmetry \( \Delta_{1s} \) .

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \Delta_{1c} \)

Definition at line 1046 of file MPll.h.

1047 {
1048 return Delta(0, q2);
1049 };

◆ getDelta2c()

double MPll::getDelta2c ( double  q2)
inlineprivate

The CP asymmetry \( \Delta_{2s} \) .

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

Definition at line 1056 of file MPll.h.

1057 {
1058 return Delta(2, q2);
1059 };

◆ getintegratedSigmaTree()

double MPll::getintegratedSigmaTree ( )
private

The integral of \( \Sigma_{tree} \) from 0 to \(q_{cut}\).

Returns
\( <\Sigma{tree}> \)

Definition at line 1564 of file MPll.cpp.

1565{
1566 return mySM.getMesons(meson).getLifetime() / HCUT * GF4 * VusVub_abs2 * fP2 * fM2 / (128. * M_PI2 * MM3 * Gammatau) * mtau * (mtau2 - MP2) * (mtau2 - MP2) * (MM2 - mtau2) * (MM2 - mtau2);
1567}
double getLifetime() const
A get method for the lifetime of the meson.
Definition: Meson.h:351
const Meson & getMesons(const QCD::meson m) const
A get method to access a meson as an object of the type Meson.
Definition: QCD.h:526

◆ 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]

_form#752 value of the function

Returns
\( <\Sigma_{i}> \)

Definition at line 1471 of file MPll.cpp.

1472{
1473 updateParameters();
1474
1475 switch (i) {
1476 case 0:
1477 return getSigma1c(q_2);
1478 break;
1479 case 2:
1480 return getSigma2c(q_2);
1481 break;
1482 case 8:
1483 return getSigma6c(q_2);
1484 break;
1485 default:
1486 std::stringstream out;
1487 out << i;
1488 throw std::runtime_error("MPll::getSigma: index " + out.str() + " not implemented");
1489 }
1490}
double getSigma1c(double q2)
The CP average .
Definition: MPll.h:1016
double getSigma2c(double q2)
The CP average .
Definition: MPll.h:1026
double getSigma6c(double q2)
The CP average .
Definition: MPll.h:1036

◆ getSigma1c()

double MPll::getSigma1c ( double  q2)
inlineprivate

The CP average \( \Sigma_{1s} \) .

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \Sigma_{1c} \)

Definition at line 1016 of file MPll.h.

1017 {
1018 return I_1c(q2);
1019 };

◆ getSigma2c()

double MPll::getSigma2c ( double  q2)
inlineprivate

The CP average \( \Sigma_{2s} \) .

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

Definition at line 1026 of file MPll.h.

1027 {
1028 return I_2c(q2);
1029 };

◆ getSigma6c()

double MPll::getSigma6c ( double  q2)
inlineprivate

The CP average \( \Sigma_{6s} \) .

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \Sigma_{6c} \)

Definition at line 1036 of file MPll.h.

1037 {
1038 return I_6c(q2);
1039 };

◆ getwidth()

double MPll::getwidth ( )
inline

The width of the meson M.

Returns
\( \Gamma_M \)

Definition at line 302 of file MPll.h.

302 {
303 updateParameters();
304 return width;
305 }
double width
Definition: MPll.h:343

◆ 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 1345 of file MPll.cpp.

1346{
1347 return (-C_10 + etaP * pow(-1, angmomP) * C_10p) *V_L(q2);
1348}
gslpp::complex V_L(double q2)
The helicity form factor .
Definition: MPll.cpp:954

◆ 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 1321 of file MPll.cpp.

1322{
1323 if (!dispersion) {
1324// if (q2 <= 1.) return 1.3e-4/MM2 * (1. + gslpp::complex::i()) / sqrt(2.);
1325// else {
1326// h_2 = (-sixteenM_PI2*4.9e-7/MM2 * (1. + gslpp::complex::i()) / sqrt(2.) - h_1 / MM2 * V_L(1.) - h_2 * sixteenM_PI2) / twoMboMM / T_L(1.);
1327// h_2 = 1.3e-4/MM2 * (1. + gslpp::complex::i()) / sqrt(2.) - h_1 * V_L(1.)/sixteenM_PI2/MM2;
1328 //else return 4.9e-7/MM2 + h_1 * (q2 * V_L(q2) - T_L(q2). * V_L(1)/T_L(1.)) / MM2 / sixteenM_PI2;
1329
1330 // return (twoMboMM * h_0 * T_L(q2) + h_1 * q2 / MM2 * V_L(q2)) / sixteenM_PI2 + h_2 * q2 * sqrt(q2);
1331 return (h_1 + h_2 * q2 ) * sqrt(q2) * sqrt(lambda(q2)) / MM3 / 2.;
1332// }
1333 } else {
1334 return -q2 / (MM2 * sixteenM_PI2) * V_L(q2) * DeltaC9_KD(q2);
1335 }
1336}
gslpp::complex DeltaC9_KD(double q2)
Definition: MPll.cpp:1315

◆ H_nunu()

gslpp::complex MPll::H_nunu ( double  q2,
QCD::lepton  lep 
)

The helicity amplitude \( H^{\lambda} \) for the invisible decay.

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

Definition at line 1360 of file MPll.cpp.

1361{
1362 if (lep == QCD::NEUTRINO_1)
1363 return -(C_L_nunu_e - etaP * pow(-1, angmomP) * C_R_nunu_e) * V_L(q2);
1364 else if (lep == QCD::NEUTRINO_2)
1365 return -(C_L_nunu_mu - etaP * pow(-1, angmomP) * C_R_nunu_mu) * V_L(q2);
1366 else if (lep == QCD::NEUTRINO_3)
1367 return -(C_L_nunu_tau - etaP * pow(-1, angmomP) * C_R_nunu_tau) * V_L(q2);
1368 else throw std::runtime_error("MPll::H_nunu: lepton not supported");
1369}
@ NEUTRINO_2
Definition: QCD.h:313
@ NEUTRINO_1
Definition: QCD.h:311
@ NEUTRINO_3
Definition: QCD.h:315

◆ 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 1355 of file MPll.cpp.

1356{
1357 return ( MboMW * (C_P - etaP * pow(-1, angmomP) * C_Pp) + twoMlepMb / q2 * (C_10 * (1. + etaP * pow(-1, angmomP) * MsoMb) - C_10p * (etaP * pow(-1, angmomP) + MsoMb))) * S_L(q2);
1358}
double S_L(double q2)
The helicity form factor .
Definition: MPll.cpp:964

◆ 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 1350 of file MPll.cpp.

1351{
1352 return MboMW * (C_S - etaP * pow(-1, angmomP) * C_Sp) * S_L(q2);
1353}

◆ 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 1338 of file MPll.cpp.

1339{
1340 return -((C_9 + deltaC9_QCDF(q2, SPLINE) + Y(q2) /*+ fDeltaC9(q2)*/ - etaP * pow(-1, angmomP) * C_9p) * V_L(q2)
1341 + MM2 / q2 * (twoMboMM * (C_7 + deltaC7_QCDF(q2, SPLINE) - etaP * pow(-1, angmomP) * C_7p) * T_L(q2)
1342 - sixteenM_PI2 * h_lambda(q2)));
1343}
gslpp::complex h_lambda(double q2)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MPll.cpp:1321
gslpp::complex deltaC9_QCDF(double q2, bool spline=false)
QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et....
Definition: MPll.cpp:1223
gslpp::complex T_L(double q2)
The helicity form factor .
Definition: MPll.cpp:959
gslpp::complex deltaC7_QCDF(double q2, bool spline=false)
QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et....
Definition: MPll.cpp:1198

◆ I1()

gslpp::complex MPll::I1 ( double  u,
double  q2 
)
private

The \( I_1 \) function from [Beneke:2001at] .

Parameters
[in]udummy variable to be integrated out
[in]q2\(q^2\) of the decay
Returns
\( I_1 \)

Definition at line 973 of file MPll.cpp.

974{
975 std::pair<double, double > uq2 = std::make_pair(u, q2);
976
977 if (I1Cached[uq2] == 0) {
978 ubar = 1. - u;
979 xp = .5 + sqrt(0.25 - (Mc2 - gslpp::complex::i()*1.e-10) / (ubar * MM2 + u * q2));
980 xm = .5 - sqrt(0.25 - (Mc2 - gslpp::complex::i()*1.e-10) / (ubar * MM2 + u * q2));
981 yp = .5 + sqrt(0.25 - (Mc2 - gslpp::complex::i()*1.e-10) / q2);
982 ym = .5 - sqrt(0.25 - (Mc2 - gslpp::complex::i()*1.e-10) / q2);
983 L1xp = log(1. - 1. / xp) * log(1. - xp) - M_PI2osix + dilog(xp / (xp - 1.));
984 L1xm = log(1. - 1. / xm) * log(1. - xm) - M_PI2osix + dilog(xm / (xm - 1.));
985 L1yp = log(1. - 1. / yp) * log(1. - yp) - M_PI2osix + dilog(yp / (yp - 1.));
986 L1ym = log(1. - 1. / ym) * log(1. - ym) - M_PI2osix + dilog(ym / (ym - 1.));
987
988 cacheI1[uq2] = 1. + twoMc2 / ubar / (MM2 - q2)*(L1xp + L1xm - L1yp - L1ym);
989 I1Cached[uq2] = 1;
990 }
991
992 return cacheI1[uq2];
993}

◆ imDC9fit()

double MPll::imDC9fit ( double *  x,
double *  p 
)
private

The fit function for the imaginary part of the QCDF correction \( \Delta C_9^{\lambda} \).

Parameters
[in]xfit variable
[in]pfit parameters vector
Returns
\( f_{Im \Delta C_9^{\lambda}} \)

Definition at line 1145 of file MPll.cpp.

1146{
1147 return p[0] / x[0] + p[1] + p[2] * x[0] + p[3] * x[0] * x[0] + p[4] * x[0] * x[0] * x[0] + p[5] * x[0] * x[0] * x[0] * x[0] + p[6] * x[0] * x[0] * x[0] * x[0] * x[0] + p[7] * x[0] * x[0] * x[0] * x[0] * x[0] * x[0];
1148
1149 //double thr = 4.*Mc2;
1150}

◆ 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 80 of file MPll.cpp.

81{
88
90 if (MPll_Lattice_flag) mpllParameters = make_vector<std::string>()
91 << "b_0_fplus" << "b_1_fplus" << "b_2_fplus" << "m_fit_fplus_lat"
92 << "b_0_fT" << "b_1_fT" << "b_2_fT" << "m_fit_fT_lat"
93 << "b_0_f0" << "b_1_f0" << "b_2_f0" << "m_fit_f0_lat" ;
94 else if (MPll_GRvDV_flag) mpllParameters = make_vector<std::string>()
95 << "b_0_fplus" << "b_1_fplus" << "b_2_fplus" << "m_fit_fplus_lat"
96 << "b_0_fT" << "b_1_fT" << "b_2_fT" << "m_fit_fT_lat"
97 << "b_1_f0" << "b_2_f0" << "m_fit_f0_lat" ;
98 else if (MPll_DM_flag) mpllParameters = make_vector<std::string>()
99 << "b_0_fplus" << "b_1_fplus" << "b_2_fplus" << "m_fit_fplus_lat"
100 << "b_0_fT" << "b_1_fT" << "b_2_fT" << "m_fit_fT_lat"
101 << "b_1_f0" << "b_2_f0" ;
102 else mpllParameters = make_vector<std::string>()
103 << "r_1_fplus" << "r_2_fplus" << "m_fit2_fplus" << "r_1_fT" << "r_2_fT" << "m_fit2_fT" << "r_2_f0" << "m_fit2_f0";
104 } else {
105 std::stringstream out;
106 out << pseudoscalar;
107 throw std::runtime_error("MPll: pseudoscalar " + out.str() + " not implemented");
108 }
109
110 if (lep != QCD::NEUTRINO_1){
111 if (dispersion) {
112 mpllParameters.insert(mpllParameters.end(), { "r1_BK", "r2_BK", "deltaC9_BK", "phDC9_BK" });
113 } else {
114#if NFPOLARBASIS_MPLL
115 mpllParameters.insert(mpllParameters.end(), { "absh_0_MP", "argh_0_MP", "absh_1_MP", "argh_1_MP", "absh_2_MP", "argh_2_MP", "Delta_C7_U", "Delta_C9_U" });
116#else
117 mpllParameters.insert(mpllParameters.end(), { "reh_0_MP", "imh_0_MP", "reh_1_MP", "imh_1_MP", "reh_2_MP", "imh_2_MP", "Delta_C7_U", "Delta_C9_U" });
118#endif
119 }
120 }
121
122 if (FixedWCbtos)
123 if (lep != QCD::NEUTRINO_1) mpllParameters.insert(mpllParameters.end(), { "C7_SM", "C9_SM", "C10_SM" });
124 else mpllParameters.insert(mpllParameters.end(), { "CLnunu_SM" });
127 return mpllParameters;
128}
bool getFlagUseDispersionRelation() const
Definition: Flavour.h:343
bool getFlagMPll_DM() const
Definition: Flavour.h:375
bool getFlagFixedWCbtos() const
Definition: Flavour.h:363
bool getFlagNeutrinoTree() const
Definition: Flavour.h:383
bool getFlagMPll_FNALMILC() const
Definition: Flavour.h:367
bool getFlagMPll_GRvDV() const
Definition: Flavour.h:371
std::vector< std::string > mpllParameters
Definition: MPll.h:319
@ K_P
Definition: QCD.h:340
@ K_0
Definition: QCD.h:339
void initializeMeson(QCD::meson meson_i) const
A method to initialize a meson.
Definition: QCD.cpp:280
const Flavour & getFlavour() const

◆ Integrand_ImTpar_pm()

double MPll::Integrand_ImTpar_pm ( double  up)
private

The sum of Integrand_ImTparplus() and Integrand_ImTparminus().

Parameters
[in]updummy variable to be integrated out
Returns
\( Im T^{\parallel}_+ \Phi^{\parallel} + Im T^{\parallel}_- \Phi^{\parallel}\)

Definition at line 1053 of file MPll.cpp.

1054{
1056}
double Integrand_ImTparminus(double up)
The imaginary part of the integral involving at fixed , according to .
Definition: MPll.cpp:1042
double Integrand_ImTparplus(double up)
The imaginary part of the integral involving at fixed , according to .
Definition: MPll.cpp:1023

◆ Integrand_ImTparminus()

double MPll::Integrand_ImTparminus ( double  up)
private

The imaginary part of the integral involving \( T^{\parallel}_- \) at fixed \( q^2 \), according to [Beneke:2001at] .

Parameters
[in]updummy variable to be integrated out
Returns
\( Im T^{\parallel}_- \Phi^{\parallel}\)

Definition at line 1042 of file MPll.cpp.

1043{
1044 double Lambdaplus = mySM.getMesons(meson).getLambdaM();
1045 gslpp::complex Lambdamin = exp(-tmpq2 / MM / Lambdaplus) / Lambdaplus * (-gsl_sf_expint_Ei(tmpq2 / MM / Lambdaplus) + gslpp::complex::i() * M_PI);
1046
1047 double u = up;
1048 return ((Tparminus(u, tmpq2)*6. * u * (1. - u)*
1049 (1 + threeGegen0 * (2. * u - 1)
1050 + threeGegen1otwo * ((10. * u - 5.)*(2. * u - 1.) - 1.))) / Lambdamin).imag();
1051}
gslpp::complex Tparminus(double u, double q2)
The function from .
Definition: MPll.cpp:1008
const double & getLambdaM() const
Definition: Meson.h:402

◆ Integrand_ImTparplus()

double MPll::Integrand_ImTparplus ( double  up)
private

The imaginary part of the integral involving \( T^{\parallel}_+ \) at fixed \( q^2 \), according to [Beneke:2001at] .

Parameters
[in]updummy variable to be integrated out
Returns
\( Im T^{\parallel}_+ \Phi^{\parallel}\)

Definition at line 1023 of file MPll.cpp.

1024{
1025 double u = up;
1026 return ((Tparplus(u, tmpq2)*6. * u * (1. - u)*
1027 (1 + threeGegen0 * (2. * u - 1)
1028 + threeGegen1otwo * ((10. * u - 5.)*(2. * u - 1.) - 1.))) / mySM.getMesons(meson).getLambdaM()).imag();
1029}
gslpp::complex Tparplus(double u, double q2)
The function from .
Definition: MPll.cpp:995

◆ Integrand_ReTpar_pm()

double MPll::Integrand_ReTpar_pm ( double  up)
private

The sum of Integrand_ReTparplus() and Integrand_ReTparminus().

Parameters
[in]updummy variable to be integrated out
Returns
\( Re T^{\parallel}_+ \Phi^{\parallel} + Re T^{\parallel}_- \Phi^{\parallel}\)

Definition at line 1058 of file MPll.cpp.

1059{
1061}
double Integrand_ReTparminus(double up)
The real part of the integral involving at fixed , according to .
Definition: MPll.cpp:1031
double Integrand_ReTparplus(double up)
The real part of the integral involving at fixed , according to .
Definition: MPll.cpp:1015

◆ Integrand_ReTparminus()

double MPll::Integrand_ReTparminus ( double  up)
private

The real part of the integral involving \( T^{\parallel}_- \) at fixed \( q^2 \), according to [Beneke:2001at] .

Parameters
[in]updummy variable to be integrated out
Returns
\( Re T^{\parallel}_- \Phi^{\parallel}\)

Definition at line 1031 of file MPll.cpp.

1032{
1033 double Lambdaplus = mySM.getMesons(meson).getLambdaM();
1034 gslpp::complex Lambdamin = exp(-tmpq2 / MM / Lambdaplus) / Lambdaplus * (-gsl_sf_expint_Ei(tmpq2 / MM / Lambdaplus) + gslpp::complex::i() * M_PI);
1035
1036 double u = up;
1037 return ((Tparminus(u, tmpq2)*6. * u * (1. - u)*
1038 (1 + threeGegen0 * (2. * u - 1)
1039 + threeGegen1otwo * ((10. * u - 5.)*(2. * u - 1.) - 1.))) / Lambdamin).real();
1040}

◆ Integrand_ReTparplus()

double MPll::Integrand_ReTparplus ( double  up)
private

The real part of the integral involving \( T^{\parallel}_+ \) at fixed \( q^2 \), according to [Beneke:2001at] .

Parameters
[in]updummy variable to be integrated out
Returns
\( Re T^{\parallel}_+ \Phi^{\parallel}\)

Definition at line 1015 of file MPll.cpp.

1016{
1017 double u = up;
1018 return ((Tparplus(u, tmpq2)*6. * u * (1. - u)*
1019 (1 + threeGegen0 * (2. * u - 1)
1020 + threeGegen1otwo * ((10. * u - 5.)*(2. * u - 1.) - 1.))) / mySM.getMesons(meson).getLambdaM()).real();
1021}

◆ 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 1492 of file MPll.cpp.

1493{
1494 updateParameters();
1495
1496 std::pair<double, double > qbin = std::make_pair(q_min, q_max);
1497
1498 old_handler = gsl_set_error_handler_off();
1499
1500 switch (i) {
1501 case 0:
1502 if (delta0Cached[qbin] == 0) {
1503 FD = convertToGslFunction(bind(&MPll::getDelta1c, &(*this), _1));
1504 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();
1505 cacheDelta0[qbin] = NN*avaDelta;
1506 delta0Cached[qbin] = 1;
1507 }
1508 return cacheDelta0[qbin];
1509 break;
1510 case 2:
1511 if (delta2Cached[qbin] == 0) {
1512 FD = convertToGslFunction(bind(&MPll::getDelta2c, &(*this), _1));
1513 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();
1514 cacheDelta2[qbin] = NN*avaDelta;
1515 delta2Cached[qbin] = 1;
1516 }
1517 return cacheDelta2[qbin];
1518 break;
1519 default:
1520 std::stringstream out;
1521 out << i;
1522 throw std::runtime_error("MPll::integrateDelta: index " + out.str() + " not implemented");
1523 }
1524
1525 gsl_set_error_handler(old_handler);
1526
1527}
double getDelta2c(double q2)
The CP asymmetry .
Definition: MPll.h:1056
double getDelta1c(double q2)
The CP asymmetry .
Definition: MPll.h:1046

◆ 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 1425 of file MPll.cpp.

1426{
1427 updateParameters();
1428
1429 std::pair<double, double > qbin = std::make_pair(q_min, q_max);
1430
1431 old_handler = gsl_set_error_handler_off();
1432
1433 switch (i) {
1434 case 0:
1435 if (sigma0Cached[qbin] == 0) {
1436 FS = convertToGslFunction(bind(&MPll::getSigma1c, &(*this), _1));
1437 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();
1438 cacheSigma0[qbin] = NN*avaSigma;
1439 sigma0Cached[qbin] = 1;
1440 }
1441 return cacheSigma0[qbin];
1442 break;
1443 case 2:
1444 if (sigma2Cached[qbin] == 0) {
1445 FS = convertToGslFunction(bind(&MPll::getSigma2c, &(*this), _1));
1446 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();
1447 cacheSigma2[qbin] = NN*avaSigma;
1448 sigma2Cached[qbin] = 1;
1449 }
1450 return cacheSigma2[qbin];
1451 break;
1452 case 8:
1453 if (sigma8Cached[qbin] == 0) {
1454 FS = convertToGslFunction(bind(&MPll::getSigma6c, &(*this), _1));
1455 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();
1456 cacheSigma8[qbin] = NN*avaSigma;
1457 sigma8Cached[qbin] = 1;
1458 }
1459 return cacheSigma8[qbin];
1460 break;
1461 default:
1462 std::stringstream out;
1463 out << i;
1464 throw std::runtime_error("MPll::integrateSigma: index " + out.str() + " not implemented");
1465 }
1466
1467 gsl_set_error_handler(old_handler);
1468
1469}

◆ integrateSigmaTree()

double MPll::integrateSigmaTree ( double  q_min,
double  q_max 
)

The integral of \( \Sigma_{tree} \) from \(q_{min}\) to \(q_{max}\) (arxiv/2301.06990)

Parameters
[in]q_minminimum q^2 of the integral
[in]q_maxmaximum q^2 of the integral
Returns
\( <\Sigma_{tree}> \)

Definition at line 1529 of file MPll.cpp.

1530{
1531 if (lep != QCD::NEUTRINO_1 or meson != QCD::B_P or !NeutrinoTree_flag) return 0.;
1532
1533 updateParameters();
1534
1535 //phase space limit where tree-level contribution is relevant (0908.1174)
1536 double q_cut = (mtau2 - MP2) * (MM2 - mtau2) / mtau2;
1537 if (q_max >= q_cut) {
1538 if (q_min == 0.) return getintegratedSigmaTree();
1539 q_max = q_cut;
1540 }
1541
1542 double prefactor = mySM.getMesons(meson).getLifetime() / HCUT * GF4 * VusVub_abs2 * fP2 * fM2 / (64. * M_PI2 * MM3 * Gammatau) * mtau2 * mtau;
1543
1544 std::pair<double, double > qbin = std::make_pair(q_min, q_max);
1545
1546 old_handler = gsl_set_error_handler_off();
1547
1548 if (sigmaTreeCached[qbin] == 0) {
1549 FD = convertToGslFunction(bind(&MPll::SigmaTree, &(*this), _1));
1550 if (gsl_integration_cquad(&FD, q_min, q_max, 1.e-2, 1.e-1, w_sigmaTree, &avaSigmaTree, &errSigmaTree, NULL) != 0) return std::numeric_limits<double>::quiet_NaN();
1551 cacheSigmaTree[qbin] = avaSigmaTree;
1552 sigmaTreeCached[qbin] = 1;
1553 }
1554 return prefactor * cacheSigmaTree[qbin];
1555
1556 gsl_set_error_handler(old_handler);
1557}
double SigmaTree(double q2)
Definition: MPll.cpp:1559
double getintegratedSigmaTree()
The integral of from 0 to .
Definition: MPll.cpp:1564
@ B_P
Definition: QCD.h:345

◆ reDC9fit()

double MPll::reDC9fit ( double *  x,
double *  p 
)
private

The fit function for the real part of the QCDF correction \( \Delta C_9^{\lambda} \).

Parameters
[in]xfit variable
[in]pfit parameters vector
Returns
\( f_{Re \Delta C_9^{\lambda}} \)

Definition at line 1140 of file MPll.cpp.

1141{
1142 return p[0] / x[0] + p[1] + p[2] * x[0] + p[3] * x[0] * x[0] + p[4] * x[0] * x[0] * x[0] + p[5] * x[0] * x[0] * x[0] * x[0] + p[6] * x[0] * x[0] * x[0] * x[0] * x[0];
1143}

◆ S_L()

double MPll::S_L ( double  q2)

The helicity form factor \( S_L \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( S_L \)

Definition at line 964 of file MPll.cpp.

965{
966 return S_L_pre * f_0(q2);
967}

◆ SigmaTree()

double MPll::SigmaTree ( double  q2)
private

Definition at line 1559 of file MPll.cpp.

1560{
1561 return MM2 * (mtau2 - MP2) - mtau2 * (mtau2 + q2 - MP2);
1562}

◆ spline_QCDF_func()

void MPll::spline_QCDF_func ( )
private

Definition at line 1248 of file MPll.cpp.

1249{
1250 int dim_DC = GSL_INTERP_DIM_DC;
1251 double min = 0.001;
1252 double interval_DC = (9.9 - min) / ((double) dim_DC);
1253 double q2_spline_DC[dim_DC];
1254 double fq2_Re_deltaC7_QCDF[dim_DC], fq2_Im_deltaC7_QCDF[dim_DC], fq2_Re_deltaC9_QCDF[dim_DC], fq2_Im_deltaC9_QCDF[dim_DC];
1255
1256 for (int i = 0; i < dim_DC; i++) {
1257 q2_spline_DC[i] = min + (double) i*interval_DC;
1258 fq2_Re_deltaC7_QCDF[i] = deltaC7_QCDF(q2_spline_DC[i], false).real();
1259 fq2_Im_deltaC7_QCDF[i] = deltaC7_QCDF(q2_spline_DC[i], false).imag();
1260 fq2_Re_deltaC9_QCDF[i] = deltaC9_QCDF(q2_spline_DC[i], false).real();
1261 fq2_Im_deltaC9_QCDF[i] = deltaC9_QCDF(q2_spline_DC[i], false).imag();
1262
1263 }
1264
1265 gsl_spline_init(spline_Re_deltaC7_QCDF, q2_spline_DC, fq2_Re_deltaC7_QCDF, dim_DC);
1266 gsl_spline_init(spline_Im_deltaC7_QCDF, q2_spline_DC, fq2_Im_deltaC7_QCDF, dim_DC);
1267 gsl_spline_init(spline_Re_deltaC9_QCDF, q2_spline_DC, fq2_Re_deltaC9_QCDF, dim_DC);
1268 gsl_spline_init(spline_Im_deltaC9_QCDF, q2_spline_DC, fq2_Im_deltaC9_QCDF, dim_DC);
1269
1270
1271}

◆ T_L()

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

The helicity form factor \( T_L^0 \).

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

Definition at line 959 of file MPll.cpp.

960{
961 return /*gslpp::complex::i() */ sqrt(lambda(q2) * q2) / (twoMM2_MMpMP) * f_T(q2);
962}

◆ Tparminus()

gslpp::complex MPll::Tparminus ( double  u,
double  q2 
)
private

The \( T^{\parallel}_- \) function from [Khodjamirian:2012rm] .

Parameters
[in]udummy variable to be integrated out
[in]q2\(q^2\) of the decay
Returns
\( T^{\parallel}_- \)

Definition at line 1008 of file MPll.cpp.

1009{
1010 double ubar = 1. - u;
1011 return -spectator_charge * (8. * C_8Lh / (ubar + u * q2 / MM2)
1012 + sixMMoMb * H_c(ubar * MM2 + u * q2, mu_h * mu_h) * C_2Lh_bar);
1013}
double mu_h
Definition: MPll.h:337
double spectator_charge
Definition: MPll.h:342

◆ Tparplus()

gslpp::complex MPll::Tparplus ( double  u,
double  q2 
)
private

The \( T^{\parallel}_+ \) function from [Khodjamirian:2012rm] .

Parameters
[in]udummy variable to be integrated out
[in]q2\(q^2\) of the decay
Returns
\( T^{\parallel}_+ \)

Definition at line 995 of file MPll.cpp.

996{
997 Ee = (MM2 - q2) / twoMM;
998 ubar = 1. - u;
999 arg1 = (fourMc2 - gslpp::complex::i()*1.e-10) / (ubar * MM2 + u * q2) - 1.;
1000 B01 = -2. * sqrt(arg1) * arctan(1. / sqrt(arg1));
1001 arg1 = (fourMc2 - gslpp::complex::i()*1.e-10) / q2 - 1.;
1002 B00 = -2. * sqrt(arg1) * arctan(1. / sqrt(arg1));
1003
1004 gslpp::complex tpar = twoMM / Ee / ubar * I1(u, q2) + (ubar * MM2 + u * q2) / Ee / Ee / ubar / ubar * (B01 - B00);
1005 return -MM / Mb * mySM.getQuarks(QCD::CHARM).getCharge() * tpar*C_2Lh_bar;
1006}
gslpp::complex I1(double u, double q2)
The function from .
Definition: MPll.cpp:973
double getCharge() const
A get method to access the particle charge.
Definition: Particle.h:97
@ CHARM
Definition: QCD.h:326
const Particle & getQuarks(const QCD::quark q) const
A get method to access a quark as an object of the type Particle.
Definition: QCD.h:536

◆ V_L()

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

The helicity form factor \( V_L^0 \).

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

Definition at line 954 of file MPll.cpp.

955{
956 return /*gslpp::complex::i() */ sqrt(lambda(q2)) / (twoMM * sqrt(q2)) * f_plus(q2);
957}

Member Data Documentation

◆ ale

double MPll::ale
private

alpha electromagnetic

Definition at line 331 of file MPll.h.

◆ alpha_s_mub

double MPll::alpha_s_mub
private

@f\aplha_s(\mu_b)$ \( */ double angmomP; /**<angular momentum of meson P; for a resonance, it's replaced by its spin */ int etaP; /**<parity of meson P */ double MW; /**<W boson mass */ gslpp::complex lambda_t; /**<Vckm factor */ gslpp::complex h_0; /**<parameter that contains the contribution from the hadronic hamiltonian */ gslpp::complex h_1; /**<parameter that contains the contribution from the hadronic hamiltonian */ gslpp::complex h_2; /**<parameter that contains the contribution from the hadronic hamiltonian */ gslpp::complex r_1; gslpp::complex r_2; gslpp::complex Delta_C9; gslpp::complex exp_Phase; double Delta_C7_U; double Delta_C9_U; /*LATTICE fit parameters*/ double b_0_fplus;/**<LATTICE fit parameter */ double b_1_fplus;/**<LATTICE fit parameter */ double b_2_fplus;/**<LATTICE fit parameter */ double m_fit2_fplus_lat;/**<LATTICE fit parameter */ double b_0_fT;/**<LATTICE fit parameter */ double b_1_fT;/**<LATTICE fit parameter */ double b_2_fT;/**<LATTICE fit parameter */ double m_fit2_fT_lat;/**<LATTICE fit parameter */ double b_0_f0;/**<LATTICE fit parameter */ double b_1_f0;/**<LATTICE fit parameter */ double b_2_f0;/**<LATTICE fit parameter */ double m_fit2_f0_lat;/**<LATTICE fit parameter */ /*LCSR fit parameters*/ double r_1_fplus;/**<LCSR fit parameter */ double r_2_fplus;/**<LCSR fit parameter */ double m_fit2_fplus;/**<LCSR fit parameter */ double r_1_fT;/**<LCSR fit parameter */ double r_2_fT;/**<LCSR fit parameter */ double m_fit2_fT;/**<LCSR fit parameter */ double r_2_f0;/**<LCSR fit parameter */ double m_fit2_f0;/**<LCSR fit parameter */ gslpp::vector<gslpp::complex> ** allcoeff;/**<vector that contains the Wilson coeffients */ gslpp::vector<gslpp::complex> ** allcoeffh;/**<Vector that contains the Wilson coeffients at scale \)\mu_h \( */ gslpp::vector<gslpp::complex> ** allcoeffprime;/**<vector that contains the primed Wilson coeffients */ gslpp::vector<gslpp::complex> ** allcoeff_noSM;/**<vector that contains the Wilson coeffients without the SM contributions.*/ gslpp::vector<gslpp::complex> ** allcoeff_nu;/**<Vector that contains the Wilson coeffients */ gslpp::vector<gslpp::complex> ** allcoeff_noSM_nu;/**<Vector that contains the Wilson coeffients without the SM contribution.*/ gslpp::complex C_1;/**<Wilson coeffients \)@_fakenlC_1 \(*/ gslpp::complex C_1L_bar;/**<Wilson coeffients \)@_fakenlC_1 \(*/ gslpp::complex C_1Lh_bar;/**<Wilson coeffients \)@_fakenlC_1 \(*/ gslpp::complex C_2;/**<Wilson coeffients \)@_fakenlC_2 \(*/ gslpp::complex C_2L_bar;/**<Leading order Wilson coeffients \)@_fakenlC_2 \(*/ gslpp::complex C_2Lh_bar;/**<Leading order Wilson coeffients \)@_fakenlC_2 \( at scale \)\mu_h \(*/ gslpp::complex C_3;/**<Wilson coeffients \)@_fakenlC_3 \(*/ gslpp::complex C_4;/**<Wilson coeffients \)@_fakenlC_4 \(*/ gslpp::complex C_5;/**<Wilson coeffients \)@_fakenlC_5 \(*/ gslpp::complex C_6;/**<Wilson coeffients \)@_fakenlC_6 \(*/ gslpp::complex C_7;/**<Wilson coeffients \)@_fakenlC_7 \(*/ gslpp::complex C_8;/**<Wilson coeffients \)@_fakenlC_8 \(*/ gslpp::complex C_8L;/**<Leading order Wilson coeffients \)@_fakenlC_8 \(*/ gslpp::complex C_8Lh;/**<Leading order Wilson coeffients \)@_fakenlC_8 \( at scale \)\mu_h \(*/ gslpp::complex C_9;/**<Wilson coeffients \)@_fakenlC_9 \(*/ gslpp::complex C_10;/**<Wilson coeffients \)@_fakenlC_{10} \(*/ gslpp::complex C_S;/**<Wilson coeffients \)@_fakenlC_S \(*/ gslpp::complex C_P;/**<Wilson coeffients \)@_fakenlC_P \(*/ gslpp::complex C_7p;/**<Wilson coeffients \)@_fakenlC_7' \(*/ gslpp::complex C_9p;/**<Wilson coeffients \)@_fakenlC_9' \(*/ gslpp::complex C_10p;/**<Wilson coeffients \)@_fakenlC_{10}' \(*/ gslpp::complex C_Sp;/**<Wilson coeffients \)@_fakenlC_S' \(*/ gslpp::complex C_Pp;/**<Wilson coeffients \)@_fakenlC_P' \(*/ gslpp::complex C_L_nunu_e;/**<Wilson coeffients \)@_fakenlC_L^{\nu_1\nu_1}' \(*/ gslpp::complex C_R_nunu_e;/**<Wilson coeffients \)@_fakenlC_R^{\nu_1\nu_1}' \(*/ gslpp::complex C_L_nunu_mu;/**<Wilson coeffients \)@_fakenlC_L^{\nu_2\nu_2}' \(*/ gslpp::complex C_R_nunu_mu;/**<Wilson coeffients \)@_fakenlC_R^{\nu_2\nu_2}' \(*/ gslpp::complex C_L_nunu_tau;/**<Wilson coeffients \)@_fakenlC_L^{\nu_3\nu_3}' \(*/ gslpp::complex C_R_nunu_tau;/**<Wilson coeffients \)@_fakenlC_R^{\nu_3\nu_3}' \(*/ gsl_interp_accel *acc_Re_deltaC7_QCDF; gsl_interp_accel *acc_Im_deltaC7_QCDF; gsl_interp_accel *acc_Re_deltaC9_QCDF; gsl_interp_accel *acc_Im_deltaC9_QCDF; gsl_spline *spline_Re_deltaC7_QCDF; gsl_spline *spline_Im_deltaC7_QCDF; gsl_spline *spline_Re_deltaC9_QCDF; gsl_spline *spline_Im_deltaC9_QCDF; std::vector<double> ReDeltaC9;/**<Vector that samples the QCDF \)\Delta C_9 \( */ std::vector<double> ImDeltaC9;/**<Vector that samples the QCDF \)\Delta C_9 \( */ std::vector<double> myq2;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ TFitResultPtr refres;/**<Vector that contains the fitted QCDF \)\Delta C_9 \( */ TFitResultPtr imfres;/**<Vector that contains the fitted QCDF \)\Delta C_9 \( */ TGraph gr1;/**<Tgraph to be used for fitting the QCDF \)\Delta C_9 \( */ TGraph gr2;/**<Tgraph to be used for fitting the QCDF \)\Delta C_9 \( */ TF1 reffit;/**<TF1 to be used for fitting the QCDF \)\Delta C_9 \( */ TF1 imffit;/**<TF1 to be used for fitting the QCDF \)\Delta C_9 \( */ double tmpq2;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex H_0_pre;/**< Cache variable */ gslpp::complex H_0_WC;/**< Cache variable */ gslpp::complex H_c_WC;/**< Cache variable */ gslpp::complex H_b_WC;/**< Cache variable */ gslpp::complex ihalfMPI;/**< Cache variable */ double fournineth;/**< Cache variable */ double half;/**< Cache variable */ double twothird;/**< Cache variable */ double Mc2;/**< Cache variable */ double Mb2;/**< Cache variable */ double logMc;/**< Cache variable */ double logMb;/**< Cache variable */ double mu_b2;/**< Cache variable */ double fourMc2;/**< Cache variable */ double fourMb2;/**< Cache variable */ double Mlep2;/**< Cache variable */ double NN;/**< Cache variable */ double MM2;/**< Cache variable */ double MM4;/**< Cache variable */ double MP2;/**< Cache variable */ double MP4;/**< Cache variable */ double MM2mMP2;/**< Cache variable */ double twoMP2;/**< Cache variable */ double twoMM;/**< Cache variable */ double twoMM2;/**< Cache variable */ double twoMM2_MMpMP;/**< Cache variable */ double twoMM_MbpMs;/**< Cache variable */ double S_L_pre;/**< Cache variable */ double fourMM2;/**< Cache variable */ double twoMboMM;/**< Cache variable */ double sixteenM_PI2;/**< Cache variable */ double ninetysixM_PI3MM3;/**< Cache variable */ double MboMW;/**< Cache variable */ double MboMM;/**< Cache variable */ double MsoMb;/**< Cache variable */ double twoMlepMb;/**< Cache variable */ double threeGegen0;/**< Cache variable */ double threeGegen1otwo;/**< Cache variable */ double M_PI2osix;/**< Cache variable */ double twoMc2;/**< Cache variable */ double sixMMoMb;/**< Cache variable */ double CF;/**< Cache variable */ double deltaT_0;/**< Cache variable */ double deltaT_1par;/**< Cache variable */ gslpp::complex ubar;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex arg1;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex B01;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex B00;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex xp;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex xm;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex yp;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex ym;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex L1xp;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex L1xm;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex L1yp;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex L1ym;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F87_0;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F87_1;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F87_2;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F87_3;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F29_0;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F29_L1;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F29_1;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F29_2;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F29_3;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F29_L1_1;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F29_L1_2;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F29_L1_3;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F27_0;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F27_1;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F27_2;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F27_3;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F27_L1_1;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F27_L1_2;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ gslpp::complex F27_L1_3;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ double F89_0;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ double F89_1;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ double F89_2;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ double F89_3;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ double Ee;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ //additional variables for B to K nu nu const double M_PI2 = M_PI * M_PI; double GF4; double MM3; double fM2; double fP2; double mtau; double mtau2; double Gammatau; double VusVub_abs2; unsigned int fplus_lat_updated;/**< Cache variable */ gslpp::vector<double> fplus_lat_cache;/**< Cache variable */ unsigned int fT_lat_updated;/**< Cache variable */ gslpp::vector<double> fT_lat_cache;/**< Cache variable */ unsigned int f0_lat_updated;/**< Cache variable */ gslpp::vector<double> f0_lat_cache;/**< Cache variable */ unsigned int fplus_updated;/**< Cache variable */ gslpp::vector<double> fplus_cache;/**< Cache variable */ unsigned int fT_updated;/**< Cache variable */ gslpp::vector<double> fT_cache;/**< Cache variable */ unsigned int f0_updated;/**< Cache variable */ double f0_cache;/**< Cache variable */ unsigned int k2_updated;/**< Cache variable */ gslpp::vector<double> k2_cache;/**< Cache variable */ unsigned int beta_updated;/**< Cache variable */ double beta_cache;/**< Cache variable */ unsigned int lambda_updated;/**< Cache variable */ unsigned int F_updated;/**< Cache variable */ unsigned int VL_updated;/**< Cache variable */ unsigned int TL_updated;/**< Cache variable */ unsigned int SL_updated;/**< Cache variable */ gslpp::vector<double> SL_cache;/**< Cache variable */ unsigned int N_updated;/**< Cache variable */ gslpp::vector<double> N_cache;/**< Cache variable */ gslpp::complex Nc_cache;/**< Cache variable */ unsigned int C_1_updated;/**< Cache variable */ gslpp::complex C_1_cache;/**< Cache variable */ unsigned int C_2_updated;/**< Cache variable */ gslpp::complex C_2_cache;/**< Cache variable */ unsigned int C_3_updated;/**< Cache variable */ gslpp::complex C_3_cache;/**< Cache variable */ unsigned int C_4_updated;/**< Cache variable */ gslpp::complex C_4_cache;/**< Cache variable */ unsigned int C_5_updated;/**< Cache variable */ gslpp::complex C_5_cache;/**< Cache variable */ unsigned int C_6_updated;/**< Cache variable */ gslpp::complex C_6_cache;/**< Cache variable */ unsigned int C_7_updated;/**< Cache variable */ gslpp::complex C_7_cache;/**< Cache variable */ unsigned int C_9_updated;/**< Cache variable */ gslpp::complex C_9_cache;/**< Cache variable */ unsigned int C_10_updated;/**< Cache variable */ gslpp::complex C_10_cache;/**< Cache variable */ unsigned int C_7p_updated;/**< Cache variable */ gslpp::complex C_7p_cache;/**< Cache variable */ unsigned int C_9p_updated;/**< Cache variable */ gslpp::complex C_9p_cache;/**< Cache variable */ unsigned int C_10p_updated;/**< Cache variable */ gslpp::complex C_10p_cache;/**< Cache variable */ unsigned int C_S_updated;/**< Cache variable */ gslpp::complex C_S_cache;/**< Cache variable */ unsigned int C_P_updated;/**< Cache variable */ gslpp::complex C_P_cache;/**< Cache variable */ unsigned int C_Sp_updated;/**< Cache variable */ gslpp::complex C_Sp_cache;/**< Cache variable */ unsigned int C_Pp_updated;/**< Cache variable */ gslpp::complex C_Pp_cache;/**< Cache variable */ unsigned int C_2Lh_updated;/**< Cache variable */ gslpp::complex C_2Lh_cache;/**< Cache variable */ unsigned int C_8Lh_updated;/**< Cache variable */ gslpp::complex C_8Lh_cache;/**< Cache variable */ unsigned int C_L_nunu_e_updated;/**< Cache variable */ unsigned int C_L_nunu_mu_updated;/**< Cache variable */ unsigned int C_L_nunu_tau_updated;/**< Cache variable */ gslpp::complex C_L_nunu_e_cache;/**< Cache variable */ gslpp::complex C_L_nunu_mu_cache;/**< Cache variable */ gslpp::complex C_L_nunu_tau_cache;/**< Cache variable */ unsigned int C_R_nunu_e_updated;/**< Cache variable */ unsigned int C_R_nunu_mu_updated;/**< Cache variable */ unsigned int C_R_nunu_tau_updated;/**< Cache variable */ gslpp::complex C_R_nunu_e_cache;/**< Cache variable */ gslpp::complex C_R_nunu_mu_cache;/**< Cache variable */ gslpp::complex C_R_nunu_tau_cache;/**< Cache variable */ unsigned int Yupdated;/**< Cache variable */ gslpp::vector<double> Ycache;/**< Cache variable */ unsigned int H_V0updated;/**< Cache variable */ gslpp::vector<double> H_V0cache;/**< Cache variable */ gslpp::complex H_V0Ccache[3];/**< Cache variable */ gslpp::complex H_V0Ccache_dispersion[4];/**< Cache variable */ unsigned int H_A0updated;/**< Cache variable */ unsigned int H_Supdated;/**< Cache variable */ gslpp::vector<double> H_Scache;/**< Cache variable */ unsigned int H_P_updated;/**< Cache variable */ gslpp::vector<double> H_P_cache;/**< Cache variable */ unsigned int I0_updated;/**< Cache variable */ unsigned int I2_updated;/**< Cache variable */ unsigned int I8_updated;/**< Cache variable */ unsigned int Itree_updated;/**< Cache variable */ gslpp::vector<double> Itree_cache;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > I1Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma0Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma2Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma8Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta0Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta2Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigmaTreeCached;/**< Cache variable */ double avaSigma;/**< Gsl integral variable */ double avaDelta;/**< Gsl integral variable */ double avaSigmaTree;/**< Gsl integral variable */ double avaDTPPR;/**< Gsl integral variable */ double errSigma;/**< Gsl integral variable */ double errDelta;/**< Gsl integral variable */ double errSigmaTree;/**< Gsl integral variable */ double errDTPPR;/**< Gsl integral variable */ gsl_function FS;/**< Gsl integral variable */ gsl_function FD;/**< Gsl integral variable */ gsl_function DTPPR;/**< Gsl integral variable */ gsl_integration_cquad_workspace * w_sigma;/**< Gsl integral variable */ gsl_integration_cquad_workspace * w_delta;/**< Gsl integral variable */ gsl_integration_cquad_workspace * w_sigmaTree;/**< Gsl integral variable */ gsl_integration_cquad_workspace * w_DTPPR;/**< Gsl integral variable */ gsl_error_handler_t * old_handler; /**< GSL error handler store */ std::map<std::pair<double, double>, gslpp::complex > cacheI1;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma0;/**< Gsl integral variable */ std::map<std::pair<double, double>, double > cacheSigma2;/**< Gsl integral variable */ std::map<std::pair<double, double>, double > cacheSigma8;/**< Gsl integral variable */ std::map<std::pair<double, double>, double > cacheDelta0;/**< Gsl integral variable */ std::map<std::pair<double, double>, double > cacheDelta2;/**< Gsl integral variable */ std::map<std::pair<double, double>, double > cacheSigmaTree;/**< Gsl integral variable */ std::map<double, unsigned int> deltaTparpCached;/**< Cache variable */ std::map<double, unsigned int> deltaTparmCached;/**< Cache variable */ std::map<double, gslpp::complex> cacheDeltaTparp;/**< Cache variable */ std::map<double, gslpp::complex> cacheDeltaTparm;/**< Cache variable */ unsigned int deltaTparpupdated;/**< Cache variable */ unsigned int deltaTparmupdated;/**< Cache variable */ unsigned int T_updated;/**< Cache variable */ gslpp::vector<double> T_cache;/**< Cache variable */ /** @brief The update parameter method for MPll. */ void updateParameters(); /** @brief The caching method for MPll. */ void checkCache(); /** @brief The first fit function from arXiv:hep-ph/0412079v1,\) f_1^{LATTICE} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] r_1 fit parameter @param[in] r_2 fit parameter @param[in] m_fit2 fit parameter @return \) f_1^{LATTICE} \( */ double LCSR_fit1(double q2, double r_1, double r_2, double m_fit2); /** @brief The second fit function from arXiv:hep-ph/0412079v1, \) f_2^{LATTICE} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] r_2 fit parameter @param[in] m_fit2 fit parameter @return \) f_2^{LATTICE} \( */ double LCSR_fit2(double q2, double r_2, double m_fit2); /** @brief The fit function from @cite Straub:2015ica, \) FF^{\rm fit} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] b_0 fit parameter @param[in] b_1 fit parameter @param[in] b_2 fit parameter @param[in] m_fit2 square of the nearest resonance mass @return \) FF^{\rm fit} \( */ double LCSR_fit3(double q2, double b_0, double b_1, double b_2, double m_fit2); /** @brief The z function from arXiv:hep-ph/0412079v1, \) z \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) z \( */ double zeta(double q2); /** @brief The z function used for the DM Form Factors, \) z \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) z \( */ double zeta_DM(double q2); /** @brief The first fit function from arXiv:hep-ph/0412079v1,\) f_+^{DM} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] m_fit2 fit parameter @return \) f_+^{DM} \( */ double phiplus_DM(double q2, double m_fit2); /** @brief The first fit function from arXiv:hep-ph/0412079v1,\) f_+^{DM} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] m_fit2 fit parameter @return \) f_+^{DM} \( */ double phi0_DM(double q2); /** @brief The first fit function from arXiv:hep-ph/0412079v1,\) f_+^{DM} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] m_fit2 fit parameter @return \) f_+^{DM} \( */ double phiT_DM(double q2, double m_fit2); /** @brief The first fit function from arXiv:hep-ph/0412079v1,\) f_+^{DM} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] b_0 fit parameter @param[in] b_1 fit parameter @param[in] b_2 fit parameter @param[in] m_fit2 fit parameter @return \) f_+^{DM} \( */ double fplus_DM(double q2, double b_0, double b_1, double b_2, double m_fit2); /** @brief The first fit function from arXiv:hep-ph/0412079v1,\) f_0^{DM} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] b_0 fit parameter @param[in] b_1 fit parameter @param[in] b_2 fit parameter @param[in] m_fit2 fit parameter @return \) f_0^{DM} \( */ double f0_DM(double q2, double b_0, double b_1, double b_2); /** @brief The first fit function from arXiv:hep-ph/0412079v1,\) f_T^{DM} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] b_0 fit parameter @param[in] b_1 fit parameter @param[in] b_2 fit parameter @param[in] m_fit2 fit parameter @return \) f_T^{DM} \( */ double fT_DM(double q2, double b_0, double b_1, double b_2, double m_fit2); /** @brief The first fit function from arXiv:1509.06235v2,\) f_1^{LATTICE} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] b_0 fit parameter @param[in] b_1 fit parameter @param[in] b_2 fit parameter @param[in] m_fit2 fit parameter @return \) f_1^{LATTICE} \( */ double LATTICE_fit1(double q2, double b_0, double b_1, double b_2, double m_fit2); /** @brief The second fit function from arXiv:1509.06235v2, \) f_2^{LATTICE} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] b_0 fit parameter @param[in] b_1 fit parameter @param[in] b_2 fit parameter @param[in] m_fit2 fit parameter @return \) f_2^{LATTICE} \( */ double LATTICE_fit2(double q2, double b_0, double b_1, double b_2, double m_fit2); /** @brief The form factor \) f_+ \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) f_+ \( */ double f_plus(double q2); /** @brief The form factor \) f_T \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) f_T \( */ double f_T(double q2); /** @brief The form factor \) f_0 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) f_0 \( */ double f_0(double q2); /** @brief The \) h(q^2,0) \( function involved into \) C_9^{eff} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) h(q^2,0) \( */ gslpp::complex H_0(double q2); /** @brief The \) h(q^2,m_c) \( function involved into \) C_9^{eff} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] mu2 mass scale @return \) h(q^2,m_c) \( */ gslpp::complex H_c(double q2, double mu2); /** @brief The \) h(q^2,m_b) \( function involved into \) C_9^{eff} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] mu2 mass scale @return \) h(q^2,m_b) \( */ gslpp::complex H_b(double q2, double mu2); /** @brief The \) Y(q^2) \( function involved into \) C_9^{eff} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Y(q^2) \( */ gslpp::complex Y(double q2); /** @brief The square of the 3-momentum of the recoiling meson in the M rest frame, \) k^2 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) k^2 \( */ double k2 (double q2); /** @brief The factor \) \beta \( used in the angular coefficients \)I_i \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \beta \( */ double beta (double q2); /** @brief The factor \) \beta^2 \( used in the angular coefficients \)I_i \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \beta^2 \( */ double beta2 (double q2); /** @brief The factor \) \lambda \( used in the angular coefficients \)I_i \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \lambda \( */ double lambda(double q2); /** @brief The factor \) F \( used in the angular coefficients \)I_i \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) F \( */ double F(double q2); /** @brief The angular coefficient \) I_{1c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) I_{1c} \( */ double I_1c(double q2); /** @brief The angular coefficient \) I_{2c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) I_{2c} \( */ double I_2c(double q2); /** @brief The angular coefficient \) I_{6c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) I_{6c} \( */ double I_6c(double q2); /** @brief The CP asymmetry \) \Delta_{i} \( . @param[in] i index of the angular coefficient \) I_{i} \( @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{i} \( */ double Delta(int i, double q2); /** @brief The value of \) \Sigma_{tree}: contains the full q2-dependence but neglects a "prefactor"

Parameters
[in]q2\(q^2\) of the decay
Returns
\( <\Sigma{tree}> \)

Definition at line 344 of file MPll.h.

◆ dispersion

bool MPll::dispersion = false
private

Definition at line 322 of file MPll.h.

◆ FixedWCbtos

bool MPll::FixedWCbtos = false
private

Definition at line 323 of file MPll.h.

◆ GF

double MPll::GF
private

Fermi constant

Definition at line 330 of file MPll.h.

◆ lep

QCD::lepton MPll::lep
private

Final leptons type

Definition at line 316 of file MPll.h.

◆ Mb

double MPll::Mb
private

b quark mass

Definition at line 335 of file MPll.h.

◆ mb_pole

double MPll::mb_pole
private

b quark pole mass

Definition at line 339 of file MPll.h.

◆ Mc

double MPll::Mc
private

c quark mass

Definition at line 338 of file MPll.h.

◆ mc_pole

double MPll::mc_pole
private

c quark pole mass

Definition at line 340 of file MPll.h.

◆ meson

QCD::meson MPll::meson
private

Initial meson type

Definition at line 317 of file MPll.h.

◆ mJ2

double MPll::mJ2
private

Definition at line 328 of file MPll.h.

◆ Mlep

double MPll::Mlep
private

muon mass

Definition at line 332 of file MPll.h.

◆ MM

double MPll::MM
private

initial meson mass

Definition at line 333 of file MPll.h.

◆ MP

double MPll::MP
private

final pseudoscalar meson mass

Definition at line 334 of file MPll.h.

◆ MPll_DM_flag

bool MPll::MPll_DM_flag
private

A flag for switching to DM FF parameterization

Definition at line 326 of file MPll.h.

◆ MPll_GRvDV_flag

bool MPll::MPll_GRvDV_flag
private

A flag for switching to GRvDV parameterization

Definition at line 325 of file MPll.h.

◆ MPll_Lattice_flag

bool MPll::MPll_Lattice_flag
private

A flag for switching to LATTICE FF parameterization

Definition at line 324 of file MPll.h.

◆ mpllParameters

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

The string of mandatory MPll parameters

Definition at line 319 of file MPll.h.

◆ Ms

double MPll::Ms
private

s quark mass

Definition at line 341 of file MPll.h.

◆ mu_b

double MPll::mu_b
private

b mass scale

Definition at line 336 of file MPll.h.

◆ mu_h

double MPll::mu_h
private

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

Definition at line 337 of file MPll.h.

◆ myF_1

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

Definition at line 320 of file MPll.h.

◆ myF_2

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

Definition at line 321 of file MPll.h.

◆ mySM

const StandardModel& MPll::mySM
private

Model type

Definition at line 315 of file MPll.h.

◆ NeutrinoTree_flag

bool MPll::NeutrinoTree_flag
private

Definition at line 327 of file MPll.h.

◆ pseudoscalar

QCD::meson MPll::pseudoscalar
private

Final pseudoscalar meson type

Definition at line 318 of file MPll.h.

◆ spectator_charge

double MPll::spectator_charge
private

spectator quark charge

Definition at line 342 of file MPll.h.

◆ width

double MPll::width
private

initial meson width

Definition at line 343 of file MPll.h.


The documentation for this class was generated from the following files: