CMS 3D CMS Logo

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>
31  public:
32  LessByMatchDistance(const edm::ParameterSet& cfg, 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, const std::pair<size_t, size_t>& p2) const {
35  return distance_(c1_[p1.first], c2_[p1.second]) < distance_(c1_[p2.first], c2_[p2.second]);
36  }
37 
38  private:
40  const C1& c1_;
41  const C2& c2_;
42  };
43  } // namespace helper
44 
45  // Template arguments:
46  // C1 .. candidate collection to be matched
47  // C2 .. target of the match (typically MC)
48  // S ... match (pre-)selector
49  // D ... match (typically cut on some distance)
50  // default: deltaR cut
51  // Q ... ranking of matches
52  // default: by smaller deltaR
53  template <typename C1,
54  typename C2,
55  typename S,
59  public:
61  ~PhysObjectMatcher() override;
62 
63  private:
64  typedef typename C1::value_type T1;
65  typedef typename C2::value_type T2;
67  typedef std::pair<size_t, size_t> IndexPair;
68  typedef std::vector<IndexPair> MatchContainer;
69  void produce(edm::Event&, const edm::EventSetup&) override;
73  bool resolveAmbiguities_; // resolve ambiguities after
74  // first pass?
75  bool resolveByMatchQuality_; // resolve by (global) quality
76  // of match (otherwise: by order
77  // of test candidates)
78  bool select(const T1& c1, const T2& c2) const { return select_(c1, c2); }
81  // DeltaR<typename C1::value_type, typename C2::value_type> testDR_;
82  };
83 
84  template <typename C1, typename C2, typename S, typename D, typename Q>
86  : config_(cfg),
87  srcToken_(consumes<C1>(cfg.template getParameter<edm::InputTag>("src"))),
88  matchedToken_(consumes<C2>(cfg.template getParameter<edm::InputTag>("matched"))),
89  resolveAmbiguities_(cfg.template getParameter<bool>("resolveAmbiguities")),
90  resolveByMatchQuality_(cfg.template getParameter<bool>("resolveByMatchQuality")),
91  select_(reco::modules::make<S>(cfg)),
93  // definition of the product
94  produces<MatchMap>();
95  // set resolveByMatchQuality only if ambiguities are to be resolved
97  }
98 
99  template <typename C1, typename C2, typename S, typename D, typename Q>
101 
102  template <typename C1, typename C2, typename S, typename D, typename Q>
104  using namespace edm;
105  using namespace std;
106  typedef std::pair<size_t, size_t> IndexPair;
107  typedef std::vector<IndexPair> MatchContainer;
108  // get collections from event
110  evt.getByToken(matchedToken_, matched);
112  evt.getByToken(srcToken_, cands);
113  // create product
114  unique_ptr<MatchMap> matchMap(new MatchMap(matched));
115  size_t size = cands->size();
116  if (size != 0) {
117  //
118  // create helpers
119  //
120  Q comparator(config_, *cands, *matched);
121  typename MatchMap::Filler filler(*matchMap);
122  ::helper::MasterCollection<C1> master(cands, evt);
123  vector<int> indices(master.size(), -1); // result: indices in target collection
124  vector<bool> mLock(matched->size(), false); // locks in target collection
125  MatchContainer matchPairs; // container of matched pairs
126  // loop over candidates
127  for (size_t c = 0; c != size; ++c) {
128  const T1& cand = (*cands)[c];
129  // no global comparison of match quality -> reset the container for each candidate
130  if (!resolveByMatchQuality_)
131  matchPairs.clear();
132  // loop over target collection
133  for (size_t m = 0; m != matched->size(); ++m) {
134  const T2& match = (*matched)[m];
135  // check lock and preselection
136  if (!mLock[m] && select(cand, match)) {
137  // double dist = testDR_(cand,match);
138  // cout << "dist between c = " << c << " and m = "
139  // << m << " is " << dist << " at pts of "
140  // << cand.pt() << " " << match.pt() << endl;
141  // matching requirement fulfilled -> store pair of indices
142  if (distance_(cand, match))
143  matchPairs.push_back(make_pair(c, m));
144  }
145  }
146  // if match(es) found and no global ambiguity resolution requested
147  if (!matchPairs.empty() && !resolveByMatchQuality_) {
148  // look for and store best match
149  size_t idx = master.index(c);
150  assert(idx < indices.size());
151  size_t index = min_element(matchPairs.begin(), matchPairs.end(), comparator)->second;
152  indices[idx] = index;
153  // if ambiguity resolution by order of (reco) candidates:
154  // lock element in target collection
155  if (resolveAmbiguities_)
156  mLock[index] = true;
157  // {
158  // MatchContainer::const_iterator i = min_element(matchPairs.begin(), matchPairs.end(), comparator);
159  // cout << "smallest distance for c = " << c << " is "
160  // << testDR_((*cands)[(*i).first],
161  // (*matched)[(*i).second]) << endl;
162  // }
163  }
164  }
165  // ambiguity resolution by global match quality (if requested)
166  if (resolveByMatchQuality_) {
167  // sort container of all matches by quality
168  sort(matchPairs.begin(), matchPairs.end(), comparator);
169  vector<bool> cLock(master.size(), false);
170  // loop over sorted container
171  for (MatchContainer::const_iterator i = matchPairs.begin(); i != matchPairs.end(); ++i) {
172  size_t c = (*i).first;
173  size_t m = (*i).second;
174  // cout << "rel dp = " << ((*cands)[c].pt()-(*matched)[m].pt())/(*matched)[m].pt() << endl;
175  // accept only pairs without any lock
176  if (mLock[m] || cLock[c])
177  continue;
178  // store index to target collection and lock the two items
179  size_t idx = master.index(c);
180  assert(idx < indices.size());
181  indices[idx] = m;
182  mLock[m] = true;
183  cLock[c] = true;
184  }
185  }
186  filler.insert(master.get(), indices.begin(), indices.end());
187  filler.fill();
188  }
189  evt.put(std::move(matchMap));
190  }
191 
192 } // namespace reco
193 
194 #endif
reco::PhysObjectMatcher::MatchMap
edm::Association< C2 > MatchMap
Definition: PhysObjectMatcher.h:66
muonTagProbeFilters_cff.matched
matched
Definition: muonTagProbeFilters_cff.py:62
bTagCombinedSVVariables_cff.indices
indices
Definition: bTagCombinedSVVariables_cff.py:67
reco::PhysObjectMatcher::select_
S select_
Definition: PhysObjectMatcher.h:79
reco::MatchByDR
Definition: MatchByDR.h:12
Handle.h
electrons_cff.bool
bool
Definition: electrons_cff.py:366
mps_fire.i
i
Definition: mps_fire.py:428
edm::Association::Filler
Definition: Association.h:78
reco::PhysObjectMatcher
Definition: PhysObjectMatcher.h:58
modules
Definition: MuonCleanerBySegments.cc:35
edm::EDGetTokenT< C1 >
reco::PhysObjectMatcher::~PhysObjectMatcher
~PhysObjectMatcher() override
Definition: PhysObjectMatcher.h:100
edm
HLT enums.
Definition: AlignableModifier.h:19
reco::PhysObjectMatcher::T2
C2::value_type T2
Definition: PhysObjectMatcher.h:65
singleTopDQM_cfi.select
select
Definition: singleTopDQM_cfi.py:50
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89281
reco::PhysObjectMatcher::T1
C1::value_type T1
Definition: PhysObjectMatcher.h:64
reco::PhysObjectMatcher::matchedToken_
edm::EDGetTokenT< C2 > matchedToken_
Definition: PhysObjectMatcher.h:72
cms::cuda::assert
assert(be >=bs)
reco::PhysObjectMatcher::distance_
D distance_
Definition: PhysObjectMatcher.h:80
EDProducer.h
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:222
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
Association.h
reco::PhysObjectMatcher::select
bool select(const T1 &c1, const T2 &c2) const
Definition: PhysObjectMatcher.h:78
edm::Handle
Definition: AssociativeIterator.h:50
class-composition.Q
Q
Definition: class-composition.py:82
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
p2
double p2[4]
Definition: TauolaWrapper.h:90
reco::helper::LessByMatchDistance::c1_
const C1 & c1_
Definition: PhysObjectMatcher.h:40
MatchByDR.h
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:535
reco::helper::LessByMatchDistance::distance_
D distance_
Definition: PhysObjectMatcher.h:39
reco::PhysObjectMatcher::MatchContainer
std::vector< IndexPair > MatchContainer
Definition: PhysObjectMatcher.h:68
S
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:97
DeltaR.h
distance_
double distance_
Definition: PFRecoTauChargedHadronFromGenericTrackPlugin.cc:198
reco::PhysObjectMatcher::srcToken_
edm::EDGetTokenT< C1 > srcToken_
Definition: PhysObjectMatcher.h:71
HLT_FULL_cff.cands
cands
Definition: HLT_FULL_cff.py:15144
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
MasterCollectionHelper.h
jetUpdater_cfi.sort
sort
Definition: jetUpdater_cfi.py:29
trigObjTnPSource_cfi.filler
filler
Definition: trigObjTnPSource_cfi.py:21
reco::PhysObjectMatcher::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: PhysObjectMatcher.h:103
match
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
helper
Definition: helper.py:1
cand
Definition: decayParser.h:32
alignmentValidation.c1
c1
do drawing
Definition: alignmentValidation.py:1025
edm::Event::put
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
edm::stream::EDProducer
Definition: EDProducer.h:38
p1
double p1[4]
Definition: TauolaWrapper.h:89
reco::PhysObjectMatcher::resolveAmbiguities_
bool resolveAmbiguities_
Definition: PhysObjectMatcher.h:73
edm::Association
Definition: Association.h:18
edm::EventSetup
Definition: EventSetup.h:58
reco::JetExtendedAssociation::value_type
Container::value_type value_type
Definition: JetExtendedAssociation.h:30
svgfig.template
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:521
reco::modules::make
S make(const edm::ParameterSet &cfg)
Definition: ParameterAdapter.h:21
InputTag.h
looper.cfg
cfg
Definition: looper.py:297
funct::D
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
reco::helper::LessByMatchDistance
Default class for ranking matches: sorting by smaller distance.
Definition: PhysObjectMatcher.h:30
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
reco::PhysObjectMatcher::IndexPair
std::pair< size_t, size_t > IndexPair
Definition: PhysObjectMatcher.h:67
reco::PhysObjectMatcher::config_
edm::ParameterSet config_
Definition: PhysObjectMatcher.h:70
reco::helper::LessByMatchDistance::LessByMatchDistance
LessByMatchDistance(const edm::ParameterSet &cfg, const C1 &c1, const C2 &c2)
Definition: PhysObjectMatcher.h:32
S
Definition: CSCDBL1TPParametersExtended.h:16
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
reco::helper::LessByMatchDistance::operator()
bool operator()(const std::pair< size_t, size_t > &p1, const std::pair< size_t, size_t > &p2) const
Definition: PhysObjectMatcher.h:34
reco::PhysObjectMatcher::PhysObjectMatcher
PhysObjectMatcher(const edm::ParameterSet &cfg)
Definition: PhysObjectMatcher.h:85
ParameterSet.h
c
auto & c
Definition: CAHitNtupletGeneratorKernelsImpl.h:46
volumeBasedMagneticField_160812_cfi.master
master
Definition: volumeBasedMagneticField_160812_cfi.py:60
edm::Event
Definition: Event.h:73
reco::PhysObjectMatcher::resolveByMatchQuality_
bool resolveByMatchQuality_
Definition: PhysObjectMatcher.h:75
reco::helper::LessByMatchDistance::c2_
const C2 & c2_
Definition: PhysObjectMatcher.h:41
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443