rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2022_I2079374

Measurement of the mass dependence of the transverse momentum of lepton pairs in Drell--Yan production in proton-proton collisions at 13 TeV
Experiment: CMS (LHC)
Inspire ID: 2079374
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Louis Moureaux
  • Buğra Bilin
  • Itana Bubanja
  • 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 double differential cross sections of the Drell-Yan lepton pair ($\ell^+\ell^-$, dielectron or dimuon) production are measured, as functions of the invariant mass $m_{\ell\ell}$, transverse momentum $p_T^{\ell\ell}$, and $\varphi^*_\eta$. The $\varphi^*_\eta$ observable, derived from angular measurements of the leptons and highly correlated with $p_T^{\ell\ell}$, is used to probe the low-$p_T^{\ell\ell}$ region in a complementary way. Dilepton masses up to 1 TeV are investigated. Additionally, a measurement is performed requiring at least one jet in the final state. To benefit from partial cancellation of the systematic uncertainty, the ratios of the differential cross sections for various $m_{\ell\ell}$ ranges to those in the Z mass peak interval are presented. The collected data correspond to an integrated luminosity of 36.3 fb$^{-1}$ of proton-proton collisions recorded with the CMS detector at the LHC at a centre-of-mass energy of 13 TeV. Measurements are compared with predictions based on perturbative quantum chromodynamics, including soft-gluon resummation.

