CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Attributes
pf2pat::IPCutPFCandidateSelectorDefinition Class Reference

Selects PFCandidates basing on their compatibility with vertex. More...

#include "CommonTools/ParticleFlow/interface/IPCutPFCandidateSelectorDefinition.h"

Inheritance diagram for pf2pat::IPCutPFCandidateSelectorDefinition:
pf2pat::PFCandidateSelectorDefinition

Public Member Functions

 IPCutPFCandidateSelectorDefinition (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
 
void select (const HandleToCollection &hc, const edm::Event &e, const edm::EventSetup &s)
 
- Public Member Functions inherited from pf2pat::PFCandidateSelectorDefinition
const_iterator begin () const
 
const_iterator end () const
 
 PFCandidateSelectorDefinition ()
 
const containerselected () const
 
size_t size () const
 

Private Attributes

double d0Cut_
 
double d0SigCut_
 
double dtCut_
 
double dtSigCut_
 
double dzCut_
 
double dzSigCut_
 
edm::EDGetTokenT
< reco::VertexCollection
verticesToken_
 

Additional Inherited Members

- Public Types inherited from pf2pat::PFCandidateSelectorDefinition
typedef reco::PFCandidateCollection collection
 
typedef
boost::transform_iterator
< Pointer,
container::const_iterator > 
const_iterator
 
typedef std::vector
< reco::PFCandidate
container
 
typedef edm::Handle< collectionHandleToCollection
 
- Protected Attributes inherited from pf2pat::PFCandidateSelectorDefinition
container selected_
 

Detailed Description

Selects PFCandidates basing on their compatibility with vertex.

Author
Giovanni Petrucciani
Version
Id:
IPCutPFCandidateSelectorDefinition.h,v 1.2 2011/04/06 12:12:38 rwolf Exp

Definition at line 24 of file IPCutPFCandidateSelectorDefinition.h.

Constructor & Destructor Documentation

pf2pat::IPCutPFCandidateSelectorDefinition::IPCutPFCandidateSelectorDefinition ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  iC 
)
inline

Definition at line 25 of file IPCutPFCandidateSelectorDefinition.h.

27  d0Cut_(cfg.getParameter<double>("d0Cut")),
28  dzCut_(cfg.getParameter<double>("dzCut")),
29  dtCut_(cfg.getParameter<double>("dtCut")),
30  d0SigCut_(cfg.getParameter<double>("d0SigCut")),
31  dzSigCut_(cfg.getParameter<double>("dzSigCut")),
32  dtSigCut_(cfg.getParameter<double>("dtSigCut")) {}
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::VertexCollection > verticesToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303

Member Function Documentation

void pf2pat::IPCutPFCandidateSelectorDefinition::select ( const HandleToCollection hc,
const edm::Event e,
const edm::EventSetup s 
)
inline

Definition at line 34 of file IPCutPFCandidateSelectorDefinition.h.

References d0, d0Cut_, d0SigCut_, dt, dtCut_, dtSigCut_, reco::TrackBase::dxy(), reco::TrackBase::dxyError(), PVValHelper::dz, reco::TrackBase::dz(), dzCut_, reco::TrackBase::dzError(), dzSigCut_, edm::Event::getByToken(), submitPVResolutionJobs::key, compareTotals::passing, reco::Vertex::position(), pf2pat::PFCandidateSelectorDefinition::selected_, mathSSE::sqrt(), reco::Vertex::t(), reco::Vertex::tError(), beam_dqm_sourceclient-live_cfg::vertices, verticesToken_, reco::Vertex::xError(), reco::Vertex::yError(), and reco::Vertex::zError().

