CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
TrackFilterForPVFinding Class Reference

#include <TrackFilterForPVFinding.h>

Inheritance diagram for TrackFilterForPVFinding:
TrackFilterForPVFindingBase HITrackFilterForPVFinding

Public Member Functions

bool operator() (const reco::TransientTrack &tracks) const
 
std::vector< reco::TransientTrackselect (const std::vector< reco::TransientTrack > &tracks) const override
 
 TrackFilterForPVFinding (const edm::ParameterSet &conf)
 
- Public Member Functions inherited from TrackFilterForPVFindingBase
 TrackFilterForPVFindingBase ()
 
 TrackFilterForPVFindingBase (const edm::ParameterSet &conf)
 
virtual ~TrackFilterForPVFindingBase ()
 

Private Attributes

float maxD0Sig_
 
float maxEta_
 
float maxNormChi2_
 
float minPt_
 
int minPxLayers_
 
int minSiLayers_
 
reco::TrackBase::TrackQuality quality_
 

Detailed Description

Description: track selection for PV finding

Definition at line 14 of file TrackFilterForPVFinding.h.

Constructor & Destructor Documentation

TrackFilterForPVFinding::TrackFilterForPVFinding ( const edm::ParameterSet conf)

Definition at line 4 of file TrackFilterForPVFinding.cc.

References edm::ParameterSet::getParameter(), maxD0Sig_, maxEta_, maxNormChi2_, minPt_, minPxLayers_, minSiLayers_, quality_, reco::TrackBase::qualityByName(), trackPseudoSelection_cff::qualityClass, AlCaHLTBitMon_QueryRunRegistry::string, and reco::TrackBase::undefQuality.

5 {
6  maxD0Sig_ = conf.getParameter<double>("maxD0Significance");
7  minPt_ = conf.getParameter<double>("minPt");
8  maxEta_ = conf.getParameter<double>("maxEta");
9  maxNormChi2_ = conf.getParameter<double>("maxNormalizedChi2");
10  minSiLayers_ = conf.getParameter<int>("minSiliconLayersWithHits");
11  minPxLayers_ = conf.getParameter<int>("minPixelLayersWithHits");
12 
13  // the next few lines are taken from RecoBTag/SecondaryVertex/interface/TrackSelector.h"
15  conf.getParameter<std::string>("trackQuality");
16  if (qualityClass == "any" || qualityClass == "Any" ||
17  qualityClass == "ANY" || qualityClass == "") {
19  } else {
21  }
22 
23 }
T getParameter(std::string const &) const
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:134
reco::TrackBase::TrackQuality quality_

Member Function Documentation

bool TrackFilterForPVFinding::operator() ( const reco::TransientTrack tracks) const

Definition at line 27 of file TrackFilterForPVFinding.cc.

References reco::TrackBase::dzError(), Measurement1D::error(), PV3DBase< T, PVType, FrameType >::eta(), fftjetproducer_cfi::etaCut, TrajectoryStateOnSurface::globalMomentum(), reco::TransientTrack::hitPattern(), reco::TransientTrack::impactPointState(), TrajectoryStateClosestToBeamLine::isValid(), maxD0Sig_, maxEta_, maxNormChi2_, minPt_, minPxLayers_, minSiLayers_, reco::TransientTrack::normalizedChi2(), reco::HitPattern::pixelLayersWithMeasurement(), reco::TrackBase::quality(), quality_, Measurement1D::significance(), reco::TransientTrack::stateAtBeamLine(), reco::TransientTrack::track(), reco::HitPattern::trackerLayersWithMeasurement(), PV3DBase< T, PVType, FrameType >::transverse(), TrajectoryStateClosestToBeamLine::transverseImpactParameter(), and reco::TrackBase::undefQuality.

28 {
29  if (!tk.stateAtBeamLine().isValid()) return false;
30  bool IPSigCut = ( tk.stateAtBeamLine().transverseImpactParameter().significance() < maxD0Sig_ )
31  && (tk.stateAtBeamLine().transverseImpactParameter().error() < 1.0)
32  && (tk.track().dzError() < 1.0);
33  bool pTCut = tk.impactPointState().globalMomentum().transverse() > minPt_;
34  bool etaCut = std::fabs( tk.impactPointState().globalMomentum().eta() ) < maxEta_;
35  bool normChi2Cut = tk.normalizedChi2() < maxNormChi2_;
36  bool nPxLayCut = tk.hitPattern().pixelLayersWithMeasurement() >= minPxLayers_;
37  bool nSiLayCut = tk.hitPattern().trackerLayersWithMeasurement() >= minSiLayers_;
38  bool trackQualityCut = (quality_==reco::TrackBase::undefQuality)|| tk.track().quality(quality_);
39 
40  return IPSigCut && pTCut && etaCut && normChi2Cut && nPxLayCut && nSiLayCut && trackQualityCut;
41 }
reco::TrackBase::TrackQuality quality_
std::vector< reco::TransientTrack > TrackFilterForPVFinding::select ( const std::vector< reco::TransientTrack > &  tracks) const
overridevirtual

Implements TrackFilterForPVFindingBase.

Definition at line 46 of file TrackFilterForPVFinding.cc.

Referenced by HITrackFilterForPVFinding::select().

47 {
48  std::vector <reco::TransientTrack> seltks;
49  for (std::vector<reco::TransientTrack>::const_iterator itk = tracks.begin();
50  itk != tracks.end(); itk++) {
51  if ( operator()(*itk) ) seltks.push_back(*itk); // calls the filter function for single tracks
52  }
53  return seltks;
54 }

Member Data Documentation

float TrackFilterForPVFinding::maxD0Sig_
private

Definition at line 24 of file TrackFilterForPVFinding.h.

Referenced by operator()(), and TrackFilterForPVFinding().

float TrackFilterForPVFinding::maxEta_
private

Definition at line 24 of file TrackFilterForPVFinding.h.

Referenced by operator()(), and TrackFilterForPVFinding().

float TrackFilterForPVFinding::maxNormChi2_
private

Definition at line 26 of file TrackFilterForPVFinding.h.

Referenced by operator()(), and TrackFilterForPVFinding().

float TrackFilterForPVFinding::minPt_
private

Definition at line 24 of file TrackFilterForPVFinding.h.

Referenced by operator()(), and TrackFilterForPVFinding().

int TrackFilterForPVFinding::minPxLayers_
private

Definition at line 25 of file TrackFilterForPVFinding.h.

Referenced by operator()(), and TrackFilterForPVFinding().

int TrackFilterForPVFinding::minSiLayers_
private

Definition at line 25 of file TrackFilterForPVFinding.h.

Referenced by operator()(), and TrackFilterForPVFinding().

reco::TrackBase::TrackQuality TrackFilterForPVFinding::quality_
private

Definition at line 27 of file TrackFilterForPVFinding.h.

Referenced by operator()(), and TrackFilterForPVFinding().