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
 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
// -*- 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_2021_I1972986 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2021_I1972986);


    /// 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
      Histo1DPtr tmp;
      for(int y = 0; y < 4; ++y) {
         _hist_sigmaAK4.add(0.5*y, 0.5*(y+1), book(tmp,y+1, 1, 1)); // d0?-x01-y01
         _hist_sigmaAK7.add(0.5*y, 0.5*(y+1), book(tmp,20+y+1, 1, 1)); // d2?-x01-y01
      }
    }


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

      // AK4 jets
      const FastJets& fjAK4 = apply<FastJets>(event, "JetsAK4");
      const Jets& jetsAK4 = fjAK4.jets(Cuts::ptIn(97*GeV, 3103*GeV) && Cuts::absrap < 2.0);
      for (const Jet& j : jetsAK4) {
        _hist_sigmaAK4.fill(j.absrap(), j.pT());
      }

      // AK7 jets
      const FastJets& fjAK7 = apply<FastJets>(event, "JetsAK7");
      const Jets& jetsAK7 = fjAK7.jets(Cuts::ptIn(97*GeV, 3103*GeV) && Cuts::absrap < 2.0);
      for (const Jet& j : jetsAK7) {
        _hist_sigmaAK7.fill(j.absrap(), j.pT());
      }

    }


    // 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);
    }


    /// @name Histograms
    //@{
    BinnedHistogram _hist_sigmaAK4;
    BinnedHistogram _hist_sigmaAK7;
    //@}

  };


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

}