rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2022_I2170533

Measurements of jet multiplicity and jet transverse momentum in multijet events in proton-proton collisions at \sqrt{s} = 13 TeV
Experiment: CMS (LHC)
Inspire ID: 2170533
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Luis Ignacio Estevez Banos
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.

Multijet events at large transverse momentum (pT) are measured at $\sqrt{s}$ = 13 TeV using data recorded with the CMS detector at the LHC, corresponding to an integrated luminosity of 36.3 fb-1. The multiplicity of jets with pT > 50 GeV that are produced in association with a high-pT dijet system is measured in various ranges of the pT of the jet with the highest transverse momentum and as a function of the azimuthal angle difference between the two highest pT jets in the dijet system. The differential production cross sections are measured as a function of the transverse momenta of the four highest pT jets. The measurements are compared with leading and next-to-leading order matrix element calculations supplemented with simulations of parton shower, hadronization, and multiparton interactions. In addition, the measurements are compared with next-to-leading order matrix element calculations combined with transverse-momentum dependent parton densities and transverse-momentum dependent parton shower.

Source code: CMS_2022_I2170533.cc
 1// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Projections/FastJets.hh"
 4#include "Rivet/Projections/FinalState.hh"
 5
 6namespace Rivet {
 7
 8  /// @brief Multijet events at 13 TeV
 9  class CMS_2022_I2170533 : public Analysis {
10
11     public:
12
13      /// Constructor
14      RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2022_I2170533);
15
16      /// @name Analysis methods
17      ///@{
18
19      /// Book histograms and initialise projections before the run
20      void init() {
21          FinalState fs;
22          FastJets akt(fs, JetAlg::ANTIKT, 0.4);
23          declare(akt, "antikT");
24
25          vector<double> edges = {0., 150., 170., 180.};
26          book(_h_Mult_ptmax200, edges);
27          book(_h_Mult_ptmax400, edges);
28          book(_h_Mult_ptmax800, edges);
29          for (size_t i = 0; i < _h_Mult_ptmax200->numBins(); ++i) {
30            book(_h_Mult_ptmax200->bin(i+1), i+1, 1, 1);
31            book(_h_Mult_ptmax400->bin(i+1), i+4, 1, 1);
32            book(_h_Mult_ptmax800->bin(i+1), i+7, 1, 1);
33          }
34
35          // Pt of the first 4 jets
36          book(_h_pT1_n1, 10, 1, 1);
37          book(_h_pT2_n2, 11, 1, 1);
38          book(_h_pT3_n3, 12, 1, 1);
39          book(_h_pT4_n4, 13, 1, 1);
40      }
41
42      /// Perform the per-event analysis
43      void analyze(const Event& event) {
44          // Preselection cuts for the jets
45          const Jets& jets = apply<FastJets>(event, "antikT").jetsByPt(Cuts::pT > 20 * GeV && Cuts::absrap < 3.2);  // rapidity and pt preselection |y| < 3.2 pT > 20 GeV
46
47          int njet = 0;  // jet counting
48
49          // cuts after preselection
50          if (jets.size() < 2) vetoEvent;                                             // dijet cut
51          if ((jets[0].pT() < 200.) or (jets[1].pT() < 100.)) vetoEvent;              // pt cut on 2 leading jets
52          if ((fabs(jets[0].rap()) > 2.5) or (fabs(jets[1].rap()) > 2.5)) vetoEvent;  // |y| < 2.5 cut on leading jets
53
54          double dphi = deltaPhi(jets[0].phi(), jets[1].phi()) / degree;
55
56          for (const Jet& j : jets) {
57              if (j.pT() > 50. && fabs(j.rap()) < 2.5) njet = njet + 1;  // Cuts on the extrajets
58          }
59
60          if (njet > 7) njet = 7;  // Last bin in Multiplicity is inclusive
61
62          if (jets[0].pT() > 200 && jets[0].pT() <= 400) _h_Mult_ptmax200->fill(dphi, njet);
63          if (jets[0].pT() > 400 && jets[0].pT() <= 800) _h_Mult_ptmax400->fill(dphi, njet);
64          if (jets[0].pT() > 800) _h_Mult_ptmax800->fill(dphi, njet);
65
66          // Filling dijet 1,2 pT (The events are already selected with leadin jetpT > 200 GeV and 2nd jet pT > 100 GeV)
67          _h_pT1_n1->fill(jets[0].pT());
68          _h_pT2_n2->fill(jets[1].pT());
69
70          // Filling extra jets pT (jet3 and jet4)
71          if (njet > 2) _h_pT3_n3->fill(jets[2].pT());
72          if (njet > 3) _h_pT4_n4->fill(jets[3].pT());
73      }
74
75      /// Normalise histograms etc., after the run
76      void finalize() {
77          const double sf = crossSection() / picobarn / sumW();
78          scale(_h_Mult_ptmax200, sf);
79          scale(_h_Mult_ptmax400, sf);
80          scale(_h_Mult_ptmax800, sf);
81          scale({_h_pT1_n1, _h_pT2_n2, _h_pT3_n3, _h_pT4_n4}, sf);
82      }
83
84      ///@}
85
86      /// @name Histograms
87      ///@{
88
89      Histo1DGroupPtr _h_Mult_ptmax200, _h_Mult_ptmax400, _h_Mult_ptmax800;
90      Histo1DPtr _h_pT1_n1, _h_pT2_n2, _h_pT3_n3, _h_pT4_n4;
91
92      ///@}
93  };
94
95  RIVET_DECLARE_PLUGIN(CMS_2022_I2170533);
96
97}  // namespace Rivet