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.Karen.Pflitsch@cern.ch
  • Katerina.Lipka@cern.ch
References:
  • Eur.Phys.J. C79 (2019) no.3, 269
  • DOI:10.1140/epjc/s10052-019-6752-1
  • arXiv: 1811.10021
  • CMS-SMP-17-014
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
  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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/WFinder.hh"
#include "Rivet/Projections/UnstableParticles.hh"

namespace Rivet {


  // Measurement of associated production of a W boson and a charm quark in proton-proton collisions at 13 TeV
  class CMS_2019_I1705068 : public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2019_I1705068);


    void init() {

      // Projections
      FinalState fs;
      WFinder wfinder_mu(fs, Cuts::abseta < 2.4 && Cuts::pT > 0*GeV, PID::MUON,
                         0*GeV, 1000000*GeV, 0*GeV, 0.1, WFinder::ChargedLeptons::PROMPT,
                         WFinder::ClusterPhotons::NODECAY, WFinder::AddPhotons::NO, WFinder::MassWindow::MT);
      declare(wfinder_mu, "WFinder_mu");

      UnstableParticles dst(Cuts::pT > 5*GeV && Cuts::abseta < 2.4);
      declare(dst, "Dstar");

      // Particle-Level histograms from the paper:
      book(_hist_WplusMinus_MuAbseta, "d04-x01-y01");
      book(_hist_Wplus_MuAbseta, "d05-x01-y01");
      book(_hist_Wminus_MuAbseta, "d06-x01-y01");

    }


    void analyze(const Event& event) {

      // Get the reconstructed W
      const WFinder& wfinder_mu = applyProjection<WFinder>(event, "WFinder_mu");
      if (wfinder_mu.bosons().size() != 1) vetoEvent;

      // No Missing Energy or MT cut at generator level:
      const FourMomentum& lepton0 = wfinder_mu.constituentLeptons()[0].momentum();
      double pt0 = lepton0.pT();
      double eta0 = fabs( lepton0.eta() );
      if ( (eta0 > 2.4) || (pt0 < 26.0*GeV) ) vetoEvent;

      int muID = wfinder_mu.constituentLeptons()[0].pid();


      // D* selection:
      // OS = W boson and D* Meson have Opposite (charge) Signs
      // SS = W Boson and D* Meson have Same (charge) Signs
      // Associated W+c only has OS contributions,
      // W+ccbar (ccbar from gluon splitting) has equal probability to be OS or SS
      // OS-SS to remove the gluon splitting background

      const UnstableParticles& dst = applyProjection<UnstableParticles>(event, "Dstar");
      for (auto p: dst.particles()) {
        if (muID == -13 && p.pid() == -413) { // OS
          _hist_Wplus_MuAbseta->fill(eta0);
          _hist_WplusMinus_MuAbseta->fill(eta0);
        }
        else if (muID == 13 && p.pid() == 413) { // OS
          _hist_Wminus_MuAbseta->fill(eta0);
          _hist_WplusMinus_MuAbseta->fill(eta0);
        }
        else if (muID == -13 && p.pid() == 413) { // SS
          _hist_Wplus_MuAbseta->fill(eta0*-1);
          _hist_WplusMinus_MuAbseta->fill(eta0*-1);
        }
        else if (muID == 13 && p.pid() == -413) { // SS
          _hist_Wminus_MuAbseta->fill(eta0*-1);
          _hist_WplusMinus_MuAbseta->fill(eta0*-1);
        }
      }

    }


    void finalize() {
      scale(_hist_Wplus_MuAbseta, crossSection()/picobarn/sumOfWeights());
      scale(_hist_Wminus_MuAbseta, crossSection()/picobarn/sumOfWeights());
      scale(_hist_WplusMinus_MuAbseta, crossSection()/picobarn/sumOfWeights());
    }


    Histo1DPtr _hist_Wplus_MuAbseta;
    Histo1DPtr _hist_Wminus_MuAbseta;
    Histo1DPtr _hist_WplusMinus_MuAbseta;

  };



  RIVET_DECLARE_PLUGIN(CMS_2019_I1705068);

}