A class for the \(M \to P l^+ l^-\) decay.
- Author
- HEPfit Collaboration
- Copyright
- GNU General Public License
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.
|
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...
|
|
@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; /*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;/**<Wilson coeffients \)@_fakenlC_L^{\nu\nu}' \(*/ gslpp::complex C_R_nunu;/**<Wilson coeffients \)@_fakenlC_R^{\nu\nu}' \(*/ 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_updated;/**< Cache variable */ gslpp::complex C_L_nunu_cache;/**< Cache variable */ unsigned int C_R_nunu_updated;/**< Cache variable */ gslpp::complex C_R_nunu_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 > 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 > 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 336 of file MPll.h.