CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
IPCutPFCandidateSelectorDefinition.h
Go to the documentation of this file.
1 #ifndef CommonTools_ParticleFlow_IPCutPFCandidateSelectorDefinition
2 #define CommonTools_ParticleFlow_IPCutPFCandidateSelectorDefinition
3 
19 
20 namespace pf2pat {
21 
23 
25  verticesToken_( iC.consumes<reco::VertexCollection>( cfg.getParameter<edm::InputTag> ( "vertices" ) ) ),
26  d0Cut_( cfg.getParameter<double>("d0Cut") ),
27  dzCut_( cfg.getParameter<double>("dzCut") ),
28  d0SigCut_( cfg.getParameter<double>("d0SigCut") ),
29  dzSigCut_( cfg.getParameter<double>("dzSigCut") ) {}
30 
31  void select( const HandleToCollection & hc,
32  const edm::Event & e,
33  const edm::EventSetup& s) {
34  selected_.clear();
35 
37  e.getByToken(verticesToken_, vertices);
38  if (vertices->empty()) return;
39  const reco::Vertex &vtx = (*vertices)[0];
40 
41  unsigned key=0;
42  for( collection::const_iterator pfc = hc->begin();
43  pfc != hc->end(); ++pfc, ++key) {
44 
45  bool passing = true;
46  const reco::Track *tk = 0;
47  if (pfc->gsfTrackRef().isNonnull()) tk = pfc->gsfTrackRef().get();
48  else if (pfc->trackRef().isNonnull()) tk = pfc->trackRef().get();
49 
50  if (tk != 0) {
51  double d0 = fabs(tk->dxy(vtx.position()));
52  double dz = fabs(tk->dz(vtx.position()));
53  double d0e = hypot(tk->dxyError(), hypot(vtx.xError(), vtx.yError()));
54  double dze = hypot(tk->dzError(), vtx.zError());
55  if (d0Cut_ > 0 && d0 > d0Cut_) passing = false;
56  if (dzCut_ > 0 && dz > dzCut_) passing = false;
57  if (d0SigCut_ > 0 && d0e > 0 && d0/d0e > d0SigCut_) passing = false;
58  if (dzSigCut_ > 0 && dze > 0 && dz/dze > dzSigCut_) passing = false;
59  }
60 
61  if( passing ) {
62  selected_.push_back( reco::PFCandidate(*pfc) );
63  reco::PFCandidatePtr ptrToMother( hc, key );
64 
65  if ( pfc->numberOfSourceCandidatePtrs() > 0 ) {
66  selected_.back().setSourceCandidatePtr( edm::Ptr<reco::PFCandidate>( pfc->sourceCandidatePtr(0)
67  ) );
68  }
69  else {
70  selected_.back().setSourceCandidatePtr( ptrToMother );
71  }
72  }
73  }
74  }
75 
76  private:
78  double d0Cut_;
79  double dzCut_;
80  double d0SigCut_;
81  double dzSigCut_;
82  };
83 }
84 
85 #endif
Selects PFCandidates basing on their compatibility with vertex.
void select(const HandleToCollection &hc, const edm::Event &e, const edm::EventSetup &s)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
double zError() const
error on z
Definition: Vertex.h:104
double dxyError() const
error on dxy
Definition: TrackBase.h:207
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const Point & position() const
position
Definition: Vertex.h:92
IPCutPFCandidateSelectorDefinition(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
edm::EDGetTokenT< reco::VertexCollection > verticesToken_
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:125
double dzError() const
error on dz
Definition: TrackBase.h:213
double xError() const
error on x
Definition: Vertex.h:100
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:38
list key
Definition: combine.py:13
susybsm::HSCParticleCollection hc
Definition: classes.h:25
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:119
double yError() const
error on y
Definition: Vertex.h:102