Rivet analyses referenceATLAS_2014_I13266413-jet cross section with 7 TeV dataExperiment: ATLAS (LHC) Inspire ID: 1326641 Status: VALIDATED Authors:
Beam energies: (3500.0, 3500.0) GeV Run details:
Double-differential three-jet production cross-sections have been measured in proton-proton collisions at a centre-of-mass energy of $\sqrt{s} = 7\,$TeV using the ATLAS detector at the Large Hadron Collider. The measurements are presented as a function of the three-jet mass ($m_{jjj}$), in bins of the sum of the absolute rapidity separations between the three leading jets ($\left|Y^\ast\right|$). Invariant masses extending up to 5\,TeV are reached for $8<\left|Y^\ast\right|<10$. These measurements use a sample of data recorded using the ATLAS detector in 2011, which corresponds to an integrated luminosity of $4.51\,\text{fb}^{-1}$. Jets are identified using the anti-$k_\text{t}$ algorithm with two different jet radius parameters, $R=0.4$ and $R=0.6$. Source code: ATLAS_2014_I1326641.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 ATLAS_2014_I1326641 : public Analysis {
10 public:
11
12 /// @name Constructors etc.
13 /// @{
14
15 /// Constructor
16 RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2014_I1326641);
17
18 /// @}
19
20
21 /// @name Analysis methods
22 /// @{
23
24 /// Book histograms and initialise projections before the run
25 void init() {
26
27 const FinalState fs;
28
29 FastJets fj04(fs, JetAlg::ANTIKT, 0.4);
30 fj04.useInvisibles();
31 declare(fj04, "AntiKT04");
32
33 FastJets fj06(fs, JetAlg::ANTIKT, 0.6);
34 fj06.useInvisibles();
35 declare(fj06, "AntiKT06");
36
37 const vector<double> ystarBins{ 0.0, 2.0, 4.0, 6.0, 8.0, 10.0 };
38
39 size_t massDsOffset(0);
40 for (size_t alg = 0; alg < 2; ++alg) {
41 book(h_trijet_Mass[alg], ystarBins);
42 for (auto& b : h_trijet_Mass[alg]->bins()) {
43 book(b, 1+massDsOffset, 1, 1);
44 ++massDsOffset;
45 }
46 }
47 }
48
49 /// Perform the per-event analysis
50 void analyze(const Event& event) {
51
52 Jets jetAr[2];
53 jetAr[AKT4] = apply<FastJets>(event, "AntiKT04").jetsByPt(Cuts::pT > 50*GeV);
54 jetAr[AKT6] = apply<FastJets>(event, "AntiKT06").jetsByPt(Cuts::pT > 50*GeV);
55
56 const size_t nJets = 3;
57 double ptCut[nJets] = { 150., 100., 50.};
58
59 // Loop over jet "radii" used in analysis
60 for (size_t alg = 0; alg < 2; ++alg) {
61 // Identify 3jets
62 vector<FourMomentum> leadJets;
63 for (const Jet& jet : jetAr[alg]) {
64 if (jet.absrap() < 3.0 && leadJets.size() < nJets) {
65 int filledJets = leadJets.size();
66 if (jet.pT() < ptCut[filledJets]) continue;
67 leadJets.push_back(jet.momentum());
68 }
69 }
70
71 if (leadJets.size() < nJets) {
72 MSG_DEBUG("Could not find three suitable leading jets");
73 continue;
74 }
75
76 const double y1 = leadJets[0].rapidity();
77 const double y2 = leadJets[1].rapidity();
78 const double y3 = leadJets[2].rapidity();
79
80 const double yStar = fabs(y1-y2) + fabs(y2-y3) + fabs(y1-y3);
81 const double m = (leadJets[0] + leadJets[1] + leadJets[2]).mass();
82 h_trijet_Mass[alg]->fill(yStar, m);
83 }
84 }
85
86
87 /// Normalise histograms etc., after the run
88 void finalize() {
89
90 const double sf( 2.0 * crossSection()/picobarn / sumOfWeights());
91 scale(h_trijet_Mass, sf);
92
93 }
94
95 /// @}
96
97
98 private:
99
100 // Data members like post-cuts event weight counters go here
101 enum Alg { AKT4=0, AKT6=1 };
102
103 private:
104
105 // The 3 jets mass spectrum for anti-kt 4 and anti-kt 6 jets (array index is jet type from enum above)
106 Histo1DGroupPtr h_trijet_Mass[2];
107
108 };
109
110 RIVET_DECLARE_PLUGIN(ATLAS_2014_I1326641);
111}
|