File indexing completed on 2024-04-06 12:24:32
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <iostream>
0021 #include <memory>
0022 #include <map>
0023 #include <vector>
0024 #include <string>
0025
0026
0027 #include "FWCore/Framework/interface/Frameworkfwd.h"
0028 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0029
0030 #include "FWCore/Framework/interface/Event.h"
0031 #include "FWCore/Framework/interface/MakerMacros.h"
0032 #include "FWCore/Framework/interface/EventSetup.h"
0033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0034 #include "FWCore/Framework/interface/ESHandle.h"
0035
0036
0037
0038
0039 #include "RecoBTag/Records/interface/BTagPerformanceRecord.h"
0040 #include "CondFormats/PhysicsToolsObjects/interface/BinningPointByMap.h"
0041 #include "RecoBTag/PerformanceDB/interface/BtagPerformance.h"
0042
0043 using namespace std;
0044
0045 class DumpBtagTable : public edm::one::EDAnalyzer<> {
0046 public:
0047 explicit DumpBtagTable(const edm::ParameterSet&);
0048
0049 private:
0050 string name;
0051 vector<string> measureName;
0052 vector<string> measureType;
0053 vector<edm::ESGetToken<BtagPerformance, BTagPerformanceRecord>> measureToken;
0054 void analyze(const edm::Event&, const edm::EventSetup&) final;
0055
0056
0057 };
0058
0059 DumpBtagTable::DumpBtagTable(const edm::ParameterSet& iConfig)
0060
0061 {
0062 measureName = iConfig.getParameter<vector<string>>("measureName");
0063 measureType = iConfig.getParameter<vector<string>>("measureType");
0064
0065 measureToken.reserve(measureName.size());
0066
0067 for (auto const& n : measureName) {
0068 measureToken.push_back(esConsumes<BtagPerformance, BTagPerformanceRecord>(edm::ESInputTag("", n)));
0069 }
0070 }
0071
0072
0073
0074
0075
0076
0077 void DumpBtagTable::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0078 map<string, PerformanceResult::ResultType> measureMap;
0079 measureMap["BTAGBEFF"] = PerformanceResult::BTAGBEFF;
0080 measureMap["BTAGBERR"] = PerformanceResult::BTAGBERR;
0081 measureMap["BTAGCEFF"] = PerformanceResult::BTAGCEFF;
0082 measureMap["BTAGCERR"] = PerformanceResult::BTAGCERR;
0083 measureMap["BTAGLEFF"] = PerformanceResult::BTAGLEFF;
0084 measureMap["BTAGLERR"] = PerformanceResult::BTAGLERR;
0085 measureMap["BTAGNBEFF"] = PerformanceResult::BTAGNBEFF;
0086 measureMap["BTAGNBERR"] = PerformanceResult::BTAGNBERR;
0087 measureMap["BTAGBEFFCORR"] = PerformanceResult::BTAGBEFFCORR;
0088 measureMap["BTAGBERRCORR"] = PerformanceResult::BTAGBERRCORR;
0089 measureMap["BTAGCEFFCORR"] = PerformanceResult::BTAGCEFFCORR;
0090 measureMap["BTAGCERRCORR"] = PerformanceResult::BTAGCERRCORR;
0091 measureMap["BTAGLEFFCORR"] = PerformanceResult::BTAGLEFFCORR;
0092 measureMap["BTAGLERRCORR"] = PerformanceResult::BTAGLERRCORR;
0093 measureMap["BTAGNBEFFCORR"] = PerformanceResult::BTAGNBEFFCORR;
0094 measureMap["BTAGNBERRCORR"] = PerformanceResult::BTAGNBERRCORR;
0095 measureMap["BTAGNBERRCORR"] = PerformanceResult::BTAGNBERRCORR;
0096 measureMap["MUEFF"] = PerformanceResult::MUEFF;
0097 measureMap["MUERR"] = PerformanceResult::MUERR;
0098 measureMap["MUFAKE"] = PerformanceResult::MUFAKE;
0099 measureMap["MUEFAKE"] = PerformanceResult::MUEFAKE;
0100
0101 if (measureName.size() != measureType.size()) {
0102 std::cout << "measureName, measureType size mismatch!" << std::endl;
0103 exit(-1);
0104 }
0105
0106 for (size_t iMeasure = 0; iMeasure < measureName.size(); iMeasure++) {
0107 std::cout << "# Dump table: " << measureName[iMeasure] << " of type " << measureType[iMeasure] << std::endl;
0108
0109
0110 const BtagPerformance& perf = iSetup.getData(measureToken[iMeasure]);
0111
0112
0113 std::cout << "# Working point: " << perf.workingPoint().cut() << std::endl;
0114 std::cout << "# [pt] [eta] [SF]" << std::endl;
0115
0116
0117
0118 BinningPointByMap measurePoint;
0119
0120 string sep = " ";
0121
0122 float bin_pt = 10.;
0123 float bin_eta = 1.2;
0124
0125 for (int i = 2; i < 50; i++) {
0126 for (int j = 0; j < 2; j++) {
0127 float pt = float(i * bin_pt);
0128 float eta = float(j * bin_eta);
0129
0130 measurePoint.insert(BinningVariables::JetEt, pt);
0131 measurePoint.insert(BinningVariables::JetAbsEta, eta);
0132
0133
0134
0135
0136 std::cout << pt << sep << pt + bin_pt << sep << eta << sep << eta + bin_eta << sep
0137 << perf.getResult(measureMap[measureType[iMeasure]], measurePoint) << std::endl;
0138
0139 measurePoint.reset();
0140 }
0141 }
0142
0143 measurePoint.reset();
0144
0145
0146
0147
0148
0149
0150
0151
0152
0153
0154
0155
0156
0157
0158
0159
0160
0161
0162
0163
0164
0165
0166
0167
0168
0169
0170 }
0171 }
0172
0173
0174 DEFINE_FWK_MODULE(DumpBtagTable);