Rivet Analyses Reference

CMS_2020_I1814328

W$^+$W$^-$ boson pair production in proton-proton collisions at $\sqrt{s} =$ 13 TeV
Experiment: CMS (LHC)
Inspire ID: 1814328
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Guillelmo Gomez-Ceballos
References:
  • arXiv: 2009.00119
  • CMS-SMP-18-004
  • Phys. Rev. D 102 (2020) 092001
Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • pp to W+W- interactions at $\sqrt{s} = 13$ TeV. Data collected by CMS during the year 2016.

A measurement of the W W boson pair production cross section in proton-proton collisions at sqrt(s) = 13 TeV is presented. The data used in this study were collected with the CMS detector at the LHC and correspond to an integrated luminosity of 35.9 fb-1. The W+W- candidate events are selected by first requiring two oppositely-charged leptons (electrons or muons). Two methods for reducing background contributions are employed. In the first one, a sequence of requirements on kinematic quantities is applied allowing a measurement of the total cross section: 117.6 +/- 6.8 pb which agrees well with the theoretical cross section. Fiducial cross sections are also reported for events with zero jets or one jet, and the change in the zero-jet fiducial cross section with the jet pT threshold is measured. Normalized differential cross sections are reported within the fiducial region; these are compared to theoretical predictions based on quantum chromodynamics at next-to-next-to-leading-order accuracy. A second method for suppressing background contributions employs two random forest classifiers. The analysis based on this method includes a measurement of the total cross section and also a measurement of the normalized jet multiplicity distribution in W+W- events. Finally, a dilepton invariant mass distribution is used to probe for physics beyond the standard model in the context of an effective field theory and limits on the presence of dimension-6 operators are derived.

Source code: CMS_2020_I1814328.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
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/DressedLeptons.hh"
#include "Rivet/Projections/MissingMomentum.hh"
#include "Rivet/Projections/PromptFinalState.hh"

namespace Rivet {


