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

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

#include <MVll.h>

Detailed Description

A class for the \(M \to V 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 V l^+ l^-\) decays, where \(M\) is a generic meson and \(V\) is a vector meson.

MVll parameters

The mandatory parameters of MVll are summarized below:

Label LaTeX symbol Description
a_0V, a_1V, a_2V, MRV \(a_0^V, a_1^V, a_2^V, \Delta m^V\) The fit parameters for the form factor \(V\) of the \(B\to K^*\).
a_0A0, a_1A0, a_2A0, MRA0 \(a_0^{A_0}, a_1^{A_0}, a_2^{A_0}, \Delta m^{A_0}\) The fit parameters for the form factor \(A_0\) of the \(B\to K^*\).
a_0A1, a_1A1, a_2A1, MRA1 \(a_0^{A_1}, a_1^{A_1}, a_2^{A_1}, \Delta m^{A_1}\) The fit parameters for the form factor \(A_1\) of the \(B\to K^*\).
a_0A12, a_1A12, a_2A12, MRA12 \(a_0^{A_{12}}, a_1^{A_{12}}, a_2^{A_{12}}, \Delta m^{A_{12}}\) The fit parameters for the form factor \(A_{12}\) of the \(B\to K^*\).
a_0T1, a_1T1, a_2T1, MRA0 \(a_0^{T_1}, a_1^{T_1}, a_2^{T_1}, \Delta m^{T_1}\) The fit parameters for the form factor \(T_1\) of the \(B\to K^*\).
a_0T2, a_1T2, a_2T2, MRA1 \(a_0^{T_2}, a_1^{T_2}, a_2^{T_2}, \Delta m^{T_2}\) The fit parameters for the form factor \(T_2\) of the \(B\to K^*\).
a_0T23, a_1T23, a_2T23, MRA1 \(a_0^{T_{23}}, a_1^{T_{23}}, a_2^{T_{23}}, \Delta m^{T_{23}}\) The fit parameters for the form factor \(T_{23}\) of the \(B\to K^*\).
a_0Vphi, a_1Vphi, a_2Vphi, MRVphi \(a_0^V, a_1^V, a_2^V, \Delta m^V\) The fit parameters for the form factor \(V\) of the \(B\to\phi\).
a_0A0phi, a_1A0phi, a_2A0phi, MRA0phi \(a_0^{A_0}, a_1^{A_0}, a_2^{A_0}, \Delta m^{A_0}\) The fit parameters for the form factor \(A_0\) of the \(B\to\phi\).
a_0A1phi, a_1A1phi, a_2A1phi, MRA1phi \(a_0^{A_1}, a_1^{A_1}, a_2^{A_1}, \Delta m^{A_1}\) The fit parameters for the form factor \(A_1\) of the \(B\to\phi\).
a_0A1phi, a_1A1phi, a_2A1phi, MRA1phi \(a_0^{A_{12}}, a_1^{A_{12}}, a_2^{A_{12}}, \Delta m^{A_{12}}\) The fit parameters for the form factor \(A_{12}\) of the \(B\to\phi\).
a_0T1phi, a_1T1phi, a_2T1phi, MRA0phi \(a_0^{T_1}, a_1^{T_1}, a_2^{T_1}, \Delta m^{T_1}\) The fit parameters for the form factor \(T_1\) of the \(B\to\phi\).
a_0T2phi, a_1T2phi, a_2T2phi, MRA1phi \(a_0^{T_2}, a_1^{T_2}, a_2^{T_2}, \Delta m^{T_2}\) The fit parameters for the form factor \(T_2\) of the \(B\to\phi\).
a_0T23phi, a_1T23phi, a_2T23phi, MRA1phi \(a_0^{T_{23}}, a_1^{T_{23}}, a_2^{T_{23}}, \Delta m^{T_{23}}\) The fit parameters for the form factor \(T_{23}\) of the \(B\to\phi\).
absh_0, absh_0_1, absh_0_2 \(\mathrm{Abs}h_0^{(0)}, \mathrm{Abs}h_0^{(1)}, \mathrm{Abs}h_0^{(2)}\) The constant, linear and quadratic terms of the absolute value of the hadronic parameter \(h_0\) of the \(B\to K^*\).
argh_0, argh_0_1, argh_0_2 \(\mathrm{Arg}h_0^{(0)}, \mathrm{Arg}h_0^{(1)}, \mathrm{Arg}h_0^{(2)}\) The constant, linear and quadratic terms of the argument of the hadronic parameter \(h_0\) of the \(B\to K^*\).
absh_p, absh_p_1, absh_p_2 \(\mathrm{Abs}h_+^{(0)}, \mathrm{Abs}h_+^{(1)}, \mathrm{Abs}h_+^{(2)}\) The constant, linear and quadratic terms of the absolute value of the hadronic parameter \(h_+\) of the \(B\to K^*\).
argh_p, argh_p_1, argh_p_2 \(\mathrm{Arg}h_+^{(0)}, \mathrm{Arg}h_+^{(1)}, \mathrm{Arg}h_+^{(2)}\) The constant, linear and quadratic terms of the argument of the hadronic parameter \(h_+\) of the \(B\to K^*\).
absh_m, absh_m_1, absh_m_2 \(\mathrm{Abs}h_-^{(0)}, \mathrm{Abs}h_-^{(1)}, \mathrm{Abs}h_-^{(2)}\) The constant, linear and quadratic terms of the absolute value of the hadronic parameter \(h_-\) of the \(B\to K^*\).
argh_m, argh_m_1, argh_m_2 \(\mathrm{Arg}h_-^{(0)}, \mathrm{Arg}h_-^{(1)}, \mathrm{Arg}h_-^{(2)}\) The constant, linear and quadratic terms of the argument of the hadronic parameter \(h_-\) 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 \(V 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 [143], we implemented the amplitude in the helicity basis; hence we made use of the helicity form factors \( \tilde{V}_\lambda(q^2), \tilde{T}_\lambda(q^2)\) and \(\tilde{S}(q^2) \) (where \(\lambda=+,-,0\) represents the helicity), which are related to the ones in the transverse basis through the following relations :

\[ \tilde{V}_0(q^2) = \frac{4m_V}{\sqrt{q^2}}A_{12}(q^2)\,,\\ \tilde{V}_{\pm}\left( q^{2}\right) = \frac{1}{2} \bigg[ \Big( 1 + \frac{m_V}{m_M} \Big) A_1\left( q^{2}\right) \mp \frac{\lambda^{1/2}(q^2)}{m_M(m_M + m_V)} V\left( q^{2}\right) \bigg]\,, \\ \tilde{T}_0(q^2)=\frac{2\sqrt{q^2}m_V}{m_M(m_M + m_V)}T_{23}(q^2)\,,\\ \tilde{T}_{\pm}\left( q^{2}\right) = \frac{m_M^2 - m_V^2}{2m_M^2}T_2\left( q^{2}\right) \mp \frac{\lambda^{1/2}(q^2)}{2m_M^2}T_1\left( q^{2}\right)\,,\\ \tilde{S}\left( q^{2}\right) = -\frac{\lambda^{1/2}(q^2)}{2m_M(m_b+m_s)}A_0\left( q^{2}\right)\,, \]

where \(\lambda(q^2) = 4m_M^2|\vec{k}|^2\), with \(\vec{k}\) as the 3-momentum of the meson \(V\) 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_\lambda(q^2) = \frac{\epsilon^*_\mu(\lambda)}{m_M^2} \int d^4x e^{iqx} \langle \bar V \vert T\{j^{\mu}_\mathrm{em} (x) \mathcal{H}_\mathrm{eff}^\mathrm{had} (0)\} \vert \bar M \rangle = h_\lambda^{(0)} + \frac{q^2}{1\,\mathrm{GeV}^2} h_\lambda^{(1)} + \frac{q^4}{1\, \mathrm{GeV}^4} h_\lambda^{(2)} \,, \]

while the effect due to exchange of hard gluons can be parametrized following the prescription of [42] as a shift to the Wilson coefficient \(C_9\) : one first have to define the corrections

\[ \Delta \mathcal{T}_a = \frac{\alpha_sC_F}{4\pi} C_a + \frac{\alpha_sC_F}{4}\frac{\pi}{N_c}\frac{f_Mf_{V,a}}{m_V F_a(q^2)}\Xi_a \sum_{\pm}\int \frac{d\omega}{\omega}\Phi_{V,\pm}(\omega)\int_0^1du\Phi_{M,a}(u)T_{a,\pm}(u,\omega)\,, \]

