Rivet Analyses Reference

ATLAS_2018_I1634970

ATLAS Inclusive jet and dijet cross section measurement at sqrt(s)=13TeV
Experiment: ATLAS (LHC)
Inspire ID: 1634970
Status: VALIDATED
Authors:
  • Jonathan Bossio
  • Zdenek Hubacek
References:Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • LHC 13TeV, ATLAS, QCD jet production, anti-kt R=0.4 jets

Inclusive jet and dijet cross-sections are measured in proton-proton collisions at a centre-of-mass energy of 13 TeV. The measurement uses a dataset with an integrated luminosity of 3.2 fb$^{-1}$ recorded in 2015 with the ATLAS detector at the Large Hadron Collider. Jets are identified using the anti-$k_\text{t}$ algorithm with a radius parameter value of $R = 0.4$. The inclusive jet cross-sections are measured double-differentially as a function of the jet transverse momentum, covering the range from 100 GeV to 3.5 TeV, and the absolute jet rapidity up to $|y| = 3$. The double-differential dijet production cross-sections are presented as a function of the dijet mass, covering the range from 300 GeV to 9 TeV, and the half absolute rapidity separation between the two leading jets within $|y| < 3$, $y^\ast$, up to $y^\ast = 3$. Next-to-leading-order, and next-to-next-to-leading-order for the inclusive jet measurement, perturbative QCD calculations corrected for non-perturbative and electroweak effects are compared to the measured cross-sections.

Source code: ATLAS_2018_I1634970.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"

namespace Rivet {


  /// @brief inclusive jet and dijet at 13 TeV
  class ATLAS_2018_I1634970 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2018_I1634970);

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

      const FinalState fs;
      declare(fs,"FinalState");
      FastJets fj04(fs, FastJets::ANTIKT, 0.4, JetAlg::Muons::ALL, JetAlg::Invisibles::DECAY);
      declare(fj04, "AntiKT04");

      // |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};
      double ystarbins[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; _pThistograms.add(ybins[i], ybins[i+1], book(tmp, i+1,1,1));}
      }
      // mjj histograms
      for(size_t i=0;i<nybins;++i){// loop over ystar bins
        {Histo1DPtr tmp; _mjjhistograms.add(ystarbins[i], ystarbins[i+1], book(tmp, i+7,1,1));}
      }

    }


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

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

      int nJets = kt4Jets.size();

      // Inclusive jet selection
      for(int ijet=0;ijet<nJets;++ijet){ // loop over jets
        FourMomentum jet = kt4Jets[ijet].momentum();
        // pT selection
        if(jet.pt()>100.0*GeV){
          // Fill distribution
          const double absy = jet.absrap();
          _pThistograms.fill(absy,jet.pt()/GeV);
        }
      }

      // Dijet selection
      if(nJets > 1){ // skip events with less than 2 jets passing pT>75GeV and |y|<3.0 cuts
        FourMomentum jet0  = kt4Jets[0].momentum(); 
        FourMomentum jet1  = kt4Jets[1].momentum();
        const double rap0  = jet0.rapidity();
        const double rap1  = jet1.rapidity();
        const double ystar = fabs(rap0-rap1)/2;
        const double mass  = (jet0 + jet1).mass(); 
        const double HT2   = jet0.pt()+jet1.pt();
        if(HT2>200*GeV && ystar<3.0){
          // Fill distribution
          _mjjhistograms.fill(ystar,mass/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() );
      _pThistograms.scale(xs_norm_factor, this);
      _mjjhistograms.scale(crossSectionPerEvent()/picobarn, this);

    }

  private:

    // The inclusive pT spectrum for akt4 jets
    BinnedHistogram _pThistograms;
    // The dijet mass spectrum for akt4 jets
    BinnedHistogram _mjjhistograms;

  };

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

}