Source code: CMS_2022_I2079374.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/LeptonFinder.hh"
  7#include "Rivet/Projections/PromptFinalState.hh"
  8
  9namespace Rivet {
 10
 11
 12  /// @brief Z pT over a wide mass range
 13  class CMS_2022_I2079374 : public Analysis {
 14  public:
 15
 16      /// Constructor
 17      RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2022_I2079374);
 18
 19      /// Book histograms and initialise projections before the run
 20      void init() {
 21        // Get options from the new option system
 22        // default to combined.
 23        _mode = 2;
 24        if (getOption("LMODE") == "EL") {
 25          _mode = 0;
 26        } else if (getOption("LMODE") == "MU") {
 27          _mode = 1;
 28        } else if (getOption("LMODE") == "EMU") {
 29          _mode = 2;
 30        } else {
 31          MSG_INFO("Assuming a mixed electron/muon sample. Set LMODE=EL/MU if your sample contains a single flavor.");
 32        }
 33
 34        FinalState fs;
 35        PromptFinalState pfs(fs);
 36        PromptFinalState bareMuons(Cuts::abspid == PID::MUON);
 37        declare(LeptonFinder(bareMuons, pfs, 0.1, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV), "muons");
 38
 39        PromptFinalState bareElectrons(Cuts::abspid == PID::ELECTRON);
 40        declare(LeptonFinder(bareElectrons, pfs, 0.1, Cuts::abseta < 2.4 && Cuts::pT > 20*GeV), "electrons");
 41
 42        FastJets jets(fs, JetAlg::ANTIKT, 0.4);
 43        declare(jets, "jets");
 44
 45        // Histograms
 46        book(_pt_0jet_mass50_76, "_0jetmass_111", refData(1, 1, 1));
 47        book(_pt_0jet_mass76_106, "_0jetmass_311", refData(3, 1, 1));
 48        book(_pt_0jet_mass106_170, "_0jetmass_511", refData(5, 1, 1));
 49        book(_pt_0jet_mass170_350, "_0jetmass_711", refData(7, 1, 1));
 50        book(_pt_0jet_mass350_1000, "_0jetmass_911", refData(9, 1, 1));
 51
 52        book(_pt_1jet_mass50_76, "_1jetmass_1111", refData(11, 1, 1));
 53        book(_pt_1jet_mass76_106, "_1jetmass_1311", refData(13, 1, 1));
 54        book(_pt_1jet_mass106_170, "_1jetmass_1511", refData(15, 1, 1));
 55        book(_pt_1jet_mass170_350, "_1jetmass_1711", refData(17, 1, 1));
 56
 57        book(_phistar_mass50_76, "_phistarmass_1911", refData(19, 1, 1));
 58        book(_phistar_mass76_106, "_phistarmass_2111", refData(21, 1, 1));
 59        book(_phistar_mass106_170, "_phistarmass_2311", refData(23, 1, 1));
 60        book(_phistar_mass170_350, "_phistarmass_2511", refData(25, 1, 1));
 61        book(_phistar_mass350_1000, "_phistarmass_2711", refData(27, 1, 1));
 62
 63      }
 64
 65      /// Z boson finder.
 66      /// Note: we don't use the standard DileptonFinder class in order to stick to
 67      /// the definition of the publication that is simpler than the DileptonFinder
 68      /// algorithm
 69      /// @param leptons pt-ordered of electron or muon collection to use to build
 70      /// the Z boson
 71      std::unique_ptr<Particle> zfinder(const DressedLeptons &leptons)
 72      {
 73        if (leptons.size() < 2) {
 74          return nullptr;
 75        }
 76        // Leading lepton pT cut
 77        if (leptons[0].pT() < 25.*GeV) {
 78          return nullptr;
 79        }
 80        if (leptons[0].charge() * leptons[1].charge() > 0) {
 81          return nullptr;
 82        }
 83
 84        std::unique_ptr<Particle> candidate(new Particle(PID::ZBOSON, leptons[0].mom() + leptons[1].mom()));
 85        if (candidate->mass() < 50*GeV || candidate->mass() > 1000*GeV) {
 86          return nullptr;
 87        }
 88        return candidate;
 89      }
 90
 91      /// Perform the per-event analysis
 92      void analyze(const Event& event)
 93      {
 94        // Fetch dressed leptons
 95        auto muons = apply<LeptonFinder>(event, "muons").dressedLeptons();
 96        auto electrons = apply<LeptonFinder>(event, "electrons").dressedLeptons();
 97
 98        const DressedLeptons* dressedLeptons = nullptr;
 99
100        //Look for Z->ee
101        std::unique_ptr<Particle> z = zfinder(electrons);
102        if (z != nullptr) {
103          dressedLeptons = &electrons;
104          if (_mode == 1) {
105            // Z->ee found, but running with LMODE=MU
106            vetoEvent;
107          }
108        } else if (_mode == 0) {
109          // LMODE=EL, but no Z->ee found
110          vetoEvent;
111        } else { //look for Z->mumu
112          z = zfinder(muons);
113          if (z.get() != nullptr) {
114            dressedLeptons = &muons;
115          } else { //no Z boson found
116            vetoEvent;
117          }
118        }
119
120        // Cluster jets
121        Jets jets = apply<FastJets>(event, "jets").jetsByPt(Cuts::absrap < 2.4 && Cuts::pT > 30*GeV);
122        idiscardIfAnyDeltaRLess(jets, *dressedLeptons, 0.4);
123
124        //------------------------------------------------
125
126        // Compute phi*
127        const Particle& lminus = dressedLeptons->at(0).charge() < 0 ? dressedLeptons->at(0) : dressedLeptons->at(1);
128        const Particle& lplus  = dressedLeptons->at(0).charge() < 0 ? dressedLeptons->at(1) : dressedLeptons->at(0);
129        double phi_acop = M_PI - deltaPhi(lminus, lplus);
130        double costhetastar = tanh( 0.5 * (lminus.eta() - lplus.eta()) );
131        double sin2thetastar = (costhetastar > 1) ? 0.0 : (1.0 - sqr(costhetastar));
132        double phistar = tan(0.5 * phi_acop) * sqrt(sin2thetastar);
133
134        if (z->mass() > 50. && z->mass() < 76.) {
135          _pt_0jet_mass50_76->fill(z->pT());
136          _phistar_mass50_76->fill(phistar);
137        } else if (z->mass() > 76. && z->mass() < 106.) {
138          _pt_0jet_mass76_106->fill(z->pT());
139          _phistar_mass76_106->fill(phistar);
140        } else if (z->mass() > 106. && z->mass() < 170.) {
141          _pt_0jet_mass106_170->fill(z->pT());
142          _phistar_mass106_170->fill(phistar);
143        } else if (z->mass() > 170. && z->mass() < 350.) {
144          _pt_0jet_mass170_350->fill(z->pT());
145          _phistar_mass170_350->fill(phistar);
146        } else if (z->mass() > 350. && z->mass() < 1000.) {
147          _pt_0jet_mass350_1000->fill(z->pT());
148          _phistar_mass350_1000->fill(phistar);
149        }
150
151        //---------------------- at least one jet --------------------------------
152        if (jets.size() < 1) {
153          return;
154        }
155
156        if (z->mass() > 50. && z->mass() < 76.) {
157          _pt_1jet_mass50_76->fill(z->pT());
158        } else if (z->mass() > 76. && z->mass() < 106.) {
159          _pt_1jet_mass76_106->fill(z->pT());
160        } else if (z->mass() > 106. && z->mass() < 170.) {
161          _pt_1jet_mass106_170->fill(z->pT());
162        } else if (z->mass() > 170. && z->mass() < 350.) {
163          _pt_1jet_mass170_350->fill(z->pT());
164        }
165      }
166
167      /// Calculates the ratio between two histograms
168      void calculateRatio(int d, const Histo1DPtr &numerator, const Histo1DPtr &denominator) {
169        Estimate1DPtr ratio;
170        book(ratio, d, 1, 1);
171
172        // The denominator has a finer binning, so rebin it
173        auto rebinned_den = denominator->clone();
174        rebinned_den.rebinX(numerator->xEdges());
175
176        divide(*numerator, rebinned_den, ratio);
177      }
178
179      /// Normalise histograms etc., after the run
180      void finalize() {
181        double norm = (sumOfWeights() != 0) ? crossSection()/picobarn/sumOfWeights() : 1.0;
182
183        if (_mode == 2)  {
184          norm /= 2.;
185        }
186
187        scale(_pt_0jet_mass50_76, norm);
188        scale(_pt_0jet_mass76_106, norm);
189        scale(_pt_0jet_mass106_170, norm);
190        scale(_pt_0jet_mass170_350, norm);
191        scale(_pt_0jet_mass350_1000, norm);
192
193        scale(_pt_1jet_mass50_76, norm);
194        scale(_pt_1jet_mass76_106, norm);
195        scale(_pt_1jet_mass106_170, norm);
196        scale(_pt_1jet_mass170_350, norm);
197
198        scale(_phistar_mass50_76, norm);
199        scale(_phistar_mass76_106, norm);
200        scale(_phistar_mass106_170, norm);
201        scale(_phistar_mass170_350, norm);
202        scale(_phistar_mass350_1000, norm);
203
204        // Calculate the ratios
205        calculateRatio(29, _pt_0jet_mass50_76, _pt_0jet_mass76_106);
206        calculateRatio(31, _pt_0jet_mass106_170, _pt_0jet_mass76_106);
207        calculateRatio(33, _pt_0jet_mass170_350, _pt_0jet_mass76_106);
208        calculateRatio(35, _pt_0jet_mass350_1000, _pt_0jet_mass76_106);
209
210        calculateRatio(37, _pt_1jet_mass50_76, _pt_1jet_mass76_106);
211        calculateRatio(39, _pt_1jet_mass106_170, _pt_1jet_mass76_106);
212        calculateRatio(41, _pt_1jet_mass170_350, _pt_1jet_mass76_106);
213
214        calculateRatio(43, _phistar_mass50_76, _phistar_mass76_106);
215        calculateRatio(45, _phistar_mass106_170, _phistar_mass76_106);
216        calculateRatio(47, _phistar_mass170_350, _phistar_mass76_106);
217        calculateRatio(49, _phistar_mass350_1000, _phistar_mass76_106);
218      }
219
220    protected:
221      size_t _mode;
222
223    private:
224      /// @name Histograms
225      Histo1DPtr _pt_0jet_mass50_76,
226                 _pt_0jet_mass76_106,
227                 _pt_0jet_mass106_170,
228                 _pt_0jet_mass170_350,
229                 _pt_0jet_mass350_1000;
230      Histo1DPtr _pt_1jet_mass50_76,
231                 _pt_1jet_mass76_106,
232                 _pt_1jet_mass106_170,
233                 _pt_1jet_mass170_350;
234      Histo1DPtr _phistar_mass50_76,
235                 _phistar_mass76_106,
236                 _phistar_mass106_170,
237                 _phistar_mass170_350,
238                 _phistar_mass350_1000;
239    };
240
241  RIVET_DECLARE_PLUGIN(CMS_2022_I2079374);
242}