rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2017_I1589844

$k_T$ splittings in $Z$ events at 8 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1589844
Status: VALIDATED
Authors:
  • Christian Gutschow
  • Frank Siegert
References: Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • $pp \to Z(\to ee/\mu\mu) +$ jets at 8 TeV

A measurement of the splitting scales occurring in the $k_\text{t}$ jet-clustering algorithm is presented for final states containing a $Z$ boson. The measurement is done using 20.2 fb$^{-1}$ of proton-proton collision data collected at a centre-of-mass energy of $\sqrt{s} = 8$ TeV by the ATLAS experiment at the LHC in 2012. The measurement is based on charged-particle track information, which is measured with excellent precision in the $p_\text{T}$ region relevant for the transition between the perturbative and the non-perturbative regimes. The data distributions are corrected for detector effects, and are found to deviate from state-of-the-art predictions in various regions of the observables. If no mode specified, will try to fill both electron and muon plots. If EL or MU is specified, only the relevant plots will be filled.

Source code: ATLAS_2017_I1589844.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/IdentifiedFinalState.hh"
  5#include "Rivet/Projections/LeptonFinder.hh"
  6#include "Rivet/Projections/FastJets.hh"
  7#include "Rivet/Projections/VetoedFinalState.hh"
  8#include "Rivet/Projections/ChargedFinalState.hh"
  9
 10namespace Rivet {
 11
 12
 13  /// kT splittings in Z events at 8 TeV
 14  class ATLAS_2017_I1589844 : public Analysis {
 15  public:
 16
 17    /// @name Constructors etc.
 18    /// @{
 19
 20    /// Constructors
 21    ATLAS_2017_I1589844(const string name="ATLAS_2017_I1589844",
 22                        const string ref_data="ATLAS_2017_I1589844") : Analysis(name) {
 23      setRefDataName(ref_data);
 24    }
 25
 26    /// @}
 27
 28
 29    /// @name Analysis methods
 30    /// @{
 31
 32    /// Book histograms and initialise projections before the run
 33    void init() {
 34
 35      // Get options from the new option system
 36      _mode = 0;
 37      if ( getOption("LMODE") == "EL" ) _mode = 1;
 38      if ( getOption("LMODE") == "MU" ) _mode = 2;
 39
 40      Cut cuts_mu = (Cuts::pT > 25*GeV) && (Cuts::abseta < 2.4);
 41      Cut cuts_el = Cuts::pT > 25*GeV && (Cuts::abseta <= 1.37 || (Cuts::abseta >= 1.52 && Cuts::abseta < 2.47));
 42
 43      FinalState fs;
 44      IdentifiedFinalState bare_mu(fs);
 45      bare_mu.acceptIdPair(PID::MUON);
 46      IdentifiedFinalState bare_el(fs);
 47      bare_el.acceptIdPair(PID::ELECTRON);
 48      LeptonFinder muons(bare_mu, fs, 0.1, cuts_mu);
 49      LeptonFinder elecs(bare_el, fs, 0.1, cuts_el);
 50      declare(muons, "muons");
 51      declare(elecs, "elecs");
 52
 53      ChargedFinalState cfs(Cuts::abseta < 2.5 && Cuts::pT > 0.4*GeV);
 54      VetoedFinalState jet_fs(cfs);
 55      jet_fs.addVetoOnThisFinalState(muons);
 56      jet_fs.addVetoOnThisFinalState(elecs);
 57      declare(FastJets(jet_fs, JetAlg::KT, 0.4), "Kt04Jets");
 58      declare(FastJets(jet_fs, JetAlg::KT, 1.0), "Kt10Jets");
 59
 60      VetoedFinalState jet_fs_all(FinalState(Cuts::abseta < 2.5 && Cuts::pT > 0.4*GeV));
 61      jet_fs_all.addVetoOnThisFinalState(muons);
 62      jet_fs_all.addVetoOnThisFinalState(elecs);
 63      FastJets jetpro04_all(jet_fs_all, JetAlg::KT, 0.4);
 64      jetpro04_all.useInvisibles();
 65      declare(jetpro04_all, "Kt04Jets_all");
 66      FastJets jetpro10_all(jet_fs_all, JetAlg::KT, 1.0);
 67      jetpro10_all.useInvisibles();
 68      declare(jetpro10_all, "Kt10Jets_all");
 69
 70      // Histograms with data binning
 71      _ndij = 8;
 72      for (size_t i = 0; i < _ndij; ++i) {
 73        if (_mode == 0 || _mode == 1) {
 74          string label = "el_d" + to_str(i) + "_kT4";
 75          book(_h[label], i + 1, 1, 1);
 76          book(_h[label + "_all"], i + 1, 1, 5);
 77          label = "el_d" + to_str(i) + "_kT10";
 78          book(_h[label], i + 1, 1, 3);
 79          book(_h[label + "_all"], i + 1, 1, 7);
 80        }
 81        if (_mode == 0 || _mode == 2) {
 82          string label = "mu_d" + to_str(i) + "_kT4";
 83          book(_h[label], i + 1, 1, 2);
 84          book(_h[label + "_all"], i + 1, 1, 6);
 85          label = "mu_d" + to_str(i) + "_kT10";
 86          book(_h[label], i + 1, 1, 4);
 87          book(_h[label + "_all"], i + 1, 1, 8);
 88        }
 89      }
 90    }
 91
 92
 93    /// Perform the per-event analysis
 94    void analyze(const Event& e) {
 95
 96      // Check we have a Z candidate:
 97      const DressedLeptons& muons = apply<LeptonFinder>(e, "muons").dressedLeptons();
 98      const DressedLeptons& elecs = apply<LeptonFinder>(e, "elecs").dressedLeptons();
 99
100      bool e_ok = elecs.size() == 2 && muons.empty();
101      bool m_ok = elecs.empty() && muons.size() == 2;
102      if (_mode == 0 && !e_ok && !m_ok) vetoEvent;
103      if (_mode == 1 && !e_ok )  vetoEvent;
104      if (_mode == 2 && !m_ok )  vetoEvent;
105
106      string lep_type = elecs.size()? "el_" : "mu_";
107
108      const DressedLeptons& leptons = elecs.size()? elecs : muons;
109      if (leptons[0].charge()*leptons[1].charge() > 0) vetoEvent;
110
111      const double dilepton_mass = (leptons[0].momentum() + leptons[1].momentum()).mass();
112      if (!inRange(dilepton_mass, 71*GeV, 111*GeV)) vetoEvent;
113
114      // Get kT splitting scales (charged particles only)
115      const FastJets& jetpro04 = apply<FastJets>(e, "Kt04Jets");
116      const shared_ptr<fastjet::ClusterSequence> seq04 = jetpro04.clusterSeq();
117      for (size_t i = 0; i < min(_ndij, (size_t)seq04->n_particles()); ++i) {
118        const double dij = sqrt(seq04->exclusive_dmerge_max(i))/GeV;
119        if (dij <= 0.0) continue;
120        const string label = lep_type + "d" + to_str(i) + "_kT4";
121        _h[label]->fill(dij);
122      }
123      const FastJets& jetpro10 = apply<FastJets>(e, "Kt10Jets");
124      const shared_ptr<fastjet::ClusterSequence> seq10 = jetpro10.clusterSeq();
125      for (size_t i = 0; i < min(_ndij, (size_t)seq10->n_particles()); ++i) {
126        const double dij = sqrt(seq10->exclusive_dmerge_max(i))/GeV;
127        if (dij <= 0.0) continue;
128        const string label = lep_type + "d" + to_str(i) + "_kT10";
129        _h[label]->fill(dij);
130      }
131
132      // Get kT splitting scales (all particles)
133      const FastJets& jetpro04_all = apply<FastJets>(e, "Kt04Jets_all");
134      const shared_ptr<fastjet::ClusterSequence> seq04_all = jetpro04_all.clusterSeq();
135      for (size_t i = 0; i < min(_ndij, (size_t)seq04_all->n_particles()); ++i) {
136        const double dij = sqrt(seq04_all->exclusive_dmerge_max(i))/GeV;
137        if (dij <= 0.0) continue;
138        const string label = lep_type + "d" + to_str(i) + "_kT4_all";
139        _h[label]->fill(dij);
140      }
141      const FastJets& jetpro10_all = apply<FastJets>(e, "Kt10Jets_all");
142      const shared_ptr<fastjet::ClusterSequence> seq10_all = jetpro10_all.clusterSeq();
143      for (size_t i = 0; i < min(_ndij, (size_t)seq10_all->n_particles()); ++i) {
144        const double dij = sqrt(seq10_all->exclusive_dmerge_max(i))/GeV;
145        if (dij <= 0.0) continue;
146        const string label = lep_type + "d" + to_str(i) + "_kT10_all";
147        _h[label]->fill(dij);
148      }
149
150    }
151
152
153    /// Normalise histograms etc., after the run
154    void finalize() {
155      const double sf = crossSectionPerEvent();
156      for (auto& kv : _h) scale(kv.second, sf);
157    }
158
159    /// @}
160
161
162  protected:
163
164    // Data members like post-cuts event weight counters go here
165    size_t _mode, _ndij;
166
167
168  private:
169
170    // Histograms
171    map<string, Histo1DPtr> _h;
172
173  };
174
175
176  RIVET_DECLARE_PLUGIN(ATLAS_2017_I1589844);
177
178}