rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2012_I1087342

Measurement of forward and forward+central jets at $\sqrt{s} = 7$ TeV
Experiment: CMS (LHC)
Inspire ID: 1087342
Status: VALIDATED
Authors:
  • Albert Knutsson
  • Rasmus Sloth Hansen
  • Bo Zhu
References:
  • JHEP 1206 (2012) 036
  • Public page: CMS-FWD-11-002
  • CERN-PH-EP-2011-179
  • doi 10.1007/JHEP06(2012)036
  • arXiv: 1202.0704
Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • pp QCD interactions at 7 TeV.

Inclusive forward jets and forward+central jets measured by CMS at $\sqrt{s} = 7$ TeV.

Source code: CMS_2012_I1087342.cc
 1// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Projections/FinalState.hh"
 4#include "Rivet/Projections/FastJets.hh"
 5
 6
 7namespace Rivet {
 8
 9  // This analysis is a derived from the class Analysis:
10  class CMS_2012_I1087342 : public Analysis {
11
12  public:
13
14    // Constructor
15    CMS_2012_I1087342() : Analysis("CMS_2012_I1087342") {
16    }
17
18    void init() {
19      const FinalState fs;
20      declare(FastJets(fs, JetAlg::ANTIKT, 0.5),"Jets");
21
22      book(_hist_jetpt_fwdincl ,1, 1, 1);
23      book(_hist_jetpt_forward ,2, 1, 1);
24      book(_hist_jetpt_central ,3, 1, 1);
25    }
26
27    void analyze(const Event& event) {
28
29      const FastJets& fj = apply<FastJets>(event,"Jets");
30      const Jets jets = fj.jets(Cuts::ptIn(35*GeV, 150*GeV) && Cuts::abseta < 4.7);
31
32      double cjet_pt = 0.0;
33      double fjet_pt = 0.0;
34
35      for(const Jet& j : jets) {
36        double pT = j.pT();
37        if (j.abseta() > 3.2) {
38          _hist_jetpt_fwdincl->fill(j.pT()/GeV);
39        }
40        if (j.abseta() < 2.8) {
41          if (cjet_pt < pT) cjet_pt = pT;
42        }
43        if (inRange(j.abseta(), 3.2, 4.7)) {
44          if (fjet_pt < pT) fjet_pt = pT;
45        }
46      }
47
48      if (cjet_pt > 35*GeV && fjet_pt > 35*GeV) {
49        _hist_jetpt_forward->fill(fjet_pt/GeV);
50        _hist_jetpt_central->fill(cjet_pt/GeV);
51      }
52
53    }
54
55
56    void finalize() {
57      scale(_hist_jetpt_fwdincl, crossSection() / picobarn / sumOfWeights() / 3.0);
58      scale(_hist_jetpt_forward, crossSection() / picobarn / sumOfWeights() / 3.0);
59      scale(_hist_jetpt_central, crossSection() / picobarn / sumOfWeights() / 5.6);
60    }
61
62
63  private:
64
65    Histo1DPtr _hist_jetpt_fwdincl;
66    Histo1DPtr _hist_jetpt_forward;
67    Histo1DPtr _hist_jetpt_central;
68
69  };
70
71
72  RIVET_DECLARE_PLUGIN(CMS_2012_I1087342);
73
74}