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
  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
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/IdentifiedFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/DressedLeptons.hh"

namespace Rivet {

  /// Colinear Z + Jets in pp at 13 TeV
  class ATLAS_2022_I2077570 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2022_I2077570);

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

      // Get options
      _mode = 0;
      if ( getOption("LMODE") == "ELEL") _mode = 1;
      if ( getOption("LMODE") == "MUMU") _mode = 2;

      // AntiKt4TruthWZJets prescription
      // Photons
      FinalState all_photons(Cuts::abspid == PID::PHOTON);
      
      // Muons
      PromptFinalState bare_mu(Cuts::abspid == PID::MUON, true); // true = use muons from prompt tau decays
      DressedLeptons all_dressed_mu(all_photons, bare_mu, 0.1, Cuts::abseta < 2.5, true);
      
      // Electrons
      PromptFinalState bare_el(Cuts::abspid == PID::ELECTRON, true); // true = use electrons from prompt tau decays
      DressedLeptons all_dressed_el(all_photons, bare_el, 0.1, Cuts::abseta < 2.5, true);
      
      //Jet forming
      VetoedFinalState vfs(FinalState(Cuts::abseta < 4.5));
      vfs.addVetoOnThisFinalState(all_dressed_el);
      vfs.addVetoOnThisFinalState(all_dressed_mu);
            
      FastJets jet(vfs, FastJets::ANTIKT, 0.4, JetAlg::Muons::ALL, JetAlg::Invisibles::NONE);
      declare(jet, "Jets");


      // Current definition of leptons + jets
      PromptFinalState photons(Cuts::abspid == PID::PHOTON);
      PromptFinalState electrons(Cuts::abspid == PID::ELECTRON);
      PromptFinalState muons(Cuts::abspid == PID::MUON);

      // Kinematic cuts for leptons
      const Cut cuts_lep = Cuts::pT > 25*GeV && Cuts::abseta < 2.5;

      DressedLeptons dressed_electrons(photons, electrons, 0.1, cuts_lep);
      declare(dressed_electrons, "DressedElectrons");

      DressedLeptons dressed_muons(photons, muons, 0.1, cuts_lep);
      declare(dressed_muons, "DressedMuons");


      book(_h["ZpT"],        1, 1, 1);
      book(_h["jetpT"],      2, 1, 1);
      book(_h["NJets"],      3, 1, 1);
      book(_h["NJets500"],   4, 1, 1);
      book(_h["minDR"],      5, 1, 1);
      book(_h["rZJ"],        6, 1, 1);
      book(_h["rZJ_coll"],   7, 1, 1);
      book(_h["rZJ_b2b"],    8, 1, 1);
      book(_h["NJets_coll"], 9, 1, 1);
      book(_h["NJets_b2b"], 10, 1, 1);
      book(_h["HT"],        11, 1, 1);
      book(_h["minDR600"],  12, 1, 1);
      book(_h["NJets600"],  13, 1, 1);

    }

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


      // access fiducial electrons and muons
      const Particle *l1 = nullptr, *l2 = nullptr;
      auto muons = apply<DressedLeptons>(event, "DressedMuons").dressedLeptons();
      auto elecs = apply<DressedLeptons>(event, "DressedElectrons").dressedLeptons();

      // lepton-jet overlap removal (note: muons are not included in jet finding)
      // Jets eta < 2.5, pT > 30GeV for overlap removal
      Jets jets = apply<FastJets>(event, "Jets").jetsByPt(Cuts::pT > 30*GeV && Cuts::absrap < 2.5);
      // Remove all jets within dR < 0.2 of a dressed lepton
      idiscardIfAnyDeltaRLess(jets, muons, 0.2);
      idiscardIfAnyDeltaRLess(jets, elecs, 0.2);
      
      // Remove leptons within dR < 0.4 of a jet
      idiscardIfAnyDeltaRLess(muons, jets, 0.4);
      idiscardIfAnyDeltaRLess(elecs, jets, 0.4);

      int nElecs = elecs.size();
      int nMuons = muons.size();

      if (nElecs + nMuons != 2) vetoEvent; // dilepton cut
      if (_mode == 2 && nMuons !=2) vetoEvent;
      if (_mode == 1 && nElecs !=2) vetoEvent;

      if (elecs.size()==2){
        l1=&elecs[0];
        l2=&elecs[1];
      }
      else if (muons.size()==2){
        l1=&muons[0];
        l2=&muons[1];
      }
      else vetoEvent;

      // if not e+e- or mu+mu- pair, veto
      if (l1->pid() + l2->pid() !=0) vetoEvent;

      // Dilepton selection :: Z mass peak.
      FourMomentum ll = (l1->mom() + l2->mom());
      double Zm  = ll.mass(); 
      if ( !inRange(Zm/GeV, 71.0, 111.0) ) vetoEvent;
      double ZpT = ll.pT();
      
      // Calculate the observables
      double jet0pT = 0.;
      double cljetpT = 0.;
      double minDR = 99.;

      // Require jets to be above 100GeV
      ifilter_select(jets, Cuts::pT > 100*GeV);
      double HTjet = sum(jets, Kin::pT, 0.);
      for (const Jet& j : jets) {
        // find minDR and closest jet to Z boson, only with 100GeV+ jets
        double dr = deltaR(j, ll ,RAPIDITY);
        if (dr < minDR) {
          minDR = dr;
          cljetpT = j.pT();
        }
      }
      const size_t Njets = jets.size();

      // NJets >= 1      
      if (Njets < 1) vetoEvent;
      // Exclusive NJets, jet pT > 100 GeV

      _h["NJets"]->fill(Njets);

      // Z pT
      _h["ZpT"]->fill(ZpT/GeV);
      
      //Leading jet pT 
      jet0pT = jets[0].pT();
      _h["jetpT"]->fill(jet0pT/GeV);

      // HT
      double HT = HTjet + l1->pT() + l2->pT();
      _h["HT"]->fill(HT/GeV);

      // HTJet > 600 GeV selection
      if (HTjet >= 600.) {
        double minDR_HTjet = 99.;
        // closest jet to Z
        for (const Jet& j : jets) {
          const double dr = deltaR(j, ll, RAPIDITY);
          if (dr < minDR_HTjet)  minDR_HTjet = dr;
        }
        // Fill histograms of HTjet > 600 GeV
        _h["NJets600"]->fill(Njets);
        _h["minDR600"]->fill(minDR_HTjet);
      } // end of HTjet > 600 GeV

      // Our high pT phase-space
      if (jet0pT/GeV < 500.0) vetoEvent;
      
      // Exclusive NJets, jet pT > 500 GeV
      _h["NJets500"]->fill(Njets);
      
      // Z/J pT ratio
      _h["rZJ"]->fill(ZpT/cljetpT);
      _h["minDR"]->fill(minDR);
      
      // Phase space with DR<1.4       
      if (minDR < 1.4) {
        _h["NJets_coll"]->fill(Njets);
        _h["rZJ_coll"]->fill(ZpT/cljetpT);
      // Phase space with DR>2.0
      } else if (minDR >2.0) {
        _h["NJets_b2b"]->fill(Njets);
        _h["rZJ_b2b"]->fill(ZpT/cljetpT);
      }

    }

    void finalize() {
      double xsec = crossSectionPerEvent()/picobarn;
      // Analysis measures Z->ll(ee or mumu)
      // If file contains both Z->ee and Z->mumu, divide xs by 2
      if (_mode == 0) xsec *= 0.5;
      // Normalize to cross section.
      scale(_h, xsec);
    } // end of finalize

    //@}

    // define histograms
    size_t _mode;
    map<string, Histo1DPtr> _h;
 
  };

  RIVET_DECLARE_PLUGIN(ATLAS_2022_I2077570);
}