PhaseSpace

https://joss.theoj.org/papers/10.21105/joss.01570/status.svg https://zenodo.org/badge/DOI/10.5281/zenodo.2591993.svg https://img.shields.io/pypi/status/phasespace.svg https://img.shields.io/pypi/pyversions/phasespace.svg https://github.com/zfit/phasespace/workflows/tests/badge.svg https://codecov.io/gh/zfit/phasespace/branch/master/graph/badge.svg Documentation Status Gitter chat https://badge.fury.io/py/phasespace.svg

Python implementation of the Raubold and Lynch method for n-body events using TensorFlow as a backend.

The code is based on the GENBOD function (W515 from CERNLIB), documented in [1] and tries to follow it as closely as possible.

Detailed documentation, including the API, can be found in https://phasespace.readthedocs.io. Don’t hesitate to join our gitter channel for questions and comments.

If you use phasespace in a scientific publication we would appreciate citations to the JOSS publication:

@article{puig_eschle_phasespace-2019,
  title = {phasespace: n-body phase space generation in Python},
  doi = {10.21105/joss.01570},
  url = {https://doi.org/10.21105/joss.01570},
  year = {2019},
  month = {oct},
  publisher = {The Open Journal},
  author = {Albert Puig and Jonas Eschle},
  journal = {Journal of Open Source Software}
}

Free software: BSD-3-Clause.

[1] F. James, Monte Carlo Phase Space, CERN 68-15 (1968)

Why?

Lately, data analysis in High Energy Physics (HEP), traditionally performed within the ROOT ecosystem, has been moving more and more towards Python. The possibility of carrying out purely Python-based analyses has become real thanks to the development of many open source Python packages, which have allowed to replace most ROOT functionality with Python-based packages.

One of the aspects where this is still not possible is in the random generation of n-body phase space events, which are widely used in the field, for example to study kinematics of the particle decays of interest, or to perform importance sampling in the case of complex amplitude models. This has been traditionally done with the TGenPhaseSpace class, which is based of the GENBOD function of the CERNLIB FORTRAN libraries and which requires a full working ROOT installation.

This package aims to address this issue by providing a TensorFlow-based implementation of such a function to generate n-body decays without requiring a ROOT installation. Additionally, an oft-needed functionality to generate complex decay chains, not included in TGenPhaseSpace, is also offered, leaving room for decaying resonances (which don’t have a fixed mass, but can be seen as a broad peak).

Installing

To install with pip:

$ pip install phasespace

This is the preferred method to install phasespace, as it will always install the most recent stable release. To install the necessary dependencies to be used with DecayLanguage, use

$ pip install "phasespace[fromdecay]"

How to use

Phasespace can directly be used to generate from a DecayChain using the DecayLanguage package as explained in the tutorial.

The generation of simple n-body decays can be done using the nbody_decay shortcut to create a decay chain with a very simple interface: one needs to pass the mass of the top particle and the masses of the children particle as a list, optionally giving the names of the particles. Then, the generate method can be used to produce the desired sample. For example, to generate \(B^0\to K\pi\), we would do:

import phasespace

B0_MASS = 5279.65
PION_MASS = 139.57018
KAON_MASS = 493.677

weights, particles = phasespace.nbody_decay(B0_MASS,
                                            [PION_MASS, KAON_MASS]).generate(n_events=1000)

Behind the scenes, this function runs the TensorFlow graph. It returns tf.Tensor, which, as TensorFlow 2.x is in eager mode, is basically a numpy array. Any tf.Tensor can be explicitly converted to a numpy array by calling tf.Tensor.numpy() on it. The generate function returns a tf.Tensor of 1000 elements in the case of weights and a list of n particles (2) arrays of (1000, 4) shape, where each of the 4-dimensions corresponds to one of the components of the generated Lorentz 4-vector. All particles are generated in the rest frame of the top particle; boosting to a certain momentum (or list of momenta) can be achieved by passing the momenta to the boost_to argument.

Sequential decays can be handled with the GenParticle class (used internally by generate) and its set_children method. As an example, to build the \(B^{0}\to K^{*}\gamma\) decay in which \(K^*\to K\pi\), we would write:

from phasespace import GenParticle

B0_MASS = 5279.65
KSTARZ_MASS = 895.55
PION_MASS = 139.57018
KAON_MASS = 493.677

kaon = GenParticle('K+', KAON_MASS)
pion = GenParticle('pi-', PION_MASS)
kstar = GenParticle('K*', KSTARZ_MASS).set_children(kaon, pion)
gamma = GenParticle('gamma', 0)
bz = GenParticle('B0', B0_MASS).set_children(kstar, gamma)

weights, particles = bz.generate(n_events=1000)

Where we have used the fact that set_children returns the parent particle. In this case, particles is a dict with the particle names as keys:

>>> particles
{'K*': array([[ 1732.79325872, -1632.88873127,   950.85807735,  2715.78804872],
       [-1633.95329448,   239.88921123, -1961.0402768 ,  2715.78804872],
       [  407.15613764, -2236.6569286 , -1185.16616251,  2715.78804872],
       ...,
       [ 1091.64603395, -1301.78721269,  1920.07503991,  2715.78804872],
       [ -517.3125083 ,  1901.39296899,  1640.15905194,  2715.78804872],
       [  656.56413668,  -804.76922982,  2343.99214816,  2715.78804872]]),
 'K+': array([[  750.08077976,  -547.22569019,   224.6920906 ,  1075.30490935],
       [-1499.90049089,   289.19714633, -1935.27960292,  2514.43047106],
       [   97.64746732, -1236.68112923,  -381.09526192,  1388.47607911],
       ...,
       [  508.66157459,  -917.93523639,  1474.7064148 ,  1876.11771642],
       [ -212.28646168,   540.26381432,   610.86656669,   976.63988936],
       [  177.16656666,  -535.98777569,   946.12636904,  1207.28744488]]),
 'gamma': array([[-1732.79325872,  1632.88873127,  -950.85807735,  2563.79195128],
       [ 1633.95329448,  -239.88921123,  1961.0402768 ,  2563.79195128],
       [ -407.15613764,  2236.6569286 ,  1185.16616251,  2563.79195128],
       ...,
       [-1091.64603395,  1301.78721269, -1920.07503991,  2563.79195128],
       [  517.3125083 , -1901.39296899, -1640.15905194,  2563.79195128],
       [ -656.56413668,   804.76922982, -2343.99214816,  2563.79195128]]),
 'pi-': array([[  982.71247896, -1085.66304109,   726.16598675,  1640.48313937],
       [ -134.0528036 ,   -49.3079351 ,   -25.76067389,   201.35757766],
       [  309.50867032,  -999.97579937,  -804.0709006 ,  1327.31196961],
       ...,
       [  582.98445936,  -383.85197629,   445.36862511,   839.6703323 ],
       [ -305.02604662,  1361.12915468,  1029.29248526,  1739.14815935],
       [  479.39757002,  -268.78145413,  1397.86577911,  1508.50060384]])}

The GenParticle class is able to cache the graphs so it is possible to generate in a loop without overhead:

for i in range(10):
    weights, particles = bz.generate(n_events=1000)
    ...
    (do something with weights and particles)
    ...

This way of generating is recommended in the case of large samples, as it allows to benefit from parallelisation while at the same time keep the memory usage low.

If we want to operate with the TensorFlow graph instead, we can use the generate_tensor method of GenParticle, which has the same signature as generate.

More examples can be found in the tests folder and in the documentation.

Physics validation

Physics validation is performed continuously in the included tests (tests/test_physics.py), run through GitHub Actions. This validation is performed at two levels:

  • In simple n-body decays, the results of phasespace are checked against TGenPhaseSpace.

  • For sequential decays, the results of phasespace are checked against RapidSim, a “fast Monte Carlo generator for simulation of heavy-quark hadron decays”. In the case of resonances, differences are expected because our tests don’t include proper modelling of their mass shape, as it would require the introduction of further dependencies. However, the results of the comparison can be expected visually.

The results of all physics validation performed by the tests_physics.py test are written in tests/plots.

Contributing

Contributions are always welcome, please have a look at the Contributing guide.

Credits

Development Lead

Developers

Contributors

None yet. Why not be the first?

Table of Contents

Usage

The base of phasespace is the GenParticle object. This object, which represents a particle, either stable or decaying, has only one mandatory argument, its name.

In most cases (except for the top particle of a decay), one wants to also specify its mass, which can be either be array-like, or a function. Functions are used to specify the mass of particles such as resonances, which are not fixed but vary according to a broad distribution. These mass functions get three arguments, and must return an array-like object of the same shape as the input arguments.

  • The minimum mass allowed by the decay chain, which will be of shape (n_events,).

  • The maximum mass available, which will be of shape (n_events,).

  • The number of events to generate.

This function signature allows to handle threshold effects cleanly, giving enough information to produce kinematically allowed decays (NB: phasespace will throw an error if a kinematically forbidden decay is requested).

A simple example

With these considerations in mind, one can build a decay chain by using the set_children method of the :py:class:~`phasespace.GenParticle` class. As an example, to build the \(B^{0}\to K^{*}\gamma\) decay in which \(K^*\to K\pi\) with a fixed mass, one would write:

from phasespace import GenParticle
import numpy as np
from particle import literals as lp
import vector

B0_MASS = lp.B_0.mass
KSTARZ_MASS = lp.Kst_892_0.mass
PION_MASS = lp.pi_plus.mass
KAON_MASS = lp.K_plus.mass

pion = GenParticle('pi-', PION_MASS)
kaon = GenParticle('K+', KAON_MASS)
kstar = GenParticle('K*', KSTARZ_MASS).set_children(pion, kaon)
gamma = GenParticle('gamma', 0)
bz = GenParticle('B0', B0_MASS).set_children(kstar, gamma)

Phasespace events can be generated using the generate method, which gets the number of events to generate as input. The method returns:

  • The normalized weights of each event, as an array of dimension (n_events,).

  • The 4-momenta of the generated particles as values of a dictionary with the particle name as key. These momenta are either expressed as arrays of dimension (n_events, 4) or :py:class:~`vector.Vector4D` objects, depending on the as_vectors flag given to generate.

N_EVENTS = 1000

weights, particles = bz.generate(n_events=N_EVENTS, as_vectors=True)
# or
weights, particles = bz.generate(n_events=N_EVENTS)

(The array-like objects can always directly be converted to a numpy array (if really needed) through np.asarray(obj).)

Boosting the particles

The particles are generated in the rest frame of the top particle. To produce them at a given momentum of the top particle, one can pass these momenta with the boost_to argument in generate. This latter approach can be useful if the momentum of the top particle is generated according to some distribution, for example the kinematics of the LHC (see test_kstargamma_kstarnonresonant_lhc and test_k1gamma_kstarnonresonant_lhc in tests/test_physics.py to see how this could be done).

The boost_to argument can be a 4-momentum array of shape (n_events, 4) with (px, py, yz, energy) or a :py:class:~`vector.Vector4D` (both a momentum and a Lorentz vector).

N_EVENTS = 1000

# Generate the top particle with a momentum of 100 GeV
top_momentum = np.array([0, 0, 100, np.sqrt(100**2 + B0_MASS**2)])
# or
top_momentum = vector.array({'px': [0], 'py': [0], 'pz': [100], 'E': [np.sqrt(100**2 + B0_MASS**2)]})
weights, particles = bz.generate(n_events=N_EVENTS, boost_to=top_momentum)
Weights

Additionally, it is possible to obtain the unnormalized weights by using the generate_unnormalized flag in generate. In this case, the method returns the unnormalized weights, the per-event maximum weight and the particle dictionary.

print(particles)
{'K*': <tf.Tensor: shape=(1000, 4), dtype=float64, numpy=
array([[ -947.16415007, -1122.95423297, -2050.20462169,  2676.47061068],
       [ 1308.52150432, -2195.01437824,   259.20092213,  2720.20432477],
       [ 1239.03247656,  1512.3898589 ,  1710.33240524,  2747.68471666],
       ...,
       [  478.28501777,   669.02354668, -2377.39504983,  2670.27453452],
       [-1191.00266947, -1248.08736695, -1845.5451777 ,  2680.346291  ],
       [ 1637.83254451, -1112.00955775,  1680.96549096,  2747.12858903]])>, 'gamma': <tf.Tensor: shape=(1000, 4), dtype=float64, numpy=
array([[  947.16415007,  1122.95423297,  2150.20462169,  2604.13633508],
       [-1308.52150432,  2195.01437824,  -159.20092213,  2560.40262099],
       [-1239.03247656, -1512.3898589 , -1610.33240524,  2532.92222909],
       ...,
       [ -478.28501777,  -669.02354668,  2477.39504983,  2610.33241124],
       [ 1191.00266947,  1248.08736695,  1945.5451777 ,  2600.26065476],
       [-1637.83254451,  1112.00955775, -1580.96549096,  2533.47835673]])>, 'pi-': <tf.Tensor: shape=(1000, 4), dtype=float64, numpy=
array([[ -82.02408852, -436.44256702, -868.28061109,  985.19093519],
       [ 775.60899228, -988.36550687,  341.35701015, 1309.36633002],
       [ 568.33555061, 1028.44597851,  859.43153282, 1462.46671148],
       ...,
       [ 130.9868867 ,   16.02328127,  -82.36496396,  208.99328   ],
       [-802.00447657, -464.83852598, -953.1837641 , 1336.90882919],
       [1102.84004854, -748.43589975, 1218.11217852, 1810.99134231]])>, 'K+': <tf.Tensor: shape=(1000, 4), dtype=float64, numpy=
array([[ -865.14006155,  -686.51166595, -1181.9240106 ,  1691.27967549],
       [  532.91251205, -1206.64887138,   -82.15608802,  1410.83799474],
       [  670.69692595,   483.94388038,   850.90087242,  1285.21800519],
       ...,
       [  347.29813106,   653.00026542, -2295.03008587,  2461.28125452],
       [ -388.99819291,  -783.24884097,  -892.3614136 ,  1343.4374618 ],
       [  534.99249597,  -363.573658  ,   462.85331244,   936.13724672]])>}

Iterative generation can be performed using normal python loops without loss in performance:

for i in range(5):
   weights, particles = bz.generate(n_events=100)
   # ...
   # (do something with weights and particles)
   # ...
Resonances with variable mass

To generate the mass of a resonance, we need to give a function as its mass instead of a floating number. This function should take as input the per-event lower mass allowed, per-event upper mass allowed and the number of events, and should return an array-like object with the generated masses and shape (nevents,). Well suited for this task are the TensorFlow Probability distributions or, for more customized mass shapes, the zfit pdfs (currently an experimental feature is needed, contact the zfit developers to learn more).

Following with the same example as above, and approximating the resonance shape by a gaussian, we could write the \(B^{0}\to K^{*}\gamma\) decay chain as (more details can be found in tests/helpers/decays.py):

from phasespace import numpy as tnp
import tensorflow_probability as tfp
from phasespace import GenParticle

KSTARZ_MASS = 895.81
KSTARZ_WIDTH = 47.4

def kstar_mass(min_mass, max_mass, n_events):
   min_mass = tnp.asarray(min_mass, tnp.float64)
   max_mass = tnp.asarray(max_mass, tnp.float64)
   kstar_width_cast = tnp.asarray(KSTARZ_WIDTH, tnp.float64)
   kstar_mass_cast = tnp.asarray(KSTARZ_MASS, tnp.float64)

   kstar_mass = tnp.broadcast_to(kstar_mass_cast, shape=(n_events,))
   if KSTARZ_WIDTH > 0:
       kstar_mass = tfp.distributions.TruncatedNormal(loc=kstar_mass,
                                                      scale=kstar_width_cast,
                                                      low=min_mass,
                                                      high=max_mass).sample()
   return kstar_mass

bz = GenParticle('B0', B0_MASS).set_children(GenParticle('K*0', mass=kstar_mass)
                                            .set_children(GenParticle('K+', mass=KAON_MASS),
                                                          GenParticle('pi-', mass=PION_MASS)),
                                            GenParticle('gamma', mass=0.0))

bz.generate(n_events=500)
Shortcut for simple decays

The generation of simple n-body decay chains can be done using the nbody_decay function of phasespace, which takes

  • The mass of the top particle.

  • The mass of children particles as a list.

  • The name of the top particle (optional).

  • The names of the children particles (optional).

If the names are not given, top and p_{i} are assigned. For example, to generate \(B^0\to K\pi\), one would do:

import phasespace
from particle import literals as lp

N_EVENTS = 1000

B0_MASS = lp.B_0.mass
PION_MASS = lp.pi_plus.mass
KAON_MASS = lp.K_plus.mass

decay = phasespace.nbody_decay(B0_MASS, [PION_MASS, KAON_MASS],
                               top_name="B0", names=["pi", "K"])
weights, particles = decay.generate(n_events=N_EVENTS)

In this example, decay is simply a GenParticle with the corresponding children.

Eager execution

By default, phasespace uses JIT (just-in-time) compilation of TensorFlow to greatly speed up the generation of events. Simplified, this means that the first time a decay is generated, a symbolic array without a concrete value is used and the computation is remembered. As a user calling the function, you will not notice this, the output will be the same as if the function was executed eagerly. The consequence is two-fold: on one hand the initial overhead is higher with a significant speedup for subsequent generations, on the other hand, the values of the generated particles inside the function are not available in pure Python (e.g. for debugging basically).

If you need to debug the internals, using tf.config.run_functions_eagerly(True) (or the environment variable "PHASESPACE_EAGER=1") will make everything run numpy-like.

Random numbers

The random number generation inside phasespace is transparent in order to allow for deterministic behavior if desired. A function that uses random number generation inside always takes a seed (or rng) argument. The behavior is as follows

  • if no seed is given, the global random number generator of TensorFlow will be used. Setting this instance explicitly or by setting the seed via tf.random.set_seed allows for a deterministic execution of a whole _script_.

  • if the seed is a number it will be used to create a random number generator from this. Using the same seed again will result in the same output.

  • if the seed is an instance of tf.random.Generator, this instance will directly be used and advances an undefined number of steps.

Tutorial for GenMultiDecay class

This tutorial shows how phasespace.fromdecay.GenMultiDecay can be used.

In order to use this functionality, you need to install the extra dependencies, for example through pip install phasespace[fromdecay].

This submodule makes it possible for phasespace and DecayLanguage to work together. More generally, GenMultiDecay can also be used as a high-level interface for simulating particles that can decay in multiple different ways.

# Import libraries
from pprint import pprint

import zfit
from particle import Particle
from decaylanguage import DecFileParser, DecayChainViewer, DecayChain, DecayMode
import tensorflow as tf

from phasespace.fromdecay import GenMultiDecay
/home/docs/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/zfit/__init__.py:59: UserWarning: TensorFlow warnings are by default suppressed by zfit. In order to show them, set the environment variable ZFIT_DISABLE_TF_WARNINGS=0. In order to suppress the TensorFlow warnings AND this warning, set ZFIT_DISABLE_TF_WARNINGS=1.
  warnings.warn(
Quick Intro to DecayLanguage

DecayLanguage can be used to parse and view .dec files. These files contain information about how a particle decays and with which probability. For more information about DecayLanguage and .dec files, see the DecayLanguage documentation.

We will begin by parsing a .dec file using DecayLanguage:

parser = DecFileParser("../tests/fromdecay/example_decays.dec")
parser.parse()

From the parser variable, one can access a certain decay for a particle using parser.build_decay_chains. This will be a dict that contains all information about how the mother particle, daughter particles etc. decay.

pi0_chain = parser.build_decay_chains("pi0")
pprint(pi0_chain)
{
'pi0'
: 
[
{
'bf'
: 
0.988228297
,
          
'fs'
: 
['gamma', 'gamma']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
0.011738247
,
          
'fs'
: 
['e+', 'e-', 'gamma']
,
          
'model'
: 
'PI0_DALITZ'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
3.3392e-05
,
          
'fs'
: 
['e+', 'e+', 'e-', 'e-']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
6.5e-08
,
          
'fs'
: 
['e+', 'e-']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
]
}

This dict can also be displayed in a more human-readable way using DecayChainViewer:

DecayChainViewer(pi0_chain)
---------------------------------------------------------------------------
FileNotFoundError                         Traceback (most recent call last)
File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/execute.py:76, in run_check(cmd, input_lines, encoding, quiet, **kwargs)
     75         kwargs['stdout'] = kwargs['stderr'] = subprocess.PIPE
---> 76     proc = _run_input_lines(cmd, input_lines, kwargs=kwargs)
     77 else:

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/execute.py:96, in _run_input_lines(cmd, input_lines, kwargs)
     95 def _run_input_lines(cmd, input_lines, *, kwargs):
---> 96     popen = subprocess.Popen(cmd, stdin=subprocess.PIPE, **kwargs)
     98     stdin_write = popen.stdin.write

File ~/.asdf/installs/python/3.11.6/lib/python3.11/subprocess.py:1026, in Popen.__init__(self, args, bufsize, executable, stdin, stdout, stderr, preexec_fn, close_fds, shell, cwd, env, universal_newlines, startupinfo, creationflags, restore_signals, start_new_session, pass_fds, user, group, extra_groups, encoding, errors, text, umask, pipesize, process_group)
   1023             self.stderr = io.TextIOWrapper(self.stderr,
   1024                     encoding=encoding, errors=errors)
-> 1026     self._execute_child(args, executable, preexec_fn, close_fds,
   1027                         pass_fds, cwd, env,
   1028                         startupinfo, creationflags, shell,
   1029                         p2cread, p2cwrite,
   1030                         c2pread, c2pwrite,
   1031                         errread, errwrite,
   1032                         restore_signals,
   1033                         gid, gids, uid, umask,
   1034                         start_new_session, process_group)
   1035 except:
   1036     # Cleanup if the child failed starting.

File ~/.asdf/installs/python/3.11.6/lib/python3.11/subprocess.py:1950, in Popen._execute_child(self, args, executable, preexec_fn, close_fds, pass_fds, cwd, env, startupinfo, creationflags, shell, p2cread, p2cwrite, c2pread, c2pwrite, errread, errwrite, restore_signals, gid, gids, uid, umask, start_new_session, process_group)
   1949         err_msg = os.strerror(errno_num)
-> 1950     raise child_exception_type(errno_num, err_msg, err_filename)
   1951 raise child_exception_type(err_msg)

FileNotFoundError: [Errno 2] No such file or directory: PosixPath('dot')

The above exception was the direct cause of the following exception:

ExecutableNotFound                        Traceback (most recent call last)
File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/IPython/core/formatters.py:977, in MimeBundleFormatter.__call__(self, obj, include, exclude)
    974     method = get_real_method(obj, self.print_method)
    976     if method is not None:
--> 977         return method(include=include, exclude=exclude)
    978     return None
    979 else:

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/decaylanguage/decay/viewer.py:259, in DecayChainViewer._repr_mimebundle_(self, include, exclude, **kwargs)
    255 """
    256 IPython display helper.
    257 """
    258 try:
--> 259     return self._graph._repr_mimebundle_(
    260         include=include, exclude=exclude, **kwargs
    261     )
    262 except AttributeError:
    263     return {"image/svg+xml": self._graph._repr_svg_()}

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/jupyter_integration.py:98, in JupyterIntegration._repr_mimebundle_(self, include, exclude, **_)
     96 include = set(include) if include is not None else {self._jupyter_mimetype}
     97 include -= set(exclude or [])
---> 98 return {mimetype: getattr(self, method_name)()
     99         for mimetype, method_name in MIME_TYPES.items()
    100         if mimetype in include}

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/jupyter_integration.py:98, in <dictcomp>(.0)
     96 include = set(include) if include is not None else {self._jupyter_mimetype}
     97 include -= set(exclude or [])
---> 98 return {mimetype: getattr(self, method_name)()
     99         for mimetype, method_name in MIME_TYPES.items()
    100         if mimetype in include}

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/jupyter_integration.py:112, in JupyterIntegration._repr_image_svg_xml(self)
    110 def _repr_image_svg_xml(self) -> str:
    111     """Return the rendered graph as SVG string."""
--> 112     return self.pipe(format='svg', encoding=SVG_ENCODING)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/piping.py:104, in Pipe.pipe(self, format, renderer, formatter, neato_no_op, quiet, engine, encoding)
     55 def pipe(self,
     56          format: typing.Optional[str] = None,
     57          renderer: typing.Optional[str] = None,
   (...)
     61          engine: typing.Optional[str] = None,
     62          encoding: typing.Optional[str] = None) -> typing.Union[bytes, str]:
     63     """Return the source piped through the Graphviz layout command.
     64 
     65     Args:
   (...)
    102         '<?xml version='
    103     """
--> 104     return self._pipe_legacy(format,
    105                              renderer=renderer,
    106                              formatter=formatter,
    107                              neato_no_op=neato_no_op,
    108                              quiet=quiet,
    109                              engine=engine,
    110                              encoding=encoding)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/_tools.py:171, in deprecate_positional_args.<locals>.decorator.<locals>.wrapper(*args, **kwargs)
    162     wanted = ', '.join(f'{name}={value!r}'
    163                        for name, value in deprecated.items())
    164     warnings.warn(f'The signature of {func.__name__} will be reduced'
    165                   f' to {supported_number} positional args'
    166                   f' {list(supported)}: pass {wanted}'
    167                   ' as keyword arg(s)',
    168                   stacklevel=stacklevel,
    169                   category=category)
--> 171 return func(*args, **kwargs)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/piping.py:121, in Pipe._pipe_legacy(self, format, renderer, formatter, neato_no_op, quiet, engine, encoding)
    112 @_tools.deprecate_positional_args(supported_number=2)
    113 def _pipe_legacy(self,
    114                  format: typing.Optional[str] = None,
   (...)
    119                  engine: typing.Optional[str] = None,
    120                  encoding: typing.Optional[str] = None) -> typing.Union[bytes, str]:
--> 121     return self._pipe_future(format,
    122                              renderer=renderer,
    123                              formatter=formatter,
    124                              neato_no_op=neato_no_op,
    125                              quiet=quiet,
    126                              engine=engine,
    127                              encoding=encoding)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/piping.py:149, in Pipe._pipe_future(self, format, renderer, formatter, neato_no_op, quiet, engine, encoding)
    146 if encoding is not None:
    147     if codecs.lookup(encoding) is codecs.lookup(self.encoding):
    148         # common case: both stdin and stdout need the same encoding
--> 149         return self._pipe_lines_string(*args, encoding=encoding, **kwargs)
    150     try:
    151         raw = self._pipe_lines(*args, input_encoding=self.encoding, **kwargs)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/piping.py:212, in pipe_lines_string(engine, format, input_lines, encoding, renderer, formatter, neato_no_op, quiet)
    206 cmd = dot_command.command(engine, format,
    207                           renderer=renderer,
    208                           formatter=formatter,
    209                           neato_no_op=neato_no_op)
    210 kwargs = {'input_lines': input_lines, 'encoding': encoding}
--> 212 proc = execute.run_check(cmd, capture_output=True, quiet=quiet, **kwargs)
    213 return proc.stdout

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/execute.py:81, in run_check(cmd, input_lines, encoding, quiet, **kwargs)
     79 except OSError as e:
     80     if e.errno == errno.ENOENT:
---> 81         raise ExecutableNotFound(cmd) from e
     82     raise
     84 if not quiet and proc.stderr:

ExecutableNotFound: failed to execute PosixPath('dot'), make sure the Graphviz executables are on your systems' PATH
<decaylanguage.decay.viewer.DecayChainViewer at 0x7f4b7c681bc0>

You can also create a decay using the DecayChain and DecayMode classes. However, a DecayChain can only contain one chain, i.e., a particle cannot decay in multiple ways.

dplus_decay = DecayMode(1, "K- pi+ pi+ pi0", model="PHSP")
pi0_decay = DecayMode(1, "gamma gamma")
dplus_single = DecayChain("D+", {"D+": dplus_decay, "pi0": pi0_decay})
DecayChainViewer(dplus_single.to_dict())
---------------------------------------------------------------------------
FileNotFoundError                         Traceback (most recent call last)
File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/execute.py:76, in run_check(cmd, input_lines, encoding, quiet, **kwargs)
     75         kwargs['stdout'] = kwargs['stderr'] = subprocess.PIPE
---> 76     proc = _run_input_lines(cmd, input_lines, kwargs=kwargs)
     77 else:

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/execute.py:96, in _run_input_lines(cmd, input_lines, kwargs)
     95 def _run_input_lines(cmd, input_lines, *, kwargs):
---> 96     popen = subprocess.Popen(cmd, stdin=subprocess.PIPE, **kwargs)
     98     stdin_write = popen.stdin.write

File ~/.asdf/installs/python/3.11.6/lib/python3.11/subprocess.py:1026, in Popen.__init__(self, args, bufsize, executable, stdin, stdout, stderr, preexec_fn, close_fds, shell, cwd, env, universal_newlines, startupinfo, creationflags, restore_signals, start_new_session, pass_fds, user, group, extra_groups, encoding, errors, text, umask, pipesize, process_group)
   1023             self.stderr = io.TextIOWrapper(self.stderr,
   1024                     encoding=encoding, errors=errors)
-> 1026     self._execute_child(args, executable, preexec_fn, close_fds,
   1027                         pass_fds, cwd, env,
   1028                         startupinfo, creationflags, shell,
   1029                         p2cread, p2cwrite,
   1030                         c2pread, c2pwrite,
   1031                         errread, errwrite,
   1032                         restore_signals,
   1033                         gid, gids, uid, umask,
   1034                         start_new_session, process_group)
   1035 except:
   1036     # Cleanup if the child failed starting.

File ~/.asdf/installs/python/3.11.6/lib/python3.11/subprocess.py:1950, in Popen._execute_child(self, args, executable, preexec_fn, close_fds, pass_fds, cwd, env, startupinfo, creationflags, shell, p2cread, p2cwrite, c2pread, c2pwrite, errread, errwrite, restore_signals, gid, gids, uid, umask, start_new_session, process_group)
   1949         err_msg = os.strerror(errno_num)
-> 1950     raise child_exception_type(errno_num, err_msg, err_filename)
   1951 raise child_exception_type(err_msg)

FileNotFoundError: [Errno 2] No such file or directory: PosixPath('dot')

The above exception was the direct cause of the following exception:

ExecutableNotFound                        Traceback (most recent call last)
File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/IPython/core/formatters.py:977, in MimeBundleFormatter.__call__(self, obj, include, exclude)
    974     method = get_real_method(obj, self.print_method)
    976     if method is not None:
--> 977         return method(include=include, exclude=exclude)
    978     return None
    979 else:

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/decaylanguage/decay/viewer.py:259, in DecayChainViewer._repr_mimebundle_(self, include, exclude, **kwargs)
    255 """
    256 IPython display helper.
    257 """
    258 try:
--> 259     return self._graph._repr_mimebundle_(
    260         include=include, exclude=exclude, **kwargs
    261     )
    262 except AttributeError:
    263     return {"image/svg+xml": self._graph._repr_svg_()}

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/jupyter_integration.py:98, in JupyterIntegration._repr_mimebundle_(self, include, exclude, **_)
     96 include = set(include) if include is not None else {self._jupyter_mimetype}
     97 include -= set(exclude or [])
---> 98 return {mimetype: getattr(self, method_name)()
     99         for mimetype, method_name in MIME_TYPES.items()
    100         if mimetype in include}

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/jupyter_integration.py:98, in <dictcomp>(.0)
     96 include = set(include) if include is not None else {self._jupyter_mimetype}
     97 include -= set(exclude or [])
---> 98 return {mimetype: getattr(self, method_name)()
     99         for mimetype, method_name in MIME_TYPES.items()
    100         if mimetype in include}

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/jupyter_integration.py:112, in JupyterIntegration._repr_image_svg_xml(self)
    110 def _repr_image_svg_xml(self) -> str:
    111     """Return the rendered graph as SVG string."""
--> 112     return self.pipe(format='svg', encoding=SVG_ENCODING)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/piping.py:104, in Pipe.pipe(self, format, renderer, formatter, neato_no_op, quiet, engine, encoding)
     55 def pipe(self,
     56          format: typing.Optional[str] = None,
     57          renderer: typing.Optional[str] = None,
   (...)
     61          engine: typing.Optional[str] = None,
     62          encoding: typing.Optional[str] = None) -> typing.Union[bytes, str]:
     63     """Return the source piped through the Graphviz layout command.
     64 
     65     Args:
   (...)
    102         '<?xml version='
    103     """
--> 104     return self._pipe_legacy(format,
    105                              renderer=renderer,
    106                              formatter=formatter,
    107                              neato_no_op=neato_no_op,
    108                              quiet=quiet,
    109                              engine=engine,
    110                              encoding=encoding)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/_tools.py:171, in deprecate_positional_args.<locals>.decorator.<locals>.wrapper(*args, **kwargs)
    162     wanted = ', '.join(f'{name}={value!r}'
    163                        for name, value in deprecated.items())
    164     warnings.warn(f'The signature of {func.__name__} will be reduced'
    165                   f' to {supported_number} positional args'
    166                   f' {list(supported)}: pass {wanted}'
    167                   ' as keyword arg(s)',
    168                   stacklevel=stacklevel,
    169                   category=category)
--> 171 return func(*args, **kwargs)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/piping.py:121, in Pipe._pipe_legacy(self, format, renderer, formatter, neato_no_op, quiet, engine, encoding)
    112 @_tools.deprecate_positional_args(supported_number=2)
    113 def _pipe_legacy(self,
    114                  format: typing.Optional[str] = None,
   (...)
    119                  engine: typing.Optional[str] = None,
    120                  encoding: typing.Optional[str] = None) -> typing.Union[bytes, str]:
--> 121     return self._pipe_future(format,
    122                              renderer=renderer,
    123                              formatter=formatter,
    124                              neato_no_op=neato_no_op,
    125                              quiet=quiet,
    126                              engine=engine,
    127                              encoding=encoding)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/piping.py:149, in Pipe._pipe_future(self, format, renderer, formatter, neato_no_op, quiet, engine, encoding)
    146 if encoding is not None:
    147     if codecs.lookup(encoding) is codecs.lookup(self.encoding):
    148         # common case: both stdin and stdout need the same encoding
--> 149         return self._pipe_lines_string(*args, encoding=encoding, **kwargs)
    150     try:
    151         raw = self._pipe_lines(*args, input_encoding=self.encoding, **kwargs)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/piping.py:212, in pipe_lines_string(engine, format, input_lines, encoding, renderer, formatter, neato_no_op, quiet)
    206 cmd = dot_command.command(engine, format,
    207                           renderer=renderer,
    208                           formatter=formatter,
    209                           neato_no_op=neato_no_op)
    210 kwargs = {'input_lines': input_lines, 'encoding': encoding}
--> 212 proc = execute.run_check(cmd, capture_output=True, quiet=quiet, **kwargs)
    213 return proc.stdout

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/execute.py:81, in run_check(cmd, input_lines, encoding, quiet, **kwargs)
     79 except OSError as e:
     80     if e.errno == errno.ENOENT:
---> 81         raise ExecutableNotFound(cmd) from e
     82     raise
     84 if not quiet and proc.stderr:

ExecutableNotFound: failed to execute PosixPath('dot'), make sure the Graphviz executables are on your systems' PATH
<decaylanguage.decay.viewer.DecayChainViewer at 0x7f4b7c1e2100>
Creating a GenMultiDecay object

A regular phasespace.GenParticle instance would not be able to simulate this decay, since the \(\pi^0\) particle can decay in four different ways. However, a GenMultiDecay object can be created directly from a DecayLanguage dict:

pi0_decay = GenMultiDecay.from_dict(pi0_chain)

When creating a GenMultiDecay object, the DecayLanguage dict is “unpacked” into separate GenParticle instances, where each GenParticle instance corresponds to one way that the particle can decay.

These GenParticle instances and the probabilities of that decay mode can be accessed via GenMultiDecay.gen_particles. This is a list of tuples, where the first element in the tuple is the probability and the second element is the GenParticle.

for probability, particle in pi0_decay.gen_particles:
    print(
        f"There is a probability of {probability} "
        f"that pi0 decays into {', '.join(child.name for child in particle.children)}"
    )
There is a probability of 0.988228297 that pi0 decays into gamma, gamma [0]

There is a probability of 0.011738247 that pi0 decays into e+, e-, gamma [1]

There is a probability of 3.3392e-05 that pi0 decays into e+ [0], e+ [1], e- [0], e- [1]

There is a probability of 6.5e-08 that pi0 decays into e+ [2], e- [2]

One can simulate this decay using the .generate method, which works the same as the GenParticle.generate method.

When calling the GenMultiDecay.generate method, it internally calls the generate method on the of the GenParticle instances in GenMultiDecay.gen_particles. The outputs are placed in a list, which is returned.

weights, events = pi0_decay.generate(n_events=10_000)
print("Number of events for each decay mode:", ", ".join(str(len(w)) for w in weights))
WARNING: All log messages before absl::InitializeLog() is called are written to STDERR
I0000 00:00:1713898691.162394     759 service.cc:145] XLA service 0x7f4b640a15f0 initialized for platform Host (this does not guarantee that XLA will be used). Devices:
I0000 00:00:1713898691.162592     759 service.cc:153]   StreamExecutor device (0): Host, Default Version
I0000 00:00:1713898691.187431     758 device_compiler.h:188] Compiled cluster using XLA!  This line is logged at most once for the lifetime of the process.
Number of events for each decay mode:
 
