CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/Validation/RecoTau/plugins/ElectronFromPVSelector.cc

Go to the documentation of this file.
00001 
00002 // Includes
00004 #include "FWCore/Framework/interface/EDProducer.h"
00005 #include "FWCore/Framework/interface/Event.h"
00006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00007 #include "FWCore/Utilities/interface/InputTag.h"
00008 
00009 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
00010 
00011 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00012 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00013 
00014 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00015 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00016 
00017 #include "DataFormats/GsfTrackReco/interface/GsfTrackExtra.h"
00018 #include "DataFormats/GsfTrackReco/interface/GsfTrackExtraFwd.h"
00019 
00020 #include "DataFormats/VertexReco/interface/Vertex.h"
00021 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00022 
00023 #include <memory>
00024 #include <vector>
00025 #include <sstream>
00026 
00027 
00029 // class definition
00031 class GsfElectronFromPVSelector : public edm::EDProducer
00032 {
00033 public:
00034   // construction/destruction
00035   GsfElectronFromPVSelector(const edm::ParameterSet& iConfig);
00036   virtual ~GsfElectronFromPVSelector();
00037   
00038   // member functions
00039   void produce(edm::Event& iEvent,const edm::EventSetup& iSetup);
00040   void endJob();
00041 
00042 private:  
00043   // member data
00044   edm::InputTag     srcPart_ ;  
00045   edm::InputTag     srcPV_   ;  
00046   double            max_dxy_ ;
00047   double            max_dz_  ;
00048 };
00049 
00050 
00051 
00053 // construction/destruction
00055 
00056 //______________________________________________________________________________
00057 GsfElectronFromPVSelector::GsfElectronFromPVSelector(const edm::ParameterSet& iConfig)
00058   : srcPart_(iConfig.getParameter<edm::InputTag>("srcElectron"))
00059   , srcPV_  (iConfig.getParameter<edm::InputTag>("srcVertex"))
00060   , max_dxy_(iConfig.getParameter<double>("max_dxy"))
00061   , max_dz_ (iConfig.getParameter<double>("max_dz"))
00062 {
00063   produces<std::vector<reco::GsfElectron> >();
00064 }
00065 
00066 
00067 //______________________________________________________________________________
00068 GsfElectronFromPVSelector::~GsfElectronFromPVSelector(){}
00069 
00071 // implementation of member functions
00073 
00074 //______________________________________________________________________________
00075 void GsfElectronFromPVSelector::produce(edm::Event& iEvent,const edm::EventSetup& iSetup)
00076 {  
00077   std::auto_ptr<std::vector<reco::GsfElectron> > goodGsfElectrons(new std::vector<reco::GsfElectron >);
00078   
00079   edm::Handle< std::vector<reco::Vertex> > VertexHandle;
00080   iEvent.getByLabel(srcPV_,VertexHandle);
00081 
00082   edm::Handle< std::vector<reco::GsfElectron> > GsfElectronHandle;
00083   iEvent.getByLabel(srcPart_,GsfElectronHandle);
00084   
00085   if( (VertexHandle->size() == 0) || (GsfElectronHandle->size() == 0) ) 
00086   {
00087     iEvent.put(goodGsfElectrons);
00088     return ;
00089   }
00090   
00091   
00092   reco::Vertex PV = VertexHandle->front();   
00093   std::vector<reco::GsfElectron>::const_iterator GsfElectronIt ;
00094 //  typename std::vector<reco::GsfElectron>::const_iterator GsfElectronIt ;
00095 
00096   for (GsfElectronIt = GsfElectronHandle->begin(); GsfElectronIt != GsfElectronHandle->end(); ++GsfElectronIt) {
00097     
00098     //int q = GsfElectronIt->gsfTrack()->charge() ;
00099     
00100     if ( fabs(GsfElectronIt->gsfTrack()->dxy(PV.position())) < max_dxy_ && 
00101          fabs(GsfElectronIt->gsfTrack()->dz(PV.position()))  < max_dz_  ) {
00102         goodGsfElectrons -> push_back(*GsfElectronIt) ;
00103     }
00104   }  
00105   
00106   iEvent.put(goodGsfElectrons);
00107   
00108 }
00109 
00110 void GsfElectronFromPVSelector::endJob()
00111 {
00112 }
00113 
00114 #include "FWCore/Framework/interface/MakerMacros.h"
00115 DEFINE_FWK_MODULE(GsfElectronFromPVSelector);