Rivet analyses referenceD0_2009_I824127Dijet angular distributionsExperiment: D0 (Tevatron Run 2) Inspire ID: 824127 Status: VALIDATED Authors:
Beam energies: (980.0, 980.0) GeV Run details:
Dijet angular distributions in different bins of dijet mass from 0.25 TeV to above 1.1 TeV in $p \bar{p}$ collisions at $\sqrt{s} = 1.96$ TeV, based on an integrated luminosity of 0.7 fb$^{-1}$. Source code: D0_2009_I824127.cc 1// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/FinalState.hh"
4#include "Rivet/Projections/ChargedFinalState.hh"
5#include "Rivet/Projections/FastJets.hh"
6
7namespace Rivet {
8
9
10 /// D0 dijet angular distributions
11 class D0_2009_I824127 : public Analysis {
12 public:
13
14 RIVET_DEFAULT_ANALYSIS_CTOR(D0_2009_I824127);
15
16
17 /// @name Analysis methods
18 /// @{
19
20 // Book histograms
21 void init() {
22 FinalState fs;
23 FastJets conefinder(fs, JetAlg::D0ILCONE, 0.7);
24 declare(conefinder, "ConeFinder");
25
26 book(_h_chi_dijet, {250., 300., 400., 500., 600., 700., 800., 900., 1000., 1100., 1960.});
27 for (auto& b : _h_chi_dijet->bins()) {
28 book(b, b.index(), 1, 1);
29 }
30 }
31
32
33 /// Do the analysis
34 void analyze(const Event & e) {
35
36 const Jets& jets = apply<JetFinder>(e, "ConeFinder").jetsByPt();
37 if (jets.size() < 2) vetoEvent;
38
39 FourMomentum j0(jets[0].momentum());
40 FourMomentum j1(jets[1].momentum());
41 double y0 = j0.rapidity();
42 double y1 = j1.rapidity();
43
44 if (fabs(y0+y1)>2) vetoEvent;
45
46 double mjj = FourMomentum(j0+j1).mass();
47 double chi = exp(fabs(y0-y1));
48 if(chi<16.) _h_chi_dijet->fill(mjj, chi);
49 }
50
51
52 /// Finalize
53 void finalize() {
54 normalize(_h_chi_dijet);
55 }
56
57 /// @}
58
59
60 private:
61
62 /// @name Histograms
63 /// @{
64 Histo1DGroupPtr _h_chi_dijet;
65 /// @}
66
67 };
68
69
70
71 RIVET_DECLARE_ALIASED_PLUGIN(D0_2009_I824127, D0_2009_S8320160);
72
73}
|