CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/CommonTools/UtilAlgos/interface/PhysObjectMatcher.h

Go to the documentation of this file.
00001 #ifndef UtilAlgos_PhysObjectMatcher_h
00002 #define UtilAlgos_PhysObjectMatcher_h
00003 /* \class PhysObjectMatcher.
00004  *
00005  * Extended version of reco::CandMatcher.
00006  * Tries to match elements from collection 1 to collection 2 with optional 
00007  * resolution of ambiguities. Uses three helper classes for 
00008  *  (1) the basic selection of the match (e.g. pdgId, charge, ..);
00009  *  (2) a distance measure (e.g. deltaR);
00010  *  (3) the ranking of several matches.
00011  *
00012  */
00013 #include "FWCore/Framework/interface/EDProducer.h"
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015 #include "FWCore/Utilities/interface/InputTag.h"
00016 #include "CommonTools/UtilAlgos/interface/DeltaR.h"
00017 #include "CommonTools/UtilAlgos/interface/MasterCollectionHelper.h"
00018 #include "FWCore/Framework/interface/Event.h"
00019 #include "DataFormats/Common/interface/Handle.h"
00020 #include "DataFormats/Common/interface/Association.h"
00021 #include "CommonTools/UtilAlgos/interface/MatchByDR.h"
00022 
00023 // #include <iostream>
00024 
00025 namespace reco {
00026 
00027   namespace helper {
00029     template<typename D, typename C1, typename C2> class LessByMatchDistance {
00030     public:
00031       LessByMatchDistance (const edm::ParameterSet& cfg,
00032                            const C1& c1, const C2& c2) :
00033         distance_(reco::modules::make<D>(cfg)), c1_(c1), c2_(c2) {}
00034       bool operator() (const std::pair<size_t,size_t>& p1,
00035                        const std::pair<size_t,size_t>& p2) const {
00036         return 
00037           distance_(c1_[p1.first],c2_[p1.second])<
00038           distance_(c1_[p2.first],c2_[p2.second]);
00039       }
00040     private:
00041       D distance_;
00042       const C1& c1_;
00043       const C2& c2_;
00044     };
00045   }
00046   
00047   // Template arguments:
00048   // C1 .. candidate collection to be matched
00049   // C2 .. target of the match (typically MC)
00050   // S ... match (pre-)selector
00051   // D ... match (typically cut on some distance)
00052   //         default: deltaR cut
00053   // Q ... ranking of matches 
00054   //         default: by smaller deltaR
00055   template<typename C1, typename C2, typename S, 
00056            typename D = reco::MatchByDR<typename C1::value_type, 
00057                                         typename C2::value_type>,
00058            typename Q = 
00059            helper::LessByMatchDistance< DeltaR<typename C1::value_type, 
00060                                                typename C2::value_type>,
00061                                         C1, C2 >
00062   >
00063   class PhysObjectMatcher : public edm::EDProducer {
00064   public:
00065     PhysObjectMatcher(const edm::ParameterSet & cfg);
00066     ~PhysObjectMatcher();
00067   private:
00068     typedef typename C1::value_type T1;
00069     typedef typename C2::value_type T2;
00070     typedef edm::Association<C2> MatchMap;
00071     typedef std::pair<size_t, size_t> IndexPair;
00072     typedef std::vector<IndexPair> MatchContainer;
00073     void produce(edm::Event&, const edm::EventSetup&);
00074     edm::ParameterSet config_;
00075     edm::InputTag src_;
00076     edm::InputTag matched_;
00077     bool resolveAmbiguities_;            // resolve ambiguities after
00078                                          //   first pass?
00079     bool resolveByMatchQuality_;         // resolve by (global) quality
00080                                          //   of match (otherwise: by order
00081                                          //   of test candidates)
00082     bool select(const T1 & c1, const T2 & c2) const { 
00083       return select_(c1, c2); 
00084     }
00085     S select_;
00086     D distance_;
00087 //     DeltaR<typename C1::value_type, typename C2::value_type> testDR_;
00088   };
00089 
00090   template<typename C1, typename C2, typename S, typename D, typename Q>
00091   PhysObjectMatcher<C1, C2, S, D, Q>::PhysObjectMatcher(const edm::ParameterSet & cfg) :
00092     config_(cfg),
00093     src_(cfg.template getParameter<edm::InputTag>("src")),
00094     matched_(cfg.template getParameter<edm::InputTag>("matched")), 
00095     resolveAmbiguities_(cfg.template getParameter<bool>("resolveAmbiguities")),
00096     resolveByMatchQuality_(cfg.template getParameter<bool>("resolveByMatchQuality")),
00097     select_(reco::modules::make<S>(cfg)), 
00098     distance_(reco::modules::make<D>(cfg)) {
00099     // definition of the product
00100     produces<MatchMap>();
00101     // set resolveByMatchQuality only if ambiguities are to be resolved
00102     resolveByMatchQuality_ = resolveByMatchQuality_ && resolveAmbiguities_;
00103   }
00104 
00105   template<typename C1, typename C2, typename S, typename D, typename Q>
00106   PhysObjectMatcher<C1, C2, S, D, Q>::~PhysObjectMatcher() { }
00107 
00108   template<typename C1, typename C2, typename S, typename D, typename Q>    
00109   void PhysObjectMatcher<C1, C2, S, D, Q>::produce(edm::Event& evt, const edm::EventSetup&) {
00110     using namespace edm;
00111     using namespace std;
00112     typedef std::pair<size_t, size_t> IndexPair;
00113     typedef std::vector<IndexPair> MatchContainer;
00114     // get collections from event
00115     Handle<C2> matched;  
00116     evt.getByLabel(matched_, matched);
00117     Handle<C1> cands;  
00118     evt.getByLabel(src_, cands);
00119     // create product
00120     auto_ptr<MatchMap> matchMap(new MatchMap(matched));
00121     size_t size = cands->size();
00122     if( size != 0 ) {
00123       //
00124       // create helpers
00125       //
00126       Q comparator(config_,*cands,*matched);
00127       typename MatchMap::Filler filler(*matchMap);
00128       ::helper::MasterCollection<C1> master(cands);
00129       vector<int> indices(master.size(), -1);      // result: indices in target collection
00130       vector<bool> mLock(matched->size(),false);   // locks in target collection
00131       MatchContainer matchPairs;                   // container of matched pairs
00132       // loop over candidates
00133       for(size_t c = 0; c != size; ++ c) {
00134         const T1 & cand = (*cands)[c];
00135         // no global comparison of match quality -> reset the container for each candidate
00136         if ( !resolveByMatchQuality_ )  matchPairs.clear();
00137         // loop over target collection
00138         for(size_t m = 0; m != matched->size(); ++m) {
00139           const T2 & match = (* matched)[m];
00140           // check lock and preselection
00141           if ( !mLock[m] && select(cand, match)) {
00142 //          double dist = testDR_(cand,match);
00143 //          cout << "dist between c = " << c << " and m = "
00144 //               << m << " is " << dist << " at pts of "
00145 //               << cand.pt() << " " << match.pt() << endl;
00146             // matching requirement fulfilled -> store pair of indices
00147             if ( distance_(cand,match) )  matchPairs.push_back(make_pair(c,m));
00148           }
00149         }
00150         // if match(es) found and no global ambiguity resolution requested
00151         if ( matchPairs.size()>0 && !resolveByMatchQuality_ ) {
00152           // look for and store best match
00153           size_t idx = master.index(c);
00154           assert(idx < indices.size());
00155           size_t index = min_element(matchPairs.begin(), matchPairs.end(), comparator)->second;
00156           indices[idx] = index;
00157           // if ambiguity resolution by order of (reco) candidates:
00158           //   lock element in target collection
00159           if ( resolveAmbiguities_ )  mLock[index] = true;
00160 //        {
00161 //          MatchContainer::const_iterator i = min_element(matchPairs.begin(), matchPairs.end(), comparator);
00162 //          cout << "smallest distance for c = " << c << " is " 
00163 //               << testDR_((*cands)[(*i).first],
00164 //                          (*matched)[(*i).second]) << endl;
00165 //        }
00166         }
00167       }
00168       // ambiguity resolution by global match quality (if requested)
00169       if ( resolveByMatchQuality_ ) {
00170         // sort container of all matches by quality
00171         sort(matchPairs.begin(),matchPairs.end(),comparator);
00172         vector<bool> cLock(master.size(),false);
00173         // loop over sorted container
00174         for ( MatchContainer::const_iterator i=matchPairs.begin();
00175               i!=matchPairs.end(); ++i ) {
00176           size_t c = (*i).first;
00177           size_t m = (*i).second;
00178 //        cout << "rel dp = " << ((*cands)[c].pt()-(*matched)[m].pt())/(*matched)[m].pt() << endl;
00179           // accept only pairs without any lock
00180           if ( mLock[m] || cLock[c] )  continue;
00181           // store index to target collection and lock the two items
00182           size_t idx = master.index(c);
00183           assert(idx < indices.size());
00184           indices[idx] = m;
00185           mLock[m] = true;
00186           cLock[c] = true;
00187         }
00188       }
00189       filler.insert(master.get(), indices.begin(), indices.end());
00190       filler.fill();
00191     }
00192     evt.put(matchMap);
00193   }    
00194      
00195 }
00196 
00197 #endif