CMS 3D CMS Logo

PATPrimaryVertexSelector.cc
Go to the documentation of this file.
5 
7  : multiplicityCut_(cfg.getParameter<unsigned int>("minMultiplicity")),
8  ptSumCut_(cfg.getParameter<double>("minPtSum")),
9  trackEtaCut_(cfg.getParameter<double>("maxTrackEta")),
10  chi2Cut_(cfg.getParameter<double>("maxNormChi2")),
11  dr2Cut_(cfg.getParameter<double>("maxDeltaR")),
12  dzCut_(cfg.getParameter<double>("maxDeltaZ")) {
13  // store squared cut (avoid using sqrt() for each vertex)
14  dr2Cut_ *= dr2Cut_;
15 }
16 
18  const edm::Event& event,
19  const edm::EventSetup& setup) {
20  selected_.clear();
21  // FIXME: the fixed reference point should be replaced with the measured beamspot
22  const math::XYZPoint beamSpot(0., 0., 0.);
23  unsigned int multiplicity;
24  double ptSum;
25  for (collection::const_iterator iv = handle->begin(); iv != handle->end(); ++iv) {
26  math::XYZVector displacement(iv->position() - beamSpot);
27  if (iv->normalizedChi2() < chi2Cut_ && fabs(displacement.z()) < dzCut_ && displacement.perp2() < dr2Cut_) {
28  getVertexVariables(*iv, multiplicity, ptSum);
29  if (multiplicity >= multiplicityCut_ && ptSum > ptSumCut_)
30  selected_.push_back(&*iv);
31  }
32  }
33  sort(selected_.begin(), selected_.end(), *this);
34 }
35 
37  unsigned int mult1;
38  double ptSum1;
39  getVertexVariables(*v1, mult1, ptSum1);
40  unsigned int mult2;
41  double ptSum2;
42  getVertexVariables(*v2, mult2, ptSum2);
43  return ptSum1 > ptSum2;
44 }
45 
47  unsigned int& multiplicity,
48  double& ptSum) const {
49  multiplicity = 0;
50  ptSum = 0.;
51  for (reco::Vertex::trackRef_iterator it = vertex.tracks_begin(); it != vertex.tracks_end(); ++it) {
52  if (fabs((**it).eta()) < trackEtaCut_) {
53  ++multiplicity;
54  ptSum += (**it).pt();
55  }
56  }
57 }
float dr2Cut_
cut on the normalized chi2
unsigned int multiplicityCut_
container of selected vertices
PATPrimaryVertexSelector(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
void select(const edm::Handle< collection > &, const edm::Event &, const edm::EventSetup &)
needed for use with an ObjectSelector
float dzCut_
cut on the (squared) transverse position
float trackEtaCut_
minimum pt sum o (selected) associated tracks
void getVertexVariables(const reco::Vertex &, unsigned int &, double &) const
access to track-related vertex quantities (multiplicity and pt-sum)
float ptSumCut_
minimum multiplicity of (selected) associated tracks
float chi2Cut_
eta cut used for the track selection
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
bool operator()(const reco::Vertex *, const reco::Vertex *) const
operator used in sorting the selected vertices
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:38
Definition: event.py:1