rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2021_I1866118

Study of Z boson plus jets events using variables sensitive to double-parton scattering in pp collisions at 13 TeV
Experiment: CMS (LHC)
Inspire ID: 1866118
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Rajat Gupta
  • Meena
  • Sunil Bansal
  • J.B. Singh
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • Drell-Yan events with $Z/\gamma^* \to \mu\mu$.

Double parton scattering is investigated using events with a Z boson and jets. The measurements are performed with proton-proton collision data recorded by the CMS experiment at the LHC at $\sqrt{s} = 13$ TeV, corresponding to an integrated luminosity of $35.8 \mathrm{fb}^{-1}$ collected in the year 2016. Differential cross sections of Z$+ \geq$ 1 jet and Z$+ \geq$ 2 jets are measured with transverse momentum of the jets above 20 GeV and pseudorapidity $|\eta| <$ 2.4. Several distributions with sensitivity to double parton scattering effects are measured as functions of the angle and the transverse momentum imbalance between the Z boson and the jets. The measured distributions are compared with predictions from several event generators with different hadronization models and different parameter settings for multiparton interactions. The measured distributions show a dependence on the hadronization and multiparton interaction simulation parameters, and are important input for future improvements of the simulations.

Source code: CMS_2021_I1866118.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/DileptonFinder.hh"
  5#include "Rivet/Projections/FinalState.hh"
  6#include "Rivet/Projections/ChargedFinalState.hh"
  7#include "Rivet/Projections/VetoedFinalState.hh"
  8#include "Rivet/Projections/FastJets.hh"
  9
 10namespace Rivet {
 11
 12
 13  /// Z boson + jets sensitive to double-parton scattering at 13 TeV
 14  class CMS_2021_I1866118 : public Analysis {
 15  public:
 16
 17    /// Constructor
 18    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2021_I1866118);
 19
 20
 21    /// Initialization
 22    void init() {
 23
 24      Cut cut = Cuts::abseta < 2.4 && Cuts::pT > 27 * GeV;
 25      DileptonFinder zmumufinder(91.2*GeV, 0.1, cut && Cuts::abspid == PID::MUON, Cuts::massIn(70*GeV, 110*GeV));
 26      declare(zmumufinder, "zmumufinder");
 27
 28      // Define veto FS in order to prevent Z-decay products entering the jet algorithm
 29
 30      VetoedFinalState had_fs;
 31      had_fs.addVetoOnThisFinalState(zmumufinder);
 32      FastJets jets(had_fs, JetAlg::ANTIKT, 0.4);
 33      jets.useInvisibles();
 34      declare(jets, "jets");
 35
 36      book(h_dphi_Z1J_cn, 1, 1, 1);
 37      book(h_reldpt_Z1J_cn, 2, 1, 1);
 38      book(h_dphi_Zdijet_Z2J_cn, 3, 1, 1);
 39      book(h_reldpt_Zdijet_Z2J_cn, 4, 1, 1);
 40      book(h_reldpt_j1j2_Z2J_cn, 5, 1, 1);
 41
 42      book(_h["dphi_Z1J"], 6, 1, 1);
 43      book(_h["reldpt_Z1J"], 7, 1, 1);
 44      book(_h["dphi_Zdijet_Z2J"], 8, 1, 1);
 45      book(_h["reldpt_Zdijet_Z2J"], 9, 1, 1);
 46      book(_h["reldpt_j1j2_Z2J"], 10, 1, 1);
 47    }
 48
 49
 50    /// Perform the per-event analysis
 51    void analyze(const Event& event) {
 52      const DileptonFinder& zmumufinder = apply<DileptonFinder>(event, "zmumufinder");
 53      const Particles& zmumus = zmumufinder.bosons();
 54      if (zmumus.size() != 1) {
 55        vetoEvent;
 56      }
 57
 58      // Find the (dressed!) leptons
 59      const Particles& leptons = zmumufinder.constituents();
 60      if (leptons.size() != 2)
 61        vetoEvent;
 62
 63      Jets jets = apply<FastJets>(event, "jets").jetsByPt(Cuts::pT > 20 * GeV && Cuts::abseta < 2.4);
 64      idiscardIfAnyDeltaRLess(jets, leptons, 0.4);
 65
 66      const size_t Njets = jets.size();
 67
 68      if (Njets < 1)
 69        vetoEvent;
 70
 71      if (Njets >= 1) {
 72        h_dphi_Z1J_cn->fill(deltaPhi(zmumus[0], jets[0]));
 73        h_reldpt_Z1J_cn->fill((zmumus[0] + jets[0].momentum()).pt() / (zmumus[0].pt() + jets[0].pt()));
 74
 75        _h["dphi_Z1J"]->fill(deltaPhi(zmumus[0], jets[0]));
 76        _h["reldpt_Z1J"]->fill((zmumus[0] + jets[0].momentum()).pt() / (zmumus[0].pt() + jets[0].pt()));
 77      }
 78
 79      if (Njets >= 2) {
 80        FourMomentum dij = jets[0].momentum() + jets[1].momentum();
 81
 82        h_dphi_Zdijet_Z2J_cn->fill(deltaPhi(zmumus[0], dij));
 83        h_reldpt_Zdijet_Z2J_cn->fill((zmumus[0] + dij).pt() / (zmumus[0].pt() + dij.pt()));
 84        h_reldpt_j1j2_Z2J_cn->fill(dij.pt() / (jets[0].pt() + jets[1].pt()));
 85
 86        _h["dphi_Zdijet_Z2J"]->fill(deltaPhi(zmumus[0], dij));
 87        _h["reldpt_Zdijet_Z2J"]->fill((zmumus[0] + dij).pt() / (zmumus[0].pt() + dij.pt()));
 88        _h["reldpt_j1j2_Z2J"]->fill(dij.pt() / (jets[0].pt() + jets[1].pt()));
 89      }
 90    }
 91
 92
 93    /// Normalise histograms etc., after the run
 94    void finalize() {
 95      double norm = (sumOfWeights() != 0) ? crossSection()/picobarn/sumOfWeights() : 1.0;
 96
 97      scale(h_dphi_Z1J_cn, norm);
 98      scale(h_reldpt_Z1J_cn, norm);
 99      scale(h_dphi_Zdijet_Z2J_cn, norm);
100      scale(h_reldpt_Zdijet_Z2J_cn, norm);
101      scale(h_reldpt_j1j2_Z2J_cn, norm);
102
103      for (auto& item : _h) {
104        double rho = item.second->densitySum(false);
105        if (rho)  scale(item.second, 1.0/rho);
106      }
107    }
108
109
110  private:
111
112    /// @name Histogram objects
113    /// @{
114    Histo1DPtr h_dphi_Z1J_cn;
115    Histo1DPtr h_reldpt_Z1J_cn;
116    Histo1DPtr h_dphi_Zdijet_Z2J_cn;
117    Histo1DPtr h_reldpt_Zdijet_Z2J_cn;
118    Histo1DPtr h_reldpt_j1j2_Z2J_cn;
119    map<string,Histo1DPtr> _h;
120    /// @}
121
122  };
123
124
125  RIVET_DECLARE_PLUGIN(CMS_2021_I1866118);
126
127}