rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

ATLAS_2014_I1325553

Measurement of the inclusive jet cross-section at 7 TeV
Experiment: ATLAS (LHC)
Inspire ID: 1325553
Status: VALIDATED
Authors:
  • Vojtech Pleskot
References: Beams: p+ p+
Beam energies: (3500.0, 3500.0) GeV
Run details:
  • p p -> jet + X. $\sqrt{s} = 7$ TeV.

Measurement of the inclusive jet cross-section in proton--proton collisions at a centre-of-mass energy of 7 TeV using a data set corresponding to an integrated luminosity of 4.5/fb collected with the ATLAS detector at the Large Hadron Collider in 2011. Jets are identified using the anti-$k_t$ algorithm with radius parameter values of 0.4 and 0.6. The double-differential cross-sections are presented as a function of the jet transverse momentum and the jet rapidity, covering jet transverse momenta from 100 GeV to 2 TeV.

Source code: ATLAS_2014_I1325553.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  /// Jet mass as a function of ystar
10  class ATLAS_2014_I1325553 : public Analysis {
11  public:
12
13    /// Constructor
14    RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2014_I1325553);
15
16    /// @name Analysis methods
17    /// @{
18
19    /// Book histograms and initialise projections before the run
20    void init() {
21
22      const FinalState fs;
23      declare(fs,"FinalState");
24
25      FastJets fj04(fs,  JetAlg::ANTIKT, 0.4);
26      fj04.useInvisibles();
27      declare(fj04, "AntiKT04");
28
29      FastJets fj06(fs,  JetAlg::ANTIKT, 0.6);
30      fj06.useInvisibles();
31      declare(fj06, "AntiKT06");
32
33      const vector<double> ybins{0.0, 0.5, 1.0, 1.5, 2.0, 2.5, 3.0};
34
35      size_t ptDsOffset = 0;
36      for (size_t alg = 0; alg < 2; ++alg) {
37        book(_pt[alg], ybins);
38        for (auto& b : _pt[alg]->bins()) {
39          book(b, 1+ptDsOffset, 1, 1);
40          ++ptDsOffset;
41        }
42      }
43    }
44
45
46    /// Perform the per-event analysis
47    void analyze(const Event& event) {
48      Jets jetAr[2];
49      jetAr[AKT4] = apply<FastJets>(event, "AntiKT04").jetsByPt(Cuts::pT > 100*GeV && Cuts::absrap < 3.0);
50      jetAr[AKT6] = apply<FastJets>(event, "AntiKT06").jetsByPt(Cuts::pT > 100*GeV && Cuts::absrap < 3.0);
51
52      // Loop over jet "radii" used in analysis
53      for (size_t alg = 0; alg < 2; ++alg) {
54
55        // fill the 1D pt histograms with all the jets passing the cuts
56        for (const Jet& jet : jetAr[alg]) {
57          const double absrap = jet.absrap();
58          if (absrap < 3.0) {
59            const double pt = jet.pT();
60            if (pt/GeV > 100*GeV) {
61              _pt[alg]->fill(absrap, pt/GeV);
62            }
63          }
64        }
65      }
66    }
67
68
69    /// Normalise histograms etc., after the run
70    void finalize() {
71
72      /// Print summary info
73      const double xs_norm_factor = 0.5*crossSection()/picobarn/sumOfWeights();
74      scale(_pt, xs_norm_factor);
75      divByGroupWidth(_pt);
76    }
77
78    /// @}
79
80
81  private:
82
83    // Data members like post-cuts event weight counters go here
84    enum Alg { AKT4=0, AKT6=1 };
85
86    /// The inclusive jet spectrum binned in rapidity for akt6 and akt4 jets (array index is jet type from enum above)
87    Histo1DGroupPtr _pt[2];
88
89  };
90
91
92  RIVET_DECLARE_PLUGIN(ATLAS_2014_I1325553);
93
94}