File indexing completed on 2024-04-06 12:23:39
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include <iostream>
0010 #include <vector>
0011 #include <iterator>
0012 #include <algorithm>
0013
0014 #include "FWCore/Framework/interface/Frameworkfwd.h"
0015 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0016 #include "FWCore/Framework/interface/Event.h"
0017 #include "FWCore/Framework/interface/MakerMacros.h"
0018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0019 #include "FWCore/Utilities/interface/InputTag.h"
0020
0021 #include "FWCore/Framework/interface/IOVSyncValue.h"
0022 #include "FWCore/ServiceRegistry/interface/Service.h"
0023 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0024
0025 #include "CondFormats/PhysicsToolsObjects/interface/Histogram.h"
0026 #include "CondFormats/PhysicsToolsObjects/interface/MVAComputer.h"
0027
0028 #include "PhysicsTools/MVAComputer/interface/BitSet.h"
0029 #include "PhysicsTools/MVAComputer/interface/MVAComputer.h"
0030
0031 using namespace PhysicsTools::Calibration;
0032
0033 class testWriteMVAComputerCondDB : public edm::one::EDAnalyzer<> {
0034 public:
0035 explicit testWriteMVAComputerCondDB(const edm::ParameterSet& params);
0036
0037 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0038
0039 void endJob() override;
0040
0041 private:
0042 std::string record;
0043 };
0044
0045 testWriteMVAComputerCondDB::testWriteMVAComputerCondDB(const edm::ParameterSet& params)
0046 : record(params.getUntrackedParameter<std::string>("record")) {}
0047
0048 void testWriteMVAComputerCondDB::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {}
0049
0050 void testWriteMVAComputerCondDB::endJob() {
0051
0052
0053
0054
0055
0056
0057
0058
0059
0060
0061 MVAComputerContainer container;
0062 MVAComputer* computer = &container.add("test");
0063
0064
0065
0066 Variable var;
0067 var.name = "test";
0068 computer->inputSet.push_back(var);
0069
0070 var.name = "normal";
0071 computer->inputSet.push_back(var);
0072
0073 var.name = "toast";
0074 computer->inputSet.push_back(var);
0075
0076
0077
0078 ProcNormalize norm;
0079
0080 PhysicsTools::BitSet testSet(3);
0081 testSet[0] = testSet[1] = true;
0082 norm.inputVars = convert(testSet);
0083
0084 HistogramF pdf(3, 4.0, 5.5);
0085 pdf.setBinContent(1, 1.0);
0086 pdf.setBinContent(2, 1.5);
0087 pdf.setBinContent(3, 1.0);
0088 norm.categoryIdx = -1;
0089 norm.distr.push_back(pdf);
0090 norm.distr.push_back(pdf);
0091
0092 computer->addProcessor(&norm);
0093
0094
0095
0096 ProcLikelihood lkh;
0097
0098 testSet = PhysicsTools::BitSet(5);
0099 testSet[2] = true;
0100 lkh.inputVars = convert(testSet);
0101
0102 pdf = HistogramF(6, 0.0, 1.0);
0103 pdf.setBinContent(1, 1.0);
0104 pdf.setBinContent(2, 1.5);
0105 pdf.setBinContent(3, 1.0);
0106 pdf.setBinContent(4, 1.0);
0107 pdf.setBinContent(5, 1.5);
0108 pdf.setBinContent(6, 1.0);
0109 ProcLikelihood::SigBkg sigBkg;
0110 sigBkg.signal = pdf;
0111 pdf = HistogramF(9, 0.0, 1.0);
0112 pdf.setBinContent(1, 1.0);
0113 pdf.setBinContent(2, 1.5);
0114 pdf.setBinContent(3, 1.0);
0115 pdf.setBinContent(4, 1.0);
0116 pdf.setBinContent(5, 1.5);
0117 pdf.setBinContent(6, 1.0);
0118 pdf.setBinContent(7, 1.5);
0119 pdf.setBinContent(8, 1.0);
0120 pdf.setBinContent(9, 1.7);
0121 sigBkg.background = pdf;
0122 sigBkg.useSplines = true;
0123 lkh.categoryIdx = -1;
0124 lkh.neverUndefined = true;
0125 lkh.individual = false;
0126 lkh.logOutput = false;
0127 lkh.keepEmpty = true;
0128 lkh.pdfs.push_back(sigBkg);
0129
0130 computer->addProcessor(&lkh);
0131
0132
0133
0134 testSet = PhysicsTools::BitSet(6);
0135 testSet[2] = testSet[3] = true;
0136 lkh.inputVars = convert(testSet);
0137 sigBkg.useSplines = true;
0138 lkh.pdfs.push_back(sigBkg);
0139
0140 computer->addProcessor(&lkh);
0141
0142
0143
0144 ProcOptional opt;
0145
0146 testSet = PhysicsTools::BitSet(7);
0147 testSet[5] = testSet[6] = true;
0148 opt.inputVars = convert(testSet);
0149
0150 opt.neutralPos.push_back(0.6);
0151 opt.neutralPos.push_back(0.7);
0152
0153 computer->addProcessor(&opt);
0154
0155
0156
0157 ProcMatrix pca;
0158
0159 testSet = PhysicsTools::BitSet(9);
0160 testSet[4] = testSet[7] = testSet[8] = true;
0161 pca.inputVars = convert(testSet);
0162
0163 pca.matrix.rows = 2;
0164 pca.matrix.columns = 3;
0165 double elements[] = {0.2, 0.3, 0.4, 0.8, 0.7, 0.6};
0166 std::copy(elements, elements + sizeof elements / sizeof elements[0], std::back_inserter(pca.matrix.elements));
0167
0168 computer->addProcessor(&pca);
0169
0170
0171
0172 ProcLinear lin;
0173
0174 testSet = PhysicsTools::BitSet(11);
0175 testSet[9] = testSet[10] = true;
0176 lin.inputVars = convert(testSet);
0177
0178 lin.coeffs.push_back(0.3);
0179 lin.coeffs.push_back(0.7);
0180 lin.offset = 0.0;
0181
0182 computer->addProcessor(&lin);
0183
0184
0185
0186 computer->output = 11;
0187
0188
0189
0190 PhysicsTools::MVAComputer comp(computer);
0191
0192 PhysicsTools::Variable::Value values[] = {PhysicsTools::Variable::Value("toast", 4.4),
0193 PhysicsTools::Variable::Value("toast", 4.5),
0194 PhysicsTools::Variable::Value("test", 4.6),
0195 PhysicsTools::Variable::Value("toast", 4.7),
0196 PhysicsTools::Variable::Value("test", 4.8),
0197 PhysicsTools::Variable::Value("normal", 4.9)};
0198
0199 std::cout << comp.eval(values, values + sizeof values / sizeof values[0]) << std::endl;
0200
0201
0202
0203 edm::Service<cond::service::PoolDBOutputService> dbService;
0204 if (!dbService.isAvailable())
0205 return;
0206
0207 dbService->createOneIOV(container, dbService->beginOfTime(), "BTauGenericMVAJetTagComputerRcd");
0208 }
0209
0210
0211 DEFINE_FWK_MODULE(testWriteMVAComputerCondDB);