Rivet analyses referenceCDF_2000_I511377Inclusive jet cross section differential in dijet massExperiment: CDF (Tevatron Run 1) Inspire ID: 511377 Status: VALIDATED Authors:
Beam energies: (900.0, 900.0) GeV Run details:
Measurement of the cross section for production of two or more jets as a function of dijet mass in the range 180 to 1000 GeV. It is based on an integrated luminosity of $86 \mathrm{pb}^{-1}$. Source code: CDF_2000_I511377.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 /// @brief CDF dijet cross-section, differential in dijet mass
10 class CDF_2000_I511377 : public Analysis {
11 public:
12
13 RIVET_DEFAULT_ANALYSIS_CTOR(CDF_2000_I511377);
14
15
16 /// @name Analysis methods
17 /// @{
18
19 /// Book histograms and initialise projections before the run
20 void init() {
21 FinalState fs((Cuts::etaIn(-4.2, 4.2)));
22 declare(FastJets(fs, JetAlg::CDFJETCLU, 0.7), "Jets");
23 book(_h_mjj ,1, 1, 1);
24 }
25
26
27 /// Perform the per-event analysis
28 void analyze(const Event& event) {
29 Jets jets = apply<FastJets>(event, "Jets").jets(cmpMomByEt);
30 if (jets.size() < 2) vetoEvent;
31 FourMomentum jet1 = jets[0].momentum();
32 FourMomentum jet2 = jets[1].momentum();
33 double eta1 = jet1.eta();
34 double eta2 = jet2.eta();
35 if (fabs(eta1) > 2.0 || fabs(eta2) > 2.0) vetoEvent;
36 if (fabs(tanh((eta1-eta2)/2)) > 2.0/3.0) vetoEvent;
37 double mjj = FourMomentum(jet1+jet2).mass()/GeV;
38 if (mjj < 180) vetoEvent;
39 _h_mjj->fill(mjj);
40 }
41
42
43 /// Normalise histograms etc., after the run
44 void finalize() {
45 scale(_h_mjj, crossSection()/picobarn/sumOfWeights());
46 }
47
48 /// @}
49
50
51 private:
52
53 /// Histogram
54 Histo1DPtr _h_mjj;
55
56 };
57
58
59
60 RIVET_DECLARE_ALIASED_PLUGIN(CDF_2000_I511377, CDF_2000_S4266730);
61
62}
|