rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2019_I1705068

Measurement of associated production of a W boson and a charm quark in proton-proton collisions at 13 TeV
Experiment: CMS (LHC)
Inspire ID: 1705068
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Svenja Pflitsch
  • Katerina Lipka
References: Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
  • 13 TeV $pp \to W+c$.

Measurements are presented of associated production of a W boson and a charm quark (W+c) in proton-proton collisions at a center-of-mass energy of 13 TeV. The data correspond to an integrated luminosity of $35.7 fb^{-1}$ collected by the CMS experiment at the CERN LHC. The W bosons are identified by their decay into a muon and a neutrino. The charm quarks are tagged via the full reconstruction of $D*(2010)^{\pm}$ mesons that decay via $D*(2010)^{\pm} \to D^{0} + \pi^{pm} \to K^{\mp} + \pi^{pm} + \pi^{pm}$. A parton-level cross section is measured in the fiducial region defined by the muon transverse momentum $\p_{T}^{\mu} > 26 GeV$, muon pseudorapidity $\abs{\eta^{\mu}} < 2.4$, and charm quark transverse momentum $p_{T}^{c} > 5\text{GeV}$. The cross section is also measured differentially as a function of the pseudorapidity of the muon from the W boson decay. A particle-level cross section is measured in the fiducial region defined by the muon transverse momentum $\p_{T}^{\mu} > 26 GeV$, muon pseudorapidity $\abs{\eta^{\mu}} < 2.4$, D*(2010)^{\pm} transverse momentum $p_{T}^{D*} > 5\text{GeV}$ and D*(2010)^{\pm} pseudorapidity $\abs{\eta^{D*}} < 2.4$. The particle-level cross section has been implemented in the Rivet plugin.

Source code: CMS_2019_I1705068.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/FastJets.hh"
  5#include "Rivet/Projections/UnstableParticles.hh"
  6#include "Rivet/Projections/LeptonFinder.hh"
  7#include "Rivet/Projections/MissingMomentum.hh"
  8
  9namespace Rivet {
 10
 11
 12  /// W + charm quark in proton-proton collisions at 13 TeV
 13  class CMS_2019_I1705068 : public Analysis {
 14  public:
 15
 16    /// Constructor
 17    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2019_I1705068);
 18
 19
 20    /// Initialise analyis
 21    void init() {
 22
 23      // Projections
 24      declare(MissingMom(), "MET");
 25      LeptonFinder lf(0.1, Cuts::abseta < 2.4 && Cuts::abspid == PID::MUON);
 26      declare(lf, "Leptons");
 27
 28      UnstableParticles dst(Cuts::pT > 5*GeV && Cuts::abseta < 2.4);
 29      declare(dst, "Dstar");
 30
 31      // Particle-Level histograms from the paper:
 32      book(_hist_WplusMinus_MuAbseta, "d04-x01-y01");
 33      book(_hist_Wplus_MuAbseta, "d05-x01-y01");
 34      book(_hist_Wminus_MuAbseta, "d06-x01-y01");
 35    }
 36
 37
 38    void analyze(const Event& event) {
 39
 40      // Identify the closest-matching mu+MET to m == mW
 41      const P4& pmiss = apply<MissingMom>(event, "MET").missingMom();
 42      const Particles& mus = apply<LeptonFinder>(event, "Leptons").particles();
 43      const Particles mus_mtfilt = select(mus, [&](const Particle& m){ return mT(m, pmiss) > 0*GeV; });
 44      const int imfound = closestMatchIndex(mus_mtfilt, pmiss, Kin::mass, 80.4*GeV);
 45
 46      // No MET or MT cut at generator level
 47      if (imfound < 0) vetoEvent;
 48      const Particle& lepton0 = mus_mtfilt[imfound];
 49      const double pt0 = lepton0.pT();
 50      const double eta0 = lepton0.abseta();
 51      const int muID = lepton0.pid();
 52      if (eta0 > 2.4 || pt0 < 26*GeV) vetoEvent;
 53
 54      // D* selection:
 55      // OS = W boson and D* Meson have Opposite (charge) Signs
 56      // SS = W Boson and D* Meson have Same (charge) Signs
 57      // Associated W+c only has OS contributions,
 58      // W+ccbar (ccbar from gluon splitting) has equal probability to be OS or SS
 59      // OS-SS to remove the gluon splitting background
 60
 61      const UnstableParticles& dst = apply<UnstableParticles>(event, "Dstar");
 62      for (const Particle& p : dst.particles()) {
 63        if (muID == -13 && p.pid() == -413) { // OS
 64          _hist_Wplus_MuAbseta->fill(eta0);
 65          _hist_WplusMinus_MuAbseta->fill(eta0);
 66        }
 67        else if (muID == 13 && p.pid() == 413) { // OS
 68          _hist_Wminus_MuAbseta->fill(eta0);
 69          _hist_WplusMinus_MuAbseta->fill(eta0);
 70        }
 71        else if (muID == -13 && p.pid() == 413) { // SS
 72          _hist_Wplus_MuAbseta->fill(eta0*-1);
 73          _hist_WplusMinus_MuAbseta->fill(eta0*-1);
 74        }
 75        else if (muID == 13 && p.pid() == -413) { // SS
 76          _hist_Wminus_MuAbseta->fill(eta0*-1);
 77          _hist_WplusMinus_MuAbseta->fill(eta0*-1);
 78        }
 79      }
 80
 81    }
 82
 83
 84    void finalize() {
 85      scale(_hist_Wplus_MuAbseta, crossSection()/picobarn/sumOfWeights());
 86      scale(_hist_Wminus_MuAbseta, crossSection()/picobarn/sumOfWeights());
 87      scale(_hist_WplusMinus_MuAbseta, crossSection()/picobarn/sumOfWeights());
 88    }
 89
 90
 91    Histo1DPtr _hist_Wplus_MuAbseta;
 92    Histo1DPtr _hist_Wminus_MuAbseta;
 93    Histo1DPtr _hist_WplusMinus_MuAbseta;
 94
 95  };
 96
 97
 98  RIVET_DECLARE_PLUGIN(CMS_2019_I1705068);
 99
100}