Rivet analyses referenceCDF_2008_I768579Jet pT and multiplicity distributions in W + jets eventsExperiment: CDF (Tevatron Run 2) Inspire ID: 768579 Status: VALIDATED Authors:
Beam energies: (980.0, 980.0) GeV Run details:
Measurement of the cross section for W boson production in association with jets in $p\bar{p}$ collisions at $\sqrt{s}=1.96$ TeV. The analysis uses 320 pb$^{-1}$ of data collected with the CDF II detector. W bosons are identified in their $e\nu$ decay channel and jets are reconstructed using an $R < 0.4$ cone algorithm. For each $W + \geq$ n-jet sample (where n = 1--4) a measurement of d$\sigma({p}\bar{p} \rightarrow W + \geq$ n jet)/d$E_T(n^{th}$-jet) $\times$ BR($W \rightarrow{e}\nu$) is made, where d$E_T(n^{th}$-jet) is the Et of the n$^{th}$-highest Et jet above 20 GeV. A measurement of the total cross section, $\sigma(p\bar{p} \rightarrow W + \geq$ $n$-jet) $\times$ BR($W \rightarrow{e}\nu)$ with $E_T(n^{th}-jet) > 25$ GeV is also made. Both measurements are made for jets with $|\eta| < 2$ and for a limited region of the $W \rightarrow{e}\nu$ decay phase space; $|\eta^{e}| < 1.1$, $p_{T}^{e} > 20$ GeV, $p_{T}^{\nu} > 30$ GeV and $M_{T} > 20$ GeV. The cross sections are corrected for all detector effects and can be directly compared to particle level $W$ + jet(s) predictions. These measurements can be used to test and tune QCD predictions for the number of jets in and kinematics of $W$ + jets events. Source code: CDF_2008_I768579.cc 1// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/FinalState.hh"
4#include "Rivet/Projections/VetoedFinalState.hh"
5#include "Rivet/Projections/InvMassFinalState.hh"
6#include "Rivet/Projections/FastJets.hh"
7#include <algorithm>
8
9namespace Rivet {
10
11
12 /// @brief CDF jet pT and multiplicity distributions in W + jets events
13 ///
14 /// This CDF analysis provides jet pT distributions for 4 jet multiplicity bins
15 /// as well as the jet multiplicity distribution in W + jets events.
16 class CDF_2008_I768579 : public Analysis {
17 public:
18
19 RIVET_DEFAULT_ANALYSIS_CTOR(CDF_2008_I768579);
20
21
22 /// @name Analysis methods
23 /// @{
24
25 void init() {
26 // Set up projections
27 // Basic FS
28 FinalState fs((Cuts::etaIn(-3.6, 3.6)));
29 declare(fs, "FS");
30
31 // Create a final state with any e-nu pair with invariant mass 65 -> 95 GeV and ET > 20 (W decay products)
32 vector<pair<PdgId,PdgId> > vids;
33 vids += make_pair(PID::ELECTRON, PID::NU_EBAR);
34 vids += make_pair(PID::POSITRON, PID::NU_E);
35 FinalState fs2(Cuts::abseta < 3.6 && Cuts::pT >= 20*GeV);
36 InvMassFinalState invfs(fs2, vids, 65*GeV, 95*GeV);
37 declare(invfs, "INVFS");
38
39 // Make a final state without the W decay products for jet clustering
40 VetoedFinalState vfs(fs);
41 vfs.addVetoOnThisFinalState(invfs);
42 declare(vfs, "VFS");
43 declare(FastJets(vfs, JetAlg::CDFJETCLU, 0.4), "Jets");
44
45 // Book histograms
46 for (int i = 0 ; i < 4 ; ++i) {
47 book(_histJetEt[i], i+1, 1, 1);
48 book(_histJetMultRatio[i], 5, 1, i+1);
49 book(_histJetMult[i], i+6, 1, 1);
50 book(_numer[i],"/TMP/numer"+std::to_string(i));
51 book(_denom[i],"/TMP/denom"+std::to_string(i));
52 }
53
54 }
55
56
57 /// Do the analysis
58 void analyze(const Event& event) {
59 // Get the W decay products (electron and neutrino)
60 const InvMassFinalState& invMassFinalState = apply<InvMassFinalState>(event, "INVFS");
61 const Particles& wDecayProducts = invMassFinalState.particles();
62
63 FourMomentum electronP, neutrinoP;
64 bool gotElectron(false), gotNeutrino(false);
65 for (const Particle& p : wDecayProducts) {
66 FourMomentum p4 = p.momentum();
67 if (p4.Et() > 20*GeV && p4.abseta() < 1.1 && p.abspid() == PID::ELECTRON) {
68 electronP = p4;
69 gotElectron = true;
70 }
71 else if (p4.Et() > 30*GeV && p.abspid() == PID::NU_E) {
72 neutrinoP = p4;
73 gotNeutrino = true;
74 }
75 }
76
77 // Veto event if the electron or MET cuts fail
78 if (!gotElectron || !gotNeutrino) vetoEvent;
79
80 // Veto event if the MTR cut fails
81 double mT2 = 2.0 * ( electronP.pT()*neutrinoP.pT() - electronP.px()*neutrinoP.px() - electronP.py()*neutrinoP.py() );
82 if (sqrt(mT2) < 20*GeV ) vetoEvent;
83
84 // Get the jets
85 const JetFinder& jetProj = apply<FastJets>(event, "Jets");
86 Jets theJets = jetProj.jets(cmpMomByEt, Cuts::Et > 20*GeV);
87 size_t njetsA(0), njetsB(0);
88 for (const Jet& j : theJets) {
89 if (j.absrap() < 2.0) {
90 // Fill differential histograms for top 4 jets with Et > 20
91 if (njetsA < 4 && j.Et() > 20*GeV) {
92 ++njetsA;
93 _histJetEt[njetsA-1]->fill(j.Et()/GeV);
94 }
95 // Count number of jets with Et > 25 (for multiplicity histograms)
96 if (j.Et() > 25*GeV) ++njetsB;
97 }
98 }
99
100 // Increment event counter
101 _denom[0]->fill();
102
103 // Jet multiplicity
104 for (size_t i = 1; i <= njetsB; ++i) {
105 _histJetMult[i-1]->fill(1960);
106 _numer[i-1]->fill();
107 if (i == 4) break;
108 _denom[i]->fill();
109 }
110 }
111
112
113 /// Finalize
114 void finalize() {
115
116 // Loop over the non-zero multiplicities
117 for (size_t i = 0; i < 4; ++i) {
118 if (!_denom[i]->val()) continue;
119 const YODA::Estimate0D ratio = (*_numer[i]) / (*_denom[i]);
120 _histJetMultRatio[i]->bin(1).set(ratio.val(), ratio.err());
121 }
122
123 // Normalize the non-ratio histograms
124 for (size_t i = 0; i < 4; ++i) {
125 scale(_histJetEt[i], crossSection()/picobarn/sumOfWeights());
126 scale(_histJetMult[i], crossSection()/picobarn/sumOfWeights());
127 }
128
129 }
130
131 /// @}
132
133
134 private:
135
136 /// @name Histograms
137 /// @{
138 Histo1DPtr _histJetEt[4];
139 BinnedEstimatePtr<int> _histJetMultRatio[4];
140 BinnedHistoPtr<int> _histJetMult[4];
141 CounterPtr _numer[4], _denom[4];
142 /// @}
143
144 };
145
146
147
148 RIVET_DECLARE_ALIASED_PLUGIN(CDF_2008_I768579, CDF_2008_S7541902);
149
150}
|