rivet is hosted by Hepforge, IPPP Durham

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// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/FastJets.hh"
  5
  6namespace Rivet {
  7
  8
  9  /// @brief inclusive jet and dijet at 13 TeV
 10  class ATLAS_2018_I1634970 : public Analysis {
 11  public:
 12
 13    /// Constructor
 14    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2018_I1634970);
 15
 16    /// Book histograms and initialise projections before the run
 17    void init() {
 18
 19      const FinalState fs;
 20      declare(fs,"FinalState");
 21      FastJets fj04(fs, FastJets::ANTIKT, 0.4, JetAlg::Muons::ALL, JetAlg::Invisibles::DECAY);
 22      declare(fj04, "AntiKT04");
 23
 24      // |y| and ystar bins
 25      const int nybins          = 6;
 26      double ybins[nybins+1]     = { 0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0};
 27      double ystarbins[nybins+1] = { 0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0};
 28
 29      // Book histograms
 30      // pT histograms
 31      for(size_t i=0;i<nybins;++i){// loop over |y| bins
 32        {Histo1DPtr tmp; _pThistograms.add(ybins[i], ybins[i+1], book(tmp, i+1,1,1));}
 33      }
 34      // mjj histograms
 35      for(size_t i=0;i<nybins;++i){// loop over ystar bins
 36        {Histo1DPtr tmp; _mjjhistograms.add(ystarbins[i], ystarbins[i+1], book(tmp, i+7,1,1));}
 37      }
 38
 39    }
 40
 41
 42    /// Perform the per-event analysis
 43    void analyze(const Event& event) {
 44
 45      const Jets& kt4Jets = apply<FastJets>(event, "AntiKT04").jetsByPt(Cuts::pT > 75*GeV && Cuts::absrap < 3.0);
 46
 47      int nJets = kt4Jets.size();
 48
 49      // Inclusive jet selection
 50      for(int ijet=0;ijet<nJets;++ijet){ // loop over jets
 51        FourMomentum jet = kt4Jets[ijet].momentum();
 52        // pT selection
 53        if(jet.pt()>100.0*GeV){
 54          // Fill distribution
 55          const double absy = jet.absrap();
 56          _pThistograms.fill(absy,jet.pt()/GeV);
 57        }
 58      }
 59
 60      // Dijet selection
 61      if(nJets > 1){ // skip events with less than 2 jets passing pT>75GeV and |y|<3.0 cuts
 62        FourMomentum jet0  = kt4Jets[0].momentum(); 
 63        FourMomentum jet1  = kt4Jets[1].momentum();
 64        const double rap0  = jet0.rapidity();
 65        const double rap1  = jet1.rapidity();
 66        const double ystar = fabs(rap0-rap1)/2;
 67        const double mass  = (jet0 + jet1).mass(); 
 68        const double HT2   = jet0.pt()+jet1.pt();
 69        if(HT2>200*GeV && ystar<3.0){
 70          // Fill distribution
 71          _mjjhistograms.fill(ystar,mass/GeV);
 72        }
 73      }
 74
 75    }
 76
 77
 78    /// Normalise histograms etc., after the run
 79    void finalize() {
 80
 81      const double xs_pb( crossSection() / picobarn );
 82      const double sumW( sumOfWeights() );
 83      const double xs_norm_factor( 0.5*xs_pb / sumW );
 84     
 85      MSG_DEBUG( "Cross-Section/pb     : " << xs_pb       );
 86      MSG_DEBUG( "ZH                   : " << crossSectionPerEvent()/ picobarn);
 87      MSG_DEBUG( "Sum of weights       : " << sumW        );
 88      MSG_DEBUG( "nEvents              : " << numEvents() );
 89      _pThistograms.scale(xs_norm_factor, this);
 90      _mjjhistograms.scale(crossSectionPerEvent()/picobarn, this);
 91
 92    }
 93
 94  private:
 95
 96    // The inclusive pT spectrum for akt4 jets
 97    BinnedHistogram _pThistograms;
 98    // The dijet mass spectrum for akt4 jets
 99    BinnedHistogram _mjjhistograms;
100
101  };
102
103  // The hook for the plugin system
104  RIVET_DECLARE_PLUGIN(ATLAS_2018_I1634970);
105
106}