1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
|
#ifndef CandUtils_CandMatcherNew_h
#define CandUtils_CandMatcherNew_h
/* class CandMatcher
*
* \author Luca Lista, INFN
*
*/
#include "DataFormats/Candidate/interface/Candidate.h"
#include "DataFormats/Common/interface/Association.h"
#include "FWCore/Utilities/interface/EDMException.h"
#include <algorithm>
#include <iterator>
#include <set>
namespace reco {
namespace utilsNew {
template <typename C>
class CandMatcher {
public:
/// map type
typedef edm::Association<C> map_type;
/// ref type
typedef typename edm::Association<C>::reference_type reference_type;
typedef std::vector<const map_type *> map_vector;
/// constructor
explicit CandMatcher(const map_vector &maps);
/// constructor
explicit CandMatcher(const map_type &map);
/// destructor
virtual ~CandMatcher();
/// get match from transient reference
reference_type operator[](const reco::Candidate &) const;
/// reference to matched collection
typename map_type::refprod_type ref() const { return map_.ref(); }
protected:
/// match map at leaf level
map_type map_;
private:
};
template <typename C>
CandMatcher<C>::CandMatcher(const typename CandMatcher<C>::map_vector &maps) : map_() {
for (typename map_vector::const_iterator i = maps.begin(); i != maps.end(); ++i)
map_ += **i;
}
template <typename C>
CandMatcher<C>::CandMatcher(const typename CandMatcher<C>::map_type &map) : map_(map) {}
template <typename C>
CandMatcher<C>::~CandMatcher() {}
template <typename C>
typename CandMatcher<C>::reference_type CandMatcher<C>::operator[](const reco::Candidate &c) const {
using namespace reco;
using namespace std;
if (c.hasMasterClone()) {
const CandidateBaseRef &master = c.masterClone();
return master->numberOfDaughters() == 0 ? map_[master] : (*this)[*master];
}
size_t nDau = c.numberOfDaughters();
if (nDau == 0)
return reference_type();
set<size_t> momIdx, common, tmp;
for (size_t i = 0; i < nDau; ++i) {
const Candidate &d = *c.daughter(i);
reference_type m = (*this)[d];
if (m.isNull())
return reference_type();
momIdx.clear();
while (m->numberOfMothers() == 1) {
m = m->motherRef();
momIdx.insert(m.key());
}
if (momIdx.empty())
return reference_type();
if (common.empty())
common = momIdx;
else {
tmp.clear();
set_intersection(common.begin(), common.end(), momIdx.begin(), momIdx.end(), inserter(tmp, tmp.begin()));
swap(common, tmp);
}
if (common.empty())
return reference_type();
}
size_t idx = *max_element(common.begin(), common.end());
return reference_type(map_.ref(), idx);
}
} // namespace utilsNew
} // namespace reco
#endif
|