Rivet analyses referenceD0_2011_I8956623-jet invariant massExperiment: D0 (Tevatron Run 2) Inspire ID: 895662 Status: VALIDATED Authors:
Beam energies: (980.0, 980.0) GeV Run details:
Inclusive three-jet differential cross-section as a function of invariant mass of the three jets with the largest transverse momenta. The measurement is made in three rapidity regions ($|y|<0.8, 1.6, 2.4$) and with jets above 40, 70, and 100 GeV. Source code: D0_2011_I895662.cc 1// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/FinalState.hh"
4#include "Rivet/Projections/FastJets.hh"
5
6namespace Rivet {
7
8
9 class D0_2011_I895662 : public Analysis {
10 public:
11
12 RIVET_DEFAULT_ANALYSIS_CTOR(D0_2011_I895662);
13
14
15 public:
16
17 void init() {
18 FastJets jets(FinalState((Cuts::etaIn(-3.6, 3.6))), JetAlg::D0ILCONE, 0.7);
19 jets.useInvisibles();
20 declare(jets, "Jets");
21
22 book(_h_m3j_08_40 ,1, 1, 1);
23 book(_h_m3j_16_40 ,2, 1, 1);
24 book(_h_m3j_24_40 ,3, 1, 1);
25 book(_h_m3j_24_70 ,4, 1, 1);
26 book(_h_m3j_24_100 ,5, 1, 1);
27 }
28
29
30 void analyze(const Event& event) {
31 Jets jets = apply<FastJets>(event, "Jets").jetsByPt(Cuts::pT > 40*GeV);
32
33 // Need three jets, leading jet above 150 GeV
34 if (jets.size() < 3 || jets[0].pT() <= 150.*GeV) vetoEvent;
35
36 std::vector<FourMomentum> p;
37 for (size_t i=0; i<3; i++) {
38 p.push_back(jets[i].momentum());
39 }
40
41 // Jets need to be separated by 2*Rcone
42 if (deltaR(p[0], p[1], RAPIDITY) < 1.4 ||
43 deltaR(p[0], p[2], RAPIDITY) < 1.4 ||
44 deltaR(p[1], p[2], RAPIDITY) < 1.4)
45 vetoEvent;
46
47 // Leading three jets need to be within |y|<2.4
48 double ymax = fabs(p[0].rapidity());
49 for (size_t i=1; i<3; i++) {
50 if (ymax < fabs(p[i].rapidity())) ymax = fabs(p[i].rapidity());
51 }
52 if (ymax >= 2.4) vetoEvent;
53
54 double m3jet = (p[0]+p[1]+p[2]).mass()/GeV;
55
56 if (ymax < 0.8) _h_m3j_08_40->fill(m3jet);
57 if (ymax < 1.6) _h_m3j_16_40->fill(m3jet);
58 if (ymax < 2.4) {
59 _h_m3j_24_40->fill(m3jet);
60 if (p[2].pT() > 70.*GeV) _h_m3j_24_70->fill(m3jet);
61 if (p[2].pT() > 100.*GeV) _h_m3j_24_100->fill(m3jet);
62 }
63
64 }
65
66
67 void finalize() {
68 // Factor of 1000 is based on GeV <-> TeV mismatch between paper and Hepdata table
69 scale(_h_m3j_08_40, 1000*crossSection()/picobarn/sumOfWeights());
70 scale(_h_m3j_16_40, 1000*crossSection()/picobarn/sumOfWeights());
71 scale(_h_m3j_24_40, 1000*crossSection()/picobarn/sumOfWeights());
72 scale(_h_m3j_24_70, 1000*crossSection()/picobarn/sumOfWeights());
73 scale(_h_m3j_24_100, 1000*crossSection()/picobarn/sumOfWeights());
74 }
75
76
77 private:
78
79 Histo1DPtr _h_m3j_08_40;
80 Histo1DPtr _h_m3j_16_40;
81 Histo1DPtr _h_m3j_24_40;
82 Histo1DPtr _h_m3j_24_70;
83 Histo1DPtr _h_m3j_24_100;
84
85 };
86
87
88 RIVET_DECLARE_PLUGIN(D0_2011_I895662);
89
90}
|