Rivet Analyses Reference

ATLAS_2010_S8817804

Inclusive jet cross section and di-jet mass and chi spectra at 7 TeV in ATLAS
Experiment: ATLAS (LHC 7TeV)
Inspire ID: 871366
Status: VALIDATED
Authors:
  • James Monk
References:Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • pp QCD jet production with a minimum jet pT of 60 GeV (inclusive) or 30 GeV (di-jets) at 7 TeV.

The first jet cross section measurement made with the ATLAS detector at the LHC. Anti-kt jets with $R=0.4$ and $R=0.6$ are resconstructed within $|y|<2.8$ and above 60 GeV for the inclusive jet cross section plots. For the di-jet plots the second jet must have pT>30 GeV. Jet pT and di-jet mass spectra are plotted in bins of rapidity between $|y| = $ 0.3, 0.8, 1.2, 2.1, and 2.8. Di-jet $\chi$ spectra are plotted in bins of di-jet mass between 340 GeV, 520 GeV, 800 GeV and 1200 GeV.

Source code: ATLAS_2010_S8817804.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
// -*- C++ -*-

#include "Rivet/Analysis.hh"
#include "Rivet/Tools/BinnedHistogram.hh"
#include "Rivet/Projections/FastJets.hh"

namespace Rivet {


  /// @brief ATLAS inclusive jet pT spectrum, di-jet Mass and di-jet chi
  class ATLAS_2010_S8817804: public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2010_S8817804);

    /// Choice of jet algorithm
    enum Alg { AKT4=0, AKT6=1 };


    void init() {
      FinalState fs;
      declare(fs, "FinalState");
      declare(FastJets(fs, FastJets::ANTIKT, 0.6), "AntiKT06");
      declare(FastJets(fs, FastJets::ANTIKT, 0.4), "AntiKT04");

      double ybins[] = { 0.0, 0.3, 0.8, 1.2, 2.1, 2.8 };
      double massBinsForChi[] = { 340, 520, 800, 1200 };


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

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

        for (size_t i = 0; i < 3; ++i) {
          Histo1DPtr tmp;
          book(tmp, i + 1 + chiDsOffset, 1, 1);
          _chiVsMass[alg].add(massBinsForChi[i], massBinsForChi[i+1], tmp);
        }
        chiDsOffset += 3;
      }
    }


    void analyze(const Event& evt) {
      Jets jetAr[2];
      jetAr[AKT6] = apply<FastJets>(evt, "AntiKT06").jetsByPt(30*GeV);
      jetAr[AKT4] = apply<FastJets>(evt, "AntiKT04").jetsByPt(30*GeV);

      // Identify the dijets
      for (size_t alg = 0; alg < 2; ++alg) {
        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, 1.0);

          if (absy < 2.8 && leadjets.size() < 2) {
            if (leadjets.empty() && pT < 60*GeV) continue;
            leadjets.push_back(jet.momentum());
          }
        }

        // Veto if no acceptable dijet found
        if (leadjets.size() < 2) {
          MSG_DEBUG("Could not find two suitable leading jets");
          continue;
        }

        const double rap1 = leadjets[0].rapidity();
        const double rap2 = leadjets[1].rapidity();
        const double mass = (leadjets[0] + leadjets[1]).mass();
        const double ymax = max(fabs(rap1), fabs(rap2));
        const double chi = exp(fabs(rap1 - rap2));
        if (fabs(rap1 + rap2) < 2.2) {
          _chiVsMass[alg].fill(mass/GeV, chi, 1.0);
        }
        _massVsY[alg].fill(ymax, mass/GeV, 1.0);

      }
    }


    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);
        _massVsY[alg].scale(crossSectionPerEvent()/picobarn, this);
        _chiVsMass[alg].scale(crossSectionPerEvent()/picobarn, this);
      }
    }


  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 _massVsY[2];

    /// The di-jet chi distribution binned in mass for akt6 and akt4 jets (array index is jet type from enum above)
    BinnedHistogram _chiVsMass[2];

  };



  RIVET_DECLARE_ALIASED_PLUGIN(ATLAS_2010_S8817804, ATLAS_2010_I871366);

}