Rivet Analyses Reference

CMS_2020_I1794169

Measurements of production cross sections of WZ and same-sign WW boson pairs in association with two jets in proton-proton collisions at $\sqrt{s} =$ 13 TeV
Experiment: CMS (LHC)
Inspire ID: 1794169
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Guillelmo Gomez-Ceballos
  • Kenneth Long
References:
  • arXiv: 2005.01173
  • CMS-SMP-19-012
  • Phys. Lett. B 809 (2020) 135710
Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • VBS WZ and W+W+ interactions at $\sqrt{s} = 13$ TeV. Data collected by CMS during the year 2016-18.

Measurements of production cross sections of WZ and same-sign W+W+ boson pairs in association with two jets in proton-proton collisions at sqrt{s}=13TeV at the LHC are reported. The data sample corresponds to an integrated luminosity of 137fb-1, collected with the CMS detector during 2016--2018. The measurements are performed in the leptonic decay modes. Differential fiducial cross sections as functions of the invariant masses of the jet and charged lepton pairs, as well as of the leading-lepton transverse momentum, are measured for WW production and are consistent with the standard model predictions. The dependence of differential cross sections on the invariant mass of the jet pair is also measured for WZ production. An observation of electroweak production of WZ boson pairs is reported with an observed (expected) significance of 6.8 (5.3) standard deviations. Constraints are obtained on the structure of quartic vector boson interactions in the framework of effective field theory.

Source code: CMS_2020_I1794169.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
162
163
164
165
166
167
168
169
170
171
172
173
174
// -*- 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 Production cross-sections of WZ and same-sign WW with two jets in pp collisions at 13 TeV
  class CMS_2020_I1794169 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2020_I1794169);


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

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

      // Initialise and register projections
      _mode = 0;
      if ( getOption("LMODE") == "WZ" ) _mode = 1;

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

      // 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 jet4p7fs(fsjet4p7, FastJets::ANTIKT, 0.4);
      declare(jet4p7fs, "jets4p7");

      // 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 > 20*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_mjj   ,  9, 1, 1);
      book(_h_WW_mll   , 11, 1, 1);
      book(_h_WW_ptlmax, 13, 1, 1);
      book(_h_WZ_mjj   , 15, 1, 1);

    }

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

      // Retrieve dressed leptons, sorted by pT
      Particles leptons = apply<DressedLeptons>(event, "leptons").particles();

      // Apply a #leptons requirement
      if (leptons.size() <= 1 || leptons.size() >= 4) return;

      // Retrieve clustered jets, sorted by pT, with a minimum pT cut
      Jets jets50 = apply<FastJets>(event, "jets4p7").jetsByPt(Cuts::pT > 50*GeV);

      // Remove all jets within dR < 0.4 of a dressed lepton
      idiscardIfAnyDeltaRLess(jets50, leptons, 0.4);

      // Apply a njets >= 2 cut
      if (jets50.size() < 2) return;

      FourMomentum dijetCand = jets50[0].momentum() + jets50[1].momentum();
      double deltaEtaJJ = std::abs(jets50[0].eta() - jets50[1].eta());

      // Apply a mjj > 500 and detajj > 2.5 cuts
      if (dijetCand.mass() <= 500*GeV || deltaEtaJJ <= 2.5) return;

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

          _h_WW_mjj   ->fill(min(dijetCand.mass()/GeV, 2999.999));
          _h_WW_mll   ->fill(min(dilCand.mass()/GeV, 499.999));
          _h_WW_ptlmax->fill(min(ptlmax/GeV, 299.999));

        }
      }

      // WZ selection
      else if (leptons.size() == 3 && _mode == 1) {
        double mllZ = 10000; int iW = -1;
        if (leptons[0].pid() * leptons[1].pid() < 0 && leptons[0].abspid() == leptons[1].abspid() &&
            fabs((leptons[0].momentum() + leptons[1].momentum()).mass() - 91.1876*GeV) < fabs(mllZ - 91.1876*GeV)) {
          mllZ = (leptons[0].momentum() + leptons[1].momentum()).mass(); iW = 2;
        }

        if (leptons[0].pid() * leptons[2].pid() < 0 && leptons[0].abspid() == leptons[2].abspid() &&
            fabs((leptons[0].momentum() + leptons[2].momentum()).mass() - 91.1876*GeV) < fabs(mllZ - 91.1876*GeV)) {
          mllZ = (leptons[0].momentum() + leptons[2].momentum()).mass(); iW = 1;
        }

        if (leptons[1].pid() * leptons[2].pid() < 0 && leptons[1].abspid() == leptons[2].abspid() &&
            fabs((leptons[1].momentum() + leptons[2].momentum()).mass() - 91.1876*GeV) < fabs(mllZ - 91.1876*GeV)) {
          mllZ = (leptons[1].momentum() + leptons[2].momentum()).mass(); iW = 0;
        }

        // Plot
        if (iW >= 0 && fabs(mllZ - 91.1876*GeV) < 15*GeV) {
          _h_WZ_mjj->fill(min(dijetCand.mass()/GeV, 2999.999));
        }
      }

    }


    /// @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();
        float width = hist->bin(i).width();
        hist->bin(i).scaleW(width != 0 ? width : 1.);
      }
      if(hist->integral() > 0) scale(hist, 1./hist->integral());
    }


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

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

      scale(_h_WW_mjj   , norm);
      scale(_h_WW_mll   , norm);
      scale(_h_WW_ptlmax, norm);
      scale(_h_WZ_mjj   , norm);

    }

    //@}


  private:

    /// Lepton-mode flag
    size_t _mode;

    /// @name Histograms
    /// @{
    Histo1DPtr _h_WW_mjj, _h_WW_mll, _h_WW_ptlmax, _h_WZ_mjj;
    /// @}

  };



  RIVET_DECLARE_PLUGIN(CMS_2020_I1794169);

}