rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

D0_2010_S8566488

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_S8566488.cc
 1// -*- C++ -*-
 2#include "Rivet/Analysis.hh"
 3#include "Rivet/Tools/BinnedHistogram.hh"
 4#include "Rivet/Projections/FinalState.hh"
 5#include "Rivet/Projections/FastJets.hh"
 6
 7namespace Rivet {
 8
 9
10  /// D0 dijet invariant mass measurement
11  class D0_2010_S8566488 : public Analysis {
12  public:
13
14    RIVET_DEFAULT_ANALYSIS_CTOR(D0_2010_S8566488);
15
16
17    /// @name Analysis methods
18    /// @{
19
20    /// Book histograms and initialise projections before the run
21    void init() {
22
23      FinalState fs;
24      FastJets conefinder(fs, FastJets::D0ILCONE, 0.7);
25      declare(conefinder, "ConeFinder");
26
27      {Histo1DPtr tmp; _h_m_dijet.add(0.0, 0.4, book(tmp, 1, 1, 1));}
28      {Histo1DPtr tmp; _h_m_dijet.add(0.4, 0.8, book(tmp, 2, 1, 1));}
29      {Histo1DPtr tmp; _h_m_dijet.add(0.8, 1.2, book(tmp, 3, 1, 1));}
30      {Histo1DPtr tmp; _h_m_dijet.add(1.2, 1.6, book(tmp, 4, 1, 1));}
31      {Histo1DPtr tmp; _h_m_dijet.add(1.6, 2.0, book(tmp, 5, 1, 1));}
32      {Histo1DPtr tmp; _h_m_dijet.add(2.0, 2.4, book(tmp, 6, 1, 1));}
33    }
34
35
36    /// Perform the per-event analysis
37    void analyze(const Event& e) {
38      const Jets& jets = apply<JetAlg>(e, "ConeFinder").jetsByPt(40.0*GeV);
39      if (jets.size() < 2) vetoEvent;
40
41      FourMomentum j0(jets[0].momentum());
42      FourMomentum j1(jets[1].momentum());
43      double ymax = std::max(j0.absrap(), j1.absrap());
44      double mjj = FourMomentum(j0+j1).mass();
45
46      _h_m_dijet.fill(ymax, mjj/TeV);
47    }
48
49
50    /// Normalise histograms etc., after the run
51    void finalize() {
52      _h_m_dijet.scale(crossSection()/sumOfWeights(), this);
53    }
54
55    /// @}
56
57
58  private:
59
60    /// @name Histograms
61    /// @{
62    BinnedHistogram _h_m_dijet;
63    /// @}
64
65  };
66
67
68
69  RIVET_DECLARE_ALIASED_PLUGIN(D0_2010_S8566488, D0_2010_I846483);
70
71}