CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/PhysicsTools/JetMCAlgos/plugins/CandOneToManyDeltaRMatcher.cc

Go to the documentation of this file.
00001 /* \class CandOneToManyDeltaRMatcher
00002  *
00003  * Producer for simple match map:
00004  * class to match two collections of candidate
00005  * with one-to-Many matching 
00006  * All elements of class "matched" are matched to each element
00007  * of class "source" orderd in DeltaR
00008  *
00009  */
00010 
00011 #include "FWCore/Framework/interface/EDProducer.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSetfwd.h"
00013 #include "FWCore/Utilities/interface/InputTag.h"
00014 
00015 #include<vector>
00016 #include<iostream>
00017 
00018 class CandOneToManyDeltaRMatcher : public edm::EDProducer {
00019  public:
00020   CandOneToManyDeltaRMatcher( const edm::ParameterSet & );
00021   ~CandOneToManyDeltaRMatcher();
00022  private:
00023   void produce( edm::Event&, const edm::EventSetup& );
00024   
00025   edm::InputTag source_;
00026   edm::InputTag matched_;
00027   bool printdebug_;
00028 };
00029 
00030 #include "FWCore/Framework/interface/ESHandle.h"
00031 #include "FWCore/Framework/interface/Event.h"
00032 #include "FWCore/Framework/interface/EventSetup.h"
00033 #include "FWCore/Utilities/interface/EDMException.h"
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00036 
00037 #include "DataFormats/Common/interface/Handle.h"
00038 #include "DataFormats/Candidate/interface/Candidate.h"
00039 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00040 #include "DataFormats/Candidate/interface/CandMatchMap.h"
00041 #include "DataFormats/Candidate/interface/CandMatchMapMany.h"
00042 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00043 
00044 
00045 #include <Math/VectorUtil.h>
00046 #include <TMath.h>
00047 
00048 using namespace edm;
00049 using namespace std;
00050 using namespace reco;
00051 using namespace ROOT::Math::VectorUtil;
00052 
00053 namespace reco {
00054   namespace helper {
00055     typedef pair<size_t, double> MatchPair;
00056     
00057     struct SortBySecond {
00058       bool operator()( const MatchPair & p1, const MatchPair & p2 ) const {
00059         return p1.second < p2.second;
00060       } 
00061     };
00062   }
00063 }
00064 
00065 CandOneToManyDeltaRMatcher::CandOneToManyDeltaRMatcher( const ParameterSet & cfg ) :
00066   source_( cfg.getParameter<InputTag>( "src" ) ),
00067   matched_( cfg.getParameter<InputTag>( "matched" ) ),
00068   printdebug_( cfg.getUntrackedParameter<bool>("printDebug", false) ) {
00069   produces<CandMatchMapMany>();
00070 }
00071 
00072 CandOneToManyDeltaRMatcher::~CandOneToManyDeltaRMatcher() {
00073 }
00074                 
00075 void CandOneToManyDeltaRMatcher::produce( Event& evt, const EventSetup& es ) {
00076   
00077   Handle<CandidateCollection> source;  
00078   Handle<CandidateCollection> matched;  
00079   evt.getByLabel( source_, source ) ;
00080   evt.getByLabel( matched_, matched ) ;
00081  
00082   if (printdebug_) {
00083     for( CandidateCollection::const_iterator c = source->begin(); c != source->end(); ++c ) {
00084       cout << "[CandOneToManyDeltaRMatcher] Et source  " << c->et() << endl;
00085     }    
00086     for( CandidateCollection::const_iterator c = matched->begin(); c != matched->end(); ++c ) {
00087       cout << "[CandOneToManyDeltaRMatcher] Et matched " << c->et() << endl;
00088     } 
00089   }
00090  
00091 
00092   auto_ptr<CandMatchMapMany> matchMap( new CandMatchMapMany( CandMatchMapMany::ref_type( CandidateRefProd( source  ),
00093                                                                                          CandidateRefProd( matched )
00094                                                                                          ) ) );
00095   for( size_t c = 0; c != source->size(); ++ c ) {
00096     const Candidate & src = (*source)[ c ];
00097     if (printdebug_) cout << "[CandOneToManyDeltaRMatcher] source (Et,Eta,Phi) =(" << src.et() << "," << 
00098                                                                                       src.eta() << "," << 
00099                                                                                       src.phi() << ")" << endl;
00100     vector<reco::helper::MatchPair> v;
00101     for( size_t m = 0; m != matched->size(); ++ m ) {
00102       const Candidate & match = ( * matched )[ m ];
00103       double dist = DeltaR( src.p4() , match.p4() );
00104       v.push_back( make_pair( m, dist ) );      
00105     }
00106     if ( v.size() > 0 ) {
00107       sort( v.begin(), v.end(), reco::helper::SortBySecond() );
00108       for( size_t m = 0; m != v.size(); ++ m ) {
00109         if (printdebug_) cout << "[CandOneToManyDeltaRMatcher]       match (Et,Eta,Phi) =(" << ( * matched )[ v[m].first ].et() << "," << 
00110                                                                                                ( * matched )[ v[m].first ].eta() << "," << 
00111                                                                                                ( * matched )[ v[m].first ].phi() << ") DeltaR=" << 
00112                                                                                                v[m].second  << endl;
00113         matchMap->insert( CandidateRef( source, c ), make_pair( CandidateRef( matched, v[m].first ), v[m].second  )  );
00114       }    
00115     } 
00116   }
00117   
00118   evt.put( matchMap );
00119 
00120 }
00121 
00122 #include "FWCore/PluginManager/interface/ModuleDef.h"
00123 #include "FWCore/Framework/interface/MakerMacros.h"
00124 
00125 DEFINE_FWK_MODULE( CandOneToManyDeltaRMatcher );