File indexing completed on 2024-04-06 12:24:35
0001 #include <functional>
0002 #include <algorithm>
0003 #include <iterator>
0004 #include <cmath>
0005 #include <set>
0006
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008
0009 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0010 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
0011 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0012 #include "DataFormats/VertexReco/interface/Vertex.h"
0013
0014 #include "RecoBTag/SecondaryVertex/interface/SecondaryVertex.h"
0015 #include "RecoBTag/SecondaryVertex/interface/TrackKinematics.h"
0016 #include "RecoBTag/SecondaryVertex/interface/V0Filter.h"
0017 #include "RecoBTag/SecondaryVertex/interface/VertexFilter.h"
0018 #include "RecoVertex/VertexTools/interface/SharedTracks.h"
0019 using namespace reco;
0020
0021 VertexFilter::VertexFilter(const edm::ParameterSet ¶ms)
0022 : useTrackWeights(params.getParameter<bool>("useTrackWeights")),
0023 minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
0024 massMax(params.getParameter<double>("massMax")),
0025 fracPV(params.getParameter<double>("fracPV")),
0026 multiplicityMin(params.getParameter<unsigned int>("multiplicityMin")),
0027 distVal2dMin(params.getParameter<double>("distVal2dMin")),
0028 distVal2dMax(params.getParameter<double>("distVal2dMax")),
0029 distVal3dMin(params.getParameter<double>("distVal3dMin")),
0030 distVal3dMax(params.getParameter<double>("distVal3dMax")),
0031 distSig2dMin(params.getParameter<double>("distSig2dMin")),
0032 distSig2dMax(params.getParameter<double>("distSig2dMax")),
0033 distSig3dMin(params.getParameter<double>("distSig3dMin")),
0034 distSig3dMax(params.getParameter<double>("distSig3dMax")),
0035 maxDeltaRToJetAxis(params.getParameter<double>("maxDeltaRToJetAxis")),
0036 v0Filter(params.getParameter<edm::ParameterSet>("v0Filter")) {}
0037
0038 bool VertexFilter::operator()(const Vertex &pv,
0039 const TemplatedSecondaryVertex<reco::Vertex> &sv,
0040 const GlobalVector &direction) const {
0041 std::vector<TrackRef> svTracks;
0042 for (std::vector<TrackBaseRef>::const_iterator iter = sv.tracks_begin(); iter != sv.tracks_end(); iter++)
0043 if (sv.trackWeight(*iter) >= minTrackWeight)
0044 svTracks.push_back(iter->castTo<TrackRef>());
0045
0046
0047
0048 if (svTracks.size() < multiplicityMin)
0049 return false;
0050
0051
0052
0053 if (sv.dist2d().error() < 0 || sv.dist3d().error() < 0)
0054 return false;
0055
0056
0057
0058 if (sv.dist2d().value() < distVal2dMin || sv.dist2d().value() > distVal2dMax || sv.dist3d().value() < distVal3dMin ||
0059 sv.dist3d().value() > distVal3dMax || sv.dist2d().significance() < distSig2dMin ||
0060 sv.dist2d().significance() > distSig2dMax || sv.dist3d().significance() < distSig3dMin ||
0061 sv.dist3d().significance() > distSig3dMax)
0062 return false;
0063
0064
0065
0066 if (Geom::deltaR(sv.position() - pv.position(), (maxDeltaRToJetAxis > 0) ? direction : -direction) >
0067 std::abs(maxDeltaRToJetAxis))
0068 return false;
0069
0070
0071
0072 TrackKinematics kin(sv);
0073
0074 double mass = useTrackWeights ? kin.weightedVectorSum().M() : kin.vectorSum().M();
0075
0076 if (mass > massMax)
0077 return false;
0078
0079
0080
0081 if (fracPV < 1.0) {
0082 double fractionSharedTracks = vertexTools::computeSharedTracks(pv, svTracks, minTrackWeight);
0083 if (fractionSharedTracks > fracPV)
0084 return false;
0085 }
0086
0087
0088
0089 if (sv.hasRefittedTracks())
0090 return v0Filter(sv.refittedTracks());
0091 else
0092 return v0Filter(svTracks);
0093 }
0094 bool VertexFilter::operator()(const reco::Vertex &pv,
0095 const TemplatedSecondaryVertex<reco::VertexCompositePtrCandidate> &sv,
0096 const GlobalVector &direction) const {
0097 const std::vector<CandidatePtr> &svTracks = sv.daughterPtrVector();
0098
0099
0100
0101 if (svTracks.size() < multiplicityMin)
0102 return false;
0103
0104
0105
0106 if (sv.dist2d().error() < 0 || sv.dist3d().error() < 0)
0107 return false;
0108
0109
0110
0111 if (sv.dist2d().value() < distVal2dMin || sv.dist2d().value() > distVal2dMax || sv.dist3d().value() < distVal3dMin ||
0112 sv.dist3d().value() > distVal3dMax || sv.dist2d().significance() < distSig2dMin ||
0113 sv.dist2d().significance() > distSig2dMax || sv.dist3d().significance() < distSig3dMin ||
0114 sv.dist3d().significance() > distSig3dMax)
0115 return false;
0116
0117
0118
0119 if (Geom::deltaR(sv.vertex() - pv.position(), (maxDeltaRToJetAxis > 0) ? direction : -direction) >
0120 std::abs(maxDeltaRToJetAxis))
0121 return false;
0122
0123
0124 TrackKinematics kin(sv);
0125
0126 double mass = useTrackWeights ? kin.weightedVectorSum().M() : kin.vectorSum().M();
0127
0128 if (mass > massMax)
0129 return false;
0130
0131
0132
0133 if (fracPV < 1.0) {
0134 double fractionSharedTracks = vertexTools::computeSharedTracks(pv, svTracks, minTrackWeight);
0135 if (fractionSharedTracks > fracPV)
0136 return false;
0137 }
0138
0139
0140
0141 return v0Filter(svTracks);
0142 }