Rivet Analyses Reference

CMS_2016_I1487288

$WZ$ production cross-section in $pp$ collisions at 8 TeV
Experiment: CMS (LHC)
Inspire ID: 1487288
Status: VALIDATED
Authors:
  • Shehu AbdusSalam
  • Andy Buckley
References:
  • Eur.Phys.J. C77 (2017) no.4, 236
  • DOI:10.1140/epjc/s10052-017-4730-z
  • CMS-SMP-14-014
  • CERN-EP-2016-205
  • arXiv: 1609.05721
Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • WZ events at 8 TeV

Fully leptonic $WZ$ diboson production at the CMS experiment, measured inclusively at 7 and 8 TeV, and differentially at 8 TeV only (this analysis). The differential cross-sections are measured with respect to $Z$ $p_T$, leading jet $p_T$, and the jet multiplicity. This data is useful for constraining MC models of QCD and EW activity in diboson events, and for constraining anomalous couplings. NOTE: total cross-section is fixed to the measured value in this analysis, since the correction factors from fiducial to inclusive $WZ$ production were not published.

Source code: CMS_2016_I1487288.cc
  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/ZFinder.hh"
#include "Rivet/Projections/WFinder.hh"

namespace Rivet {


  /// @brief WZ production cross section in pp collisions at 7 and 8 TeV
  class CMS_2016_I1487288 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2016_I1487288);


    /// @name Analysis methods
    //@{

    /// Book histograms and initialise projections before the run
    void init() {

      FinalState fs(Cuts::abseta < 4.9);

      FastJets fj(fs, FastJets::ANTIKT, 0.5, JetAlg::Muons::ALL, JetAlg::Invisibles::DECAY);
      declare(fj, "Jets");

      ZFinder zeeFinder(fs, Cuts::abseta < 2.5 && Cuts::pT > 20*GeV, PID::ELECTRON, 71*GeV, 111*GeV);
      declare(zeeFinder, "Zee");

      ZFinder zmumuFinder(fs, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV, PID::MUON, 71*GeV, 111*GeV);
      declare(zmumuFinder, "Zmumu");

      WFinder weFinder(fs, Cuts::abseta < 2.5 && Cuts::pT > 20*GeV, PID::ELECTRON, 60*GeV, 100*GeV, 30*GeV);
      declare(weFinder, "We");

      WFinder wmuFinder(fs, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV, PID::MUON, 60*GeV, 100*GeV, 30*GeV);
      declare(wmuFinder, "Wmu");

      book(_h_ZpT, "d03-x01-y01");
      book(_h_Njet, "d04-x01-y01", {-0.5, 0.5, 1.5, 2.5, 3.5}); ///< @todo Ref data has null bin widths
      book(_h_JpT, "d05-x01-y01");

      MSG_WARNING("\033[91;1mLIMITED VALIDITY - check info file for details!\033[m");
    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {

      // Find Z -> l+ l-
      const ZFinder& zeeFS = apply<ZFinder>(event, "Zee");
      const ZFinder& zmumuFS = apply<ZFinder>(event, "Zmumu");
      const Particles zlls = zeeFS.bosons() + zmumuFS.bosons();
      if (zlls.empty()) vetoEvent;

      // Next find the W
      const WFinder& weFS = apply<WFinder>(event, "We");
      const WFinder& wmuFS = apply<WFinder>(event, "Wmu");
      const Particles wls = weFS.bosons() + wmuFS.bosons();
      if (wls.empty()) vetoEvent;


      // If more than one Z candidate, use the one with Mll nearest to MZ
      const Particles zlls_mz = sortBy(zlls, [](const Particle& a, const Particle& b){
          return fabs(a.mass() - 91.2*GeV) < fabs(b.mass() - 91.2*GeV); });
      const Particle& Z = zlls_mz.front();
      // const bool isZee = any(Z.constituents(), hasAbsPID(PID::ELECTRON));

      // If more than one Z candidate, use the one with Mll nearest to MZ
      const Particles wls_mw = sortBy(wls, [](const Particle& a, const Particle& b){
          return fabs(a.mass() - 80.4*GeV) < fabs(b.mass() - 80.4*GeV); });
      const Particle& W = wls_mw.front();
      // const bool isWe = any(W.constituents(), hasAbsPID(PID::ELECTRON));

      // Isolate W and Z charged leptons from each other
      for (const Particle& lw : W.constituents()) {
        if (lw.charge3() == 0) continue;
        for (const Particle& lz : Z.constituents()) {
          if (deltaR(lw, lz) < 0.1) vetoEvent;
        }
      }

      // Fill Z pT histogram
      _h_ZpT->fill(Z.pT()/GeV);


      // Isolate jets from W and Z charged leptons
      const Particles wzleps = filter_select(W.constituents()+Z.constituents(), isChLepton);
      const Jets& jets = apply<FastJets>("Jets", event).jetsByPt(Cuts::pT > 30*GeV and Cuts::abseta < 2.5);
      const Jets isojets = discardIfAnyDeltaRLess(jets, wzleps, 0.5);

      // Fill jet histograms
      _h_Njet->fill(isojets.size());
      if (!isojets.empty()) _h_JpT->fill(isojets[0].pT()/GeV);

    }


    /// Normalise histograms etc., after the run
    void finalize() {
      // Total cross-section is corrected for BR(W->lnu), BR(Z->ll), leptonic-tau fraction f_tau = 6.5-7.6%,
      // and unpublished detector/acceptance signal efficiencies epsilon_sig. Fix to published value: valid for shape comparison only
      const double xsec8tev = 24.09; // picobarn;
      normalize(_h_ZpT,  xsec8tev);
      normalize(_h_Njet, xsec8tev);
      normalize(_h_JpT,  xsec8tev);
    }


  private:

    /// Histogram
    Histo1DPtr _h_ZpT, _h_Njet, _h_JpT;

  };



  RIVET_DECLARE_PLUGIN(CMS_2016_I1487288);

}