9884, 116

We can confirm that the counts above are close to the expected counts based on the probabilities.

Changing mass settings

Since DecayLanguage dicts do not contain any information about the mass of a particle, the fromdecay submodule uses the particle package to find the mass of a particle based on its name. The mass can either be a constant value or a function (besides the top particle, which is always a constant). These settings can be modified by passing in additional parameters to GenMultiDecay.from_dict. There are two optional parameters that can be passed to GenMultiDecay.from_dict: tolerance and mass_converter.

Constant vs variable mass

If a particle has a width less than tolerance, its mass is set to a constant value. This will be demonsttrated with the decay below:

dsplus_chain = parser.build_decay_chains("D*+", stable_particles=["D+"])
DecayChainViewer(dsplus_chain)
---------------------------------------------------------------------------
FileNotFoundError                         Traceback (most recent call last)
File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/execute.py:76, in run_check(cmd, input_lines, encoding, quiet, **kwargs)
     75         kwargs['stdout'] = kwargs['stderr'] = subprocess.PIPE
---> 76     proc = _run_input_lines(cmd, input_lines, kwargs=kwargs)
     77 else:

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/execute.py:96, in _run_input_lines(cmd, input_lines, kwargs)
     95 def _run_input_lines(cmd, input_lines, *, kwargs):
