-
Notifications
You must be signed in to change notification settings - Fork 31
/
measure.py
343 lines (268 loc) · 12.8 KB
/
measure.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
r"""### Energy and Angular Measures
The appropriate notions of energy and angle depend on the collider context. Typically, one wants
to work with observables that respect the appropriate Lorentz subgroup for the collision type
of interest. EnergyFlow is capable of handling two broad classes of measures: $e^+e^-$ and
hadronic, which are selected using the required `measure` argument. For substructure applications,
it is often convenient to normalize the energies so that $\sum_iz_i=1$. The `normed` keyword
argument is provided to control normalization of the energies (default is `True`).
Each measure comes with a parameter $\beta>0$ which controls the relative weighting between
smaller and larger anglular structures. This can be set using the `beta` keyword argument
(default is `1`). There is also a $\kappa$ parameter to control the relative weighting
between soft and hard energies. This can be set using the `kappa` keyword argument
(default is `1`). Only `kappa=1` yields collinear-safe observables.
Beyond the measures implemented here, the user can implement their own custom measure by
passing in $\{z_i\}$ and $\{\theta_{ij}\}$ directly to the EFP classes.
#### Hadronic Measures
For hadronic collisions, observables are typically desired to be invariant under boosts along
the beam direction and rotations about the beam direction. Thus, particle transverse momentum
$p_T$ and rapidity-azimuth coordinates $(y,\phi)$ are used.
There are two hadronic measures implemented in EnergyFlow: `'hadr'` and `'hadrdot'`.
These are listed explicitly below.
`'hadr'`:
$$z_i=p_{T,i}^{\kappa},\quad\quad \theta_{ij}=(\Delta y_{ij}^2 + \Delta\phi_{ij}^2)^{\beta/2}.$$
`'hadrdot'`:
$$z_i=p_{T,i}^{\kappa},\quad\quad \theta_{ij}=\left(\frac{2p^\mu_ip_{j\mu}}{p_{T,i}p_{T,j}}
\right)^{\beta/2}.$$
#### *e+e-* Measures
For $e^+e^-$ collisions, observables are typically desired to be invariant under the full
group of rotations about the interaction point. Since the center of momentum energy is known,
the particle energy $E$ is typically used. For the angular measure, pairwise Lorentz contractions
of the normalized particle four-momenta are used.
There is one $e^+e^-$ measure implemented.
`'ee'`:
$$z_i = E_{i}^{\kappa},
\quad\quad \theta_{ij} = \left(\frac{2p_i^\mu p_{j \mu}}{E_i E_j}\right)^{\beta/2}.$$
"""
from __future__ import absolute_import, division, print_function
from abc import ABCMeta, abstractmethod
import warnings
import numpy as np
from six import with_metaclass
from energyflow.utils import transfer
from energyflow.utils.particle_utils import *
__all__ = ['Measure']
# special value of kappa indicating "particle flow"
pf_marker = 'pf'
# form theta_ij**2 matrix from array of (rapidity,phi) values
# theta_ij**2 = (y_i - y_j)**2 + (phi_i - phi_j)**2
def thetas2_from_yphis(yphis):
X = yphis[:,np.newaxis] - yphis[np.newaxis,:]
X[...,0] **= 2
X[...,1] = (np.pi - np.abs(np.abs(X[...,1]) - np.pi))**2
return X[...,0] + X[...,1]
# get theta_ij**2 matrix from four-vectors using combination of above functions
def thetas2_from_p4s(p4s):
return thetas2_from_yphis(np.vstack([ys_from_p4s(p4s), phis_from_p4s(p4s)]).T)
# kappa is a number, so raise energies to that number and form phats
def kappa_func(Es, ps, kappa):
return Es**kappa, ps/Es[:,np.newaxis]
# kappa indicates particle flow, so make energies 1 and leave ps alone
def pf_func(Es, ps, kappa):
return np.ones(Es.shape), ps
###############################################################################
# Measure
###############################################################################
class Measure(with_metaclass(ABCMeta, object)):
"""Class for dealing with any kind of measure."""
def __new__(cls, *args, **kwargs):
if cls is Measure:
measure = args[0]
if 'hadr' in measure:
return super(Measure, cls).__new__(HadronicMeasure.factory(measure))
if 'ee' in measure:
return super(Measure, cls).__new__(EEMeasure.factory(measure))
raise NotImplementedError('measure {} is unknown'.format(measure))
else:
return super(Measure, cls).__new__(cls)
def __init__(self, measure, beta=1, kappa=1, normed=True, coords=None, check_input=True):
r"""Processes inputs according to the measure choice.
**Arguments**
- **measure** : _string_
- The string specifying the energy and angular measures to use.
- **beta** : _float_
- The angular weighting exponent $\beta$. Must be positive.
- **kappa** : {_float_, `'pf'`}
- If a number, the energy weighting exponent $\kappa$. If `'pf'`,
use $\kappa=v$ where $v$ is the valency of the vertex. `'pf'`
cannot be used with measure `'hadr'`. Only IRC-safe for `kappa=1`.
- **normed** : bool
- Whether or not to use normalized energies.
- **coords** : {`'ptyphim'`, `'epxpypz'`, `None`}
- Controls which coordinates are assumed for the input. If `'ptyphim'`, the
fourth column (the masses) is optional and massless particles are assumed
if it is not present. If `None`, coords with be `'ptyphim'` if using a
hadronic measure and `'epxpypz'` if using the e+e- measure.
- **check_input** : bool
- Whether to check the type of input each time or assume the first input type.
"""
transfer(self, locals(), ['measure', 'kappa', 'normed', 'coords', 'check_input'])
self.beta = float(beta)
self.half_beta = self.beta/2
assert self.beta > 0
if self.coords not in [None, 'epxpypz', 'ptyphim']:
raise ValueError('coords must be one of epxpypz, ptyphim, or None')
self.need_meas_func = True
def evaluate(self, arg):
"""Evaluate the measure on a set of particles.
**Arguments**
- **arg** : _2-d numpy.ndarray_
- A two-dimensional array of the particles with each row being a
particle and the columns specified by the `coords` attribute.
**Returns**
- (_ 1-d numpy.ndarray_, _2-d numpy.ndarray_)
- (`zs`, `thetas`) where `zs` is a vector of the energy fractions for
each particle and `thetas` is the distance matrix between the particles.
"""
# check type only if needed
if self.need_meas_func or self.check_input:
self.set_meas_func(arg)
# get zs and angles
zs, angles = self.meas_func(arg)
return (zs/np.sum(zs) if self.normed else zs), angles
def set_meas_func(self, arg):
# support arg as numpy.ndarray
if isinstance(arg, np.ndarray):
self.meas_func = self.array_handler(arg.shape[1])
# support arg as list (of lists)
elif isinstance(arg, list):
array_meas = self.array_handler(len(arg[0]))
def wrapped_meas(arg):
return array_meas(np.asarray(arg))
self.meas_func = wrapped_meas
# support arg as fastjet pseudojet
elif hasattr(arg, 'constituents'):
self.meas_func = self.pseudojet
# raise error if not one of these types
else:
raise TypeError('arg is not one of numpy.ndarray, list, or fastjet.PseudoJet')
self.need_meas_func = False
@abstractmethod
def array_handler(self, dim):
pass
@abstractmethod
def pseudojet(self, arg):
pass
def _ps_dot(self, ps):
return np.abs(2*np.dot(ps[:,np.newaxis]*ps[np.newaxis,:], self.metric))
def _set_k_func(self):
self._k_func = kappa_func
if self.kappa == pf_marker:
if self.normed:
warnings.warn('Normalization not supported when kappa=\'' + pf_marker + '\'.')
self.normed = False
self._k_func = pf_func
###############################################################################
# HadronicMeasure
###############################################################################
class HadronicMeasure(Measure):
@staticmethod
def factory(measure):
if 'dot' in measure:
return HadronicDotMeasure
return HadronicDefaultMeasure
def __init__(self, *args, **kwargs):
super(HadronicMeasure, self).__init__(*args, **kwargs)
self._set_k_func()
if self.coords is None:
self.coords = 'ptyphim'
self.epxpypz = (self.coords == 'epxpypz')
def array_handler(self, dim):
if dim == 3:
if self.epxpypz:
raise ValueError('epxpypz coords require second dimension of arg to be 4')
return self.ndarray_dim3
elif dim == 4:
return self.ndarray_dim4
else:
raise ValueError('second dimension of arg must be either 3 or 4')
@abstractmethod
def ndarray_dim3(self, arg):
pass
@abstractmethod
def ndarray_dim4(self, arg):
pass
@abstractmethod
def pseudojet(self, arg):
constituents = arg.constituents()
pts = np.asarray([c.pt() for c in constituents])
return pts, constituents
###############################################################################
# EEMeasure
###############################################################################
class EEMeasure(Measure):
@staticmethod
def factory(measure):
return EEDefaultMeasure
def __init__(self, *args, **kwargs):
super(EEMeasure, self).__init__(*args, **kwargs)
self._set_k_func()
if self.coords is None:
self.coords = 'epxpypz'
self.epxpypz = (self.coords == 'epxpypz')
def array_handler(self, dim):
if dim < 2:
raise ValueError('second dimension of arg must be >= 2')
if not self.epxpypz and dim not in [3, 4]:
raise ValueError('ptyphim coords only work with inputs of dimension 3 or 4')
self.metric = flat_metric(dim)
return self.ndarray_dim_arb
@abstractmethod
def ndarray_dim_arb(self, arg):
pass
@abstractmethod
def pseudojet(self, arg):
constituents = arg.constituents()
Es = np.asarray([c.e() for c in constituents])
return Es, constituents
###############################################################################
# HadronicDefaultMeasure
###############################################################################
class HadronicDefaultMeasure(HadronicMeasure):
def __init__(self, *args, **kwargs):
super(HadronicDefaultMeasure, self).__init__(*args, **kwargs)
if self.kappa == pf_marker:
raise ValueError('particle flow not available for HadronicDefaultMeasure')
def ndarray_dim3(self, arg):
return arg[:,0]**self.kappa, thetas2_from_yphis(arg[:,(1,2)])**self.half_beta
def ndarray_dim4(self, arg):
if self.epxpypz:
return pts_from_p4s(arg)**self.kappa, thetas2_from_p4s(arg)**self.half_beta
else:
return self.ndarray_dim3(arg[:,:3])
def pseudojet(self, arg):
pts, constituents = super(HadronicDefaultMeasure, self).pseudojet(arg)
thetas = np.asarray([[c1.delta_R(c2) for c2 in constituents] for c1 in constituents])
return pts**self.kappa, thetas**self.beta
###############################################################################
# HadronicDotMeasure
###############################################################################
class HadronicDotMeasure(HadronicMeasure):
metric = flat_metric(4)
def ndarray_dim3(self, arg):
pts, p4s = self._k_func(arg[:,0], p4s_from_ptyphims(arg), self.kappa)
return pts, self._ps_dot(p4s)**self.half_beta
def ndarray_dim4(self, arg):
if self.epxpypz:
pts, p4s = self._k_func(pts_from_p4s(arg), arg, self.kappa)
else:
pts, p4s = self._k_func(arg[:,0], p4s_from_ptyphims(arg), self.kappa)
return pts, self._ps_dot(p4s)**self.half_beta
def pseudojet(self, arg):
pts, constituents = super(HadronicDotMeasure, self).pseudojet(arg)
p4s = np.asarray([[c.e(), c.px(), c.py(), c.pz()] for c in constituents])
pts, p4s = self._k_func(pts, p4s, self.kappa)
return pts, self._ps_dot(p4s)**self.half_beta
###############################################################################
# EEDefaultMeasure
###############################################################################
class EEDefaultMeasure(EEMeasure):
def ndarray_dim_arb(self, arg):
if not self.epxpypz:
arg = p4s_from_ptyphims(arg)
Es, ps = self._k_func(arg[:,0], arg, self.kappa)
return Es, self._ps_dot(ps)**self.half_beta
def pseudojet(self, arg):
Es, constituents = super(EEDefaultMeasure, self).pseudojet(arg)
p4s = np.asarray([[c.e(), c.px(), c.py(), c.pz()] for c in constituents])
Es, p4s = self._k_func(Es, p4s, self.kappa)
return Es, self._ps_dot(p4s)**self.half_beta