Rivet Analyses Reference

ATLAS_2016_I1467454

High-mass Drell-Yan at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1467454
Status: VALIDATED
Authors:
  • Markus Zinser
References:Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • Drell-Yan production in $pp$ collisions at 8 TeV, leptonic $Z$ decays

This paper presents a measurement of the double-differential cross section for the Drell-Yan $Z/\gamma^\ast \rightarrow \ell^+\ell^-$ and photon-induced $\gamma\gamma \rightarrow \ell^+\ell^-$ processes where $\ell$ is an electron or muon. The measurement is performed for invariant masses of the lepton pairs, $m_{\ell\ell}$, between 116 GeV and 1500 GeV using a sample of 20.3 fb$^{-1}$ of pp collisions data at centre-of-mass energy of $\sqrt{s} = 8$ TeV collected by the ATLAS detector at the LHC in 2012. The data are presented double differentially in invariant mass and absolute dilepton rapidity as well as in invariant mass and absolute pseudorapidity separation of the lepton pair. The single-differential cross section as a function of $m_{\ell\ell}$ is also reported. The electron and muon channel measurements are combined and a total experimental precision of better than 1\% is achieved at low $m_{\ell\ell}$. A comparison to next-to-next-to-leading order perturbative QCD predictions using several recent parton distribution functions and including next-to-leading order electroweak effects indicates the potential of the data to constrain parton distribution functions. In particular, a large impact of the data on the photon PDF is demonstrated. Decay channel (e or mu) is selected by LMODE=EL,MU

Source code: ATLAS_2016_I1467454.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ZFinder.hh"
#include "Rivet/Tools/BinnedHistogram.hh"

namespace Rivet {


  /// @brief High-mass Drell-Yan at 8 TeV
  class ATLAS_2016_I1467454 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_I1467454);
    //@}


    /// @name Analysis methods
    //@{

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

      // Get options from the new option system
      _mode = 0;
      if ( getOption("LMODE") == "EL" ) _mode = 0;
      if ( getOption("LMODE") == "MU" ) _mode = 1;

      const FinalState fs;
      Cut cuts = Cuts::abseta < 2.5 && Cuts::pT > 30*GeV;
      ZFinder zfinder(fs, cuts, _mode? PID::MUON : PID::ELECTRON, 116*GeV, 1500*GeV, 0.1);
      declare(zfinder, "ZFinder");

      size_t ch = _mode? 11 : 0; // offset
      book(_hist_mll, 18 + ch, 1, 1);

      vector<double> mll_bins = { 116., 150., 200., 300., 500., 1500. };
      for (size_t i = 0; i < (mll_bins.size() - 1); ++i) {
        {Histo1DPtr tmp; _hist_rap.add( mll_bins[i], mll_bins[i+1], book(tmp, 19 + ch + i, 1, 1));}
        {Histo1DPtr tmp; _hist_deta.add(mll_bins[i], mll_bins[i+1], book(tmp, 24 + ch + i, 1, 1));}
      }

    }


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

      const ZFinder& zfinder = apply<ZFinder>(event, "ZFinder");
      if (zfinder.bosons().size() != 1)  vetoEvent;

      const Particle z0  = zfinder.bosons()[0];
      /// @todo Could use z0.constituents()
      const Particle el1 = zfinder.constituentLeptons()[0];
      const Particle el2 = zfinder.constituentLeptons()[1];

      if (el1.pT() > 40*GeV || el2.pT() > 40*GeV) {
        const double mass = z0.mass();
        _hist_mll->fill(mass/GeV);
        _hist_rap. fill(mass/GeV, z0.absrap());
        _hist_deta.fill(mass/GeV, deltaEta(el1,el2));
      }
    }


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

      const double sf = crossSection()/sumOfWeights();
      scale(_hist_mll, sf);
      _hist_rap.scale(sf*0.5,  this);
      _hist_deta.scale(sf*0.5, this);

    }

    //@}


    /// Choose to work in electron or muon mode
    size_t _mode;


    /// @name Histograms
    //@{
    Histo1DPtr _hist_mll;
    BinnedHistogram _hist_rap, _hist_deta;
    //@}

  };

  RIVET_DECLARE_PLUGIN(ATLAS_2016_I1467454);

}