#include <PhysicsTools/CandUtils/interface/CandMatcher.h>
Public Types | |
typedef reco::helper::CandMapTrait< C1, C2 >::type | map_type |
map type | |
typedef std::vector< const map_type * > | map_vector |
map vector | |
typedef reco::helper::CandRefTrait< C2 > ::ref_type | ref_type |
ref type | |
typedef map_type::key_type | reference_type |
concrete candidate reference type | |
typedef reco::helper::CandRefTrait< C2 > ::refProd_type | refProd_type |
refProd type | |
typedef reference_type::value_type | value_type |
concrete candidate reference type | |
Public Member Functions | |
CandMatcherBase (const map_type &map) | |
constructor | |
CandMatcherBase (const map_vector &maps) | |
constructor | |
ref_type | operator() (const reco::Candidate &) const |
get match from transient reference | |
virtual | ~CandMatcherBase () |
destructor | |
Protected Member Functions | |
virtual bool | compositePreselect (const reco::Candidate &c, const reco::Candidate &m) const =0 |
composite candidate preselection | |
virtual std::vector< const reco::Candidate * > | getDaughters (const reco::Candidate *) const =0 |
get ultimate daughter (can skip status = 3 in MC) | |
void | initMaps () |
init maps | |
const std::vector< const map_type * > & | maps () const |
Private Types | |
typedef std::map< const reco::Candidate *, reference_type > | CandRefMap |
pointer map type | |
typedef std::map< const reco::Candidate *, ref_type > | MatchedRefMap |
pointer map type | |
Private Member Functions | |
void | init () |
init at constructor | |
Private Attributes | |
CandRefMap | candRefs_ |
pointer map of candidates (e.g.: reco) | |
std::vector< const map_type * > | maps_ |
pointers to stored maps | |
refProd_type | matched_ |
reference to matched collectino | |
std::vector< std::set< size_t > > | matchedMothers_ |
mother + n.daughters indices from matched | |
MatchedRefMap | matchedRefs_ |
pointer map of matched candidates (e.g.: MC truth) |
Definition at line 16 of file CandMatcher.h.
typedef std::map<const reco::Candidate *, reference_type> CandMatcherBase< C1, C2 >::CandRefMap [private] |
typedef reco::helper::CandMapTrait<C1, C2>::type CandMatcherBase< C1, C2 >::map_type |
typedef std::vector<const map_type *> CandMatcherBase< C1, C2 >::map_vector |
typedef std::map<const reco::Candidate *, ref_type> CandMatcherBase< C1, C2 >::MatchedRefMap [private] |
typedef reco::helper::CandRefTrait<C2>::ref_type CandMatcherBase< C1, C2 >::ref_type |
typedef map_type::key_type CandMatcherBase< C1, C2 >::reference_type |
typedef reco::helper::CandRefTrait<C2>::refProd_type CandMatcherBase< C1, C2 >::refProd_type |
typedef reference_type::value_type CandMatcherBase< C1, C2 >::value_type |
CandMatcherBase< C1, C2 >::CandMatcherBase | ( | const map_vector & | maps | ) | [explicit] |
constructor
CandMatcherBase< C1, C2 >::CandMatcherBase | ( | const map_type & | map | ) | [explicit] |
constructor
CandMatcherBase< C1, C2 >::~CandMatcherBase | ( | ) | [inline, virtual] |
virtual bool CandMatcherBase< C1, C2 >::compositePreselect | ( | const reco::Candidate & | c, | |
const reco::Candidate & | m | |||
) | const [protected, pure virtual] |
composite candidate preselection
Implemented in CandMatcher< C1, C2 >.
Referenced by CandMatcherBase< C1, C2 >::operator()().
virtual std::vector<const reco::Candidate *> CandMatcherBase< C1, C2 >::getDaughters | ( | const reco::Candidate * | ) | const [protected, pure virtual] |
get ultimate daughter (can skip status = 3 in MC)
Implemented in CandMatcher< C1, C2 >.
Referenced by CandMatcherBase< C1, C2 >::initMaps().
void CandMatcherBase< C1, C2 >::init | ( | void | ) | [inline, private] |
init at constructor
Definition at line 89 of file CandMatcher.h.
References edm::errors::InvalidReference, m, CandMatcherBase< C1, C2 >::maps_, and CandMatcherBase< C1, C2 >::matched_.
00089 { 00090 matched_ = maps_.front()->refProd().val; 00091 for( typename map_vector::const_iterator m = maps_.begin() + 1; 00092 m != maps_.end(); ++ m ) { 00093 if( (*m)->refProd().val != matched_ ) 00094 throw edm::Exception( edm::errors::InvalidReference ) 00095 << "Multiple match maps specified matching different MC truth collections.\n" 00096 << "Please, specify maps all matching to the same MC truth collection.\n" 00097 << "In most of the cases you may want to match to genParticleCandidate."; 00098 } 00099 }
void CandMatcherBase< C1, C2 >::initMaps | ( | ) | [inline, protected] |
init maps
Definition at line 114 of file CandMatcher.h.
References c, CandMatcherBase< C1, C2 >::candRefs_, d, f, CandMatcherBase< C1, C2 >::getDaughters(), i, j, k, m, CandMatcherBase< C1, C2 >::maps_, CandMatcherBase< C1, C2 >::matched_, CandMatcherBase< C1, C2 >::matchedMothers_, CandMatcherBase< C1, C2 >::matchedRefs_, HcalSimpleRecAlgoImpl::reco(), and std.
Referenced by CandMatcher< C1, C2 >::CandMatcher().
00114 { 00115 using namespace reco; 00116 using namespace std; 00117 for( typename map_vector::const_iterator m = maps_.begin(); 00118 m != maps_.end(); ++ m ) { 00119 typename CandMatcherBase<C1, C2>::map_type::ref_type::key_type cands = (*m)->refProd().key; 00120 for( size_t i = 0; i < cands->size(); ++ i ) { 00121 candRefs_[ & (*cands)[ i ] ] = reference_type( cands, i ); 00122 } 00123 const C2 & matched = * matched_; 00124 size_t matchedSize = matched.size(); 00125 for( size_t i = 0; i < matchedSize; ++ i ) 00126 matchedRefs_[ & matched[ i ] ] = ref_type( matched_, i ); 00127 matchedMothers_.resize( matchedSize ); 00128 for( size_t i = 0; i < matchedSize; ++ i ) { 00129 const Candidate & c = matched[ i ]; 00130 for( Candidate::const_iterator d = c.begin(); d != c.end(); ++ d ) { 00131 vector<const Candidate *> daus = getDaughters( & * d ); 00132 for( size_t j = 0; j < daus.size(); ++ j ) { 00133 const Candidate * daughter = daus[ j ]; 00134 typename MatchedRefMap::const_iterator f = matchedRefs_.find( daughter ); 00135 if ( f == matchedRefs_.end() ) continue; 00136 size_t k = f->second.key(); 00137 assert( k < matchedMothers_.size() ); 00138 matchedMothers_[ k ].insert( i ); 00139 } 00140 } 00141 } 00142 } 00143 }
const std::vector<const map_type *>& CandMatcherBase< C1, C2 >::maps | ( | ) | const [inline, protected] |
Definition at line 48 of file CandMatcher.h.
References CandMatcherBase< C1, C2 >::maps_.
00048 { return maps_; }
CandMatcherBase< C1, C2 >::ref_type CandMatcherBase< C1, C2 >::operator() | ( | const reco::Candidate & | c | ) | const [inline] |
get match from transient reference
Definition at line 150 of file CandMatcher.h.
References reco::Candidate::begin(), CandMatcherBase< C1, C2 >::candRefs_, CandMatcherBase< C1, C2 >::compositePreselect(), d, reco::Candidate::end(), f, reco::Candidate::hasMasterClone(), edm::Ref< C, T, F >::isNull(), k, edm::Ref< C, T, F >::key(), m, CandMatcherBase< C1, C2 >::maps_, reco::Candidate::masterClone(), CandMatcherBase< C1, C2 >::matched_, CandMatcherBase< C1, C2 >::matchedMothers_, reco::Candidate::numberOfDaughters(), HcalSimpleRecAlgoImpl::reco(), std, swap(), and tmp.
00150 { 00151 using namespace reco; 00152 using namespace std; 00153 if ( c.hasMasterClone() ) 00154 return (*this)( * c.masterClone() ); 00155 unsigned int nDau = c.numberOfDaughters(); 00156 const C2 & matched = * matched_; 00157 if ( nDau > 0 ) { 00158 // check for composite candidate c 00159 // navigate to daughters and find parent matches 00160 set<size_t> momsIntersection, momDaughters, tmp; 00161 for( Candidate::const_iterator d = c.begin(); d != c.end(); ++ d ) { 00162 // check here generically if status == 3, then descend down to one more level 00163 ref_type m = (*this)( * d ); 00164 // if a daughter does not match, return a null ref. 00165 if ( m.isNull() ) return ref_type(); 00166 // get matched mother indices (fetched previously) 00167 const set<size_t> & allMomDaughters = matchedMothers_[ m.key() ]; 00168 momDaughters.clear(); 00169 for( set<size_t>::const_iterator k = allMomDaughters.begin(); 00170 k != allMomDaughters.end(); ++ k ) { 00171 size_t m = * k; 00172 if( compositePreselect( c, matched[ m ] ) ) 00173 momDaughters.insert( m ); 00174 } 00175 // if no mother was found return null reference 00176 if ( momDaughters.size() == 0 ) return ref_type(); 00177 // the first time, momsIntersection is set to momDaughters 00178 if ( momsIntersection.size() == 0 ) momsIntersection = momDaughters; 00179 else { 00180 tmp.clear(); 00181 set_intersection( momsIntersection.begin(), momsIntersection.end(), 00182 momDaughters.begin(), momDaughters.end(), 00183 inserter( tmp, tmp.begin() ) ); 00184 swap( momsIntersection, tmp ); 00185 } 00186 if ( momsIntersection.size() == 0 ) return ref_type(); 00187 } 00188 // if multiple mothers are found, return a null reference 00189 if ( momsIntersection.size() > 1 ) return ref_type(); 00190 // return a reference to the unique mother 00191 return ref_type( matched_, * momsIntersection.begin() ); 00192 } else { 00193 // check for non-composite (leaf) candidate 00194 // if one of the maps contains the candidate c 00195 for( typename std::vector<const map_type *>::const_iterator m = maps_.begin(); 00196 m != maps_.end(); ++ m ) { 00197 typename CandRefMap::const_iterator f = candRefs_.find( & c ); 00198 if ( f != candRefs_.end() ) { 00199 reference_type ref = f->second; 00200 typename map_type::const_iterator f = (*m)->find( ref ); 00201 if ( f != (*m)->end() ) { 00202 return f->val; 00203 } 00204 } 00205 } 00206 return ref_type(); 00207 } 00208 }
CandRefMap CandMatcherBase< C1, C2 >::candRefs_ [private] |
pointer map of candidates (e.g.: reco)
Definition at line 59 of file CandMatcher.h.
Referenced by CandMatcherBase< C1, C2 >::initMaps(), and CandMatcherBase< C1, C2 >::operator()().
std::vector<const map_type *> CandMatcherBase< C1, C2 >::maps_ [private] |
pointers to stored maps
Definition at line 51 of file CandMatcher.h.
Referenced by CandMatcherBase< C1, C2 >::init(), CandMatcherBase< C1, C2 >::initMaps(), CandMatcherBase< C1, C2 >::maps(), and CandMatcherBase< C1, C2 >::operator()().
refProd_type CandMatcherBase< C1, C2 >::matched_ [private] |
reference to matched collectino
Definition at line 53 of file CandMatcher.h.
Referenced by CandMatcherBase< C1, C2 >::init(), CandMatcherBase< C1, C2 >::initMaps(), and CandMatcherBase< C1, C2 >::operator()().
std::vector<std::set<size_t> > CandMatcherBase< C1, C2 >::matchedMothers_ [private] |
mother + n.daughters indices from matched
Definition at line 63 of file CandMatcher.h.
Referenced by CandMatcherBase< C1, C2 >::initMaps(), and CandMatcherBase< C1, C2 >::operator()().
MatchedRefMap CandMatcherBase< C1, C2 >::matchedRefs_ [private] |
pointer map of matched candidates (e.g.: MC truth)
Definition at line 61 of file CandMatcher.h.
Referenced by CandMatcherBase< C1, C2 >::initMaps().