Lecture 12 – Phase space simulation#
How to perform a phase space simulation with Python?
Prepare the notebook by importing numpy, matplotlib.pyplot, and pylorentz and download the data file (CSV format) from Google Drive.
Two-body decay#
Let’s start with a simple two body decay at rest: \(B^0\rightarrow K^+\pi^-\).
B0_MASS = 5279.65
PION_MASS = 139.57018
KAON_MASS = 493.677
n_events = 100_000
decay = phasespace.nbody_decay(B0_MASS, [PION_MASS, KAON_MASS])
weights, four_momenta = decay.generate(n_events=n_events)
WARNING: All log messages before absl::InitializeLog() is called are written to STDERR
I0000 00:00:1771522795.176880 3133 device_compiler.h:196] Compiled cluster using XLA! This line is logged at most once for the lifetime of the process.
The simulation produces a dictionary (four_momenta) of tf.Tensor objects. Each object can be addressed with particles['p_i'], where i is the number of the \(i\)-th generated particle.
four_momenta
{'p_0': <tf.Tensor: shape=(100000, 4), dtype=float64, numpy=
array([[ -803.67993678, 467.06666165, 2444.06953313, 2618.58901417],
[-1672.1817886 , -918.54150558, 1788.1885131 , 2618.58901417],
[ 347.48298548, 1948.04287955, 1709.36043718, 2618.58901417],
...,
[ 270.40729597, 1131.51605322, -2341.81124465, 2618.58901417],
[ 1091.26615979, 358.17813685, -2349.12221495, 2618.58901417],
[ 1374.26290117, -1795.61596184, 1313.27582253, 2618.58901417]],
shape=(100000, 4))>,
'p_1': <tf.Tensor: shape=(100000, 4), dtype=float64, numpy=
array([[ 803.67993678, -467.06666165, -2444.06953313, 2661.06098583],
[ 1672.1817886 , 918.54150558, -1788.1885131 , 2661.06098583],
[ -347.48298548, -1948.04287955, -1709.36043718, 2661.06098583],
...,
[ -270.40729597, -1131.51605322, 2341.81124465, 2661.06098583],
[-1091.26615979, -358.17813685, 2349.12221495, 2661.06098583],
[-1374.26290117, 1795.61596184, -1313.27582253, 2661.06098583]],
shape=(100000, 4))>}
Each tf.Tensor can be converted to a NumPy array, which can then be converted to a pylorentz.
def to_lorentz(p: tf.Tensor) -> Momentum4:
p = p.numpy().T
return Momentum4(p[3], *p[:3])
pion = to_lorentz(four_momenta["p_0"])
kaon = to_lorentz(four_momenta["p_1"])
These objects can be used to do kinematic computations. Let’s first verify that the invariant mass of the kaon+pion system corresponds to the mass of the mother \(B^0\):
B0 = pion + kaon
np.testing.assert_almost_equal(B0.m.mean(), B0_MASS)
Let’s also plot the momentum components of the two daugther particles.
But it’s monochromatic!! of course it is… it’s a decay at rest. The momentum components are uniformly distributed in the available phase space.
Three-body decay#
Let’s consider now a three body decay like \(B^0\rightarrow K^+\pi^-\pi^0\) and repeat the plot of the relevant kinematic variables. We can also make Dalitz plots this time.
n_events = 50_000
PION0_MASS = 134.9766
decay = phasespace.nbody_decay(B0_MASS, [PION_MASS, PION0_MASS, KAON_MASS])
weights, four_momenta = decay.generate(n_events=n_events)
pim = to_lorentz(four_momenta["p_0"])
pi0 = to_lorentz(four_momenta["p_1"])
kaon = to_lorentz(four_momenta["p_2"])
s1 = (kaon + pim).m2
s2 = (kaon + pi0).m2
s3 = (pim + pi0).m2
Decay chain#
The phasespace package allows to treat also multiple decays. Let’s consider the \(B^0\rightarrow K^{\ast 0}\gamma\) decay, followed by \(K^{\ast 0}\rightarrow \pi^-K^+\). It can be simulated using the following procedure:
from phasespace import GenParticle
B0_MASS = 5279.65
K0STAR_MASS = 895.55
PION_MASS = 139.57018
KAON_MASS = 493.677
GAMMA_MASS = 0.0
Kp = GenParticle("K+", KAON_MASS)
pim = GenParticle("pi-", PION_MASS)
Kstar = GenParticle("KStar", K0STAR_MASS).set_children(Kp, pim)
gamma = GenParticle("gamma", GAMMA_MASS)
B0 = GenParticle("B0", B0_MASS).set_children(Kstar, gamma)
weights, four_momenta = B0.generate(n_events=100_000)
four_momenta
{'KStar': <tf.Tensor: shape=(100000, 4), dtype=float64, numpy=
array([[ 94.75778934, -1393.8028554 , -2149.83128152, 2715.77793272],
[-2360.02600944, 46.79944535, 1000.76322076, 2715.77793272],
[ 100.52124732, -984.24522485, 2365.290002 , 2715.77793272],
...,
[-1012.95910278, 1606.48478225, 1722.37059829, 2715.77793272],
[ 1565.82882328, -2012.13627529, -270.05126032, 2715.77793272],
[ -840.61741583, 1974.59350643, 1402.77682542, 2715.77793272]],
shape=(100000, 4))>,
'gamma': <tf.Tensor: shape=(100000, 4), dtype=float64, numpy=
array([[ -94.75778934, 1393.8028554 , 2149.83128152, 2563.87206728],
[ 2360.02600944, -46.79944535, -1000.76322076, 2563.87206728],
[ -100.52124732, 984.24522485, -2365.290002 , 2563.87206728],
...,
[ 1012.95910278, -1606.48478225, -1722.37059829, 2563.87206728],
[-1565.82882328, 2012.13627529, 270.05126032, 2563.87206728],
[ 840.61741583, -1974.59350643, -1402.77682542, 2563.87206728]],
shape=(100000, 4))>,
'K+': <tf.Tensor: shape=(100000, 4), dtype=float64, numpy=
array([[ 142.93198093, -1416.56649873, -2000.38478234, 2504.46530306],
[-2143.14065442, -4.570941 , 682.10819271, 2302.6205343 ],
[ 309.55113446, -865.42786578, 1698.34406967, 1993.20767985],
...,
[-1079.12932539, 1395.72886108, 1433.02674139, 2325.91095626],
[ 1388.10371883, -1977.5008459 , -304.02391684, 2484.65069016],
[ -238.60697106, 680.53305591, 372.56884341, 950.04370964]],
shape=(100000, 4))>,
'pi-': <tf.Tensor: shape=(100000, 4), dtype=float64, numpy=
array([[ -48.17419159, 22.76364333, -149.44649918, 211.31262966],
[-216.88535503, 51.37038635, 318.65502805, 413.15739842],
[-209.02988714, -118.81735907, 666.94593233, 722.57025287],
...,
[ 66.17022261, 210.75592117, 289.3438569 , 389.86697647],
[ 177.72510445, -34.6354294 , 33.97265652, 231.12724257],
[-602.01044477, 1294.06045052, 1030.20798201, 1765.73422308]],
shape=(100000, 4))>}
gamma = to_lorentz(four_momenta["gamma"])
pion = to_lorentz(four_momenta["pi-"])
kaon = to_lorentz(four_momenta["K+"])
Kstar = to_lorentz(four_momenta["KStar"])
Let’s build the Dalitz plots matching particle pairs. The particles measured in the final state are \(K^-,\; \pi^-\) and \(\gamma\).
s1 = (pion + kaon).m2
s2 = (gamma + kaon).m2
s3 = (gamma + pion).m2
Width distribution#
These distributions aren’t so interesting, because the masses of each particle are one fixed value. So let’s simulate a more realistic \(K^\ast\) particle; not monochromatic, but with a width of 47 MeV.[1] The mass is extracted from a Gaussian distribution centered at the B0_MASS value and with \(\sigma = 47/2.36 \sim 20\) MeV. See more info on how to do this with the phasespace package here.
import tensorflow as tf
import tensorflow_probability as tfp
K0STAR_WIDTH = 47 / 2.36
def kstar_mass(min_mass, max_mass, n_events):
min_mass = tf.cast(min_mass, tf.float64)
max_mass = tf.cast(max_mass, tf.float64)
kstar_mass_cast = tf.cast(K0STAR_MASS, dtype=tf.float64)
tf.cast(K0STAR_WIDTH, tf.float64)
tf.broadcast_to(kstar_mass_cast, shape=(n_events,))
return tfp.distributions.TruncatedNormal(
loc=K0STAR_MASS,
scale=K0STAR_WIDTH,
low=min_mass,
high=max_mass,
).sample()
K = GenParticle("K+", KAON_MASS)
pion = GenParticle("pi-", PION_MASS)
Kstar = GenParticle("KStar", kstar_mass).set_children(K, pion)
gamma = GenParticle("gamma", GAMMA_MASS)
B0 = GenParticle("B0", B0_MASS).set_children(Kstar, gamma)
weights, four_momenta = B0.generate(n_events=100_000)
gamma = to_lorentz(four_momenta["gamma"])
pion = to_lorentz(four_momenta["pi-"])
kaon = to_lorentz(four_momenta["K+"])
Kstar = to_lorentz(four_momenta["KStar"])
Now you have all the 4-vectors to plot the invariant mass distributions for the different steps of the decay chains.
s1 = (pion + kaon).m2
s2 = (gamma + kaon).m2
s3 = (gamma + pion).m2