rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2014_I1300647

Measurement of $Z/\gamma^*$ boson $p_T$ at $\sqrt{s} = 7\text{TeV}$
Experiment: ATLAS (LHC)
Inspire ID: 1300647
Status: VALIDATED
Authors:
  • Elena Yatsenko
References: Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • $Z/\gamma^*$ production with decays to electrons and/or muons.

A measurement of the $Z/\gamma^*$ transverse momentum spectrum using ATLAS proton-proton collision data at a center of mass energy of $\sqrt{s} = 7{\text{TeV}}$ at the LHC. The measurement is performed in both $Z/\gamma^* \rightarrow ee$ and $Z/\gamma^* \rightarrow \mu \mu$ channels.

Source code: ATLAS_2014_I1300647.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/DileptonFinder.hh"
  5
  6/// @todo Include more projections as required, e.g. ChargedFinalState, FastJets, DileptonFinder...
  7
  8namespace Rivet {
  9
 10
 11  class ATLAS_2014_I1300647 : public Analysis {
 12  public:
 13
 14    /// Constructor
 15    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2014_I1300647);
 16
 17
 18  public:
 19
 20    /// @name Analysis methods
 21    /// @{
 22
 23    /// Book histograms and initialise projections before the run
 24    void init() {
 25
 26      DileptonFinder zfinder_dressed_el(91.2*GeV, 0.1, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV && Cuts::abspid == PID::ELECTRON, Cuts::massIn(66.0*GeV, 116.0*GeV));
 27      declare(zfinder_dressed_el, "DileptonFinder_dressed_el");
 28
 29      DileptonFinder zfinder_bare_el(91.2*GeV, 0.0, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV && Cuts::abspid == PID::ELECTRON, Cuts::massIn(66.0*GeV, 116.0*GeV));
 30      declare(zfinder_bare_el,	"DileptonFinder_bare_el");
 31
 32      DileptonFinder zfinder_dressed_mu(91.2*GeV, 0.1, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV && Cuts::abspid == PID::MUON, Cuts::massIn(66.0*GeV, 116.0*GeV));
 33      declare(zfinder_dressed_mu, "DileptonFinder_dressed_mu");
 34
 35      DileptonFinder zfinder_bare_mu(91.2*GeV, 0.0, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV && Cuts::abspid == PID::MUON, Cuts::massIn(66.0*GeV, 116.0*GeV));
 36      declare(zfinder_bare_mu,	"DileptonFinder_bare_mu");
 37
 38      // Book histograms
 39      book(_hist_zpt_el_dressed ,1, 1, 1);  // electron "dressed"
 40      book(_hist_zpt_mu_dressed ,1, 1, 2);  // muon "dressed"
 41      book(_hist_zpt_el_bare    ,1, 2, 1);  // electron "bare"
 42      book(_hist_zpt_mu_bare    ,1, 2, 2);  // muon "bare"
 43
 44      //double-differential plots
 45      book(_h_zpt_el_mu_dressed, {0., 1., 2., 2.4}, {"d03-x01-y02", "d03-x01-y04", "d03-x01-y06"});
 46
 47    }
 48
 49
 50    /// Perform the per-event analysis
 51    void analyze(const Event& event) {
 52
 53      const DileptonFinder& zfinder_dressed_el = apply<DileptonFinder>(event, "DileptonFinder_dressed_el");
 54      const DileptonFinder& zfinder_bare_el    = apply<DileptonFinder>(event, "DileptonFinder_bare_el");
 55      const DileptonFinder& zfinder_dressed_mu = apply<DileptonFinder>(event, "DileptonFinder_dressed_mu");
 56      const DileptonFinder& zfinder_bare_mu    = apply<DileptonFinder>(event, "DileptonFinder_bare_mu");
 57
 58      FillPlots1d(zfinder_dressed_el, _hist_zpt_el_dressed);
 59
 60      FillPlots1d(zfinder_bare_el,    _hist_zpt_el_bare);
 61
 62      FillPlots1d(zfinder_dressed_mu, _hist_zpt_mu_dressed);
 63
 64      FillPlots1d(zfinder_bare_mu,    _hist_zpt_mu_bare);
 65
 66      FillPlots3d(zfinder_dressed_el, _h_zpt_el_mu_dressed);
 67      FillPlots3d(zfinder_dressed_mu, _h_zpt_el_mu_dressed);
 68
 69    }
 70
 71    void FillPlots1d(const DileptonFinder& zfinder, Histo1DPtr hist) {
 72      if(zfinder.bosons().size() != 1) return;
 73      const FourMomentum pZ = zfinder.bosons()[0].momentum();
 74      hist->fill(pZ.pT()/GeV);
 75      return;
 76    }
 77
 78    void FillPlots3d(const DileptonFinder& zfinder, Histo1DGroupPtr& binnedHist) {
 79      if(zfinder.bosons().size() != 1) return;
 80      const FourMomentum pZ = zfinder.bosons()[0].momentum();
 81      binnedHist->fill(pZ.rapidity(), pZ.pT()/GeV);
 82      return;
 83    }
 84
 85    /// Normalise histograms etc., after the run
 86    void finalize() {
 87
 88      normalize(_hist_zpt_el_dressed);
 89      normalize(_hist_zpt_el_bare);
 90
 91      normalize(_hist_zpt_mu_dressed);
 92      normalize(_hist_zpt_mu_bare);
 93
 94      normalize(_h_zpt_el_mu_dressed);
 95
 96    }
 97
 98    /// @}
 99
100
101  private:
102
103    // Data members like post-cuts event weight counters go here
104
105
106  private:
107
108    /// @name Histograms
109    /// @{
110    Histo1DGroupPtr _h_zpt_el_mu_dressed;
111
112
113    Histo1DPtr _hist_zpt_el_dressed;
114    Histo1DPtr _hist_zpt_el_bare;
115    Histo1DPtr _hist_zpt_mu_dressed;
116    Histo1DPtr _hist_zpt_mu_bare;
117
118    /// @}
119
120  };
121
122
123  RIVET_DECLARE_PLUGIN(ATLAS_2014_I1300647);
124}