Rivet Analyses Reference

CMS_2013_I1208923

Jet-pT and dijet mass at sqrt(s) = 7 TeV
Experiment: CMS (LHC)
Inspire ID: 1208923
Status: VALIDATED
Authors:
  • Markus Radziej
  • Samatha Dooling
No references listed
Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • Hard QCD process with a flat pT spectrum is recommended.

The single jet differential cross section has been measured as a function of the jet momentum and the dijet differential cross section as a function of the dijet mass. The data has been recorded by the CMS detector at the center of mass energy of 7 TeV. To reconstruct the jets, the anti-$k_T$ algorithm with a cone radius of $R = 0.7$ (AK7) has been used. The results are split into five rapidity bins, ranging from 0.0 to 2.5 with a width of 0.5 each.

Source code: CMS_2013_I1208923.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
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Tools/BinnedHistogram.hh"

namespace Rivet {

  // This analysis is a derived from the class Analysis:
  class CMS_2013_I1208923 : public Analysis {

  public:
    // Constructor
    CMS_2013_I1208923()
      : Analysis("CMS_2013_I1208923") {

    }

    // Book histograms and initialize projections:
    void init() {
      const FinalState fs;
      declare(fs, "FS");

      // Initialize the projections
      declare(FastJets(fs, FastJets::ANTIKT, 0.7), "Jets");

      // Book histograms
      {Histo1DPtr tmp; _h_sigma.add(0.0, 0.5,   book(tmp, 1, 1, 1));}
      {Histo1DPtr tmp; _h_sigma.add(0.5, 1.0,   book(tmp, 1, 1, 2));}
      {Histo1DPtr tmp; _h_sigma.add(1.0, 1.5,   book(tmp, 1, 1, 3));}
      {Histo1DPtr tmp; _h_sigma.add(1.5, 2.0,   book(tmp, 1, 1, 4));}
      {Histo1DPtr tmp; _h_sigma.add(2.0, 2.5,   book(tmp, 1, 1, 5));}
      
      {Histo1DPtr tmp; _h_invMass.add(0.0, 0.5, book(tmp, 2, 1, 1));}
      {Histo1DPtr tmp; _h_invMass.add(0.5, 1.0, book(tmp, 2, 1, 2));}
      {Histo1DPtr tmp; _h_invMass.add(1.0, 1.5, book(tmp, 2, 1, 3));}
      {Histo1DPtr tmp; _h_invMass.add(1.5, 2.0, book(tmp, 2, 1, 4));}
      {Histo1DPtr tmp; _h_invMass.add(2.0, 2.5, book(tmp, 2, 1, 5));}
    }

    // Analysis
    void analyze(const Event &event) {
      const double weight = 1.0;
      const FastJets &fJets = apply<FastJets>(event, "Jets");
      
      // Fill the jet pT spectra
      const Jets& jets = fJets.jetsByPt(Cuts::pt>100.*GeV && Cuts::absrap <2.5);
      for (const Jet &j : jets) {
        _h_sigma.fill(fabs(j.momentum().rapidity()), j.momentum().pT() / GeV, weight);
      }

      // Require two jets
      const Jets& dijets = fJets.jetsByPt(Cuts::pt>30.*GeV && Cuts::absrap < 2.5);
      if (dijets.size() > 1) {
        if (dijets[0].momentum().pT() / GeV > 60.) {
          // Fill the invariant mass histogram
          double ymax = max(dijets[0].momentum().absrapidity(), dijets[1].momentum().absrapidity());
          double invMass = FourMomentum(dijets[0].momentum() + dijets[1].momentum()).mass();
          _h_invMass.fill(fabs(ymax), invMass, weight);
        }
      } 

    }


    // Scale histograms by the production cross section
    void finalize() {
      _h_sigma.scale(  crossSection() / sumOfWeights() / 2.0, this);
      _h_invMass.scale(crossSection() / sumOfWeights() / 2.0, this);
    }

  private:
    BinnedHistogram _h_sigma;
    BinnedHistogram _h_invMass;
  };

  // This global object acts as a hook for the plugin system.
  RIVET_DECLARE_PLUGIN(CMS_2013_I1208923);
}