rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2021_I1972986

Measurement and QCD analysis of double-differential inclusive jet cross sections in proton-proton collisions at 13 TeV
Experiment: CMS (LHC)
Inspire ID: 1972986
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Patrick L.S. Connor
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • pp to jets at $\sqrt{s} = 13$ TeV. Data collected by CMS during the year 2016.

A measurement of the inclusive jet production in proton-proton collisions at the LHC at √s = 13 TeV is presented. The double-differential cross sections are measured as a function of the jet transverse momentum pT and the absolute jet rapidity |y|. The anti-kT clustering algorithm is used with distance parameter of 0.4 (0.7) in a phase space region with jet pT from 97 GeV up to 3.1 TeV and |y|< 2.0. Data collected with the CMS detector are used, corresponding to an integrated luminosity of 36.3/fb (33.5/fb). The measurement is used in a comprehensive QCD analysis at next-to-next-to-leading order, which results in significant improvement in the accuracy of the parton distributions in the proton. Simultaneously, the value of the strong coupling constant at the Z boson mass is extracted as alpha_S(Z) = 0.1170 ± 0.0019. For the first time, these data are used in a standard model effective field theory analysis at next-to-leading order, where parton distributions and the QCD parameters are extracted simultaneously with imposed constraints on the Wilson coefficient c1 of 4-quark contact interactions.

Source code: CMS_2021_I1972986.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 pT at 13 TeV
10  class CMS_2021_I1972986 : public Analysis {
11  public:
12
13    /// Constructor
14    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2021_I1972986);
15
16
17    /// Book histograms and initialize projections:
18    void init() {
19
20      // Initialize the projections
21      const FinalState fs;
22      declare(FastJets(fs, JetAlg::ANTIKT, 0.4), "JetsAK4");
23      declare(FastJets(fs, JetAlg::ANTIKT, 0.7), "JetsAK7");
24
25      // Book sets of histograms, binned in absolute rapidity
26      book(_hist_sigmaAK4, {0., 0.5, 1.0, 1.5, 2.0},
27                           {"d01-x01-y01", "d02-x01-y01", "d03-x01-y01", "d04-x01-y01"});
28      book(_hist_sigmaAK7, {0., 0.5, 1.0, 1.5, 2.0},
29                           {"d21-x01-y01", "d22-x01-y01", "d23-x01-y01", "d24-x01-y01"});
30    }
31
32
33    /// Per-event analysis
34    void analyze(const Event &event) {
35
36      // AK4 jets
37      const FastJets& fjAK4 = apply<FastJets>(event, "JetsAK4");
38      const Jets& jetsAK4 = fjAK4.jets(Cuts::ptIn(97*GeV, 3103*GeV) && Cuts::absrap < 2.0);
39      for (const Jet& j : jetsAK4) {
40        _hist_sigmaAK4->fill(j.absrap(), j.pT());
41      }
42
43      // AK7 jets
44      const FastJets& fjAK7 = apply<FastJets>(event, "JetsAK7");
45      const Jets& jetsAK7 = fjAK7.jets(Cuts::ptIn(97*GeV, 3103*GeV) && Cuts::absrap < 2.0);
46      for (const Jet& j : jetsAK7) {
47        _hist_sigmaAK7->fill(j.absrap(), j.pT());
48      }
49
50    }
51
52
53    // Finalize
54    void finalize() {
55      /// @note Extra division factor is the *signed* dy, i.e. 2*d|y|
56      scale(_hist_sigmaAK4, 0.5*crossSection()/picobarn/sumOfWeights());
57      scale(_hist_sigmaAK7, 0.5*crossSection()/picobarn/sumOfWeights());
58      divByGroupWidth({_hist_sigmaAK4, _hist_sigmaAK7});
59    }
60
61
62    /// @name Histograms
63    /// @{
64    Histo1DGroupPtr _hist_sigmaAK4;
65    Histo1DGroupPtr _hist_sigmaAK7;
66    /// @}
67
68  };
69
70
71  RIVET_DECLARE_PLUGIN(CMS_2021_I1972986);
72
73}