Processing math: 100%
rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2019_I1719955

Azimuthal separation in nearly back-to-back jet topologies in inclusive 2- and 3-jet events in pp collisions at s=13 TeV
Experiment: CMS (LHC)
Status: VALIDATED
Authors:
  • Armando Bermudez Martinez
  • SMP conveners
References:
  • 10.1140/epjc/s10052-019-7276-4
  • arXiv: 1902.04374
  • CERN-EP-2018-344
Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • QCD at s=13 TeV, ptHat (or equivalent) greater than 100 GeV for the 2-jet observables and greater than 10 GeV for the 3-jet observables

A measurement for inclusive 2- and 3-jet events of the azimuthal correlation between the two jets with the largest transverse momenta, Δϕ12, is presented. The measurement considers events where the two leading jets are nearly collinear ("back-to-back") in the transverse plane and is performed for several ranges of the leading jet transverse momentum. Proton-proton collision data collected with the CMS experiment at a center-of-mass energy of 13 TeV and corresponding to an integrated luminosity of 35.9 fb1 are used.

Source code: CMS_2019_I1719955.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  /// CMS azimuthal decorrelations in back-to-back dijet events at 13 TeV
 9  class CMS_2019_I1719955 : public Analysis {
10  public:
11    /// Constructor
12    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2019_I1719955);
13
14
15    /// Book projections and histograms
16    void init() {
17      const FinalState fs;
18      declare(FastJets(fs, JetAlg::ANTIKT, 0.4), "ANTIKT");
19
20      book(_h_deltaPhi_2J, {200., 300., 400., 500., 600., 700., 800., 1000., 1200., 4000.});
21      book(_h_deltaPhi_3J, {200., 300., 400., 500., 600., 700., 800., 1000., 1200., 4000.});
22      for (size_t i=1; i<_h_deltaPhi_2J->numBins()+1; ++i) {
23        book(_h_deltaPhi_2J->bin(i), i, 1, 1);
24        book(_h_deltaPhi_3J->bin(i), i+9, 1, 1);
25      }
26    }
27
28    /// Per-event analysis
29    void analyze(const Event& event) {
30      const Jets& jets = apply<JetFinder>(event, "ANTIKT").jetsByPt(Cuts::absrap < 5. && Cuts::pT > 100*GeV);
31      const Jets& lowjets = apply<JetFinder>(event, "ANTIKT").jetsByPt(Cuts::absrap < 2.5 && Cuts::pT > 30*GeV);
32      if (jets.size() < 2) vetoEvent;
33      if (jets[0].absrap() > 2.5 || jets[1].absrap() > 2.5) vetoEvent;
34
35      const double dphi = 180./M_PI*deltaPhi(jets[0].phi(), jets[1].phi());
36      _h_deltaPhi_2J->fill(jets[0].pT(), dphi);
37      if (lowjets.size() > 2) _h_deltaPhi_3J->fill(jets[0].pT(), dphi);
38    }
39
40    /// Scale histograms
41    void finalize() {
42      int region_ptmax_2J = 0;
43      double norm_finalize[9];
44      for (auto& histo_2J : _h_deltaPhi_2J->bins()) {
45        norm_finalize[region_ptmax_2J] = histo_2J->integral();
46        if (norm_finalize[region_ptmax_2J] != 0) scale(histo_2J, 1.0/norm_finalize[region_ptmax_2J]);
47        region_ptmax_2J++;
48      }
49      int region_ptmax_3J = 0;
50      for (auto& histo_3J : _h_deltaPhi_3J->bins()) {
51        if (norm_finalize[region_ptmax_3J] != 0) scale(histo_3J, 1.0/norm_finalize[region_ptmax_3J]);
52        region_ptmax_3J++;
53      }
54    }
55
56
57  private:
58
59    Histo1DGroupPtr _h_deltaPhi_2J, _h_deltaPhi_3J;
60
61  };
62
63
64  RIVET_DECLARE_PLUGIN(CMS_2019_I1719955);
65
66}