where \(a=\perp,\parallel\), \(F_\perp(q^2) = T_1(q^2) \), \(F_\parallel(q^2) = T_1(q^2) - T_3(q^2)\), \(\Xi_\perp(q^2) = 1 \), \(\Xi_\parallel(q^2) = \frac{2m_Vm_M}{m_M^2-q^2}\), and \(\Phi_X\) are leading twist light-cone distributions; the term proportional to \(C_a\) is the one describing the corrections where the spectator quark is connected to the hard process only through soft interactions, while the one proportional to \(T_{a,\pm}\) is the one describing the corrections where the spectactor quark is involved in the hard process. Therefore, it is possible to define the correction to the Wilson coefficient in the following way:

\[ \Delta C_{9,\pm} = \frac{1}{q^2}\frac{m_b}{m_M} \left((m_M^2-m_V^2) \frac{m_M^2 - q^2}{m_M^2} \mp \sqrt{\lambda(q^2)}\right) \Delta T_{\perp}(q^2)\,,\\ \Delta C_{9,0} = \frac{1}{ 2 m_V m_M \sqrt{q^2} } \left(\left[(m_M^2-m_V^2) ( m_M^2-m_V^2 - q^2) - \lambda(q^2)\right] (m_M^2 - q^2) \frac{m_b}{m_M^2q^2} \Delta T_{\perp}(q^2) - \lambda(q^2) \frac{m_b}{m_M^2-m_V^2}\left(\Delta T_{\parallel}(q^2) + \Delta T_\perp(q^2)\right)\right)\,. \]

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

\[ H_V(\lambda) = -i\, N \Big\{C_{9} \tilde{V}_{L\lambda} +C_{9}' \tilde{V}_{R\lambda} + \frac{m_M^2}{q^2} \Big[\frac{2\, m_b}{m_M} (C_{7} \tilde{T}_{L\lambda} + C_{7}' \tilde{T}_{R\lambda}) - 16 \pi^2 h_\lambda \Big] \Big\} \,, \\ H_A(\lambda) = -i\, N (C_{10} \tilde{V}_{L\lambda} + C_{10}'\tilde{V}_{R\lambda}) \,, \\ 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\pm}(q^2) = -\tilde{V}_{R\mp}(q^2)=\tilde{V}_\pm(q^2)\,,\\ \tilde{T}_{L\pm}(q^2) = -\tilde{T}_{R\mp}(q^2)=\tilde{T}_\pm(q^2)\,,\\ \tilde{S}_L(q^2) = -\tilde{S}_R(q^2)=\tilde{S}(q^2)\,. \]

The hadronic non-factorizable contribution has been parameterized in the following way:

\begin{eqnarray} h_-(q^2) &=& -\frac{m_b}{8\pi^2 m_B} \tilde T_{L -}(q^2) h_-^{(0)} -\frac{\tilde V_{L -}(q^2)}{16\pi^2 m_B^2} h_-^{(1)} q^2 + h_-^{(2)} q^4+{\cal O}(q^6)\,,\\ h_+(q^2) &=& -\frac{m_b}{8\pi^2 m_B} \tilde T_{L +}(q^2) h_-^{(0)} -\frac{\tilde V_{L +}(q^2)}{16\pi^2 m_B^2} h_-^{(1)} q^2 + h_+^{(0)} + h_+^{(1)}q^2 + h_+^{(2)} q^4+{\cal O}(q^6)\,,\\ h_0(q^2) &=& -\frac{m_b}{8\pi^2 m_B} \tilde T_{L 0}(q^2) h_-^{(0)} -\frac{\tilde V_{L 0}(q^2)}{16\pi^2 m_B^2} h_-^{(1)} q^2 + h_0^{(0)}\sqrt{q^2} + h_0^{(1)}(q^2)^\frac{3}{2} +{\cal O}((q^2)^\frac{5}{2})\,. \end{eqnarray}

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)d(\cos\theta_k)d\phi} = \frac{9}{32\,\pi} \Big( I^s_1\sin^2\theta_k+I^c_1\cos^2\theta_k +(I^s_2\sin^2\theta_k+I^c_2\cos^2\theta_k)\cos2\theta_l \\ + I_3\sin^2\theta_k\sin^2\theta_l\cos2\phi +I_4\sin2\theta_k\sin2\theta_l\cos\phi + I_5\sin2\theta_k\sin\theta_l\cos\phi \\ +(I_6^s\sin^2\theta_k + I_6^c \cos^2\theta_K) \cos\theta_l + I_7\sin2\theta_k\sin\theta_l\sin\phi+I_8\sin2\theta_k\sin2\theta_l\sin\phi +I_9\sin^2\theta_k\sin^2\theta_l\sin2\phi \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_1^s = F \left\{\frac{\beta^2\!+\!2}{8}\left(|H_V^+|^2+|H_V^-|^2+(V\rightarrow A)\right) +\frac{m_\ell^2}{q^2}\left(|H_V^+|^2+|H_V^-|^2-(V\rightarrow A)\right)\right\}\,,\\ I_2^c = -F\, \frac{\beta^2}{2}\left(|H_V^0|^2+|H_A^0|^2\right)\,,\\ I_2^s = F\, \frac{\beta^2}{8}\left(|H_V^+|^2+|H_V^-|^2\right)+(V\rightarrow A)\,,\\ I_3 = -\frac{F}{2}{\rm Re} \left[H_V^+(H_V^-)^*\right]+(V\rightarrow A)\,,\\ I_4 = F\, \frac{\beta^2}{4}{\rm Re}\left[(H_V^-+H_V^+)\left(H_V^0\right)^*\right]+(V\rightarrow A)\,,\\ I_5 = F\left\{ \frac{\beta}{2}{\rm Re}\left[(H_V^--H_V^+)\left(H_A^0\right)^*\right] +(V\leftrightarrow A) - \frac{\beta\,m_\ell}{\sqrt{q^2}} {\rm Re} \left[H_S^* (H_V^+ + H_V^-)\right]\right\}\,,\\ I_6^s = F \beta\,{\rm Re}\left[H_V^-(H_A^-)^*-H_V^+(H_A^+)^*\right]\,,\\ I_6^c = 2 F \frac{\beta\, m_\ell}{\sqrt{q^2}} {\rm Re} \left[ H_S^* H_V^0 \right]\,,\\ I_7 = F \left\{ \frac{\beta}{2}\,{\rm Im}\left[\left(H_A^++H_A^-\right) (H_V^0)^* \, +(V\leftrightarrow A) \right] - \frac{\beta\, m_\ell}{\sqrt{q^2} }\, {\rm Im} \left[ H_S^*(H_V^{-} - H_V^{+}) \right] \right\}\,,\\ I_8 = F\, \frac{\beta^2}{4}{\rm Im}\left[(H_V^--H_V^+)(H_V^0)^*\right]+(V\rightarrow A)\,,\\ I_9 = F\, \frac{\beta^2}{2}{\rm Im}\left[H_V^+(H_V^-)^*\right]+(V\rightarrow A)\,, \]

where

