rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2018_I1667854

Differential cross section of Z boson production in association with jets in proton-proton collisions at $\sqrt{s} = 13\,$TeV
Experiment: CMS (LHC)
Inspire ID: 1667854
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Philippe Gras
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • Run MC generators with Z decaying leptonically into both electrons and muons at 13 TeV CoM energy. If only one of the two decay channels is included, set LMODE accordingly.

The production of a Z boson, decaying to two charged leptons, in association with jets in proton-proton collisions at a centre-of-mass energy of 13 TeV is measured. Data recorded with the CMS detector at the LHC are used that correspond to an integrated luminosity of 2.19$\,$fb$^{-1}$. The cross section is measured as a function of the jet multiplicity and its dependence on the transverse momentum of the Z boson, the jet kinematic variables (transverse momentum and rapidity), the scalar sum of the jet momenta, which quantifies the hadronic activity, and the balance in transverse momentum between the reconstructed jet recoil and the Z boson. The measurements are compared with predictions from four different calculations. The first two merge matrix elements with different parton multiplicities in the final state and parton showering, one of which includes one-loop corrections. The third is a fixed-order calculation with next-to-next-to-leading order accuracy for the process with a Z boson and one parton in the final state. The fourth combines the fully differential next-to-next-to-leading order calculation with next-to-next-to-leading logarithm resummation and parton showering.

