This file is indexed.

/usr/include/Rivet/Jet.hh is in librivet-dev 1.8.3-1.3.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  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
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
// -*- C++ -*-
#ifndef RIVET_Jet_HH
#define RIVET_Jet_HH

#include "Rivet/Rivet.hh"
#include <numeric>

namespace Rivet {


  /// @brief Representation of a clustered jet of particles.
  class Jet : public ParticleBase {
  public:

    /// @name Constructors
    //@{

    Jet() : ParticleBase() { clear(); }

    /// Set all the jet data, with full particle information.
    Jet(const vector<Particle>& particles, const FourMomentum& pjet)
      : ParticleBase() {
      setState(particles, pjet);
    }

    // /// Set all the jet data, without particle ID information.
    // Jet(const vector<FourMomentum>& momenta, const FourMomentum& pjet)
    //   : ParticleBase() {
    //   setState(momenta, pjet);
    // }

    //@}


    /// @name Access jet constituents
    //@{

    /// Number of particles in this jet.
    size_t size() const { return _particles.size(); }

    // /// Define a Jet::iterator via a typedef.
    // typedef vector<FourMomentum>::iterator iterator;

    // /// Define a Jet::const_iterator via a typedef.
    // typedef vector<FourMomentum>::const_iterator const_iterator;

    // /// Get a begin iterator over the particle/track four-momenta in this jet.
    // iterator begin() {
    //   return _momenta.begin();
    // }

    // /// Get an end iterator over the particle/track four-momenta in this jet.
    // iterator end() {
    //   return _momenta.end();
    // }

    // /// Get a const begin iterator over the particle/track four-momenta in this jet.
    // const_iterator begin() const {
    //   return _momenta.begin();
    // }

    // /// Get a const end iterator over the particle/track four-momenta in this jet.
    // const_iterator end() const {
    //   return _momenta.end();
    // }

    // /// Get the track momenta in this jet.
    // vector<FourMomentum>& momenta() {
    //   return _momenta;
    // }

    // /// Get the track momenta in this jet (const version).
    // const vector<FourMomentum>& momenta() const {
    //   return _momenta;
    // }

    /// Get the particles in this jet.
    vector<Particle>& particles() { return _particles; }

    /// Get the particles in this jet (const version)
    const vector<Particle>& particles() const { return _particles; }

    /// Check whether this jet contains a particular particle.
    bool containsParticle(const Particle& particle) const;

    /// Check whether this jet contains a certain particle type.
    bool containsParticleId(PdgId pid) const;

    /// Check whether this jet contains at least one of certain particle types.
    bool containsParticleId(const vector<PdgId>& pids) const;

    /// Check whether this jet contains a charm-flavoured hadron (or decay products from one).
    bool containsCharm() const;

    /// Check whether this jet contains a bottom-flavoured hadron (or decay products from one).
    bool containsBottom() const;

    //@}


    /// @name Access the effective jet 4-vector properties
    //@{

    /// Get equivalent single momentum four-vector.
    const FourMomentum& momentum() const { return _momentum; }

    /// Get the unweighted average \f$ \eta \f$ for this jet. (caches)
    double eta() const { return momentum().eta(); }

    /// Get the unweighted average \f$ \phi \f$ for this jet. (caches)
    double phi() const { return momentum().phi(); }

    /// Get the total energy of this jet.
    double totalEnergy() const { return momentum().E(); }

    /// Get the energy carried in this jet by neutral particles.
    double neutralEnergy() const;

    /// Get the energy carried in this jet by hadrons.
    double hadronicEnergy() const;

    /// Get the sum of the \f$ p_T \f$ values of the constituent tracks. (caches)
    double ptSum() const { return momentum().pT(); }

    /// Get the sum of the \f$ E_T \f$ values of the constituent tracks. (caches)
    double EtSum() const { return momentum().Et(); }

    //@}


    /// @name Set the jet constituents and properties
    //@{

    /// Set all the jet data, with full particle information.
    Jet& setState(const vector<Particle>& particles, const FourMomentum& pjet);

    // /// Set all the jet data, without particle ID information.
    // Jet& setState(const vector<FourMomentum>& momenta, const FourMomentum& pjet);

    /// Set the effective 4-momentum of the jet.
    Jet& setMomentum(const FourMomentum& momentum);

    /// Set the particles collection with full particle information.
    Jet& setParticles(const vector<Particle>& particles);

    // /// Set the particles collection with momentum information only.
    // Jet& setParticles(const vector<FourMomentum>& momenta);

    // /// Add a particle/track to this jet.
    // Jet& addParticle(const FourMomentum& particle);

    // /// Add a particle/track to this jet.
    // Jet& addParticle(const Particle& particle);

    /// Reset this jet as empty.
    Jet& clear();

    //@}


  private:

    /// Full particle information including tracks, ID etc
    ParticleVector _particles;