---> 96     popen = subprocess.Popen(cmd, stdin=subprocess.PIPE, **kwargs)
     98     stdin_write = popen.stdin.write

File ~/.asdf/installs/python/3.11.6/lib/python3.11/subprocess.py:1026, in Popen.__init__(self, args, bufsize, executable, stdin, stdout, stderr, preexec_fn, close_fds, shell, cwd, env, universal_newlines, startupinfo, creationflags, restore_signals, start_new_session, pass_fds, user, group, extra_groups, encoding, errors, text, umask, pipesize, process_group)
   1023             self.stderr = io.TextIOWrapper(self.stderr,
   1024                     encoding=encoding, errors=errors)
-> 1026     self._execute_child(args, executable, preexec_fn, close_fds,
   1027                         pass_fds, cwd, env,
   1028                         startupinfo, creationflags, shell,
   1029                         p2cread, p2cwrite,
   1030                         c2pread, c2pwrite,
   1031                         errread, errwrite,
   1032                         restore_signals,
   1033                         gid, gids, uid, umask,
   1034                         start_new_session, process_group)
   1035 except:
   1036     # Cleanup if the child failed starting.

File ~/.asdf/installs/python/3.11.6/lib/python3.11/subprocess.py:1950, in Popen._execute_child(self, args, executable, preexec_fn, close_fds, pass_fds, cwd, env, startupinfo, creationflags, shell, p2cread, p2cwrite, c2pread, c2pwrite, errread, errwrite, restore_signals, gid, gids, uid, umask, start_new_session, process_group)
   1949         err_msg = os.strerror(errno_num)
