rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

CMS_2011_I902309

Measurement of the inclusive jet cross-section in $pp$ collisions at $\sqrt{s} = 7$ TeV
Experiment: CMS (LHC)
Inspire ID: 902309
Status: VALIDATED
Authors:
  • Rasmus Sloth Hansen<rsh07@phys.au.dk>
References:
  • http://cdsweb.cern.ch/record/1355680
Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Inclusive QCD at 7TeV comEnergy, ptHat (or equivalent) greater than 10 GeV

The inclusive jet cross section is measured in pp collisions with a center-of-mass energy of 7 TeV at the LHC using the CMS experiment. The data sample corresponds to an integrated luminosity of 34 inverse picobarns. The measurement is made for jet transverse momenta in the range 18-1100 GeV and for absolute values of rapidity less than 3. Jets are anti-kt with $R=0.5$, $p_\perp>18$ GeV and $|y|<3.0$.

Source code: CMS_2011_I902309.cc
 1// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Projections/FinalState.hh"
 4#include "Rivet/Projections/FastJets.hh"
 5
 6namespace Rivet {
 7
 8
 9  // Inclusive jet pT
10  class CMS_2011_I902309 : public Analysis {
11  public:
12
13    RIVET_DEFAULT_ANALYSIS_CTOR(CMS_2011_I902309);
14
15
16    /// @{
17
18    /// Book histograms and initialize projections:
19    void init() {
20      const FinalState fs;
21
22      // Initialize the projectors:
23      declare(FastJets(fs, JetAlg::ANTIKT, 0.5),"Jets");
24
25      // Book histograms:
26      book(_hist_sigma, {0., 0.5, 1., 1.5, 2., 2.5, 3.});
27      for (auto& b : _hist_sigma->bins()) {
28        book(b, b.index(), 1, 1);
29      }
30    }
31
32    /// Analysis
33    void analyze(const Event &event) {
34
35      const FastJets& fj = apply<FastJets>(event,"Jets");
36      const Jets& jets = fj.jets(Cuts::ptIn(18*GeV, 1100.0*GeV) && Cuts::absrap < 4.7);
37
38      // Fill the relevant histograms:
39      for(const Jet& j : jets) {
40        _hist_sigma->fill(j.absrap(), j.pT());
41      }
42    }
43
44    /// Finalize
45    void finalize() {
46      scale(_hist_sigma, crossSection()/picobarn/sumOfWeights()/2.0);
47      divByGroupWidth(_hist_sigma);
48    }
49
50    /// @}
51
52
53  private:
54
55    Histo1DGroupPtr _hist_sigma;
56
57  };
58
59
60
61  RIVET_DECLARE_ALIASED_PLUGIN(CMS_2011_I902309, CMS_2011_S9086218);
62
63}