Rivet analyses referenceCMS_2016_I1487277Measurement 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 TeVExperiment: CMS (LHC) Inspire ID: 1487277 Status: VALIDATED Authors:
Beam energies: (4000.0, 4000.0) GeV Run details:
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#include "Rivet/Tools/BinnedHistogram.hh"
6
7namespace Rivet {
8
9 // Inclusive jet pT
10 class CMS_2016_I1487277 : public Analysis {
11 public:
12
13 // Constructor
14 CMS_2016_I1487277() : Analysis("CMS_2016_I1487277") {}
15
16
17 // Book histograms and initialize projections:
18 void init() {
19 const FinalState fs;
20
21 // Initialize the projectors:
22 declare(FastJets(fs, FastJets::ANTIKT, 0.7),"Jets");
23
24 // Book histograms:
25
26
27 {Histo1DPtr tmp; _hist_sigma.add(0.0, 0.5, book(tmp, 1, 1, 1));}
28 {Histo1DPtr tmp; _hist_sigma.add(0.5, 1.0, book(tmp, 2, 1, 1));}
29 {Histo1DPtr tmp; _hist_sigma.add(1.0, 1.5, book(tmp, 3, 1, 1));}
30 {Histo1DPtr tmp; _hist_sigma.add(1.5, 2.0, book(tmp, 4, 1, 1));}
31 {Histo1DPtr tmp; _hist_sigma.add(2.0, 2.5, book(tmp, 5, 1, 1));}
32 {Histo1DPtr tmp; _hist_sigma.add(2.5, 3.0, book(tmp, 6, 1, 1));}
33 {Histo1DPtr tmp; _hist_sigma.add(3.2, 4.7, book(tmp, 7, 1, 1));}
34
35 }
36
37 // Analysis
38 void analyze(const Event &event) {
39 const FastJets &fj = applyProjection<FastJets>(event,"Jets");
40 const Jets& jets = fj.jets(Cuts::ptIn(18*GeV, 5000.0*GeV) && Cuts::absrap < 5.2);
41
42 // Fill the relevant histograms:
43 for(const Jet &j : jets) {
44 _hist_sigma.fill(j.absrap(), j.pT());
45 }
46 }
47
48 // Finalize
49 void finalize() {
50 _hist_sigma.scale(crossSection()/sumOfWeights()/2.0, this);
51 }
52
53 private:
54 BinnedHistogram _hist_sigma;
55 Histo1DPtr _hist_ptbins_y1;
56 Histo1DPtr _hist_ptbins_y2;
57 Histo1DPtr _hist_ptbins_y3;
58 Histo1DPtr _hist_ptbins_y4;
59 Histo1DPtr _hist_ptbins_y5;
60 Histo1DPtr _hist_ptbins_y6;
61 Histo1DPtr _hist_ptbins_y7;
62
63 };
64
65 // This global object acts as a hook for the plugin system.
66 RIVET_DECLARE_PLUGIN(CMS_2016_I1487277);
67
68}
|