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 declare(DISLepton(), "Lepton");
28 const DISFinalState disfs(DISFrame::LAB);
29 FastJets jets(disfs, JetAlg::KT, 1.0);
30 declare(jets, "Jets");
31
32 // Table 11
33 book(_h_eta_incl[0], 11, 1, 1);
34 book(_h_eta_incl[1], 11, 1, 2);
35 // Table 12
36 book(_h_eta_di[0], 12, 1, 1);
37 book(_h_eta_di[1], 12, 1, 2);
38 // Table 13
39 book(_h_eta_tri[0], 13, 1, 1);
40 book(_h_eta_tri[1], 13, 1, 2);
41 // Table 14
42 book(_h_et_incl[0], 14, 1, 1);
43 book(_h_et_incl[1], 14, 1, 2);
44 // Table 15
45 book(_h_et_di[0], 15, 1, 1);
46 book(_h_et_di[1], 15, 1, 2);
47 // Table 16
48 book(_h_et_tri[0], 16, 1, 1);
49 book(_h_et_tri[1], 16, 1, 2);
50 // Table 17
51 book(_h_q2_incl[0], 17, 1, 1);
52 book(_h_q2_incl[1], 17, 1, 2);
53 // Table 18
54 book(_h_q2_di[0], 18, 1, 1);
55 book(_h_q2_di[1], 18, 1, 2);
56 // Table 19
57 book(_h_q2_tri[0], 19, 1, 1);
58 book(_h_q2_tri[1], 19, 1, 2);
59
60 // Table 20
61 book(_h_x_incl[0], 20, 1, 1);
62 book(_h_x_incl[1], 20, 1, 2);
63 // Table 22
64 book(_h_m_di[0], 22, 1, 1);
65 book(_h_m_di[1], 22, 1, 2);
66 // Table 23
67 book(_h_m_tri[0], 23, 1, 1);
68 book(_h_m_tri[1], 23, 1, 2);
69
70 }
71
72
73 /// Perform the per-event analysis
74 void analyze(const Event& event) {
75
76 int fLepton = 0;
77
78 const ParticlePair bs = event.beams();
79 if (bs.first.pid() == PID::POSITRON || bs.second.pid() == PID::POSITRON) fLepton = 1;
80 const Particle& bproton = (bs.first.pid() == PID::PROTON) ? bs.first : bs.second;
81 const int orientation = sign(bproton.momentum().pz());
82 // DIS kinematics
83 const DISKinematics& dk = apply<DISKinematics>(event, "Kinematics");
84 double q2 = dk.Q2();
85 double x = dk.x();
86 double y = dk.y();
87
88 if (q2 < 200) vetoEvent;
89 if (y > 0.9) vetoEvent;
90 // make sure charged current
91 const DISLepton& dl = apply<DISLepton>(event,"Lepton");
92 if ( ! PID::isNeutrino(dl.out().abspid()) ) vetoEvent;
93 // Jet selection
94 const Jets jets = apply<FastJets>(event, "Jets").jets(Cuts::Et > 5*GeV && Cuts::etaIn(-1*orientation, 2.5*orientation), cmpMomByEt);
95 MSG_DEBUG("Jet multiplicity = " << jets.size());
96 if (jets.size() < 1) vetoEvent;
97 if (jets[0].Et() < 14*GeV) vetoEvent;
98
99 double eta12 = 0;
100 double et12 = 0;
101
102 double et123 = 0;
103 double eta123 =0;
104
105 for (size_t i = 0; i < jets.size(); ++i) {
106 if (jets[i].Et() < 14*GeV) continue;
107 _h_eta_incl[fLepton]->fill(orientation*jets[i].eta());
108 _h_et_incl[fLepton]->fill(jets[i].Et());
109 _h_q2_incl[fLepton]->fill(q2);
110 _h_x_incl[fLepton]->fill(x);
111 }
112
113 if (jets.size() > 1) {
114 eta12 = orientation*(jets[0].eta() + jets[1].eta())/2;
115 et12 = (jets[0].Et() + jets[1].Et())/2;
116 _h_eta_di[fLepton]->fill(eta12);
117 _h_et_di[fLepton]->fill(et12);
118 _h_q2_di[fLepton]->fill(q2);
119 _h_m_di[fLepton]->fill( (jets[0].momentum()+jets[1].momentum()).mass());
120 }
121
122 if (jets.size() > 2) {
123 eta123 = orientation*(jets[0].eta() + jets[1].eta()+jets[2].eta())/3;
124 et123 = (jets[0].Et() + jets[1].Et()+jets[2].Et())/3;
125
126 _h_eta_tri[fLepton]->fill(eta123);
127 _h_et_tri[fLepton]->fill(et123);
128 _h_q2_tri[fLepton]->fill(q2);
129 _h_m_tri[fLepton]->fill( (jets[0].momentum()+jets[1].momentum()+jets[2].momentum()).mass());
130 }
131 }
132
133
134 /// Normalise histograms etc., after the run
135 void finalize() {
136 const double sf = crossSection()/picobarn/sumOfWeights();
137
138
139 scale(_h_eta_incl, sf);
140 scale(_h_eta_di, sf);
141 scale(_h_eta_tri, sf);
142 scale(_h_et_incl, sf);
143 scale(_h_et_di, sf);
144 scale(_h_et_tri, sf);
145 scale(_h_q2_incl, sf);
146 scale(_h_q2_di, sf);
147 scale(_h_q2_tri, sf);
148 scale(_h_x_incl, sf);
149 scale(_h_m_di, sf);
150 scale(_h_m_tri, sf);
151
152 }
153
154 /// @}
155
156
157 /// @name Histograms
158 /// @{
159 Histo1DPtr _h_eta_incl[2], _h_eta_di[2], _h_eta_tri[2];
160 Histo1DPtr _h_et_incl[2], _h_et_di[2], _h_et_tri[2];
161 Histo1DPtr _h_q2_incl[2], _h_q2_di[2], _h_q2_tri[2];
162 Histo1DPtr _h_x_incl[2], _h_x_di[2], _h_x_tri[2];
163 Histo1DPtr _h_m_di[2], _h_m_tri[2];
164 /// @}
165 };
166
167
168 RIVET_DECLARE_PLUGIN(ZEUS_2008_I780108);
169
170}
|