Rivet Analyses Reference

ATLAS_2012_I1082936

Inclusive jet and dijet cross sections at 7 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1082936
Status: VALIDATED
Authors:
  • Holger Schulz
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • QCD jet production with a minimum leading jet pT of 30 GeV and minimum second jet pT of 20 GeV at 7 TeV.

Inclusive jet and dijet cross sections have been measured in proton-proton collisions at a centre-of-mass energy of 7 TeV using the ATLAS detector at the Large Hadron Collider. The cross sections were measured using jets clustered with the anti-kT algorithm with parameters R=0.4 and R=0.6. These measurements are based on the 2010 data sample, consisting of a total integrated luminosity of 37 inverse picobarns. Inclusive jet double-differential cross sections are presented as a function of jet transverse momentum, in bins of jet rapidity. Dijet double-differential cross sections are studied as a function of the dijet invariant mass, in bins of half the rapidity separation of the two leading jets. The measurements are performed in the jet rapidity range $|y|<4.4$, covering jet transverse momenta from 20 GeV to 1.5 TeV and dijet invariant masses from 70 GeV to 5 TeV. This is the successor analysis of ATLAS_2010_S8817804

Source code: ATLAS_2012_I1082936.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
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Tools/BinnedHistogram.hh"
#include "Rivet/Projections/FinalState.hh"

namespace Rivet {


  class ATLAS_2012_I1082936 : public Analysis {
  public:

    /// @name Constructors etc.
    //@{

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

    //@}


  public:

    /// @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");


      // Histogram booking copied from the previous analysis
      double ybins[] = { 0.0, 0.3, 0.8, 1.2, 2.1, 2.8, 3.6, 4.4 };
      double ystarbins[] = { 0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0, 3.5, 4.0, 4.4};

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

        for (size_t i = 0; i < 9; ++i) {
          Histo1DPtr tmp;
          _mass[alg].add(ystarbins[i], ystarbins[i+1], book(tmp, 1 + massDsOffset, 1, i+1));
        }
        massDsOffset += 1;
      }
    }

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

      // Loop over jet "radii" used in analysis
      for (size_t alg = 0; alg < 2; ++alg) {
        // Identify dijets
        vector<FourMomentum> leadjets;
        for (const Jet& jet : jetAr[alg]) {
          const double pT = jet.pT();
          const double absy = jet.absrap();
          _pThistos[alg].fill(absy, pT/GeV);

          if (absy < 4.4 && leadjets.size() < 2) {
            if (leadjets.empty() && pT < 30*GeV) continue;
            leadjets.push_back(jet.momentum());
          }
        }
        // Make sure we have a leading jet with pT >30 GeV and a second to leading jet with pT>20 GeV
        if (leadjets.size() < 2) {
          MSG_DEBUG("Could not find two suitable leading jets");
          continue;
        }

        const double y1 = leadjets[0].rapidity();
        const double y2 = leadjets[1].rapidity();
        const double ystar = fabs(y1-y2)/2.;
        const double m = (leadjets[0] + leadjets[1]).mass();
        // Fill mass histogram
        _mass[alg].fill(ystar, m/TeV);
      }
    }


    /// Normalise histograms etc., after the run
    void finalize() {
      for (size_t alg = 0; alg < 2; ++alg) {
        // factor 0.5 needed because it is differential in dy and not d|y|
        _pThistos[alg].scale(0.5*crossSectionPerEvent()/picobarn, this);
        _mass[alg].scale(crossSectionPerEvent()/picobarn, this);
}
    }

    //@}


  private:

    // Data members like post-cuts event weight counters go here

    enum Alg { AKT4=0, AKT6=1 };

  private:

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

    /// The di-jet mass spectrum binned in rapidity for akt6 and akt4 jets (array index is jet type from enum above)
    BinnedHistogram _mass[2];
  };

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

}