Rivet Analyses Reference

ATLAS_2014_I1325553

Measurement of the inclusive jet cross-section at 7 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1325553
Status: VALIDATED
Authors:
  • Vojtech Pleskot
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • p p -> jet + X. $\sqrt{s} = 7$ TeV.

Measurement of the inclusive jet cross-section in proton--proton collisions at a centre-of-mass energy of 7 TeV using a data set corresponding to an integrated luminosity of 4.5/fb collected with the ATLAS detector at the Large Hadron Collider in 2011. Jets are identified using the anti-$k_t$ algorithm with radius parameter values of 0.4 and 0.6. The double-differential cross-sections are presented as a function of the jet transverse momentum and the jet rapidity, covering jet transverse momenta from 100 GeV to 2 TeV.

Source code: ATLAS_2014_I1325553.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
 99
100
101
102
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Tools/BinnedHistogram.hh"

namespace Rivet {


  /// Jet mass as a function of ystar
  class ATLAS_2014_I1325553 : public Analysis {
  public:

    /// Constructor
    ATLAS_2014_I1325553()
      : Analysis("ATLAS_2014_I1325553")
    {    }


    /// @name Analysis methods
    //@{

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

      const FinalState fs;
      declare(fs,"FinalState");

      FastJets fj04(fs,  FastJets::ANTIKT, 0.4);
      fj04.useInvisibles();
      declare(fj04, "AntiKT04");

      FastJets fj06(fs,  FastJets::ANTIKT, 0.6);
      fj06.useInvisibles();
      declare(fj06, "AntiKT06");

      double ybins[] = {0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0};

      size_t ptDsOffset(0);
      for (size_t alg = 0; alg < 2; ++alg) {
        for (size_t i = 0; i < 6; ++i) {
          Histo1DPtr tmp;
          _pt[alg].add(ybins[i], ybins[i + 1], book(tmp, 1 + ptDsOffset, 1, 1));
          ptDsOffset += 1;
        }
      }
    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {
      Jets jetAr[2];
      jetAr[AKT4] = apply<FastJets>(event, "AntiKT04").jetsByPt(Cuts::pT > 100*GeV && Cuts::absrap < 3.0);
      jetAr[AKT6] = apply<FastJets>(event, "AntiKT06").jetsByPt(Cuts::pT > 100*GeV && Cuts::absrap < 3.0);

      // Loop over jet "radii" used in analysis
      for (size_t alg = 0; alg < 2; ++alg) {

        // fill the 1D pt histograms with all the jets passing the cuts
        for (const Jet& jet : jetAr[alg]) {
          const double absrap = jet.absrap();
          if (absrap < 3.0) {
            const double pt = jet.pT();
            if (pt/GeV > 100*GeV) {
              _pt[alg].fill(absrap, pt/GeV);
            }
          }
        }
      }
    }


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

      /// Print summary info
      const double xs_pb( crossSection() / picobarn );
      const double sumW( sumOfWeights() );
      const double xs_norm_factor( 0.5*xs_pb / sumW );

      for (size_t alg = 0; alg < 2; ++alg) {
        _pt[alg].scale(xs_norm_factor, this);
      }
    }

    //@}


  private:

    // Data members like post-cuts event weight counters go here
    enum Alg { AKT4=0, AKT6=1 };

    /// The inclusive jet spectrum binned in rapidity for akt6 and akt4 jets (array index is jet type from enum above)
    BinnedHistogram _pt[2];

  };


  RIVET_DECLARE_PLUGIN(ATLAS_2014_I1325553);

}