rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CDF_2007_I743342

CDF Run II inclusive jet cross-section using the kT algorithm
Experiment: CDF (Tevatron Run 2)
Inspire ID: 743342
Status: VALIDATED
Authors:
  • David Voong
  • Frank Siegert
References:
  • Phys.Rev.D75:092006,2007
  • Erratum-ibid.D75:119901,2007
  • FNAL-PUB 07/026-E
  • hep-ex/0701051
Beams: p- p+
Beam energies: (980.0, 980.0) GeV
Run details:
  • p-pbar collisions at 1960 GeV. Jet pT bins from 54 GeV to 700 GeV. Jet rapidity $< |2.1|$.

CDF Run II measurement of inclusive jet cross sections at a p-pbar collision energy of 1.96 TeV. Jets are reconstructed in the central part of the detector ($|y|<2.1$) using the kT algorithm with an $R$ parameter of 0.7. The minimum jet pT considered is 54 GeV, with a maximum around 700 GeV. The inclusive jet pT is plotted in bins of rapidity $|y|<0.1$, $0.1<|y|<0.7$, $0.7<|y|<1.1$, $1.1<|y|<1.6$ and $1.6<|y|<2.1$.

Source code: CDF_2007_I743342.cc
 1// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Projections/FastJets.hh"
 4
 5namespace Rivet {
 6
 7
 8  /// @brief CDF inclusive jet cross-section using the \f$ k_\perp \f$ algorithm
 9  class CDF_2007_I743342 : public Analysis {
10  public:
11
12    RIVET_DEFAULT_ANALYSIS_CTOR(CDF_2007_I743342);
13
14
15    void init() {
16      // Set up projections
17      const FinalState fs;
18      declare(FastJets(fs, JetAlg::KT, 0.5), "JetsD05");
19      declare(FastJets(fs, JetAlg::KT, 0.7), "JetsD07");
20      declare(FastJets(fs, JetAlg::KT, 1.0), "JetsD10");
21
22      // Book histos
23      book(_binnedHistosD07, {0., 0.1, 0.7, 1.1, 1.6, 2.1},
24                             {"d01-x01-y01", "d02-x01-y01", "d03-x01-y01", "d04-x01-y01", "d05-x01-y01"});
25      book(_histoD05, 6, 1, 1);
26      book(_histoD10, 7, 1, 1);
27    }
28
29
30    void analyze(const Event& event) {
31
32      for (const Jet& jet : apply<JetFinder>(event, "JetsD07").jets(Cuts::pT > 54*GeV)) {
33        _binnedHistosD07->fill(jet.absrap(), jet.pT());
34      }
35
36      for (const Jet& jet : apply<JetFinder>(event, "JetsD05").jets(Cuts::pT > 54*GeV)) {
37        if (inRange(jet.absrap(), 0.1, 0.7))  _histoD05->fill(jet.pT());
38      }
39
40      for (const Jet& jet : apply<JetFinder>(event, "JetsD10").jets(Cuts::pT > 54*GeV)) {
41        if (inRange(jet.absrap(), 0.1, 0.7))  _histoD10->fill(jet.pT());
42      }
43    }
44
45
46    // Normalise histograms to cross-section
47    void finalize() {
48      const double xSec = crossSectionPerEvent()/nanobarn;
49
50      scale(_histoD05, xSec);
51      scale(_histoD10, xSec);
52      // scale to xSec/yBinWidth and take into account the double yBinWidth due
53      // to the absolute value of y
54      scale(_binnedHistosD07, xSec/2.0);
55      divByGroupWidth(_binnedHistosD07);
56    }
57
58
59  private:
60
61    Histo1DGroupPtr _binnedHistosD07;
62
63    /// Single histogram for the \f$R=0.5\f$ \f$k_\perp\f$ jets
64    Histo1DPtr _histoD05;
65
66    /// Single histogram for the \f$R=1.0\f$ \f$k_\perp\f$ jets
67    Histo1DPtr _histoD10;
68
69  };
70
71
72
73  RIVET_DECLARE_ALIASED_PLUGIN(CDF_2007_I743342, CDF_2007_S7057202);
74
75}