    // /// The particle momenta.
    // /// @todo Eliminate this to ensure consistency.
    // std::vector<FourMomentum> _momenta;

    /// Effective jet 4-vector
    FourMomentum _momentum;

  };


  /// Typedef for a collection of Jet objects.
  typedef std::vector<Jet> Jets;


  /// @name Jet comparison functions for STL sorting
  //@{

  /// @brief Compare jets by \f$ p_\perp \f$ (descending - usual sorting for HEP)
  /// Use this so that highest \f$ p_\perp \f$ is at the front of the list
  inline bool cmpJetsByPt(const Jet& a, const Jet& b) {
    return a.ptSum() > b.ptSum();
  }
  /// @brief Compare jets by \f$ p_\perp \f$ (ascending)
  /// Use this so that lowest \f$ p_\perp \f$ is at the front of the list
  inline bool cmpJetsByAscPt(const Jet& a, const Jet& b) {
    return a.ptSum() < b.ptSum();
  }

  /// @brief Compare jets by descending momentum, \f$ p \f$
  inline bool cmpJetsByP(const Jet& a, const Jet& b) {
    return a.momentum().vector3().mod() > b.momentum().vector3().mod();
  }
  /// @brief Compare jets by ascending momentum, \f$ p \f$
  inline bool cmpJetsByAscP(const Jet& a, const Jet& b) {
    return a.momentum().vector3().mod() < b.momentum().vector3().mod();
  }

  // @brief Compare jets by \f$ E_\perp \f$ (descending - usual sorting for HEP)
  /// Use this so that highest \f$ E_\perp \f$ is at the front of the list
  inline bool cmpJetsByEt(const Jet& a, const Jet& b) {
    return a.EtSum() > b.EtSum();
  }
  // @brief Compare jets by \f$ E_\perp \f$ (ascending)
  /// Use this so that lowest \f$ E_\perp \f$ is at the front of the list
  inline bool cmpJetsByEtDesc(const Jet& a, const Jet& b) {
    return a.EtSum() < b.EtSum();
  }

  /// @brief Compare jets by \f$ E \f$ (descending - usual sorting for HEP)
  /// Use this so that highest \f$ E \f$ is at the front of the list
  inline bool cmpJetsByE(const Jet& a, const Jet& b) {
    return a.momentum().E() > b.momentum().E();
  }
  /// @brief Compare jets by \f$ E \f$ (ascending)
  /// Use this so that lowest \f$ E \f$ is at the front of the list
  inline bool cmpJetsByAscE(const Jet& a, const Jet& b) {
    return a.momentum().E() < b.momentum().E();
  }

  /// @brief Compare jets by \f$ \eta \f$ (descending)
  /// Use this so that highest \f$ \eta \f$ is at the front of the list
  inline bool cmpJetsByDescPseudorapidity(const Jet& a, const Jet& b) {
    return a.momentum().pseudorapidity() > b.momentum().pseudorapidity();
  }
  /// @brief Compare jets by \f$ \eta \f$ (ascending)
  /// Use this so that lowest \f$ \eta \f$ is at the front of the list
  inline bool cmpJetsByAscPseudorapidity(const Jet& a, const Jet& b) {
    return a.momentum().pseudorapidity() < b.momentum().pseudorapidity();
  }

  /// @brief Compare jets by \f$ |\eta| \f$ (descending)
  /// Use this so that highest \f$ |\eta| \f$ is at the front of the list
  inline bool cmpJetsByDescAbsPseudorapidity(const Jet& a, const Jet& b) {
    return fabs(a.momentum().pseudorapidity()) > fabs(b.momentum().pseudorapidity());
  }
  /// @brief Compare jets by \f$ |\eta| \f$ (ascending)
  /// Use this so that lowest \f$ |\eta| \f$ is at the front of the list
  inline bool cmpJetsByAscAbsPseudorapidity(const Jet& a, const Jet& b) {
    return fabs(a.momentum().pseudorapidity()) < fabs(b.momentum().pseudorapidity());
  }

  /// @brief Compare jets by \f$ y \f$ (descending)
  /// Use this so that highest \f$ y \f$ is at the front of the list
  inline bool cmpJetsByDescRapidity(const Jet& a, const Jet& b) {
    return a.momentum().rapidity() > b.momentum().rapidity();
  }
  /// @brief Compare jets by \f$ y \f$ (ascending)
  /// Use this so that lowest \f$ y \f$ is at the front of the list
  inline bool cmpJetsByAscRapidity(const Jet& a, const Jet& b) {
    return a.momentum().rapidity() < b.momentum().rapidity();
  }

  /// @brief Compare jets by \f$ |y| \f$ (descending)
  /// Use this so that highest \f$ |y| \f$ is at the front of the list
  inline bool cmpJetsByDescAbsRapidity(const Jet& a, const Jet& b) {
    return fabs(a.momentum().rapidity()) > fabs(b.momentum().rapidity());
  }
  /// @brief Compare jets by \f$ |y| \f$ (ascending)
  /// Use this so that lowest \f$ |y| \f$ is at the front of the list
  inline bool cmpJetsByAscAbsRapidity(const Jet& a, const Jet& b) {
    return fabs(a.momentum().rapidity()) < fabs(b.momentum().rapidity());
  }

