CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PhysObjectMatcher.h
Go to the documentation of this file.
1 #ifndef UtilAlgos_PhysObjectMatcher_h
2 #define UtilAlgos_PhysObjectMatcher_h
3 /* \class PhysObjectMatcher.
4  *
5  * Extended version of reco::CandMatcher.
6  * Tries to match elements from collection 1 to collection 2 with optional
7  * resolution of ambiguities. Uses three helper classes for
8  * (1) the basic selection of the match (e.g. pdgId, charge, ..);
9  * (2) a distance measure (e.g. deltaR);
10  * (3) the ranking of several matches.
11  *
12  */
22 
23 // #include <iostream>
24 
25 namespace reco {
26 
27  namespace helper {
29  template<typename D, typename C1, typename C2> class LessByMatchDistance {
30  public:
32  const C1& c1, const C2& c2) :
33  distance_(reco::modules::make<D>(cfg)), c1_(c1), c2_(c2) {}
34  bool operator() (const std::pair<size_t,size_t>& p1,
35  const std::pair<size_t,size_t>& p2) const {
36  return
37  distance_(c1_[p1.first],c2_[p1.second])<
38  distance_(c1_[p2.first],c2_[p2.second]);
39  }
40  private:
42  const C1& c1_;
43  const C2& c2_;
44  };
45  }
46 
47  // Template arguments:
48  // C1 .. candidate collection to be matched
49  // C2 .. target of the match (typically MC)
50  // S ... match (pre-)selector
51  // D ... match (typically cut on some distance)
52  // default: deltaR cut
53  // Q ... ranking of matches
54  // default: by smaller deltaR
55  template<typename C1, typename C2, typename S,
56  typename D = reco::MatchByDR<typename C1::value_type,
57  typename C2::value_type>,
58  typename Q =
60  typename C2::value_type>,
61  C1, C2 >
62  >
64  public:
67  private:
68  typedef typename C1::value_type T1;
69  typedef typename C2::value_type T2;
71  typedef std::pair<size_t, size_t> IndexPair;
72  typedef std::vector<IndexPair> MatchContainer;
73  void produce(edm::Event&, const edm::EventSetup&) override;
77  bool resolveAmbiguities_; // resolve ambiguities after
78  // first pass?
79  bool resolveByMatchQuality_; // resolve by (global) quality
80  // of match (otherwise: by order
81  // of test candidates)
82  bool select(const T1 & c1, const T2 & c2) const {
83  return select_(c1, c2);
84  }
87 // DeltaR<typename C1::value_type, typename C2::value_type> testDR_;
88  };
89 
90  template<typename C1, typename C2, typename S, typename D, typename Q>
92  config_(cfg),
93  srcToken_(consumes<C1>(cfg.template getParameter<edm::InputTag>("src"))),
94  matchedToken_(consumes<C2>(cfg.template getParameter<edm::InputTag>("matched"))),
95  resolveAmbiguities_(cfg.template getParameter<bool>("resolveAmbiguities")),
96  resolveByMatchQuality_(cfg.template getParameter<bool>("resolveByMatchQuality")),
97  select_(reco::modules::make<S>(cfg)),
98  distance_(reco::modules::make<D>(cfg)) {
99  // definition of the product
100  produces<MatchMap>();
101  // set resolveByMatchQuality only if ambiguities are to be resolved
102  resolveByMatchQuality_ = resolveByMatchQuality_ && resolveAmbiguities_;
103  }
104 
105  template<typename C1, typename C2, typename S, typename D, typename Q>
107 
108  template<typename C1, typename C2, typename S, typename D, typename Q>
110  using namespace edm;
111  using namespace std;
112  typedef std::pair<size_t, size_t> IndexPair;
113  typedef std::vector<IndexPair> MatchContainer;
114  // get collections from event
115  Handle<C2> matched;
116  evt.getByToken(matchedToken_, matched);
117  Handle<C1> cands;
118  evt.getByToken(srcToken_, cands);
119  // create product
120  auto_ptr<MatchMap> matchMap(new MatchMap(matched));
121  size_t size = cands->size();
122  if( size != 0 ) {
123  //
124  // create helpers
125  //
126  Q comparator(config_,*cands,*matched);
127  typename MatchMap::Filler filler(*matchMap);
129  vector<int> indices(master.size(), -1); // result: indices in target collection
130  vector<bool> mLock(matched->size(),false); // locks in target collection
131  MatchContainer matchPairs; // container of matched pairs
132  // loop over candidates
133  for(size_t c = 0; c != size; ++ c) {
134  const T1 & cand = (*cands)[c];
135  // no global comparison of match quality -> reset the container for each candidate
136  if ( !resolveByMatchQuality_ ) matchPairs.clear();
137  // loop over target collection
138  for(size_t m = 0; m != matched->size(); ++m) {
139  const T2 & match = (* matched)[m];
140  // check lock and preselection
141  if ( !mLock[m] && select(cand, match)) {
142 // double dist = testDR_(cand,match);
143 // cout << "dist between c = " << c << " and m = "
144 // << m << " is " << dist << " at pts of "
145 // << cand.pt() << " " << match.pt() << endl;
146  // matching requirement fulfilled -> store pair of indices
147  if ( distance_(cand,match) ) matchPairs.push_back(make_pair(c,m));
148  }
149  }
150  // if match(es) found and no global ambiguity resolution requested
151  if ( matchPairs.size()>0 && !resolveByMatchQuality_ ) {
152  // look for and store best match
153  size_t idx = master.index(c);
154  assert(idx < indices.size());
155  size_t index = min_element(matchPairs.begin(), matchPairs.end(), comparator)->second;
156  indices[idx] = index;
157  // if ambiguity resolution by order of (reco) candidates:
158  // lock element in target collection
159  if ( resolveAmbiguities_ ) mLock[index] = true;
160 // {
161 // MatchContainer::const_iterator i = min_element(matchPairs.begin(), matchPairs.end(), comparator);
162 // cout << "smallest distance for c = " << c << " is "
163 // << testDR_((*cands)[(*i).first],
164 // (*matched)[(*i).second]) << endl;
165 // }
166  }
167  }
168  // ambiguity resolution by global match quality (if requested)
169  if ( resolveByMatchQuality_ ) {
170  // sort container of all matches by quality
171  sort(matchPairs.begin(),matchPairs.end(),comparator);
172  vector<bool> cLock(master.size(),false);
173  // loop over sorted container
174  for ( MatchContainer::const_iterator i=matchPairs.begin();
175  i!=matchPairs.end(); ++i ) {
176  size_t c = (*i).first;
177  size_t m = (*i).second;
178 // cout << "rel dp = " << ((*cands)[c].pt()-(*matched)[m].pt())/(*matched)[m].pt() << endl;
179  // accept only pairs without any lock
180  if ( mLock[m] || cLock[c] ) continue;
181  // store index to target collection and lock the two items
182  size_t idx = master.index(c);
183  assert(idx < indices.size());
184  indices[idx] = m;
185  mLock[m] = true;
186  cLock[c] = true;
187  }
188  }
189  filler.insert(master.get(), indices.begin(), indices.end());
190  filler.fill();
191  }
192  evt.put(matchMap);
193  }
194 
195 }
196 
197 #endif
const edm::Handle< C1 > & get() const
std::vector< IndexPair > MatchContainer
int i
Definition: DBlmapReader.cc:9
bool select(const T1 &c1, const T2 &c2) const
tuple cfg
Definition: looper.py:259
edm::EDGetTokenT< C2 > matchedToken_
edm::ParameterSet config_
Definition: deltaR.h:79
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
assert(m_qm.get())
S make(const edm::ParameterSet &cfg)
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
edm::EDGetTokenT< C1 > srcToken_
Default class for ranking matches: sorting by smaller distance.
void produce(edm::Event &, const edm::EventSetup &) override
edm::Association< C2 > MatchMap
U second(std::pair< T, U > const &p)
PhysObjectMatcher(const edm::ParameterSet &cfg)
tuple c2
Definition: counter.py:145
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
size_t index(size_t i) const
LessByMatchDistance(const edm::ParameterSet &cfg, const C1 &c1, const C2 &c2)
Container::value_type value_type
double p2[4]
Definition: TauolaWrapper.h:90
std::pair< size_t, size_t > IndexPair
bool operator()(const std::pair< size_t, size_t > &p1, const std::pair< size_t, size_t > &p2) const
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:150
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:99
double p1[4]
Definition: TauolaWrapper.h:89
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
tuple size
Write out results.
def template
Definition: svgfig.py:520