CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackFromPVSelector.cc
Go to the documentation of this file.
1 // Includes
8 
10 
13 
16 
17 #include <memory>
18 #include <vector>
19 #include <sstream>
20 
21 
23 // class definition
26 {
27 public:
28  // construction/destruction
29  TrackFromPVSelector(const edm::ParameterSet& iConfig);
30  virtual ~TrackFromPVSelector();
31 
32  // member functions
33  void produce(edm::Event& iEvent,const edm::EventSetup& iSetup) override;
34 
35 private:
36  // member data
37  double max_dxy_ ;
38  double max_dz_ ;
41 };
42 
43 
44 
46 // construction/destruction
48 
49 //______________________________________________________________________________
51  : max_dxy_ ( iConfig.getParameter<double>( "max_dxy" ) )
52  , max_dz_ ( iConfig.getParameter<double>( "max_dz" ) )
53  , v_recoVertexToken_( consumes< std::vector<reco::Vertex> >( iConfig.getParameter<edm::InputTag>( "srcVertex" ) ) )
54  , v_recoTrackToken_ ( consumes< std::vector<reco::Track> >( iConfig.getParameter<edm::InputTag>( "srcTrack" ) ) )
55 {
56  produces<std::vector<reco::Track> >();
57 }
58 
59 
60 //______________________________________________________________________________
62 
64 // implementation of member functions
66 
67 //______________________________________________________________________________
69 {
70  std::auto_ptr<std::vector<reco::Track> > goodTracks(new std::vector<reco::Track >);
71 
73  iEvent.getByToken( v_recoVertexToken_, VertexHandle );
74 
76  iEvent.getByToken( v_recoTrackToken_, TrackHandle );
77 
78  if( (VertexHandle->size() == 0) || (TrackHandle->size() == 0) )
79  {
80  iEvent.put(goodTracks);
81  return ;
82  }
83 
84  reco::Vertex PV = VertexHandle->front();
85  //typename std::vector<reco::Track>::const_iterator TrackIt ;
86  std::vector<reco::Track>::const_iterator TrackIt ;
87 
88  for (TrackIt = TrackHandle->begin(); TrackIt != TrackHandle->end(); ++TrackIt) {
89  if ( fabs(TrackIt->dxy(PV.position())) < max_dxy_ &&
90  fabs(TrackIt->dz(PV.position())) < max_dz_ ){
91  goodTracks -> push_back(*TrackIt) ;
92  }
93  }
94 
95  iEvent.put(goodTracks);
96 
97 }
98 
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< std::vector< reco::Vertex > > v_recoVertexToken_
const Point & position() const
position
Definition: Vertex.h:106
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
return((rh^lh)&mask)
edm::EDGetTokenT< std::vector< reco::Track > > v_recoTrackToken_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
TrackFromPVSelector(const edm::ParameterSet &iConfig)