  //@}

  inline double deltaR(const Jet& j1, const Jet& j2,
                       RapScheme scheme = PSEUDORAPIDITY) {
    return deltaR(j1.momentum(), j2.momentum(), scheme);
  }

  inline double deltaR(const Jet& j, const Particle& p,
                       RapScheme scheme = PSEUDORAPIDITY) {
    return deltaR(j.momentum(), p.momentum(), scheme);
  }

  inline double deltaR(const Particle& p, const Jet& j,
                       RapScheme scheme = PSEUDORAPIDITY) {
    return deltaR(p.momentum(), j.momentum(), scheme);
  }

  inline double deltaR(const Jet& j, const FourMomentum& v,
                       RapScheme scheme = PSEUDORAPIDITY) {
    return deltaR(j.momentum(), v, scheme);
  }

  inline double deltaR(const Jet& j, const FourVector& v,
                       RapScheme scheme = PSEUDORAPIDITY) {
    return deltaR(j.momentum(), v, scheme);
  }

  inline double deltaR(const Jet& j, const Vector3& v) {
    return deltaR(j.momentum(), v);
  }

  inline double deltaR(const Jet& j, double eta, double phi) {
    return deltaR(j.momentum(), eta, phi);
  }

  inline double deltaR(const FourMomentum& v, const Jet& j,
                       RapScheme scheme = PSEUDORAPIDITY) {
    return deltaR(v, j.momentum(), scheme);
  }

  inline double deltaR(const FourVector& v, const Jet& j,
                       RapScheme scheme = PSEUDORAPIDITY) {
    return deltaR(v, j.momentum(), scheme);
  }

  inline double deltaR(const Vector3& v, const Jet& j) {
    return deltaR(v, j.momentum());
  }

  inline double deltaR(double eta, double phi, const Jet& j) {
    return deltaR(eta, phi, j.momentum());
  }


  inline double deltaPhi(const Jet& j1, const Jet& j2) {
    return deltaPhi(j1.momentum(), j2.momentum());
  }

  inline double deltaPhi(const Jet& j, const Particle& p) {
    return deltaPhi(j.momentum(), p.momentum());
  }

  inline double deltaPhi(const Particle& p, const Jet& j) {
    return deltaPhi(p.momentum(), j.momentum());
  }

  inline double deltaPhi(const Jet& j, const FourMomentum& v) {
    return deltaPhi(j.momentum(), v);
  }

  inline double deltaPhi(const Jet& j, const FourVector& v) {
    return deltaPhi(j.momentum(), v);
  }

  inline double deltaPhi(const Jet& j, const Vector3& v) {
    return deltaPhi(j.momentum(), v);
  }

  inline double deltaPhi(const Jet& j, double phi) {
    return deltaPhi(j.momentum(), phi);
  }

  inline double deltaPhi(const FourMomentum& v, const Jet& j) {
    return deltaPhi(v, j.momentum());
  }

  inline double deltaPhi(const FourVector& v, const Jet& j) {
    return deltaPhi(v, j.momentum());
  }

  inline double deltaPhi(const Vector3& v, const Jet& j) {
    return deltaPhi(v, j.momentum());
  }

  inline double deltaPhi(double phi, const Jet& j) {
    return deltaPhi(phi, j.momentum());
  }


  inline double deltaEta(const Jet& j1, const Jet& j2) {
    return deltaEta(j1.momentum(), j2.momentum());
  }

  inline double deltaEta(const Jet& j, const Particle& p) {
    return deltaEta(j.momentum(), p.momentum());
  }

  inline double deltaEta(const Particle& p, const Jet& j) {
    return deltaEta(p.momentum(), j.momentum());
  }

  inline double deltaEta(const Jet& j, const FourMomentum& v) {
    return deltaEta(j.momentum(), v);
  }

  inline double deltaEta(const Jet& j, const FourVector& v) {
    return deltaEta(j.momentum(), v);
  }

  inline double deltaEta(const Jet& j, const Vector3& v) {
    return deltaEta(j.momentum(), v);
  }

  inline double deltaEta(const Jet& j, double eta) {
    return deltaEta(j.momentum(), eta);
  }

  inline double deltaEta(const FourMomentum& v, const Jet& j) {
    return deltaEta(v, j.momentum());
  }

  inline double deltaEta(const FourVector& v, const Jet& j) {
    return deltaEta(v, j.momentum());
  }

  inline double deltaEta(const Vector3& v, const Jet& j) {
    return deltaEta(v, j.momentum());
  }

  inline double deltaEta(double eta, const Jet& j) {
    return deltaEta(eta, j.momentum());
  }

}

#endif