rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2021_I1847230

Measurements of angular distance and momentum ratio distributions in three-jet and Z + two-jet final states in pp collisions at 8 and 13 TeV
Experiment: CMS (LHC)
Inspire ID: 1847230
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Hannes Jung
  • Hyunyong Kim
  • Olga Kodolova
References: Beams: p+ p+
Beam energies: (4000.0, 4000.0); (6500.0, 6500.0) GeV
Run details:
  • QCD with pt>300 or Z -> mu+ mu-

Collinear (small-angle) and large-angle, as well as soft and hard radiations are investigated in three-jet and Z + two-jet events collected in proton-proton collisions at the LHC. The normalized production cross sections are measured as a function of the ratio of transverse momenta of two jets and their angular separation. The measurements in the three-jet and Z + two-jet events are based on data collected at a center-of-mass energy of 8 TeV, corresponding to an integrated luminosity of 19.8 fb^{-1}. The Z + two-jet events are reconstructed in the dimuon decay channel of the Z boson. The three-jet measurement is extended to include \sqrt{s} = 13 TeV data corresponding to an integrated luminosity of 2.3 fb^{-1}. The results are compared to predictions from event generators that include parton showers, multiple parton interactions, and hadronization. The collinear and soft regions are in general well described by parton showers, whereas the regions of large angular separation are often best described by calculations using higher-order matrix elements.

