rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2022_I2077570

Z + high transverse momentum jets at ATLAS
Experiment: ATLAS (LHC)
Inspire ID: 2077570
Status: VALIDATED
Authors:
  • Alexandre Laurier
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • pp -> Z + jets at 13 TeV with pT(jet) > 100 GeV

Cross-section measurements for a $Z$ boson produced in association with high-transverse-momentum jets ($p_\text{T}\geq 100$ GeV) and decaying into a charged-lepton pair ($e^+e^-$, $\mu^+\mu^-$) are presented. The measurements are performed using proton-proton collisions at $\sqrt{s} = 13$ TeV corresponding to an integrated luminosity of 139 fb${}^{-1}$ collected by the ATLAS experiment at the LHC. Measurements of angular correlations between the $Z$ boson and the closest jet are performed in events with at least one jet with $p_\text{T}\geq 500$ GeV. Event topologies of particular interest are the collinear emission of a $Z$ boson in dijet events and a boosted $Z$ boson recoiling against a jet. Fiducial cross sections are compared with state-of-the-art theoretical predictions. The data are found to agree with next-to-next-to-leading-order predictions by NNLOjet and with the next-to-leading-order multi-leg generators MadGraph5_aMC@NLO and Sherpa.

Source code: ATLAS_2022_I2077570.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FastJets.hh"
  4#include "Rivet/Projections/FinalState.hh"
  5#include "Rivet/Projections/IdentifiedFinalState.hh"
  6#include "Rivet/Projections/VetoedFinalState.hh"
  7#include "Rivet/Projections/PromptFinalState.hh"
  8#include "Rivet/Projections/LeptonFinder.hh"
  9
 10namespace Rivet {
 11
 12
 13  /// Collinear Z + Jets in pp at 13 TeV
 14  class ATLAS_2022_I2077570 : public Analysis {
 15  public:
 16
 17    /// Constructor
 18    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2022_I2077570);
 19
 20    /// Book histograms and initialise projections before the run
 21    void init() {
 22
 23      // Get options
 24      _mode = 0;
 25      if ( getOption("LMODE") == "ELEL") _mode = 1;
 26      if ( getOption("LMODE") == "MUMU") _mode = 2;
 27
 28      // AntiKt4TruthWZJets prescription
 29      // Photons
 30      FinalState all_photons(Cuts::abspid == PID::PHOTON);
 31
 32      // Muons
 33      PromptFinalState bare_mu(Cuts::abspid == PID::MUON, TauDecaysAs::PROMPT);
 34      LeptonFinder all_dressed_mu(bare_mu, all_photons, 0.1, Cuts::abseta < 2.5);
 35
 36      // Electrons
 37      PromptFinalState bare_el(Cuts::abspid == PID::ELECTRON, TauDecaysAs::PROMPT);
 38      LeptonFinder all_dressed_el(bare_el, all_photons, 0.1, Cuts::abseta < 2.5);
 39
 40      //Jet forming
 41      VetoedFinalState vfs(FinalState(Cuts::abseta < 4.5));
 42      vfs.addVetoOnThisFinalState(all_dressed_el);
 43      vfs.addVetoOnThisFinalState(all_dressed_mu);
 44
 45      FastJets jet(vfs, JetAlg::ANTIKT, 0.4, JetMuons::ALL, JetInvisibles::NONE);
 46      declare(jet, "Jets");
 47
 48
 49      // Current definition of leptons + jets
 50      PromptFinalState photons(Cuts::abspid == PID::PHOTON);
 51      PromptFinalState electrons(Cuts::abspid == PID::ELECTRON);
 52      PromptFinalState muons(Cuts::abspid == PID::MUON);
 53
 54      // Kinematic cuts for leptons
 55      const Cut cuts_lep = Cuts::pT > 25*GeV && Cuts::abseta < 2.5;
 56
 57      LeptonFinder dressed_electrons(electrons, photons, 0.1, cuts_lep);
 58      declare(dressed_electrons, "DressedElectrons");
 59
 60      LeptonFinder dressed_muons(muons, photons, 0.1, cuts_lep);
 61      declare(dressed_muons, "DressedMuons");
 62
 63
 64      book(_h["ZpT"],        1, 1, 1);
 65      book(_h["jetpT"],      2, 1, 1);
 66      book(_d["NJets"],      3, 1, 1);
 67      book(_d["NJets500"],   4, 1, 1);
 68      book(_h["minDR"],      5, 1, 1);
 69      book(_h["rZJ"],        6, 1, 1);
 70      book(_h["rZJ_coll"],   7, 1, 1);
 71      book(_h["rZJ_b2b"],    8, 1, 1);
 72      book(_d["NJets_coll"], 9, 1, 1);
 73      book(_d["NJets_b2b"], 10, 1, 1);
 74      book(_h["HT"],        11, 1, 1);
 75      book(_h["minDR600"],  12, 1, 1);
 76      book(_d["NJets600"],  13, 1, 1);
 77
 78    }
 79
 80    /// Perform the per-event analysis
 81    void analyze(const Event& event) {
 82
 83
 84      // access fiducial electrons and muons
 85      const Particle *l1 = nullptr, *l2 = nullptr;
 86      auto muons = apply<LeptonFinder>(event, "DressedMuons").dressedLeptons();
 87      auto elecs = apply<LeptonFinder>(event, "DressedElectrons").dressedLeptons();
 88
 89      // lepton-jet overlap removal (note: muons are not included in jet finding)
 90      // Jets eta < 2.5, pT > 30GeV for overlap removal
 91      Jets jets = apply<FastJets>(event, "Jets").jetsByPt(Cuts::pT > 30*GeV && Cuts::absrap < 2.5);
 92      // Remove all jets within dR < 0.2 of a dressed lepton
 93      idiscardIfAnyDeltaRLess(jets, muons, 0.2);
 94      idiscardIfAnyDeltaRLess(jets, elecs, 0.2);
 95
 96      // Remove leptons within dR < 0.4 of a jet
 97      idiscardIfAnyDeltaRLess(muons, jets, 0.4);
 98      idiscardIfAnyDeltaRLess(elecs, jets, 0.4);
 99
100      int nElecs = elecs.size();
101      int nMuons = muons.size();
102
103      if (nElecs + nMuons != 2) vetoEvent; // dilepton cut
104      if (_mode == 2 && nMuons !=2) vetoEvent;
105      if (_mode == 1 && nElecs !=2) vetoEvent;
106
107      if (elecs.size()==2){
108        l1=&elecs[0];
109        l2=&elecs[1];
110      }
111      else if (muons.size()==2){
112        l1=&muons[0];
113        l2=&muons[1];
114      }
115      else vetoEvent;
116
117      // if not e+e- or mu+mu- pair, veto
118      if (l1->pid() + l2->pid() !=0) vetoEvent;
119
120      // Dilepton selection :: Z mass peak.
121      FourMomentum ll = (l1->mom() + l2->mom());
122      double Zm  = ll.mass();
123      if ( !inRange(Zm/GeV, 71.0, 111.0) ) vetoEvent;
124      double ZpT = ll.pT();
125
126      // Calculate the observables
127      double jet0pT = 0.;
128      double cljetpT = 0.;
129      double minDR = 99.;
130
131      // Require jets to be above 100GeV
132      iselect(jets, Cuts::pT > 100*GeV);
133      double HTjet = sum(jets, Kin::pT, 0.);
134      for (const Jet& j : jets) {
135        // find minDR and closest jet to Z boson, only with 100GeV+ jets
136        double dr = deltaR(j, ll ,RAPIDITY);
137        if (dr < minDR) {
138          minDR = dr;
139          cljetpT = j.pT();
140        }
141      }
142      const size_t Njets = jets.size();
143
144      // NJets >= 1
145      if (Njets < 1) vetoEvent;
146      // Exclusive NJets, jet pT > 100 GeV
147
148      _d["NJets"]->fill(Njets);
149
150      // Z pT
151      _h["ZpT"]->fill(ZpT/GeV);
152
153      //Leading jet pT
154      jet0pT = jets[0].pT();
155      _h["jetpT"]->fill(jet0pT/GeV);
156
157      // HT
158      double HT = HTjet + l1->pT() + l2->pT();
159      _h["HT"]->fill(HT/GeV);
160
161      // HTJet > 600 GeV selection
162      if (HTjet >= 600.) {
163        double minDR_HTjet = 99.;
164        // closest jet to Z
165        for (const Jet& j : jets) {
166          const double dr = deltaR(j, ll, RAPIDITY);
167          if (dr < minDR_HTjet)  minDR_HTjet = dr;
168        }
169        // Fill histograms of HTjet > 600 GeV
170        _d["NJets600"]->fill(Njets);
171        _h["minDR600"]->fill(minDR_HTjet);
172      } // end of HTjet > 600 GeV
173
174      // Our high pT phase-space
175      if (jet0pT/GeV < 500.0) vetoEvent;
176
177      // Exclusive NJets, jet pT > 500 GeV
178      _d["NJets500"]->fill(Njets);
179
180      // Z/J pT ratio
181      _h["rZJ"]->fill(ZpT/cljetpT);
182      _h["minDR"]->fill(minDR);
183
184      // Phase space with DR<1.4
185      if (minDR < 1.4) {
186        _d["NJets_coll"]->fill(Njets);
187        _h["rZJ_coll"]->fill(ZpT/cljetpT);
188      // Phase space with DR>2.0
189      } else if (minDR >2.0) {
190        _d["NJets_b2b"]->fill(Njets);
191        _h["rZJ_b2b"]->fill(ZpT/cljetpT);
192      }
193
194    }
195
196    void finalize() {
197      double xsec = crossSectionPerEvent()/picobarn;
198      // Analysis measures Z->ll(ee or mumu)
199      // If file contains both Z->ee and Z->mumu, divide xs by 2
200      if (_mode == 0) xsec *= 0.5;
201      // Normalize to cross section.
202      scale(_h, xsec);
203      scale(_d, xsec);
204    } // end of finalize
205
206
207    // Define histograms
208    size_t _mode;
209    map<string, Histo1DPtr> _h;
210    map<string, BinnedHistoPtr<int>> _d;
211
212  };
213
214
215  RIVET_DECLARE_PLUGIN(ATLAS_2022_I2077570);
216
217}