rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2016_CONF_2016_092

Inclusive jet cross sections using early 13 TeV data
Experiment: ATLAS (LHC)
Status: OBSOLETE
Authors:
  • Stefan von Buddenbrock
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • Soft and hard QCD events

This note presents the measurement of inclusive-jet cross-sections in proton--proton collisions at a centre-of-mass energy of 13TeV. The measurement uses data corresponding to an integrated luminosity of 3.2 fb$^{-1}$ collected with the ATLAS detector at the Large Hadron Collider in 2015. Jets are clustered using the anti-$k_\text{t}$ algorithm with a radius parameter value of $R = 0.4$. Double differential inclusive-jet cross-sections are presented as a function of the jet transverse momentum, covering the range from 100 GeV to about 3.2 TeV, and of the absolute jet rapidity, up to $|y| = 3$. The predictions of next-to-leading-order QCD calculations using several parton distribution function (PDF) sets and corrected for non-perturbative and electroweak effects are compared to the measured cross-sections. The predictions are consistent with the measured cross-sections within uncertainties.

Source code: ATLAS_2016_CONF_2016_092.cc
 1// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Projections/FinalState.hh"
 4#include "Rivet/Projections/FastJets.hh"
 5
 6namespace Rivet {
 7
 8
 9  /// Inclusive jet cross sections using early 13 TeV data
10  /// @author Stefan von Buddenbrock <stef.von.b@cern.ch>
11  class ATLAS_2016_CONF_2016_092 : public Analysis {
12  public:
13
14    /// Constructor
15    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_CONF_2016_092);
16
17
18    /// Bookings
19    void init() {
20
21      // Define the jets
22      FastJets antiKT04Jets(FinalState(Cuts::open()),  JetAlg::ANTIKT, 0.4);
23      antiKT04Jets.useInvisibles();
24      declare(antiKT04Jets, "antiKT04Jets");
25
26      // Book histograms
27      const vector<double> y_bins{0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0};
28      book(_h_pT, y_bins);
29      for (auto& b : _h_pT->bins()) {
30        book(b, b.index(), 1, 1);
31      }
32    }
33
34
35    /// Per-event analysis
36    void analyze(const Event& event) {
37      const Jets& jets = apply<FastJets>(event, "antiKT04Jets").jetsByPt(Cuts::pT > 100*GeV && Cuts::absrap < 3.0);
38      for (const Jet& j : jets) {
39        _h_pT->fill(j.absrap(), j.pT()/GeV);
40      }
41    }
42
43
44    /// Post-run scaling
45    void finalize() {
46      // Divide by 2 to only get positive rapidity values
47      scale(_h_pT, 0.5*crossSection()/picobarn/sumOfWeights());
48    }
49
50
51    /// Histograms
52    Histo1DGroupPtr _h_pT;
53
54  };
55
56
57  RIVET_DECLARE_PLUGIN(ATLAS_2016_CONF_2016_092);
58
59}