rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2017_I1467451

Measurement of the transverse momentum spectrum of the Higgs boson produced in $pp$collisions at $\sqrt{s} = 8$ TeV using H to WW decays
Experiment: CMS (LHC)
Inspire ID: 1467451
Status: VALIDATED
Authors:
  • Lorenzo Viliani
References: Beams: p+ p+
Beam energies: (4000.0, 4000.0) GeV
Run details:
  • H->WW->lvlv, both ggH and XH (VBF anf VH) productions modes are included

The cross section for Higgs boson production in pp collisions is studied using the $H\rightarrow W^+ W^-$ decay mode, followed by leptonic decays of the W bosons to an oppositely charged electron-muon pair in the final state. The measurements are performed using data collected by the CMS experiment at the LHC at a centre-of-mass energy of 8 TeV, corresponding to an integrated luminosity of $19.4\,\text{fb}^{-1}$. The Higgs boson transverse momentum ($p_\text{T}$) is reconstructed using the lepton pair $p_\text{T}$ and missing $p_\text{T}$. The differential cross section times branching fraction is measured as a function of the Higgs boson $p_\text{T}$ in a fiducial phase space defined to match the experimental acceptance in terms of the lepton kinematics and event topology. The production cross section times branching fraction in the fiducial phase space is measured to be $39 \pm 8\,\text{(stat)} \pm 9\,\text{(syst)}\,\text{fb}$. The measurements are found to agree, within experimental uncertainties, with theoretical calculations based on the standard model. The fiducial region is defined at particle level (using Born level leptons) as: - emununu, leptons are required not to originate from leptonic tau decays - leading lepton $p_T > 20$ GeV - subleading lepton $p_T > 10$ GeV - psudorapidity of electrons and muons $|\eta| < 2.5$ - Invariant mass of the two charged leptons $m_{ll} > 12$ GeV - Charged lepton pair $p_T > 30$ GeV - Invariant mass of the leptonic system in the transverse plane $m_T^{e\mu\nu\nu} > 50$ GeV Notes on reinterpretation: Major- Despite the existance of a "fiducial" cross section definition, it is not clear this can really be used to constrain any BSM signal other than an excess/deficit in SM Higgs->WW events. There is a b-jet veto is applied to suppress ttbar background, which is however not implemented in the fiducial cross section defintion, and so is extrapolated out during the unfolding. Even after this, the background remains very large and is subtracted using a fit to data. How all this might affect a hypothetical non-SM-H source of electrons and muons (on which the "fiducial" region is based) is very difficult to evaluate. Minor- The rivet routine uses dressed leptons even though the above (and the paper) say they are born level.

Source code: CMS_2017_I1467451.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/FinalState.hh"
  4#include "Rivet/Projections/ChargedLeptons.hh"
  5#include "Rivet/Projections/LeptonFinder.hh"
  6#include "Rivet/Projections/IdentifiedFinalState.hh"
  7#include "Rivet/Projections/PromptFinalState.hh"
  8#include "Rivet/Projections/MissingMomentum.hh"
  9
 10namespace Rivet {
 11
 12
 13  /// Higgs -> WW -> emu + MET in 8 TeV pp collisions
 14  class CMS_2017_I1467451 : public Analysis {
 15  public:
 16
 17    /// Constructor
 18    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2017_I1467451);
 19
 20
 21    /// Book histograms and initialise projections before the run
 22    void init() {
 23
 24      const double lepConeSize = 0.1;
 25      const double lepMaxEta = 2.5;
 26      const Cut lepton_cut = (Cuts::abseta < lepMaxEta);
 27
 28      MSG_WARNING("\033[91;1mLIMITED VALIDITY - check info file for details!\033[m");
 29
 30
 31      // Initialise and register projections
 32      FinalState fs((Cuts::etaIn(-2.5,2.5)));
 33      FinalState fsm((Cuts::etaIn(-5,5)));
 34      declare(fs, "FS");
 35      declare(fsm, "FSM");
 36
 37      ChargedLeptons charged_leptons(fs);
 38      IdentifiedFinalState photons(fs);
 39      photons.acceptIdPair(PID::PHOTON);
 40
 41      PromptFinalState prompt_leptons(charged_leptons);
 42      prompt_leptons.acceptMuonDecays(true);
 43      prompt_leptons.acceptTauDecays(false);
 44
 45      PromptFinalState prompt_photons(photons);
 46      prompt_photons.acceptMuonDecays(true);
 47      prompt_photons.acceptTauDecays(false);
 48
 49      LeptonFinder dressed_leptons(prompt_leptons, prompt_photons, lepConeSize, lepton_cut);
 50      declare(dressed_leptons, "LeptonFinder");
 51
 52      MissingMomentum Met(fsm);
 53      declare(Met, "MET");
 54
 55
 56      // Book histograms
 57      book(histoPtH , 1,1,1);
 58      book(histoXsec, 2,1,1);
 59    }
 60
 61
 62    /// Perform the per-event analysis
 63    void analyze(const Event& event) {
 64
 65      Particles leptons = apply<LeptonFinder>(event, "LeptonFinder").particlesByPt(10.0*GeV);
 66      if (leptons.size() < 2) vetoEvent;
 67      if (leptons[0].pT() < 20*GeV || leptons[1].pT() < 10*GeV) vetoEvent;
 68      if (leptons[0].charge() == leptons[1].charge()) vetoEvent;
 69      if (leptons[0].abspid() == leptons[1].abspid()) vetoEvent;
 70
 71      FourMomentum LL = (leptons[0].momentum() + leptons[1].momentum());
 72      if (LL.mass() < 12*GeV) vetoEvent;
 73      if (LL.pT() < 30*GeV) vetoEvent;
 74
 75      FourMomentum EtMiss = apply<MissingMomentum>(event,"MET").missingMomentum();
 76      FourMomentum P4H = LL + EtMiss;
 77
 78      double dphi = deltaPhi(LL, EtMiss);
 79
 80      double mT = sqrt(2*LL.pT()*EtMiss.pT()*(1-cos(dphi)));
 81      if (mT < 50*GeV) vetoEvent;
 82
 83      histoPtH->fill(min(P4H.pT()/GeV, 199.));
 84      histoXsec->fill(8000); ///< @todo Should probably be a Counter
 85    }
 86
 87
 88    /// Normalise histograms etc., after the run
 89    void finalize() {
 90      scale(histoPtH, crossSection()/sumOfWeights()/femtobarn);
 91      scale(histoXsec, (histoXsec->xMax()-histoXsec->xMin())*crossSection()/sumOfWeights()/femtobarn);
 92    }
 93
 94
 95  private:
 96
 97    Histo1DPtr histoPtH, histoXsec;
 98
 99  };
100
101
102  RIVET_DECLARE_PLUGIN(CMS_2017_I1467451);
103
104}