  /// @brief Measurements of W+W- boson pair production in pp collisions at 13 TeV
  class CMS_2020_I1814328 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2020_I1814328);


    /// @name Analysis methods
    /// @{

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

      // Initialise and register projections

      // The basic final-state projection:
      // all final-state particles within
      // the given eta acceptance
      const FinalState fs(Cuts::abseta < 4.9);
      const FinalState fsjet4p5(Cuts::abseta < 4.5);
      const FinalState fsjet2p5(Cuts::abseta < 2.5);

      // The final-state particles declared above are clustered using FastJet with
      // the anti-kT algorithm and a jet-radius parameter 0.4
      // muons and neutrinos are excluded from the clustering
      FastJets jet4p5fs(fsjet4p5, FastJets::ANTIKT, 0.4);
      declare(jet4p5fs, "jets4p5");
      FastJets jet2p5fs(fsjet2p5, FastJets::ANTIKT, 0.4);
      declare(jet2p5fs, "jets2p5");

      // FinalState of prompt photons and bare muons and electrons in the event
      PromptFinalState photons(Cuts::abspid == PID::PHOTON);
      PromptFinalState bare_leps(Cuts::abspid == PID::MUON || Cuts::abspid == PID::ELECTRON);
      bare_leps.acceptTauDecays(false);

      // Dress the prompt bare leptons with prompt photons within dR < 0.1,
      // and apply some fiducial cuts on the dressed leptons
      Cut lepton_cuts = Cuts::abseta < 2.5 && Cuts::pT > 25*GeV;
      DressedLeptons dressed_leps(photons, bare_leps, 0.1, lepton_cuts);
      declare(dressed_leps, "leptons");

      // Missing momentum
      declare(MissingMomentum(fs), "MET");

      // Book histograms
      book(_h_WW_njets_norm , 2, 1, 1);
      book(_h_WW_mll_norm   , 4, 1, 1);
      book(_h_WW_ptlmax_norm, 5, 1, 1);
      book(_h_WW_ptlmin_norm, 6, 1, 1);
      book(_h_WW_dphill_norm, 7, 1, 1);
      book(_h_WW_njet0      , 8, 1, 1);

    }

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

      // Apply a missing-momentum cut
      if (apply<MissingMomentum>(event, "MET").missingPt() < 20*GeV) return;

      // Retrieve dressed leptons, sorted by pT
      vector<DressedLepton> leptons = apply<DressedLeptons>(event, "leptons").dressedLeptons();

      // Retrieve clustered jets, sorted by pT, with a minimum pT cut
      Jets jets25 = apply<FastJets>(event, "jets4p5").jetsByPt(Cuts::pT > 25*GeV);
      Jets jets30 = apply<FastJets>(event, "jets4p5").jetsByPt(Cuts::pT > 30*GeV);
      Jets jets35 = apply<FastJets>(event, "jets4p5").jetsByPt(Cuts::pT > 35*GeV);
      Jets jets45 = apply<FastJets>(event, "jets4p5").jetsByPt(Cuts::pT > 45*GeV);
      Jets jets60 = apply<FastJets>(event, "jets4p5").jetsByPt(Cuts::pT > 60*GeV);
      Jets jetsNj = apply<FastJets>(event, "jets2p5").jetsByPt(Cuts::pT > 30*GeV);

      // Remove all jets within dR < 0.4 of a dressed lepton
      idiscardIfAnyDeltaRLess(jets25, leptons, 0.4);
      idiscardIfAnyDeltaRLess(jets30, leptons, 0.4);
      idiscardIfAnyDeltaRLess(jets35, leptons, 0.4);
      idiscardIfAnyDeltaRLess(jets45, leptons, 0.4);
      idiscardIfAnyDeltaRLess(jets60, leptons, 0.4);
      idiscardIfAnyDeltaRLess(jetsNj, leptons, 0.4);

      if (leptons.size() == 2 && leptons[0].pid() * leptons[1].pid() < 0) {
        FourMomentum dilCand = leptons[0].momentum() + leptons[1].momentum();
        if (dilCand.mass() > 20*GeV && dilCand.pT() > 30*GeV){
          double ptlmax = leptons[0].pT(); double ptlmin = leptons[1].pT();
          if (ptlmax < ptlmin) {
            ptlmax = leptons[1].pT(); ptlmin = leptons[0].pT();
          }

          if (leptons[0].abspid() != leptons[1].abspid()) {
            _h_WW_njets_norm ->fill(min((double)jetsNj.size()+1, 2.999));
            _h_WW_mll_norm   ->fill(min(dilCand.mass()/GeV, 1499.999));
            _h_WW_ptlmax_norm->fill(min(ptlmax/GeV, 399.999));
            _h_WW_ptlmin_norm->fill(min(ptlmin/GeV, 149.999));
            _h_WW_dphill_norm->fill(deltaPhi(leptons[0], leptons[1]));
          }

          if (jets25.size() == 0) _h_WW_njet0->fill(1.0);
          if (jets30.size() == 0) _h_WW_njet0->fill(2.0);
          if (jets35.size() == 0) _h_WW_njet0->fill(3.0);
          if (jets45.size() == 0) _h_WW_njet0->fill(4.0);
          if (jets60.size() == 0) _h_WW_njet0->fill(5.0);

        }
      }

    }


    /// @todo Replace with barchart()
    void normalizeToSum(Histo1DPtr hist) {
      double sum = 0.;
      for (size_t i = 0; i < hist->numBins(); ++i) {
        sum += hist->bin(i).height();
      }
      scale(hist, 1./sum);
    }


    /// Normalise histograms etc., after the run
    void finalize() {

      double norm = (sumOfWeights() != 0) ? crossSection()/picobarn/sumOfWeights() : 1.0;

      normalizeToSum(_h_WW_njets_norm );
      normalizeToSum(_h_WW_mll_norm   );
      normalizeToSum(_h_WW_ptlmax_norm);
      normalizeToSum(_h_WW_ptlmin_norm);
      normalizeToSum(_h_WW_dphill_norm);

      scale(_h_WW_njet0, norm);

    }

    /// @}


    /// @name Histograms
    /// @{
    Histo1DPtr _h_WW_njets_norm;
    Histo1DPtr _h_WW_mll_norm, _h_WW_ptlmax_norm, _h_WW_ptlmin_norm, _h_WW_dphill_norm;
    Histo1DPtr _h_WW_njet0;
    /// @}

  };



  RIVET_DECLARE_PLUGIN(CMS_2020_I1814328);

}