-> 1950     raise child_exception_type(errno_num, err_msg, err_filename)
   1951 raise child_exception_type(err_msg)

FileNotFoundError: [Errno 2] No such file or directory: PosixPath('dot')

The above exception was the direct cause of the following exception:

ExecutableNotFound                        Traceback (most recent call last)
File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/IPython/core/formatters.py:977, in MimeBundleFormatter.__call__(self, obj, include, exclude)
    974     method = get_real_method(obj, self.print_method)
    976     if method is not None:
--> 977         return method(include=include, exclude=exclude)
    978     return None
    979 else:

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/decaylanguage/decay/viewer.py:259, in DecayChainViewer._repr_mimebundle_(self, include, exclude, **kwargs)
    255 """
    256 IPython display helper.
    257 """
    258 try:
--> 259     return self._graph._repr_mimebundle_(
    260         include=include, exclude=exclude, **kwargs
    261     )
    262 except AttributeError:
    263     return {"image/svg+xml": self._graph._repr_svg_()}

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/jupyter_integration.py:98, in JupyterIntegration._repr_mimebundle_(self, include, exclude, **_)
     96 include = set(include) if include is not None else {self._jupyter_mimetype}
     97 include -= set(exclude or [])
---> 98 return {mimetype: getattr(self, method_name)()
     99         for mimetype, method_name in MIME_TYPES.items()
    100         if mimetype in include}

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/jupyter_integration.py:98, in <dictcomp>(.0)
     96 include = set(include) if include is not None else {self._jupyter_mimetype}
     97 include -= set(exclude or [])
