CMS 3D CMS Logo

NewMatcher.h
Go to the documentation of this file.
1 #ifndef UtilAlgos_NewMatcher_h
2 #define UtilAlgos_NewMatcher_h
3 /* \class Matcher
4  *
5  * \author Luca Lista, INFN
6  *
7  */
18 
19 namespace reco {
20  namespace modulesNew {
21 
22  template <typename C1, typename C2, typename S, typename D = DeltaR<typename C1::value_type, typename C2::value_type> >
23  class Matcher : public edm::global::EDProducer<> {
24  public:
26  ~Matcher() override;
27 
28  private:
29  typedef typename C1::value_type T1;
30  typedef typename C2::value_type T2;
32  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
35  double distMin_;
36  double matchDistance(const T1& c1, const T2& c2) const { return distance_(c1, c2); }
37  bool select(const T1& c1, const T2& c2) const { return select_(c1, c2); }
40  };
41 
42  namespace helper {
43  typedef std::pair<size_t, double> MatchPair;
44 
45  struct SortBySecond {
46  bool operator()(const MatchPair& p1, const MatchPair& p2) const { return p1.second < p2.second; }
47  };
48  } // namespace helper
49 
50  template <typename C1, typename C2, typename S, typename D>
52  : srcToken_(consumes<C1>(cfg.template getParameter<edm::InputTag>("src"))),
53  matchedToken_(consumes<C2>(cfg.template getParameter<edm::InputTag>("matched"))),
54  distMin_(cfg.template getParameter<double>("distMin")),
55  select_(reco::modules::make<S>(cfg)),
57  produces<MatchMap>();
58  }
59 
60  template <typename C1, typename C2, typename S, typename D>
62 
63  template <typename C1, typename C2, typename S, typename D>
65  using namespace edm;
66  using namespace std;
68  evt.getByToken(matchedToken_, matched);
70  evt.getByToken(srcToken_, cands);
71  unique_ptr<MatchMap> matchMap(new MatchMap(matched));
72  size_t size = cands->size();
73  if (size != 0) {
74  typename MatchMap::Filler filler(*matchMap);
75  ::helper::MasterCollection<C1> master(cands, evt);
76  std::vector<int> indices(master.size(), -1);
77  for (size_t c = 0; c != size; ++c) {
78  const T1& cand = (*cands)[c];
79  vector<helper::MatchPair> v;
80  for (size_t m = 0; m != matched->size(); ++m) {
81  const T2& match = (*matched)[m];
82  if (select(cand, match)) {
83  double dist = matchDistance(cand, match);
84  if (dist < distMin_)
85  v.push_back(make_pair(m, dist));
86  }
87  }
88  if (!v.empty()) {
89  size_t idx = master.index(c);
90  assert(idx < indices.size());
91  indices[idx] = min_element(v.begin(), v.end(), helper::SortBySecond())->first;
92  }
93  }
94  filler.insert(master.get(), indices.begin(), indices.end());
95  filler.fill();
96  }
97  evt.put(std::move(matchMap));
98  }
99 
100  } // namespace modulesNew
101 } // namespace reco
102 
103 #endif
size
Write out results.
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
Definition: helper.py:1
edm::EDGetTokenT< C1 > srcToken_
Definition: NewMatcher.h:33
edm::Association< C2 > MatchMap
Definition: NewMatcher.h:31
S make(const edm::ParameterSet &cfg)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
edm::EDGetTokenT< C2 > matchedToken_
Definition: NewMatcher.h:34
assert(be >=bs)
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: NewMatcher.h:64
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:521
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
double matchDistance(const T1 &c1, const T2 &c2) const
Definition: NewMatcher.h:36
fixed size matrix
HLT enums.
std::pair< size_t, double > MatchPair
Definition: NewMatcher.h:43
bool operator()(const MatchPair &p1, const MatchPair &p2) const
Definition: NewMatcher.h:46
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
bool select(const T1 &c1, const T2 &c2) const
Definition: NewMatcher.h:37
def move(src, dest)
Definition: eostools.py:511
Matcher(const edm::ParameterSet &cfg)
Definition: NewMatcher.h:51