Rivet Analyses Reference

CMS_2015_I1310737

Jet multiplicity and differential cross-sections of $Z$+jets events in $pp$ at $\sqrt{s} = 7$ TeV
Experiment: CMS (LHC)
Inspire ID: 1310737
Status: VALIDATED
Authors:
  • Fabio Cossutti
  • Chiara La Licata
References:
  • Phys.Rev. D91 (2015) 052008
  • http://dx.doi.org/10.1103/PhysRevD.91.052008
  • arxiv:1408.3104
  • http://inspirehep.net/record/1310737
Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Run MC generators with Z decaying leptonically into both electrons and muons at 7 TeV CoM energy. Order of 5 million unweighted events can give a reasonable global comparison, but precision in the high jet multiplicity region/high jet pt may require substantially larger samples or statistical enhancement of high jet multiplicities.

Measurements of differential cross sections are presented for the production of a Z boson and at least one hadronic jet in proton-proton collisions at $\sqrt{s}=7$ TeV, recorded by the CMS detector, using a data sample corresponding to an integrated luminosity of 4.9 $\text{fb}^{-1}$. The jet multiplicity distribution is measured for up to six jets. The differential cross sections are measured as a function of jet transverse momentum and pseudorapidity for the four highest transverse momentum jets. The distribution of the scalar sum of jet transverse momenta is also measured as a function of the jet multiplicity. The measurements are compared with theoretical predictions at leading and next-to-leading order in perturbative QCD. Cuts: First two leading electrons or muons with $p_T > 20$ GeV and $|\eta| < 2.4$ Dilepton invariant mass in the [71,111] GeV range Jets $p_T > 30$ GeV and $|\eta| < 2.4$ $\Delta{R}($lepton,jet$) > 0.5$

Source code: CMS_2015_I1310737.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/ChargedFinalState.hh"
#include "Rivet/Projections/InvMassFinalState.hh"
#include "Rivet/Projections/VisibleFinalState.hh"
#include "Rivet/Projections/VetoedFinalState.hh"
#include "Rivet/Projections/ZFinder.hh"

namespace Rivet {


