rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2016_I1487726

Collinear W emissions at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1487726
Status: VALIDATED
Authors:
  • marek.schoenherr@cern.ch
  • mileswu@uchicago.edu
  • chris.g@cern.ch
References:
  • DOI:10.1016/j.physletb.2016.12.005
  • arXiv: 1609.07045
  • Phys.Lett. B765 (2017) 132-153
Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • 8 TeV $pp \to W+\text{jets} \to \mu\nu+\text{jets}$.

The $W$ boson angular distribution in events with high transverse momentum jets is measured using data collected by the ATLAS experiment from proton-proton collisions at a centre-of-mass energy $\sqrt{s}=8$ TeV at the Large Hadron Collider, corresponding to an integrated luminosity of 20.3 fb$^{-1}$. The focus is on the contributions to $W+$jets processes from real $W$ emission, which is achieved by studying events where a muon is observed close to a high transverse momentum jet. At small angular separations, these contributions are expected to be large. Various theoretical models of this process are compared to the data in terms of the absolute cross-section and the angular distributions of the muon from the leptonic $W$ decay.

Source code: ATLAS_2016_I1487726.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/DressedLeptons.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6
  7namespace Rivet {
  8
  9
 10class ATLAS_2016_I1487726 : public Analysis {
 11
 12    public:
 13
 14    /// Constructor
 15    /// @brief Collinear W emissions at 8 TeV
 16    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_I1487726);
 17
 18    public:
 19
 20        /// @name Analysis methods
 21        //@{
 22
 23        /// Book histograms and initialise projections before the run
 24        void init() {
 25
 26            _mode = 0;
 27            if ( getOption("LMODE") == "EL" )  _mode = 1;
 28
 29            // These really should include non-prompt leptons/photons
 30            FinalState mufs(Cuts::abspid == PID::MUON);
 31            FinalState elfs(Cuts::abspid == PID::ELECTRON);
 32            FinalState phs(Cuts::abspid == PID::PHOTON);
 33
 34            Cut lep_fid = (Cuts::abseta < 2.4 && Cuts::pT >= 25*GeV);
 35            DressedLeptons dlep(phs, _mode? elfs : mufs, 0.1, lep_fid, true);
 36            declare(dlep, "DressedLeptons");
 37
 38            FastJets fj(FinalState(), FastJets::ANTIKT, 0.4, JetAlg::Muons::NONE, JetAlg::Invisibles::NONE);
 39            declare(fj, "AntiKt4Jets");
 40
 41            book(h_mu_jet_dr,          2, 1, 1);
 42            book(h_mu_jet_dr_pt500600, 4, 1, 1);
 43            book(h_mu_jet_dr_pt650,    5, 1, 1);
 44        }
 45
 46
 47        /// Perform the per-event analysis
 48        void analyze(const Event& event) {
 49
 50          const vector<DressedLepton> leptons = apply<DressedLeptons>(event, "DressedLeptons").dressedLeptons();
 51          const Jets jets = apply<FastJets>(event, "AntiKt4Jets").jetsByPt(Cuts::pT >= 100*GeV && Cuts::abseta <= 2.1);
 52
 53          if (leptons.size() != 1)       vetoEvent;
 54          if (jets.size() < 1)           vetoEvent;
 55          if (jets[0].pt() < 500.0*GeV)  vetoEvent;
 56
 57          // find closest jet to the lepton.
 58          Jet jet;
 59          double drmin = 999;
 60          for (const Jet &j : jets) {
 61            double dr = deltaR(leptons[0], j);
 62            if (dr < drmin) {
 63                drmin = dr;
 64                jet = j;
 65            }
 66          }
 67
 68          h_mu_jet_dr->fill(drmin);
 69          if (jets[0].pT() > 650*GeV)  h_mu_jet_dr_pt650->fill(drmin);
 70          else if (jets[0].pT() > 500*GeV && jets[0].pT() < 600*GeV) {
 71            h_mu_jet_dr_pt500600->fill(drmin);
 72          }
 73        }
 74
 75        /// Normalise histograms etc., after the run
 76        void finalize() {
 77          const double sf = crossSection() / femtobarn / sumOfWeights();
 78          scale(h_mu_jet_dr, sf);
 79          scale(h_mu_jet_dr_pt500600, sf);
 80          scale(h_mu_jet_dr_pt650, sf);
 81        }
 82
 83        //@}
 84
 85    protected:
 86
 87        size_t _mode;
 88
 89    private:
 90
 91        /// @name Histograms
 92        //@{
 93        Histo1DPtr h_mu_jet_dr;
 94        Histo1DPtr h_mu_jet_dr_pt500600;
 95        Histo1DPtr h_mu_jet_dr_pt650;
 96        //@}
 97  };
 98
 99
100  // The hook for the plugin system
101  RIVET_DECLARE_PLUGIN(ATLAS_2016_I1487726);
102}