Dynamics

By default, the dynamic terms in an amplitude model are set to \(1\) by the HelicityAmplitudeBuilder. The method set_dynamics() can then be used to set dynamics lineshapes for specific resonances. The dynamics.builder module provides some tools to set standard lineshapes (see below), but it is also possible to set custom dynamics.

The standard lineshapes provided by AmpForm are illustrated below. For more info, have a look at the following pages:

%matplotlib widget
import logging
import warnings

import matplotlib.pyplot as plt
import mpl_interactions.ipyplot as iplt
import numpy as np
import sympy as sp
from IPython.display import Math

import symplot

logging.basicConfig()
logging.getLogger().setLevel(logging.ERROR)

warnings.filterwarnings("ignore")

Form factor

AmpForm uses Blatt-Weisskopf functions \(B_L\) as barrier factors (also called form factors, see BlattWeisskopfSquared):

from ampform.dynamics import BlattWeisskopfSquared

L = sp.Symbol("L", integer=True)
z = sp.Symbol("z", real=True)
ff2 = BlattWeisskopfSquared(L, z)
Math(sp.multiline_latex(ff2, ff2.doit(), environment="eqnarray"))
\[\begin{split}\displaystyle \begin{eqnarray} B_{L}^2\left(z\right) & = & \begin{cases} 1 & \text{for}\: L = 0 \\\frac{2 z}{z + 1} & \text{for}\: L = 1 \\\frac{13 z^{2}}{9 z + \left(z - 3\right)^{2}} & \text{for}\: L = 2 \\\frac{277 z^{3}}{z \left(z - 15\right)^{2} + \left(2 z - 5\right) \left(18 z - 45\right)} & \text{for}\: L = 3 \\\frac{12746 z^{4}}{25 z \left(2 z - 21\right)^{2} + \left(z^{2} - 45 z + 105\right)^{2}} & \text{for}\: L = 4 \\\frac{998881 z^{5}}{z^{5} + 15 z^{4} + 315 z^{3} + 6300 z^{2} + 99225 z + 893025} & \text{for}\: L = 5 \\\frac{118394977 z^{6}}{z^{6} + 21 z^{5} + 630 z^{4} + 18900 z^{3} + 496125 z^{2} + 9823275 z + 108056025} & \text{for}\: L = 6 \\\frac{19727003738 z^{7}}{z^{7} + 28 z^{6} + 1134 z^{5} + 47250 z^{4} + 1819125 z^{3} + 58939650 z^{2} + 1404728325 z + 18261468225} & \text{for}\: L = 7 \\\frac{4392846440677 z^{8}}{z^{8} + 36 z^{7} + 1890 z^{6} + 103950 z^{5} + 5457375 z^{4} + 255405150 z^{3} + 9833098275 z^{2} + 273922023375 z + 4108830350625} & \text{for}\: L = 8 \end{cases} \end{eqnarray}\end{split}\]

The Blatt-Weisskopf form factor is used to ‘dampen’ the breakup-momentum at threshold and when going to infinity. A usual choice for \(z\) is therefore \(z=q^2d^2\) with \(q^2\) the BreakupMomentumSquared and \(d\) the impact parameter (also called meson radius):

from ampform.dynamics import BreakupMomentumSquared

m, m_a, m_b, d = sp.symbols("m, m_a, m_b, d")
s = m ** 2
q_squared = BreakupMomentumSquared(s, m_a, m_b)
ff2 = BlattWeisskopfSquared(L, z=q_squared * d ** 2)
np_blatt_weisskopf, sliders = symplot.prepare_sliders(
    plot_symbol=m,
    expression=ff2.doit(),
)
np_breakup_momentum = sp.lambdify((m, L, d, m_a, m_b), q_squared.doit())
../_images/dynamics_13_0.svg

Relativistic Breit-Wigner

AmpForm has two types of relativistic Breit-Wigner functions. Both are compared below ― for more info, see the links to the API.

Without form factor

The ‘normal’ relativistic_breit_wigner() looks as follows:

from ampform.dynamics import relativistic_breit_wigner

m, m0, w0 = sp.symbols("m, m0, Gamma0")
rel_bw = relativistic_breit_wigner(s=m ** 2, mass0=m0, gamma0=w0)
rel_bw
\[\displaystyle \frac{\Gamma_{0} m_{0}}{- i \Gamma_{0} m_{0} - m^{2} + m_{0}^{2}}\]

With form factor

The relativistic Breit-Wigner can be adapted slightly, so that its amplitude goes to zero at threshold (\(m_0 = m_a + m_b\)) and that it becomes normalizable. This is done with form factors and can be obtained with the function relativistic_breit_wigner_with_ff():

from ampform.dynamics import PhaseSpaceFactor  # noqa: F401
from ampform.dynamics import BreakupMomentumSquared, ComplexSqrt


def breakup_momentum(s: sp.Symbol, m_a: sp.Symbol, m_b: sp.Symbol) -> sp.Expr:
    q_squared = BreakupMomentumSquared(s, m_a, m_b)
    return ComplexSqrt(q_squared)
from ampform.dynamics import (
    PhaseSpaceFactorAnalytic,
    relativistic_breit_wigner_with_ff,
)

rel_bw_with_ff = relativistic_breit_wigner_with_ff(
    s=s,
    mass0=m0,
    gamma0=w0,
    m_a=m_a,
    m_b=m_b,
    angular_momentum=L,
    meson_radius=1,
    phsp_factor=PhaseSpaceFactorAnalytic,
)
rel_bw_with_ff
\[\displaystyle \frac{\Gamma_{0} m_{0} \sqrt{B_{L}^2\left(q^2\left(m^{2}\right)\right)}}{- m^{2} + m_{0}^{2} - i m_{0} \Gamma_{0}\left(m^{2}\right)}\]

