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:
  • ATLAS-CONF-2016-092
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
 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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Tools/BinnedHistogram.hh"

namespace Rivet {


  /// Inclusive jet cross sections using early 13 TeV data
  /// @author Stefan von Buddenbrock <stef.von.b@cern.ch>
  class ATLAS_2016_CONF_2016_092 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_CONF_2016_092);


    /// Bookings
    void init() {

      // Define the jets
      FastJets antiKT04Jets(FinalState(Cuts::open()),  FastJets::ANTIKT, 0.4);
      antiKT04Jets.useInvisibles();
      declare(antiKT04Jets, "antiKT04Jets");

      // Book histograms
      const double y_bins[] = {0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0};
      for (size_t i = 0; i < 6; i++) {
        Histo1DPtr tmp;
        _h_pT.add(y_bins[i], y_bins[i+1], book(tmp, i+1, 1, 1));
      }
    }


    /// Per-event analysis
    void analyze(const Event& event) {
      const Jets& jets = apply<FastJets>(event, "antiKT04Jets")
        .jetsByPt(Cuts::pT > 100*GeV && Cuts::absrap < 3.0);
      for (const Jet& j : jets)
        _h_pT.fill(j.absrap(), j.pT()/GeV, 1.0);
    }


    /// Post-run scaling
    void finalize() {
      // Divide by 2 to only get positive rapidity values
      _h_pT.scale(0.5*crossSection()/picobarn/sumOfWeights(), this);
    }


    /// Histograms
    BinnedHistogram _h_pT;

  };


  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2016_CONF_2016_092);

}