Rivet analyses referenceATLAS_2017_I1604271ATLAS Inclusive jet cross section measurement at sqrt(s)=8TeVExperiment: ATLAS (LHC) Inspire ID: 1604271 Status: VALIDATED Authors:
Beam energies: (4000.0, 4000.0) GeV Run details:
Inclusive jet production cross-sections are measured in proton-proton collisions at a centre-of-mass energy of $\sqrt{s}=8$ TeV recorded by the ATLAS experiment at the Large Hadron Collider at CERN. The total integrated luminosity of the analysed data set amounts to 20.2 fb$^{-1}$. Double-differential cross-sections are measured for jets defined by the anti-$k_\text{t}$ jet clustering algorithm with radius parameters of $R=0.4$ and $R=0.6$ and are presented as a function of the jet transverse momentum, in the range between 70 GeV and 2.5 TeV and in six bins of the absolute jet rapidity, between 0 and 3.0. The measured cross-sections are compared to predictions of quantum chromodynamics, calculated at next-to-leading order in perturbation theory, and corrected for non-perturbative and electroweak effects. The level of agreement with predictions, using a selection of different parton distribution functions for the proton, is quantified. Tensions between the data and the theory predictions are observed. Source code: ATLAS_2017_I1604271.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 /// @brief Inclusive jets at 8 TeV
10 class ATLAS_2017_I1604271 : public Analysis {
11 public:
12
13 /// Constructor
14 RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2017_I1604271);
15
16 /// Book histograms and initialise projections before the run
17 void init() {
18
19 const FinalState fs;
20 declare(fs,"FinalState");
21 FastJets fj04(fs, FastJets::ANTIKT, 0.4);
22 FastJets fj06(fs, FastJets::ANTIKT, 0.6);
23 fj04.useInvisibles();
24 declare(fj04, "AntiKT04");
25 fj06.useInvisibles();
26 declare(fj06, "AntiKT06");
27
28 // |y| and ystar bins
29 const int nybins = 6;
30 double ybins[nybins+1] = { 0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0};
31
32 // Book histograms
33 // pT histograms
34 for (size_t i = 0; i < nybins; ++i){ // loop over |y| bins
35 {Histo1DPtr tmp; _pThistograms6.add(ybins[i], ybins[i+1], book(tmp, i+1,1,1));}
36 {Histo1DPtr tmp; _pThistograms4.add(ybins[i], ybins[i+1], book(tmp, i+7,1,1));}
37 }
38
39 }
40
41
42 /// Perform the per-event analysis
43 void analyze(const Event& event) {
44
45 const Jets& kt4Jets = applyProjection<FastJets>(event, "AntiKT04").jetsByPt(Cuts::pT>70*GeV && Cuts::absrap < 3.0);
46 const Jets& kt6Jets = applyProjection<FastJets>(event, "AntiKT06").jetsByPt(Cuts::pT>70*GeV && Cuts::absrap < 3.0);
47
48 int nJets4 = kt4Jets.size();
49 int nJets6 = kt6Jets.size();
50
51 // Inclusive jet selection
52 for(int ijet=0;ijet<nJets4;++ijet){ // loop over jets
53 FourMomentum jet = kt4Jets[ijet].momentum();
54 // pT selection
55 const double absy = jet.absrap();
56 _pThistograms4.fill(absy,jet.pt()/GeV);
57 }
58
59 for(int ijet=0;ijet<nJets6;++ijet){ // loop over jets
60 FourMomentum jet = kt6Jets[ijet].momentum();
61 // pT selection
62 const double absy = jet.absrap();
63 _pThistograms6.fill(absy,jet.pt()/GeV);
64 }
65
66
67 }
68
69
70 /// Normalise histograms etc., after the run
71 void finalize() {
72
73 const double xs_pb( crossSection() / picobarn );
74 const double sumW( sumOfWeights() );
75 const double xs_norm_factor( 0.5*xs_pb / sumW );
76
77 MSG_DEBUG( "Cross-Section/pb : " << xs_pb );
78 MSG_DEBUG( "ZH : " << crossSectionPerEvent()/ picobarn);
79 MSG_DEBUG( "Sum of weights : " << sumW );
80 MSG_DEBUG( "nEvents : " << numEvents() );
81 _pThistograms4.scale(xs_norm_factor, this);
82 _pThistograms6.scale(xs_norm_factor, this);
83
84 }
85
86 private:
87
88 // The inclusive pT spectrum for akt4 jets
89 BinnedHistogram _pThistograms4;
90 // The inclusive pT spectrum for akt6 jets
91 BinnedHistogram _pThistograms6;
92
93 };
94
95 // The hook for the plugin system
96 RIVET_DECLARE_PLUGIN(ATLAS_2017_I1604271);
97
98}
|