---> 98 return {mimetype: getattr(self, method_name)()
     99         for mimetype, method_name in MIME_TYPES.items()
    100         if mimetype in include}

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/jupyter_integration.py:112, in JupyterIntegration._repr_image_svg_xml(self)
    110 def _repr_image_svg_xml(self) -> str:
    111     """Return the rendered graph as SVG string."""
--> 112     return self.pipe(format='svg', encoding=SVG_ENCODING)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/piping.py:104, in Pipe.pipe(self, format, renderer, formatter, neato_no_op, quiet, engine, encoding)
     55 def pipe(self,
     56          format: typing.Optional[str] = None,
     57          renderer: typing.Optional[str] = None,
   (...)
     61          engine: typing.Optional[str] = None,
     62          encoding: typing.Optional[str] = None) -> typing.Union[bytes, str]:
     63     """Return the source piped through the Graphviz layout command.
     64 
     65     Args:
   (...)
    102         '<?xml version='
    103     """
--> 104     return self._pipe_legacy(format,
    105                              renderer=renderer,
    106                              formatter=formatter,
    107                              neato_no_op=neato_no_op,
    108                              quiet=quiet,
    109                              engine=engine,
    110                              encoding=encoding)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/_tools.py:171, in deprecate_positional_args.<locals>.decorator.<locals>.wrapper(*args, **kwargs)
    162     wanted = ', '.join(f'{name}={value!r}'
    163                        for name, value in deprecated.items())
    164     warnings.warn(f'The signature of {func.__name__} will be reduced'
    165                   f' to {supported_number} positional args'
    166                   f' {list(supported)}: pass {wanted}'
    167                   ' as keyword arg(s)',
    168                   stacklevel=stacklevel,
    169                   category=category)
--> 171 return func(*args, **kwargs)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/piping.py:121, in Pipe._pipe_legacy(self, format, renderer, formatter, neato_no_op, quiet, engine, encoding)
    112 @_tools.deprecate_positional_args(supported_number=2)
    113 def _pipe_legacy(self,
    114                  format: typing.Optional[str] = None,
   (...)
    119                  engine: typing.Optional[str] = None,
    120                  encoding: typing.Optional[str] = None) -> typing.Union[bytes, str]:
--> 121     return self._pipe_future(format,
    122                              renderer=renderer,
    123                              formatter=formatter,
    124                              neato_no_op=neato_no_op,
    125                              quiet=quiet,
    126                              engine=engine,
    127                              encoding=encoding)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/piping.py:149, in Pipe._pipe_future(self, format, renderer, formatter, neato_no_op, quiet, engine, encoding)
    146 if encoding is not None:
    147     if codecs.lookup(encoding) is codecs.lookup(self.encoding):
    148         # common case: both stdin and stdout need the same encoding
--> 149         return self._pipe_lines_string(*args, encoding=encoding, **kwargs)
    150     try:
    151         raw = self._pipe_lines(*args, input_encoding=self.encoding, **kwargs)

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/piping.py:212, in pipe_lines_string(engine, format, input_lines, encoding, renderer, formatter, neato_no_op, quiet)
    206 cmd = dot_command.command(engine, format,
    207                           renderer=renderer,
    208                           formatter=formatter,
    209                           neato_no_op=neato_no_op)
    210 kwargs = {'input_lines': input_lines, 'encoding': encoding}
--> 212 proc = execute.run_check(cmd, capture_output=True, quiet=quiet, **kwargs)
    213 return proc.stdout

File ~/checkouts/readthedocs.org/user_builds/phasespace/envs/stable/lib/python3.11/site-packages/graphviz/backend/execute.py:81, in run_check(cmd, input_lines, encoding, quiet, **kwargs)
     79 except OSError as e:
     80     if e.errno == errno.ENOENT:
---> 81         raise ExecutableNotFound(cmd) from e
     82     raise
     84 if not quiet and proc.stderr:

ExecutableNotFound: failed to execute PosixPath('dot'), make sure the Graphviz executables are on your systems' PATH
<decaylanguage.decay.viewer.DecayChainViewer at 0x7f4b7c3c63c0>
print(
    f"pi0 width = {Particle.from_evtgen_name('pi0').width}\n"
    f"D0 width = {Particle.from_evtgen_name('D0').width}"
)
pi0 width = 7.81e-06
D0 width = 1.604e-09

\(\pi^0\) has a greater width than \(D^0\). If the tolerance is set to a value between their widths, the \(D^0\) particle will have a constant mass while \(\pi^0\) will not.

dstar_decay = GenMultiDecay.from_dict(dsplus_chain, tolerance=1e-8)
# Loop over D0 and pi+ particles, see graph above
for particle in dstar_decay.gen_particles[0][1].children:
    # If a particle width is less than tolerance or if it does not have any children, its mass will be fixed.
    assert particle.has_fixed_mass

# Loop over D+ and pi0. See above.
for particle in dstar_decay.gen_particles[1][1].children:
    if particle.name == "pi0":
        assert not particle.has_fixed_mass
Configuring mass functions

By default, the mass function used for variable mass is the relativistic Breit-Wigner distribution. This can however be changed. If you want the mother particle to have a specific mass function for a specific decay, you can add a zfit parameter to the DecayLanguage dict. Consider for example the previous \(D^{*+}\) example:

dsplus_custom_mass_func = dsplus_chain.copy()
dsplus_chain_subset = dsplus_custom_mass_func["D*+"][1]["fs"][1]
print("Before:")
pprint(dsplus_chain_subset)
# Set the mass function of pi0 to a gaussian distribution when it decays into two photons (gamma)
dsplus_chain_subset["pi0"][0]["zfit"] = "gauss"
print("After:")
pprint(dsplus_chain_subset)
Before:

