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
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/Beam.hh"
#include "Rivet/Projections/DISKinematics.hh"
#include "Rivet/Projections/FastJets.hh"
#include "fastjet/SISConePlugin.hh"

namespace Rivet {


  /// @brief ZEUS inclusive jet photoproduction study used to measure alpha_s
  ///
  /// @author Jon Butterworth
  class ZEUS_2012_I1116258 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ZEUS_2012_I1116258);

    /// @name Analysis methods
    //@{

    // Book projections and histograms
    void init() {

      // Projections

      // Jet schemes checked with original code, M.Wing, A.Geiser
      FinalState fs;
      double jet_radius = 1.0;
      declare(FastJets(fs, fastjet::JetAlgorithm::kt_algorithm, fastjet::RecombinationScheme::Et_scheme, jet_radius), "Jets");
      declare(FastJets(fs, fastjet::JetAlgorithm::antikt_algorithm, fastjet::RecombinationScheme::Et_scheme, jet_radius), "Jets_akt");

      // bit of messing about to use the correct recombnation scheme for SISCone.
      double overlap_threshold = 0.75;
      fastjet::SISConePlugin * plugin = new fastjet::SISConePlugin(jet_radius, overlap_threshold);
      plugin->set_use_jet_def_recombiner(true);
      JetDefinition siscone(plugin);
      siscone.set_recombination_scheme(fastjet::RecombinationScheme::Et_scheme);
      declare(FastJets(fs, siscone), "Jets_sis");


      declare(DISKinematics(), "Kinematics");

      // all eta
      book(_h_etjet[0], 1, 1, 1);

      // two ET cuts.
      book(_h_etajet[0], 2, 1, 1);
      book(_h_etajet[1], 3, 1, 1);

      // in eta regions
      book(_h_etjet[1], 4, 1, 1);
      book(_h_etjet[2], 5, 1, 1);
      book(_h_etjet[3], 6, 1, 1);
      book(_h_etjet[4], 7, 1, 1);
      book(_h_etjet[5], 8, 1, 1);

      // antiKT
      book(_h_etjet[6], 9, 1, 1);
      book(_h_etajet[2], 11, 1, 1);

      // SiSCone
      book(_h_etjet[7], 10, 1, 1);
      book(_h_etajet[3], 12, 1, 1);

    }


    // Do the analysis
    void analyze(const Event& event) {

      // Determine kinematics, including event orientation since ZEUS coord system is for +z = proton direction
      const DISKinematics& kin = apply<DISKinematics>(event, "Kinematics");
      const int orientation = kin.orientation();

      // Q2 and inelasticity cuts
      if (kin.Q2() > 1*GeV2) vetoEvent;
      if (!inRange(sqrt(kin.W2()), 142.0, 293.0)) vetoEvent;

      // Jet selection
      // @TODO check the recombination scheme
      const Jets jets = apply<FastJets>(event, "Jets")          \
        .jets(Cuts::Et > 17*GeV && Cuts::etaIn(-1*orientation, 2.5*orientation), cmpMomByEt);
      MSG_DEBUG("kT Jet multiplicity = " << jets.size());

      const Jets jets_akt = apply<FastJets>(event, "Jets_akt")      \
        .jets(Cuts::Et > 17*GeV && Cuts::etaIn(-1*orientation, 2.5*orientation), cmpMomByEt);

      const Jets jets_sis = apply<FastJets>(event, "Jets_sis")  \
        .jets(Cuts::Et > 17*GeV && Cuts::etaIn(-1*orientation, 2.5*orientation), cmpMomByEt);


      // Fill histograms
      for (const Jet& jet : jets ){
        _h_etjet[0]->fill(jet.pt());
        _h_etajet[0]->fill(orientation*jet.eta());
        if (jet.pt()>21*GeV) {
          _h_etajet[1]->fill(orientation*jet.eta());
        }
        if (orientation*jet.eta() < 0) {
          _h_etjet[1]->fill(jet.pt());
        } else if (orientation*jet.eta() < 1) {
          _h_etjet[2]->fill(jet.pt());
        } else if (orientation*jet.eta() < 1.5) {
          _h_etjet[3]->fill(jet.pt());
        } else if (orientation*jet.eta() < 2) {
          _h_etjet[4]->fill(jet.pt());
        } else {
          _h_etjet[5]->fill(jet.pt());
        }
      }

      for (const Jet& jet : jets_akt ){
        _h_etjet[6]->fill(jet.pt());
        _h_etajet[2]->fill(orientation*jet.eta());
      }
      for (const Jet& jet : jets_sis ){
        _h_etjet[7]->fill(jet.pt());
        _h_etajet[3]->fill(orientation*jet.eta());
      }

    }


    // Finalize
    void finalize() {
      const double sf = crossSection()/picobarn/sumOfWeights();
      for( int i = 0; i < 8; i++ ) {
        scale(_h_etjet[i], sf);
      }
      for( int i = 0; i < 4; i++ ) {
        scale(_h_etajet[i], sf);
      }
    }

    //@}


  private:

    /// @name Histograms
    //@{
    Histo1DPtr _h_etjet[8], _h_etajet[4];
    //@}

  };


  RIVET_DECLARE_PLUGIN(ZEUS_2012_I1116258);

}