\[ F=\frac{ \lambda^{1/2}\beta\, q^2}{3 \times 2^{5} \,\pi^3\, m_M^3} BF(V \to {\rm final \, state})\,, \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 V(), A_0(), A_1(), A_2(), T_1() and T_2() using the fit function FF_fit() from [45] . The form factor are consequentely translated in the helicity basis through the functions V_0t(), V_p(), V_m(), T_0t(), T_p(), T_m() and S_L() . The basic elements required to compute the hard gluon corrections to the Wilson coefficient \(C_9\) are build in the functions Tperpplus(), Tparplus(), Tparminus(), Cperp() and Cpar(); these corrections have to be integrated to be computed, so the final correction is either obtaind through direct integration in the functions DeltaC9_p(), DeltaC9_m() and DeltaC9_0(), or obtained through fitting in the functions fDeltaC9_p(), fDeltaC9_m() and fDeltaC9_0(). Form factors, Wilson coefficients and parameters are combined together in the functions H_V_0(), H_V_p(), H_V_m(), H_A_0(), H_A_p(), H_A_m(), H_S() and H_P() in order to build the helicity aplitudes, which are consequentely combined to create the angular coefficients in the function I_1c(), I_1s(), I_2c(), I_2s(), I_3(), I_4(), I_5(), I_6c(), I_6s(), I_7(), I_8(), I_9(). Those coefficients are used to create the CP averaged coefficients in the functions getSigma1c(), getSigma1s(), getSigma2c(), getSigma2s(), getSigma3(), getSigma4(), getSigma5(), getSigma6c(), getSigma6s(), getSigma7(), getSigma8(), getSigma9(), and the CP asymmetric coefficients in the function Delta(). The CP averaged and asymmetric coefficients are integrated over the \(q^2\) bin in the functions integrateSigma() and integrateDelta(), in order to be further used to build the observables.

Definition at line 308 of file MVll.h.

Public Member Functions

double beta (double q2)
 The factor \( \beta \) used in the angular coefficients \(I_i\). More...
 
double getgtilde_1_im (double q2)
 The immaginary part of \( \tilde{g}^1 \). More...
 
double getgtilde_1_re (double q2)
 The real part of \( \tilde{g}^1 \). More...
 
double getgtilde_2_im (double q2)
 The immaginary part of \( \tilde{g}^2 \). More...
 
double getgtilde_2_re (double q2)
 The real part of \( \tilde{g}^2 \). More...
 
double getgtilde_3_im (double q2)
 The imaginary part of \( \tilde{g}^3 \). More...
 
double getgtilde_3_re (double q2)
 The real part of \( \tilde{g}^3 \). More...
 
double geth_0_im (double q2)
 The imaginary part of \( h_0 \).
More...
 
double geth_0_re (double q2)
 The real part of \( h_0 \).
More...
 
gslpp::complex geth_m_0 ()
 \( h_-(0) \). More...
 
double geth_m_im (double q2)
 The imaginary part of \( h_- \).
More...
 
double geth_m_re (double q2)
 The real part of \( h_- \).
More...
 
gslpp::complex geth_p_0 ()
 \( h_+(0) \). More...
 
double geth_p_im (double q2)
 The imaginary part of \( h_+ \).
More...
 
double geth_p_re (double q2)
 The real part of \( h_+ \).
More...
 
double getMlep ()
 The mass of the lepton l. More...
 
gslpp::complex getQCDf_1 (double q2)
 
gslpp::complex getQCDf_2 (double q2)
 
gslpp::complex getQCDf_3 (double q2)
 
double getQCDfC9_1 (double q2, double cutoff)
 
double getQCDfC9_2 (double q2, double cutoff)
 
double getQCDfC9_3 (double q2, double cutoff)
 
double getQCDfC9p_1 (double cutoff)
 
double getQCDfC9p_2 (double cutoff)
 
double getQCDfC9p_3 (double cutoff)
 
double getS (double q2)
 The form factor \( S \) . More...
 
double getSigma (int i, double q_2)
 The value of \( \Sigma_{i} \) from \(q_{min}\) to \(q_{max}\). More...
 
double getT0 (double q2)
 The form factor \( T_0 \) . More...
 
double getTm (double q2)
 The form factor \( T_- \) . More...
 
double getTp (double q2)
 The form factor \( T_+ \) . More...
 
double getV0 (double q2)
 The form factor \( V_0 \) . More...
 
double getVm (double q2)
 The form factor \( V_- \) . More...
 
double getVp (double q2)
 The form factor \( V_+ \) . More...
 
double getwidth ()
 The width of the meson M. More...
 
gslpp::complex H_A_0 (double q2, bool bar)
 The helicity amplitude \( H_A^0 \) . More...
 
gslpp::complex H_A_m (double q2, bool bar)
 The helicity amplitude \( H_A^- \) . More...
 
gslpp::complex H_A_p (double q2, bool bar)
 The helicity amplitude \( H_A^+ \) . More...
 
gslpp::complex h_lambda (int hel, double q2)
 The non-pertubative ccbar contributions to the helicity amplitudes. More...
 
gslpp::complex H_P (double q2, bool bar)
 The helicity amplitude \( H_P \) . More...
 
gslpp::complex H_S (double q2, bool bar)
 The helicity amplitude \( H_S \) . More...
 
gslpp::complex H_V_0 (double q2, bool bar)
 The helicity amplitude \( H_V^0 \) . More...
 
gslpp::complex H_V_m (double q2, bool bar)
 The helicity amplitude \( H_V^- \) . More...
 
gslpp::complex H_V_p (double q2, bool bar)
 The helicity amplitude \( H_V^+ \) . More...
 
std::vector< std::string > initializeMVllParameters ()
 A method for initializing the parameters necessary for MVll. 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...
 
 MVll (const StandardModel &SM_i, QCD::meson meson_i, QCD::meson vector_i, QCD::lepton lep_i)
 Constructor. More...
 
virtual ~MVll ()
 Destructor. More...
 

Private Attributes

double ale
 
double angmomV
 
bool dispersion
 
int etaV
 
gslpp::complex exp_Phase [3]
 
bool FixedWCbtos
 
double GF
 
QCD::lepton lep
 
double Mb
 
double mb_pole
 
double Mc
 
double mc_pole
 
QCD::meson meson
 
double mJ2
 
double Mlep
 
double MM
 
double Ms
 
double mu_b
 
double mu_h
 
double MV
 
std::vector< std::string > mvllParameters
 
std::unique_ptr< F_1myF_1
 
std::unique_ptr< F_2myF_2
 
const StandardModelmySM
 
double spectator_charge
 
QCD::meson vectorM
 
double width
 
double xs
 
double ys
 

Constructor & Destructor Documentation

◆ MVll()

MVll::MVll ( const StandardModel SM_i,
QCD::meson  meson_i,
QCD::meson  vector_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]vector_ifinal vector meson of the decay
[in]lep_ifinal leptons of the decay

Definition at line 21 of file MVll.cpp.

22 : mySM(SM_i), myF_1(new F_1()), myF_2(new F_2()),
23 N_cache(3, 0.),
24 V_cache(3, 0.),
25 A0_cache(3, 0.),
26 A1_cache(3, 0.),
27 T1_cache(3, 0.),
28 T2_cache(3, 0.),
29 k2_cache(2, 0.),
30 VL0_cache(3, 0.),
31 TL0_cache(3, 0.),
32 SL_cache(2, 0.),
33 Ycache(2, 0.),
34 H_V0cache(2, 0.),
35 H_V1cache(2, 0.),
36 H_V2cache(2, 0.),
37 H_Scache(2, 0.),
38 H_Pcache(4, 0.),
39 T_cache(5, 0.)
40 {
41  lep = lep_i;
42  meson = meson_i;
43  vectorM = vector_i;
44  dispersion = false;
45  FixedWCbtos = false;
46  mJ2 = 3.096 * 3.096;
47 
48  I0_updated = 0;
49  I1_updated = 0;
50  I2_updated = 0;
51  I3_updated = 0;
52  I4_updated = 0;
53  I5_updated = 0;
54  I6_updated = 0;
55  I7_updated = 0;
56  I8_updated = 0;
57  I9_updated = 0;
58  I10_updated = 0;
59  I11_updated = 0;
60 
61  VL1_updated = 0;
62  VL2_updated = 0;
63  TL1_updated = 0;
64  TL2_updated = 0;
65  VR1_updated = 0;
66  VR2_updated = 0;
67  TR1_updated = 0;
68  TR2_updated = 0;
69  VL0_updated = 0;
70  TL0_updated = 0;
71  VR0_updated = 0;
72  TR0_updated = 0;
73  SL_updated = 0;
74  SR_updated = 0;
75 
76  deltaTparpupdated = 0;
77  deltaTparmupdated = 0;
78  deltaTperpupdated = 0;
79 
80  w_sigma = gsl_integration_cquad_workspace_alloc(100);
81  // w_DTPPR = gsl_integration_cquad_workspace_alloc (100);
82  w_delta = gsl_integration_cquad_workspace_alloc(100);
83 
84  acc_Re_T_perp = gsl_interp_accel_alloc();
85  acc_Im_T_perp = gsl_interp_accel_alloc();
86  acc_Re_T_para = gsl_interp_accel_alloc();
87  acc_Im_T_para = gsl_interp_accel_alloc();
88 
89  spline_Re_T_perp = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM);
90  spline_Im_T_perp = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM);
91  spline_Re_T_para = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM);
92  spline_Im_T_para = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM);
93 
94 #if COMPUTECP
95  acc_Re_T_perp_conj = gsl_interp_accel_alloc();
96  acc_Im_T_perp_conj = gsl_interp_accel_alloc();
97  acc_Re_T_para_conj = gsl_interp_accel_alloc();
98  acc_Im_T_para_conj = gsl_interp_accel_alloc();
99 
100  spline_Re_T_perp_conj = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM);
101  spline_Im_T_perp_conj = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM);
102  spline_Re_T_para_conj = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM);
103  spline_Im_T_para_conj = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM);
104 #endif
105 
106  acc_Re_deltaC7_QCDF = gsl_interp_accel_alloc();
107  acc_Im_deltaC7_QCDF = gsl_interp_accel_alloc();
108  acc_Re_deltaC9_QCDF = gsl_interp_accel_alloc();
109  acc_Im_deltaC9_QCDF = gsl_interp_accel_alloc();
110 
111  spline_Re_deltaC7_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
112  spline_Im_deltaC7_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
113  spline_Re_deltaC9_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
114  spline_Im_deltaC9_QCDF = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
115 
116 #if COMPUTECP
117  acc_Re_deltaC7_QCDF_conj = gsl_interp_accel_alloc();
118  acc_Im_deltaC7_QCDF_conj = gsl_interp_accel_alloc();
119  acc_Re_deltaC9_QCDF_conj = gsl_interp_accel_alloc();
120  acc_Im_deltaC9_QCDF_conj = gsl_interp_accel_alloc();
121 
122  spline_Re_deltaC7_QCDF_conj = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
123  spline_Im_deltaC7_QCDF_conj = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
124  spline_Re_deltaC9_QCDF_conj = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
125  spline_Im_deltaC9_QCDF_conj = gsl_spline_alloc(gsl_interp_cspline, GSL_INTERP_DIM_DC);
126 #endif
127 
128  h_pole = false;
129 
130  M_PI2 = M_PI*M_PI;
131 
132  F87_1 = (4. / 3. * M_PI2 - 40. / 3.);
133  F87_2 = (32. / 9. * M_PI2 - 316. / 9.);
134  F87_3 = (200. / 27. * M_PI2 - 658. / 9.);
135 
136  F89_0 = (104. / 9. - 32. / 27. * M_PI2);
137  F89_1 = (1184. / 27. - 40. / 9. * M_PI2);
138  F89_2 = (-32. / 3. * M_PI2 + 14212. / 135.);
139  F89_3 = (-560. / 27. * M_PI2 + 193444. / 945.);
140 
141  CF = 4. / 3.;
142 
143 }

