rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

HERA_2015_I1377206

Measurement of sigma_red (F2) of H1 and ZEUS at different beam energies
Experiment: H1 (HERA)
Inspire ID: 1377206
Status: VALIDATED
Authors:
  • Hannes Jung
  • Tymoteusz Strozniak
References: Beams: e+ p+, p+ e+, e- p+, p+ e-
Beam energies: (27.5, 920.0); (920.0, 27.5); (27.5, 820.0); (820.0, 27.5); (27.5, 575.0); (575.0, 27.5); (27.5, 460.0); (460.0, 27.5) GeV
Run details:
  • Cover full range in Bjorken-y

A combination is presented of all inclusive deep inelastic cross sections previously published by the H1 and ZEUS collaborations at HERA for neutral and charged current ep scattering for zero beam polarisation. The data were taken at proton beam energies of 920, 820, 575 and 460 GeV and an electron beam energy of 27.5GeV. The data correspond to an integrated luminosity of about 1 fb^-1 and span six orders of magnitude in negative four-momentum-transfer squared, Q2, and Bjorken x . The correlations of the systematic uncertainties were evaluated and taken into account for the combination. The combined cross sections were input to QCD analyses at leading order, next-to-leading order and at next-to-next-to-leading order, providing a new set of parton distribution functions, called HERAPDF2.0.. The energy and beam options need to be specified when using rivet-merge.

