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

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

phasespace is available on conda-forge and pip.

To install phasespace with conda, run:

$ conda install phasespace -c conda-forge

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.

For the newest development version, which may be unstable, you can install the version from git with

$ pip install git+https://github.com/zfit/phasespace

How to use

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.58
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.58
KSTARZ_MASS = 895.81
PION_MASS = 139.57018
KAON_MASS = 493.677

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)

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

Core 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 a number or tf.constant, 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 a TensorFlow Tensor:

  • 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 GenParticle class (which returns the class itself). As an example, to build the \(B^{0}\to K^{*}\gamma\) decay in which \(K^*\to K\pi\) with a fixed mass, we would write:

from phasespace import GenParticle

B0_MASS = 5279.58
KSTARZ_MASS = 895.81
PION_MASS = 139.57018
KAON_MASS = 493.677

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 expressed as arrays of dimension (n_events, 4).

N_EVENTS = 1000

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

The generate method return an eager Tensor: this is basically a wrapped numpy array. With Tensor.numpy(), it can always directly be converted to a numpy array (if really needed).

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 both generate and ~`tf.Tensor`. 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).

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([[ -603.44477608, -1945.08421179,  1557.41162855,  2715.78804872],
       [  296.88919618,  1337.65953131,  2166.92246122,  2715.78804872],
       [ 2323.74734923,    57.66407503, -1081.620211  ,  2715.78804872],
       ...,
       [-1862.62572576,  -581.33360391, -1663.04113484,  2715.78804872],
       [ 1158.02015945, -1341.81849422,  1852.44206612,  2715.78804872],
       [ -616.0590558 ,  -396.55169591, -2456.87752273,  2715.78804872]])>, 'gamma': <tf.Tensor: shape=(1000, 4), dtype=float64, numpy=
array([[  603.44477608,  1945.08421179, -1557.41162855,  2563.79195128],
       [ -296.88919618, -1337.65953131, -2166.92246122,  2563.79195128],
       [-2323.74734923,   -57.66407503,  1081.620211  ,  2563.79195128],
       ...,
       [ 1862.62572576,   581.33360391,  1663.04113484,  2563.79195128],
       [-1158.02015945,  1341.81849422, -1852.44206612,  2563.79195128],
       [  616.0590558 ,   396.55169591,  2456.87752273,  2563.79195128]])>, 'pi-': <tf.Tensor: shape=(1000, 4), dtype=float64, numpy=
array([[-341.4054273 , -711.31308624,  265.76212775,  844.17611678],
       [ 194.41001606,  585.73173112,  455.42086878,  779.59278989],
       [1201.8176238 , -208.97447493, -420.46318868, 1297.80779017],
       ...,
       [-282.09369454,  -75.01062151, -606.43007722,  687.34323271],
       [ 336.51407506, -698.33259367, 1037.94177656, 1302.96320008],
       [   6.63377537,  -64.15941644,  -53.62883711,  162.83834011]])>, 'K+': <tf.Tensor: shape=(1000, 4), dtype=float64, numpy=
array([[ -262.03934878, -1233.77112554,  1291.64950079,  1871.61193194],
       [  102.47918012,   751.9278002 ,  1711.50159243,  1936.19525883],
       [ 1121.92972544,   266.63854996,  -661.15702232,  1417.98025855],
       ...,
       [-1580.53203121,  -506.32298239, -1056.61105762,  2028.44481601],
       [  821.50608439,  -643.48590055,   814.50028956,  1412.82484863],
       [ -622.69283117,  -332.39227947, -2403.24868562,  2552.94970861]])>}

It is worth noting that the graph generation is cached even when using generate, so 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 a ~`tf.Tensor` 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):

import tensorflow as tf
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 = tf.cast(min_mass, tf.float64)
   max_mass = tf.cast(max_mass, tf.float64)
   kstar_width_cast = tf.cast(KSTARZ_WIDTH, tf.float64)
   kstar_mass_cast = tf.cast(KSTARZ_MASS, dtype=tf.float64)

   kstar_mass = tf.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

N_EVENTS = 1000

B0_MASS = 5279.58
PION_MASS = 139.57018
KAON_MASS = 493.677

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 TensorFlow to build a graph of the computations. This is usually more performant, especially if used multiple times. However, this has the disadvantage that _inside_ phasespac, the actual values are not computed on Python runtime, e.g. if a breakpoint is set the values of a tf.Tensor won’t be available.

TensorFlow (since version 2.0) however can easily switch to so called “eager execution”: in this mode, it behaves the same as Numpy; values are computed instantly and the Python code is not only executed once but every time.

To switch this on or off, the global flag in TensorFlow tf.config.run_functions_eagerly(True) or the enviroment variable “PHASESPACE_EAGER” (which switches this flag) can be used.

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.

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: Union[Callable, int, float])[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, Tensor`, callable) – Mass of the particle. If it’s a float, it get converted to a tf.constant.

generate(n_events: Union[int, tensorflow.python.framework.ops.Tensor, tensorflow.python.ops.variables.Variable], boost_to: Optional[tensorflow.python.framework.ops.Tensor] = None, normalize_weights: bool = True, seed: Optional[Union[int, tensorflow.python.ops.stateful_random_ops.Generator]] = None)Tuple[tensorflow.python.framework.ops.Tensor, Dict[str, tensorflow.python.framework.ops.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. 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.

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: tensorflow.python.framework.ops.Tensor = None, max_mass: tensorflow.python.framework.ops.Tensor = None, n_events: Union[tensorflow.python.framework.ops.Tensor, tensorflow.python.ops.variables.Variable] = None, seed: Optional[Union[int, tensorflow.python.ops.stateful_random_ops.Generator]] = None)tensorflow.python.framework.ops.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: Optional[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.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.

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 Python 2.7, 3.4, 3.5 and 3.6, and for PyPy. 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
Core 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.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