Rivet analyses referenceCMS_2016_I1459051Measurement of the inclusive jet cross-section in $pp$ collisions at $\sqrt{s} = 13 \text{TeV}$Experiment: CMS (LHC) Inspire ID: 1459051 Status: VALIDATED Authors:
Beam energies: (6500.0, 6500.0) GeV Run details:
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// -*- 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
10 /// Inclusive jet pT at 13 TeV
11 class CMS_2016_I1459051 : public Analysis {
12 public:
13
14 /// Constructor
15 RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2016_I1459051);
16
17
18 /// Book histograms and initialize projections:
19 void init() {
20
21 // Initialize the projections
22 const FinalState fs;
23 declare(FastJets(fs, FastJets::ANTIKT, 0.4), "JetsAK4");
24 declare(FastJets(fs, FastJets::ANTIKT, 0.7), "JetsAK7");
25
26 // Book sets of histograms, binned in absolute rapidity
27 // AK7
28 {Histo1DPtr tmp; _hist_sigmaAK7.add(0.0, 0.5, book(tmp, 1, 1, 1));}
29 {Histo1DPtr tmp; _hist_sigmaAK7.add(0.5, 1.0, book(tmp, 2, 1, 1));}
30 {Histo1DPtr tmp; _hist_sigmaAK7.add(1.0, 1.5, book(tmp, 3, 1, 1));}
31 {Histo1DPtr tmp; _hist_sigmaAK7.add(1.5, 2.0, book(tmp, 4, 1, 1));}
32 {Histo1DPtr tmp; _hist_sigmaAK7.add(2.0, 2.5, book(tmp, 5, 1, 1));}
33 {Histo1DPtr tmp; _hist_sigmaAK7.add(2.5, 3.0, book(tmp, 6, 1, 1));}
34 book(_hist_sigmaAK7Forward, 7, 1, 1);
35 // AK4
36 {Histo1DPtr tmp; _hist_sigmaAK4.add(0.0, 0.5, book(tmp, 8, 1, 1));}
37 {Histo1DPtr tmp; _hist_sigmaAK4.add(0.5, 1.0, book(tmp, 9, 1, 1));}
38 {Histo1DPtr tmp; _hist_sigmaAK4.add(1.0, 1.5, book(tmp, 10, 1, 1));}
39 {Histo1DPtr tmp; _hist_sigmaAK4.add(1.5, 2.0, book(tmp, 11, 1, 1));}
40 {Histo1DPtr tmp; _hist_sigmaAK4.add(2.0, 2.5, book(tmp, 12, 1, 1));}
41 {Histo1DPtr tmp; _hist_sigmaAK4.add(2.5, 3.0, book(tmp, 13, 1, 1));}
42 book(_hist_sigmaAK4Forward, 14, 1, 1);
43
44 }
45
46
47 /// Per-event analysis
48 void analyze(const Event &event) {
49
50 const double weight = 1.0;
51
52 // AK4 jets
53 const FastJets& fjAK4 = applyProjection<FastJets>(event, "JetsAK4");
54 const Jets& jetsAK4 = fjAK4.jets(Cuts::ptIn(114*GeV, 2200.0*GeV) && Cuts::absrap < 4.7);
55 for (const Jet& j : jetsAK4) {
56 _hist_sigmaAK4.fill(j.absrap(), j.pT(), weight);
57 if (inRange(j.absrap(), 3.2, 4.7)) _hist_sigmaAK4Forward->fill(j.pT(), weight);
58 }
59
60 // AK7 jets
61 const FastJets& fjAK7 = applyProjection<FastJets>(event, "JetsAK7");
62 const Jets& jetsAK7 = fjAK7.jets(Cuts::ptIn(114*GeV, 2200.0*GeV) && Cuts::absrap < 4.7);
63 for (const Jet& j : jetsAK7) {
64 _hist_sigmaAK7.fill(j.absrap(), j.pT(), weight);
65 if (inRange(j.absrap(), 3.2, 4.7)) _hist_sigmaAK7Forward->fill(j.pT(), weight);
66 }
67
68 }
69
70
71 // Finalize
72 void finalize() {
73 /// @note Extra division factor is the *signed* dy, i.e. 2*d|y|
74 _hist_sigmaAK4.scale(crossSection()/picobarn/sumOfWeights()/2.0, this);
75 _hist_sigmaAK7.scale(crossSection()/picobarn/sumOfWeights()/2.0, this);
76 scale(_hist_sigmaAK4Forward,crossSection()/picobarn/sumOfWeights()/3.0);
77 scale(_hist_sigmaAK7Forward,crossSection()/picobarn/sumOfWeights()/3.0);
78 }
79
80
81 /// @name Histograms
82 //@{
83 BinnedHistogram _hist_sigmaAK4;
84 BinnedHistogram _hist_sigmaAK7;
85 Histo1DPtr _hist_sigmaAK4Forward;
86 Histo1DPtr _hist_sigmaAK7Forward;
87 //@}
88
89 };
90
91
92 // This global object acts as a hook for the plugin system.
93 RIVET_DECLARE_PLUGIN(CMS_2016_I1459051);
94
95}
|