Source code: CMS_2018_I1667854.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FastJets.hh"
  4#include "Rivet/Projections/VisibleFinalState.hh"
  5#include "Rivet/Projections/FinalState.hh"
  6#include "Rivet/Projections/VetoedFinalState.hh"
  7#include "Rivet/Projections/LeptonFinder.hh"
  8
  9namespace Rivet {
 10
 11
 12  /// @brief Differential cross section of Z boson production in association with jets at 13 TeV
 13  ///
 14  /// @note Code copied from CMS_2015_I1410737 and adapted by P. Gras to CMS-SMP-16-015,
 15  class CMS_2018_I1667854 : public Analysis {
 16  public:
 17
 18    /// Constructor
 19    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2018_I1667854);
 20
 21
 22    /// Book histograms and initialise projections before the run
 23    void init() {
 24
 25      // Get options from the new option system; defaults to combined e+mu
 26      _mode = 2;
 27      if ( getOption("LMODE") == "EL" ) _mode = 0;
 28      if ( getOption("LMODE") == "MU" ) _mode = 1;
 29      if ( getOption("LMODE") == "EMU" ) _mode = 2;
 30
 31      // Projections
 32      FinalState fs;
 33      VisibleFinalState visfs(fs);
 34      VetoedFinalState fs_notaudecay(fs);
 35      fs_notaudecay.addDecayProductsVeto(PID::TAU);
 36      fs_notaudecay.addDecayProductsVeto(-PID::TAU);
 37
 38      IdentifiedFinalState bareMuons(fs_notaudecay);
 39      bareMuons.acceptIdPair(PID::MUON);
 40      declare(LeptonFinder(bareMuons, fs, /*dRmax = */0.1,
 41                           Cuts::abseta < 2.4 && Cuts::pT > 20*GeV), "muons");
 42
 43      IdentifiedFinalState bareElectrons(fs_notaudecay);
 44      bareElectrons.acceptIdPair(PID::ELECTRON);
 45      declare(LeptonFinder(bareElectrons, fs, /*dRmax =*/ 0.1,
 46                           Cuts::abseta < 2.4 && Cuts::pT > 20*GeV), "electrons");
 47
 48      FastJets jets(visfs, JetAlg::ANTIKT, 0.4);
 49      declare(jets, "jets");
 50
 51      // Histograms
 52      book(_h_excmult_jets_tot, 1, 1, 1);
 53      book(_h_incmult_jets_tot, 2, 1, 1);
 54      book(_h_zpt1, 3, 1, 1);
 55      book(_h_leading_jet_pt_tot, 4, 1, 1);
 56      book(_h_second_jet_pt_tot, 5, 1, 1);
 57      book(_h_third_jet_pt_tot, 6, 1, 1);
 58      book(_h_leading_jet_y_tot, 7, 1, 1);
 59      book(_h_second_jet_y_tot, 8, 1, 1);
 60      book(_h_third_jet_y_tot, 9, 1, 1);
 61      book(_h_ht1_tot, 10, 1, 1);
 62      book(_h_ht2_tot, 11, 1, 1);
 63      book(_h_ht3_tot, 12, 1, 1);
 64      book(_h_ptbal1, 13, 1, 1);
 65      book(_h_ptbal2, 14, 1, 1);
 66      book(_h_ptbal3, 15, 1, 1);
 67      book(_h_jzb, 16, 1, 1);
 68      book(_h_jzb_ptHigh, 17, 1, 1);
 69      book(_h_jzb_ptLow, 18, 1, 1);
 70    }
 71
 72
 73    /// @brief Z boson finder
 74    ///
 75    /// @note We don't use the standard DileptonFinder class in order to stick to
 76    /// the definition of the publication that is simpler than the DileptonFinder algorithm.
 77    ///
 78    /// @param leptons pt-ordered list of electrons or muons from which to build the Z boson
 79    std::unique_ptr<Particle> zfinder(const Particles& leptons) {
 80      if (leptons.size() < 2) return 0;
 81      if (leptons[0].charge()*leptons[1].charge() > 0) return 0;
 82      std::unique_ptr<Particle> cand(new Particle(PID::ZBOSON, leptons[0].mom() + leptons[1].mom()));
 83      if (cand->mass() < 71.*GeV || cand->mass() > 111.*GeV) return 0;
 84      return cand;
 85    }
 86
 87
 88    /// Perform the per-event analysis
 89    void analyze(const Event& event) {
 90
 91      // Get leptons
 92      const Particles& muons = apply<LeptonFinder>(event, "muons").particlesByPt();
 93      const Particles& electrons = apply<LeptonFinder>(event, "electrons").particlesByPt();
 94
 95      // Look for Z->ee
 96      std::unique_ptr<Particle> z = zfinder(electrons);
 97      const Particles* dressedLeptons = 0;
 98      if (z.get() != nullptr) {
 99        dressedLeptons = &electrons;
100        if (_mode == 1)
101          vetoEvent;
102      } else { // look for Z->mumu
103        z = zfinder(muons);
104        if (z.get() != nullptr) {
105          dressedLeptons = &muons;
106          if (_mode == 0)
107            vetoEvent;
108        } else { // no Z boson found
109          vetoEvent;
110        }
111      }
112
113      // Cluster jets
114      const FastJets& fj = apply<FastJets>(event, "jets");
115      const Jets& jets = fj.jetsByPt(Cuts::absrap < 2.4 && Cuts::pT > 30*GeV);
116
117      // Remove jets overlapping with any of the two selected leptons
118      Jets goodjets = discard(jets, [dressedLeptons](const ParticleBase& j) {
119          return deltaR(j, (*dressedLeptons)[0]) < 0.4 ||  deltaR(j, (*dressedLeptons)[1]) < 0.4;
120        });
121
122      // Compute jet pt scalar sum, H_T:
123      double ht = sum(goodjets, [](const ParticleBase& j){return j.pT();}, 0.);
124
125      // Fill jet number integral histograms
126      _h_excmult_jets_tot->fill(goodjets.size());
127      /// @todo Could be better computed by toIntegral transform on exclusive histo
128      for (size_t iJet = 0; iJet <= goodjets.size(); iJet++ )
129        _h_incmult_jets_tot->fill(iJet);
130
131      if (goodjets.size() < 1) return;
132
133      // Hadronic recoil:
134      FourMomentum recoil;
135      for (const auto& j: goodjets) {
136        recoil += j.momentum();
137      }
138
139      // Jet-Z balance = |recoil_T| - |pt(Z)|
140      double jzb  = recoil.pT() - z->pT();
141
142      // pT balance:
143      double ptbal = (recoil + z->momentum()).pT();
144
145      // Fill leading-jet histograms
146      _h_zpt1->fill(z->pT());
147      const Jet& j1 = goodjets[0];
148      _h_leading_jet_pt_tot->fill(j1.pT()/GeV);
149      _h_leading_jet_y_tot->fill(j1.absrapidity());
150      _h_ht1_tot->fill(ht/GeV);
151      _h_jzb->fill(jzb/GeV);
152      if (z->pT() > 50*GeV) {
153        _h_jzb_ptHigh->fill(jzb/GeV);
154      } else {
155        _h_jzb_ptLow->fill(jzb/GeV);
156      }
157      _h_ptbal1->fill(ptbal/GeV);
158
159      // Fill 2nd jet histograms
160      if (goodjets.size() < 2) return;
161      const Jet& j2 = goodjets[1];
162      _h_second_jet_pt_tot->fill(j2.pT()/GeV);
163      _h_second_jet_y_tot->fill(j2.absrapidity());
164      _h_ht2_tot->fill(ht/GeV);
165      _h_ptbal2->fill(ptbal/GeV);
166
167      // Fill 3rd jet histograms
168      if (goodjets.size() < 3) return;
169      const Jet& j3 = goodjets[2];
170      _h_third_jet_pt_tot->fill(j3.pT()/GeV);
171      _h_third_jet_y_tot->fill(j3.absrapidity());
172      _h_ht3_tot->fill(ht/GeV);
173      _h_ptbal3->fill(ptbal/GeV);
174    }
175
176
177    /// Normalise histograms etc., after the run
178    void finalize() {
179
180      // Normalisation factor
181      double norm = (sumOfWeights() != 0) ? crossSection()/picobarn/sumOfWeights() : 1.0;
182      // When running in combined mode, need to average to get lepton xsec
183      if (_mode == 2) norm /= 2.;
184
185      // MSG_INFO("Cross section = " << std::setfill(' ') << std::setw(14)
186      //          << std::fixed << std::setprecision(3) << crossSection()/picobarn << " pb");
187      // MSG_INFO("# Events      = " << std::setfill(' ') << std::setw(14)
188      //          << std::fixed << std::setprecision(3) << numEvents() );
189      // MSG_INFO("SumW          = " << std::setfill(' ') << std::setw(14)
190      //          << std::fixed << std::setprecision(3) << sumOfWeights());
191      // MSG_INFO("Norm factor   = " << std::setfill(' ') << std::setw(14)
192      //          << std::fixed << std::setprecision(6) << norm);
193
194      scale(_h_excmult_jets_tot, norm);
195      scale(_h_incmult_jets_tot, norm);
196      scale(_h_zpt1, norm);
197      scale(_h_leading_jet_pt_tot, norm);
198      scale(_h_second_jet_pt_tot, norm);
199      scale(_h_third_jet_pt_tot, norm);
200      scale(_h_leading_jet_y_tot, norm);
201      scale(_h_second_jet_y_tot, norm);
202      scale(_h_third_jet_y_tot, norm);
203      scale(_h_ht1_tot, norm);
204      scale(_h_ht2_tot, norm);
205      scale(_h_ht3_tot, norm);
206      scale(_h_ptbal1, norm);
207      scale(_h_ptbal2, norm);
208      scale(_h_ptbal3, norm);
209      scale(_h_jzb, norm);
210      scale(_h_jzb_ptHigh, norm);
211      scale(_h_jzb_ptLow, norm);
212    }
213
214
215  protected:
216
217    size_t _mode;
218
219
220  private:
221
222    /// @name Histograms
223    /// @{
224    Histo1DPtr _h_excmult_jets_tot,  _h_incmult_jets_tot;
225    Histo1DPtr _h_leading_jet_pt_tot, _h_second_jet_pt_tot;
226    Histo1DPtr _h_third_jet_pt_tot, _h_fourth_jet_pt_tot;
227    Histo1DPtr _h_leading_jet_y_tot, _h_second_jet_y_tot;
228    Histo1DPtr _h_third_jet_y_tot, _h_fourth_jet_y_tot;
229    Histo1DPtr _h_ht1_tot, _h_ht2_tot, _h_ht3_tot, _h_ht4_tot;
230    Histo1DPtr _h_ptbal1, _h_ptbal2, _h_ptbal3;
231    Histo1DPtr _h_jzb, _h_jzb_ptHigh, _h_jzb_ptLow;
232    Histo1DPtr _h_zpt1;
233    /// @}
234
235  };
236
237
238
239  RIVET_DECLARE_PLUGIN(CMS_2018_I1667854);
240
241}