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

namespace Rivet {


  /// D0 dijet invariant mass measurement
  class D0_2010_S8566488 : public Analysis {
  public:

    RIVET_DEFAULT_ANALYSIS_CTOR(D0_2010_S8566488);


    /// @name Analysis methods
    /// @{

    /// Book histograms and initialise projections before the run
    void init() {

      FinalState fs;
      FastJets conefinder(fs, FastJets::D0ILCONE, 0.7);
      declare(conefinder, "ConeFinder");

      {Histo1DPtr tmp; _h_m_dijet.add(0.0, 0.4, book(tmp, 1, 1, 1));}
      {Histo1DPtr tmp; _h_m_dijet.add(0.4, 0.8, book(tmp, 2, 1, 1));}
      {Histo1DPtr tmp; _h_m_dijet.add(0.8, 1.2, book(tmp, 3, 1, 1));}
      {Histo1DPtr tmp; _h_m_dijet.add(1.2, 1.6, book(tmp, 4, 1, 1));}
      {Histo1DPtr tmp; _h_m_dijet.add(1.6, 2.0, book(tmp, 5, 1, 1));}
      {Histo1DPtr tmp; _h_m_dijet.add(2.0, 2.4, book(tmp, 6, 1, 1));}
    }


    /// Perform the per-event analysis
    void analyze(const Event& e) {
      const Jets& jets = apply<JetAlg>(e, "ConeFinder").jetsByPt(40.0*GeV);
      if (jets.size() < 2) vetoEvent;

      FourMomentum j0(jets[0].momentum());
      FourMomentum j1(jets[1].momentum());
      double ymax = std::max(j0.absrap(), j1.absrap());
      double mjj = FourMomentum(j0+j1).mass();

      _h_m_dijet.fill(ymax, mjj/TeV);
    }


    /// Normalise histograms etc., after the run
    void finalize() {
      _h_m_dijet.scale(crossSection()/sumOfWeights(), this);
    }

    /// @}


  private:

    /// @name Histograms
    /// @{
    BinnedHistogram _h_m_dijet;
    /// @}

  };



  RIVET_DECLARE_ALIASED_PLUGIN(D0_2010_S8566488, D0_2010_I846483);

}