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
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
| // -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/HeavyHadrons.hh"
#include "Rivet/Tools/BinnedHistogram.hh"
namespace Rivet {
/// @brief ATLAS inclusive b-jet pT spectrum, di-jet mass and di-jet chi
class ATLAS_2011_I930220: public Analysis {
public:
ATLAS_2011_I930220()
: Analysis("ATLAS_2011_I930220")
{ }
void init() {
FinalState fs((Cuts::etaIn(-3.5, 3.5)));
declare(fs, "FinalState");
FastJets fj(fs, FastJets::ANTIKT, 0.4);
fj.useInvisibles();
declare(fj, "Jets");
declare(HeavyHadrons(Cuts::abseta < 3.5 && Cuts::pT > 5*GeV), "BHadrons");
double ybins[] = { 0.0, 0.3, 0.8, 1.2, 2.1 };
for (size_t i = 0; i < 4; ++i) {
Histo1DPtr tmp;
_bjetpT_SV0.add(ybins[i], ybins[i+1], book(tmp, i+1, 1, 1));
}
book(_bjetpT_SV0_All ,5, 1, 1);
book(_bjetpT_pTRel ,6, 1, 1);
book(_dijet_mass ,7, 1, 1);
book(_dijet_phi ,8, 1, 1);
book(_dijet_chi_110_370 ,9, 1, 1);
book(_dijet_chi_370_850 ,10, 1, 1);
book(_chiCounter1, "_chiCounter1");
book(_chiCounter2, "_chiCounter2");
book(_phiCounter, "_phiCounter1");
}
void analyze(const Event& evt) {
const Particles& bHadrons = apply<HeavyHadrons>(evt, "BHadrons").bHadrons();
const Jets& jets = apply<JetAlg>(evt, "Jets").jetsByPt(15*GeV);
FourMomentum leadingJet, subleadingJet;
int leadJet = 0, subJet = 0;
for (const Jet& j : jets) {
bool hasB = false;
for (const Particle& b : bHadrons)
if (deltaR(j, b) < 0.3) { hasB = true; break; }
// Identify and classify the leading and subleading jets
if (j.absrap() < 2.1) { ///< Move this into the jets defn
if (!leadJet) {
leadingJet = j.momentum();
leadJet = (hasB && j.pT() > 40*GeV) ? 2 : 1;
}
else if (leadJet && !subJet) {
subleadingJet = j.momentum();
subJet = (hasB && j.pT() > 40*GeV) ? 2 : 1;
}
if (hasB) {
_bjetpT_SV0.fill(j.absrap(), j.pT()/GeV);
_bjetpT_SV0_All->fill(j.pT()/GeV);
_bjetpT_pTRel->fill(j.pT()/GeV);
}
}
}
// Di-b-jet plots require both the leading and subleading jets to be b-tagged and have pT > 40 GeV
if (leadJet == 2 && subJet == 2) {
const double mass = FourMomentum( leadingJet + subleadingJet ).mass();
_dijet_mass->fill(mass/GeV);
// Plot dphi for high-mass di-b-jets
if (mass > 110*GeV) {
_phiCounter->fill();
const double d_phi = deltaPhi( leadingJet.phi(), subleadingJet.phi() );
_dijet_phi->fill(fabs(d_phi));
}
// Plot chi for low y_boost di-b-jets (in two high-mass bins)
const double y_boost = 0.5 * (leadingJet.rapidity() + subleadingJet.rapidity());
const double chi = exp( fabs( leadingJet.rapidity() - subleadingJet.rapidity() ) );
if ( fabs(y_boost) < 1.1 ) {
if (inRange(mass/GeV, 110, 370)) {
_chiCounter1->fill();
_dijet_chi_110_370->fill(chi);
} else if (inRange(mass/GeV, 370, 850)) {
_chiCounter2->fill();
_dijet_chi_370_850->fill(chi);
}
}
}
}
void finalize() {
// Normalizing to cross-section and mass
// Additional factors represent the division by rapidity
const double xsec = crossSectionPerEvent()/(picobarn);
const double chiScale1 = 1 / dbl(*_chiCounter1) / 260.0;
const double chiScale2 = 1 / dbl(*_chiCounter2) / 480.0;
const double phiScale = 1 / dbl(*_phiCounter);
_bjetpT_SV0.scale(xsec/2, this);
scale(_bjetpT_SV0_All, xsec);
scale(_bjetpT_pTRel, xsec);
scale(_dijet_mass, xsec);
scale(_dijet_phi, phiScale );
scale(_dijet_chi_110_370, chiScale1);
scale(_dijet_chi_370_850, chiScale2);
}
private:
BinnedHistogram _bjetpT_SV0;
Histo1DPtr _bjetpT_SV0_All;
Histo1DPtr _bjetpT_pTRel;
Histo1DPtr _dijet_mass;
Histo1DPtr _dijet_phi;
Histo1DPtr _dijet_chi_110_370;
Histo1DPtr _dijet_chi_370_850;
CounterPtr _chiCounter1;
CounterPtr _chiCounter2;
CounterPtr _phiCounter;
};
// The hook for the plugin system
RIVET_DECLARE_PLUGIN(ATLAS_2011_I930220);
}
|