Rivet Analyses Reference

CMS_2016_I1459051

Measurement of the inclusive jet cross-section in $pp$ collisions at $\sqrt{s} = 13 \text{TeV}$
Experiment: CMS (LHC)
Inspire ID: 1459051
Status: VALIDATED
Authors:
  • Paolo Gunnellini
References:
  • Eur.Phys.J. C76 (2016) no.8, 451
  • CERN-EP-2016-104
  • CMS-SMP-15-007
Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • Inclusive hard QCD at 13 TeV center-of-mass energy, and ptHat (or equivalent) greater than 90 GeV

A measurement of the double-differential inclusive jet cross section as a function of jet transverse momentum pT and absolute jet rapidity |y| is presented. The analysis is based on proton-proton collisions collected by the CMS experiment at the LHC at a centre-of-mass energy of 13 TeV. The data samples correspond to integrated luminosities of 71 and 44 pb-1 for |y| < 3 and 3.2 < |y| < 4.7, respectively. Jets are reconstructed with the anti-kt clustering algorithm for two jet sizes, R, of 0.7 and 0.4, in a phase space region covering jet pT up to 2 TeV and jet rapidity up to |y| = 4.7. Predictions of perturbative quantum chromodynamics at next-to-leading order precision, complemented with electroweak and nonperturbative corrections, are used to compute the absolute scale and the shape of the inclusive jet cross section. The cross-section difference in $R$, when going to a smaller jet size of 0.4, is best described by Monte Carlo event generators with next-to-leading order predictions matched to parton showering, hadronisation, and multiparton interactions. In the phase space accessible with the new data, this measurement provides a first indication that jet physics is as well understood at $\sqrt(s) = 13 \text{TeV}$ as at smaller centre-of-mass energies.

Source code: CMS_2016_I1459051.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
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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Tools/BinnedHistogram.hh"

namespace Rivet {


  /// Inclusive jet pT at 13 TeV
  class CMS_2016_I1459051 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2016_I1459051);


    /// Book histograms and initialize projections:
    void init() {

      // Initialize the projections
      const FinalState fs;
      declare(FastJets(fs, FastJets::ANTIKT, 0.4), "JetsAK4");
      declare(FastJets(fs, FastJets::ANTIKT, 0.7), "JetsAK7");

      // Book sets of histograms, binned in absolute rapidity
      // AK7
      {Histo1DPtr tmp; _hist_sigmaAK7.add(0.0, 0.5, book(tmp, 1, 1, 1));}
      {Histo1DPtr tmp; _hist_sigmaAK7.add(0.5, 1.0, book(tmp, 2, 1, 1));}
      {Histo1DPtr tmp; _hist_sigmaAK7.add(1.0, 1.5, book(tmp, 3, 1, 1));}
      {Histo1DPtr tmp; _hist_sigmaAK7.add(1.5, 2.0, book(tmp, 4, 1, 1));}
      {Histo1DPtr tmp; _hist_sigmaAK7.add(2.0, 2.5, book(tmp, 5, 1, 1));}
      {Histo1DPtr tmp; _hist_sigmaAK7.add(2.5, 3.0, book(tmp, 6, 1, 1));}
      book(_hist_sigmaAK7Forward, 7, 1, 1);
      // AK4
      {Histo1DPtr tmp; _hist_sigmaAK4.add(0.0, 0.5, book(tmp, 8, 1, 1));}
      {Histo1DPtr tmp; _hist_sigmaAK4.add(0.5, 1.0, book(tmp, 9, 1, 1));}
      {Histo1DPtr tmp; _hist_sigmaAK4.add(1.0, 1.5, book(tmp, 10, 1, 1));}
      {Histo1DPtr tmp; _hist_sigmaAK4.add(1.5, 2.0, book(tmp, 11, 1, 1));}
      {Histo1DPtr tmp; _hist_sigmaAK4.add(2.0, 2.5, book(tmp, 12, 1, 1));}
      {Histo1DPtr tmp; _hist_sigmaAK4.add(2.5, 3.0, book(tmp, 13, 1, 1));}
      book(_hist_sigmaAK4Forward, 14, 1, 1);

    }


    /// Per-event analysis
    void analyze(const Event &event) {

      const double weight = 1.0;

      // AK4 jets
      const FastJets& fjAK4 = applyProjection<FastJets>(event, "JetsAK4");
      const Jets& jetsAK4 = fjAK4.jets(Cuts::ptIn(114*GeV, 2200.0*GeV) && Cuts::absrap < 4.7);
      for (const Jet& j : jetsAK4) {
        _hist_sigmaAK4.fill(j.absrap(), j.pT(), weight);
        if (inRange(j.absrap(), 3.2, 4.7)) _hist_sigmaAK4Forward->fill(j.pT(), weight);
      }

      // AK7 jets
      const FastJets& fjAK7 = applyProjection<FastJets>(event, "JetsAK7");
      const Jets& jetsAK7 = fjAK7.jets(Cuts::ptIn(114*GeV, 2200.0*GeV) && Cuts::absrap < 4.7);
      for (const Jet& j : jetsAK7) {
        _hist_sigmaAK7.fill(j.absrap(), j.pT(), weight);
        if (inRange(j.absrap(), 3.2, 4.7)) _hist_sigmaAK7Forward->fill(j.pT(), weight);
      }

    }


    // Finalize
    void finalize() {
      /// @note Extra division factor is the *signed* dy, i.e. 2*d|y|
      _hist_sigmaAK4.scale(crossSection()/picobarn/sumOfWeights()/2.0, this);
      _hist_sigmaAK7.scale(crossSection()/picobarn/sumOfWeights()/2.0, this);
      scale(_hist_sigmaAK4Forward,crossSection()/picobarn/sumOfWeights()/3.0);
      scale(_hist_sigmaAK7Forward,crossSection()/picobarn/sumOfWeights()/3.0);
    }


    /// @name Histograms
    //@{
    BinnedHistogram _hist_sigmaAK4;
    BinnedHistogram _hist_sigmaAK7;
    Histo1DPtr _hist_sigmaAK4Forward;
    Histo1DPtr _hist_sigmaAK7Forward;
    //@}

  };


  // This global object acts as a hook for the plugin system.
  RIVET_DECLARE_PLUGIN(CMS_2016_I1459051);

}