rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CDF_2012_I1124333

CDF measurement of the Z pT in the electron channel using 2.1 fb-1
Experiment: CDF (Tevatron)
Inspire ID: 1124333
Status: VALIDATED
Authors:
  • Simone Amoroso
References: Beams: p- p+
Beam energies: (980.0, 980.0) GeV
Run details:
  • $Z\to ee$ event, a cut on the dilepton invariant mass between 60 and 120 GeV can be added to increase statistics

Measurement of the Z boson pT in the electron channel performed by the CDF experiment, using ppbar at 1960 GeV. A cut on the dielectron invariant mass of $66 GeV < m_{ee} < 116 GeV$ is applied.

Source code: CDF_2012_I1124333.cc
 1// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Projections/DileptonFinder.hh"
 4
 5namespace Rivet {
 6
 7
 8  /// @ CDF Run II Z \f$ p_\perp \f$ in Drell-Yan events
 9  /// @author Simone Amoroso
10  class CDF_2012_I1124333 : public Analysis {
11  public:
12
13    /// Constructor
14    RIVET_DEFAULT_ANALYSIS_CTOR(CDF_2012_I1124333);
15
16
17    /// @name Analysis methods
18    /// @{
19
20    /// Book histograms and initialise projections before the run
21    void init() {
22
23      ///  Initialise and register projections
24      DileptonFinder zfinder(91.2*GeV, 0.0, Cuts::abspid == PID::ELECTRON, Cuts::massIn(66*GeV, 116*GeV));
25      declare(zfinder, "DileptonFinder");
26
27      ///  Book histograms
28      book(_hist_zpt, 2, 1, 1);
29    }
30
31
32    /// Perform the per-event analysis
33    void analyze(const Event& event) {
34      const DileptonFinder& zfinder = apply<DileptonFinder>(event, "DileptonFinder");
35      if (zfinder.bosons().size() != 1) {
36        MSG_DEBUG("Num e+ e- pairs found = " << zfinder.bosons().size());
37        vetoEvent;
38      }
39
40      const FourMomentum& pZ = zfinder.bosons()[0].momentum();
41      if (pZ.mass2() < 0) {
42        MSG_DEBUG("Negative Z mass**2 = " << pZ.mass2()/GeV2 << "!");
43        vetoEvent;
44      }
45
46      MSG_DEBUG("Dilepton mass = " << pZ.mass()/GeV << " GeV");
47      _hist_zpt->fill(pZ.pT());
48      //      _hist_z_xs->fill(1);
49    }
50
51
52    /// Normalise histograms etc., after the run
53    void finalize() {
54      scale(_hist_zpt, crossSection()/picobarn/sumOfWeights());
55    }
56
57    /// @}
58
59
60  private:
61
62    /// @name Histograms
63    /// @{
64    Histo1DPtr _hist_zpt;
65    //    Histo1DPtr _hist_z_xs;
66    /// @}
67
68  };
69
70
71  RIVET_DECLARE_PLUGIN(CDF_2012_I1124333);
72
73}