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

namespace Rivet {


  /// Differential cross-sections of $\mathrm{Z}/\gamma^* \to e^{+}e^{-}$ vs rapidity and $\phi^*$
  class LHCB_2012_I1208102 : public Analysis {
  public:


    /// Constructor
    LHCB_2012_I1208102()
      : Analysis("LHCB_2012_I1208102")
    {    }


    /// @name Analysis methods
    //@{

    /// Book histograms
    void init() {
      ZFinder zeefinder(FinalState(), Cuts::etaIn(2.0, 4.5) && Cuts::pT > 20*GeV, PID::ELECTRON, 60*GeV, 120*GeV);
      declare(zeefinder, "ZeeFinder");

      book(_h_sigma_vs_y ,2, 1, 1);
      book(_h_sigma_vs_phi ,3, 1, 1);
    }


    /// Do the analysis
    void analyze(const Event& e) {
      const ZFinder& zeefinder = apply<ZFinder>(e, "ZeeFinder");
      if (zeefinder.empty()) vetoEvent;
      if (zeefinder.bosons().size() > 1)
        MSG_WARNING("Found multiple (" << zeefinder.bosons().size() << ") Z -> e+ e- decays!");

      // Z momenta
      const FourMomentum zee = zeefinder.bosons()[0].momentum();

      if (zeefinder.constituents().size() < 2) vetoEvent;

      const Particle pozitron = zeefinder.constituents()[0];
      const Particle electron = zeefinder.constituents()[1];

      // Calculation of the angular variable
      const double diffphi = deltaPhi(pozitron, electron);
      const double diffpsd = deltaEta(pozitron, electron);
      const double accphi = M_PI - diffphi;
      const double angular = tan(accphi/2) / cosh(diffpsd/2);

      // Fill histograms
      _h_sigma_vs_y->fill(zee.rapidity());
      _h_sigma_vs_phi->fill(angular);
    }


    /// Finalize
    void finalize() {
      const double xs = crossSection()/picobarn;
      scale(_h_sigma_vs_y, xs/sumOfWeights());
      scale(_h_sigma_vs_phi, xs/sumOfWeights());
    }

    //@}


  private:

    /// @name Histograms
    //@{
    Histo1DPtr _h_sigma_vs_y, _h_sigma_vs_phi;
    //@}

  };


  // The hook for the plugin system
  RIVET_DECLARE_PLUGIN(LHCB_2012_I1208102);

}