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
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
namespace Rivet {
/// CMS 4-jet production at 7 TeV
class CMS_2013_I1273574 : public Analysis {
public :
/// Constructor
CMS_2013_I1273574()
: Analysis("CMS_2013_I1273574" )
{ }
/// Book histograms and initialise projections before the run
void init () {
const FinalState cnfs((Cuts:: etaIn(- 4.7 , 4.7 )));
declare(FastJets(cnfs, FastJets:: ANTIKT, 0.5 ), "Jets" );
// Modified to match the HEPDATA record.
// eta of highest pT jet
//book(_h_jetetas[0] ,1,1,1);
book(_h_jetetas[0 ] ,6 ,1 ,1 );
// pt of the highest pT jet
book(_h_jetpts[0 ] ,2 ,1 ,1 );
//book(_h_DeltaS ,3,1,1);
book(_h_DeltaS ,12 ,1 ,1 );
//book(_h_DeltaPhiSoft ,4,1,1);
book(_h_DeltaPhiSoft ,10 ,1 ,1 );
//book(_h_DeltaPtRelSoft ,5,1,1);
book(_h_DeltaPtRelSoft ,11 ,1 ,1 );
// eta and pT of 3rd highest pT jet
//book(_h_jetetas[2] ,6,1,1);
//book(_h_jetpts[2] ,7,1,1);
book(_h_jetetas[2 ] ,8 ,1 ,1 );
book(_h_jetpts[2 ] ,4 ,1 ,1 );
// eta and pT of 4th highest pT jet
//book(_h_jetetas[3] ,8,1,1);
//book(_h_jetpts[3] ,9,1,1);
book(_h_jetetas[3 ] ,9 ,1 ,1 );
book(_h_jetpts[3 ] ,5 ,1 ,1 );
// eta and pT of 2nd highest pT jet
//book(_h_jetetas[1] ,10,1,1);
//book(_h_jetpts[1] ,11,1,1);
book(_h_jetetas[1 ] ,7 ,1 ,1 );
book(_h_jetpts[1 ] ,3 ,1 ,1 );
}
/// Perform the per-event analysis
void analyze (const Event& event) {
/// @todo Use jetsByPt(ptGtr(20*GeV) & absetaIn(4.7)), then no need for the lower loop;
const Jets jets = apply< FastJets> (event, "Jets" ).jetsByPt(20 * GeV);
if (jets.size() < 4 ) vetoEvent;
// Ensure that there are exactly 4 jets > 20 GeV, with two above 50 GeV
Jets hardjets, alljets;
for (const Jet& j : jets) {
if (j.abseta() > 4.7 ) continue ;
if (j.pT() > 50 * GeV) hardjets.push_back(j);
if (j.pT() > 20 * GeV) alljets.push_back(j);
}
if (hardjets.size() < 2 || alljets.size() != 4 ) vetoEvent;
const double weight = 1.0 ;
// Histogram pT and eta of all 4 jets
for (size_t i = 0 ; i < 4 ; ++ i) {
_h_jetpts[i]-> fill(alljets[i].pT()/ GeV, weight);
_h_jetetas[i]-> fill(alljets[i].eta(), weight);
}
// Create vector sums of the hard and soft pairs of jets
const FourMomentum p12 = alljets[0 ].momentum() + alljets[1 ].momentum();
const FourMomentum p34 = alljets[2 ].momentum() + alljets[3 ].momentum();
// Fill the delta(phi) between the soft jets
const double dphisoft = deltaPhi(alljets[2 ], alljets[3 ]);
_h_DeltaPhiSoft-> fill(dphisoft, weight);
// Fill the pT balance between the soft jets
const double ptbalanceSoft = p34.pT() / (alljets[2 ].pT() + alljets[3 ].pT());
_h_DeltaPtRelSoft-> fill(ptbalanceSoft, weight);
// Fill the azimuthal angle difference between the two jet pairs
const double p12p34_trans = p12.px()* p34.px() + p12.py()* p34.py();
const double DeltaS = acos( p12p34_trans / p12.pT() / p34.pT() );
_h_DeltaS-> fill(DeltaS, weight);
}
/// Normalise histograms (mostly to cross-section)
void finalize () {
const double invlumi = crossSection()/ picobarn/ sumOfWeights();
for (size_t i = 0 ; i < 4 ; ++ i) {
scale(_h_jetpts[i], invlumi);
scale(_h_jetetas[i], invlumi);
}
normalize(_h_DeltaPtRelSoft);
normalize(_h_DeltaPhiSoft);
normalize(_h_DeltaS);
}
private :
Histo1DPtr _h_jetpts[4 ], _h_jetetas[4 ];
Histo1DPtr _h_DeltaS, _h_DeltaPhiSoft, _h_DeltaPtRelSoft;
};
RIVET_DECLARE_PLUGIN(CMS_2013_I1273574);
}