CMS 3D CMS Logo

ClassifierMerger.cc
Go to the documentation of this file.
3 
8 
9 
10 #include<vector>
11 #include<memory>
12 
13 namespace {
14  class ClassifierMerger final : public edm::global::EDProducer<> {
15  public:
16  explicit ClassifierMerger(const edm::ParameterSet& conf) {
17  for (auto const & it : conf.getParameter<std::vector<std::string> >("inputClassifiers")) {
18  srcMVAs.push_back(consumes<MVACollection>(edm::InputTag(it,"MVAValues")));
19  srcQuals.push_back(consumes<QualityMaskCollection>(edm::InputTag(it,"QualityMasks")));
20  }
21 
22  produces<MVACollection>("MVAValues");
23  produces<QualityMaskCollection>("QualityMasks");
24 
25  }
26 
27  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
29  desc.add<std::vector<std::string> >("inputClassifiers",std::vector<std::string>());
30  descriptions.add("ClassifierMerger", desc);
31  }
32 
33  private:
34 
35  using MVACollection = std::vector<float>;
36  using QualityMaskCollection = std::vector<unsigned char>;
37 
38  void produce(edm::StreamID, edm::Event& evt, const edm::EventSetup&) const override {
39 
40  // get Master
42  evt.getByToken(srcMVAs[0], hmva);
43  auto size = (*hmva).size();
44 
46  evt.getByToken(srcQuals[0], hqual);
47 
48 
49  // products
50  auto mvas = std::make_unique<MVACollection>(*hmva);
51  auto quals = std::make_unique<QualityMaskCollection>(*hqual);
52 
53  for (auto i=1U; i<srcQuals.size(); ++i) {
54  evt.getByToken(srcQuals[i], hqual);
55  auto const & iq = *hqual;
56  assert(iq.size()==size);
57  for (auto j=0U; j!=size; ++j) (*quals)[j] |= iq[j];
58  }
59 
60  evt.put(std::move(mvas),"MVAValues");
61  evt.put(std::move(quals),"QualityMasks");
62 
63  }
64 
65  std::vector<edm::EDGetTokenT<MVACollection>> srcMVAs;
66  std::vector<edm::EDGetTokenT<QualityMaskCollection>> srcQuals;
67 
68 
69  };
70 }
71 
72 
75 
76 DEFINE_FWK_MODULE(ClassifierMerger);
77 
size
Write out results.
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void produce(StreamID, Event &, EventSetup const &) const =0
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(ConfigurationDescriptions &descriptions)
def move(src, dest)
Definition: eostools.py:510