EnergyFlow Demo

EnergyFlow website

Energy Flow Polynomials

In this tutorial, we introduce the EnergyFlow package for computing arbitrary Energy Flow Polynomials (EFPs) and collections of EFPs. The package is built with both usability, flexibility, and computational efficiency in mind.

For a collection of $M$ particles with energy measure $z_i$ and pairwise angular distance measures $\theta_{ij}$, the EFPs are multiparticle energy correlator observables, indexed by multigraphs $G$, and defined as: $$ \mathrm{EFP}_G = \sum_{i_1=1}^M \sum_{i_2=1}^M\cdots \sum_{i_N=1}^M z_{i_1}z_{i_2}\cdots z_{i_N} \prod_{(k,\ell)\in G} \theta_{i_ki_\ell},$$ where $(k,\ell)\in G$ iterates over the edges in the multigraph.

Choices of Measure

The specific choice of energy and angular measures depends on the collider context. We provide the following measure options (default is hadrdot):

hadr: $$ z_i = p_{T,i}^\kappa,\,\,\,\,\,\, \theta_{ij} = ((y_i - y_j)^2 + (\phi_i - \phi_j)^2)^{\beta/2} $$ hadrdot: $$ z_i = p_{T,i}^\kappa,\,\,\,\,\,\, \theta_{ij} = (2\, \hat p_i^\mu \hat p_{i\,\mu})^{\beta/2},\,\,\,\,\,\, \mathrm{where }\,\,\,\,\,\,\hat p_i^\mu \equiv p_i^\mu/p_{T,i} $$ ee: $$ z_i = E_i^\kappa,\,\,\,\,\,\, \theta_{ij} = (2\, \hat p_i^\mu \hat p_{i\,\mu})^{\beta/2},\,\,\,\,\,\, \mathrm{where }\,\,\,\,\,\,\hat p_i^\mu \equiv p_i^\mu/E_i $$

The energy and angular weighting parameters kappa and beta default to $\kappa=1$ and $\beta = 2$. The choice of $\kappa = 1$ is required for infrared and collinear (IRC) safety of the observables. Any choice of $\beta > 0$ guarantees IRC safety. We also provide the normed option (default is True) to use a normalized and dimensionless energy measure $z_i/\sum_{j=1}^M z_j$ rather than $z_i$.

With this refresher, we have enough to begin using EnergyFlow to compute arbitrary EFPs! Ensure you have EnergyFlow installed. It is easily pip installable by executing: pip install energyflow

We start by importing the EnergyFlow package as well as some other helpful Python libraries.

In [1]:
import energyflow as ef
import numpy as np

Obtaining Events

Let's get some events. Typically these would be read in from your favorite event generator for your physics case of interest.

For the purposes of this tutorial, we will uniformly sample massless $M$-body phase space with the use of our implementation of the RAMBO algorithm via the function ef.gen_massless_phase_space. It returns nevents events consisting of nparticles massless four-momenta with center of mass energy of energy in the center of momentum frame. In general, EnergyFlow supports events as arrays of four-momenta [E,px,py,pz] or arrays of hadronic coordinates [pT,y,phi] for hadronic measures.

Let's generate 50 events with 20 particles each at center of mass energy 100 GeV.

In [2]:
events = ef.gen_massless_phase_space(nevents=50, nparticles=20, energy=100)

Defining Graphs

To specify a particular EFP to be computed for our events, we must simply specify the corresponding multigraph.

In EnergyFlow, multigraphs are specified as lists of edges, where edges are pairs of vertices. Here are several examples of graphs given as edge lists, where we label the vertices with integers from $0$ to $N-1$.

In [3]:
dot      = []
line     = [(0,1)]
wedge    = [(0,1), (0,2)]
triangle = [(0,1), (0,2), (1,2)]
square   = [(0,1), (1,2), (2,3), (3,1)]

Multigraphs can have multiple edges per pair of vertices. Here are several examples of multigraphs as edgelists, with one or more doubled edges.

In [4]:
multiline      = [(0,1), (0,1)]
multiwedge1    = [(0,1), (0,1), (0,2)]
multiwedge2    = [(0,1), (0,1), (0,2), (0,2)]
multitriangle1 = [(0,1), (0,2), (0,2), (1,2)]

In fact, arbitrary objects can be used to label the vertices. We typically use integer-labeling for simplicity and readibility, though this is not required. The following two edge lists both define the same graph from the perspective of EnergyFlow.

