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