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/ZFinder.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
13    /// Constructor
14    LHCB_2012_I1208102()
15      : Analysis("LHCB_2012_I1208102")
16    {    }
17
18
19    /// @name Analysis methods
20    //@{
21
22    /// Book histograms
23    void init() {
24      ZFinder zeefinder(FinalState(), Cuts::etaIn(2.0, 4.5) && Cuts::pT > 20*GeV, PID::ELECTRON, 60*GeV, 120*GeV);
25      declare(zeefinder, "ZeeFinder");
26
27      book(_h_sigma_vs_y ,2, 1, 1);
28      book(_h_sigma_vs_phi ,3, 1, 1);
29    }
30
31
32    /// Do the analysis
33    void analyze(const Event& e) {
34      const ZFinder& zeefinder = apply<ZFinder>(e, "ZeeFinder");
35      if (zeefinder.empty()) vetoEvent;
36      if (zeefinder.bosons().size() > 1)
37        MSG_WARNING("Found multiple (" << zeefinder.bosons().size() << ") Z -> e+ e- decays!");
38
39      // Z momenta
40      const FourMomentum zee = zeefinder.bosons()[0].momentum();
41
42      if (zeefinder.constituents().size() < 2) vetoEvent;
43
44      const Particle pozitron = zeefinder.constituents()[0];
45      const Particle electron = zeefinder.constituents()[1];
46
47      // Calculation of the angular variable
48      const double diffphi = deltaPhi(pozitron, electron);
49      const double diffpsd = deltaEta(pozitron, electron);
50      const double accphi = M_PI - diffphi;
51      const double angular = tan(accphi/2) / cosh(diffpsd/2);
52
53      // Fill histograms
54      _h_sigma_vs_y->fill(zee.rapidity());
55      _h_sigma_vs_phi->fill(angular);
56    }
57
58
59    /// Finalize
60    void finalize() {
61      const double xs = crossSection()/picobarn;
62      scale(_h_sigma_vs_y, xs/sumOfWeights());
63      scale(_h_sigma_vs_phi, xs/sumOfWeights());
64    }
65
66    //@}
67
68
69  private:
70
71    /// @name Histograms
72    //@{
73    Histo1DPtr _h_sigma_vs_y, _h_sigma_vs_phi;
74    //@}
75
76  };
77
78
79  // The hook for the plugin system
80  RIVET_DECLARE_PLUGIN(LHCB_2012_I1208102);
81
82}