rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2019_I1753680

Measurements of differential Z boson production cross sections in proton-proton collisions at 13 TeV
Experiment: CMS (LHC)
Inspire ID: 1753680
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Jay Lawhorn
  • Markus Seidel
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • pp to Z interactions at $\sqrt{s} = 13$ TeV. Data collected by CMS during the year 2016.

Measurements are presented of the differential cross sections for Z bosons produced in proton-proton collisions at $\sqrt{s} = 13 \text{TeV}$ and decaying to muons and electrons. The data analyzed were collected in 2016 with the CMS detector at the LHC and correspond to an integrated luminosity of $35.9/fb$. The measured fiducial inclusive product of cross section and branching fraction agrees with next-to-next-to-leading order quantum chromodynamics calculations. Differential cross sections of the transverse momentum \pt, the optimized angular variable $\phi^{*}_{\eta}$, and the rapidity of lepton pairs are measured. The data are corrected for detector effects and compared to theoretical predictions using fixed order, resummed, and parton shower calculations. The uncertainties of the measured normalized cross sections are smaller than 0.5% for $\phi^{*}_{\eta} < 0.5$ and for $p_{T}^{Z} < 50 \text{GeV}$.

Source code: CMS_2019_I1753680.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/FastJets.hh"
  5#include "Rivet/Projections/LeptonFinder.hh"
  6#include "Rivet/Projections/DileptonFinder.hh"
  7
  8namespace Rivet {
  9
 10
 11  /// @brief Measurements of differential Z boson production cross sections in proton-proton collisions at 13 TeV
 12  class CMS_2019_I1753680 : public Analysis {
 13  public:
 14
 15    /// Constructor
 16    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2019_I1753680);
 17
 18
 19    /// @name Analysis methods
 20    /// @{
 21
 22    /// Book histograms and initialise projections before the run
 23    void init() {
 24
 25      // Get options from the new option system
 26      // default to combined.
 27      _mode = 2;
 28      if ( getOption("LMODE") == "EL" ) _mode = 0;
 29      if ( getOption("LMODE") == "MU" ) _mode = 1;
 30      if ( getOption("LMODE") == "EMU" ) _mode = 2;
 31
 32      // Initialise and register projections
 33      Cut cut = Cuts::abseta < 2.4 && Cuts::pT > 25*GeV;
 34      DileptonFinder zeeFind(91.2*GeV, 0.1, cut && Cuts::abspid == PID::ELECTRON, Cuts::massIn(76.2*GeV, 106.2*GeV));
 35      declare(zeeFind, "ZeeFind");
 36      DileptonFinder zmmFind(91.2*GeV, 0.1, cut && Cuts::abspid == PID::MUON    , Cuts::massIn(76.2*GeV, 106.2*GeV));
 37      declare(zmmFind, "ZmmFind");
 38
 39      // Book histograms
 40      book(_h_Zmm_absY,    26, 1, 1);
 41      book(_h_Zee_absY,    26, 1, 2);
 42      book(_h["absY"],     26, 1, 3);
 43      book(_h_Zmm_pt,      27, 1, 1);
 44      book(_h_Zee_pt,      27, 1, 2);
 45      book(_h["pt"],       27, 1, 3);
 46      book(_h_Zmm_phiStar, 28, 1, 1);
 47      book(_h_Zee_phiStar, 28, 1, 2);
 48      book(_h["phiStar"],  28, 1, 3);
 49      book(_h["pt_Y0"],    29, 1, 1);
 50      book(_h["pt_Y1"],    29, 1, 2);
 51      book(_h["pt_Y2"],    29, 1, 3);
 52      book(_h["pt_Y3"],    29, 1, 4);
 53      book(_h["pt_Y4"],    29, 1, 5);
 54
 55      book(_h_norm["pt"],      30, 1, 1);
 56      book(_h_norm["phiStar"], 31, 1, 1);
 57      book(_h_norm["absY"],    32, 1, 1);
 58      book(_h_norm["pt_Y0"],   33, 1, 1);
 59      book(_h_norm["pt_Y1"],   33, 1, 2);
 60      book(_h_norm["pt_Y2"],   33, 1, 3);
 61      book(_h_norm["pt_Y3"],   33, 1, 4);
 62      book(_h_norm["pt_Y4"],   33, 1, 5);
 63
 64    }
 65
 66
 67    /// Perform the per-event analysis
 68    void analyze(const Event& event) {
 69
 70      const DileptonFinder& zeeFS = apply<DileptonFinder>(event, "ZeeFind");
 71      const DileptonFinder& zmumuFS = apply<DileptonFinder>(event, "ZmmFind");
 72
 73      const Particles& zees = zeeFS.bosons();
 74      const Particles& zmumus = zmumuFS.bosons();
 75
 76      if (zees.size() + zmumus.size() != 1) {
 77        MSG_DEBUG("Did not find exactly one good Z candidate");
 78        vetoEvent;
 79      }
 80
 81      //event identification depending on mass window
 82      bool ee_event=false;
 83      bool mm_event=false;
 84
 85      if (zees.size() == 1) {
 86        ee_event = true;
 87      }
 88      if (zmumus.size() == 1) {
 89        mm_event = true;
 90      }
 91
 92      if (ee_event && _mode == 1)
 93        vetoEvent;
 94      if (mm_event && _mode == 0)
 95        vetoEvent;
 96
 97      const Particles& theLeptons = ee_event ? zeeFS.constituents() : zmumuFS.constituents();
 98      const Particle& lminus = theLeptons[0].charge() < 0 ? theLeptons[0] : theLeptons[1];
 99      const Particle& lplus = theLeptons[0].charge() < 0 ? theLeptons[1] : theLeptons[0];
100
101      //calculate phi*
102      const double thetaStar = acos(tanh( 0.5 * (lminus.eta() - lplus.eta()) ));
103      const double dPhi = M_PI - deltaPhi(lminus, lplus);
104      const double phiStar = tan(0.5 * dPhi) * sin(thetaStar);
105
106      const Particle& zcand = ee_event ? zees[0] : zmumus[0];
107
108      if (ee_event) {
109        _h_Zee_absY->fill(zcand.absrap());
110        _h_Zee_pt->fill(zcand.pt()/GeV);
111        _h_Zee_phiStar->fill(phiStar);
112      }
113      else if (mm_event) {
114        _h_Zmm_absY->fill(zcand.absrap());
115        _h_Zmm_pt->fill(zcand.pt()/GeV);
116        _h_Zmm_phiStar->fill(phiStar);
117      }
118
119      _h["pt"]->fill(zcand.pt()/GeV);
120      _h_norm["pt"]->fill(zcand.pt()/GeV);
121      _h["phiStar"]->fill(phiStar);
122      _h_norm["phiStar"]->fill(phiStar);
123      _h["absY"]->fill(zcand.absrap());
124      _h_norm["absY"]->fill(zcand.absrap());
125
126      if (zcand.absrap() < 0.4) {
127        _h["pt_Y0"]->fill(zcand.pt()/GeV);
128        _h_norm["pt_Y0"]->fill(zcand.pt()/GeV);
129      }
130      else if (zcand.absrap() < 0.8) {
131        _h["pt_Y1"]->fill(zcand.pt()/GeV);
132        _h_norm["pt_Y1"]->fill(zcand.pt()/GeV);
133      }
134      else if (zcand.absrap() < 1.2) {
135        _h["pt_Y2"]->fill(zcand.pt()/GeV);
136        _h_norm["pt_Y2"]->fill(zcand.pt()/GeV);
137      }
138      else if (zcand.absrap() < 1.6) {
139        _h["pt_Y3"]->fill(zcand.pt()/GeV);
140        _h_norm["pt_Y3"]->fill(zcand.pt()/GeV);
141      }
142      else if (zcand.absrap() < 2.4) {
143        _h["pt_Y4"]->fill(zcand.pt()/GeV);
144        _h_norm["pt_Y4"]->fill(zcand.pt()/GeV);
145      }
146
147    }
148
149    /// Normalise histograms etc., after the run
150    void finalize() {
151
152      double norm = (sumOfWeights() != 0) ? crossSection()/picobarn/sumOfWeights() : 1.0;
153
154      scale(_h_Zmm_pt,      norm);
155      scale(_h_Zmm_absY,    norm);
156      scale(_h_Zmm_phiStar, norm);
157
158      scale(_h_Zee_pt,      norm);
159      scale(_h_Zee_absY,    norm);
160      scale(_h_Zee_phiStar, norm);
161
162      // when running in combined mode, need to average to get lepton xsec
163      if (_mode == 2) norm /= 2.;
164      scale(_h, norm);
165
166      for (auto& item : _h_norm) {
167        const double rho = item.second->densitySum(false);
168        if (rho)  scale(item.second, 1.0/rho);
169      }
170
171    }
172
173    /// @}
174
175  protected:
176
177    size_t _mode;
178
179    /// @name Histograms
180
181  private:
182
183    Histo1DPtr   _h_Zmm_pt, _h_Zmm_phiStar, _h_Zmm_absY;
184    Histo1DPtr   _h_Zee_pt, _h_Zee_phiStar, _h_Zee_absY;
185
186    map<string,Histo1DPtr> _h;
187    map<string,Histo1DPtr> _h_norm;
188
189  };
190
191
192  RIVET_DECLARE_PLUGIN(CMS_2019_I1753680);
193
194
195}