rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

D0_2010_I846483

Dijet invariant mass
Experiment: D0 (Tevatron Run 2)
Inspire ID: 846483
Status: VALIDATED
Authors:
  • Frank Siegert
References: Beams: p- p+
Beam energies: (980.0, 980.0) GeV
Run details:
  • $p \bar{p} \to$ jets at 1960 GeV. Analysis needs two hard jets above 40 GeV.

The inclusive dijet production double differential cross section as a function of the dijet invariant mass and of the largest absolute rapidity ($|y|_\text{max}$) of the two jets with the largest transverse momentum in an event is measured using 0.7 fb$^{-1}$ of data. The measurement is performed in six rapidity regions up to $|y|_\text{max}=2.4$.

Source code: D0_2010_I846483.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  /// D0 dijet invariant mass measurement
10  class D0_2010_I846483 : public Analysis {
11  public:
12
13    RIVET_DEFAULT_ANALYSIS_CTOR(D0_2010_I846483);
14
15
16    /// @name Analysis methods
17    /// @{
18
19    /// Book histograms and initialise projections before the run
20    void init() {
21
22      FinalState fs;
23      FastJets conefinder(fs, JetAlg::D0ILCONE, 0.7);
24      declare(conefinder, "ConeFinder");
25
26      book(_h_m_dijet, {0., 0.4, 0.8, 1.2, 1.6, 2., 2.4});
27      for (auto& b : _h_m_dijet->bins()) {
28        book(b, b.index(), 1, 1);
29      }
30    }
31
32
33    /// Perform the per-event analysis
34    void analyze(const Event& e) {
35      const Jets& jets = apply<JetFinder>(e, "ConeFinder").jetsByPt(Cuts::pT > 40.0*GeV);
36      if (jets.size() < 2) vetoEvent;
37
38      const double ymax = std::max(jets[0].absrap(), jets[1].absrap());
39      const double mjj = FourMomentum(jets[0].mom() + jets[1].mom()).mass();
40
41      _h_m_dijet->fill(ymax, mjj/TeV);
42    }
43
44
45    /// Normalise histograms etc., after the run
46    void finalize() {
47      scale(_h_m_dijet, crossSection()/picobarn/sumOfWeights());
48      divByGroupWidth(_h_m_dijet);
49    }
50
51    /// @}
52
53
54  private:
55
56    /// @name Histograms
57    /// @{
58    Histo1DGroupPtr _h_m_dijet;
59    /// @}
60
61  };
62
63
64
65  RIVET_DECLARE_ALIASED_PLUGIN(D0_2010_I846483, D0_2010_S8566488);
66
67}