  class CMS_2015_I1310737 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2015_I1310737);


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

      FinalState fs; ///< @todo No cuts?
      VisibleFinalState visfs(fs);
      // Prompt leptons only
      PromptFinalState pfs(fs);

      ZFinder zeeFinder(pfs, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV, PID::ELECTRON, 71.0*GeV, 111.0*GeV);
      declare(zeeFinder, "ZeeFinder");

      ZFinder zmumuFinder(pfs, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV, PID::MUON, 71.0*GeV, 111.0*GeV);
      declare(zmumuFinder, "ZmumuFinder");

      VetoedFinalState jetConstits(visfs);
      jetConstits.addVetoOnThisFinalState(zeeFinder);
      jetConstits.addVetoOnThisFinalState(zmumuFinder);

      FastJets akt05Jets(jetConstits, FastJets::ANTIKT, 0.5);
      declare(akt05Jets, "AntiKt05Jets");


      book(_h_excmult_jets_tot ,1, 1, 1);
      book(_h_incmult_jets_tot ,2, 1, 1);
      book(_h_leading_jet_pt_tot ,3, 1, 1);
      book(_h_second_jet_pt_tot ,4, 1, 1);
      book(_h_third_jet_pt_tot ,5, 1, 1);
      book(_h_fourth_jet_pt_tot ,6, 1, 1);
      book(_h_leading_jet_eta_tot ,7, 1, 1);
      book(_h_second_jet_eta_tot ,8, 1, 1);
      book(_h_third_jet_eta_tot ,9, 1, 1);
      book(_h_fourth_jet_eta_tot ,10, 1, 1);
      book(_h_ht1_tot ,11, 1, 1);
      book(_h_ht2_tot ,12, 1, 1);
      book(_h_ht3_tot ,13, 1, 1);
      book(_h_ht4_tot ,14, 1, 1);
    }


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

      const ZFinder& zeeFS = apply<ZFinder>(event, "ZeeFinder");
      const ZFinder& zmumuFS = apply<ZFinder>(event, "ZmumuFinder");

      const Particles& zees = zeeFS.bosons();
      const Particles& zmumus = zmumuFS.bosons();

      // We did not find exactly one Z. No good.
      if (zees.size() + zmumus.size() != 1) {
        MSG_DEBUG("Did not find exactly one good Z candidate");
        vetoEvent;
      }

      // Find the (dressed!) leptons
      const Particles& dressedLeptons = zees.size() ? zeeFS.constituents() : zmumuFS.constituents();

      // Cluster jets
      // NB. Veto has already been applied on leptons and photons used for dressing
      const FastJets& fj = apply<FastJets>(event, "AntiKt05Jets");
      const Jets& jets = fj.jetsByPt(Cuts::abseta < 2.4 && Cuts::pT > 30*GeV);

      // Perform lepton-jet overlap and HT calculation
      double ht = 0;
      Jets goodjets;
      for (const Jet& j : jets) {
        // Decide if this jet is "good", i.e. isolated from the leptons
        /// @todo Nice use-case for any() and a C++11 lambda
        bool overlap = false;
        for (const Particle& l : dressedLeptons) {
          if (deltaR(j, l) < 0.5) {
            overlap = true;
            break;
          }
        }

        // Fill HT and good-jets collection
        if (overlap) continue;
        goodjets.push_back(j);
        ht += j.pT();
      }

      // We don't care about events with no isolated jets
      if (goodjets.empty()) {
        MSG_DEBUG("No jets in event");
        vetoEvent;
      }


      /////////////////


      // Weight to be used for histo filling
      const double w = 0.5 * 1.0;

      // Fill jet number integral histograms
      _h_excmult_jets_tot->fill(goodjets.size(), w);
      /// @todo Could be better computed by toIntegral transform on exclusive histo
      for (size_t iJet = 1; iJet <= goodjets.size(); iJet++ )
        _h_incmult_jets_tot->fill(iJet, w);

      // Fill leading jet histograms
      const Jet& j1 = goodjets[0];
      _h_leading_jet_pt_tot->fill(j1.pT()/GeV, w);
      _h_leading_jet_eta_tot->fill(j1.abseta(), w);
      _h_ht1_tot->fill(ht/GeV, w);

      // Fill 2nd jet histograms
      if (goodjets.size() < 2) return;
      const Jet& j2 = goodjets[1];
      _h_second_jet_pt_tot->fill(j2.pT()/GeV, w);
      _h_second_jet_eta_tot->fill(j2.abseta(), w);
      _h_ht2_tot->fill(ht/GeV, w);

      // Fill 3rd jet histograms
      if (goodjets.size() < 3) return;
      const Jet& j3 = goodjets[2];
      _h_third_jet_pt_tot->fill(j3.pT()/GeV, w);
      _h_third_jet_eta_tot->fill(j3.abseta(), w);
      _h_ht3_tot->fill(ht/GeV, w);

      // Fill 4th jet histograms
      if (goodjets.size() < 4) return;
      const Jet& j4 = goodjets[3];
      _h_fourth_jet_pt_tot->fill(j4.pT()/GeV, w);
      _h_fourth_jet_eta_tot->fill(j4.abseta(), w);
      _h_ht4_tot->fill(ht/GeV, w);
    }


    /// Normalise histograms etc., after the run
    void finalize() {

      const double norm = (sumOfWeights() != 0) ? crossSection()/sumOfWeights() : 1.0;

      scale(_h_excmult_jets_tot, norm );
      scale(_h_incmult_jets_tot, norm );
      scale(_h_leading_jet_pt_tot, norm );
      scale(_h_second_jet_pt_tot, norm );
      scale(_h_third_jet_pt_tot, norm );
      scale(_h_fourth_jet_pt_tot, norm );
      scale(_h_leading_jet_eta_tot, norm );
      scale(_h_second_jet_eta_tot, norm );
      scale(_h_third_jet_eta_tot, norm );
      scale(_h_fourth_jet_eta_tot, norm );
      scale(_h_ht1_tot, norm );
      scale(_h_ht2_tot, norm );
      scale(_h_ht3_tot, norm );
      scale(_h_ht4_tot, norm );
    }


  private:

    /// @name Histograms

    Histo1DPtr _h_excmult_jets_tot,  _h_incmult_jets_tot;
    Histo1DPtr _h_leading_jet_pt_tot, _h_second_jet_pt_tot, _h_third_jet_pt_tot, _h_fourth_jet_pt_tot;
    Histo1DPtr _h_leading_jet_eta_tot, _h_second_jet_eta_tot, _h_third_jet_eta_tot, _h_fourth_jet_eta_tot;
    Histo1DPtr _h_ht1_tot, _h_ht2_tot, _h_ht3_tot, _h_ht4_tot;

  };


  RIVET_DECLARE_PLUGIN(CMS_2015_I1310737);


}