Rivet Analyses Reference

ATLAS_2014_I1288706

Measurement of the low-mass Drell-Yan differential cross section at 7 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1288706
Status: VALIDATED
Authors:
  • Elena Yatsenko
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Run with inclusive $Z$ events, with $Z/\gamma^*$ decays to electrons and/or muons at $\sqrt{s} = 7$ TeV.

Measurements of the differential cross section for the process $Z/\gamma^\ast \rightarrow l$ $(l = e, \mu)$ as a function of dilepton invariant mass in pp collisions at $\sqrt{s} = 7$ TeV. The measurement is performed in the $e$ and $\mu$ channels for invariant masses between 26 GeV and 66 GeV in the fiducial region $p_\text{T}^\text{leading} > 15$ GeV, $p_\text{T}^\text{subleading} > 12$ GeV, $|\eta| < 2.4$ using an integrated luminosity of 1.6 $\text{fb}^{-1}$. The analysis is extended to invariant masses as low as 12 GeV in the muon channel within a fiducal region of $p_\text{T}^\text{leading} > 9$ GeV, $p_\text{T}^\text{subleading} > 6$ GeV, $|\eta| < 2.4$ with 35 $\text{pb}^{-1}$.

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

namespace Rivet {


  class ATLAS_2014_I1288706 : public Analysis {
  public:

    /// Constructor
    ATLAS_2014_I1288706()
      : Analysis("ATLAS_2014_I1288706")
    {
     }


    /// @name Analysis methods
    //@{

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

       // Set up projections
       FinalState fs;

       ZFinder zfinder_ext_dressed_mu(fs, Cuts::abseta<2.4 && Cuts::pT>6.0*GeV, PID::MUON, 12.0*GeV, 66.0*GeV, 0.1);
       declare(zfinder_ext_dressed_mu, "ZFinder_ext_dressed_mu");

       ZFinder zfinder_dressed_mu(fs, Cuts::abseta<2.4 && Cuts::pT>12*GeV, PID::MUON, 26.0*GeV, 66.0*GeV, 0.1);
       declare(zfinder_dressed_mu, "ZFinder_dressed_mu");

       ZFinder zfinder_dressed_el(fs, Cuts::abseta<2.4 && Cuts::pT>12*GeV, PID::ELECTRON, 26.0*GeV, 66.0*GeV, 0.1);
       declare(zfinder_dressed_el, "ZFinder_dressed_el");

       book(_hist_ext_mu_dressed ,1, 1, 1); // muon, dressed level, extended phase space
       book(_hist_mu_dressed	    ,2, 1, 1); // muon, dressed level, nominal phase space
       book(_hist_el_dressed	    ,2, 1, 2); // electron, dressed level, nominal phase space
    }


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


      const ZFinder& zfinder_ext_dressed_mu = apply<ZFinder>(event, "ZFinder_ext_dressed_mu");
      const ZFinder& zfinder_dressed_mu     = apply<ZFinder>(event, "ZFinder_dressed_mu"    );
      const ZFinder& zfinder_dressed_el     = apply<ZFinder>(event, "ZFinder_dressed_el"    );

      fillPlots(zfinder_ext_dressed_mu, _hist_ext_mu_dressed, 9*GeV);
      fillPlots(zfinder_dressed_mu,     _hist_mu_dressed,    15*GeV);
      fillPlots(zfinder_dressed_el,     _hist_el_dressed,    15*GeV);
    }


    /// Helper function to fill the histogram if a Z is found with the required lepton cuts
    void fillPlots(const ZFinder& zfinder, Histo1DPtr hist, double leading_pT) {
      if (zfinder.bosons().size() != 1) return;
      const FourMomentum el1 = zfinder.rawParticles()[0];
      const FourMomentum el2 = zfinder.rawParticles()[1];
      if (el1.pT() < leading_pT && el2.pT() < leading_pT) return;
      const double mass = zfinder.bosons()[0].mass();
      hist->fill(mass/GeV);
    }


    /// Normalise histograms etc., after the run
    void finalize() {
      scale(_hist_ext_mu_dressed, crossSection()/sumOfWeights());
      scale(_hist_mu_dressed,     crossSection()/sumOfWeights());
      scale(_hist_el_dressed,     crossSection()/sumOfWeights());
    }

    //@}


  private:

    /// Histograms
    Histo1DPtr _hist_ext_mu_dressed, _hist_mu_dressed, _hist_el_dressed;

  };


  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2014_I1288706);

}