rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2014_I1268975

High-mass dijet cross section
Experiment: ATLAS (LHC)
Inspire ID: 1268975
Status: VALIDATED
Authors:
  • Christopher Meyer
References: Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • QCD jet production with a minimum leading jet pT of 100 GeV and minimum second jet pT of 50 GeV at 7 TeV.

Double-differential dijet cross sections measured in pp collisions at the LHC with a 7 TeV centre-of-mass energy are presented as functions of dijet mass and rapidity separation of the two highest-pT jets. These measurements are obtained using data corresponding to an integrated luminosity of 4.5/fb, recorded by the ATLAS detector in 2011. The data are corrected for detector effects so that cross sections are presented at the particle level. Cross sections are measured up to 5 TeV dijet mass using jets reconstructed with the anti-kt algorithm for values of the jet radius parameter of 0.4 and 0.6. The cross sections are compared with next-to-leading-order perturbative QCD calculations by NLOJET++ corrected to account for non-perturbative effects. Comparisons with POWHEG predictions, using a next-to-leading-order matrix element calculation interfaced to a parton-shower Monte Carlo simulation, are also shown. Electroweak effects are accounted for in both cases. The quantitative comparison of data and theoretical predictions obtained using various parameterizations of the parton distribution functions is performed using a frequentist method. An example setting a lower limit on the compositeness scale for a model of contact interactions is presented, showing that the unfolded results can be used to constrain contributions to dijet production beyond that predicted by the Standard Model.

Source code: ATLAS_2014_I1268975.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  /// Jet mass as a function of ystar
 10  class ATLAS_2014_I1268975 : public Analysis {
 11  public:
 12
 13    /// Constructor
 14    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2014_I1268975);
 15
 16
 17    /// @name Analysis methods
 18    /// @{
 19
 20    /// Book histograms and initialise projections before the run
 21    void init() {
 22
 23      const FinalState fs;
 24      declare(fs,"FinalState");
 25
 26      FastJets fj04(fs,  JetAlg::ANTIKT, 0.4);
 27      fj04.useInvisibles();
 28      declare(fj04, "AntiKT04");
 29
 30      FastJets fj06(fs,  JetAlg::ANTIKT, 0.6);
 31      fj06.useInvisibles();
 32      declare(fj06, "AntiKT06");
 33
 34      vector<double> ystarbins{ 0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0};
 35
 36      size_t massDsOffset(0);
 37      for (size_t alg = 0; alg < 2; ++alg) {
 38        book(_mass[alg], ystarbins);
 39        for (auto& b : _mass[alg]->bins()) {
 40          book(b, 1+massDsOffset, 1, b.index());
 41        }
 42        ++massDsOffset;
 43      }
 44    }
 45
 46
 47    /// Perform the per-event analysis
 48    void analyze(const Event& event) {
 49      Jets jetAr[2];
 50      jetAr[AKT4] = apply<FastJets>(event, "AntiKT04").jetsByPt(Cuts::pT > 50*GeV);
 51      jetAr[AKT6] = apply<FastJets>(event, "AntiKT06").jetsByPt(Cuts::pT > 50*GeV);
 52
 53      // Loop over jet "radii" used in analysis
 54      for (size_t alg = 0; alg < 2; ++alg) {
 55
 56        // Identify dijets
 57        vector<FourMomentum> leadjets;
 58        for (const Jet& jet : jetAr[alg]) {
 59          if (jet.absrap() < 3.0 && leadjets.size() < 2) {
 60            if (leadjets.empty() && jet.pT() < 100*GeV) continue;
 61            leadjets.push_back(jet.momentum());
 62          }
 63        }
 64
 65        // Make sure we have a leading jet with pT > 100 GeV and a second to leading jet with pT > 50 GeV
 66        if (leadjets.size() < 2) {
 67          MSG_DEBUG("Could not find two suitable leading jets");
 68          continue;
 69        }
 70
 71        const double y1    = leadjets[0].rapidity();
 72        const double y2    = leadjets[1].rapidity();
 73        const double ystar = fabs(y1-y2) / 2.;
 74        const double m     = (leadjets[0] + leadjets[1]).mass();
 75
 76        // Fill mass histogram
 77        _mass[alg]->fill(ystar, m/TeV);
 78      }
 79    }
 80
 81
 82    /// Normalise histograms etc., after the run
 83    void finalize() {
 84      scale(_mass, crossSectionPerEvent()/picobarn);
 85      divByGroupWidth(_mass);
 86    }
 87
 88    /// @}
 89
 90
 91  private:
 92
 93    // Data members like post-cuts event weight counters go here
 94    enum Alg { AKT4=0, AKT6=1 };
 95
 96    /// The di-jet mass spectrum binned in rapidity for akt6 and akt4 jets (array index is jet type from enum above)
 97    Histo1DGroupPtr _mass[2];
 98
 99  };
100
101
102  RIVET_DECLARE_PLUGIN(ATLAS_2014_I1268975);
103
104}