Here, \(\Gamma(m)\) is the CoupledWidth (also called running width or mass-dependent width), defined as:

from ampform.dynamics import CoupledWidth

L = sp.Symbol("L", integer=True)
width = CoupledWidth(
    s=s,
    mass0=m0,
    gamma0=w0,
    m_a=m_a,
    m_b=m_b,
    angular_momentum=L,
    meson_radius=1,
    phsp_factor=PhaseSpaceFactorAnalytic,
)
Math(sp.multiline_latex(width, width.evaluate(), environment="eqnarray"))
\[\displaystyle \begin{eqnarray} \Gamma_{0}\left(m^{2}\right) & = & \frac{\Gamma_{0} B_{L}^2\left(q^2\left(m^{2}\right)\right) \rho^\mathrm{ac}\left(m^{2}\right)}{B_{L}^2\left(q^2\left(m_{0}^{2}\right)\right) \rho^\mathrm{ac}\left(m_{0}^{2}\right)} \end{eqnarray}\]

It is possible to choose different formulations for the phase space factor \(\rho\), see Analytic continuation.

Analytic continuation

The following shows the effect of Analytic continuation a on relativistic Breit-Wigner:

from ampform.dynamics import PhaseSpaceFactorComplex

# Two types of relativistic Breit-Wigners
rel_bw_with_ff = relativistic_breit_wigner_with_ff(
    s=s,
    mass0=m0,
    gamma0=w0,
    m_a=m_a,
    m_b=m_b,
    angular_momentum=L,
    meson_radius=d,
    phsp_factor=PhaseSpaceFactorComplex,
)
rel_bw_with_ff_ac = relativistic_breit_wigner_with_ff(
    s=s,
    mass0=m0,
    gamma0=w0,
    m_a=m_a,
    m_b=m_b,
    angular_momentum=L,
    meson_radius=d,
    phsp_factor=PhaseSpaceFactorAnalytic,
)

# Lambdify
np_rel_bw_with_ff, sliders = symplot.prepare_sliders(
    plot_symbol=m,
    expression=rel_bw_with_ff.doit(),
)
np_rel_bw_with_ff_ac = sp.lambdify(
    args=(m, w0, L, d, m0, m_a, m_b),
    expr=rel_bw_with_ff_ac.doit(),
)
np_rel_bw = sp.lambdify(
    args=(m, w0, L, d, m0, m_a, m_b),
    expr=rel_bw.doit(),
)

# Set sliders
plot_domain = np.linspace(start=0, stop=4, num=500)
sliders.set_ranges(
    m0=(0, 4, 200),
    Gamma0=(0, 1, 100),
    L=(0, 8),
    m_a=(0, 2, 200),
    m_b=(0, 2, 200),
    d=(0, 5),
)
sliders.set_values(
    m0=1.8,
    Gamma0=0.6,
    L=1,
    m_a=0.6,
    m_b=0.6,
    d=1,
)

fig, axes = plt.subplots(
    nrows=2,
    figsize=(8, 6),
    sharex=True,
    #     gridspec_kw={"height_ratios": [1, 1.8]},
)
ax_ff, ax_ac = axes
ax_ac.set_xlabel("$m$")
for ax in axes:
    ax.axhline(0, c="gray", linewidth=0.5)

rho_c = PhaseSpaceFactorComplex(m, m_a, m_b)
ax_ff.set_title(
    f"'Complex' phase space factor: ${sp.latex(rho_c.evaluate())}$"
)
ax_ac.set_title("Analytic continuation of the phase space factor")

ylim = "auto"  # (-0.6, 1.2)
controls = iplt.plot(
    plot_domain,
    lambda *args, **kwargs: np_rel_bw_with_ff(*args, **kwargs).real,
    label="real",
    **sliders,
    ylim=ylim,
    ax=ax_ff,
)
iplt.plot(
    plot_domain.real,
    lambda *args, **kwargs: np_rel_bw_with_ff(*args, **kwargs).imag,
    label="imaginary",
    controls=controls,
    ylim=ylim,
    ax=ax_ff,
)
iplt.plot(
    plot_domain.real,
    lambda *args, **kwargs: np.abs(np_rel_bw_with_ff(*args, **kwargs)) ** 2,
    label="absolute",
    controls=controls,
    ylim=ylim,
    ax=ax_ff,
    c="black",
    linestyle="dotted",
)


iplt.plot(
    plot_domain.real,
    lambda *args, **kwargs: np_rel_bw_with_ff_ac(*args, **kwargs).real,
    label="real",
    controls=controls,
    ylim=ylim,
    ax=ax_ac,
)
iplt.plot(
    plot_domain.real,
    lambda *args, **kwargs: np_rel_bw_with_ff_ac(*args, **kwargs).imag,
    label="imaginary",
    controls=controls,
    ylim=ylim,
    ax=ax_ac,
)
iplt.plot(
    plot_domain.real,
    lambda *args, **kwargs: np.abs(np_rel_bw_with_ff_ac(*args, **kwargs)) ** 2,
    label="absolute",
    controls=controls,
    ylim=ylim,
    ax=ax_ac,
    c="black",
    linestyle="dotted",
)

for ax in axes:
    iplt.axvline(
        controls["m0"],
        ax=ax,
        c="red",
        label=f"${sp.latex(m0)}$",
        alpha=0.3,
    )
    iplt.axvline(
        lambda m_a, m_b, **kwargs: m_a + m_b,
        controls=controls,
        ax=ax,
        c="black",
        alpha=0.3,
        label=f"${sp.latex(m_a)} + {sp.latex(m_b)}$",
    )
ax_ac.legend(loc="upper right")
fig.tight_layout()
plt.show()
../_images/dynamics_29_0.svg