◆ ~MVll()

MVll::~MVll ( )
virtual

Destructor.

Definition at line 145 of file MVll.cpp.

146 {
147 }

Member Function Documentation

◆ beta()

double MVll::beta ( double  q2)

The factor \( \beta \) used in the angular coefficients \(I_i\).

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

Definition at line 1681 of file MVll.cpp.

1681 {
1682  return (MM4 + q2 * q2 + MV4 - twoMV2 * q2 - twoMM2 * (q2 + MV2)) / fourMM2;
1683 }
1684 

◆ getgtilde_1_im()

double MVll::getgtilde_1_im ( double  q2)
inline

The immaginary part of \( \tilde{g}^1 \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{IM}(\tilde{g}^1) \)

Definition at line 589 of file MVll.h.

590  {
591  updateParameters();
592  return C2_inv * (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * (h_lambda(2,q2)-h_lambda(1,q2))).imag()/q2;
593  }

◆ getgtilde_1_re()

double MVll::getgtilde_1_re ( double  q2)
inline

The real part of \( \tilde{g}^1 \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{RE}(\tilde{g}^1) \)

Definition at line 578 of file MVll.h.

579  {
580  updateParameters();
581  return C2_inv * (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * (h_lambda(2,q2)-h_lambda(1,q2))).real()/q2;
582  }

◆ getgtilde_2_im()

double MVll::getgtilde_2_im ( double  q2)
inline

The immaginary part of \( \tilde{g}^2 \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{IM}(\tilde{g}^2) \)

Definition at line 611 of file MVll.h.

612  {
613  updateParameters();
614  return C2_inv * (gtilde_2_pre/A_1(q2) * (h_lambda(1,q2)+h_lambda(2,q2))).imag()/q2;
615  }

◆ getgtilde_2_re()

double MVll::getgtilde_2_re ( double  q2)
inline

The real part of \( \tilde{g}^2 \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{RE}(\tilde{g}^2) \)

Definition at line 600 of file MVll.h.

601  {
602  updateParameters();
603  return C2_inv * (gtilde_2_pre/A_1(q2) * (h_lambda(1,q2)+h_lambda(2,q2))).real()/q2;
604  }

◆ getgtilde_3_im()

double MVll::getgtilde_3_im ( double  q2)
inline

The imaginary part of \( \tilde{g}^3 \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{IM}(\tilde{g}^3) \)

Definition at line 634 of file MVll.h.

635  {
636  updateParameters();
637  return C2_inv * (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2)*h_lambda(0,q2)/q2-
638  (MM2mMV2 - q2)/(4.*MV) * (h_lambda(1,q2)+h_lambda(2,q2))/q2)).imag();
639  }

◆ getgtilde_3_re()

double MVll::getgtilde_3_re ( double  q2)
inline

The real part of \( \tilde{g}^3 \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{RE}(\tilde{g}^3) \)

Definition at line 622 of file MVll.h.

623  {
624  updateParameters();
625  return C2_inv * (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2)*h_lambda(0,q2)/q2-
626  (MM2mMV2 - q2)/(4.*MV) * (h_lambda(1,q2)+h_lambda(2,q2))/q2)).real();
627  }

◆ geth_0_im()

double MVll::geth_0_im ( double  q2)
inline

The imaginary part of \( h_0 \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{IM}(h_0) \)

Definition at line 656 of file MVll.h.

657  {
658  return (sixteenM_PI2MM2 * h_lambda(0,q2)/q2).imag();
659  }

◆ geth_0_re()

double MVll::geth_0_re ( double  q2)
inline

The real part of \( h_0 \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{RE}(h_0) \)

Definition at line 646 of file MVll.h.

647  {
648  return (sixteenM_PI2MM2 * h_lambda(0,q2)/q2).real();
649  }

◆ geth_m_0()

gslpp::complex MVll::geth_m_0 ( )
inline

\( h_-(0) \).

Returns
\( h_-(0) \)

Definition at line 694 of file MVll.h.

695  {
696  return h_lambda(2,0.);
697  }

◆ geth_m_im()

double MVll::geth_m_im ( double  q2)
inline

The imaginary part of \( h_- \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{IM}(h_-) \)

Definition at line 714 of file MVll.h.

715  {
716  return (sixteenM_PI2MM2 * h_lambda(2,q2)/q2).imag();
717  }

◆ geth_m_re()

double MVll::geth_m_re ( double  q2)
inline

The real part of \( h_- \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{RE}(h_-) \)

Definition at line 704 of file MVll.h.

705  {
706  return (sixteenM_PI2MM2 * h_lambda(2,q2)/q2).real();
707  }

◆ geth_p_0()

gslpp::complex MVll::geth_p_0 ( )
inline

\( h_+(0) \).

Returns
\( h_+(0) \)

Definition at line 665 of file MVll.h.

666  {
667  return h_lambda(1,0.);
668  }

◆ geth_p_im()

double MVll::geth_p_im ( double  q2)
inline

The imaginary part of \( h_+ \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{IM}(h_+) \)

Definition at line 685 of file MVll.h.

686  {
687  return (sixteenM_PI2MM2 * h_lambda(1,q2)/q2).imag();
688  }

◆ geth_p_re()

double MVll::geth_p_re ( double  q2)
inline

The real part of \( h_+ \).

Parameters
[in]q2\(q^2\) of the decay
Returns
\( \mathrm{RE}(h_+) \)

Definition at line 675 of file MVll.h.

676  {
677  return (sixteenM_PI2MM2 * h_lambda(1,q2)/q2).real();
678  }

◆ getMlep()

double MVll::getMlep ( )
inline

The mass of the lepton l.

Returns
\( m_l \)

Definition at line 364 of file MVll.h.

364  {
365  updateParameters();
366  return Mlep;
367  }

◆ getQCDf_1()

gslpp::complex MVll::getQCDf_1 ( double  q2)
inline

Definition at line 516 of file MVll.h.

