CMS 3D CMS Logo

HGC3DClusterTMVASelector.cc
Go to the documentation of this file.
6 
11 
12 #include "TMVA/Factory.h"
13 #include "TMVA/Reader.h"
14 
15 namespace l1t {
17  public:
20 
21  private:
22  class Var {
23  public:
25  void declare(TMVA::Reader &r) { r.AddVariable(name_, &val_); }
26  void fill(const l1t::HGCalMulticluster &c) { val_ = expr_(c); }
27 
28  private:
31  float val_;
32  };
33 
36  std::vector<Var> variables_;
38  std::unique_ptr<TMVA::Reader> reader_;
40 
41  void produce(edm::Event &, const edm::EventSetup &) override;
42 
43  }; // class
44 } // namespace l1t
45 
47  : src_(consumes<l1t::HGCalMulticlusterBxCollection>(iConfig.getParameter<edm::InputTag>("src"))),
48  preselection_(iConfig.getParameter<std::string>("preselection")),
49  method_(iConfig.getParameter<std::string>("method")),
50  weightsFile_(iConfig.getParameter<std::string>("weightsFile")),
51  reader_(new TMVA::Reader()),
52  wp_(iConfig.getParameter<std::string>("wp")) {
53  // first create all the variables
54  for (const auto &psvar : iConfig.getParameter<std::vector<edm::ParameterSet>>("variables")) {
55  variables_.emplace_back(psvar.getParameter<std::string>("name"), psvar.getParameter<std::string>("value"));
56  }
57  // then declare them
58  for (auto &var : variables_)
59  var.declare(*reader_);
60  // then read the weights
61  if (weightsFile_[0] != '/' && weightsFile_[0] != '.') {
63  }
65  // finally, declare outputs
66  produces<l1t::HGCalMulticlusterBxCollection>();
67  produces<l1t::HGCalMulticlusterBxCollection>("fail");
68 }
69 
71  std::unique_ptr<l1t::HGCalMulticlusterBxCollection> out = std::make_unique<l1t::HGCalMulticlusterBxCollection>();
72  std::unique_ptr<l1t::HGCalMulticlusterBxCollection> fail = std::make_unique<l1t::HGCalMulticlusterBxCollection>();
73 
75  iEvent.getByToken(src_, multiclusters);
76 
77  for (int bx = multiclusters->getFirstBX(); bx <= multiclusters->getLastBX(); ++bx) {
78  for (auto it = multiclusters->begin(bx), ed = multiclusters->end(bx); it != ed; ++it) {
79  const auto &c = *it;
80  if (preselection_(c)) {
81  for (auto &var : variables_)
82  var.fill(c);
83  float mvaOut = reader_->EvaluateMVA(method_);
84  if (mvaOut > wp_(c)) {
85  out->push_back(bx, c);
86  } else {
87  fail->push_back(bx, c);
88  }
89  }
90  }
91  }
92 
93  iEvent.put(std::move(out));
94  iEvent.put(std::move(fail), "fail");
95 }
Var(const std::string &name, const std::string &expr)
void produce(edm::Event &, const edm::EventSetup &) override
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::unique_ptr< TMVA::Reader > reader_
delete x;
Definition: CaloConfig.h:22
StringCutObjectSelector< l1t::HGCalMulticluster > preselection_
int iEvent
Definition: GenABIO.cc:224
StringObjectFunction< l1t::HGCalMulticluster > expr_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
StringObjectFunction< l1t::HGCalMulticluster > wp_
edm::EDGetTokenT< l1t::HGCalMulticlusterBxCollection > src_
TMVA::IMethod * loadTMVAWeights(TMVA::Reader *reader, const std::string &method, const std::string &weightFile, bool verbose=false)
HLT enums.
const std::string & fullPath() const
Definition: FileInPath.cc:144
HGC3DClusterTMVASelector(const edm::ParameterSet &)
void fill(const l1t::HGCalMulticluster &c)
def move(src, dest)
Definition: eostools.py:511