Rivet analyses referenceZEUS_2008_I780108Multi-jet cross-sections in charged current $e^{\pm} p$ scattering at HERAExperiment: ZEUS (HERA) Inspire ID: 780108 Status: VALIDATED Authors:
Beam energies: ANY Run details:
Jet cross sections were measured in charged current deep inelastic $e^{\pm}p$ scattering at high boson virtualities $Q^2$ with the ZEUS detector at HERA II using an integrated luminosity of $0.36 fb^{-1}$. Differential cross sections are presented for inclusive-jet production as functions of $Q^2$, Bjorken x and the jet transverse energy and pseudorapidity. The dijet invariant mass cross section is also presented. Observation of three- and four-jet events in charged-current $e^{\pm}p$ processes is reported for the first time. The predictions of next-to-leading-order (NLO) QCD calculations are compared to the measurements. The measured inclusive-jet cross sections are well described in shape and normalization by the NLO predictions. The data have the potential to constrain the u and d valence quark distributions in the proton if included as input to global fits. Source code: ZEUS_2008_I780108.cc 1// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/FinalState.hh"
4#include "Rivet/Projections/FastJets.hh"
5#include "Rivet/Projections/DISFinalState.hh"
6#include "Rivet/Projections/DISKinematics.hh"
7
8namespace Rivet {
9
10
11 /// @brief Multi-jet cross-sections in charged current $e^{\pm} p$ scattering at HERA
12 class ZEUS_2008_I780108 : public Analysis {
13 public:
14
15 /// Constructor
16 RIVET_DEFAULT_ANALYSIS_CTOR(ZEUS_2008_I780108);
17
18 /// @name Analysis methods
19 /// @{
20
21 /// Book histograms and initialise projections before the run
22 void init() {
23
24 // Projections
25 const DISKinematics diskin;
26 declare(diskin, "Kinematics");
27 const DISFinalState disfs(DISFinalState::BoostFrame::LAB);
28 FastJets jets(disfs, FastJets::KT, 1.0);
29 declare(jets, "Jets");
30
31 // Table 11
32 book(_h_eta_incl[0], 11, 1, 1);
33 book(_h_eta_incl[1], 11, 1, 2);
34 // Table 12
35 book(_h_eta_di[0], 12, 1, 1);
36 book(_h_eta_di[1], 12, 1, 2);
37 // Table 13
38 book(_h_eta_tri[0], 13, 1, 1);
39 book(_h_eta_tri[1], 13, 1, 2);
40 // Table 14
41 book(_h_et_incl[0], 14, 1, 1);
42 book(_h_et_incl[1], 14, 1, 2);
43 // Table 15
44 book(_h_et_di[0], 15, 1, 1);
45 book(_h_et_di[1], 15, 1, 2);
46 // Table 16
47 book(_h_et_tri[0], 16, 1, 1);
48 book(_h_et_tri[1], 16, 1, 2);
49 // Table 17
50 book(_h_q2_incl[0], 17, 1, 1);
51 book(_h_q2_incl[1], 17, 1, 2);
52 // Table 18
53 book(_h_q2_di[0], 18, 1, 1);
54 book(_h_q2_di[1], 18, 1, 2);
55 // Table 19
56 book(_h_q2_tri[0], 19, 1, 1);
57 book(_h_q2_tri[1], 19, 1, 2);
58
59 // Table 20
60 book(_h_x_incl[0], 20, 1, 1);
61 book(_h_x_incl[1], 20, 1, 2);
62 // Table 22
63 book(_h_m_di[0], 22, 1, 1);
64 book(_h_m_di[1], 22, 1, 2);
65 // Table 23
66 book(_h_m_tri[0], 23, 1, 1);
67 book(_h_m_tri[1], 23, 1, 2);
68
69 }
70
71
72 /// Perform the per-event analysis
73 void analyze(const Event& event) {
74
75 int fLepton = 0;
76
77 const ParticlePair bs = event.beams();
78 if (bs.first.pid() == PID::POSITRON || bs.second.pid() == PID::POSITRON) fLepton = 1;
79 const Particle& bproton = (bs.first.pid() == PID::PROTON) ? bs.first : bs.second;
80 const int orientation = sign(bproton.momentum().pz());
81 // DIS kinematics
82 const DISKinematics& dk = apply<DISKinematics>(event, "Kinematics");
83 double q2 = dk.Q2();
84 double x = dk.x();
85 double y = dk.y();
86
87 if (q2 < 200) vetoEvent;
88 if (y > 0.9) vetoEvent;
89 // Jet selection
90 const Jets jets = apply<FastJets>(event, "Jets").jets(Cuts::Et > 5*GeV && Cuts::etaIn(-1*orientation, 2.5*orientation), cmpMomByEt);
91 MSG_DEBUG("Jet multiplicity = " << jets.size());
92 if (jets.size() < 1) vetoEvent;
93 if (jets[0].Et() < 14*GeV) vetoEvent;
94
95 double eta12 = 0;
96 double et12 = 0;
97
98 double et123 = 0;
99 double eta123 =0;
100
101 for (size_t i = 0; i < jets.size(); i++)
102 {
103 if (jets[i].Et() < 14*GeV) continue;
104 _h_eta_incl[fLepton]->fill(orientation*jets[i].eta());
105 _h_et_incl[fLepton]->fill(jets[i].Et());
106 _h_q2_incl[fLepton]->fill(q2);
107 _h_x_incl[fLepton]->fill(x);
108 }
109
110 if (jets.size() > 1)
111 {
112 eta12 = orientation*(jets[0].eta() + jets[1].eta())/2;
113 et12 = (jets[0].Et() + jets[1].Et())/2;
114 _h_eta_di[fLepton]->fill(eta12);
115 _h_et_di[fLepton]->fill(et12);
116 _h_q2_di[fLepton]->fill(q2);
117 _h_m_di[fLepton]->fill( (jets[0].momentum()+jets[1].momentum()).mass());
118 }
119
120 if (jets.size() > 2)
121 {
122 eta123 = orientation*(jets[0].eta() + jets[1].eta()+jets[2].eta())/3;
123 et123 = (jets[0].Et() + jets[1].Et()+jets[2].Et())/3;
124
125 _h_eta_tri[fLepton]->fill(eta123);
126 _h_et_tri[fLepton]->fill(et123);
127 _h_q2_tri[fLepton]->fill(q2);
128 _h_m_tri[fLepton]->fill( (jets[0].momentum()+jets[1].momentum()+jets[2].momentum()).mass());
129 }
130 }
131
132
133 /// Normalise histograms etc., after the run
134 void finalize() {
135 const double sf = crossSection()/picobarn/sumOfWeights();
136
137
138 scale(_h_eta_incl, sf);
139 scale(_h_eta_di, sf);
140 scale(_h_eta_tri, sf);
141 scale(_h_et_incl, sf);
142 scale(_h_et_di, sf);
143 scale(_h_et_tri, sf);
144 scale(_h_q2_incl, sf);
145 scale(_h_q2_di, sf);
146 scale(_h_q2_tri, sf);
147 scale(_h_x_incl, sf);
148 scale(_h_m_di, sf);
149 scale(_h_m_tri, sf);
150
151 }
152
153 /// @}
154
155
156 /// @name Histograms
157 /// @{
158 Histo1DPtr _h_eta_incl[2], _h_eta_di[2], _h_eta_tri[2];
159 Histo1DPtr _h_et_incl[2], _h_et_di[2], _h_et_tri[2];
160 Histo1DPtr _h_q2_incl[2], _h_q2_di[2], _h_q2_tri[2];
161 Histo1DPtr _h_x_incl[2], _h_x_di[2], _h_x_tri[2];
162 Histo1DPtr _h_m_di[2], _h_m_tri[2];
163 /// @}
164 };
165
166
167 RIVET_DECLARE_PLUGIN(ZEUS_2008_I780108);
168
169}
|