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
|
#include "DQMServices/Core/interface/DQMEDAnalyzer.h"
#include "DataFormats/PatCandidates/interface/Jet.h"
#include "DQMOffline/RecoB/interface/JetTagPlotter.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/Event.h"
/** \class MiniAODTaggerAnalyzer
*
* Tagger analyzer to run on MiniAOD
*
*/
class MiniAODTaggerAnalyzer : public DQMEDAnalyzer {
public:
explicit MiniAODTaggerAnalyzer(const edm::ParameterSet& pSet);
~MiniAODTaggerAnalyzer() override = default;
void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
private:
void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
typedef std::vector<std::string> vstring;
// using JetTagPlotter object for all the hard work ;)
std::unique_ptr<JetTagPlotter> jetTagPlotter_;
const edm::EDGetTokenT<std::vector<pat::Jet> > jetToken_;
const edm::ParameterSet discrParameters_;
const std::string folder_;
const vstring discrNumerator_;
const vstring discrDenominator_;
const int mclevel_;
const bool doCTagPlots_;
const bool dodifferentialPlots_;
const double discrCut_;
const bool etaActive_;
const double etaMin_;
const double etaMax_;
const bool ptActive_;
const double ptMin_;
const double ptMax_;
};
MiniAODTaggerAnalyzer::MiniAODTaggerAnalyzer(const edm::ParameterSet& pSet)
: jetToken_(consumes<std::vector<pat::Jet> >(pSet.getParameter<edm::InputTag>("JetTag"))),
discrParameters_(pSet.getParameter<edm::ParameterSet>("parameters")),
folder_(pSet.getParameter<std::string>("folder")),
discrNumerator_(pSet.getParameter<vstring>("numerator")),
discrDenominator_(pSet.getParameter<vstring>("denominator")),
mclevel_(pSet.getParameter<int>("MClevel")),
doCTagPlots_(pSet.getParameter<bool>("CTagPlots")),
dodifferentialPlots_(pSet.getParameter<bool>("differentialPlots")),
discrCut_(pSet.getParameter<double>("discrCut")),
etaActive_(pSet.getParameter<bool>("etaActive")),
etaMin_(pSet.getParameter<double>("etaMin")),
etaMax_(pSet.getParameter<double>("etaMax")),
ptActive_(pSet.getParameter<bool>("ptActive")),
ptMin_(pSet.getParameter<double>("ptMin")),
ptMax_(pSet.getParameter<double>("ptMax")) {}
void MiniAODTaggerAnalyzer::bookHistograms(DQMStore::IBooker& ibook, edm::Run const& run, edm::EventSetup const& es) {
jetTagPlotter_ = std::make_unique<JetTagPlotter>(folder_,
EtaPtBin(etaActive_, etaMin_, etaMax_, ptActive_, ptMin_, ptMax_),
discrParameters_,
mclevel_,
false,
ibook,
doCTagPlots_,
dodifferentialPlots_,
discrCut_);
}
void MiniAODTaggerAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
edm::Handle<std::vector<pat::Jet> > jetCollection;
iEvent.getByToken(jetToken_, jetCollection);
// Loop over the pat::Jets
for (std::vector<pat::Jet>::const_iterator jet = jetCollection->begin(); jet != jetCollection->end(); ++jet) {
// apply basic jet cuts
if (jet->pt() > ptMin_ && std::abs(jet->eta()) < etaMax_) {
// fill numerator
float numerator = 0;
for (const auto& discrLabel : discrNumerator_) {
numerator += jet->bDiscriminator(discrLabel);
}
// fill denominator
float denominator;
if (discrDenominator_.empty()) {
denominator = 1; // no division performed
} else {
denominator = 0;
for (const auto& discrLabel : discrDenominator_) {
denominator += jet->bDiscriminator(discrLabel);
}
}
const float jec = 1.; // JEC not implemented!
// only add to histograms when discriminator values are valid
if (numerator >= 0 && denominator > 0) {
const reco::Jet& recoJet = *jet;
if (jetTagPlotter_->etaPtBin().inBin(recoJet, jec)) {
jetTagPlotter_->analyzeTag(recoJet, jec, numerator / denominator, jet->partonFlavour());
}
}
}
}
// fill JetMultiplicity (once per event)
if (mclevel_ > 0) {
jetTagPlotter_->analyzeTag(1.);
} else {
jetTagPlotter_->analyzeTag();
}
}
//define this as a plug-in
DEFINE_FWK_MODULE(MiniAODTaggerAnalyzer);
|