34  {
35  selected_.clear();
36 
38  e.getByToken(verticesToken_, vertices);
39  if (vertices->empty())
40  return;
41  const reco::Vertex &vtx = (*vertices)[0];
42  double vt = vtx.t();
43  double vte = vtx.tError();
44 
45  unsigned key = 0;
46  for (collection::const_iterator pfc = hc->begin(); pfc != hc->end(); ++pfc, ++key) {
47  bool passing = true;
48  const reco::Track *tk = nullptr;
49  if (pfc->gsfTrackRef().isNonnull())
50  tk = pfc->gsfTrackRef().get();
51  else if (pfc->trackRef().isNonnull())
52  tk = pfc->trackRef().get();
53 
54  if (tk != nullptr) {
55  double d0 = fabs(tk->dxy(vtx.position()));
56  double dz = fabs(tk->dz(vtx.position()));
57  double d0e = hypot(tk->dxyError(), hypot(vtx.xError(), vtx.yError()));
58  double dze = hypot(tk->dzError(), vtx.zError());
59  if (d0Cut_ > 0 && d0 > d0Cut_)
60  passing = false;
61  if (dzCut_ > 0 && dz > dzCut_)
62  passing = false;
63  if (d0SigCut_ > 0 && d0e > 0 && d0 / d0e > d0SigCut_)
64  passing = false;
65  if (dzSigCut_ > 0 && dze > 0 && dz / dze > dzSigCut_)
66  passing = false;
67  }
68  double pfct = pfc->time();
69  double pfcte = pfc->timeError();
70  double dt = fabs(pfct - vt);
71  double dte = std::sqrt(pfcte * pfcte + vte * vte);
72  if (dtCut_ > 0 && pfcte > 0 && vte > 0 && dt > dtCut_)
73  passing = false;
74  if (dtSigCut_ > 0 && pfcte > 0 && vte > 0 && dt / dte > dtSigCut_)
75  passing = false;
76 
77  if (passing) {
78  selected_.push_back(reco::PFCandidate(*pfc));
79  reco::PFCandidatePtr ptrToMother(hc, key);
80 
81  if (pfc->numberOfSourceCandidatePtrs() > 0) {
82  selected_.back().setSourceCandidatePtr(edm::Ptr<reco::PFCandidate>(pfc->sourceCandidatePtr(0)));
83  } else {
84  selected_.back().setSourceCandidatePtr(ptrToMother);
85  }
86  }
87  }
88  }
float dt
Definition: AMPTWrapper.h:136
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
double zError() const
error on z
Definition: Vertex.h:141
double dxyError() const
error on dxy
Definition: TrackBase.h:769
const Point & position() const
position
Definition: Vertex.h:127
edm::EDGetTokenT< reco::VertexCollection > verticesToken_
T sqrt(T t)
Definition: SSEVec.h:19
tuple key
prepare the HTCondor submission files and eventually submit them
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:622
double dzError() const
error on dz
Definition: TrackBase.h:778
static constexpr float d0
double xError() const
error on x
Definition: Vertex.h:137
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
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:608
double yError() const
error on y
Definition: Vertex.h:139
double tError() const
error on t
Definition: Vertex.h:143
double t() const
t coordinate
Definition: Vertex.h:135

Member Data Documentation

double pf2pat::IPCutPFCandidateSelectorDefinition::d0Cut_
private

Definition at line 92 of file IPCutPFCandidateSelectorDefinition.h.

Referenced by select().

double pf2pat::IPCutPFCandidateSelectorDefinition::d0SigCut_
private

Definition at line 95 of file IPCutPFCandidateSelectorDefinition.h.

Referenced by select().

double pf2pat::IPCutPFCandidateSelectorDefinition::dtCut_
private

Definition at line 94 of file IPCutPFCandidateSelectorDefinition.h.

Referenced by select().

double pf2pat::IPCutPFCandidateSelectorDefinition::dtSigCut_
private

Definition at line 97 of file IPCutPFCandidateSelectorDefinition.h.

Referenced by select().

double pf2pat::IPCutPFCandidateSelectorDefinition::dzCut_
private

Definition at line 93 of file IPCutPFCandidateSelectorDefinition.h.

Referenced by select().

double pf2pat::IPCutPFCandidateSelectorDefinition::dzSigCut_
private

Definition at line 96 of file IPCutPFCandidateSelectorDefinition.h.

Referenced by select().

edm::EDGetTokenT<reco::VertexCollection> pf2pat::IPCutPFCandidateSelectorDefinition::verticesToken_
private

Definition at line 91 of file IPCutPFCandidateSelectorDefinition.h.

Referenced by select().