In [5]:
graph_LHC  = [('atlas','cms'),('atlas','cms'),('atlas','lhcb'),('cms','lhcb'),('lhcb','alice')]
graph_ints = [(0,1),(0,1),(0,2),(1,2),(2,3)]

Computing an Energy Flow Polynomial

EFPs are defined by passing a graph and measure choies to ef.EFP. The compute method of ef.EFP can then be called on an event in order to compute the EFP on that event.

For concreteness, let's begin by defining the EFP corresponding to the line graph [(0,1)] (one edge connecting two vertices) which should be equal to twice the squared center of mass energy of the event with a suitable choice of measure parameters.

In [6]:
# Computing a single EFP on a single event

# specify a graph and event
graph = [(0,1)]
event = events[0]

# define the EFP corresponding to the specified graph
EFP_graph = ef.EFP(graph, measure='hadrdot', beta=2, normed=False, coords='epxpypz')

# compute the EFP on the specified event
result = EFP_graph.compute(event)

print("EFP Value:", result)
EFP Value: 20000.000000000007

We can see that the value of this EFP is indeed $\mathrm{EFP}_{[(0,1)]} = 2\times E_\mathrm{CM}^2 = 2\times (100\,\mathrm{GeV})^2 = 20\,000\, \mathrm{GeV}^2$, as expected for our events.

The framework above can be immediately extended to computing the value of an EFP on a collection of many events. This can either be done with list comprehension or even more simply using the batch_compute method, which tries to use as many processes as there are CPUs in the machine. The number of worker processes can be controlled by passing in n_jobs.

In [7]:
# Computing a single EFP on many events

# specify a graph
graph = [(0,1), (0,2), (0,2), (1,2)]

# define the EFP corresponding to the specified graph
efp_graph = ef.EFP(graph, measure='hadr', beta=1, normed=True)

# compute the EFP on the collection of events with list comprehension
results = np.asarray([efp_graph.compute(event) for event in events])
print("EFP w. list comprehension:", results)
EFP w. list comprehension: [ 3425.95623748  7415.20417255 10639.10756502   538.54450213
  2257.7204505   1388.99529475  1804.2802471   2379.91818331
   744.76895496  2448.44326014   791.65788366  1777.23263608
  2549.01010263  1345.0581705   3173.89847025   583.10565749
  6402.76565968   588.39404401  1111.39724888  1300.48733859
  1124.37600351  2011.89562027  2421.43308899   366.13033974
 24689.53697442  1601.08307778  3088.81736771  3751.02824457
   475.10900468  3256.96085805   545.84972603  5007.8556729
   863.92331652  1750.0094433   2483.6595287  11504.25260136
  2083.06997218  1249.40659285  1745.32638252  3072.47495042
   884.73536234   628.91049857  4019.04287973  1189.00573617
  1085.78678264  1998.76301115  2675.57238823 18365.35483771
   355.40085464  1601.72095566]

Computing a set of Energy Flow Polynomials

If we are interested in using the spanning properties of EFPs, we typically want to compute large numbers of EFPs. The EnergyFlow package does the heavy lifting for you! It contains information about all multigraphs with up to 10 edges, which can be easily and intuitively accessed using EFPSet.

Relevant graph quantities which can easily be used to select a set of multigraphs are summarized below. More options are available: for a full list see the documentation.

n : Number of vertices in the multigraph.

d : Degree, or number of edges in the multigraph.

v : Maximum valency (number of edges touching a vertex).

c : Variable Elimination computational complexity $\mathcal O(M^c)$

p : Number of prime factors (or connected components of the multigraph).

As a basis, EFPs are typically organized by the number of edges d. Not only does this correspond to the degree of the polynomial, but there are also a finite number of multigraphs up to a specified degree d. Lets get all EFPs with up to five edges by passing d<=5 to EFPSet and compute them on our events!

In [8]:
# get all EFPs with d<=5 (up to d<=10 available by default)
efpset = ef.EFPSet('d<=5', measure='hadr', beta=1, normed=True, verbose=True)

# compute their values on our events
results = np.asarray([efpset.compute(event) for event in events])

print("Results:", results)
Originally Available EFPs:
  Prime: 23691
  Composite: 21540
  Total:  45231
Current Stored EFPs:
  Prime: 54
  Composite: 48
  Total:  102
