Form Factors (hazma.form_factors)

Overview

hazma has several form factors for computing matrix elements into mesonic final states. Currently, hazma has all the important vector form factors relevant for center of mass energies below the \(\tau\) mass. In a future release, we may include scalar form factors or form factors for other tensor interactions.

Below, we describe how to use the form factors to compute cross-sections, widths and distributions.

Vector Form Factors

Each of the form factor classes listed in the class table follow a particular structure. They all have the following functions:

Method

Description

form_factor(q,..., **kwargs)

Compute the form factor (without Lorentz structure).

integrated_form_factor(q, **kwargs)

Compute the squared current integrated over phase space.

width(mv, **kwargs)

Compute the partial decay width of a massive vector.

cross_section(q, mx, mv, gvxx, wv, **kwargs)

Compute the annihilation cross-section of dark matter.

q represents the center-of-mass energy in these functions. The ... is a place-holder for additional arguments that depend on the type of form factor.

  • Two-body final states: For two-body final states, q is the only needed argument (aside from additional arguments that depend on the specific final state.)

  • Three-body final states: For three-body final states, the squared invariant masses \(s=(p_2+p_3)^2\) and \(t=(p_1+p_3)^2\) are required, where \(p_1\), \(p_2\) and \(p_3\) are the momenta of particles 1, 2, and 3.

  • N>3-body final states: Currently, the maximum number of final states is 4. However, for \(N\)-body final states with \(N>3\), the form factors require the full momentum information of the final states. In these cases, the signature is form_factor(q, momenta, **kwargs), where momenta is a NumPy array containing the 4-momenta of all final state particles.

The **kwargs contains all the additional information needed for the specific final state (e.g. couplings.) All of these functions are vectorized, so you can pass an array of values to compute the quantity for each value at once.

The form factors are given as follows. Let \(J^{\mu} = \bra{\mathcal{H}}j^{\mu}\ket{0}\) be the hadronic current for the given final state \(\mathcal{H}\). We decompose the current into the form

\[J^{\mu}_{\mathcal{H}} = \sum_{a=1}^{N}F^{a}_{\mathcal{H}}j^{a,\mu}_{\mathcal{H}}\]

where \(F^{a}_{\mathcal{H}}\) is the form factor associated with current \(j^{a,\mu}_{\mathcal{H}}\). For all the two- and three-body final states, there is only one current (i.e. \(N=1\).) For these cases, the hadronic currents for final states containing two pseudo-scalars, a pseudo-scalar and photon, a pseudo-scalar and vector and three pseudo-scalars are given by:

\[\begin{split}J^{\mu}_{P_1P_2} &= -(p_{1}^{\mu}-p_{2}^{\mu})F_{P_1P_2}(q^2), & q &=p_1+p_2\\ J^{\mu}_{P\gamma} &= \epsilon_{\mu\nu\alpha\beta}q^{\nu}\epsilon^{\alpha}_{\gamma}(p_{\gamma})p^{\beta}_{\gamma}F_{P\gamma}(q^2), & q &=p_{P}+p_{\gamma}\\ J^{\mu}_{PV} &= \epsilon_{\mu\nu\alpha\beta}q^{\nu}\epsilon^{\alpha}_{V}(p_{V})p^{\beta}_{P}F_{PV}(q^2), & q &=p_{P}+p_{V}\\ J^{\mu}_{P_1P_2P_3} &= \epsilon_{\mu\nu\alpha\beta}p_{1}^{\nu}p_{2}^{\alpha}p_{3}^{\beta}F_{P_1P_2P_3}(s,t)\end{split}\]

In the two-body form factors, the momentum \(q\) is given by the sum of the final state momenta. In the three pseudo-scalar form factor, \(s=(p_2+p_3)^2\) and \(t=(p_1+p_3)^2\).

In terms of these currents, the form_factor functions compute the \(F_{\mathcal{H}}\). The integrated_form_factor functions compute the following:

\[\mathcal{J}_{\mathcal{H}}(q^2) = -\frac{1}{3q^2}g_{\mu\nu}\int\dd{\Pi}_{\mathrm{LIPS}} J_{\mathcal{H}}^{\mu}\bar{J}_{\mathcal{H}}^{\nu}\]

From the integrated current, the widths and cross-sections are easily computed. The expressions are:

\[\begin{split}\sigma_{\bar{\chi}\chi\to\mathcal{H}}(q^2) &= \frac{g_{V\chi\chi}^2(q^2+2m_{\chi}^2)} {\sqrt{q^2-4m_{\chi}^2}\qty((q^2-m_{V})^2+m_{V}^2\Gamma_{V}^2)} \frac{\sqrt{q^2}}{2}\mathcal{J}_{\mathcal{H}}(q^2)\\ \Gamma_{V\to\mathcal{H}} &= \frac{m_{V}}{2}\mathcal{J}_{\mathcal{H}}(m_{V}^2)\end{split}\]

