rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

D0_2009_I826756

Z+jets angular distributions
Experiment: D0 (Tevatron Run 2)
Inspire ID: 826756
Status: VALIDATED
Authors:
  • Frank Siegert
References: Beams: p- p+
Beam energies: (980.0, 980.0) GeV
Run details:
  • $p \bar{p} \to \mu^+ \mu^-$ + jets at 1960 GeV. Needs mass cut on lepton pair to avoid photon singularity, looser than $65 < m_{ee} < 115$ GeV.

First measurements at a hadron collider of differential cross sections for $Z (\to \mu\mu)$+jet+X production in $\Delta\phi(Z, j)$, $|\Delta y(Z, j)|$ and $|y_\mathrm{boost}(Z, j)|$. Vector boson production in association with jets is an excellent probe of QCD and constitutes the main background to many small cross section processes, such as associated Higgs production. These measurements are crucial tests of the predictions of perturbative QCD and current event generators, which have varied success in describing the data. Using these measurements as inputs in tuning event generators will increase the experimental sensitivity to rare signals.

Source code: D0_2009_I826756.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/DileptonFinder.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6
  7namespace Rivet {
  8
  9
 10  /// D0 Z+jets angular distributions
 11  class D0_2009_I826756 : public Analysis {
 12  public:
 13
 14    RIVET_DEFAULT_ANALYSIS_CTOR(D0_2009_I826756);
 15
 16
 17    /// @name Analysis methods
 18    /// @{
 19
 20    /// Book histograms
 21    void init() {
 22      Cut cut = Cuts::abseta < 1.7 && Cuts::pT > 15*GeV;
 23      DileptonFinder zfinder(91.2*GeV, 0.2, cut && Cuts::abspid == PID::MUON, Cuts::massIn(65*GeV, 115*GeV));
 24      declare(zfinder, "DileptonFinder");
 25
 26      FastJets conefinder(zfinder.remainingFinalState(), JetAlg::D0ILCONE, 0.5);
 27      declare(conefinder, "ConeFinder");
 28
 29      book(_h_dphi_jet_Z25 ,1, 1, 1);
 30      book(_h_dphi_jet_Z45 ,2, 1, 1);
 31
 32      book(_h_dy_jet_Z25 ,3, 1, 1);
 33      book(_h_dy_jet_Z45 ,4, 1, 1);
 34
 35      book(_h_yboost_jet_Z25 ,5, 1, 1);
 36      book(_h_yboost_jet_Z45 ,6, 1, 1);
 37
 38      book(_h_dphi_jet_Z25_xs ,1, 1, 2);
 39      book(_h_dphi_jet_Z45_xs ,2, 1, 2);
 40
 41      book(_h_dy_jet_Z25_xs ,3, 1, 2);
 42      book(_h_dy_jet_Z45_xs ,4, 1, 2);
 43
 44      book(_h_yboost_jet_Z25_xs ,5, 1, 2);
 45      book(_h_yboost_jet_Z45_xs ,6, 1, 2);
 46
 47      book(_inclusive_Z_sumofweights, "_inclusive_Z_sumofweights");
 48    }
 49
 50
 51    void analyze(const Event& event) {
 52      const DileptonFinder& zfinder = apply<DileptonFinder>(event, "DileptonFinder");
 53      if (zfinder.bosons().size() == 1) {
 54        // count inclusive sum of weights for histogram normalisation
 55        _inclusive_Z_sumofweights->fill();
 56
 57        const FourMomentum& zmom = zfinder.bosons()[0].momentum();
 58        if (zmom.pT() < 25*GeV) vetoEvent;
 59
 60        Jets jets = apply<JetFinder>(event, "ConeFinder").jetsByPt(Cuts::pT > 20*GeV && Cuts::abseta < 2.8);
 61
 62        // Return if there are no jets:
 63        if (jets.size() < 1) {
 64          MSG_DEBUG("Skipping event " << numEvents() << " because no jets pass cuts ");
 65          vetoEvent;
 66        }
 67
 68        const FourMomentum& jetmom = jets[0].momentum();
 69        const double yZ = zmom.rapidity();
 70        const double yjet = jetmom.rapidity();
 71        const double dphi = deltaPhi(zmom, jetmom);
 72        const double dy = deltaRap(zmom, jetmom);
 73        const double yboost = fabs(yZ+yjet)/2;
 74
 75        if (zmom.pT() > 25*GeV) {
 76          _h_dphi_jet_Z25->fill(dphi);
 77          _h_dy_jet_Z25->fill(dy);
 78          _h_yboost_jet_Z25->fill(yboost);
 79          _h_dphi_jet_Z25_xs->fill(dphi);
 80          _h_dy_jet_Z25_xs->fill(dy);
 81          _h_yboost_jet_Z25_xs->fill(yboost);
 82        }
 83        if (zmom.pT() > 45*GeV) {
 84          _h_dphi_jet_Z45->fill(dphi);
 85          _h_dy_jet_Z45->fill(dy);
 86          _h_yboost_jet_Z45->fill(yboost);
 87          _h_dphi_jet_Z45_xs->fill(dphi);
 88          _h_dy_jet_Z45_xs->fill(dy);
 89          _h_yboost_jet_Z45_xs->fill(yboost);
 90        }
 91      }
 92
 93    }
 94
 95
 96    void finalize() {
 97      if (_inclusive_Z_sumofweights->val() == 0) return;
 98      scale(_h_dphi_jet_Z25, 1/ *_inclusive_Z_sumofweights);
 99      scale(_h_dphi_jet_Z45, 1/ *_inclusive_Z_sumofweights);
100      scale(_h_dy_jet_Z25, 1/ *_inclusive_Z_sumofweights);
101      scale(_h_dy_jet_Z45, 1/ *_inclusive_Z_sumofweights);
102      scale(_h_yboost_jet_Z25, 1/ *_inclusive_Z_sumofweights);
103      scale(_h_yboost_jet_Z45, 1/ *_inclusive_Z_sumofweights);
104
105      scale(_h_dphi_jet_Z25_xs, crossSectionPerEvent());
106      scale(_h_dphi_jet_Z45_xs, crossSectionPerEvent());
107      scale(_h_dy_jet_Z25_xs, crossSectionPerEvent());
108      scale(_h_dy_jet_Z45_xs, crossSectionPerEvent());
109      scale(_h_yboost_jet_Z25_xs, crossSectionPerEvent());
110      scale(_h_yboost_jet_Z45_xs, crossSectionPerEvent());
111    }
112
113    /// @}
114
115
116  private:
117
118    /// @name Histograms (normalised)
119    /// @{
120    Histo1DPtr _h_dphi_jet_Z25;
121    Histo1DPtr _h_dphi_jet_Z45;
122
123    Histo1DPtr _h_dy_jet_Z25;
124    Histo1DPtr _h_dy_jet_Z45;
125
126    Histo1DPtr _h_yboost_jet_Z25;
127    Histo1DPtr _h_yboost_jet_Z45;
128    /// @}
129
130    /// @name Histograms (absolute cross sections)
131    /// @{
132    Histo1DPtr _h_dphi_jet_Z25_xs;
133    Histo1DPtr _h_dphi_jet_Z45_xs;
134
135    Histo1DPtr _h_dy_jet_Z25_xs;
136    Histo1DPtr _h_dy_jet_Z45_xs;
137
138    Histo1DPtr _h_yboost_jet_Z25_xs;
139    Histo1DPtr _h_yboost_jet_Z45_xs;
140    /// @}
141
142    CounterPtr _inclusive_Z_sumofweights;
143
144  };
145
146
147
148  RIVET_DECLARE_ALIASED_PLUGIN(D0_2009_I826756, D0_2009_S8349509);
149
150}