Source code: HERA_2015_I1377206.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/FastJets.hh"
  5#include "Rivet/Projections/DISKinematics.hh"
  6#include "Rivet/Projections/Beam.hh"
  7
  8namespace Rivet {
  9
 10
 11  /// @brief Measurement of sigma_red (F2) of H1 and ZEUS at different beam energies
 12  class HERA_2015_I1377206 : public Analysis {
 13  public:
 14
 15    /// Constructor
 16    RIVET_DEFAULT_ANALYSIS_CTOR(HERA_2015_I1377206);
 17
 18
 19    /// @name Analysis methods
 20    /// @{
 21
 22    /// Book histograms and initialise projections before the run
 23    void init() {
 24
 25      // Initialise and register projections
 26      declare(FinalState(Cuts::abseta < 5 && Cuts::pT > 100*MeV), "FS");
 27      declare(DISLepton(), "Lepton");
 28      declare(DISKinematics(), "Kinematics");
 29
 30      string beamOpt = getOption<string>("BEAM","NONE");
 31      if (beamOpt == "NONE") {
 32        const ParticlePair& beam = beams();
 33        _positron = (beam.first.pid() == PID::POSITRON || beam.second.pid() == PID::POSITRON);
 34      }
 35      else {
 36        if (beamOpt == "EMINUS")      _positron = false;
 37        else if (beamOpt == "EPLUS")  _positron = true;
 38        else {
 39          throw BeamError("Beam species not supported\n");
 40        }
 41      }
 42
 43      // Book beams-dependent histograms
 44      const double eps = 1e-2;
 45      for (double eVal : allowedEnergies()) {
 46        const string en = toString(int(eVal)) + toString(_positron);
 47        if (isCompatibleWithSqrtS(eVal, eps))  _sqs = en;
 48
 49        if (fuzzyEquals(eVal, 318.1, eps) && _positron) {
 50          // NC e+ p at sqrts=318
 51          const vector<double> Q2edges = {
 52            0.1, 0.15, 0.2, 0.3, 0.4, 0.45, 0.6, 0.7, 1.0, 1.1, 1.3, 1.7, 2.3,
 53            3.1, 3.8, 5.3, 8.0, 9.1, 11., 13., 17.4, 19.1, 25.8, 28., 30., 42.,
 54            49., 54., 65., 75., 108., 134., 180., 225., 280., 325., 355., 455., 460.,
 55            545., 560., 765., 770., 835., 900., 1120., 1295., 1300., 1755., 1800.,
 56            2270., 2500., 3685., 4000., 6520., 7000., 9275., 10000., 15000., 17000.,
 57            24770., 25000., 42000.
 58          };
 59          book(_g[en+"sigred"], Q2edges);
 60          _g[en+"sigred"]->maskBins({9, 24, 27, 36, 38, 40, 42, 44, 47, 49, 51, 53, 55, 57, 59, 61});
 61          size_t idx = 0;
 62          for (auto& b : _g[en+"sigred"]->bins()) {
 63            book(b, 1, 1, ++idx);
 64          }
 65          // CC e+ p at sqrts=318
 66          book(_g[en+"sigred_cc"], {280., 325., 460., 545., 900., 1120., 1300., 1755., 1800., 2270.,
 67                                2500., 3685., 4000., 6520., 7000., 9275., 10000., 20000., 42000.});
 68          _g[en+"sigred_cc"]->maskBins({2, 4, 6, 8, 10, 12, 14, 16});
 69          idx = 0;
 70          for (auto& b : _g[en+"sigred_cc"]->bins()) {
 71            book(b, 6, 1, ++idx);
 72          }
 73        }
 74        else if (fuzzyEquals(eVal, 300.3, eps) && _positron) {
 75          // NC e+ p at sqrts=300
 76          const vector<double> Q2edges = {
 77            0.01, 0.05, 0.07, 0.09, 0.12, 0.18, 0.22, 0.32, 0.4, 0.45, 0.6, 0.7, 1.0, 1.1, 1.3,
 78            1.7, 2.3, 3.1, 3.8, 5.3, 8., 9.1, 11., 13., 17.4, 19.1, 25.8, 28., 30., 42., 49.,
 79            54., 65., 75., 108., 134., 180., 225., 280., 325., 355., 455., 460., 545., 560.,
 80            765., 770., 835., 900., 1120., 1295., 1300., 1755., 1800., 2270., 2500., 3685.,
 81            4000., 6520., 7000., 9275., 10000., 15000., 17000., 24770., 25000., 42000.
 82          };
 83          book(_g[en+"sigred"], Q2edges);
 84          _g[en+"sigred"]->maskBins({13, 28, 31, 40, 42, 44, 46, 48, 51, 53, 55, 57, 59, 61, 63, 65});
 85          size_t idx = 0;
 86          for (auto& b : _g[en+"sigred"]->bins()) {
 87            book(b, 2, 1, ++idx);
 88          }
 89        }
 90        else if (fuzzyEquals(eVal, 251.5, eps) && _positron) {
 91          // NC e+ p at sqrts=251
 92          const vector<double> Q2edges = {
 93            1., 1.7, 2.3, 3.1, 3.8, 5.3, 8., 9.1, 11., 13., 17.4, 22.1, 28.,
 94            30., 42., 49., 54., 65., 75., 108., 134., 180., 225., 280., 325.,
 95            355., 455., 460., 545., 560., 765., 770., 835.
 96          };
 97          book(_g[en+"sigred"], Q2edges);
 98          _g[en+"sigred"]->maskBins({8, 13, 16, 18, 25, 27, 29, 31});
 99          size_t idx = 0;
100          for (auto& b : _g[en+"sigred"]->bins()) {
101            book(b, 3, 1, ++idx);
102          }
103        }
104        else if (fuzzyEquals(eVal, 224.9, eps) && _positron) {
105          // NC e+ p at sqrts=225
106          const vector<double> Q2edges = {
107            1., 1.7, 2.3, 3.1, 3.8, 5.3, 8., 9.1, 11., 13., 17.4, 22.1, 28.,
108            30., 42., 49., 54., 65., 75., 108., 134., 180., 225., 280., 325.,
109            355., 455., 460., 545., 560., 765., 770., 835.
110          };
111          book(_g[en+"sigred"], Q2edges);
112          _g[en+"sigred"]->maskBins({8, 13, 16, 18, 25, 27, 29, 31});
113          size_t idx = 0;
114          for (auto& b : _g[en+"sigred"]->bins()) {
115            book(b, 4, 1, ++idx);
116          }
117        }
118        else if (fuzzyEquals(eVal, 318.1, eps) && !_positron) {
119          // NC e- p at sqrts=318
120          const vector<double> Q2edges = {
121            54., 65., 75., 108., 134., 180., 225., 280., 325., 355., 455.,
122            460., 545., 560., 765., 770., 835., 900., 1120., 1295., 1300.,
123            1755., 1800., 2270., 2500., 3685., 4000., 6520., 7000., 9275.,
124            10000., 15000., 17000., 24770., 25000., 42000., 70000.
125          };
126          book(_g[en+"sigred"], Q2edges);
127          _g[en+"sigred"]->maskBins({2, 9, 11, 13, 15, 17, 20, 22, 24, 26, 28, 30, 32, 34});
128          size_t idx = 0;
129          for (auto& b : _g[en+"sigred"]->bins()) {
130            book(b, 5, 1, ++idx);
131          }
132          // CC e- p at sqrts=318
133          book(_g[en+"sigred_cc"], {280., 325., 460., 545., 900., 1120., 1300., 1755., 1800., 2270.,
134                                   2500., 3685., 4000., 6520., 7000., 9275., 10000., 20000., 42000.});
135          _g[en+"sigred_cc"]->maskBins({2, 4, 6, 8, 10, 12, 14, 16});
136          idx = 0;
137          for (auto& b : _g[en+"sigred_cc"]->bins()) {
138            book(b, 7, 1, ++idx);
139          }
140        } // end of if
141      } // end of for
142      if (_sqs == "" && !merging()) {
143        throw BeamError("Invalid beam energy for " + name() + "\n");
144      }
145    } // end of init
146
147
148    void analyze(const Event& event) {
149
150      /// @todo Do the event by event analysis here
151      const DISKinematics& dk = apply<DISKinematics>(event, "Kinematics");
152      const DISLepton& dl = apply<DISLepton>(event,"Lepton");
153
154      // Get the DIS kinematics
155      double x  = dk.x();
156      double y = dk.y();
157      double Q2 = dk.Q2()/GeV;
158
159      // Flux factor
160      const double alpha = 7.29927e-3;
161      // GF = 1.16638e-5 Fermi constant
162      const double GF2 = 1.16638e-5*1.16638e-5;
163      // MW = 80.385 W-boson mass
164      const double MW2 = 80.385 * 80.385;
165
166      if (PID::isNeutrino(dl.out().abspid()) ) {
167        // fill histo for CC
168        double F = 2.0*M_PI*x/GF2 * sqr((MW2 + Q2)/MW2);
169        _g[_sqs+"sigred_cc"]->fill(Q2,x,F); // fill histogram x,Q2
170      }
171      else {
172        // fill histo for NC
173        double F = x*sqr(Q2)/(2.0*M_PI*sqr(alpha)*(1.0+sqr(1-y)));
174        _g[_sqs+"sigred"]->fill(Q2,x,F); // fill histogram x,Q2
175      }
176    }
177
178
179    /// Normalise histograms etc., after the run
180    void finalize() {
181      const double gev2nb = 0.389e6;
182      const double scalefactor=crossSection()/nanobarn/sumOfWeights()/gev2nb ;
183      // with _h_sigred.scale also q2 bin width is scaled
184      scale(_g, scalefactor);
185      divByGroupWidth(_g);
186    }
187
188    /// @}
189
190
191    /// @name Histograms
192    /// @{
193    map<string,Histo1DGroupPtr> _g;
194
195    string _sqs = "";
196
197    bool _positron;
198    /// @}
199
200  };
201
202
203  RIVET_DECLARE_PLUGIN(HERA_2015_I1377206);
204
205}