File indexing completed on 2024-09-12 04:16:49
0001 #include "RecoVertex/PrimaryVertexProducer/interface/TrackFilterForPVFinding.h"
0002 #include <cmath>
0003
0004 TrackFilterForPVFinding::TrackFilterForPVFinding(const edm::ParameterSet& conf) {
0005 maxD0Sig_ = conf.getParameter<double>("maxD0Significance");
0006 maxD0Error_ = conf.getParameter<double>("maxD0Error");
0007 maxDzError_ = conf.getParameter<double>("maxDzError");
0008 minPt_ = conf.getParameter<double>("minPt");
0009 maxEta_ = conf.getParameter<double>("maxEta");
0010 maxNormChi2_ = conf.getParameter<double>("maxNormalizedChi2");
0011 minSiLayers_ = conf.getParameter<int>("minSiliconLayersWithHits");
0012 minPxLayers_ = conf.getParameter<int>("minPixelLayersWithHits");
0013 minStripHits_ = conf.getParameter<int>("minValidStripHits");
0014
0015
0016 std::string qualityClass = conf.getParameter<std::string>("trackQuality");
0017 if (qualityClass == "any" || qualityClass == "Any" || qualityClass == "ANY" || qualityClass.empty()) {
0018 quality_ = reco::TrackBase::undefQuality;
0019 } else {
0020 quality_ = reco::TrackBase::qualityByName(qualityClass);
0021 }
0022 }
0023
0024
0025 bool TrackFilterForPVFinding::operator()(const reco::TransientTrack& tk) const {
0026 if (!tk.stateAtBeamLine().isValid())
0027 return false;
0028 bool IPSigCut = (tk.stateAtBeamLine().transverseImpactParameter().significance() < maxD0Sig_) &&
0029 (tk.stateAtBeamLine().transverseImpactParameter().error() < maxD0Error_) &&
0030 (tk.track().dzError() < maxDzError_);
0031 bool pTCut = tk.impactPointState().globalMomentum().transverse() > minPt_;
0032 bool etaCut = std::fabs(tk.impactPointState().globalMomentum().eta()) < maxEta_;
0033 bool normChi2Cut = tk.normalizedChi2() < maxNormChi2_;
0034 bool nPxLayCut = tk.hitPattern().pixelLayersWithMeasurement() >= minPxLayers_;
0035 bool nSiLayCut = tk.hitPattern().trackerLayersWithMeasurement() >= minSiLayers_;
0036 bool trackQualityCut = (quality_ == reco::TrackBase::undefQuality) || tk.track().quality(quality_);
0037 bool nStripHitsCut = tk.hitPattern().numberOfValidStripHits() >= minStripHits_;
0038
0039 return IPSigCut && pTCut && etaCut && normChi2Cut && nPxLayCut && nSiLayCut && trackQualityCut && nStripHitsCut;
0040 }
0041
0042
0043 std::vector<reco::TransientTrack> TrackFilterForPVFinding::select(
0044 const std::vector<reco::TransientTrack>& tracks) const {
0045 std::vector<reco::TransientTrack> seltks;
0046 for (std::vector<reco::TransientTrack>::const_iterator itk = tracks.begin(); itk != tracks.end(); itk++) {
0047 if (operator()(*itk))
0048 seltks.push_back(*itk);
0049 }
0050 return seltks;
0051 }
0052
0053
0054 std::vector<reco::TransientTrack> TrackFilterForPVFinding::selectTight(const std::vector<reco::TransientTrack>& tracks,
0055 double minPtTight) const {
0056 std::vector<reco::TransientTrack> seltks;
0057 for (std::vector<reco::TransientTrack>::const_iterator itk = tracks.begin(); itk != tracks.end(); itk++) {
0058 if (itk->impactPointState().globalMomentum().transverse() < minPtTight)
0059 continue;
0060 if (operator()(*itk))
0061 seltks.push_back(*itk);
0062 }
0063 return seltks;
0064 }
0065
0066 void TrackFilterForPVFinding::fillPSetDescription(edm::ParameterSetDescription& desc) {
0067 desc.add<double>("maxNormalizedChi2", 10.0);
0068 desc.add<double>("minPt", 0.0);
0069 desc.add<std::string>("algorithm", "filter");
0070 desc.add<double>("maxEta", 2.4);
0071 desc.add<double>("maxD0Significance", 4.0);
0072 desc.add<double>("maxD0Error", 1.0);
0073 desc.add<double>("maxDzError", 1.0);
0074 desc.add<std::string>("trackQuality", "any");
0075 desc.add<int>("minPixelLayersWithHits", 2);
0076 desc.add<int>("minSiliconLayersWithHits", 5);
0077 desc.add<int>("minValidStripHits", 0);
0078 }