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
6namespace Rivet {
7
8
9 /// Inclusive jet pT at 13 TeV
10 class CMS_2016_I1459051 : public Analysis {
11 public:
12
13 /// Constructor
14 RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2016_I1459051);
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_sigmaAK7, {0., 0.5, 1., 1.5, 2., 2.5, 3.});
27 book(_hist_sigmaAK4, {0., 0.5, 1., 1.5, 2., 2.5, 3.});
28 for (size_t i=1; i < _hist_sigmaAK7->numBins()+1; ++i) {
29 book(_hist_sigmaAK7->bin(i), i, 1, 1);
30 book(_hist_sigmaAK4->bin(i), 7+i, 1, 1);
31 }
32 book(_hist_sigmaAK7Forward, 7, 1, 1);
33 book(_hist_sigmaAK4Forward, 14, 1, 1);
34
35 }
36
37
38 /// Per-event analysis
39 void analyze(const Event &event) {
40
41 // AK4 jets
42 const FastJets& fjAK4 = apply<FastJets>(event, "JetsAK4");
43 const Jets& jetsAK4 = fjAK4.jets(Cuts::ptIn(114*GeV, 2200.0*GeV) && Cuts::absrap < 4.7);
44 for (const Jet& j : jetsAK4) {
45 _hist_sigmaAK4->fill(j.absrap(), j.pT());
46 if (inRange(j.absrap(), 3.2, 4.7)) _hist_sigmaAK4Forward->fill(j.pT());
47 }
48
49 // AK7 jets
50 const FastJets& fjAK7 = apply<FastJets>(event, "JetsAK7");
51 const Jets& jetsAK7 = fjAK7.jets(Cuts::ptIn(114*GeV, 2200.0*GeV) && Cuts::absrap < 4.7);
52 for (const Jet& j : jetsAK7) {
53 _hist_sigmaAK7->fill(j.absrap(), j.pT());
54 if (inRange(j.absrap(), 3.2, 4.7)) _hist_sigmaAK7Forward->fill(j.pT());
55 }
56
57 }
58
59
60 // Finalize
61 void finalize() {
62 /// @note Extra division factor is the *signed* dy, i.e. 2*d|y|
63 scale(_hist_sigmaAK4, crossSection()/picobarn/sumOfWeights()/2.0);
64 scale(_hist_sigmaAK7, crossSection()/picobarn/sumOfWeights()/2.0);
65 scale(_hist_sigmaAK4Forward,crossSection()/picobarn/sumOfWeights()/3.0);
66 scale(_hist_sigmaAK7Forward,crossSection()/picobarn/sumOfWeights()/3.0);
67 divByGroupWidth({_hist_sigmaAK4, _hist_sigmaAK7});
68 }
69
70
71 /// @name Histograms
72 /// @{
73 Histo1DGroupPtr _hist_sigmaAK4;
74 Histo1DGroupPtr _hist_sigmaAK7;
75 Histo1DPtr _hist_sigmaAK4Forward;
76 Histo1DPtr _hist_sigmaAK7Forward;
77 /// @}
78
79 };
80
81
82 RIVET_DECLARE_PLUGIN(CMS_2016_I1459051);
83
84}
|