|
Rivet analyses reference
LHCB_2021_I1990313
Forward $Z$ boson production cross-section in proton-proton collisions at $\sqrt{s}=13$ TeV
Experiment: LHCB (LHC)
Inspire ID: 1990313
Status: VALIDATED
Authors:
References:
- JHEP 07 (2022) 026
- DOI:10.1007/JHEP07(2022)026
- arXiv: 2112.07458
Beams: p+ p+
Beam energies: (6500.0, 6500.0) GeV
Run details:
- Proton-proton collisions at $\sqrt{s}=13$ TeV, measurement of the $Z$ boson production cross section in the forward region. The production cross-section is measured using the $Z\to\mu^+\mu^-$ events in the fiducial region, defined as pseudorapidity $2.0<\eta<4.5$ and $p_T>20$ GeV/$c$ for both constituent muons, and dimuon invariant mass $60<M_{\mu\mu}<120$ GeV/$c^{2}$.
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
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 | // -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/ZFinder.hh"
namespace Rivet {
/// @brief Forward Z production cross-section in pp collisions at 13 TeV
class LHCB_2021_I1990313 : public Analysis {
public:
/// Constructor
RIVET_DEFAULT_ANALYSIS_CTOR(LHCB_2021_I1990313);
/// @name Analysis methods
///@{
/// Book histograms and initialise projections before the run
void init() {
// Initialise and register projections
ZFinder zmumufinder(FinalState(), Cuts::absetaIn(2.0, 4.5) && Cuts::pT > 20 * GeV, PID::MUON, 60*GeV, 120*GeV);
declare(zmumufinder, "ZmumuFinder");
// Book histograms
// specify custom binning
book(_h_sigma_vs_y, 18, 1, 1);
book(_h_sigma_vs_pt, 19, 1, 1);
book(_h_sigma_vs_phi, 20, 1, 1);
double ylow, yhigh;
for (int i = 1; i < 6; ++i) {
ylow = 1.5 + 0.5*(double)i;
yhigh = ylow + 0.5;
{Histo1DPtr tmp; _h_sigma_vs_ypt.add (ylow, yhigh, book(tmp, 21, 1, i) );};
{Histo1DPtr tmp; _h_sigma_vs_yphi.add(ylow, yhigh, book(tmp, 22, 1, i) );};
}
}
/// Perform the per-event analysis
void analyze(const Event& event) {
// Retrieve dressed leptons, sorted by pT
const ZFinder& zmumufinder = apply<ZFinder>(event, "ZmumuFinder");
if(zmumufinder.empty()) vetoEvent;
if(zmumufinder.bosons().size() > 1)
MSG_WARNING("Found multiple (" << zmumufinder.bosons().size() << ") Z -> mu+ mu- decays!");
// Z momenta
FourMomentum zmumu = zmumufinder.bosons()[0].momentum();
if(zmumufinder.constituentLeptons().size() < 2) vetoEvent;
const Particle& muon_p = zmumufinder.constituentLeptons()[0];
const Particle& muon_m = zmumufinder.constituentLeptons()[1];
const double diffphi = deltaPhi(muon_p, muon_m);
const double diffpsd = deltaEta(muon_p, muon_m);
const double accphi = M_PI - diffphi;
const double angular = tan(accphi/2) / cosh(diffpsd/2);
_h_sigma_vs_y->fill(zmumu.rapidity());
_h_sigma_vs_pt->fill(zmumu.pT()/GeV);
_h_sigma_vs_phi->fill(angular);
_h_sigma_vs_ypt.fill(zmumu.rapidity(), zmumu.pT()/GeV);
_h_sigma_vs_yphi.fill(zmumu.rapidity(), angular);
}
/// Normalise histograms etc., after the run
void finalize() {
const double xs = crossSection()/picobarn;
double scale_f = xs/sumOfWeights()/2.;
_h_sigma_vs_y->scaleW(scale_f);
_h_sigma_vs_pt->scaleW(scale_f);
_h_sigma_vs_phi->scaleW(scale_f);
for (Histo1DPtr h : _h_sigma_vs_ypt.histos()) h->scaleW(scale_f*2.);
for (Histo1DPtr h : _h_sigma_vs_yphi.histos()) h->scaleW(scale_f*2.);
}
///@}
Histo1DPtr _h_sigma_vs_y, _h_sigma_vs_pt, _h_sigma_vs_phi;
BinnedHistogram _h_sigma_vs_ypt, _h_sigma_vs_yphi;
};
RIVET_DECLARE_PLUGIN(LHCB_2021_I1990313);
}
|
|