Examples

Partial Widths of Vector

\(V\to\pi^{+}\pi^{-}\)
import hazma.form_factors.vector as ffv

# Compute the width V -> pi-pi for mv = 1 GeV with couplings of vector to
# quarks gvuu = 2/3 and gvdd=-1/3
ffv.VectorFormFactorPiPi().width(mv=1e3, gvuu=2.0/3.0, gvdd=-1.0/3.0)
# output: [17.30309111083309]
\(V\to\pi^{+}K^{-}K^{0}\)
import hazma.form_factors.vector as ffv

# Compute the width V -> pi-k-k0 for mv = 1.3 GeV with couplings of vector to
# quarks gvuu = 2/3, gvdd=-1/3, gvss=-1/3
ffv.VectorFormFactorPiKK0().width(mv=1.3e3, gvuu=2.0/3.0, gvdd=-1.0/3.0, gvss=-1.0/3.0)
# output: 0.0004685155427290262 [MeV]
\(V\to\pi^{+}\pi^{-}\pi^{+}\pi^{-}\)
import hazma.form_factors.vector as ffv

# Compute the width V -> pi-k-k0 for mv = 1.3 GeV with couplings of vector to
# quarks gvuu = 2/3, gvdd=-1/3, gvss=-1/3
ffv.VectorFormFactorPiPiPiPi().width(mv=1.3e3, gvuu=2.0/3.0, gvdd=-1.0/3.0)
# output: 10.799920290416575 [MeV]

Dark Matter Annihilation

\(\bar{\chi}\chi\to K^{+}K^{-}\)
import hazma.form_factors.vector as ffv

mx = 300.0  # 300 MeV
q = 1.2e3   # 1.2 GeV
mv = 600.0  # 600 MeV
wv = 1.0    # 1 MeV
gvxx, gvuu, gvdd, gvss = 1.0, 2.0 / 3.0, -1.0 / 3.0, -1.0 / 3.0
ffv.VectorFormFactorKK().cross_section(
    q=q, mx=mx, mv=mv, gvxx=gvxx, wv=wv, gvuu=gvuu, gvdd=gvss, gvss=gvss
)
# output: 5.235577288150283e-09 [MeV^-2]

API

Class

Description

VectorFormFactorPiPi()

\(\pi^{+}\pi^{-}\) form factor.

VectorFormFactorPi0Pi0()

\(\pi^{0}\pi^{0}\) form factor.

VectorFormFactorKK()

\(K^{+}K^{-}\) form factor.

VectorFormFactorK0K0()

\(K^{0}\bar{K}^{0}\) form factor.

VectorFormFactorPi0Gamma()

\(\pi^{0}\gamma\) form factor.

VectorFormFactorPi0Omega()

\(\pi^{0}\omega\) form factor.

VectorFormFactorPi0Phi()

\(\pi^{0}\phi\) form factor.

VectorFormFactorEtaGamma()

\(\eta\gamma\) form factor.

VectorFormFactorEtaOmega()

\(\eta\omega\) form factor.

VectorFormFactorEtaPhi()

\(\eta\phi\) form factor.

VectorFormFactorPi0K0K0()

\(\pi^{0}K^{0}\bar{K}^{0}\) form factor.

VectorFormFactorPi0KpKm()

\(\pi^{0}K^{+}K^{-}\) form factor.

VectorFormFactorPiKK0()

\(\pi^{+}K^{-}K^{0}\) or \(\pi^{-}K^{+}\bar{K}^{0}\) form factor.

VectorFormFactorPiPiEta()

\(\pi^{+}\pi^{-}\eta\) form factor.

VectorFormFactorPiPiEtaPrime()

\(\pi^{+}\pi^{-}\eta'\) form factor.

VectorFormFactorPiPiOmega()

\(\pi^{+}\pi^{-}\omega\) form factor.

VectorFormFactorPi0Pi0Omega()

\(\pi^{0}\pi^{0}\omega\) form factor.

VectorFormFactorPiPiPi0()

\(\pi^{+}\pi^{-}\pi^{0}\) form factor.

VectorFormFactorPiPiPi0Pi0()

\(\pi^{+}\pi^{-}\pi^{0}\pi^{0}\) form factor.

VectorFormFactorPiPiPiPi()

\(\pi^{+}\pi^{-}\pi^{+}\pi^{-}\) form factor.