Rivet analyses referenceCMS_2021_I1972986Measurement and QCD analysis of double-differential inclusive jet cross sections in proton-proton collisions at 13 TeVExperiment: CMS (LHC) Inspire ID: 1972986 Status: VALIDATED Authors:
Beam energies: (6500.0, 6500.0) GeV Run details:
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#include "Rivet/Tools/BinnedHistogram.hh"
6
7
8namespace Rivet {
9
10
11 /// Inclusive jet pT at 13 TeV
12 class CMS_2021_I1972986 : public Analysis {
13 public:
14
15 /// Constructor
16 RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2021_I1972986);
17
18
19 /// Book histograms and initialize projections:
20 void init() {
21
22 // Initialize the projections
23 const FinalState fs;
24 declare(FastJets(fs, FastJets::ANTIKT, 0.4), "JetsAK4");
25 declare(FastJets(fs, FastJets::ANTIKT, 0.7), "JetsAK7");
26
27 // Book sets of histograms, binned in absolute rapidity
28 Histo1DPtr tmp;
29 for(int y = 0; y < 4; ++y) {
30 _hist_sigmaAK4.add(0.5*y, 0.5*(y+1), book(tmp,y+1, 1, 1)); // d0?-x01-y01
31 _hist_sigmaAK7.add(0.5*y, 0.5*(y+1), book(tmp,20+y+1, 1, 1)); // d2?-x01-y01
32 }
33 }
34
35
36 /// Per-event analysis
37 void analyze(const Event &event) {
38
39 // AK4 jets
40 const FastJets& fjAK4 = apply<FastJets>(event, "JetsAK4");
41 const Jets& jetsAK4 = fjAK4.jets(Cuts::ptIn(97*GeV, 3103*GeV) && Cuts::absrap < 2.0);
42 for (const Jet& j : jetsAK4) {
43 _hist_sigmaAK4.fill(j.absrap(), j.pT());
44 }
45
46 // AK7 jets
47 const FastJets& fjAK7 = apply<FastJets>(event, "JetsAK7");
48 const Jets& jetsAK7 = fjAK7.jets(Cuts::ptIn(97*GeV, 3103*GeV) && Cuts::absrap < 2.0);
49 for (const Jet& j : jetsAK7) {
50 _hist_sigmaAK7.fill(j.absrap(), j.pT());
51 }
52
53 }
54
55
56 // Finalize
57 void finalize() {
58 /// @note Extra division factor is the *signed* dy, i.e. 2*d|y|
59 _hist_sigmaAK4.scale(crossSection()/picobarn/sumOfWeights()/2.0, this);
60 _hist_sigmaAK7.scale(crossSection()/picobarn/sumOfWeights()/2.0, this);
61 }
62
63
64 /// @name Histograms
65 //@{
66 BinnedHistogram _hist_sigmaAK4;
67 BinnedHistogram _hist_sigmaAK7;
68 //@}
69
70 };
71
72
73 // This global object acts as a hook for the plugin system.
74 RIVET_DECLARE_PLUGIN(CMS_2021_I1972986);
75
76}
|