rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2020_I1837084

Measurement of the Z boson differential production cross section using its invisible decay mode in proton-proton collisions at 13 TeV
Experiment: CMS (LHC)
Inspire ID: 1837084
Status: VALIDATED
Authors:
  • cms-pag-conveners-smp@cern.ch
  • Guillelmo Gomez-Ceballos
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 of the total and differential fiducial cross sections for the Z boson decaying into two neutrinos are presented at the LHC in proton-proton collisions at a center-of-mass energy of 13 TeV. The data were collected by the CMS detector in 2016 and correspond to an integrated luminosity of 35.9 fb-1. In these measurements, events are selected containing an imbalance in transverse momentum and one or more energetic jets. The fiducial differential cross section is measured as a function of the Z boson transverse momentum. The results are combined with a previous measurement of charged-lepton decays of the Z boson.

Source code: CMS_2020_I1837084.cc
 1// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Projections/FinalState.hh"
 4#include "Rivet/Projections/FastJets.hh"
 5#include "Rivet/Projections/DressedLeptons.hh"
 6#include "Rivet/Projections/MissingMomentum.hh"
 7#include "Rivet/Projections/PromptFinalState.hh"
 8#include "Rivet/Projections/ZFinder.hh"
 9
10namespace Rivet {
11
12
13  /// Measurements of differential Z-boson -> ll and vv production cross-sections in 13 TeV pp collisions
14  class CMS_2020_I1837084 : public Analysis {
15  public:
16
17    /// Constructor
18    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2020_I1837084);
19
20
21    /// @name Analysis methods
22    /// @{
23
24    /// Book histograms and initialise projections before the run
25    void init() {
26
27      // Initialise and register projections
28      ZFinder zmmFind(FinalState(), Cuts::pT > 0*GeV, PID::MUON, 76.1876*GeV, 106.1876*GeV, 0.1,
29                      ZFinder::ChargedLeptons::PROMPT, ZFinder::ClusterPhotons::NODECAY, ZFinder::AddPhotons::YES );
30      declare(zmmFind, "ZmmFind");
31
32      // Book histograms
33      book(_h_Z_pt,      12, 1, 1);
34      book(_h_Z_pt_norm, 13, 1, 1);
35
36    }
37
38
39    /// Perform the per-event analysis
40    void analyze(const Event& event) {
41
42      const Particles& zmms = apply<ZFinder>(event, "ZmmFind").bosons();
43
44      if (zmms.size() == 1 && zmms[0].pT() > 200*GeV) {
45        _h_Z_pt     ->fill(min(zmms[0].pT()/GeV, 1499.999));
46        _h_Z_pt_norm->fill(min(zmms[0].pT()/GeV, 1499.999));
47      }
48
49    }
50
51
52    /// @todo Replace with barchart()
53    void normalizeToSum(Histo1DPtr hist) {
54      //double sum = 0.;
55      for (size_t i = 0; i < hist->numBins(); ++i) {
56        //sum += hist->bin(i).height();
57        double width = hist->bin(i).width();
58        hist->bin(i).scaleW(width != 0 ? width : 1.);
59      }
60      if (hist->integral() > 0) scale(hist, 1./hist->integral());
61    }
62
63
64    /// Normalise histograms etc., after the run
65    void finalize() {
66
67      double norm = (sumOfWeights() != 0) ? crossSection()/femtobarn/sumOfWeights() : 1.0;
68
69      scale(_h_Z_pt, norm);
70
71      normalizeToSum(_h_Z_pt_norm);
72
73    }
74
75    /// @}
76
77
78    /// @name Histograms
79    /// @{
80    Histo1DPtr _h_Z_pt, _h_Z_pt_norm;
81    /// @}
82
83
84  };
85
86
87  RIVET_DECLARE_PLUGIN(CMS_2020_I1837084);
88
89}