CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ElectronFromPVSelector.cc
Go to the documentation of this file.
1 // Includes
8 
10 
13 
16 
19 
22 
23 #include <memory>
24 #include <vector>
25 #include <sstream>
26 
27 
29 // class definition
32 {
33 public:
34  // construction/destruction
37 
38  // member functions
39  void produce(edm::Event& iEvent,const edm::EventSetup& iSetup) override;
40 
41 private:
42  // member data
43  double max_dxy_ ;
44  double max_dz_ ;
47 };
48 
49 
50 
52 // construction/destruction
54 
55 //______________________________________________________________________________
57  : max_dxy_ ( iConfig.getParameter<double>( "max_dxy" ) )
58  , max_dz_ ( iConfig.getParameter<double>( "max_dz" ) )
59  , v_recoVertexToken_ ( consumes< std::vector<reco::Vertex> >( iConfig.getParameter<edm::InputTag>( "srcVertex" ) ) )
60  , v_recoGsfElectronToken_( consumes< std::vector<reco::GsfElectron> >( iConfig.getParameter<edm::InputTag>( "srcElectron" ) ) )
61 {
62  produces<std::vector<reco::GsfElectron> >();
63 }
64 
65 
66 //______________________________________________________________________________
68 
70 // implementation of member functions
72 
73 //______________________________________________________________________________
75 {
76  std::auto_ptr<std::vector<reco::GsfElectron> > goodGsfElectrons(new std::vector<reco::GsfElectron >);
77 
79  iEvent.getByToken( v_recoVertexToken_, VertexHandle );
80 
82  iEvent.getByToken( v_recoGsfElectronToken_, GsfElectronHandle );
83 
84  if( (VertexHandle->size() == 0) || (GsfElectronHandle->size() == 0) )
85  {
86  iEvent.put(goodGsfElectrons);
87  return ;
88  }
89 
90 
91  reco::Vertex PV = VertexHandle->front();
92  std::vector<reco::GsfElectron>::const_iterator GsfElectronIt ;
93 // typename std::vector<reco::GsfElectron>::const_iterator GsfElectronIt ;
94 
95  for (GsfElectronIt = GsfElectronHandle->begin(); GsfElectronIt != GsfElectronHandle->end(); ++GsfElectronIt) {
96 
97  //int q = GsfElectronIt->gsfTrack()->charge() ;
98 
99  if ( fabs(GsfElectronIt->gsfTrack()->dxy(PV.position())) < max_dxy_ &&
100  fabs(GsfElectronIt->gsfTrack()->dz(PV.position())) < max_dz_ ) {
101  goodGsfElectrons -> push_back(*GsfElectronIt) ;
102  }
103  }
104 
105  iEvent.put(goodGsfElectrons);
106 
107 }
108 
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const Point & position() const
position
Definition: Vertex.h:106
return((rh^lh)&mask)
edm::EDGetTokenT< std::vector< reco::GsfElectron > > v_recoGsfElectronToken_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
GsfElectronFromPVSelector(const edm::ParameterSet &iConfig)
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
edm::EDGetTokenT< std::vector< reco::Vertex > > v_recoVertexToken_