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
 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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/DressedLeptons.hh"
#include "Rivet/Projections/MissingMomentum.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/ZFinder.hh"

namespace Rivet {


  /// Measurements of differential Z-boson -> ll and vv production cross-sections in 13 TeV pp collisions
  class CMS_2020_I1837084 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2020_I1837084);


    /// @name Analysis methods
    /// @{

    /// Book histograms and initialise projections before the run
    void init() {

      // Initialise and register projections
      ZFinder zmmFind(FinalState(), Cuts::pT > 0*GeV, PID::MUON, 76.1876*GeV, 106.1876*GeV, 0.1,
                      ZFinder::ChargedLeptons::PROMPT, ZFinder::ClusterPhotons::NODECAY, ZFinder::AddPhotons::YES );
      declare(zmmFind, "ZmmFind");

      // Book histograms
      book(_h_Z_pt,      12, 1, 1);
      book(_h_Z_pt_norm, 13, 1, 1);

    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {

      const Particles& zmms = apply<ZFinder>(event, "ZmmFind").bosons();

      if (zmms.size() == 1 && zmms[0].pT() > 200*GeV) {
        _h_Z_pt     ->fill(min(zmms[0].pT()/GeV, 1499.999));
        _h_Z_pt_norm->fill(min(zmms[0].pT()/GeV, 1499.999));
      }

    }


    /// @todo Replace with barchart()
    void normalizeToSum(Histo1DPtr hist) {
      double sum = 0.;
      for (size_t i = 0; i < hist->numBins(); ++i) {
        sum += hist->bin(i).height();
        double width = hist->bin(i).width();
        hist->bin(i).scaleW(width != 0 ? width : 1.);
      }
      if (hist->integral() > 0) scale(hist, 1./hist->integral());
    }


    /// Normalise histograms etc., after the run
    void finalize() {

      double norm = (sumOfWeights() != 0) ? crossSection()/femtobarn/sumOfWeights() : 1.0;

      scale(_h_Z_pt, norm);

      normalizeToSum(_h_Z_pt_norm);

    }

    /// @}


    /// @name Histograms
    /// @{
    Histo1DPtr _h_Z_pt, _h_Z_pt_norm;
    /// @}


  };


  RIVET_DECLARE_PLUGIN(CMS_2020_I1837084);

}