A class for the \(M \to V 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 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 [Jager:2012uw], 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 [Beneke:2001at] 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 [Straub:2015ica] . 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.
|
gslpp::complex | AmpMVpsi_zExpansion (double mpsi, int tran) |
| Polarization amplitudes for M to V psi, Eq. B.16 of arXiv:2206.03797. More...
|
|
double | beta (double q2) |
| The factor \( \beta \) used in the angular coefficients \(I_i\). More...
|
|
double | Delta_C9_zExp (int hel) |
| The non-pertubative ccbar contributions to the helicity amplitudes. More...
|
|
double | getDelta_C9_zExp_0 () |
| The non-pertubative ccbar contributions to the helicity amplitudes. More...
|
|
double | getDelta_C9_zExp_m () |
| The non-pertubative ccbar contributions to the helicity amplitudes. More...
|
|
double | getDelta_C9_zExp_p () |
| The non-pertubative ccbar contributions to the helicity amplitudes. 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...
|
|
double | integrateSigmaTree (double q_min, double q_max) |
| The integral of \( \Sigma_{tree} \) from \(q_{min}\) to \(q_{max}\) (arxiv/2301.06990) More...
|
|
| MVll (const StandardModel &SM_i, QCD::meson meson_i, QCD::meson vector_i, QCD::lepton lep_i) |
| Constructor. More...
|
|
virtual | ~MVll () |
| Destructor. More...
|
|
|
gslpp::complex | A_Seidel (double q2, double mb2) |
|
gslpp::complex | B0 (double s, double m2) |
|
gslpp::complex | B0diff (double q2, double u, double m2) |
|
gslpp::complex | B_Seidel (double q2, double mb2) |
|
gslpp::complex | C_Seidel (double q2) |
|
gslpp::complex | Cq34 (bool conjugate) |
| QCDF Correction from various BFS paper (hep-ph/0412400). Part of Weak Annihilation. More...
|
|
gslpp::complex | deltaC7_QCDF (double q2, bool conjugate, 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_QCDF (double q2, bool conjugate, bool spline=false) |
| QCDF Correction from various BFS papers (hep-ph/0403185, hep-ph/0412400) and Greub et. al (arXiv:0810.4077).. More...
|
|
double | FF_fit (double q2, double a_0, double a_1, double a_2, double MR2) |
| The fit function from [Straub:2015ica], \( FF^{\rm fit} \). More...
|
|
void | fit_QCDF_func () |
|
double | getintegratedSigmaTree () |
| The integral of \( \Sigma_{tree} \) from 0 to \(q_{cut}\). More...
|
|
gslpp::complex | h_func (double s, double m2) |
|
gslpp::complex | I1 (double q2, double u, double m2) |
|
gslpp::complex | lambda_B_minus (double q2) |
|
double | phi_V (double u) |
| QCDF Correction from various BFS paper (hep-ph/0106067).Vector meson distribution amplitude. More...
|
|
double | QCDF_fit_func (double *x, double *p) |
|
double | SigmaTree (double q2) |
|
void | spline_QCDF_func () |
|
gslpp::complex | T_0 (double q2, bool conjugate) |
|
gslpp::complex | T_minus (double q2, bool conjugate) |
|
gslpp::complex | t_para (double q2, double u, double m2) |
| QCDF Correction from various BFS paper (hep-ph/0106067). Part of 4 quark operator contribution. More...
|
|
double | T_para_imag (double q2, bool conjugate) |
| QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total. More...
|
|
double | T_para_imag (double q2, double u, bool conjugate) |
| QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total. More...
|
|
gslpp::complex | T_para_minus_O8 (double q2, double u) |
| QCDF Correction from various BFS paper (hep-ph/0106067). Chromomagnetic dipole contribution contribution. More...
|
|
gslpp::complex | T_para_minus_QSS (double q2, double u, bool conjugate) |
| QCDF Correction from various BFS paper (hep-ph/0106067). 4 quark operator contribution. More...
|
|
gslpp::complex | T_para_minus_WA (bool conjugate) |
| QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation. More...
|
|
gslpp::complex | T_para_plus_QSS (double q2, double u, bool conjugate) |
| QCDF Correction from various BFS paper (hep-ph/0106067). 4 quark operator contribution. More...
|
|
double | T_para_real (double q2, bool conjugate) |
| QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total. More...
|
|
double | T_para_real (double q2, double u, bool conjugate) |
| QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total. More...
|
|
gslpp::complex | t_perp (double q2, double u, double m2) |
| QCDF Correction from various BFS paper (hep-ph/0106067). Part of 4 quark operator contribution. More...
|
|
double | T_perp_imag (double q2, bool conjugate) |
| QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total. More...
|
|
double | T_perp_imag (double q2, double u, bool conjugate) |
| QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total. More...
|
|
gslpp::complex | T_perp_plus_O8 (double q2, double u) |
| QCDF Correction from various BFS paper (hep-ph/0106067). Chromomagnetic dipole contribution contribution. More...
|
|
gslpp::complex | T_perp_plus_QSS (double q2, double u, bool conjugate) |
| QCDF Correction from various BFS paper (hep-ph/0106067). 4 quark operator contribution. More...
|
|
double | T_perp_real (double q2, bool conjugate) |
| QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total. More...
|
|
double | T_perp_real (double q2, double u, bool conjugate) |
| QCDF Correction from various BFS papers (hep-ph/0106067, hep-ph/0412400). Total. More...
|
|
gslpp::complex | T_perp_WA_1 () |
| QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation. More...
|
|
gslpp::complex | T_perp_WA_2 (bool conjugate) |
| QCDF Correction from various BFS paper (hep-ph/0412400). Weak Annihilation. More...
|
|
@f\aplha_s(\mu_b)$ \( */ double fB; double fpara; double fperp; double MW; /**<W boson mass */ gslpp::complex lambda_t; /**<Vckm factor */ gslpp::complex lambda_u; /**<Vckm factor */ double b; /**<BF of the decay V -> final states */ gslpp::complex h_0[3]; /**<Parameter that contains the contribution from the hadronic hamiltonian */ gslpp::complex h_1[3]; /**<Parameter that contains the contribution from the hadronic hamiltonian */ gslpp::complex h_2[3]; /**<Parameter that contains the contribution from the hadronic hamiltonian */ gslpp::complex SU3_breaking; gslpp::complex beta_0[7]; /**<Parameter that contains the contribution from the hadronic hamiltonian */ gslpp::complex beta_1[7]; /**<Parameter that contains the contribution from the hadronic hamiltonian */ gslpp::complex beta_2[7]; /**<Parameter that contains the contribution from the hadronic hamiltonian */ double t_p;/**< Cache variable */ double t_m;/**< Cache variable */ double t_0;/**< Cache variable */ double z_0;/**< Cache variable */ double s_p; double s_0; double Q2; double chiOPE; double twoalphaBtoKst; double rho_0; double rho_1; double rho_2; double rho_3; double rho_4; double rho_5; double onemrho_0_2; double onemrho_1_2; double onemrho_2_2; double onemrho_3_2; double onemrho_4_2; double onemrho_5_2; double DeltaC9; double DeltaC10; double MMpMV;/**< Cache variable */ double MMpMV2;/**< Cache variable */ double MMmMV;/**< Cache variable */ double MMmMV2;/**< Cache variable */ double rV;/**< Cache variable */ double Chi1minus;/**< Cache variable */ double Chi1plus;/**< Cache variable */ double Chi0plus;/**< Cache variable */ double Chi0minus;/**< Cache variable */ double ChiTT;/**< Cache variable */ double ChiBB;/**< Cache variable */ double MM2;/**< Cache variable */ double MM4;/**< Cache variable */ double MV2;/**< Cache variable */ double MV4;/**< Cache variable */ double MMMV;/**< Cache variable */ double MM2mMV2;/**< Cache variable */ double MM2pMV2;/**< Cache variable */ double fourMV;/**< Cache variable */ double onepMMoMV;/**< Cache variable */ double MM_MMpMV;/**< Cache variable */ double twoMM2;/**< Cache variable */ double twoMV2;/**< Cache variable */ double twoMM_mbpms;/**< Cache variable */ double fourMM2;/**< Cache variable */ double Mlep2;/**< Cache variable */ double twoMlepMb;/**< Cache variable */ double MboMW;/**< Cache variable */ double MsoMb;/**< Cache variable */ double M_PI2osix;/**< Cache variable */ double N_QCDF; double twoMM;/**< Cache variable */ double ninetysixM_PI3MM3;/**< Cache variable */ double M_PI2; double sixteenM_PI2;/**< Cache variable */ double sixteenM_PI2MM2;/**< Cache variable */ double twoMboMM;/**< Cache variable */ gslpp::complex H_0_pre;/**< Cache variable */ double mu_b2;/**< Cache variable */ double Mc2;/**< Cache variable */ double Mb2;/**< Cache variable */ double fourMc2;/**< Cache variable */ double fourMb2;/**< Cache variable */ double logMc;/**< Cache variable */ double logMb;/**< Cache variable */ gslpp::complex H_0_WC;/**< Cache variable */ gslpp::complex H_c_WC;/**< Cache variable */ gslpp::complex H_b_WC;/**< Cache variable */ double fournineth;/**< Cache variable */ double half;/**< Cache variable */ double twothird;/**< Cache variable */ double sqrt3;/**< Cache variable */ gslpp::complex ihalfMPI;/**< Cache variable */ double twoMM3;/**< Cache variable */ double gtilde_1_pre;/**< Cache variable */ double gtilde_2_pre;/**< Cache variable */ double gtilde_3_pre;/**< Cache variable */ double C2_inv;/**< Cache variable */ double S_L_pre;/**< Cache variable */ gslpp::complex NN;/**< coupling including the CKM element */ gslpp::complex NN_conjugate;/**< conjugate of the coupling including the CKM element */ double CF;/**< Cache variable */ double deltaT_0;/**< Cache variable */ double deltaT_1par;/**< Cache variable */ double deltaT_1perp;/**< Cache variable */ bool h_pole; 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 \( */ 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 a_0V;/**<LCSR fit parameter */ double a_1V;/**<LCSR fit parameter */ double a_2V;/**<LCSR fit parameter */ double MRV_2;/**<LCSR fit parameter */ double a_0A0;/**<LCSR fit parameter */ double a_1A0;/**<LCSR fit parameter */ double a_2A0;/**<LCSR fit parameter */ double MRA0_2;/**<LCSR fit parameter */ double a_0A1;/**<LCSR fit parameter */ double a_1A1;/**<LCSR fit parameter */ double a_2A1;/**<LCSR fit parameter */ double MRA1_2;/**<LCSR fit parameter */ double a_0A12;/**<LCSR fit parameter */ double a_1A12;/**<LCSR fit parameter */ double a_2A12;/**<LCSR fit parameter */ double MRA12_2;/**<LCSR fit parameter */ double a_0T1;/**<LCSR fit parameter */ double a_1T1;/**<LCSR fit parameter */ double a_2T1;/**<LCSR fit parameter */ double MRT1_2;/**<LCSR fit parameter */ double a_0T2;/**<LCSR fit parameter */ double a_1T2;/**<LCSR fit parameter */ double a_2T2;/**<LCSR fit parameter */ double MRT2_2;/**<LCSR fit parameter */ double a_0T23;/**<LCSR fit parameter */ double a_1T23;/**<LCSR fit parameter */ double a_2T23;/**<LCSR fit parameter */ double MRT23_2;/**<LCSR fit parameter */ double a_0f;/**<DM fit parameter */ double a_1f;/**<DM fit parameter */ double a_2f;/**<DM fit parameter */ double MRf_2;/**<DM fit parameter */ double a_0g;/**<DM fit parameter */ double a_1g;/**<DM fit parameter */ double a_2g;/**<DM fit parameter */ double MRg_2;/**<DM fit parameter */ double a_0F1;/**<DM fit parameter */ double a_1F1;/**<DM fit parameter */ double a_2F1;/**<DM fit parameter */ double MRF1_2;/**<DM fit parameter */ double a_0F2;/**<DM fit parameter */ double a_1F2;/**<DM fit parameter */ double a_2F2;/**<DM fit parameter */ double MRF2_2;/**<DM fit parameter */ double a_0T0;/**<DM fit parameter */ double a_1T0;/**<DM fit parameter */ double a_2T0;/**<DM fit parameter */ double MRT0_2;/**<DM fit parameter */ //additional variables for B to K nu nu double GF4; double MM3; double fM2; double fV2; double mtau; double mtau2; double Gammatau; double VusVub_abs2; 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 contribution.*/ 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}' \(*/ std::vector<double> Re_T_perp;/**<Vector that samples the QCDF \)@_fakenlRe(T_{perp}) \( */ std::vector<double> Im_T_perp;/**<Vector that samples the QCDF \)@_fakenlIm(T_{perp}) \( */ std::vector<double> Re_T_para;/**<Vector that samples the QCDF \)@_fakenlRe(T_{para}) \( */ std::vector<double> Im_T_para;/**<Vector that samples the QCDF \)@_fakenlIm(T_{para}) \( */ gsl_interp_accel *acc_Re_T_perp; gsl_interp_accel *acc_Im_T_perp; gsl_interp_accel *acc_Re_T_para; gsl_interp_accel *acc_Im_T_para; gsl_spline *spline_Re_T_perp; gsl_spline *spline_Im_T_perp; gsl_spline *spline_Re_T_para; gsl_spline *spline_Im_T_para; 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; #if COMPUTECP std::vector<double> Re_T_perp_conj;/**<Vector that samples the QCDF \)@_fakenlRe(T_{perp}) \( */ std::vector<double> Im_T_perp_conj;/**<Vector that samples the QCDF \)@_fakenlIm(T_{perp}) \( */ std::vector<double> Re_T_para_conj;/**<Vector that samples the QCDF \)@_fakenlRe(T_{para}) \( */ std::vector<double> Im_T_para_conj;/**<Vector that samples the QCDF \)@_fakenlIm(T_{para}) \( */ gsl_interp_accel *acc_Re_T_perp_conj; gsl_interp_accel *acc_Im_T_perp_conj; gsl_interp_accel *acc_Re_T_para_conj; gsl_interp_accel *acc_Im_T_para_conj; gsl_interp_accel *acc_Re_deltaC7_QCDF_conj; gsl_interp_accel *acc_Im_deltaC7_QCDF_conj; gsl_interp_accel *acc_Re_deltaC9_QCDF_conj; gsl_interp_accel *acc_Im_deltaC9_QCDF_conj; gsl_spline *spline_Re_T_perp_conj; gsl_spline *spline_Im_T_perp_conj; gsl_spline *spline_Re_T_para_conj; gsl_spline *spline_Im_T_para_conj; gsl_spline *spline_Re_deltaC7_QCDF_conj; gsl_spline *spline_Im_deltaC7_QCDF_conj; gsl_spline *spline_Re_deltaC9_QCDF_conj; gsl_spline *spline_Im_deltaC9_QCDF_conj; #endif std::vector<double> myq2;/**<Variable used to compute the QCDF \)\Delta C_9 \( */ TFitResultPtr Re_T_perp_res;/**<Vector that contains the fitted QCDF \)@_fakenlRe(T_{perp}) \( */ TFitResultPtr Im_T_perp_res;/**<Vector that contains the fitted QCDF \)@_fakenlIm(T_{perp}) \( */ TFitResultPtr Re_T_para_res;/**<Vector that contains the fitted QCDF \)@_fakenlRe(T_{para}) \( */ TFitResultPtr Im_T_para_res;/**<Vector that contains the fitted QCDF \)@_fakenlIm(T_{para}) \( */ TFitResultPtr Re_T_perp_res_conj;/**<Vector that contains the fitted QCDF \)@_fakenlRe(T_{perp}) \( */ TFitResultPtr Im_T_perp_res_conj;/**<Vector that contains the fitted QCDF \)@_fakenlIm(T_{perp}) \( */ TFitResultPtr Re_T_para_res_conj;/**<Vector that contains the fitted QCDF \)@_fakenlRe(T_{para}) \( */ TFitResultPtr Im_T_para_res_conj;/**<Vector that contains the fitted QCDF \)@_fakenlIm(T_{para}) \( */ 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 QCDFfit;/**<TF1 to be used for fitting the QCDF. */ 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 avaSigma;/**< GSL integral variable */ double avaDelta;/**< GSL integral variable */ double avaSigmaTree;/**< Gsl integral variable */ double errSigma;/**< GSL integral variable */ double errDelta;/**< GSL integral variable */ double errSigmaTree;/**< Gsl integral variable */ gsl_function FS;/**< GSL integral variable */ gsl_function FD;/**< 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_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;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma1;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma2;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma3;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma4;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma5;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma6;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma7;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma9;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma10;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigma11;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheDelta0;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheDelta1;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheDelta2;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheDelta3;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheDelta6;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheDelta7;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheDelta8;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheDelta10;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheDelta11;/**< Cache variable */ std::map<std::pair<double, double>, double > cacheSigmaTree;/**< Gsl integral variable */ unsigned int N_updated;/**< Cache variable */ gslpp::vector<double> N_cache;/**< Cache variable */ gslpp::complex Nc_cache;/**< Cache variable */ unsigned int V_updated;/**< Cache variable */ gslpp::vector<double> V_cache;/**< Cache variable */ unsigned int A0_updated;/**< Cache variable */ gslpp::vector<double> A0_cache;/**< Cache variable */ unsigned int A1_updated;/**< Cache variable */ gslpp::vector<double> A1_cache;/**< Cache variable */ unsigned int T1_updated;/**< Cache variable */ gslpp::vector<double> T1_cache;/**< Cache variable */ unsigned int T2_updated;/**< Cache variable */ gslpp::vector<double> T2_cache;/**< Cache variable */ unsigned int k2_updated;/**< Cache variable */ gslpp::vector<double> k2_cache;/**< Cache variable */ unsigned int z_updated;/**< Cache variable */ unsigned int lambda_updated;/**< Cache variable */ unsigned int beta_updated;/**< Cache variable */ double beta_cache;/**< Cache variable */ unsigned int F_updated;/**< Cache variable */ unsigned int VL1_updated;/**< Cache variable */ unsigned int VL2_updated;/**< Cache variable */ unsigned int TL1_updated;/**< Cache variable */ unsigned int TL2_updated;/**< Cache variable */ unsigned int VR1_updated;/**< Cache variable */ unsigned int VR2_updated;/**< Cache variable */ unsigned int TR1_updated;/**< Cache variable */ unsigned int TR2_updated;/**< Cache variable */ unsigned int VL0_updated;/**< Cache variable */ gslpp::vector<double> VL0_cache;/**< Cache variable */ unsigned int TL0_updated;/**< Cache variable */ gslpp::vector<double> TL0_cache;/**< Cache variable */ unsigned int VR0_updated;/**< Cache variable */ unsigned int TR0_updated;/**< Cache variable */ unsigned int Mb_Ms_updated;/**< Cache variable */ unsigned int SL_updated;/**< Cache variable */ gslpp::vector<double> SL_cache;/**< Cache variable */ unsigned int SR_updated;/**< 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 */ gslpp::complex h0Ccache[4];/**< Cache variable */ gslpp::complex h1Ccache[4];/**< Cache variable */ gslpp::complex h2Ccache[4];/**< Cache variable */ gslpp::complex beta0Ccache[8];/**< Cache variable */ gslpp::complex beta1Ccache[8];/**< Cache variable */ gslpp::complex beta2Ccache[8];/**< Cache variable */ unsigned int h0_updated;/**< Cache variable */ unsigned int h1_updated;/**< Cache variable */ unsigned int h2_updated;/**< Cache variable */ unsigned int H_V0updated;/**< Cache variable */ gslpp::vector<double> H_V0cache;/**< Cache variable */ unsigned int H_V1updated;/**< Cache variable */ gslpp::vector<double> H_V1cache;/**< Cache variable */ unsigned int H_V2updated;/**< Cache variable */ gslpp::vector<double> H_V2cache;/**< Cache variable */ unsigned int H_A0updated;/**< Cache variable */ unsigned int H_A1updated;/**< Cache variable */ unsigned int H_A2updated;/**< Cache variable */ unsigned int H_Supdated;/**< Cache variable */ gslpp::vector<double> H_Scache;/**< Cache variable */ unsigned int H_Pupdated;/**< Cache variable */ gslpp::vector<double> H_Pcache;/**< Cache variable */ unsigned int I0_updated;/**< Cache variable */ unsigned int I1_updated;/**< Cache variable */ unsigned int I2_updated;/**< Cache variable */ unsigned int I3_updated;/**< Cache variable */ unsigned int I4_updated;/**< Cache variable */ unsigned int I5_updated;/**< Cache variable */ unsigned int I6_updated;/**< Cache variable */ unsigned int I7_updated;/**< Cache variable */ unsigned int I8_updated;/**< Cache variable */ unsigned int I9_updated;/**< Cache variable */ unsigned int I10_updated;/**< Cache variable */ unsigned int I11_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 > sigma1Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma2Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma3Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma4Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma5Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma6Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma7Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma9Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma10Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigma11Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta0Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta1Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta2Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta3Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta6Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta7Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta8Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta10Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > delta11Cached;/**< Cache variable */ std::map<std::pair<double, double>, unsigned int > sigmaTreeCached;/**< Cache variable */ std::map<double, unsigned int> deltaTparpCached;/**< Cache variable */ std::map<double, unsigned int> deltaTparmCached;/**< Cache variable */ std::map<double, unsigned int> deltaTperpCached;/**< Cache variable */ std::map<double, gslpp::complex> cacheDeltaTparp;/**< Cache variable */ std::map<double, gslpp::complex> cacheDeltaTparm;/**< Cache variable */ std::map<double, gslpp::complex> cacheDeltaTperp;/**< Cache variable */ unsigned int deltaTparpupdated;/**< Cache variable */ unsigned int deltaTparmupdated;/**< Cache variable */ unsigned int deltaTperpupdated;/**< Cache variable */ unsigned int T_updated;/**< Cache variable */ gslpp::vector<double> T_cache;/**< Cache variable */ /** @brief The update parameter method for MVll. */ void updateParameters(); /** @brief The caching method for MVll. */ void checkCache(); /** @brief The lattice parameter \) z \( from arXiv:1310.3722v3. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) z \( */ double z(double q2); /** @brief The DM parameter \) z \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) z \( */ double z_DM(double q2); /** @brief The prefactor function of the form factor \) f \(, \) \phi_f \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param MRf_2 fit parameter @return \) \phi_f \( */ double phi_f(double q2, double MRf_2); /** @brief The prefactor function of the form factor \) g \(, \) \phi_g \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param MRg_2 fit parameter @return \) \phi_g \( */ double phi_g(double q2, double MRg_2); /** @brief The prefactor function of the form factor \) F_1 \(, \) \phi_{F_1} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param MRF1_2 fit parameter @return \) \phi_{F_1} \( */ double phi_F1(double q2, double MRF1_2); /** @brief The prefactor function of the form factor \) F_2 \(, \) \phi_{F_2} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param MRF2_2 fit parameter @return \) \phi_{F_2} \( */ double phi_F2(double q2, double MRF2_2); /** @brief The prefactor function of the form factor \) T_0 \(, \) \phi_{T_0} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param MRT0_2 fit parameter @return \) \phi_{T_0} \( */ double phi_T0(double q2, double MRT0_2); /** @brief The prefactor function of the form factor \) T_1 \(, \) \phi_{T_1} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param MRT1_2 fit parameter @return \) \phi_{T_1} \( */ double phi_T1(double q2, double MRT1_2); /** @brief The prefactor function of the form factor \) T_2 \(, \) \phi_{T_2} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param MRT2_2 fit parameter @return \) \phi_{T_2} \( */ double phi_T2(double q2, double MRT2_2); /** @brief The transverse form factor \) f \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] a_0f fit parameter @param[in] a_1f fit parameter @param[in] a_2f fit parameter @param[in] MRf_2 fit parameter @return \) f \( */ double f_DM(double q2, double a_0f, double a_1f, double a_2f, double MRf_2); /** @brief The transverse form factor \) g \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] a_0g fit parameter @param[in] a_1g fit parameter @param[in] a_2g fit parameter @param[in] MRg_2 fit parameter @return \) g \( */ double g_DM(double q2, double a_0g, double a_1g, double a_2g, double MRg_2); /** @brief The transverse form factor \) F_1 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] a_0F1 fit parameter @param[in] a_1F1 fit parameter @param[in] a_2F1 fit parameter @param[in] MRF1_2 fit parameter @return \) F_1 \( */ double F1_DM(double q2, double a_0F1, double a_1F1, double a_2F1, double MRF1_2); /** @brief The transverse form factor \) F_2 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] a_0F2 fit parameter @param[in] a_1F2 fit parameter @param[in] a_2F2 fit parameter @param[in] MRF2_2 fit parameter @return \) F_2 \( */ double F2_DM(double q2, double a_0F2, double a_1F2, double a_2F2, double MRF2_2); /** @brief The transverse form factor \) T_0 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] a_0T0 fit parameter @param[in] a_1T0 fit parameter @param[in] a_2T0 fit parameter @param[in] MRT0_2 fit parameter @return \) T_0 \( */ double T0_DM(double q2, double a_0T0, double a_1T0, double a_2T0, double MRT0_2); /** @brief The transverse form factor \) T_1 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] a_0T1 fit parameter @param[in] a_1T1 fit parameter @param[in] a_2T1 fit parameter @param[in] MRT1_2 fit parameter @return \) T_1 \( */ double T1_DM(double q2, double a_0T1, double a_1T1, double a_2T1, double MRT1_2); /** @brief The transverse form factor \) T_2 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] a_0T2 fit parameter @param[in] a_1T2 fit parameter @param[in] a_2T2 fit parameter @param[in] MRT2_2 fit parameter @return \) T_2 \( */ double T2_DM(double q2, double a_0T2, double a_1T2, double a_2T2, double MRT2_2); /** @brief The transverse form factor \) V \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) V \( */ double V(double q2); /** @brief The transverse form factor \) A_0 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) A_0 \( */ double A_0(double q2); /** @brief The transverse form factor \) A_1 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) A_1 \( */ double A_1(double q2); /** @brief The transverse form factor \) A_2 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) A_2 \( */ double A_2(double q2); /** @brief The transverse form factor \) T_1 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) T_1 \( */ double T_1(double q2); /** @brief The transverse form factor \) T_2 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) T_2 \( */ double T_2(double q2); /** @brief The helicity form factor \) \tilde{V}_0 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \tilde{V}_0 \( */ double V_0t(double q2); /** @brief The helicity form factor \) V_+ \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) V_+ \( */ double V_p(double q2); /** @brief The helicity form factor \) V_- \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) V_- \( */ double V_m(double q2); /** @brief The helicity form factor \) \tilde{T}_0 \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \tilde{T}_0 \( */ double T_0t(double q2); /** @brief The helicity form factor \) T_+ \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) T_+ \( */ double T_p(double q2); /** @brief The helicity form factor \) T_- \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) T_- \( */ double T_m(double q2); /** @brief The helicity form factor \) S_L \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) S_L \( */ double S_L(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^2) \( function involved into \) C_9^{eff} \(. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] m2 squared mass @param[in] mu2 squared mass scale @return \) h(q^2,m^2) \( */ gslpp::complex H(double q2, double m2, 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); gslpp::complex funct_g(double q2); gslpp::complex DeltaC9_KD(double q2, int com); /** @brief The expansion parameter \) \hat{z} \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) zh \( */ gslpp::complex zh(double q2); /** @brief The Blaschke factor \) \mathcal{P} \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) P \( */ gslpp::complex P(double q2); /** @brief The recursive function \) Phi_1 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_1 \( */ gslpp::complex Phi_1(double q2); /** @brief The recursive function \) Phi_1^* \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_1^* \( */ gslpp::complex Phi_1_st(double q2); /** @brief The recursive function \) Phi_2 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_2 \( */ gslpp::complex Phi_2(double q2); /** @brief The recursive function \) Phi_2^* \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_2^* \( */ gslpp::complex Phi_2_st(double q2); /** @brief The recursive function \) Phi_3 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_3 \( */ gslpp::complex Phi_3(double q2); /** @brief The recursive function \) Phi_3^* \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_3^* \( */ gslpp::complex Phi_3_st(double q2); /** @brief The recursive function \) Phi_4 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_4 \( */ gslpp::complex Phi_4(double q2); /** @brief The recursive function \) Phi_4^* \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_4^* \( */ gslpp::complex Phi_4_st(double q2); /** @brief The recursive function \) Phi_5 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_5 \( */ gslpp::complex Phi_5(double q2); /** @brief The recursive function \) Phi_5^* \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_5^* \( */ gslpp::complex Phi_5_st(double q2); /** @brief The recursive function \) Phi_6 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_6 \( */ gslpp::complex Phi_6(double q2); /** @brief The recursive function \) Phi_6^* \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) Phi_6^* \( */ gslpp::complex Phi_6_st(double q2); /** @brief The 0th normalized Szego polynomial \) p_0 \( from arXiv:2206.03797. @return \) p0 \( */ gslpp::complex p0(); /** @brief The 1st normalized Szego polynomial \) p_1 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) p1 \( */ gslpp::complex p1(double q2); /** @brief The 2nd normalized Szego polynomial \) p_2 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) p2 \( */ gslpp::complex p2(double q2); /** @brief The 3rd normalized Szego polynomial \) p_3 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) p3 \( */ gslpp::complex p3(double q2); /** @brief The 4th normalized Szego polynomial \) p_4 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) p4 \( */ gslpp::complex p4(double q2); /** @brief The 5th normalized Szego polynomial \) p_5 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) p5 \( */ gslpp::complex p5(double q2); /** @brief The 6th normalized Szego polynomial \) p_6 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) p6 \( */ gslpp::complex p6(double q2); /** @brief The outer function \) phi_1 \( from arXiv:2011.09813. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) phi_1 \( */ gslpp::complex phi_1(double q2); /** @brief The outer function \) phi_2 \( from arXiv:2011.09813. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) phi_2 \( */ gslpp::complex phi_2(double q2); /** @brief The outer function \) phi_3 \( from arXiv:2011.09813. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) phi_3 \( */ gslpp::complex phi_3(double q2); /** @brief The outer function \) phi_4 \( from arXiv:2011.09813. @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) phi_4 \( */ gslpp::complex phi_4(double q2); /** @brief The correction to \) C_9 \( from arXiv:2206.03797. @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] tran transversity @return \) \Delta C_9 \( */ gslpp::complex DeltaC9_zExpansion(double q2, int tran); /** @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^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 @param[in] b_i BF of the decay \) V \to M_1 M_2 \( @return \) F \( */ double F(double q2, double b_i); /** @brief The angular coefficient \) I_{1c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_{1c} \( */ double I_1c(double q2, bool bar); /** @brief The angular coefficient \) I_{1s} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_{1s} \( */ double I_1s(double q2, bool bar); /** @brief The angular coefficient \) I_{2c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_{2c} \( */ double I_2c(double q2, bool bar); /** @brief The angular coefficient \) I_{2s} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_{2s} \( */ double I_2s(double q2, bool bar); /** @brief The angular coefficient \) I_3 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_3 \( */ double I_3(double q2, bool bar); /** @brief The angular coefficient \) I_4 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_4 \( */ double I_4(double q2, bool bar); /** @brief The angular coefficient \) I_5 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_5 \( */ double I_5(double q2, bool bar); /** @brief The angular coefficient \) I_{6c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_{6c} \( */ double I_6c(double q2, bool bar); /** @brief The angular coefficient \) I_{6s} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_{6s} \( */ double I_6s(double q2, bool bar); /** @brief The angular coefficient \) I_7 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_7 \( */ double I_7(double q2, bool bar); /** @brief The angular coefficient \) I_8 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_8 \( */ double I_8(double q2, bool bar); /** @brief The angular coefficient \) I_9 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) I_9 \( */ double I_9(double q2, bool bar); /** @brief The angular coefficient \) h_1s \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) h_1s \( */ double h_1s(double q2, bool bar); /** @brief The angular coefficient \) h_1c \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) h_1c \( */ double h_1c(double q2, bool bar); /** @brief The angular coefficient \) h_2s \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) h_2s \( */ double h_2s(double q2, bool bar); /** @brief The angular coefficient \) h_2c \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) h_2c \( */ double h_2c(double q2, bool bar); /** @brief The angular coefficient \) h_3 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) h_3 \( */ double h_3(double q2, bool bar); /** @brief The angular coefficient \) h_4 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) h_4 \( */ double h_4(double q2, bool bar); /** @brief The angular coefficient \) h_7 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) h_7 \( */ double h_7(double q2, bool bar); /** @brief The angular coefficient \) s_5 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) s_5 \( */ double s_5(double q2, bool bar); /** @brief The angular coefficient \) s_6s \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) s_6s \( */ double s_6s(double q2, bool bar); /** @brief The angular coefficient \) s_6c \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) s_6c \( */ double s_6c(double q2, bool bar); /** @brief The angular coefficient \) s_8 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) s_8 \( */ double s_8(double q2, bool bar); /** @brief The angular coefficient \) s_9 \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @param[in] bar boolean variable to distinguish the decay \) M \to V \ell \ell \( (true) from the CP-conjugate \) \bar{M} \to \bar{V} \ell \ell \( (false) @return \) s_9 \( */ double s_9(double q2, bool bar); /** @brief The CP average \) \Sigma_{1s} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{1s} \( */ double getSigma1c(double q2) { switch(vectorM){ case QCD::K_star: return (I_1c(q2, 0) + I_1c(q2, 1))/2.; break; case QCD::K_star_P: return (I_1c(q2, 0) + I_1c(q2, 1))/2.; break; case QCD::PHI: return (I_1c(q2, 0) + I_1c(q2, 1) - ys * h_1c(q2, 0) )/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getSigma1c : vector " + out.str() + " not implemented"); } }; /** @brief The CP average \) \Sigma_{1c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{1c} \( */ double getSigma1s(double q2) { switch(vectorM){ case QCD::K_star: return (I_1s(q2, 0) + I_1s(q2, 1))/2.; break; case QCD::K_star_P: return (I_1s(q2, 0) + I_1s(q2, 1))/2.; break; case QCD::PHI: return (I_1s(q2, 0) + I_1s(q2, 1) - ys * h_1s(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getSigma1s : vector " + out.str() + " not implemented"); } }; /** @brief The CP average \) \Sigma_{2s} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{2s} \( */ double getSigma2c(double q2) { switch(vectorM){ case QCD::K_star: return (I_2c(q2, 0) + I_2c(q2, 1))/2.; break; case QCD::K_star_P: return (I_2c(q2, 0) + I_2c(q2, 1))/2.; break; case QCD::PHI: return (I_2c(q2, 0) + I_2c(q2, 1) - ys * h_2c(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getSigma2c : vector " + out.str() + " not implemented"); } }; /** @brief The CP average \) \Sigma_{2c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{2c} \( */ double getSigma2s(double q2) { switch(vectorM){ case QCD::K_star: return (I_2s(q2, 0) + I_2s(q2, 1))/2.; break; case QCD::K_star_P: return (I_2s(q2, 0) + I_2s(q2, 1))/2.; break; case QCD::PHI: return (I_2s(q2, 0) + I_2s(q2, 1) - ys * h_2s(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getSigma2s : vector " + out.str() + " not implemented"); } }; /** @brief The CP average \) \Sigma_{3} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{3} \( */ double getSigma3(double q2) { switch(vectorM){ case QCD::K_star: return (I_3(q2, 0) + I_3(q2, 1))/2.; break; case QCD::K_star_P: return (I_3(q2, 0) + I_3(q2, 1))/2.; break; case QCD::PHI: return (I_3(q2, 0) + I_3(q2, 1) - ys * h_3(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getSigma3 : vector " + out.str() + " not implemented"); } }; /** @brief The CP average \) \Sigma_{4} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{4} \( */ double getSigma4(double q2) { switch(vectorM){ case QCD::K_star: return (I_4(q2, 0) + I_4(q2, 1))/2.; break; case QCD::K_star_P: return (I_4(q2, 0) + I_4(q2, 1))/2.; break; case QCD::PHI: return (I_4(q2, 0) + I_4(q2, 1) - ys * h_4(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getSigma4 : vector " + out.str() + " not implemented"); } }; /** @brief The CP average \) \Sigma_{5} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{5} \( */ double getSigma5(double q2) { return (I_5(q2, 0) + I_5(q2, 1))/2.; }; /** @brief The CP average \) \Sigma_{6s} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{6s} \( */ double getSigma6s(double q2) { return (I_6s(q2, 0) + I_6s(q2, 1))/2.; }; /** @brief The CP average \) \Sigma_{6c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{6c} \( */ double getSigma6c(double q2) { return (I_6c(q2, 0) + I_6c(q2, 1))/2.; }; /** @brief The CP average \) \Sigma_{7} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{7} \( */ double getSigma7(double q2) { switch(vectorM){ case QCD::K_star: return (I_7(q2, 0) + I_7(q2, 1))/2.; break; case QCD::K_star_P: return (I_7(q2, 0) + I_7(q2, 1))/2.; break; case QCD::PHI: return (I_7(q2, 0) + I_7(q2, 1) - ys * h_7(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getSigma7 : vector " + out.str() + " not implemented"); } }; /** @brief The CP average \) \Sigma_{8} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{8} \( */ double getSigma8(double q2) { return (I_8(q2, 0) + I_8(q2, 1))/2.; }; /** @brief The CP average \) \Sigma_{9} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Sigma_{9} \( */ double getSigma9(double q2) { return (I_9(q2, 0) + I_9(q2, 1))/2.; }; /** @brief The CP difference \) \Delta_{1s} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{1s} \( */ double getDelta1c(double q2) { return (I_1c(q2, 0) - I_1c(q2, 1)) / 2.; }; /** @brief The CP difference \) \Delta_{1c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{1c} \( */ double getDelta1s(double q2) { return (I_1s(q2, 0) - I_1s(q2, 1)) / 2.; }; /** @brief The CP difference \) \Delta_{2s} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{2s} \( */ double getDelta2c(double q2) { return (I_2c(q2, 0) - I_2c(q2, 1)) / 2.; }; /** @brief The CP difference \) \Delta_{2c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{2c} \( */ double getDelta2s(double q2) { return (I_2s(q2, 0) - I_2s(q2, 1))/2.; }; /** @brief The CP difference \) \Delta_{3} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{3} \( */ double getDelta3(double q2) { return (I_3(q2, 0) - I_3(q2, 1))/2.; }; /** @brief The CP difference \) \Delta_{4} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{4} \( */ double getDelta4(double q2) { return (I_4(q2, 0) - I_4(q2, 1))/2.; }; /** @brief The CP difference \) \Delta_{5} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{5} \( */ double getDelta5(double q2) { switch(vectorM){ case QCD::K_star: return (I_5(q2, 0) - I_5(q2, 1))/2.; break; case QCD::K_star_P: return (I_5(q2, 0) - I_5(q2, 1))/2.; break; case QCD::PHI: return (1. - ys*ys)/(1. + xs*xs) * (I_5(q2, 0) - I_5(q2, 1) - xs * s_5(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getDelta5 : vector " + out.str() + " not implemented"); } }; /** @brief The CP difference \) \Delta_{6s} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{6s} \( */ double getDelta6s(double q2) { switch(vectorM){ case QCD::K_star: return (I_6s(q2, 0) - I_6s(q2, 1))/2.; break; case QCD::K_star_P: return (I_6s(q2, 0) - I_6s(q2, 1))/2.; break; case QCD::PHI: return (1. - ys*ys)/(1. + xs*xs) * (I_6s(q2, 0) - I_6s(q2, 1) - xs * s_6s(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getDelta6s : vector " + out.str() + " not implemented"); } }; /** @brief The CP difference \) \Delta_{6c} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{6c} \( */ double getDelta6c(double q2) { switch(vectorM){ case QCD::K_star: return (I_6c(q2, 0) - I_6c(q2, 1))/2.; break; case QCD::K_star_P: return (I_6c(q2, 0) - I_6c(q2, 1))/2.; break; case QCD::PHI: return (1. - ys*ys)/(1. + xs*xs) * (I_6c(q2, 0) - I_6c(q2, 1) - xs * s_6c(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getDelta6c : vector " + out.str() + " not implemented"); } }; /** @brief The CP difference \) \Delta_{7} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{7} \( */ double getDelta7(double q2) { return (I_7(q2, 0) - I_7(q2, 1))/2.; }; /** @brief The CP difference \) \Delta_{8} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{8} \( */ double getDelta8(double q2) { switch(vectorM){ case QCD::K_star: return (I_8(q2, 0) - I_8(q2, 1))/2.; break; case QCD::K_star_P: return (I_8(q2, 0) - I_8(q2, 1))/2.; break; case QCD::PHI: return (1. - ys*ys)/(1. + xs*xs) * (I_8(q2, 0) - I_8(q2, 1) - xs * s_8(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getDelta8 : vector " + out.str() + " not implemented"); } }; /** @brief The CP difference \) \Delta_{9} \( . @param[in] q2 \)@_fakenlq^2 \( of the decay @return \) \Delta_{9} \( */ double getDelta9(double q2) { switch(vectorM){ case QCD::K_star: return (I_9(q2, 0) - I_9(q2, 1))/2.; break; case QCD::K_star_P: return (I_9(q2, 0) - I_9(q2, 1))/2.; break; case QCD::PHI: return (1. - ys*ys)/(1. + xs*xs) * (I_9(q2, 0) - I_9(q2, 1) - xs * s_9(q2, 0))/2.; break; default: std::stringstream out; out << vectorM; throw std::runtime_error("MVll::getDelta9 : vector " + out.str() + " not implemented"); } }; /** @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 814 of file MVll.h.