rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ZEUS_2012_I1116258

Inclusive-jet photoproduction at HERA and determination of the strong coupling
Experiment: ZEUS (HERA Run II)
Inspire ID: 1116258
Status: UNVALIDATED
Authors:
  • Jon Butterworth
References:
  • Nucl.Phys. B864 (2012) 1-37
  • DESY 12/045
  • hep-ex/1205.6153
Beams: p+ e+
Beam energies: (920.0, 27.5) GeV
Run details:
  • 920 GeV protons colliding with 27.5 GeV positrons; Direct and resolved photoproduction of dijets; Jet $pT > 17$ GeV Jet pseudorapidity $-1 < |\eta| < 2.5$

Inclusive-jet cross sections have been measured in the reaction ep->e+jet+X for photon virtuality Q2 < 1 GeV2 and gamma-p centre-of-mass energies in the region 142 < W(gamma-p) < 293 GeV with the ZEUS detector at HERA using an integrated luminosity of 300 pb-1. Jets were identified using the kT, anti-kT or SIScone jet algorithms in the laboratory frame. Single-differential cross sections are presented as functions of the jet transverse energy, ETjet, and pseudorapidity, etajet, for jets with ETjet > 17 GeV and -1 < etajet < 2.5. In addition, measurements of double-differential inclusive-jet cross sections are presented as functions of ETjet in different regions of etajet. Next-to-leading-order QCD calculations give a good description of the measurements, except for jets with low ETjet and high etajet. The influence of non-perturbative effects not related to hadronisation was studied. Measurements of the ratios of cross sections using different jet algorithms are also presented; the measured ratios are well described by calculations including up to O(alphas2) terms. Values of alphas(Mz) were extracted from the measurements and the energy-scale dependence of the coupling was determined. The value of alphas(Mz) extracted from the measurements based on the kT jet algorithm is alphas(Mz) = 0.1206 +0.0023 -0.0022 (exp.) +0.0042 -0.0035 (th.); the results from the anti-kT and SIScone algorithms are compatible with this value and have a similar precision.

