Rivet analyses referenceLHCB_2021_I1990313Forward $Z$ boson production cross-section in proton-proton collisions at $\sqrt{s}=13$ TeVExperiment: LHCB (LHC) Inspire ID: 1990313 Status: VALIDATED Authors:
Beam energies: (6500.0, 6500.0) GeV Run details:
A precision measurement of the $Z^0$ boson production cross-section at $\sqrt{s} = 13$ TeV in the forward region is presented, using $pp$ collision data collected by the LHCb detector, corresponding to an integrated luminosity of 5.1 fb$^{-1}$. The production cross-section is measured using $Z \to \mu^+ \mu^-$ events within the fiducial region defined as pseudorapidity $2.0<\eta<4.5$ and transverse momentum $p_{T}>20$ GeV$/c$ for both muons and dimuon invariant mass $60<M_{\mu\mu}<120$ GeV$/c^{2}$. The integrated cross-section is determined to be \begin{equation*} \sigma(Z\to\mu^+\mu^-) = 196.4 \pm0.2\pm 1.6\pm 3.9 \ \mathrm{pb}, \end{equation*} where the first uncertainty is statistical, the second is systematic, and the third is due to the luminosity determination. The measured results are in agreement with theoretical predictions within uncertainties. Source code: LHCB_2021_I1990313.cc 1// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/FinalState.hh"
4#include "Rivet/Projections/DileptonFinder.hh"
5
6namespace Rivet {
7
8
9 /// @brief Forward Z production cross-section in pp collisions at 13 TeV
10 class LHCB_2021_I1990313 : public Analysis {
11 public:
12
13 /// Constructor
14 RIVET_DEFAULT_ANALYSIS_CTOR(LHCB_2021_I1990313);
15
16
17 /// @name Analysis methods
18 ///@{
19
20 /// Book histograms and initialise projections before the run
21 void init() {
22
23 // Initialise and register projections
24 DileptonFinder zmumufinder(91.2*GeV, 0.1, Cuts::absetaIn(2.0, 4.5) && Cuts::pT > 20 * GeV &&
25 Cuts::abspid == PID::MUON, Cuts::massIn(60*GeV, 120*GeV));
26 declare(zmumufinder, "ZmumuFinder");
27
28 // Book histograms
29 // specify custom binning
30 book(_h_sigma_vs_y, 18, 1, 1);
31 book(_h_sigma_vs_pt, 19, 1, 1);
32 book(_h_sigma_vs_phi, 20, 1, 1);
33 book(_h_sigma_vs_ypt, {2., 2.5, 3., 3.5, 4., 4.5});
34 book(_h_sigma_vs_yphi, {2., 2.5, 3., 3.5, 4., 4.5});
35 for (size_t i = 1; i < _h_sigma_vs_ypt->numBins()+1; ++i) {
36 book(_h_sigma_vs_ypt->bin(i), 21, 1, i);
37 book(_h_sigma_vs_yphi->bin(i), 22, 1, i);
38 }
39 }
40
41
42 /// Perform the per-event analysis
43 void analyze(const Event& event) {
44
45 // Retrieve dressed leptons, sorted by pT
46 const DileptonFinder& zmumufinder = apply<DileptonFinder>(event, "ZmumuFinder");
47 if(zmumufinder.empty()) vetoEvent;
48 if(zmumufinder.bosons().size() > 1)
49 MSG_WARNING("Found multiple (" << zmumufinder.bosons().size() << ") Z -> mu+ mu- decays!");
50
51 // Z momenta
52 FourMomentum zmumu = zmumufinder.bosons()[0].momentum();
53 if (zmumufinder.leptons().size() < 2) vetoEvent;
54
55 const Particle& muon_p = zmumufinder.constituents()[0];
56 const Particle& muon_m = zmumufinder.constituents()[1];
57
58 const double diffphi = deltaPhi(muon_p, muon_m);
59 const double diffpsd = deltaEta(muon_p, muon_m);
60 const double accphi = M_PI - diffphi;
61 const double angular = tan(accphi/2) / cosh(diffpsd/2);
62
63 _h_sigma_vs_y->fill(zmumu.absrapidity());
64 _h_sigma_vs_pt->fill(zmumu.pT()/GeV);
65 _h_sigma_vs_phi->fill(angular);
66 _h_sigma_vs_ypt->fill(zmumu.absrapidity(), zmumu.pT()/GeV);
67 _h_sigma_vs_yphi->fill(zmumu.absrapidity(), angular);
68 }
69
70
71 /// Normalise histograms etc., after the run
72 void finalize() {
73 const double xs = crossSection()/picobarn;
74 double scale_f = xs/sumOfWeights()/2.;
75 _h_sigma_vs_y->scaleW(scale_f);
76 _h_sigma_vs_pt->scaleW(scale_f);
77 _h_sigma_vs_phi->scaleW(scale_f);
78 scale(_h_sigma_vs_ypt, scale_f*2.);
79 scale(_h_sigma_vs_yphi, scale_f*2.);
80 }
81
82 ///@}
83
84
85 Histo1DPtr _h_sigma_vs_y, _h_sigma_vs_pt, _h_sigma_vs_phi;
86 Histo1DGroupPtr _h_sigma_vs_ypt, _h_sigma_vs_yphi;
87
88 };
89
90
91 RIVET_DECLARE_PLUGIN(LHCB_2021_I1990313);
92
93}
|