517  {
518  updateParameters();
519 // return (gtilde_1_pre/(sqrt(lambda(q2)) * V(q2)) * 1./(16. * M_PI * M_PI * MM*MM) * (fDeltaC9_m(q2) * V_m(q2) - fDeltaC9_p(q2) * V_p(q2)));
520  return 0.;
521  }

◆ getQCDf_2()

gslpp::complex MVll::getQCDf_2 ( double  q2)
inline

Definition at line 523 of file MVll.h.

524  {
525  updateParameters();
526 // return (gtilde_2_pre/A_1(q2) * 1./(16. * M_PI * M_PI * MM*MM) * (fDeltaC9_m(q2) * V_m(q2) + fDeltaC9_p(q2) * V_p(q2)));
527  return 0.;
528  }

◆ getQCDf_3()

gslpp::complex MVll::getQCDf_3 ( double  q2)
inline

Definition at line 530 of file MVll.h.

531  {
532  updateParameters();
533 // return (gtilde_3_pre/(lambda(q2) * A_2(q2)) * (sqrt(q2) * 1./(16. * M_PI * M_PI * MM*MM) * fDeltaC9_0(q2) * V_0t(q2) - (MM2mMV2 - q2)/(4.*MV) * 1./(16. * M_PI * M_PI * MM*MM) * (fDeltaC9_m(q2) * V_m(q2) + fDeltaC9_p(q2) * V_p(q2))));
534  return 0.;
535  }

◆ getQCDfC9_1()

double MVll::getQCDfC9_1 ( double  q2,
double  cutoff 
)
inline

Definition at line 537 of file MVll.h.

538  {
539  updateParameters();
540  return (getQCDf_1(q2) - getQCDf_1(cutoff) * cutoff/q2).abs();
541  }

◆ getQCDfC9_2()

double MVll::getQCDfC9_2 ( double  q2,
double  cutoff 
)
inline

Definition at line 543 of file MVll.h.

544  {
545  updateParameters();
546  return (getQCDf_2(q2) - getQCDf_2(cutoff) * cutoff/q2).abs();
547  }

◆ getQCDfC9_3()

double MVll::getQCDfC9_3 ( double  q2,
double  cutoff 
)
inline

Definition at line 549 of file MVll.h.

550  {
551  updateParameters();
552  return (getQCDf_3(q2) - getQCDf_3(cutoff) * cutoff/q2).abs();
553  }

◆ getQCDfC9p_1()

double MVll::getQCDfC9p_1 ( double  cutoff)
inline

Definition at line 555 of file MVll.h.

556  {
557  updateParameters();
558  return (getQCDf_1(cutoff) * cutoff).abs();
559  }

◆ getQCDfC9p_2()

double MVll::getQCDfC9p_2 ( double  cutoff)
inline

Definition at line 561 of file MVll.h.

562  {
563  updateParameters();
564  return (getQCDf_2(cutoff) * cutoff).abs();
565  }

◆ getQCDfC9p_3()

double MVll::getQCDfC9p_3 ( double  cutoff)
inline

Definition at line 567 of file MVll.h.

568  {
569  updateParameters();
570  return (getQCDf_3(cutoff) * cutoff).abs();
571  }

◆ getS()

double MVll::getS ( double  q2)
inline

The form factor \( S \) .

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

Definition at line 446 of file MVll.h.

447  {
448  updateParameters();
449  return S_L_pre/sqrt(lambda(q2)) * S_L(q2);
450  };

◆ getSigma()

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

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

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

Definition at line 1956 of file MVll.cpp.

