rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2018_I1643640

Azimuthal correlations for inclusive 2-jet, 3-jet, and 4-jet events in pp collisions at $\sqrt{s}$ = 13 TeV
Experiment: CMS (LHC)
Inspire ID: 1643640
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Hannes Jung
  • Paolo Gunnellini
  • Panos Kokkas
References:
  • 10.1140/epjc/s10052-018-6033-4
  • arxiv:1712.05471
  • CMS-SMP-16-014
Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • QCD at $\sqrt{s} = 13 \text{TeV}$

Azimuthal correlations between the two jets with the largest transverse momenta $p_T$ in inclusive 2-, 3-, and 4-jet events are presented for several regions of the leading jet $p_T$ up to 4 TeV. For 3- and 4-jet scenarios, measurements of the minimum azimuthal angles between any two of the three or four leading $p_T$ jets are also presented. The analysis is based on data from proton-proton collisions collected by the CMS Collaboration at a centre-of-mass energy of 13 TeV, corresponding to an integrated luminosity of 35.9 fb$^{-1}$. Calculations based on leading-order matrix elements supplemented with parton showering and hadronization do not fully describe the data, so next-to-leading-order calculations matched with parton shower and hadronization models are needed to better describe the measured distributions. Furthermore, we show that azimuthal jet correlations are sensitive to details of the parton showering, hadronization, and multiparton interactions. A next-to-leading-order calculation matched with parton showers in the MC\@NLO method, as implemented in HERWIG 7, gives a better overall description of the measurements than the POWHEG method.

