rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

LHCB_2012_I1208102

Differential cross-sections of $\mathrm{Z}/\gamma^* \to e^{+}e^{-}$ vs rapidity and $\phi^*$
Experiment: LHCb (LHC)
Inspire ID: 1208102
Status: VALIDATED
Authors:
  • Ana Elena Dumitriu
  • Marek Sirendi
  • David Ward
  • Alex Grecu
References: Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • $\mathrm{Z}/\gamma^* \to e^{+}e^{-}_$ decays with di-lepton invariant mass $> 40$ GeV/$c^2$ produced in $pp$ collisions at $\sqrt{s}=7$ TeV.

Measurement of the $pp \to \mathrm{Z}^0$ cross-section in the $\mathrm{Z}/\gamma \to e^{+}e^{-}$ mode at $\sqrt{s} = 7$ TeV. Daughter electrons are required to have $p_T > 20$ GeV/$c$, $2 < \eta < 4.5$ and the dielectron invariant mass in range 60-120 GeV/$c^2$. The cross-section is given as a function of $Z$ rapidity and an angular variable ($\phi^*$) closely related to $Z$ transverse momentum (derived from the lepton pseudorapidity and azimuthal angle differences). For event generators implementing cross-section QCD corrections only at LO the distributions are normalized to the cross-section measured in data $76.0 \pm 0.8 \pm 2.0 \pm 2.6 \pm 0.4$ pb, where the first uncertainty is statistical, the second is systematic, the third is due to luminosity uncertainty and the fourth to FSR corrections.

Source code: LHCB_2012_I1208102.cc
 1// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Projections/DileptonFinder.hh"
 4
 5namespace Rivet {
 6
 7
 8  /// Differential cross-sections of $\mathrm{Z}/\gamma^* \to e^{+}e^{-}$ vs rapidity and $\phi^*$
 9  class LHCB_2012_I1208102 : public Analysis {
10  public:
11
12    /// Constructor
13    RIVET_DEFAULT_ANALYSIS_CTOR(LHCB_2012_I1208102);
14
15
16    /// @name Analysis methods
17    /// @{
18
19    /// Book histograms
20    void init() {
21      DileptonFinder zeefinder(91.2*GeV, 0.1, Cuts::etaIn(2.0, 4.5) && Cuts::pT > 20*GeV &&
22                               Cuts::abspid == PID::ELECTRON, Cuts::massIn(60*GeV, 120*GeV));
23      declare(zeefinder, "ZeeFinder");
24
25      book(_h_sigma_vs_y,   2, 1, 1);
26      book(_h_sigma_vs_phi, 3, 1, 1);
27    }
28
29
30    /// Do the analysis
31    void analyze(const Event& e) {
32      const DileptonFinder& zeefinder = apply<DileptonFinder>(e, "ZeeFinder");
33      if (zeefinder.empty()) vetoEvent;
34      if (zeefinder.bosons().size() > 1)
35        MSG_WARNING("Found multiple (" << zeefinder.bosons().size() << ") Z -> e+ e- decays!");
36
37      // Z momenta
38      const FourMomentum zee = zeefinder.bosons()[0].momentum();
39
40      if (zeefinder.constituents().size() < 2) vetoEvent;
41
42      const Particle pozitron = zeefinder.constituents()[0];
43      const Particle electron = zeefinder.constituents()[1];
44
45      // Calculation of the angular variable
46      const double diffphi = deltaPhi(pozitron, electron);
47      const double diffpsd = deltaEta(pozitron, electron);
48      const double accphi = M_PI - diffphi;
49      const double angular = tan(accphi/2) / cosh(diffpsd/2);
50
51      // Fill histograms
52      _h_sigma_vs_y->fill(zee.rapidity());
53      _h_sigma_vs_phi->fill(angular);
54    }
55
56
57    /// Finalize
58    void finalize() {
59      const double xs = crossSection()/picobarn;
60      scale(_h_sigma_vs_y, xs/sumOfWeights());
61      scale(_h_sigma_vs_phi, xs/sumOfWeights());
62    }
63
64    /// @}
65
66
67  private:
68
69    /// @name Histograms
70    /// @{
71    Histo1DPtr _h_sigma_vs_y, _h_sigma_vs_phi;
72    /// @}
73
74  };
75
76
77  RIVET_DECLARE_PLUGIN(LHCB_2012_I1208102);
78
79}