CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/CommonTools/UtilAlgos/interface/NewMatcher.h

Go to the documentation of this file.
00001 #ifndef UtilAlgos_NewMatcher_h
00002 #define UtilAlgos_NewMatcher_h
00003 /* \class Matcher
00004  *
00005  * \author Luca Lista, INFN
00006  *
00007  */
00008 #include "FWCore/Framework/interface/EDProducer.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 #include "FWCore/Utilities/interface/InputTag.h"
00011 #include "CommonTools/UtilAlgos/interface/DeltaR.h"
00012 #include "CommonTools/UtilAlgos/interface/MasterCollectionHelper.h"
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "DataFormats/Common/interface/Handle.h"
00015 #include "DataFormats/Common/interface/Association.h"
00016 #include "DataFormats/Common/interface/getRef.h"
00017 #include "DataFormats/Common/interface/View.h"
00018 
00019 namespace reco {
00020   namespace modulesNew {
00021 
00022     template<typename C1, typename C2,
00023              typename S, typename D = DeltaR<typename C1::value_type, typename C2::value_type> >
00024     class Matcher : public edm::EDProducer {
00025     public:
00026       Matcher(const edm::ParameterSet & cfg);
00027       ~Matcher();
00028     private:
00029       typedef typename C1::value_type T1;
00030       typedef typename C2::value_type T2;
00031       typedef edm::Association<C2> MatchMap;
00032       void produce(edm::Event&, const edm::EventSetup&);
00033       edm::InputTag src_;
00034       edm::InputTag matched_;
00035       double distMin_;
00036       double matchDistance(const T1 & c1, const T2 & c2) const {
00037         return distance_(c1, c2);
00038       }
00039       bool select(const T1 & c1, const T2 & c2) const { 
00040         return select_(c1, c2); 
00041       }
00042       S select_;
00043       D distance_;
00044     };
00045 
00046     namespace helper {
00047       typedef std::pair<size_t, double> MatchPair;
00048 
00049       struct SortBySecond {
00050         bool operator()(const MatchPair & p1, const MatchPair & p2) const {
00051           return p1.second < p2.second;
00052         } 
00053       };
00054     }
00055 
00056     template<typename C1, typename C2, typename S, typename D>
00057     Matcher<C1, C2, S, D>::Matcher(const edm::ParameterSet & cfg) :
00058       src_(cfg.template getParameter<edm::InputTag>("src")),
00059       matched_(cfg.template getParameter<edm::InputTag>("matched")), 
00060       distMin_(cfg.template getParameter<double>("distMin")),
00061       select_(reco::modules::make<S>(cfg)), 
00062       distance_(reco::modules::make<D>(cfg)) {
00063       produces<MatchMap>();
00064     }
00065 
00066     template<typename C1, typename C2, typename S, typename D>
00067     Matcher<C1, C2, S, D>::~Matcher() { }
00068 
00069     template<typename C1, typename C2, typename S, typename D>    
00070     void Matcher<C1, C2, S, D>::produce(edm::Event& evt, const edm::EventSetup&) {
00071       using namespace edm;
00072       using namespace std;
00073       Handle<C2> matched;  
00074       evt.getByLabel(matched_, matched);
00075       Handle<C1> cands;  
00076       evt.getByLabel(src_, cands);
00077       auto_ptr<MatchMap> matchMap(new MatchMap(matched));
00078       size_t size = cands->size();
00079       if( size != 0 ) {
00080         typename MatchMap::Filler filler(*matchMap);
00081 	::helper::MasterCollection<C1> master(cands);
00082         std::vector<int> indices(master.size(), -1);
00083         for(size_t c = 0; c != size; ++ c) {
00084           const T1 & cand = (*cands)[c];
00085           vector<helper::MatchPair> v;
00086           for(size_t m = 0; m != matched->size(); ++m) {
00087             const T2 & match = (* matched)[m];
00088             if (select(cand, match)) {
00089               double dist = matchDistance(cand, match);
00090               if (dist < distMin_) v.push_back(make_pair(m, dist));
00091             }
00092           }
00093           if(v.size() > 0) {
00094             size_t idx = master.index(c);
00095             assert(idx < indices.size());
00096             indices[idx] = min_element(v.begin(), v.end(), helper::SortBySecond())->first;
00097           }
00098         }
00099         filler.insert(master.get(), indices.begin(), indices.end());
00100         filler.fill();
00101       }
00102       evt.put(matchMap);
00103     }    
00104      
00105   }
00106 }
00107 
00108 #endif