Rivet Analyses Reference

ATLAS_2019_I1744201

Z+jet at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1744201
Status: VALIDATED
Authors:
  • Aliaksei Hrynevich
  • Deepak Kar
References:Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • pp -> Z(ee)+jets at 8 TeV

The inclusive cross-section for jet production in association with a $Z$ boson decaying into an electron-positron pair is measured as a function of the transverse momentum and the absolute rapidity of jets using 19.9 fb$^{-1}$ of $\sqrt{s}=8$ TeV proton-proton collision data collected with the ATLAS detector at the Large Hadron Collider. The measured Z+ jets cross-section is unfolded to the particle level. The cross-section is compared with state-of-the-art Standard Model calculations, including the next-to-leading-order and next-to-next-to-leading-order perturbative QCD calculations, corrected for non-perturbative and QED radiation effects. The results of the measurements cover final-state jets with transverse momenta up to 1 TeV, and show good agreement with fixed-order calculations.

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

namespace Rivet{

  /// @brief:  Z+jet at 8 TeV
  class ATLAS_2019_I1744201 : public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2019_I1744201);
  
    void init() {

      const FinalState fs(Cuts::abseta < 5.0);
      Cut cut = Cuts::abseta < 2.47 && Cuts::pT >= 20*GeV;
          
      ZFinder zfinder_el(fs, cut, PID::ELECTRON, 66*GeV, 116*GeV, 0.1, ZFinder::ChargedLeptons::PROMPT);
      declare(zfinder_el, "ZFinder_el");
      
      declare(FastJets(zfinder_el.remainingFinalState(), FastJets::ANTIKT, 0.4, 
                       JetAlg::Muons::NONE, JetAlg::Invisibles::NONE), "AKT04");

      h_jet_y_pt.resize(6);
      for (size_t iPtBin=0; iPtBin < h_jet_y_pt.size(); ++iPtBin) {
        book(h_jet_y_pt[iPtBin], iPtBin+2, 1, 1); 
      }

    }
    
    void analyze(const Event& event) {

      // electrons selection
      const ZFinder& zfinder = apply<ZFinder>(event, "ZFinder_el");
      if ( zfinder.bosons().size() != 1)  vetoEvent; 

      const Particles& leptons = zfinder.constituents();
      if ( leptons.size() != 2)  vetoEvent; 

      if (deltaR(leptons[0], leptons[1]) < 0.2)  vetoEvent; 


      // jets selection
      Jets jets = apply<FastJets>(event, "AKT04").jetsByPt( Cuts::pT > 25*GeV && Cuts::absrap < 3.4 );
      idiscardIfAnyDeltaRLess(jets, leptons, 0.4);
      if (jets.empty())  vetoEvent;  // require at least one jet in event

      for (const Jet& jet : jets) {
        const double jet_pt = jet.pT() / GeV;
        for(size_t iPtBin = 0; iPtBin < (ptBins.size() - 1); ++iPtBin) {
          if (jet_pt >= ptBins[iPtBin] && jet_pt < ptBins[iPtBin+1]) {
            h_jet_y_pt[iPtBin]->fill(jet.absrap());
          }
        }
	    }
   }
    
  void finalize() {
 
    const double norm = crossSection()/femtobarn/sumOfWeights();
    for(int iPtBin=0; iPtBin < 6; ++iPtBin){
    	scale(h_jet_y_pt[iPtBin], norm / ( ptBins[iPtBin+1] - ptBins[iPtBin] ));
    }
      
  }

   protected:

     vector<double> ptBins = { 25., 50., 100., 200., 300., 400., 1050. };
  
   private:
     vector<Histo1DPtr> h_jet_y_pt;

  };
  
  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(ATLAS_2019_I1744201);
}