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
  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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/IdentifiedFinalState.hh"
#include "Rivet/Projections/DressedLeptons.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/ChargedFinalState.hh"

namespace Rivet {


  /// kT splittings in Z events at 8 TeV
  class ATLAS_2017_I1589844 : public Analysis {
  public:

    /// @name Constructors etc.
    //@{

    /// Constructors
    ATLAS_2017_I1589844(const string name="ATLAS_2017_I1589844",
                        const string ref_data="ATLAS_2017_I1589844") : Analysis(name) {
      setRefDataName(ref_data);
    }

    //@}


    /// @name Analysis methods
    //@{

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

      // Get options from the new option system
      _mode = 0;
      if ( getOption("LMODE") == "EL" ) _mode = 1;
      if ( getOption("LMODE") == "MU" ) _mode = 2;

      const FinalState fs;

      const Cut cuts_mu = (Cuts::pT > 25*GeV) && (Cuts::abseta < 2.4);
      const Cut cuts_el = Cuts::pT > 25*GeV && (Cuts::abseta <= 1.37 || (Cuts::abseta >= 1.52 && Cuts::abseta < 2.47));

      IdentifiedFinalState bare_mu(fs);
      bare_mu.acceptIdPair(PID::MUON);
      IdentifiedFinalState bare_el(fs);
      bare_el.acceptIdPair(PID::ELECTRON);
      const DressedLeptons muons(fs, bare_mu, 0.1, cuts_mu, true);
      const DressedLeptons elecs(fs, bare_el, 0.1, cuts_el, true);
      declare(muons, "muons");
      declare(elecs, "elecs");

      const ChargedFinalState cfs(Cuts::abseta < 2.5 && Cuts::pT > 0.4*GeV);
      VetoedFinalState jet_fs(cfs);
      jet_fs.addVetoOnThisFinalState(muons);
      jet_fs.addVetoOnThisFinalState(elecs);
      declare(FastJets(jet_fs, FastJets::KT, 0.4), "Kt04Jets");
      declare(FastJets(jet_fs, FastJets::KT, 1.0), "Kt10Jets");

      VetoedFinalState jet_fs_all(FinalState(Cuts::abseta < 2.5 && Cuts::pT > 0.4*GeV));
      jet_fs_all.addVetoOnThisFinalState(muons);
      jet_fs_all.addVetoOnThisFinalState(elecs);
      FastJets jetpro04_all(jet_fs_all, FastJets::KT, 0.4);
      jetpro04_all.useInvisibles();
      declare(jetpro04_all, "Kt04Jets_all");
      FastJets jetpro10_all(jet_fs_all, FastJets::KT, 1.0);
      jetpro10_all.useInvisibles();
      declare(jetpro10_all, "Kt10Jets_all");

      // Histograms with data binning
      _ndij = 8;
      for (size_t i = 0; i < _ndij; ++i) {
        if (_mode == 0 || _mode == 1) {
          string label = "el_d" + to_str(i) + "_kT4";
          book(_h[label], i + 1, 1, 1);
          book(_h[label + "_all"], i + 1, 1, 5);
          label = "el_d" + to_str(i) + "_kT10";
          book(_h[label], i + 1, 1, 3);
          book(_h[label + "_all"], i + 1, 1, 7);
        }
        if (_mode == 0 || _mode == 2) {
          string label = "mu_d" + to_str(i) + "_kT4";
          book(_h[label], i + 1, 1, 2);
          book(_h[label + "_all"], i + 1, 1, 6);
          label = "mu_d" + to_str(i) + "_kT10";
          book(_h[label], i + 1, 1, 4);
          book(_h[label + "_all"], i + 1, 1, 8);
        }
      }
    }


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

      // Check we have a Z candidate:
      const vector<DressedLepton>& muons = apply<DressedLeptons>(e, "muons").dressedLeptons();
      const vector<DressedLepton>& elecs = apply<DressedLeptons>(e, "elecs").dressedLeptons();

      bool e_ok = elecs.size() == 2 && muons.empty();
      bool m_ok = elecs.empty() && muons.size() == 2;
      if (_mode == 0 && !e_ok && !m_ok) vetoEvent;
      if (_mode == 1 && !e_ok )  vetoEvent;
      if (_mode == 2 && !m_ok )  vetoEvent;

      string lep_type = elecs.size()? "el_" : "mu_";

      const vector<DressedLepton>& leptons = elecs.size()? elecs : muons;
      if (leptons[0].charge()*leptons[1].charge() > 0) vetoEvent;

      const double dilepton_mass = (leptons[0].momentum() + leptons[1].momentum()).mass();
      if (!inRange(dilepton_mass, 71*GeV, 111*GeV)) vetoEvent;

      // Get kT splitting scales (charged particles only)
      const FastJets& jetpro04 = applyProjection<FastJets>(e, "Kt04Jets");
      const shared_ptr<fastjet::ClusterSequence> seq04 = jetpro04.clusterSeq();
      for (size_t i = 0; i < min(_ndij, (size_t)seq04->n_particles()); ++i) {
        const double dij = sqrt(seq04->exclusive_dmerge_max(i))/GeV;
        if (dij <= 0.0) continue;
        const string label = lep_type + "d" + to_str(i) + "_kT4";
        _h[label]->fill(dij);
      }
      const FastJets& jetpro10 = applyProjection<FastJets>(e, "Kt10Jets");
      const shared_ptr<fastjet::ClusterSequence> seq10 = jetpro10.clusterSeq();
      for (size_t i = 0; i < min(_ndij, (size_t)seq10->n_particles()); ++i) {
        const double dij = sqrt(seq10->exclusive_dmerge_max(i))/GeV;
        if (dij <= 0.0) continue;
        const string label = lep_type + "d" + to_str(i) + "_kT10";
        _h[label]->fill(dij);
      }

      // Get kT splitting scales (all particles)
      const FastJets& jetpro04_all = applyProjection<FastJets>(e, "Kt04Jets_all");
      const shared_ptr<fastjet::ClusterSequence> seq04_all = jetpro04_all.clusterSeq();
      for (size_t i = 0; i < min(_ndij, (size_t)seq04_all->n_particles()); ++i) {
        const double dij = sqrt(seq04_all->exclusive_dmerge_max(i))/GeV;
        if (dij <= 0.0) continue;
        const string label = lep_type + "d" + to_str(i) + "_kT4_all";
        _h[label]->fill(dij);
      }
      const FastJets& jetpro10_all = applyProjection<FastJets>(e, "Kt10Jets_all");
      const shared_ptr<fastjet::ClusterSequence> seq10_all = jetpro10_all.clusterSeq();
      for (size_t i = 0; i < min(_ndij, (size_t)seq10_all->n_particles()); ++i) {
        const double dij = sqrt(seq10_all->exclusive_dmerge_max(i))/GeV;
        if (dij <= 0.0) continue;
        const string label = lep_type + "d" + to_str(i) + "_kT10_all";
        _h[label]->fill(dij);
      }

    }


    /// Normalise histograms etc., after the run
    void finalize() {
      const double sf = crossSectionPerEvent();
      for (auto& kv : _h) scale(kv.second, sf);
    }

    //@}


  protected:

    // Data members like post-cuts event weight counters go here
    size_t _mode, _ndij;


  private:

    // Histograms
    map<string, Histo1DPtr> _h;

  };


  RIVET_DECLARE_PLUGIN(ATLAS_2017_I1589844);

}