File indexing completed on 2023-03-17 11:17:13
0001 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0002
0003 #include "DataFormats/Candidate/interface/Candidate.h"
0004 #include "DataFormats/VertexReco/interface/Vertex.h"
0005 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
0006
0007 #include "RecoBTag/SecondaryVertex/interface/TrackSelector.h"
0008
0009 using namespace reco;
0010
0011 TrackSelector::TrackSelector(const edm::ParameterSet ¶ms)
0012 : minPixelHits(params.getParameter<unsigned int>("pixelHitsMin")),
0013 minTotalHits(params.getParameter<unsigned int>("totalHitsMin")),
0014 minPt(params.getParameter<double>("ptMin")),
0015 maxNormChi2(params.getParameter<double>("normChi2Max")),
0016 maxJetDeltaR(params.getParameter<double>("jetDeltaRMax")),
0017 maxDistToAxis(params.getParameter<double>("maxDistToAxis")),
0018 maxDecayLen(params.getParameter<double>("maxDecayLen")),
0019 sip2dValMin(params.getParameter<double>("sip2dValMin")),
0020 sip2dValMax(params.getParameter<double>("sip2dValMax")),
0021 sip2dSigMin(params.getParameter<double>("sip2dSigMin")),
0022 sip2dSigMax(params.getParameter<double>("sip2dSigMax")),
0023 sip3dValMin(params.getParameter<double>("sip3dValMin")),
0024 sip3dValMax(params.getParameter<double>("sip3dValMax")),
0025 sip3dSigMin(params.getParameter<double>("sip3dSigMin")),
0026 sip3dSigMax(params.getParameter<double>("sip3dSigMax")),
0027 useVariableJTA_(params.existsAs<bool>("useVariableJTA") ? params.getParameter<bool>("useVariableJTA") : false) {
0028 std::string qualityClass = params.getParameter<std::string>("qualityClass");
0029 if (qualityClass == "any" || qualityClass == "Any" || qualityClass == "ANY" || qualityClass.empty()) {
0030 selectQuality = false;
0031 quality = reco::TrackBase::undefQuality;
0032 } else {
0033 selectQuality = true;
0034 quality = reco::TrackBase::qualityByName(qualityClass);
0035 }
0036 if (useVariableJTA_) {
0037 varJTApars = {params.getParameter<double>("a_dR"),
0038 params.getParameter<double>("b_dR"),
0039 params.getParameter<double>("a_pT"),
0040 params.getParameter<double>("b_pT"),
0041 params.getParameter<double>("min_pT"),
0042 params.getParameter<double>("max_pT"),
0043 params.getParameter<double>("min_pT_dRcut"),
0044 params.getParameter<double>("max_pT_dRcut"),
0045 params.getParameter<double>("max_pT_trackPTcut")};
0046 }
0047 }
0048
0049 bool TrackSelector::operator()(const Track &track,
0050 const btag::TrackIPData &ipData,
0051 const Jet &jet,
0052 const GlobalPoint &pv) const {
0053 bool jtaPassed = false;
0054 if (useVariableJTA_) {
0055 jtaPassed = TrackIPTagInfo::passVariableJTA(
0056 varJTApars, jet.pt(), track.pt(), reco::deltaR(jet.momentum(), track.momentum()));
0057 } else
0058 jtaPassed = true;
0059
0060 return track.pt() >= minPt && reco::deltaR2(jet.momentum(), track.momentum()) < maxJetDeltaR * maxJetDeltaR &&
0061 jtaPassed && trackSelection(track, ipData, jet, pv);
0062 }
0063
0064 bool TrackSelector::operator()(const CandidatePtr &track,
0065 const btag::TrackIPData &ipData,
0066 const Jet &jet,
0067 const GlobalPoint &pv) const {
0068 bool jtaPassed = false;
0069 if (useVariableJTA_) {
0070 jtaPassed = TrackIPTagInfo::passVariableJTA(
0071 varJTApars, jet.pt(), track->pt(), reco::deltaR(jet.momentum(), track->momentum()));
0072 } else
0073 jtaPassed = true;
0074
0075 return track->pt() >= minPt && reco::deltaR2(jet.momentum(), track->momentum()) < maxJetDeltaR * maxJetDeltaR &&
0076 jtaPassed && trackSelection(*reco::btag::toTrack(track), ipData, jet, pv);
0077 }
0078
0079 bool TrackSelector::trackSelection(const Track &track,
0080 const btag::TrackIPData &ipData,
0081 const Jet &jet,
0082 const GlobalPoint &pv) const {
0083 return (!selectQuality || track.quality(quality)) &&
0084 (minPixelHits <= 0 || track.hitPattern().numberOfValidPixelHits() >= (int)minPixelHits) &&
0085 (minTotalHits <= 0 || track.hitPattern().numberOfValidHits() >= (int)minTotalHits) &&
0086 track.normalizedChi2() < maxNormChi2 && std::abs(ipData.distanceToJetAxis.value()) <= maxDistToAxis &&
0087 (ipData.closestToJetAxis - pv).mag() <= maxDecayLen && ipData.ip2d.value() >= sip2dValMin &&
0088 ipData.ip2d.value() <= sip2dValMax && ipData.ip2d.significance() >= sip2dSigMin &&
0089 ipData.ip2d.significance() <= sip2dSigMax && ipData.ip3d.value() >= sip3dValMin &&
0090 ipData.ip3d.value() <= sip3dValMax && ipData.ip3d.significance() >= sip3dSigMin &&
0091 ipData.ip3d.significance() <= sip3dSigMax;
0092 }