Source code: ZEUS_2012_I1116258.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/Beam.hh"
  4#include "Rivet/Projections/DISKinematics.hh"
  5#include "Rivet/Projections/FastJets.hh"
  6#include "fastjet/SISConePlugin.hh"
  7
  8namespace Rivet {
  9
 10
 11  /// @brief ZEUS inclusive jet photoproduction study used to measure alpha_s
 12  ///
 13  /// @author Jon Butterworth
 14  class ZEUS_2012_I1116258 : public Analysis {
 15  public:
 16
 17    /// Constructor
 18    RIVET_DEFAULT_ANALYSIS_CTOR(ZEUS_2012_I1116258);
 19
 20    /// @name Analysis methods
 21    /// @{
 22
 23    // Book projections and histograms
 24    void init() {
 25
 26      // Projections
 27
 28      // Jet schemes checked with original code, M.Wing, A.Geiser
 29      FinalState fs;
 30      double jet_radius = 1.0;
 31      declare(FastJets(fs, fastjet::JetAlgorithm::kt_algorithm, fastjet::RecombinationScheme::Et_scheme, jet_radius), "Jets");
 32      declare(FastJets(fs, fastjet::JetAlgorithm::antikt_algorithm, fastjet::RecombinationScheme::Et_scheme, jet_radius), "Jets_akt");
 33
 34      // bit of messing about to use the correct recombnation scheme for SISCone.
 35      double overlap_threshold = 0.75;
 36      fastjet::SISConePlugin * plugin = new fastjet::SISConePlugin(jet_radius, overlap_threshold);
 37      plugin->set_use_jet_def_recombiner(true);
 38      JetDefinition siscone(plugin);
 39      siscone.set_recombination_scheme(fastjet::RecombinationScheme::Et_scheme);
 40      declare(FastJets(fs, siscone), "Jets_sis");
 41
 42
 43      declare(DISKinematics(), "Kinematics");
 44
 45      // all eta
 46      book(_h_etjet[0], 1, 1, 1);
 47
 48      // two ET cuts.
 49      book(_h_etajet[0], 2, 1, 1);
 50      book(_h_etajet[1], 3, 1, 1);
 51
 52      // in eta regions
 53      book(_h_etjet[1], 4, 1, 1);
 54      book(_h_etjet[2], 5, 1, 1);
 55      book(_h_etjet[3], 6, 1, 1);
 56      book(_h_etjet[4], 7, 1, 1);
 57      book(_h_etjet[5], 8, 1, 1);
 58
 59      // antiKT
 60      book(_h_etjet[6], 9, 1, 1);
 61      book(_h_etajet[2], 11, 1, 1);
 62
 63      // SiSCone
 64      book(_h_etjet[7], 10, 1, 1);
 65      book(_h_etajet[3], 12, 1, 1);
 66
 67    }
 68
 69
 70    // Do the analysis
 71    void analyze(const Event& event) {
 72
 73      // Determine kinematics, including event orientation since ZEUS coord system is for +z = proton direction
 74      const DISKinematics& kin = apply<DISKinematics>(event, "Kinematics");
 75      const int orientation = kin.orientation();
 76
 77      // Q2 and inelasticity cuts
 78      if (kin.Q2() > 1*GeV2) vetoEvent;
 79      if (!inRange(sqrt(kin.W2()), 142.0, 293.0)) vetoEvent;
 80
 81      // Jet selection
 82      /// @todo check the recombination scheme
 83      const Jets jets = apply<FastJets>(event, "Jets")          \
 84        .jets(Cuts::Et > 17*GeV && Cuts::etaIn(-1*orientation, 2.5*orientation), cmpMomByEt);
 85      MSG_DEBUG("kT Jet multiplicity = " << jets.size());
 86
 87      const Jets jets_akt = apply<FastJets>(event, "Jets_akt")		\
 88        .jets(Cuts::Et > 17*GeV && Cuts::etaIn(-1*orientation, 2.5*orientation), cmpMomByEt);
 89
 90      const Jets jets_sis = apply<FastJets>(event, "Jets_sis")          \
 91        .jets(Cuts::Et > 17*GeV && Cuts::etaIn(-1*orientation, 2.5*orientation), cmpMomByEt);
 92
 93
 94      // Fill histograms
 95      for (const Jet& jet : jets ){
 96        _h_etjet[0]->fill(jet.pt());
 97        _h_etajet[0]->fill(orientation*jet.eta());
 98        if (jet.pt()>21*GeV) {
 99          _h_etajet[1]->fill(orientation*jet.eta());
100        }
101        if (orientation*jet.eta() < 0) {
102          _h_etjet[1]->fill(jet.pt());
103        } else if (orientation*jet.eta() < 1) {
104          _h_etjet[2]->fill(jet.pt());
105        } else if (orientation*jet.eta() < 1.5) {
106          _h_etjet[3]->fill(jet.pt());
107        } else if (orientation*jet.eta() < 2) {
108          _h_etjet[4]->fill(jet.pt());
109        } else {
110          _h_etjet[5]->fill(jet.pt());
111        }
112      }
113
114      for (const Jet& jet : jets_akt ){
115        _h_etjet[6]->fill(jet.pt());
116        _h_etajet[2]->fill(orientation*jet.eta());
117      }
118      for (const Jet& jet : jets_sis ){
119        _h_etjet[7]->fill(jet.pt());
120        _h_etajet[3]->fill(orientation*jet.eta());
121      }
122
123    }
124
125
126    // Finalize
127    void finalize() {
128      const double sf = crossSection()/picobarn/sumOfWeights();
129      for( int i = 0; i < 8; i++ ) {
130        scale(_h_etjet[i], sf);
131      }
132      for( int i = 0; i < 4; i++ ) {
133        scale(_h_etajet[i], sf);
134      }
135    }
136
137    /// @}
138
139
140  private:
141
142    /// @name Histograms
143    /// @{
144    Histo1DPtr _h_etjet[8], _h_etajet[4];
145    /// @}
146
147  };
148
149
150  RIVET_DECLARE_PLUGIN(ZEUS_2012_I1116258);
151
152}