CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
17 void
19  const edm::Event& event,
20  const edm::EventSetup& setup) {
21  selected_.clear();
22  // FIXME: the fixed reference point should be replaced with the measured beamspot
23  const math::XYZPoint beamSpot(0.,0.,0.);
24  unsigned int multiplicity;
25  double ptSum;
26  for ( collection::const_iterator iv=handle->begin(); iv!=handle->end(); ++iv ) {
27  math::XYZVector displacement(iv->position()-beamSpot);
28  if ( iv->normalizedChi2()<chi2Cut_ &&
29  fabs(displacement.z())<dzCut_ && displacement.perp2()<dr2Cut_ ) {
30  getVertexVariables(*iv,multiplicity,ptSum);
31  if ( multiplicity>=multiplicityCut_ && ptSum>ptSumCut_ ) selected_.push_back(&*iv);
32  }
33  }
34  sort(selected_.begin(),selected_.end(),*this);
35 }
36 
37 bool
39  unsigned int mult1;
40  double ptSum1;
41  getVertexVariables(*v1,mult1,ptSum1);
42  unsigned int mult2;
43  double ptSum2;
44  getVertexVariables(*v2,mult2,ptSum2);
45  return ptSum1>ptSum2;
46 }
47 
48 void
50  unsigned int& multiplicity, double& ptSum) const {
51  multiplicity = 0;
52  ptSum = 0.;
54  it!=vertex.tracks_end(); ++it) {
55  if(fabs((**it).eta())<trackEtaCut_) {
56  ++multiplicity;
57  ptSum += (**it).pt();
58  }
59  }
60 }
float dr2Cut_
cut on the normalized chi2
unsigned int multiplicityCut_
container of selected vertices
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
void getVertexVariables(const reco::Vertex &, unsigned int &, double &) const
access to track-related vertex quantities (multiplicity and pt-sum)
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
tuple handle
Definition: patZpeak.py:22
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
float ptSumCut_
minimum multiplicity of (selected) associated tracks
float chi2Cut_
eta cut used for the track selection
bool operator()(const reco::Vertex *, const reco::Vertex *) const
operator used in sorting the selected vertices
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")