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
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
| // -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/PromptFinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/Sphericity.hh"
#include "Rivet/Projections/SmearedParticles.hh"
#include "Rivet/Projections/SmearedJets.hh"
#include "Rivet/Projections/SmearedMET.hh"
#include "Rivet/Tools/Cutflow.hh"
namespace Rivet {
/// @brief ATLAS 2016 1-lepton + many jets SUSY search, from 14.8/fb CONF note
class ATLAS_2016_CONF_2016_094 : public Analysis {
public:
/// Constructor
RIVET_DEFAULT_ANALYSIS_CTOR(ATLAS_2016_CONF_2016_094);
/// @name Analysis methods
//@{
/// Book histograms and initialise projections before the run
void init() {
// Initialise and register projections
FinalState calofs(Cuts::abseta < 4.9);
FastJets fj(calofs, FastJets::ANTIKT, 0.4);
declare(fj, "TruthJets");
declare(SmearedJets(fj, JET_SMEAR_ATLAS_RUN2, [](const Jet& j) {
if (j.abseta() > 2.5) return 0.;
return j.bTagged(Cuts::pT > 5*GeV) ? 0.80 :
j.cTagged(Cuts::pT > 5*GeV) ? 1/6. : 1/106.; }), "Jets");
// MissingMomentum mm(calofs);
// declare(mm, "TruthMET");
// declare(SmearedMET(mm, MET_SMEAR_ATLAS_RUN2), "MET");
FinalState es(Cuts::abspid == PID::ELECTRON && Cuts::abseta < 2.47 && !Cuts::absetaIn(1.37, 1.52) && Cuts::pT > 10*GeV);
declare(es, "TruthElectrons");
declare(SmearedParticles(es, ELECTRON_RECOEFF_ATLAS_RUN2, ELECTRON_SMEAR_ATLAS_RUN2), "Electrons");
FinalState mus(Cuts::abspid == PID::MUON && Cuts::abseta < 2.4 && Cuts::pT > 10*GeV);
declare(mus, "TruthMuons");
declare(SmearedParticles(mus, MUON_EFF_ATLAS_RUN2, MUON_SMEAR_ATLAS_RUN2), "Muons");
// Book histograms/counters
book(_h_08j40_0b,"08j40_0b");
book(_h_09j40_0b,"09j40_0b");
book(_h_10j40_0b,"10j40_0b");
book(_h_08j40_3b,"08j40_3b");
book(_h_09j40_3b,"09j40_3b");
book(_h_10j40_3b,"10j40_3b");
book(_h_08j60_0b,"08j60_0b");
book(_h_09j60_0b,"09j60_0b");
book(_h_10j60_0b,"10j60_0b");
book(_h_08j60_3b,"08j60_3b");
book(_h_09j60_3b,"09j60_3b");
book(_h_10j60_3b,"10j60_3b");
}
/// Perform the per-event analysis
void analyze(const Event& event) {
// Get baseline electrons, muons, and jets
// NB. for electrons, we don't apply the loose ID here, since we don't want to double-count effs with later use of tight ID
Particles elecs = apply<ParticleFinder>(event, "Electrons").particles();
Particles muons = apply<ParticleFinder>(event, "Muons").particles();
Jets jets = apply<JetAlg>(event, "Jets").jetsByPt(Cuts::pT > 20*GeV && Cuts::abseta < 2.4);
ifilter_select(jets, JetEffFilter([](const Jet& j) { return j.pT() > 60*GeV ? 1.0 : 0.94; }));
// Jet/electron/muon overlap removal and selection
// Remove any untagged jet within dR = 0.2 of an electron
for (const Particle& e : elecs)
ifilter_discard(jets, [&](const Jet& j) { return !j.bTagged(Cuts::pT > 5*GeV) && deltaR(e, j, RAPIDITY) < 0.2; });
// Remove any untagged low-multiplicity/muon-dominated jet within dR = 0.4 of a muon
for (const Particle& m : muons)
ifilter_discard(jets, [&](const Jet& j) {
if (j.bTagged(Cuts::pT > 5*GeV)) return false; /// @note A different b-tag working point, 85%, was actually used here *sigh*
if (deltaR(m, j, RAPIDITY) > 0.4) return false;
if (j.particles(Cuts::abscharge != 0).size() < 3) return true;
return m.pT()/j.pT() > 0.5;
});
// Removing leptons within dR = 0.4 of remaining jets
for (const Jet& j : jets) {
ifilter_discard(elecs, deltaRLess(j, 0.4, RAPIDITY));
ifilter_discard(muons, deltaRLess(j, 0.4, RAPIDITY));
}
// Signal jet and lepton selection
const Jets sigjets40 = filter_select(jets, Cuts::pT > 40*GeV);
const Jets sigjets60 = filter_select(sigjets40, Cuts::pT > 60*GeV);
const Jets sigbjets40 = filter_select(sigjets40, [](const Jet& j) { return j.bTagged(Cuts::pT > 5*GeV); });
const Jets sigbjets60 = filter_select(sigjets60, [](const Jet& j) { return j.bTagged(Cuts::pT > 5*GeV); });
const Particles sigmuons = filter_select(muons, Cuts::pT > 35*GeV);
Particles sigelecs = filter_select(elecs, Cuts::pT > 35*GeV);
ifilter_select(sigelecs, ParticleEffFilter(ELECTRON_EFF_ATLAS_RUN2_TIGHT));
//////////////////
// Event selection cuts
if (sigelecs.size() + sigmuons.size() != 1) vetoEvent;
const Particle siglepton = sigelecs.empty() ? sigmuons.front() : sigelecs.front();
/// @note The note describes Nj = 5, 6, 7, 8, 9, >= 10 and Nb = 0, 1, 2, 3, >= 4 = 30 2D bins
/// for each jet cut... but only provides data for six Nj = >= 8, 9, 10, Nb = 0, >= 3 bins.
/// We just implement the latter for now.
// Fill counters
if (sigjets40.size() >= 8 && sigbjets40.empty()) _h_08j40_0b->fill();
if (sigjets40.size() >= 9 && sigbjets40.empty()) _h_09j40_0b->fill();
if (sigjets40.size() >= 10 && sigbjets40.empty()) _h_10j40_0b->fill();
if (sigjets40.size() >= 8 && sigbjets40.size() >= 3) _h_08j40_3b->fill();
if (sigjets40.size() >= 9 && sigbjets40.size() >= 3) _h_09j40_3b->fill();
if (sigjets40.size() >= 10 && sigbjets40.size() >= 3) _h_10j40_3b->fill();
if (sigjets60.size() >= 8 && sigbjets60.empty()) _h_08j60_0b->fill();
if (sigjets60.size() >= 9 && sigbjets60.empty()) _h_09j60_0b->fill();
if (sigjets60.size() >= 10 && sigbjets60.empty()) _h_10j60_0b->fill();
if (sigjets60.size() >= 8 && sigbjets60.size() >= 3) _h_08j60_3b->fill();
if (sigjets60.size() >= 9 && sigbjets60.size() >= 3) _h_09j60_3b->fill();
if (sigjets60.size() >= 10 && sigbjets60.size() >= 3) _h_10j60_3b->fill();
}
/// Normalise counters after the run
void finalize() {
const double sf = 14.8*crossSection()/femtobarn/sumOfWeights();
scale(_h_08j40_0b, sf); scale(_h_09j40_0b, sf); scale(_h_10j40_0b, sf);
scale(_h_08j40_3b, sf); scale(_h_09j40_3b, sf); scale(_h_10j40_3b, sf);
scale(_h_08j60_0b, sf); scale(_h_09j60_0b, sf); scale(_h_10j60_0b, sf);
scale(_h_08j60_3b, sf); scale(_h_09j60_3b, sf); scale(_h_10j60_3b, sf);
}
//@}
private:
/// @name Histograms
//@{
CounterPtr _h_08j40_0b, _h_09j40_0b, _h_10j40_0b, _h_08j40_3b, _h_09j40_3b, _h_10j40_3b;
CounterPtr _h_08j60_0b, _h_09j60_0b, _h_10j60_0b, _h_08j60_3b, _h_09j60_3b, _h_10j60_3b;
//@}
};
// The hook for the plugin system
RIVET_DECLARE_PLUGIN(ATLAS_2016_CONF_2016_094);
}
|