Source code: CMS_2018_I1643640.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Tools/BinnedHistogram.hh"
  4#include "Rivet/Projections/FinalState.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6
  7namespace Rivet {
  8
  9
 10  /// CMS Azimuthal corellations at 13 TeV
 11  class CMS_2018_I1643640 : public Analysis {
 12  public:
 13
 14    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2018_I1643640);
 15
 16    void init() {
 17      FinalState fs;
 18      FastJets akt(fs, FastJets::ANTIKT, 0.4);
 19      declare(akt, "antikT");
 20
 21      Histo1DPtr tmp;
 22
 23      _h_deltaPhi_2J_phi12.add( 200.,  300., book(tmp, 1, 1, 1));
 24      _h_deltaPhi_2J_phi12.add( 300.,  400., book(tmp, 2, 1, 1));
 25      _h_deltaPhi_2J_phi12.add( 400.,  500., book(tmp, 3, 1, 1));
 26      _h_deltaPhi_2J_phi12.add( 500.,  600., book(tmp, 4, 1, 1));
 27      _h_deltaPhi_2J_phi12.add( 600.,  700., book(tmp, 5, 1, 1));
 28      _h_deltaPhi_2J_phi12.add( 700.,  800., book(tmp, 6, 1, 1));
 29      _h_deltaPhi_2J_phi12.add( 800.,  1000., book(tmp, 7, 1, 1));
 30      _h_deltaPhi_2J_phi12.add( 1000., 1200., book(tmp, 8, 1, 1));
 31      _h_deltaPhi_2J_phi12.add( 1200., 7000., book(tmp, 9, 1, 1));
 32
 33      _h_deltaPhi_3J_phi12.add( 200.,  300., book(tmp, 10, 1, 1));
 34      _h_deltaPhi_3J_phi12.add( 300.,  400., book(tmp, 11, 1, 1));
 35      _h_deltaPhi_3J_phi12.add( 400.,  500., book(tmp, 12, 1, 1));
 36      _h_deltaPhi_3J_phi12.add( 500.,  600., book(tmp, 13, 1, 1));
 37      _h_deltaPhi_3J_phi12.add( 600.,  700., book(tmp, 14, 1, 1));
 38      _h_deltaPhi_3J_phi12.add( 700.,  800., book(tmp, 15, 1, 1));
 39      _h_deltaPhi_3J_phi12.add( 800.,  1000., book(tmp, 16, 1, 1));
 40      _h_deltaPhi_3J_phi12.add( 1000., 7000., book(tmp, 17, 1, 1));
 41
 42      _h_deltaPhi_4J_phi12.add( 200.,  300., book(tmp, 18, 1, 1));
 43      _h_deltaPhi_4J_phi12.add( 300.,  400., book(tmp, 19, 1, 1));
 44      _h_deltaPhi_4J_phi12.add( 400.,  500., book(tmp, 20, 1, 1));
 45      _h_deltaPhi_4J_phi12.add( 500.,  600., book(tmp, 21, 1, 1));
 46      _h_deltaPhi_4J_phi12.add( 600.,  700., book(tmp, 22, 1, 1));
 47      _h_deltaPhi_4J_phi12.add( 700.,  800., book(tmp, 23, 1, 1));
 48      _h_deltaPhi_4J_phi12.add( 800.,  1000., book(tmp, 24, 1, 1));
 49      _h_deltaPhi_4J_phi12.add( 1000., 7000., book(tmp, 25, 1, 1));
 50
 51      _h_deltaPhi_3J_phimin2J.add( 200.,  300., book(tmp, 26, 1, 1));
 52      _h_deltaPhi_3J_phimin2J.add( 300.,  400., book(tmp, 27, 1, 1));
 53      _h_deltaPhi_3J_phimin2J.add( 400.,  500., book(tmp, 28, 1, 1));
 54      _h_deltaPhi_3J_phimin2J.add( 500.,  600., book(tmp, 29, 1, 1));
 55      _h_deltaPhi_3J_phimin2J.add( 600.,  700., book(tmp, 30, 1, 1));
 56      _h_deltaPhi_3J_phimin2J.add( 700.,  800., book(tmp, 31, 1, 1));
 57      _h_deltaPhi_3J_phimin2J.add( 800.,  1000., book(tmp, 32, 1, 1));
 58      _h_deltaPhi_3J_phimin2J.add( 1000., 7000., book(tmp, 33, 1, 1));
 59
 60      _h_deltaPhi_4J_phimin2J.add( 200.,  300., book(tmp, 34, 1, 1));
 61      _h_deltaPhi_4J_phimin2J.add( 300.,  400., book(tmp, 35, 1, 1));
 62      _h_deltaPhi_4J_phimin2J.add( 400.,  500., book(tmp, 36, 1, 1));
 63      _h_deltaPhi_4J_phimin2J.add( 500.,  600., book(tmp, 37, 1, 1));
 64      _h_deltaPhi_4J_phimin2J.add( 600.,  700., book(tmp, 38, 1, 1));
 65      _h_deltaPhi_4J_phimin2J.add( 700.,  800., book(tmp, 39, 1, 1));
 66      _h_deltaPhi_4J_phimin2J.add( 800.,  1000., book(tmp, 40, 1, 1));
 67      _h_deltaPhi_4J_phimin2J.add( 1000., 7000., book(tmp, 41, 1, 1));
 68
 69    }
 70
 71
 72    void analyze(const Event & event) {
 73      const Jets& jets = apply<JetAlg>(event, "antikT").jetsByPt();
 74
 75      // 2 jet case and Delta_phi12
 76      if( jets.size() >= 2 ) {
 77        if ( (jets[0].pT() >= 200.*GeV)  &&  (jets[1].pT() >= 100.*GeV) ) {
 78          if ( (fabs(jets[0].rap()) <= 2.5)  &&  (fabs(jets[1].rap()) <= 2.5) ) {
 79            double dphi = deltaPhi(jets[0].phi(), jets[1].phi());
 80            _h_deltaPhi_2J_phi12.fill(jets[0].pT(), dphi, 1.0);
 81          }
 82        }
 83      }
 84
 85      // 3 jet case and Delta_phi12
 86      if ( jets.size() >= 3 ) {
 87        if ( (jets[0].pT() >= 200.*GeV)  &&  (jets[1].pT() >= 100.*GeV)  && (jets[2].pT() >= 100.*GeV) ) {
 88          if ( (fabs(jets[0].rap()) <= 2.5)  &&  (fabs(jets[1].rap()) <= 2.5) &&  (fabs(jets[2].rap()) <= 2.5)) {
 89            double dphi = deltaPhi(jets[0].phi(), jets[1].phi());
 90            _h_deltaPhi_3J_phi12.fill(jets[0].pT(), dphi, 1.0);
 91          }
 92        }
 93      }
 94
 95      // 4 jet case and Delta_phi12
 96      if ( jets.size() >= 4 ) {
 97        if ( (jets[0].pT() >= 200.*GeV)  &&  (jets[1].pT() >= 100.*GeV)  && (jets[2].pT() >= 100.*GeV)   && (jets[3].pT() >= 100.*GeV)) {
 98          if ( (fabs(jets[0].rap()) <= 2.5)  &&  (fabs(jets[1].rap()) <= 2.5) &&  (fabs(jets[2].rap()) <= 2.5) &&  (fabs(jets[3].rap()) <= 2.5)) {
 99            double dphi = deltaPhi(jets[0].phi(), jets[1].phi());
100            _h_deltaPhi_4J_phi12.fill(jets[0].pT(), dphi, 1.0);
101          }
102        }
103      }
104
105      // 3 jet case and Delta_Phi_min2j
106      if ( jets.size() >= 3 ) {
107        if ( (jets[0].pT() >= 200.*GeV)  &&  (jets[1].pT() >= 100.*GeV)  && (jets[2].pT() >= 100.*GeV) ) {
108          if ( (fabs(jets[0].rap()) <= 2.5)  &&  (fabs(jets[1].rap()) <= 2.5) &&  (fabs(jets[2].rap()) <= 2.5)) {
109            double dphi01 = deltaPhi(jets[0].phi(), jets[1].phi());
110            if (dphi01 >= PI/2. ){
111              double dphi02 = deltaPhi(jets[0].phi(), jets[2].phi());
112              double dphi12 = deltaPhi(jets[1].phi(), jets[2].phi());
113              // evaluate DPhi2Jmin
114              vector<double> Dphis2J{dphi01,dphi02,dphi12};
115              double DPhi2Jmin = min(Dphis2J);
116              // double Dphis2J[3] = {dphi01,dphi02,dphi12};
117              // double DPhi2Jmin = Dphis2J[0];
118              // for (int gg=1; gg<3; ++gg) { if (DPhi2Jmin>Dphis2J[gg]) DPhi2Jmin = Dphis2J[gg]; }
119              _h_deltaPhi_3J_phimin2J.fill(jets[0].pT(), DPhi2Jmin, 1.0);
120            }
121          }
122        }
123      }
124
125      // 4 jet case and Delta_Phi_min2j
126      if ( jets.size() >= 4 ) {
127        if ( (jets[0].pT() >= 200.*GeV)  &&  (jets[1].pT() >= 100.*GeV)  && (jets[2].pT() >= 100.*GeV)   && (jets[3].pT() >= 100.*GeV)) {
128          if ( (fabs(jets[0].rap()) <= 2.5)  &&  (fabs(jets[1].rap()) <= 2.5) &&  (fabs(jets[2].rap()) <= 2.5) &&  (fabs(jets[3].rap()) <= 2.5)) {
129            double dphi01 = deltaPhi(jets[0].phi(), jets[1].phi());
130            if (dphi01 >= PI/2.) {
131              double dphi02 = deltaPhi(jets[0].phi(), jets[2].phi());
132              double dphi03 = deltaPhi(jets[0].phi(), jets[3].phi());
133              double dphi12 = deltaPhi(jets[1].phi(), jets[2].phi());
134              double dphi13 = deltaPhi(jets[1].phi(), jets[3].phi());
135              double dphi23 = deltaPhi(jets[2].phi(), jets[3].phi());
136              /// evaluate DPhi2Jmin
137              // double Dphis2J[6]={dphi01,dphi02,dphi03,dphi12,dphi13,dphi23};
138              // double DPhi2Jmin=Dphis2J[0];
139              // for(int gg=1; gg<6; ++gg){ if(DPhi2Jmin>Dphis2J[gg]){DPhi2Jmin=Dphis2J[gg];} }
140              vector<double> Dphis2J{dphi01,dphi02,dphi03,dphi12,dphi13,dphi23};
141              double DPhi2Jmin = min(Dphis2J);
142              _h_deltaPhi_4J_phimin2J.fill(jets[0].pT(), DPhi2Jmin, 1.0);
143            }
144          }
145        }
146      }
147    }  // end analyze
148
149
150    void finalize() {
151      for (Histo1DPtr histo : _h_deltaPhi_2J_phi12.histos()) normalize(histo);
152      for (Histo1DPtr histo : _h_deltaPhi_3J_phi12.histos()) normalize(histo);
153      for (Histo1DPtr histo : _h_deltaPhi_4J_phi12.histos()) normalize(histo);
154      for (Histo1DPtr histo : _h_deltaPhi_3J_phimin2J.histos()) normalize(histo);
155      for (Histo1DPtr histo : _h_deltaPhi_4J_phimin2J.histos()) normalize(histo);
156    }
157
158
159  private:
160
161    BinnedHistogram _h_deltaPhi_2J_phi12;
162    BinnedHistogram _h_deltaPhi_3J_phi12;
163    BinnedHistogram _h_deltaPhi_4J_phi12;
164    BinnedHistogram _h_deltaPhi_3J_phimin2J;
165    BinnedHistogram _h_deltaPhi_4J_phimin2J;
166
167  };
168
169
170  RIVET_DECLARE_PLUGIN(CMS_2018_I1643640);
171
172}