Results: [[1.00000000e+00 6.47741037e+00 7.09165172e+01 ... 1.92731027e+04
  1.33207431e+04 1.14026833e+04]
 [1.00000000e+00 8.38195548e+00 1.02475169e+02 ... 6.03468623e+04
  4.40736076e+04 4.13739275e+04]
 [1.00000000e+00 9.14192496e+00 1.25226262e+02 ... 9.56771818e+04
  6.87190500e+04 6.38540225e+04]
 ...
 [1.00000000e+00 1.08365705e+01 1.91047936e+02 ... 2.43118459e+05
  1.62181681e+05 1.49437402e+05]
 [1.00000000e+00 3.92684387e+00 2.14485599e+01 ... 1.29876041e+03
  9.77240038e+02 9.33723246e+02]
 [1.00000000e+00 5.61951356e+00 4.45468201e+01 ... 7.90520029e+03
  5.97753065e+03 5.60394181e+03]]

We can very easily do much more sophisticated EFP selections using any of the graph quantities. For example to select all prime EFPs with at most 4 vertices and at most 5 edges, we simply use EFPSet with the following intuitive syntax:

In [9]:
# get all EFPs with n<=4, d<=5, that are prime (i.e. p==1)
efpset = ef.EFPSet(('n<=',6), ('d<=',5), ('p==',1), measure='hadr', beta=1, normed=True, verbose=True)

# compute their values on our events
results = np.asarray([efpset.compute(event) for event in events])
Originally Available EFPs:
  Prime: 23691
  Composite: 21540
  Total:  45231
Current Stored EFPs:
  Prime: 54
  Composite: 0
  Total:  54

To study the properties of an individual EFP within EFPSet we can easily do this using specs and graphs. Suppose we are interested in the 95th EFP. We can find out its graph and all of its relevant information simply with the following syntax. Here we show only a subset of all the information that specs provides about an EFP: for a full list see the documentation.

In [11]:
efpset = ef.EFPSet('d<=5', measure='hadr', beta=1, normed=True, verbose=True)

ind = 96
graph = efpset.graphs(ind)
n, _, d, v, _, c, p, _ = efpset.specs[ind]

print("Graph:", graph)
print("Number of vertices, n:", n)
print("Number of edges,    d:", d)
print("Maximum valency,    v:", v)
print("VE complexity,      c:", c)
print("Number of primes,   p:", p)
Originally Available EFPs:
  Prime: 23691
  Composite: 21540
  Total:  45231
Current Stored EFPs:
  Prime: 54
  Composite: 48
  Total:  102
Graph: [(0, 1), (2, 3), (4, 5), (4, 6), (4, 7)]
Number of vertices, n: 8
Number of edges,    d: 5
Maximum valency,    v: 3
VE complexity,      c: 2
Number of primes,   p: 3

Custom measures

If you want to specify your own custom measure for $z_i$ and $\theta_{ij}$ to be used in the EFP formula, that's also possible within the EnergyFlow package. You simply compute your own custom-defined zs and thetas on the events and pass them to the compute methods.

We demonstrate this below by using random numbers for zs and thetas, which can be replaced with any custom values.

In [12]:
# zs and thetas can be passed in explicitly if you want to use a custom measure
(zs, thetas) = (np.random.rand(100,25), np.random.rand(100,25,25))

results = np.asarray([efpset.compute(zs=z, thetas=theta) for (z,theta) in zip(zs,thetas)])
print(results)
[[1.04637461e+01 5.15570630e+01 3.20219020e+01 ... 4.38845484e+06
  3.52265327e+07 3.64284630e+08]
 [1.36306266e+01 9.25341706e+01 6.12626231e+01 ... 4.85402489e+07
  5.05375289e+08 6.78438811e+09]
 [1.11692742e+01 6.42131647e+01 4.32774631e+01 ... 1.14586649e+07
  9.95863740e+07 1.09174289e+09]
 ...
 [1.14937327e+01 6.52426006e+01 4.37276509e+01 ... 1.21436695e+07
  1.03871374e+08 1.18210572e+09]
 [1.33279983e+01 8.80756053e+01 5.82535104e+01 ... 3.98005442e+07
  4.01299768e+08 5.30002822e+09]
 [1.20689116e+01 7.29417525e+01 4.87920981e+01 ... 1.89355565e+07
  1.73257390e+08 2.06481415e+09]]

And that's it! Now you should be able to specify any EFP (i.e. multigraph) or set of EFPs that you want to comput with EFP or EFPSet, compute them on a set of events with compute or batch_compute, and study the results with specs and graphs! As always, see the documentation for a full description of the EnergyFlow package.