CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/PhysicsTools/TagAndProbe/plugins/ElectronMatchedCandidateProducer.cc

Go to the documentation of this file.
00001 #include "PhysicsTools/TagAndProbe/interface/ElectronMatchedCandidateProducer.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "FWCore/Framework/interface/MakerMacros.h"
00004 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00005 
00006 #include "DataFormats/Math/interface/deltaR.h" // reco::deltaR
00007 
00008 
00009 ElectronMatchedCandidateProducer::ElectronMatchedCandidateProducer(const edm::ParameterSet &params)
00010 {
00011 
00012   const edm::InputTag allelectrons("gsfElectrons");
00013   electronCollection_ = 
00014     params.getUntrackedParameter<edm::InputTag>("ReferenceElectronCollection", 
00015                                                 allelectrons);
00016   scCollection_ = 
00017     params.getParameter<edm::InputTag>("src");
00018 
00019   delRMatchingCut_ = params.getUntrackedParameter<double>("deltaR",
00020                                                            0.30);
00021   
00022   produces< edm::PtrVector<reco::Candidate> >();
00023   produces< edm::RefToBaseVector<reco::Candidate> >();
00024 }
00025 
00026 
00027 
00028 
00029 ElectronMatchedCandidateProducer::~ElectronMatchedCandidateProducer()
00030 {
00031 
00032 }
00033 
00034 
00035 //
00036 // member functions
00037 //
00038 
00039 
00040 // ------------ method called to produce the data  ------------
00041 
00042 void ElectronMatchedCandidateProducer::produce(edm::Event &event, 
00043                               const edm::EventSetup &eventSetup)
00044 {
00045    // Create the output collection
00046   std::auto_ptr< edm::RefToBaseVector<reco::Candidate> > 
00047     outColRef( new edm::RefToBaseVector<reco::Candidate> );
00048   std::auto_ptr< edm::PtrVector<reco::Candidate> > 
00049     outColPtr( new edm::PtrVector<reco::Candidate> );
00050 
00051 
00052   // Read electrons
00053   edm::Handle<edm::View<reco::GsfElectron> > electrons;
00054   event.getByLabel(electronCollection_, electrons);
00055    
00056 
00057 
00058   //Read candidates
00059   edm::Handle<edm::View<reco::Candidate> > recoCandColl; 
00060   event.getByLabel( scCollection_ , recoCandColl); 
00061 
00062 
00063   const edm::PtrVector<reco::Candidate>& ptrVect = recoCandColl->ptrVector();
00064   const edm::RefToBaseVector<reco::Candidate>& refs = recoCandColl->refVector();
00065   unsigned int counter=0;
00066 
00067   // Loop over candidates
00068   for(edm::View<reco::Candidate>::const_iterator scIt = recoCandColl->begin();
00069       scIt != recoCandColl->end(); ++scIt, ++counter){
00070     // Now loop over electrons
00071     for(edm::View<reco::GsfElectron>::const_iterator  elec = electrons->begin(); 
00072         elec != electrons->end();  ++elec) {
00073 
00074       reco::SuperClusterRef eSC = elec->superCluster();
00075 
00076       double dRval = reco::deltaR((float)eSC->eta(), (float)eSC->phi(), 
00077                                   scIt->eta(), scIt->phi());    
00078        
00079       if( dRval < delRMatchingCut_ ) {
00080         //outCol->push_back( *scIt );
00081         outColRef->push_back( refs[counter] );
00082         outColPtr->push_back( ptrVect[counter]  );
00083       } // end if loop
00084     } // end electron loop
00085 
00086   } // end candidate loop
00087 
00088   event.put(outColRef);
00089   event.put(outColPtr);
00090 }
00091 
00092 
00093 
00094 
00095 // ------ method called once each job just before starting event loop  ---
00096 
00097 
00098 
00099 void ElectronMatchedCandidateProducer::beginJob() {}
00100 
00101 
00102 
00103 void ElectronMatchedCandidateProducer::endJob() {}
00104 
00105 
00106 
00107 //define this as a plug-in
00108 DEFINE_FWK_MODULE( ElectronMatchedCandidateProducer );
00109