CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FastCandMatcher.h
Go to the documentation of this file.
1 #ifndef CandUtils_FastCandMatcher_h
2 #define CandUtils_FastCandMatcher_h
3 /* class CandMatcher
4  *
5  * \author Luca Lista, INFN
6  *
7  */
8 #include <set>
13 
14 template<typename C>
15 class FastCandMatcher {
16 public:
20  typedef std::vector<const map_type *> map_vector;
22  explicit FastCandMatcher( const map_vector & maps );
24  explicit FastCandMatcher( const map_type & map );
26  const reco::Candidate * operator()( const reco::Candidate & ) const;
27 
28 protected:
29  const std::vector<const map_type *> & maps() const { return maps_; }
30 private:
32  std::vector<const map_type *> maps_;
33 };
34 
35 template<typename C>
37  maps_( maps ) {
38 }
39 
40 template<typename C>
42  maps_( 1, & map ) {
43 }
44 
45 template<typename C>
47  if ( c.hasMasterClone() )
48  return (*this)( * c.masterClone() );
49  unsigned int nDau = c.numberOfDaughters();
50  if ( nDau > 0 ) {
51  // check for composite candidate c
52  // navigate to daughters and find parent matches
53  std::set<const reco::Candidate *> momsIntersection, momDaughters, tmp;
54  for( reco::Candidate::const_iterator dau = c.begin(); dau != c.end(); ++ dau ) {
55  // check here generically if status == 3, then descend down to one more level
56  const reco::Candidate * dauMatch = (*this)( * dau );
57  // if a daughter does not match, return a null ref.
58  if ( dauMatch == 0 ) return 0;
59  // get matched mothers
60  size_t mothers = dauMatch->numberOfMothers();
61  for( size_t i = 0; i < mothers; ++ i ) {
62  const reco::Candidate * mom = dauMatch->mother( i );
63  if ( mom != 0 && mom->pdgId() == dauMatch->pdgId() &&
64  mom->status() == 3 && dauMatch->status() == 1 ) {
65  // assume a single mother at this point...
66  mom = mom->mother( 0 );
67  }
68  momDaughters.insert( mom );
69  }
70  // if no mother was found return null reference
71  if ( momDaughters.size() == 0 ) return 0;
72  // the first time, momsIntersection is set to momDaughters
73  if ( momsIntersection.size() == 0 ) momsIntersection = momDaughters;
74  else {
75  tmp.clear();
76  set_intersection( momsIntersection.begin(), momsIntersection.end(),
77  momDaughters.begin(), momDaughters.end(),
78  inserter( tmp, tmp.begin() ) );
79  swap( momsIntersection, tmp );
80  }
81  if ( momsIntersection.size() == 0 ) return 0;
82  }
83  // if multiple mothers are found, return a null reference
84  if ( momsIntersection.size() > 1 ) return 0;
85  // return a reference to the unique mother
86  return * momsIntersection.begin();
87  } else {
88  // check for non-composite (leaf) candidate
89  // if one of the maps contains the candidate c
90  for( typename std::vector<const map_type *>::const_iterator m = maps_.begin();
91  m != maps_.end(); ++ m ) {
92  for( typename map_type::const_iterator i = (*m)->begin(); i != (*m)->end(); ++ i ) {
93  if ( & * i->key == & c )
94  return & * i->val;
95  }
96  }
97  return 0;
98  }
99 }
100 
101 #endif
void swap(ora::Record &rh, ora::Record &lh)
Definition: Record.h:70
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
virtual int status() const =0
status word
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
FastCandMatcher(const map_vector &maps)
constructor
virtual size_type numberOfDaughters() const =0
number of daughters
virtual bool hasMasterClone() const =0
edm::AssociationMap< edm::OneToOne< C, reco::CandidateCollection > > map_type
map type
virtual const_iterator end() const =0
last daughter const_iterator
virtual int pdgId() const =0
PDG identifier.
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
virtual const_iterator begin() const =0
first daughter const_iterator
std::vector< const map_type * > map_vector
map vector
const reco::Candidate * operator()(const reco::Candidate &) const
get match from transient reference
const std::vector< const map_type * > & maps() const
std::vector< const map_type * > maps_
pointers to stored maps
virtual const CandidateBaseRef & masterClone() const =0