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
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"

namespace Rivet {


  /// @brief Inclusive jets at 8 TeV
  class ATLAS_2017_I1604271 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2017_I1604271);

    /// Book histograms and initialise projections before the run
    void init() {

      const FinalState fs;
      declare(fs,"FinalState");
      FastJets fj04(fs, FastJets::ANTIKT, 0.4);
      FastJets fj06(fs, FastJets::ANTIKT, 0.6);
      fj04.useInvisibles();
      declare(fj04, "AntiKT04");
      fj06.useInvisibles();
      declare(fj06, "AntiKT06");

      // |y| and ystar bins
      const int nybins = 6;
      double ybins[nybins+1] = { 0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0};

      // Book histograms
      // pT histograms
      for (size_t i = 0; i < nybins; ++i){ // loop over |y| bins
        {Histo1DPtr tmp; _pThistograms6.add(ybins[i], ybins[i+1], book(tmp, i+1,1,1));}
        {Histo1DPtr tmp; _pThistograms4.add(ybins[i], ybins[i+1], book(tmp, i+7,1,1));}
      }

    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {

      const Jets& kt4Jets = applyProjection<FastJets>(event, "AntiKT04").jetsByPt(Cuts::pT>70*GeV && Cuts::absrap < 3.0);
      const Jets& kt6Jets = applyProjection<FastJets>(event, "AntiKT06").jetsByPt(Cuts::pT>70*GeV && Cuts::absrap < 3.0);

      int nJets4 = kt4Jets.size();
      int nJets6 = kt6Jets.size();

      // Inclusive jet selection
      for(int ijet=0;ijet<nJets4;++ijet){ // loop over jets
        FourMomentum jet = kt4Jets[ijet].momentum();
        // pT selection
        const double absy = jet.absrap();
        _pThistograms4.fill(absy,jet.pt()/GeV);
      }

      for(int ijet=0;ijet<nJets6;++ijet){ // loop over jets
        FourMomentum jet = kt6Jets[ijet].momentum();
        // pT selection
        const double absy = jet.absrap();
        _pThistograms6.fill(absy,jet.pt()/GeV);
      }


    }


    /// Normalise histograms etc., after the run
    void finalize() {

      const double xs_pb( crossSection() / picobarn );
      const double sumW( sumOfWeights() );
      const double xs_norm_factor( 0.5*xs_pb / sumW );
     
      MSG_DEBUG( "Cross-Section/pb     : " << xs_pb       );
      MSG_DEBUG( "ZH                   : " << crossSectionPerEvent()/ picobarn);
      MSG_DEBUG( "Sum of weights       : " << sumW        );
      MSG_DEBUG( "nEvents              : " << numEvents() );
      _pThistograms4.scale(xs_norm_factor, this);
      _pThistograms6.scale(xs_norm_factor, this);

    }

  private:

    // The inclusive pT spectrum for akt4 jets
    BinnedHistogram _pThistograms4;
    // The inclusive pT spectrum for akt6 jets
    BinnedHistogram _pThistograms6;

  };

  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2017_I1604271);

}