rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2014_I1298023

Same-sign WW cross section
Experiment: ATLAS (LHC)
Inspire ID: 1298023
Status: VALIDATED
Authors:
  • Simone Pagan Griso
  • Christian Gutschow
References: Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • p + p -> W W j j, where the W bosons decay leptonically

First measurement of $W^\pm W^\pm jj$, same-electric-charge diboson production in association with two jets, using $20.3 \mathrm{fb}^{-1}$ of proton-proton collision data at $\sqrt{s} = 8$\,TeV recorded by the ATLAS detector at the Large Hadron Collider. Events with two reconstructed same-charge leptons ($e^\pm e^\pm$, $e^\pm \mu^\pm$, and $\mu^\pm \mu^\pm$) and two or more jets are analyzed. Production cross sections are measured in two fiducial regions, with different sensitivities to the electroweak and strong production mechanisms. First evidence for $W^\pm W^\pm j j$ production and electroweak-only $W^\pm W^\pm jj$ production is observed with a significance of 4.5 and 3.6 standard deviations, respectively.

Source code: ATLAS_2014_I1298023.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/ChargedLeptons.hh"
  4#include "Rivet/Projections/VetoedFinalState.hh"
  5#include "Rivet/Projections/MissingMomentum.hh"
  6#include "Rivet/Projections/FastJets.hh"
  7#include "Rivet/Projections/LeptonFinder.hh"
  8
  9namespace Rivet {
 10
 11  /// ATLAS same-sign WW production at 8 TeV (PRL)
 12  class ATLAS_2014_I1298023 : public Analysis {
 13  public:
 14
 15    /// Constructor
 16    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2014_I1298023);
 17
 18
 19    /// @name Analysis methods
 20    /// @{
 21
 22    void init() {
 23
 24      const FinalState fs;
 25
 26      // bare leptons
 27      ChargedLeptons bare_leptons(fs);
 28
 29      // dressed leptons
 30      Cut cuts = (Cuts::abseta < 2.5) & (Cuts::pT > 25*GeV);
 31      LeptonFinder leptons(bare_leptons, fs, 0.1, cuts);
 32      declare(leptons, "leptons");
 33
 34      // MET
 35      declare(MissingMomentum(fs), "MissingET");
 36
 37      // jets
 38      VetoedFinalState vfs(fs);
 39      vfs.addVetoPairId(PID::MUON);
 40      vfs.vetoNeutrinos();
 41      declare(FastJets(vfs, JetAlg::ANTIKT, 0.4), "jets");
 42
 43      // book histogram
 44      book(_hist ,1, 1, 1);
 45    }
 46
 47
 48    /// Perform the per-event analysis
 49    void analyze(const Event& event) {
 50      const DressedLeptons& leptons = apply<LeptonFinder>(event, "leptons").dressedLeptons();
 51      if ( leptons.size() < 2 )  vetoEvent;
 52
 53      double minDR_ll = DBL_MAX, mll = -1.0;
 54      for (unsigned int i = 0; i < leptons.size(); ++i) {
 55        DressedLepton lep1 = leptons[i];
 56        for (unsigned int j = i + 1; j < leptons.size(); ++j) {
 57          DressedLepton lep2 = leptons[j];
 58          double dr = deltaR(lep1, lep2);
 59          if ( dr < minDR_ll )  minDR_ll = dr;
 60          double m = (lep1.momentum() + lep2.momentum()).mass();
 61          if ( mll < 0. || m < mll )  mll = m;
 62        }
 63      }
 64      if ( minDR_ll <= 0.3 || mll <= 20*GeV )  vetoEvent;
 65
 66      if ( leptons[0].charge() * leptons[1].charge() < 0.0 )  vetoEvent;
 67
 68      const MissingMomentum& met = apply<MissingMomentum>(event, "MissingET");
 69      if ( met.vectorEt().mod() <= 40*GeV )  vetoEvent;
 70
 71      const Jets& all_jets = apply<FastJets>(event, "jets").jetsByPt( (Cuts::abseta < 4.5) && (Cuts::pT > 30*GeV) );
 72      Jets jets;
 73      double minDR_overall = DBL_MAX;
 74      for (const Jet& jet : all_jets) {
 75        double minDR_jet = DBL_MAX, minDR_electrons = DBL_MAX;
 76        for( DressedLepton lep : leptons ) {
 77          double dr = deltaR(jet, lep);
 78          if ( dr < minDR_jet )  minDR_jet = dr;
 79          if ( lep.abspid() == 11 && dr < minDR_electrons )  minDR_electrons = dr;
 80        }
 81        if ( minDR_electrons < 0.05 )  continue; // veto jet if it overlaps with electron
 82        if ( minDR_jet < minDR_overall )  minDR_overall = minDR_jet;
 83        jets += jet;
 84      }
 85      if ( jets.size() < 2 || minDR_overall <= 0.3 )  vetoEvent;
 86
 87      FourMomentum dijet = jets[0].momentum() + jets[1].momentum();
 88      if ( dijet.mass() <= 500*GeV )  vetoEvent;
 89
 90      // inclusive region
 91      _hist->fill(1);
 92
 93      // VBS region
 94      if ( deltaRap(jets[0], jets[1]) > 2.4 )  _hist->fill(2);
 95    }
 96
 97    /// Normalise histograms etc., after the run
 98    void finalize() {
 99
100      const double scalefactor( crossSection()/picobarn / sumOfWeights() );
101      scale(_hist, scalefactor);
102
103    }
104    /// @}
105
106  private:
107
108    /// @name Histograms
109    /// @{
110    Histo1DPtr _hist;
111    /// @}
112
113  };
114
115
116  RIVET_DECLARE_PLUGIN(ATLAS_2014_I1298023);
117
118}