rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2014_I1322726

Z production in PbPb and pp collisions at 2.76 TeV in the dimuon and dielectron decay channels
Experiment: CMS (LHC)
Inspire ID: 1322726
Status: VALIDATED
Authors:
  • Hannes Jung
  • Sara Taheri Monfared
  • Muhammad Ibrahim Abdulhamid
References: Beams: p+ p+
Beam energies: (1380.0, 1380.0) GeV
Run details:
  • DY events with m(ll)>60 GeV, run in pp mode (PbPb to be added)

The differential cross-section for Z production in PbPb and pp collisions at 2.76 TeV in the dimuon and dielectron decay channels has been measured from data acquired by the CMS collaboration at the LHC during the year 2010.

Source code: CMS_2014_I1322726.cc
 1// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Projections/FinalState.hh"
 4#include "Rivet/Projections/FastJets.hh"
 5#include "Rivet/Projections/DileptonFinder.hh"
 6
 7namespace Rivet {
 8
 9
10  /// Z production in PbPb and pp collisions at 2.76 TeV in the dimuon and dielectron decay channels
11  class CMS_2014_I1322726 : public Analysis {
12  public:
13
14    /// Constructor
15    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2014_I1322726);
16
17    /// @name Analysis methods
18    ///@{
19
20    /// Book histograms and initialise projections before the run
21    void init() {
22
23      DileptonFinder zmumuFinder(91.2*GeV, 0.1, Cuts::abseta < 4.9 && Cuts::abspid == PID::MUON, Cuts::massIn(60*GeV, 120*GeV));
24      declare(zmumuFinder, "ZmumuFinder");
25
26      DileptonFinder zeeFinder(91.2*GeV, 0.1, Cuts::abseta < 4.9 && Cuts::abspid == PID::ELECTRON, Cuts::massIn(60*GeV, 120*GeV));
27      declare(zeeFinder, "ZeeFinder");
28
29      FastJets jetfs(FinalState(Cuts::abseta < 4.9), JetAlg::ANTIKT, 0.4, JetMuons::NONE, JetInvisibles::NONE);
30      declare(jetfs, "jets");
31
32      book(_h["ds/dydp-mu"], 1, 1, 1);
33      book(_h["ds/dydp-el"], 2, 1, 1);
34      book(_h["ds/dy-mu"], 3, 1, 1);
35      book(_h["ds/dy-el"], 4, 1, 1);
36    }
37
38
39    /// Perform the per-event analysis
40    void analyze(const Event& event) {
41      const DileptonFinder& zfindermu = apply<DileptonFinder>(event, "ZmumuFinder");  // muon
42      const Particles& zmumus = zfindermu.bosons();
43
44      const DileptonFinder& zfinderel = apply<DileptonFinder>(event, "ZeeFinder");  // electron
45      const Particles& zelels = zfinderel.bosons();
46
47      if ((zmumus.size() + zelels.size()) != 1)
48        vetoEvent;
49
50      if (zmumus.size() == 1) {                          //cout <<"zmumus " <<zmumus <<endl;
51        if (abs(zmumus[0].momentum().rapidity()) < 2) {  //5 //2
52
53          _h["ds/dydp-mu"]->fill(zmumus[0].momentum().pT() / GeV, 0.25);
54          _h["ds/dy-mu"]->fill(zmumus[0].momentum().rapidity());
55        }
56      }
57
58      if (zelels.size() == 1) {                             //cout <<"electron-find"<< endl;
59        if (abs(zelels[0].momentum().rapidity()) < 1.44) {  //5 //1.44
60
61          _h["ds/dydp-el"]->fill(zelels[0].momentum().pT() / GeV, 1 / 2.88);
62          _h["ds/dy-el"]->fill(zelels[0].momentum().rapidity());
63        }
64      }
65    }
66
67
68    void finalize() {
69      double norm = crossSection() / picobarn / sumW();
70
71      scale(_h["ds/dydp-mu"], norm);
72      scale(_h["ds/dydp-el"], norm);
73      scale(_h["ds/dy-mu"], norm);
74      scale(_h["ds/dy-el"], norm);
75    }
76
77    ///@}
78
79
80    /// Histograms
81    map<string, Histo1DPtr> _h;
82
83  };
84
85
86  RIVET_DECLARE_PLUGIN(CMS_2014_I1322726);
87
88}