{
'pi0'
: 
[
{
'bf'
: 
0.988228297
,
          
'fs'
: 
['gamma', 'gamma']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
0.011738247
,
          
'fs'
: 
['e+', 'e-', 'gamma']
,
          
'model'
: 
'PI0_DALITZ'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
3.3392e-05
,
          
'fs'
: 
['e+', 'e+', 'e-', 'e-']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
6.5e-08
,
          
'fs'
: 
['e+', 'e-']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
]
}

After:

{
'pi0'
: 
[
{
'bf'
: 
0.988228297
,
          
'fs'
: 
['gamma', 'gamma']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
,
          
'zfit'
: 
'gauss'
}
,
         
{
'bf'
: 
0.011738247
,
          
'fs'
: 
['e+', 'e-', 'gamma']
,
          
'model'
: 
'PI0_DALITZ'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
3.3392e-05
,
          
'fs'
: 
['e+', 'e+', 'e-', 'e-']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
6.5e-08
,
          
'fs'
: 
['e+', 'e-']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
]
}

Notice the added zfit field to the first decay mode of the \(\pi^0\) particle. This dict can then be passed to GenMultiDecay.from_dict, like before.

GenMultiDecay.from_dict(dsplus_custom_mass_func)
<phasespace.fromdecay.genmultidecay.GenMultiDecay at 0x7f4b536fe710>

If you want all \(\pi^0\) particles to decay with the same mass function, you do not need to specify the zfit parameter for each decay in the dict. Instead, one can pass the particle_model_map parameter to the constructor:

GenMultiDecay.from_dict(
    dsplus_chain, particle_model_map={"pi0": "gauss"}
)  # pi0 always decays with a gaussian mass distribution.
<phasespace.fromdecay.genmultidecay.GenMultiDecay at 0x7f4b5376de90>

When using DecayChains, the syntax for specifying the mass function becomes cleaner:

dplus_decay = DecayMode(
    1, "K- pi+ pi+ pi0", model="PHSP"
)  # The model parameter will be ignored by GenMultiDecay
pi0_decay = DecayMode(
    1, "gamma gamma", zfit="gauss"
)  # Make pi0 have a gaussian mass distribution
dplus_single = DecayChain("D+", {"D+": dplus_decay, "pi0": pi0_decay})
GenMultiDecay.from_dict(dplus_single.to_dict())
<phasespace.fromdecay.genmultidecay.GenMultiDecay at 0x7f4b79401950>
Custom mass functions

The built-in supported mass function names are gauss, bw, and relbw, with gauss being the gaussian distribution, bw being the Breit-Wigner distribution, and relbw being the relativistic Breit-Wigner distribution.

If a non-supported value for the zfit parameter is not specified, it will automatically use the relativistic Breit-Wigner distribution. This behavior can be changed by changing the value of GenMultiDecay.DEFAULT_MASS_FUNC to a different string, e.g., "gauss". If an invalid value for the zfit parameter is used, a KeyError is raised.

It is also possible to add your own mass functions besides the built-in ones. You should then create a function that takes the mass and width of a particle and returns a mass function which with the format that is used for all phasespace mass functions. Below is an example of a custom gaussian distribution (implemented in the same way as the built-in gaussian distribution), which uses zfit PDFs:

def custom_gauss(mass, width):
    particle_mass = tf.cast(mass, tf.float64)
    particle_width = tf.cast(width, tf.float64)

    # This is the actual mass function that will be returned
    def mass_func(min_mass, max_mass, n_events):
        min_mass = tf.cast(min_mass, tf.float64)
        max_mass = tf.cast(max_mass, tf.float64)
        # Use a zfit PDF
        pdf = zfit.pdf.Gauss(mu=particle_mass, sigma=particle_width, obs="")
        iterator = tf.stack([min_mass, max_mass], axis=-1)
        return tf.vectorized_map(
            lambda lim: pdf.sample(1, limits=(lim[0], lim[1])), iterator
        )

    return mass_func

This function can then be passed to GenMultiDecay.from_dict as a dict, where the key specifies the zfit parameter name. In the example below, it is set to "custom_gauss". However, this name can be chosen arbitrarily and does not need to be the same as the function name.

dsplus_chain_subset = dsplus_custom_mass_func["D*+"][1]["fs"][1]
print("Before:")
pprint(dsplus_chain_subset)

# Set the mass function of pi0 to the custom gaussian distribution
#  when it decays into an electron-positron pair and a photon (gamma)
dsplus_chain_subset["pi0"][1]["zfit"] = "custom_gauss"
print("After:")
pprint(dsplus_chain_subset)
Before:

{
'pi0'
: 
[
{
'bf'
: 
0.988228297
,
          
'fs'
: 
['gamma', 'gamma']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
,
          
'zfit'
: 
'gauss'
}
,
         
{
'bf'
: 
0.011738247
,
          
'fs'
: 
['e+', 'e-', 'gamma']
,
          
'model'
: 
'PI0_DALITZ'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
3.3392e-05
,
          
'fs'
: 
['e+', 'e+', 'e-', 'e-']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
6.5e-08
,
          
'fs'
: 
['e+', 'e-']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
]
}

After:

{
'pi0'
: 
[
{
'bf'
: 
0.988228297
,
          
'fs'
: 
['gamma', 'gamma']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
,
          
'zfit'
: 
'gauss'
}
,
         
{
'bf'
: 
0.011738247
,
          
'fs'
: 
['e+', 'e-', 'gamma']
,
          
'model'
: 
'PI0_DALITZ'
,
          
'model_params'
: 
''
,
          
'zfit'
: 
'custom_gauss'
}
,
         
{
'bf'
: 
3.3392e-05
,
          
'fs'
: 
['e+', 'e+', 'e-', 'e-']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
,
         
{
'bf'
: 
6.5e-08
,
          
'fs'
: 
['e+', 'e-']
,
          
'model'
: 
'PHSP'
,
          
'model_params'
: 
''
}
]
}

GenMultiDecay.from_dict(dsplus_custom_mass_func, {"custom_gauss": custom_gauss})
<phasespace.fromdecay.genmultidecay.GenMultiDecay at 0x7f4b537dfd10>

phasespace package

phasespace.phasespace module

Implementation of the Raubold and Lynch method to generate n-body events.

The code is based on the GENBOD function (W515 from CERNLIB), documented in:

  1. James, Monte Carlo Phase Space, CERN 68-15 (1968)

class phasespace.phasespace.GenParticle(name: str, mass: Callable | int | float | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes])[source]

Bases: object

Representation of a particle.

Instances of this class can be combined with each other to build decay chains, which can then be used to generate phase space events through the generate or generate_tensor method.

A GenParticle must have
  • a name, which is ensured not to clash with any others in

    the decay chain.

  • a mass, which can be either a number or a function to generate it according to

    a certain distribution. The returned ~`tf.Tensor` needs to have shape (nevents,). In this case, the particle is not considered as having a fixed mass and the has_fixed_mass method will return False.

It may also have:

  • Children, ie, decay products, which are also GenParticle instances.

Parameters:
  • name (str) – Name of the particle.

  • mass (float, array-like, callable) – Mass of the particle. If it’s a float, it get converted to array-like.

generate(n_events: int | tf.Tensor | tf.Variable, boost_to: tf.Tensor | vector.Momentum | None = None, normalize_weights: bool = True, seed: SeedLike = None, *, as_vectors: bool | None = None) tuple[tf.Tensor, dict[str, tf.Tensor]][source]

Generate normalized n-body phase space as tensorflow tensors.

Any TensorFlow tensor can always be converted to a numpy array with the method numpy().

Events are generated in the rest frame of the particle, unless boost_to is given.

Note

In this method, the event weights are returned normalized to their maximum.

Parameters:
  • n_events (int) – Number of events to generate.

  • boost_to (optional) – Momentum vector of shape (x, 4), where x is optional, to where the resulting events will be boosted in the (px, py, pz, E) format. Can also be a vector momentum Lorentz vector. If not specified, events are generated in the rest frame of the particle.

  • normalize_weights (bool, optional) – Normalize the event weight to its max?

  • seed (SeedLike) – The seed can be a number or a tf.random.Generator that are used as a seed to create a random number generator inside the function or directly as the random number generator instance, respectively.

  • as_vectors (bool, optional) – If True, the output momenta are returned as vector objects.

Returns:

Result of the generation, which varies with the value of normalize_weights:

  • If True, the tuple elements are the normalized event weights as a tensor of shape (n_events, ), and the momenta generated particles as a dictionary of tensors of shape (4, n_events) with particle names as keys.

  • If False, the tuple weights are the unnormalized event weights as a tensor of shape (n_events, ), the maximum per-event weights as a tensor of shape (n_events, ) and the momenta generated particles as a dictionary of tensors of shape (4, n_events) with particle names as keys.

Return type:

tuple

Raises:
  • tf.errors.InvalidArgumentError – If the the decay is kinematically forbidden.

  • ValueError – If n_events and the size of boost_to don’t match.

  • See GenParticle.generate_unnormalized.

generate_tensor(n_events: int, boost_to=None, normalize_weights: bool = True)[source]

Generate normalized n-body phase space as numpy arrays.

Events are generated in the rest frame of the particle, unless boost_to is given.

Note

In this method, the event weights are returned normalized to their maximum.

Parameters:
  • n_events (int) – Number of events to generate.

  • boost_to (optional) – Momentum vector of shape (x, 4), where x is optional, to where the resulting events will be boosted. If not specified, events are generated in the rest frame of the particle.

  • normalize_weights (bool, optional) – Normalize the event weight to its max

Returns:

Result of the generation, which varies with the value of normalize_weights:

  • If True, the tuple elements are the normalized event weights as an array of shape (n_events, ), and the momenta generated particles as a dictionary of arrays of shape (4, n_events) with particle names as keys.

  • If False, the tuple weights are the unnormalized event weights as an array of shape (n_events, ), the maximum per-event weights as an array of shape (n_events, ) and the momenta generated particles as a dictionary of arrays of shape (4, n_events) with particle names as keys.

Return type:

tuple

Raises:
  • tf.errors.InvalidArgumentError – If the the decay is kinematically forbidden.

  • ValueError – If n_events and the size of boost_to don’t match.

  • See GenParticle.generate_unnormalized.

get_mass(min_mass: tf.Tensor = None, max_mass: tf.Tensor = None, n_events: tf.Tensor | tf.Variable = None, seed: SeedLike = None) tf.Tensor[source]

Get the particle mass.

If the particle is resonant, the mass function will be called with the min_mass, max_mass, n_events and optionally a seed parameter.

Parameters:
  • min_mass (tensor) – Lower mass range. Defaults to None, which is only valid in the case of fixed mass.

  • max_mass (tensor) – Upper mass range. Defaults to None, which is only valid in the case of fixed mass.

  • () (n_events) – Number of events to produce. Has to be specified if the particle is resonant.

  • seed (SeedLike) – The seed can be a number or a tf.random.Generator that are used as a seed to create a random number generator inside the function or directly as the random number generator instance, respectively.

Returns:

Mass of the particles, either a scalar or shape (nevents,)

Return type:

Tensor`

Raises:

ValueError – If the mass is requested and has not been set.

property has_children

Does the particle have children?

Type:

bool

property has_fixed_mass

Is the mass a callable function?

Type:

bool

property has_grandchildren

Does the particle have grandchildren?

Type:

bool

set_children(*children)[source]

Assign children.

Parameters:

children (GenParticle) – Two or more children to assign to the current particle.

Returns:

self

Raises:
  • ValueError – If there is an inconsistency in the parent/children relationship, ie,

  • if children were already set, if their parent was or if less than two children were given.

  • KeyError – If there is a particle name clash.

  • RuntimeError – If generate was already called before.

class phasespace.phasespace.Particle[source]

Bases: object

Deprecated Particle class.

Renamed to GenParticle.

phasespace.phasespace.generate_decay(*args, **kwargs)[source]

Deprecated.

phasespace.phasespace.nbody_decay(mass_top: float, masses: list, top_name: str = '', names: list = None)[source]

Shortcut to build an n-body decay of a GenParticle.

If the particle names are not given, the top particle is called ‘top’ and the children ‘p_{i}’, where i corresponds to their position in the masses sequence.

Parameters:
  • mass_top (tensor, list) – Mass of the top particle. Can be a list of 4-vectors.

  • masses (list) – Masses of the child particles.

  • name_top (str, optional) – Name of the top particle. If not given, the top particle is named top.

  • names (list, optional) – Names of the child particles. If not given, they are build as ‘p_{i}’, where i is given by their ordering in the masses list.

Returns:

Particle decay.

Return type:

GenParticle

Raises:

ValueError – If the length of masses and names doesn’t match.

phasespace.phasespace.pdk(a, b, c)[source]

Calculate the PDK (2-body phase space) function.

Based on Eq. (9.17) in CERN 68-15 (1968).

Parameters:
  • a (Tensor`) – \(M_{i+1}\) in Eq. (9.17).

  • b (Tensor`) – \(M_{i}\) in Eq. (9.17).

  • c (Tensor`) – \(m_{i+1}\) in Eq. (9.17).

Returns:

~`tf.Tensor`

phasespace.phasespace.process_list_to_tensor(lst)[source]

Convert a list to a tensor.

The list is converted to a tensor and transposed to get the proper shape.

Note

If lst is a tensor, nothing is done to it other than convert it to tf.float64.

Parameters:

lst (list) – List to convert.

Returns:

~`tf.Tensor`

phasespace.phasespace.to_vectors(particles: dict[str, tf.Tensor]) dict[str, vector.Momentum][source]

Convert a dictionary of particles to a dictionary of vector.Momentum instances.

Parameters:

particles (dict) – Dictionary of particles, with the keys being the particle names and the values being the momenta.

Returns:

Dictionary of vector.Momentum instances with numpy arrays

Return type:

dict

phasespace.kinematics module

Basic kinematics.

phasespace.kinematics.beta(vector)[source]

Calculate beta of a given 4-vector.

Parameters:

vector – Input Lorentz momentum vector.

phasespace.kinematics.boost_components(vector)[source]

Get the boost components of a given 4-vector.

Parameters:

vector – Input Lorentz momentum vector.

phasespace.kinematics.lorentz_boost(vector, boostvector)[source]

Perform Lorentz boost.

Parameters:
  • vector – 4-vector to be boosted

  • boostvector – Boost vector. Can be either 3-vector or 4-vector, since only spatial components are used.

phasespace.kinematics.lorentz_vector(space, time)[source]

Make a Lorentz vector from spatial and time components.

Parameters:
  • space – 3-vector of spatial components.

  • time – Time component.

phasespace.kinematics.mass(vector)[source]

Calculate mass scalar for Lorentz 4-momentum.

Parameters:

vector – Input Lorentz momentum vector.

phasespace.kinematics.metric_tensor()[source]

Metric tensor for Lorentz space (constant).

phasespace.kinematics.scalar_product(vec1, vec2)[source]

Calculate scalar product of two 3-vectors.

Parameters:
  • vec1 – First vector.

  • vec2 – Second vector.

phasespace.kinematics.spatial_component(vector)[source]

Extract spatial components of the input Lorentz vector.

Parameters:

vector – Input Lorentz vector (where indexes 0-2 are space, index 3 is time).

phasespace.kinematics.time_component(vector)[source]

Extract time component of the input Lorentz vector.

Parameters:

vector – Input Lorentz vector (where indexes 0-2 are space, index 3 is time).

phasespace.kinematics.x_component(vector)[source]

Extract spatial X component of the input Lorentz or 3-vector.

Parameters:

vector – Input vector.

phasespace.kinematics.y_component(vector)[source]

Extract spatial Y component of the input Lorentz or 3-vector.

Parameters:

vector – Input vector.

phasespace.kinematics.z_component(vector)[source]

Extract spatial Z component of the input Lorentz or 3-vector.

Parameters:

vector – Input vector.

phasespace.fromdecay module
class phasespace.fromdecay.genmultidecay.GenMultiDecay(gen_particles: list[tuple[float, GenParticle]])[source]

Bases: object

A GenParticle-type container that can handle multiple decays.

Parameters:

gen_particles – All the GenParticles and their corresponding probabilities. The list must be of the format [[probability, GenParticle instance], [probability, …

DEFAULT_MASS_FUNC = 'relbw'
MASS_WIDTH_TOLERANCE = 0.01
classmethod from_dict(dec_dict: dict, mass_converter: dict[str, Callable] = None, tolerance: float = None, particle_model_map: dict[str, str] = None)[source]

Create a GenMultiDecay instance from a dict in the DecayLanguage package format.

Parameters:
  • dec_dict – The input dict from which the GenMultiDecay object will be created from. A dec_dict has the same structure as the dicts used in DecayLanguage, see the examples below.

  • mass_converter – A dict with mass function names and their corresponding mass functions. These functions should take the particle mass and the mass width as inputs and return a mass function that phasespace can understand. This dict will be combined with the predefined mass functions in this package. See the Example below or the tutorial for how to use this parameter.

  • tolerance – Minimum mass width of the particle to use a mass function instead of assuming the mass to be constant. If None, the default value, defined by the class variable MASS_WIDTH_TOLERANCE, will be used. This value can be customized if desired.

  • particle_model_map – A dict where the key is a particle name and the value is a mass function name. If a particle is specified in the particle_model_map, then all appearances of this particle in dec_dict will get the same mass function. This way, one does not have to manually add the zfit parameter to every place where this particle appears in dec_dict. If the zfit parameter is specified for a particle which is also included in particle_model_map, the zfit parameter mass function name will be prioritized.

Returns:

The created GenMultiDecay object.

Examples

DecayLanguage is usually used to create a dict that can be understood by GenMultiDecay. >>> from decaylanguage import DecFileParser >>> from phasespace.fromdecay import GenMultiDecay

Parse a .dec file to create a DecayLanguage dict describing a D*+ particle that can decay in 2 different ways: D*+ -> D0(->K- pi+) pi+ or D*+ -> D+ gamma.

>>> parser = DecFileParser('example_decays.dec')
>>> parser.parse()
>>> dst_chain = parser.build_decay_chains("D*+")
>>> dst_chain
{'D*+': [{'bf': 0.984,
    'fs': [{'D0': [{'bf': 1.0,
            'fs': ['K-', 'pi+']}]},
        'pi+']},
    {'bf': 0.016,
     'fs': ['D+', 'gamma']}]}

If the D0 particle should have a mass distribution of a gaussian when it decays, one can pass the particle_model_map parameter to from_dict: >>> dst_gen = GenMultiDecay.from_dict(dst_chain, particle_model_map={“D0”: “gauss”}) This will then set the mass function of D0 to a gaussian for all its decays.

If more custom control is required, e.g., if D0 can decay in multiple ways and one of the decays should have a specific mass function, a zfit parameter can be added to its decay dict: >>> dst_chain[“D*+”][0][“fs”][0][“D0”][0][“zfit”] = “gauss” >>> dst_chain {‘D*+’: [{‘bf’: 0.984,

‘fs’: [{‘D0’: [{‘bf’: 1.0,

‘fs’: [‘K-’, ‘pi+’], ‘zfit’: ‘gauss’}]},

‘pi+’]},

{‘bf’: 0.016, ‘fs’: [‘D+’, ‘gamma’]}]}

This dict can then be passed to GenMultiDecay.from_dict: >>> dst_gen = GenMultiDecay.from_dict(dst_chain) This will now convert make the D0 particle have a gaussian mass function, only when it decays into K- and pi+. In this case, there are no other ways that D0 can decay, so using particle_model_map is a cleaner and easier option.

If the decay of the D0 particle should be modelled by a mass distribution that does not come with the package, a custom distribution can be created: >>> def custom_gauss(mass, width): >>> particle_mass = tf.cast(mass, tf.float64) >>> particle_width = tf.cast(width, tf.float64) >>> def mass_func(min_mass, max_mass, n_events): >>> min_mass = tf.cast(min_mass, tf.float64) >>> max_mass = tf.cast(max_mass, tf.float64) >>> # Use a zfit PDF >>> pdf = zfit.pdf.Gauss(mu=particle_mass, sigma=particle_width, obs=””) >>> iterator = tf.stack([min_mass, max_mass], axis=-1) >>> return tf.vectorized_map( >>> lambda lim: pdf.sample(1, limits=(lim[0], lim[1])), iterator >>> ) >>> return mass_func

Once again change the distribution in the dst_chain dict. Here, it is changed to “custom_gauss” but any name can be used. >>> dst_chain[“D*+”][0][“fs”][0][“D0”][0][“zfit”] = “custom_gauss”

One can then pass the custom_gauss function and its name (in this case “custom_gauss”) as a dict`to `from_dict as the mass_converter parameter: >>> dst_gen = GenMultiDecay.from_dict(dst_chain, mass_converter={“custom_gauss”: custom_gauss})

Notes

For a more in-depth tutorial, see the tutorial on GenMultiDecay in the documentation.

generate(n_events: int, normalize_weights: bool = True, **kwargs) tuple[list[Tensor], list[Tensor]] | tuple[list[Tensor], list[Tensor], list[Tensor]][source]

Generate four-momentum vectors from the decay(s).

Parameters:
  • n_events – Total number of events combined, for all the decays.

  • normalize_weights – Normalize weights according to all events generated. This also changes the return values. See the phasespace documentation for more details.

  • kwargs – Additional parameters passed to all calls of GenParticle.generate

Returns:

The arguments returned by GenParticle.generate are returned. See the phasespace documentation for details. However, instead of being 2 or 3 tensors, it is 2 or 3 lists of tensors, each entry in the lists corresponding to the return arguments from the corresponding GenParticle instances in self.gen_particles. Note that when normalize_weights is True, the weights are normalized to the maximum of all returned events.

phasespace.fromdecay.mass_functions.breitwigner_factory(mass, width)[source]
phasespace.fromdecay.mass_functions.gauss_factory(mass, width)[source]
phasespace.fromdecay.mass_functions.relativistic_breitwigner_factory(mass, width)[source]

Contributing

Contributions are welcome, and they are greatly appreciated! Every little bit helps, and credit will always be given.

You can contribute in many ways:

Types of Contributions
Report Bugs

Report bugs at https://github.com/zfit/phasespace/issues.

If you are reporting a bug, please include:

  • Your operating system name and version.

  • Any details about your local setup that might be helpful in troubleshooting.

  • Detailed steps to reproduce the bug.

Fix Bugs

Look through the GitHub issues for bugs. Anything tagged with “bug” and “help wanted” is open to whoever wants to implement it.

Implement Features

Look through the GitHub issues for features. Anything tagged with “enhancement” and “help wanted” is open to whoever wants to implement it.

Write Documentation

TensorFlow PhaseSpace could always use more documentation, whether as part of the official TensorFlow PhaseSpace docs, in docstrings, or even on the web in blog posts, articles, and such.

Submit Feedback

The best way to send feedback is to file an issue at https://github.com/zfit/phasespace/issues.

If you are proposing a feature:

  • Explain in detail how it would work.

  • Keep the scope as narrow as possible, to make it easier to implement.

  • Remember that this is a volunteer-driven project, and that contributions are welcome :)

Get Started!

Ready to contribute? Here’s how to set up phasespace for local development.

  1. Fork the phasespace repo on GitHub.

  2. Clone your fork locally:

    $ git clone git@github.com:your_name_here/phasespace.git
    
  3. Install your local copy into a virtualenv. Assuming you have virtualenvwrapper installed, this is how you set up your fork for local development:

    $ mkvirtualenv phasespace
    $ cd phasespace/
    $ pip install -e .[dev]
    
  4. Create a branch for local development:

    $ git checkout -b name-of-your-bugfix-or-feature
    

    Now you can make your changes locally.

  5. When you’re done making changes, check that your changes pass pre-commit and the tests:

    $ pytest $ pre-commit run -a

  6. Commit your changes and push your branch to GitHub:

    $ git add .
    $ git commit -m "Your detailed description of your changes."
    $ git push origin name-of-your-bugfix-or-feature
    
  7. Submit a pull request through the GitHub website.

Pull Request Guidelines

Before you submit a pull request, check that it meets these guidelines:

  1. The pull request should include tests.

  2. If the pull request adds functionality, the docs should be updated. Put your new functionality into a function with a docstring, and add the feature to the list in README.rst.

  3. The pull request should work for all supported Python versions. Check https://github.com/zfit/phasespace/actions/workflows/ci.yml and make sure that the tests pass for all supported Python versions.

Tips

To run a subset of tests (for example those in tests/test_generate.py):

$ pytest -k test_generate
Deploying

A reminder for the maintainers on how to deploy. Make sure all your changes are committed (including an entry in HISTORY.rst). Then run:

$ bumpversion patch # possible: major / minor / patch
$ git push
$ git push --tags

GitHub Actions will then deploy to PyPI if tests pass.

Credits

Development Lead
Developers
Contributors

None yet. Why not be the first?

Changelog

Develop
Major Features and Improvements
Behavioral changes
Bug fixes and small changes
Requirement changes
Thanks
1.10.0 (16 Apr 2024)

Add support for Python 3.12, drop support for 3.8

Major Features and Improvements
  • integrating vector support for generate: boost_to can be a Momentum Lorentz vector and return the boosted particles as a vector using as_vectors=True.

Requirement changes

Upgrade to TensorFlow > 0.16

1.9.0 (20 Jul 2023)

Add support for Python 3.11, drop support for 3.7

1.8.0 (27 Jan 2023)
Requirement changes
  • upgrade to zfit >= 0.10.0 and zfit-physics >= 0.3.0

  • pinning uproot and awkward to ~4 and ~1, respectively

1.7.0 (1. Sep 2022)

Upgraded Python and TensorFlow version.

Added tf and tensorflow extra to requirements. If you intend to use phasespace with TensorFlow in the future (and not another backend like numpy or JAX), make sure to always install with phasespace[tf].

Requirement changes
  • upgrade to TensorFlow >= 2.7

  • Python from 3.7 to 3.10 is now supported

1.6.0 (14 Apr 2022)
Major Features and Improvements
  • Improved GenMultiDecay to have better control on the decay mass of non-stable particles.

  • Added a particle_model_map argument to the GenMultiDecay class. This is a dict where the key is a particle name and the value is a mass function name. The feature can be seen in the GenMultiDecay Tutorial.

1.5.0 (27 Nov 2021)
Major Features and Improvements
  • add support to generate from a DecayChain using the decaylanguage package from Scikit-HEP. This is in the new subpackage “fromdecay” and can be used by installing the extra with pip install phasespace[fromdecay].

Requirement changes
  • drop Python 3.6 support

Thanks
  • to Simon Thor for contributing the fromdecay subpackage.

1.4.2 (5.11.2021)
Requirement changes
  • Losen restriction on TensorFlow, allow version 2.7 (and 2.5, 2.6)

1.4.1 (27.08.2021)
Requirement changes
  • Losen restriction on TensorFlow, allow version 2.6 (and 2.5)

1.4.0 (11.06.2021)
Requirement changes
  • require TensorFlow 2.5 as 2.4 breaks some functionality

1.3.0 (28.05.2021)
Major Features and Improvements
  • Support Python 3.9

  • Support TensorFlow 2.5

  • improved compilation in tf.functions, use of XLA where applicable

  • developer: modernization of setup, CI and more

Thanks
  • Remco de Boer for many commits and cleanups

1.2.0 (17.12.20)
Major Features and Improvements
  • Python 3.8 support

  • Allow eager execution by setting with tf.config.run_functions_eagerly(True) or the environment variable “PHASESPACE_EAGER”

  • Deterministic random number generation via seed or tf.random.Generator instance

Behavioral changes
Bug fixes and small changes
Requirement changes
  • tighten TensorFlow to 2.3/2.4

  • tighten TensorFlow Probability to 0.11/0.12

Thanks
  • Remco de Boer and Stefan Pflüger for discussions on random number genration

1.1.0 (27.1.2020)

This release switched to TensorFlow 2.0 eager mode. Please upgrade your TensorFlow installation if possible and change your code (minimal changes) as described under “Behavioral changes”. In case this is currently impossible to do, please downgrade to < 1.1.0.

Major Features and Improvements
  • full TF2 compatibility

Behavioral changes
  • generate now returns an eager Tensor. This is basically a numpy array wrapped by TensorFlow. To explicitly convert it to a numpy array, use the numpy() method of the eager Tensor.

  • generate_tensor is now depreceated, generate can directly be used instead.

Bug fixes and small changes
Requirement changes
  • requires now TensorFlow >= 2.0.0

Thanks
1.0.4 (13-10-2019)
Major Features and Improvements

Release to conda-forge, thanks to Chris Burr