rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2017_I1604271

ATLAS Inclusive jet cross section measurement at sqrt(s)=8TeV
Experiment: ATLAS (LHC)
Inspire ID: 1604271
Status: VALIDATED
Authors:
  • Jon Butterworth
  • Jonathan Bossio
  • Zdenek Hubacek
References: Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • LHC 8TeV, ATLAS, QCD jet production, anti-kt R=0.4 jets, anti-kt R=0.6 jets

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, JetAlg::ANTIKT, 0.4);
22      FastJets fj06(fs, JetAlg::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 vector<double> ybins{ 0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0};
30
31      // Book histograms
32      // pT histograms
33      book(_pThistograms6, ybins);
34      book(_pThistograms4, ybins);
35      for (size_t i = 0; i < _pThistograms4->numBins(); ++i) {
36        book(_pThistograms6->bin(i+1), i+1, 1, 1);
37        book(_pThistograms4->bin(i+1), i+7, 1, 1);
38      }
39    }
40
41
42    /// Perform the per-event analysis
43    void analyze(const Event& event) {
44
45      const Jets& kt4Jets = apply<FastJets>(event, "AntiKT04").jetsByPt(Cuts::pT>70*GeV && Cuts::absrap < 3.0);
46      const Jets& kt6Jets = apply<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_norm_factor = 0.5*crossSection() / picobarn / sumOfWeights();
74
75      scale(_pThistograms4, xs_norm_factor);
76      scale(_pThistograms6, xs_norm_factor);
77      divByGroupWidth({_pThistograms4, _pThistograms6});
78
79    }
80
81  private:
82
83    // The inclusive pT spectrum for akt4 jets
84    Histo1DGroupPtr _pThistograms4;
85    // The inclusive pT spectrum for akt6 jets
86    Histo1DGroupPtr _pThistograms6;
87
88  };
89
90  RIVET_DECLARE_PLUGIN(ATLAS_2017_I1604271);
91
92}