Source code: CMS_2021_I1847230.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/FastJets.hh"
  5#include "Rivet/Projections/DileptonFinder.hh"
  6
  7namespace Rivet {
  8
  9
 10  /// @brief Angular distance and momentum ratios in 3-jet and Z+2-jet final states
 11  class CMS_2021_I1847230 : public Analysis {
 12  public:
 13
 14    CMS_2021_I1847230 ()
 15      : Analysis("CMS_2021_I1847230")
 16    {}
 17
 18
 19    void init() {
 20      _mode = 0;
 21      if ( getOption("MODE") == "QCD8TeV" ) _mode = 1;
 22      else if ( getOption("MODE") == "QCD13TeV" ) _mode = 2;
 23      else if ( getOption("MODE") == "ZJet" ) _mode = 3;
 24
 25      if (_mode == 1) {
 26         _jr = 0.5;
 27         book(_h1, "d01-x01-y01");
 28         book(_h2, "d02-x01-y01");
 29         book(_h3, "d03-x01-y01");
 30         book(_h4, "d04-x01-y01");
 31       }
 32       if (_mode == 2) {
 33         _jr = 0.4;
 34         book(_h1, "d05-x01-y01");
 35         book(_h2, "d06-x01-y01");
 36         book(_h3, "d07-x01-y01");
 37         book(_h4, "d08-x01-y01");
 38       }
 39       if (_mode == 1 or _mode == 2) {
 40         const FastJets jets(FinalState(), JetAlg::ANTIKT, _jr);
 41         declare(jets, "jets");
 42       }
 43       if (_mode == 3) {
 44         //PromptFinalState pfs(Cuts::abseta < 2.4 and Cuts::pT > 100*MeV);
 45         DileptonFinder zfinder(91.2*GeV, 0.2, Cuts::abseta < 5 and Cuts::pT > 30*GeV and
 46                                Cuts::abspid == PID::MUON, Cuts::massIn(70*GeV, 110*GeV));
 47         declare(zfinder, "DileptonFinder");
 48         declare(FastJets(zfinder.remainingFinalState(), JetAlg::ANTIKT, 0.5), "JetsAK5_zj");
 49
 50         book(_h1, "d09-x01-y01");
 51         book(_h2, "d10-x01-y01");
 52         book(_h3, "d11-x01-y01");
 53         book(_h4, "d12-x01-y01");
 54
 55         book(_ZJw_gen, "TMP/ZJw_gen");
 56       }
 57    }
 58
 59    void analyze(const Event& event) {
 60
 61      if (_mode == 1 or _mode ==2) {
 62        const Jets& jets = apply<JetFinder>(event, "jets").jetsByPt(Cuts::pT > 30.0*GeV);
 63        if (jets.size() < 3) vetoEvent;
 64        const FourMomentum jet1 = jets[0].momentum();
 65        const FourMomentum jet2 = jets[1].momentum();
 66        const FourMomentum jet3 = jets[2].momentum();
 67
 68        if (jet1.pT() < 510.0*GeV) vetoEvent;
 69        if (jet1.absrapidity() > 2.5 or jet2.absrapidity() > 2.5) vetoEvent;
 70        const double del_phi12 = mapAngle0ToPi(jet2.phi() - jet1.phi());
 71        if (abs(del_phi12 - M_PI) > 1.0) vetoEvent;
 72        const double jet3_pt_jet2_pt = jet3.pT()/jet2.pT();
 73        if (!inRange(jet3_pt_jet2_pt, 0.1, 0.9)) vetoEvent;
 74        const double del_r23 = deltaR(jet3.rapidity(), jet3.phi(), jet2.rapidity(), jet2.phi());
 75        if (!inRange(del_r23, _jr+0.1, 1.5)) vetoEvent;
 76
 77        if (del_r23 < 1.0) _h1->fill(jet3_pt_jet2_pt);
 78        if (del_r23 > 1.0) _h2->fill(jet3_pt_jet2_pt);
 79        if (jet3_pt_jet2_pt < 0.3) _h3->fill(del_r23);
 80        if (jet3_pt_jet2_pt > 0.6) _h4->fill(del_r23);
 81      }
 82
 83      if (_mode == 3) {
 84        const DileptonFinder& zfinder = apply<DileptonFinder>(event, "DileptonFinder");
 85        if (zfinder.bosons().size() != 1) vetoEvent;
 86        const Particle& z = zfinder.bosons()[0];
 87        const Particles leptons = sortBy(zfinder.constituents(), cmpMomByPt);
 88        if (leptons[0].pT() < 25.0*GeV or leptons[1].pT() < 10.0*GeV or z.pT() < 80.0*GeV) vetoEvent;
 89        if (leptons[0].absrapidity() > 2.1 or leptons[1].absrapidity() > 2.4) vetoEvent;
 90
 91        const PseudoJets& psjetsAK5_zj = apply<FastJets>(event, "JetsAK5_zj").pseudojetsByPt(20.0*GeV);
 92
 93        if (psjetsAK5_zj.empty()) vetoEvent;
 94
 95        const PseudoJet& j0 = psjetsAK5_zj[0];
 96        const FourMomentum jmom0(j0.e(), j0.px(), j0.py(), j0.pz());
 97
 98        if (jmom0.absrapidity() > 1.0 or jmom0.pT() < 80.0*GeV) vetoEvent;
 99        if (!(deltaPhi(z, jmom0) > 2.0 and deltaR(leptons[0], jmom0) > 0.5 and deltaR(leptons[1], jmom0) > 0.5)) vetoEvent;
100
101        _ZJw_gen ->fill();
102
103        if(psjetsAK5_zj.size() < 2) vetoEvent;
104
105        const PseudoJet& j1 = psjetsAK5_zj[1];
106        const FourMomentum jmom1(j1.e(), j1.px(), j1.py(), j1.pz());
107        if(deltaR(leptons[0], jmom1) < 0.5 or deltaR(leptons[1], jmom1) < 0.5 or jmom1.absrapidity() > 2.4) vetoEvent;
108
109        const double dR_gen_Jj = deltaR(jmom0, jmom1);
110        if (!inRange(dR_gen_Jj, 0.5, 1.5)) vetoEvent;
111        const double rPt_gen_Jj = jmom1.pT()/jmom0.pT();
112
113        if (dR_gen_Jj < 1.0) _h1->fill(rPt_gen_Jj);
114        if (dR_gen_Jj > 1.0) _h2->fill(rPt_gen_Jj);
115        if (rPt_gen_Jj < 0.3) _h3->fill(dR_gen_Jj);
116        if (rPt_gen_Jj > 0.6) _h4->fill(dR_gen_Jj);
117      }
118    }
119
120
121    void finalize() {
122      if (_mode == 1 or _mode == 2) {
123        normalize(_h1);
124        normalize(_h2);
125        normalize(_h3);
126        normalize(_h4);
127        //scale(_h1, 1.0/_h1->effNumEntries());
128        //scale(_h2, 1.0/_h2->effNumEntries());
129        //scale(_h3, 1.0/_h3->effNumEntries());
130        //scale(_h4, 1.0/_h4->effNumEntries());
131      }
132      if (_mode == 3) {
133        for (size_t i = 0; i < _h1->numBins(); i++) {
134          _h1->bin(i).scaleW(_h1->bin(i).xWidth());
135        }
136        for (size_t i = 0; i < _h2->numBins(); i++) {
137          _h2->bin(i).scaleW(_h2->bin(i).xWidth());
138        }
139        for (size_t i = 0; i < _h3->numBins(); i++) {
140          _h3->bin(i).scaleW(_h3->bin(i).xWidth());
141        }
142        for (size_t i = 0; i < _h4->numBins(); i++) {
143          _h4->bin(i).scaleW(_h4->bin(i).xWidth());
144        }
145        scale(_h1, 1.0/ *_ZJw_gen);
146        scale(_h2, 1.0/ *_ZJw_gen);
147        scale(_h3, 1.0/ *_ZJw_gen);
148        scale(_h4, 1.0/ *_ZJw_gen);
149      }
150    }
151
152
153  private:
154
155    Histo1DPtr  _h1;
156    Histo1DPtr  _h2;
157    Histo1DPtr  _h3;
158    Histo1DPtr  _h4;
159
160    CounterPtr _ZJw_gen;
161
162    double _jr;
163
164    size_t _mode;
165
166  };
167
168
169  RIVET_DECLARE_PLUGIN(CMS_2021_I1847230);
170
171}