rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2016_I1487277

Measurement and QCD analysis of double-differential inclusive jet cross sections in $pp$ collisions at $\sqrt{s} = 8 \text{TeV}$ and cross-section ratios to 2.76 and 7 TeV
Experiment: CMS (LHC)
Inspire ID: 1487277
Status: VALIDATED
Authors:
  • Hannes Jung
References: Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • QCD at $\sqrt{s} = 8 \text{TeV}$

A measurement of the double-differential inclusive jet cross section as a function of the jet transverse momentum $\pt$ and the absolute jet rapidity $\abs{y}$ is presented. Data from LHC proton-proton collisions at $\sqrt{s}=8\text{TeV}$, corresponding to an integrated luminosity of 19.7 fb^{-1}, have been collected with the CMS detector. Jets are reconstructed using the anti-$\kt$ clustering algorithm with a size parameter of 0.7 in a phase space region covering jet $\pt$ from $74\text{GeV}$ up to $2.5\text{TeV}$ and jet absolute rapidity up to $\abs{y}=3.0$. The low-$\pt$ jet range between 21 and $74\text{GeV}$ is also studied up to $\abs{y}=4.7$, using a dedicated data sample corresponding to an integrated luminosity of 5.6\pbinv. The measured jet cross section is corrected for detector effects and compared with the predictions from perturbative QCD at next-to-leading order (NLO) using various sets of parton distribution functions (PDF). Cross section ratios to the corresponding measurements performed at 2.76 and $7\text{TeV}$ are presented. From the measured double-differential jet cross section, the value of the strong coupling constant evaluated at the $Z$ mass is $\alpha_\mathrm{S}(M_{Z}) = 0.1164^{+0.0060}_{-0.0043}$, where the errors include the PDF, scale, nonperturbative effects and experimental uncertainties, using the CT10 NLO PDFs. Improved constraints on PDFs based on the inclusive jet cross section measurement are presented.

Source code: CMS_2016_I1487277.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  // Inclusive jet pT
 9  class CMS_2016_I1487277 : public Analysis {
10  public:
11
12    // Constructor
13    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2016_I1487277);
14
15
16    // Book histograms and initialize projections:
17    void init() {
18      const FinalState fs;
19
20      // Initialize the projectors:
21      declare(FastJets(fs, JetAlg::ANTIKT, 0.7),"Jets");
22
23      // Book histograms:
24      book(_hist_sigma, {0., 0.5, 1., 1.5, 2., 2.5, 3., 4.7});
25      for (auto& b : _hist_sigma->bins()) {
26        book(b, b.index(), 1, 1);
27      }
28
29    }
30
31    // Analysis
32    void analyze(const Event &event) {
33      const FastJets &fj = apply<FastJets>(event,"Jets");
34      const Jets& jets = fj.jets(Cuts::ptIn(18*GeV, 5000.0*GeV) && Cuts::absrap < 5.2);
35
36      // Fill the relevant histograms:
37      for (const Jet &j : jets) {
38        _hist_sigma->fill(j.absrap(), j.pT());
39      }
40    }
41
42    // Finalize
43    void finalize() {
44      scale(_hist_sigma, crossSection()/picobarn/sumOfWeights()/2.0);
45    }
46
47  private:
48
49    Histo1DGroupPtr _hist_sigma;
50
51  };
52
53  RIVET_DECLARE_PLUGIN(CMS_2016_I1487277);
54
55}