1961 {
1962  updateParameters();
1963 
1964  switch (i) {
1965  case 0:
1966  return getSigma1c(q_2);
1967  break;
1968  case 1:
1969  return getSigma1s(q_2);
1970  break;
1971  case 2:
1972  return getSigma2c(q_2);
1973  break;
1974  case 3:
1975  return getSigma2s(q_2);
1976  break;
1977  case 4:
1978  return getSigma3(q_2);
1979  break;
1980  case 5:
1981  return getSigma4(q_2);
1982  break;
1983  case 6:
1984  return getSigma5(q_2);
1985  break;
1986  case 7:
1987  return getSigma6s(q_2);
1988  break;
1989  case 9:
1990  return getSigma7(q_2);
1991  break;
1992  case 10:
1993  return getSigma8(q_2);
1994  break;
1995  case 11:
1996  return getSigma9(q_2);
1997  break;
1998  default:
1999  std::stringstream out;

◆ getT0()

double MVll::getT0 ( double  q2)
inline

The form factor \( T_0 \) .

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

Definition at line 413 of file MVll.h.

414  {
415  updateParameters();
416  return twoMM3/sqrt(q2 * lambda(q2)) * T_0t(q2);
417  };

◆ getTm()

double MVll::getTm ( double  q2)
inline

The form factor \( T_- \) .

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

Definition at line 435 of file MVll.h.

436  {
437  updateParameters();
438  return T_m(q2);
439  };

◆ getTp()

double MVll::getTp ( double  q2)
inline

The form factor \( T_+ \) .

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

Definition at line 424 of file MVll.h.

425  {
426  updateParameters();
427  return T_p(q2);
428  };

◆ getV0()

double MVll::getV0 ( double  q2)
inline

The form factor \( V_0 \) .

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

Definition at line 381 of file MVll.h.

382  {
383  updateParameters();
384  return (2. * MM * sqrt(q2))/sqrt(lambda(q2)) * V_0t(q2);
385  };

◆ getVm()

double MVll::getVm ( double  q2)
inline

The form factor \( V_- \) .

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

Definition at line 403 of file MVll.h.

404  {
405  return V_m(q2);
406  };

◆ getVp()

double MVll::getVp ( double  q2)
inline

The form factor \( V_+ \) .

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

Definition at line 392 of file MVll.h.

393  {
394  updateParameters();
395  return V_p(q2);
396  };

◆ getwidth()

double MVll::getwidth ( )
inline

The width of the meson M.

Returns
\( \Gamma_M \)

Definition at line 355 of file MVll.h.

355  {
356  updateParameters();
357  return width;
358  }

◆ H_A_0()

gslpp::complex MVll::H_A_0 ( double  q2,
bool  bar 
)

The helicity amplitude \( H_A^0 \) .

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

Definition at line 1644 of file MVll.cpp.

1648 {

◆ H_A_m()

gslpp::complex MVll::H_A_m ( double  q2,
bool  bar 
)

The helicity amplitude \( H_A^- \) .

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

Definition at line 1656 of file MVll.cpp.

1660 {

◆ H_A_p()

gslpp::complex MVll::H_A_p ( double  q2,
bool  bar 
)

The helicity amplitude \( H_A^+ \) .

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

Definition at line 1650 of file MVll.cpp.

1654 {

◆ h_lambda()

gslpp::complex MVll::h_lambda ( int  hel,
double  q2 
)

The non-pertubative ccbar contributions to the helicity amplitudes.

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

Definition at line 1602 of file MVll.cpp.

1606 {
1607  if (!dispersion) {
1608  if (h_pole == true) return SU3_breaking * (h_0[hel]+(1. - h_2[hel]) * q2 * (h_1[hel] - h_0[hel]) / (q2 - h_2[hel]));
1609  else if (hel == 1) return SU3_breaking * (h_0[1] + h_1[1] * q2 + h_2[1] * q2 * q2 + (twoMboMM * h_0[2] * T_p(q2) + h_1[2] * q2 / MM2 * V_p(q2)) / sixteenM_PI2);
1610  else if (hel == 2) return SU3_breaking * (twoMboMM * h_0[2] * T_m(q2) + h_1[2] * q2 / MM2 * V_m(q2)) / sixteenM_PI2 + h_2[2] * q2 * q2;
1611  else return SU3_breaking * ((h_0[hel] + h_1[hel] * q2) * sqrt(q2) + (twoMboMM * h_0[2] * T_0t(q2) + h_1[2] * q2 * V_0t(q2) / MM2) / sixteenM_PI2);
1612  } else {
1613  if (hel == 0) return SU3_breaking * (-sqrt(q2) / (MM2 * 16. * M_PI * M_PI) * ((MMpMV2 * (MM2mMV2 - q2) * A_1(q2) * DeltaC9_KD(q2, 1) - lambda(q2) * A_2(q2) * DeltaC9_KD(q2, 2)) / (4. * MV * MM * MMpMV)));
1614  else if (hel == 1) {
1615  if (q2 == 0.) return SU3_breaking * (-1. / (MM2 * 16. * M_PI * M_PI) * (
1616  (MMpMV * A_1(0.)) / (2. * MM) * ((-h_0[1] + h_2[1]) / (1. + h_1[1] / mJ2)) * exp_Phase[1]
1617  - sqrt(lambda(0.)) / (2. * MM * MMpMV) * V(0.) * ((-h_0[0] + h_2[0]) / (1. + h_1[0] / mJ2)) * exp_Phase[1]));
1618  else return SU3_breaking * (-q2 / (MM2 * 16. * M_PI * M_PI) * ((MMpMV * A_1(q2)) / (2. * MM) * DeltaC9_KD(q2, 1) - sqrt(lambda(q2)) / (2. * MM * MMpMV) * V(q2) * DeltaC9_KD(q2, 0)));
1619  } else {
1620  if (q2 == 0.) return SU3_breaking * (-1. / (MM2 * 16. * M_PI * M_PI) *
1621  ((MMpMV * A_1(0.)) / (2. * MM) * ((-h_0[1] + h_2[1]) / (1. + h_1[1] / mJ2)) * exp_Phase[1]
1622  + sqrt(lambda(0.)) / (2. * MM * MMpMV) * V(0.) * ((-h_0[0] + h_2[0]) / (1. + h_1[0] / mJ2)) * exp_Phase[1]));
1623  else return SU3_breaking * (-q2 / (MM2 * 16. * M_PI * M_PI) * ((MMpMV * A_1(q2)) / (2. * MM) * DeltaC9_KD(q2, 1) + sqrt(lambda(q2)) / (2. * MM * MMpMV) * V(q2) * DeltaC9_KD(q2, 0)));

◆ H_P()

gslpp::complex MVll::H_P ( double  q2,
bool  bar 
)

The helicity amplitude \( H_P \) .

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

Definition at line 1668 of file MVll.cpp.

1672 {

◆ H_S()

gslpp::complex MVll::H_S ( double  q2,
bool  bar 
)

The helicity amplitude \( H_S \) .

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

Definition at line 1662 of file MVll.cpp.

1666 {

◆ H_V_0()

gslpp::complex MVll::H_V_0 ( double  q2,
bool  bar 
)

The helicity amplitude \( H_V^0 \) .

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

Definition at line 1625 of file MVll.cpp.

1629 {
1630  if (!bar) return -gslpp::complex::i() * NN * (((C_9 + deltaC9_QCDF(q2, !bar, SPLINE) /*+ fDeltaC9_0(q2)*/ + Y(q2)) - etaV * pow(-1, angmomV) * C_9p) * V_0t(q2) + T_0(q2, !bar) + MM2 / q2 * (twoMboMM * (C_7 + deltaC7_QCDF(q2, !bar, SPLINE) - etaV * pow(-1, angmomV) * C_7p) * T_0t(q2) - sixteenM_PI2 * h_lambda(0, q2)));

◆ H_V_m()

gslpp::complex MVll::H_V_m ( double  q2,
bool  bar 
)

The helicity amplitude \( H_V^- \) .

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

Definition at line 1638 of file MVll.cpp.

1642 {

◆ H_V_p()

gslpp::complex MVll::H_V_p ( double  q2,
bool  bar 
)

The helicity amplitude \( H_V^+ \) .

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

Definition at line 1632 of file MVll.cpp.

1636 {

◆ initializeMVllParameters()

std::vector< std::string > MVll::initializeMVllParameters ( )

A method for initializing the parameters necessary for MVll.

Returns
the vector of MVll specific parameters

Definition at line 149 of file MVll.cpp.

150 {
153 
154 #if NFPOLARBASIS_MVLL
156  << "a_0Vphi" << "a_1Vphi" << "a_2Vphi" << "MRV" << "a_0A0phi" << "a_1A0phi" << "a_2A0phi" << "MRA0"
157  << "a_0A1phi" << "a_1A1phi" << "a_2A1phi" << "MRA1" << "a_1A12phi" << "a_2A12phi" << "MRA12" /*a_0A12 and a_0T2 are not independent*/
158  << "a_0T1phi" << "a_1T1phi" << "a_2T1phi" << "MRT1" << "a_1T2phi" << "a_2T2phi" << "MRT2"
159  << "a_0T23phi" << "a_1T23phi" << "a_2T23phi" << "MRT23"
160  << "absh_0" << "absh_p" << "absh_m" << "argh_0" << "argh_p" << "argh_m"
161  << "absh_0_1" << "absh_p_1" << "absh_m_1" << "argh_0_1" << "argh_p_1" << "argh_m_1"
162  << "absh_p_2" << "absh_m_2" << "argh_p_2" << "argh_m_2" << "xs_phi" << "SU3_breaking_abs" << "SU3_breaking_arg";
164  << "a_0V" << "a_1V" << "a_2V" << "MRV" << "a_0A0" << "a_1A0" << "a_2A0" << "MRA0"
165  << "a_0A1" << "a_1A1" << "a_2A1" << "MRA1" << "a_1A12" << "a_2A12" << "MRA12" /*a_0A12 and a_0T2 are not independent*/
166  << "a_0T1" << "a_1T1" << "a_2T1" << "MRT1" << "a_1T2" << "a_2T2" << "MRT2"
167  << "a_0T23" << "a_1T23" << "a_2T23" << "MRT23"
168  << "absh_0" << "absh_p" << "absh_m" << "argh_0" << "argh_p" << "argh_m"
169  << "absh_0_1" << "absh_p_1" << "absh_m_1" << "argh_0_1" << "argh_p_1" << "argh_m_1"
170  << "absh_p_2" << "absh_m_2" << "argh_p_2" << "argh_m_2";
171 #else
173  << "a_0Vphi" << "a_1Vphi" << "a_2Vphi" << "MRV" << "a_0A0phi" << "a_1A0phi" << "a_2A0phi" << "MRA0"
174  << "a_0A1phi" << "a_1A1phi" << "a_2A1phi" << "MRA1" << "a_1A12phi" << "a_2A12phi" << "MRA12" /*a_0A12 and a_0T2 are not independent*/
175  << "a_0T1phi" << "a_1T1phi" << "a_2T1phi" << "MRT1" << "a_1T2phi" << "a_2T2phi" << "MRT2"
176  << "a_0T23phi" << "a_1T23phi" << "a_2T23phi" << "MRT23"
177  << "reh_0" << "reh_p" << "reh_m" << "imh_0" << "imh_p" << "imh_m"
178  << "reh_0_1" << "reh_p_1" << "reh_m_1" << "imh_0_1" << "imh_p_1" << "imh_m_1"
179  << "reh_p_2" << "reh_m_2" << "imh_p_2" << "imh_m_2" << "xs_phi" << "SU3_breaking_abs" << "SU3_breaking_arg";
181  << "a_0V" << "a_1V" << "a_2V" << "MRV" << "a_0A0" << "a_1A0" << "a_2A0" << "MRA0"
182  << "a_0A1" << "a_1A1" << "a_2A1" << "MRA1" << "a_1A12" << "a_2A12" << "MRA12" /*a_0A12 and a_0T2 are not independent*/
183  << "a_0T1" << "a_1T1" << "a_2T1" << "MRT1" << "a_1T2" << "a_2T2" << "MRT2"
184  << "a_0T23" << "a_1T23" << "a_2T23" << "MRT23"
185  << "reh_0" << "reh_p" << "reh_m" << "imh_0" << "imh_p" << "imh_m"
186  << "reh_0_1" << "reh_p_1" << "reh_m_1" << "imh_0_1" << "imh_p_1" << "imh_m_1"
187  << "reh_p_2" << "reh_m_2" << "imh_p_2" << "imh_m_2";
188 #endif
189  else {
190  std::stringstream out;
191  out << vectorM;
192  throw std::runtime_error("MVll: vector " + out.str() + " not implemented");
193  }
194 
195  if (dispersion) {
196  mvllParameters.clear();
198  << "a_0Vphi" << "a_1Vphi" << "a_2Vphi" << "MRV" << "a_0A0phi" << "a_1A0phi" << "a_2A0phi" << "MRA0"
199  << "a_0A1phi" << "a_1A1phi" << "a_2A1phi" << "MRA1" << "a_1A12phi" << "a_2A12phi" << "MRA12" /*a_0A12 and a_0T2 are not independent*/
200  << "a_0T1phi" << "a_1T1phi" << "a_2T1phi" << "MRT1" << "a_1T2phi" << "a_2T2phi" << "MRT2"
201  << "a_0T23phi" << "a_1T23phi" << "a_2T23phi" << "MRT23"
202  << "r1_1" << "r2_1" << "deltaC9_1" << "phDC9_1"
203  << "r1_2" << "r2_2" << "deltaC9_2" << "phDC9_2"
204  << "r1_3" << "r2_3" << "deltaC9_3" << "phDC9_3" << "xs_phi" << "SU3_breaking_abs" << "SU3_breaking_arg";
206  << "a_0V" << "a_1V" << "a_2V" << "MRV" << "a_0A0" << "a_1A0" << "a_2A0" << "MRA0"
207  << "a_0A1" << "a_1A1" << "a_2A1" << "MRA1" << "a_1A12" << "a_2A12" << "MRA12" /*a_0A12 and a_0T2 are not independent*/
208  << "a_0T1" << "a_1T1" << "a_2T1" << "MRT1" << "a_1T2" << "a_2T2" << "MRT2"
209  << "a_0T23" << "a_1T23" << "a_2T23" << "MRT23"
210  << "r1_1" << "r2_1" << "deltaC9_1" << "phDC9_1"
211  << "r1_2" << "r2_2" << "deltaC9_2" << "phDC9_2"
212  << "r1_3" << "r2_3" << "deltaC9_3" << "phDC9_3";
213  }
214 
215  if (FixedWCbtos) mvllParameters.insert(mvllParameters.end(), { "C7_SM", "C9_SM", "C10_SM" });
216 
219  return mvllParameters;
220 }

◆ integrateDelta()

double MVll::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 2001 of file MVll.cpp.

2001  : index " + out.str() + " not implemented");
2002  }
2003 }
2004 
2005 double MVll::integrateDelta(int i, double q_min, double q_max)
2006 {
2007  updateParameters();
2008 
2009  std::pair<double, double > qbin = std::make_pair(q_min, q_max);
2010 
2011  old_handler = gsl_set_error_handler_off();
2012 
2013  switch (i) {
2014  case 0:
2015  if (delta0Cached[qbin] == 0) {
2016  FD = convertToGslFunction(boost::bind(&MVll::getDelta1c, &(*this), _1));
2017  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();
2018  cacheDelta0[qbin] = avaDelta;
2019  delta0Cached[qbin] = 1;
2020  }
2021  return cacheDelta0[qbin];
2022  break;
2023  case 1:
2024  if (delta1Cached[qbin] == 0) {
2025  FD = convertToGslFunction(boost::bind(&MVll::getDelta1s, &(*this), _1));
2026  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();
2027  cacheDelta1[qbin] = avaDelta;
2028  delta1Cached[qbin] = 1;
2029  }
2030  return cacheDelta1[qbin];
2031  break;
2032  case 2:
2033  if (delta2Cached[qbin] == 0) {
2034  FD = convertToGslFunction(boost::bind(&MVll::getDelta2c, &(*this), _1));
2035  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();
2036  cacheDelta2[qbin] = avaDelta;
2037  delta2Cached[qbin] = 1;
2038  }
2039  return cacheDelta2[qbin];
2040  break;
2041  case 3:
2042  if (delta3Cached[qbin] == 0) {
2043  FD = convertToGslFunction(boost::bind(&MVll::getDelta2s, &(*this), _1));
2044  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();
2045  cacheDelta3[qbin] = avaDelta;
2046  delta3Cached[qbin] = 1;
2047  }
2048  return cacheDelta3[qbin];
2049  break;
2050  case 6:
2051  if (delta6Cached[qbin] == 0) {
2052  FD = convertToGslFunction(boost::bind(&MVll::getDelta5, &(*this), _1));
2053  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();
2054  cacheDelta6[qbin] = avaDelta;
2055  delta6Cached[qbin] = 1;
2056  }
2057  return cacheDelta6[qbin];
2058  break;
2059  case 7:
2060  if (delta7Cached[qbin] == 0) {
2061  FD = convertToGslFunction(boost::bind(&MVll::getDelta6s, &(*this), _1));
2062  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();
2063  cacheDelta7[qbin] = avaDelta;
2064  delta7Cached[qbin] = 1;
2065  }
2066  return cacheDelta7[qbin];
2067  break;
2068  case 8:
2069  if (delta8Cached[qbin] == 0) {
2070  FD = convertToGslFunction(boost::bind(&MVll::getDelta6c, &(*this), _1));
2071  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();
2072  cacheDelta8[qbin] = avaDelta;
2073  delta8Cached[qbin] = 1;
2074  }
2075  return cacheDelta8[qbin];
2076  break;
2077  case 10:
2078  if (delta10Cached[qbin] == 0) {
2079  FD = convertToGslFunction(boost::bind(&MVll::getDelta8, &(*this), _1));
2080  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();
2081  cacheDelta10[qbin] = avaDelta;
2082  delta10Cached[qbin] = 1;
2083  }
2084  return cacheDelta10[qbin];
2085  break;
2086  case 11:
2087  if (delta11Cached[qbin] == 0) {
2088  FD = convertToGslFunction(boost::bind(&MVll::getDelta9, &(*this), _1));
2089  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();
2090  cacheDelta11[qbin] = avaDelta;
2091  delta11Cached[qbin] = 1;
2092  }
2093  return cacheDelta11[qbin];
2094  break;
2095  default:
2096  std::stringstream out;
2097  out << i;
2098  throw std::runtime_error("MVll::integrateDelta: index " + out.str() + " not implemented");
2099  }

◆ integrateSigma()

double MVll::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 1838 of file MVll.cpp.

1843 {
1844  updateParameters();
1845 
1846  std::pair<double, double > qbin = std::make_pair(q_min, q_max);
1847 
1848  old_handler = gsl_set_error_handler_off();
1849 
1850  switch (i) {
1851  case 0:
1852  if (sigma0Cached[qbin] == 0) {
1853  FS = convertToGslFunction(boost::bind(&MVll::getSigma1c, &(*this), _1));
1854  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();
1855  cacheSigma0[qbin] = avaSigma;
1856  sigma0Cached[qbin] = 1;
1857  }
1858  return cacheSigma0[qbin];
1859  break;
1860  case 1:
1861  if (sigma1Cached[qbin] == 0) {
1862  FS = convertToGslFunction(boost::bind(&MVll::getSigma1s, &(*this), _1));
1863  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();
1864  cacheSigma1[qbin] = avaSigma;
1865  sigma1Cached[qbin] = 1;
1866  }
1867  return cacheSigma1[qbin];
1868  break;
1869  case 2:
1870  if (sigma2Cached[qbin] == 0) {
1871  FS = convertToGslFunction(boost::bind(&MVll::getSigma2c, &(*this), _1));
1872  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();
1873  cacheSigma2[qbin] = avaSigma;
1874  sigma2Cached[qbin] = 1;
1875  }
1876  return cacheSigma2[qbin];
1877  break;
1878  case 3:
1879  if (sigma3Cached[qbin] == 0) {
1880  FS = convertToGslFunction(boost::bind(&MVll::getSigma2s, &(*this), _1));
1881  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();
1882  cacheSigma3[qbin] = avaSigma;
1883  sigma3Cached[qbin] = 1;
1884  }
1885  return cacheSigma3[qbin];
1886  break;
1887  case 4:
1888  if (sigma4Cached[qbin] == 0) {
1889  FS = convertToGslFunction(boost::bind(&MVll::getSigma3, &(*this), _1));
1890  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();
1891  cacheSigma4[qbin] = avaSigma;
1892  sigma4Cached[qbin] = 1;
1893  }
1894  return cacheSigma4[qbin];
1895  break;
1896  case 5:
1897  if (sigma5Cached[qbin] == 0) {
1898  FS = convertToGslFunction(boost::bind(&MVll::getSigma4, &(*this), _1));
1899  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();
1900  cacheSigma5[qbin] = avaSigma;
1901  sigma5Cached[qbin] = 1;
1902  }
1903  return cacheSigma5[qbin];
1904  break;
1905  case 6:
1906  if (sigma6Cached[qbin] == 0) {
1907  FS = convertToGslFunction(boost::bind(&MVll::getSigma5, &(*this), _1));
1908  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();
1909  cacheSigma6[qbin] = avaSigma;
1910  sigma6Cached[qbin] = 1;
1911  }
1912  return cacheSigma6[qbin];
1913  break;
1914  case 7:
1915  if (sigma7Cached[qbin] == 0) {
1916  FS = convertToGslFunction(boost::bind(&MVll::getSigma6s, &(*this), _1));
1917  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();
1918  cacheSigma7[qbin] = avaSigma;
1919  sigma7Cached[qbin] = 1;
1920  }
1921  return cacheSigma7[qbin];
1922  break;
1923  case 9:
1924  if (sigma9Cached[qbin] == 0) {
1925  FS = convertToGslFunction(boost::bind(&MVll::getSigma7, &(*this), _1));
1926  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();
1927  cacheSigma9[qbin] = avaSigma;
1928  sigma9Cached[qbin] = 1;
1929  }
1930  return cacheSigma9[qbin];
1931  break;
1932  case 10:
1933  if (sigma10Cached[qbin] == 0) {
1934  FS = convertToGslFunction(boost::bind(&MVll::getSigma8, &(*this), _1));
1935  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();
1936  cacheSigma10[qbin] = avaSigma;
1937  sigma10Cached[qbin] = 1;
1938  }
1939  return cacheSigma10[qbin];
1940  break;
1941  case 11:
1942  if (sigma11Cached[qbin] == 0) {
1943  FS = convertToGslFunction(boost::bind(&MVll::getSigma9, &(*this), _1));
1944  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();
1945  cacheSigma11[qbin] = avaSigma;
1946  sigma11Cached[qbin] = 1;
1947  }
1948  return cacheSigma11[qbin];
1949  break;
1950  default:
1951  std::stringstream out;
1952  out << i;
1953  throw std::runtime_error("MVll::integrateSigma: index " + out.str() + " not implemented");
1954  }

Member Data Documentation

◆ ale

double MVll::ale
private

Alpha electromagnetic

Definition at line 739 of file MVll.h.

◆ angmomV

double MVll::angmomV
private

angular momentum of meson V; for a resonance, it's replaced by its spin

Definition at line 754 of file MVll.h.

◆ dispersion

bool MVll::dispersion
private

Definition at line 733 of file MVll.h.

◆ etaV

int MVll::etaV
private

parity of meson V

Definition at line 755 of file MVll.h.

◆ exp_Phase

gslpp::complex MVll::exp_Phase[3]
private

Definition at line 736 of file MVll.h.

◆ FixedWCbtos

bool MVll::FixedWCbtos
private

Definition at line 734 of file MVll.h.

◆ GF

double MVll::GF
private

Fermi constant

Definition at line 738 of file MVll.h.

◆ lep

QCD::lepton MVll::lep
private

Final leptons type

Definition at line 727 of file MVll.h.

◆ Mb

double MVll::Mb
private

b quark mass

Definition at line 743 of file MVll.h.

◆ mb_pole

double MVll::mb_pole
private

b quark pole mass

Definition at line 747 of file MVll.h.

◆ Mc

double MVll::Mc
private

c quark mass

Definition at line 746 of file MVll.h.

◆ mc_pole

double MVll::mc_pole
private

c quark pole mass

Definition at line 748 of file MVll.h.

◆ meson

QCD::meson MVll::meson
private

Initial meson type

Definition at line 728 of file MVll.h.

◆ mJ2

double MVll::mJ2
private

Definition at line 735 of file MVll.h.

◆ Mlep

double MVll::Mlep
private

Muon mass

Definition at line 740 of file MVll.h.

◆ MM

double MVll::MM
private

Initial meson mass

Definition at line 741 of file MVll.h.

◆ Ms

double MVll::Ms
private

s quark mass

Definition at line 749 of file MVll.h.

◆ mu_b

double MVll::mu_b
private

b mass scale

Definition at line 744 of file MVll.h.

◆ mu_h

double MVll::mu_h
private

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

Definition at line 745 of file MVll.h.

◆ MV

double MVll::MV
private

Final vector meson mass

Definition at line 742 of file MVll.h.

◆ mvllParameters

std::vector<std::string> MVll::mvllParameters
private

The string of mandatory MVll parameters

Definition at line 730 of file MVll.h.

◆ myF_1

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

Definition at line 731 of file MVll.h.

◆ myF_2

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

Definition at line 732 of file MVll.h.

◆ mySM

const StandardModel& MVll::mySM
private

Model type

Definition at line 726 of file MVll.h.

◆ spectator_charge

double MVll::spectator_charge
private

spectator quark charge

Definition at line 750 of file MVll.h.

◆ vectorM

QCD::meson MVll::vectorM
private

Final vector meson type

Definition at line 729 of file MVll.h.

◆ width

double MVll::width
private

Initial meson width

Definition at line 751 of file MVll.h.

◆ xs

double MVll::xs
private

CP-violation factor \(\frac{\Delta m}{\Gamma}\)

Definition at line 753 of file MVll.h.

◆ ys

double MVll::ys
private

CP-violation factor \(\frac{\Delta \Gamma}{2\Gamma}\)

Definition at line 752 of file MVll.h.


The documentation for this class was generated from the following files:
MVll::vectorM
QCD::meson vectorM
Definition: MVll.h:729
MVll::MV
double MV
Definition: MVll.h:742
MVll::h_lambda
gslpp::complex h_lambda(int hel, double q2)
The non-pertubative ccbar contributions to the helicity amplitudes.
Definition: MVll.cpp:1602
MVll::getQCDf_2
gslpp::complex getQCDf_2(double q2)
Definition: MVll.h:523
QCD::initializeMeson
void initializeMeson(QCD::meson meson_i) const
A method to initialize a meson.
Definition: QCD.cpp:236
MVll::getQCDf_1
gslpp::complex getQCDf_1(double q2)
Definition: MVll.h:516
MVll::exp_Phase
gslpp::complex exp_Phase[3]
Definition: MVll.h:736
make_vector
Definition: std_make_vector.h:15
h_0
A class for the correction in .
Definition: MVllObservables.h:1587
MVll::lep
QCD::lepton lep
Definition: MVll.h:727
MVll::mySM
const StandardModel & mySM
Definition: MVll.h:726
MVll::myF_1
std::unique_ptr< F_1 > myF_1
Definition: MVll.h:731
MVll::etaV
int etaV
Definition: MVll.h:755
QCD::K_star
Definition: QCD.h:348
MVll::Mlep
double Mlep
Definition: MVll.h:740
F_2
Definition: F_2.h:15
gslpp::pow
complex pow(const complex &z1, const complex &z2)
Definition: gslpp_complex.cpp:395
MVll::mJ2
double mJ2
Definition: MVll.h:735
gslpp::sqrt
complex sqrt(const complex &z)
Definition: gslpp_complex.cpp:385
gslpp::complex::i
static const complex & i()
Definition: gslpp_complex.cpp:154
MVll::meson
QCD::meson meson
Definition: MVll.h:728
StandardModel::getFlavour
const Flavour & getFlavour() const
Definition: StandardModel.h:1006
F_1
Definition: F_1.h:15
MVll::FixedWCbtos
bool FixedWCbtos
Definition: MVll.h:734
QCD::K_star_P
Definition: QCD.h:349
MVll::width
double width
Definition: MVll.h:751
QCD::PHI
Definition: QCD.h:347
MVll::myF_2
std::unique_ptr< F_2 > myF_2
Definition: MVll.h:732
MVll::angmomV
double angmomV
Definition: MVll.h:754
convertToGslFunction
gsl_function convertToGslFunction(const F &f)
Definition: gslpp_function_adapter.h:24
MVll::integrateDelta
double integrateDelta(int i, double q_min, double q_max)
The integral of from to .
Definition: MVll.cpp:2001
MVll::MM
double MM
Definition: MVll.h:741
MVll::dispersion
bool dispersion
Definition: MVll.h:733
MVll::getQCDf_3
gslpp::complex getQCDf_3(double q2)
Definition: MVll.h:530
Flavour::getFlagUseDispersionRelation
bool getFlagUseDispersionRelation() const
Definition: Flavour.h:242
Flavour::getFlagFixedWCbtos
bool getFlagFixedWCbtos() const
Definition: Flavour.h:252
MVll::mvllParameters
std::vector< std::string > mvllParameters
Definition: MVll.h:730