File indexing completed on 2024-04-06 12:24:34
0001 #include "RecoBTag/SecondaryVertex/interface/CombinedSVComputer.h"
0002
0003 using namespace reco;
0004
0005 inline edm::ParameterSet CombinedSVComputer::dropDeltaR(const edm::ParameterSet &pset) const {
0006 edm::ParameterSet psetCopy(pset);
0007 psetCopy.addParameter<double>("jetDeltaRMax", 99999.0);
0008 return psetCopy;
0009 }
0010
0011 CombinedSVComputer::CombinedSVComputer(const edm::ParameterSet ¶ms)
0012 : trackFlip(params.getParameter<bool>("trackFlip")),
0013 vertexFlip(params.getParameter<bool>("vertexFlip")),
0014 charmCut(params.getParameter<double>("charmCut")),
0015 sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))),
0016 trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
0017 trackNoDeltaRSelector(dropDeltaR(params.getParameter<edm::ParameterSet>("trackSelection"))),
0018 trackPseudoSelector(params.getParameter<edm::ParameterSet>("trackPseudoSelection")),
0019 pseudoMultiplicityMin(params.getParameter<unsigned int>("pseudoMultiplicityMin")),
0020 trackMultiplicityMin(params.getParameter<unsigned int>("trackMultiplicityMin")),
0021 minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
0022 useTrackWeights(params.getParameter<bool>("useTrackWeights")),
0023 vertexMassCorrection(params.getParameter<bool>("correctVertexMass")),
0024 pseudoVertexV0Filter(params.getParameter<edm::ParameterSet>("pseudoVertexV0Filter")),
0025 trackPairV0Filter(params.getParameter<edm::ParameterSet>("trackPairV0Filter")) {}
0026
0027 inline double CombinedSVComputer::flipValue(double value, bool vertex) const {
0028 return (vertex ? vertexFlip : trackFlip) ? -value : value;
0029 }
0030
0031 inline CombinedSVComputer::IterationRange CombinedSVComputer::flipIterate(int size, bool vertex) const {
0032 IterationRange range;
0033 if (vertex ? vertexFlip : trackFlip) {
0034 range.begin = size - 1;
0035 range.end = -1;
0036 range.increment = -1;
0037 } else {
0038 range.begin = 0;
0039 range.end = size;
0040 range.increment = +1;
0041 }
0042
0043 return range;
0044 }
0045
0046 const btag::TrackIPData &CombinedSVComputer::threshTrack(const CandIPTagInfo &trackIPTagInfo,
0047 const btag::SortCriteria sort,
0048 const reco::Jet &jet,
0049 const GlobalPoint &pv) const {
0050 const CandIPTagInfo::input_container &tracks = trackIPTagInfo.selectedTracks();
0051 const std::vector<btag::TrackIPData> &ipData = trackIPTagInfo.impactParameterData();
0052 std::vector<std::size_t> indices = trackIPTagInfo.sortedIndexes(sort);
0053
0054 IterationRange range = flipIterate(indices.size(), false);
0055 TrackKinematics kin;
0056 range_for(i, range) {
0057 std::size_t idx = indices[i];
0058 const btag::TrackIPData &data = ipData[idx];
0059 const CandidatePtr &track = tracks[idx];
0060
0061 if (!trackNoDeltaRSelector(track, data, jet, pv))
0062 continue;
0063
0064 kin.add(track);
0065 if (kin.vectorSum().M() > charmCut)
0066 return data;
0067 }
0068 if (trackFlip) {
0069 static const btag::TrackIPData dummy = {GlobalPoint(),
0070 GlobalPoint(),
0071 Measurement1D(1.0, 1.0),
0072 Measurement1D(1.0, 1.0),
0073 Measurement1D(1.0, 1.0),
0074 Measurement1D(1.0, 1.0),
0075 0.};
0076 return dummy;
0077 } else {
0078 static const btag::TrackIPData dummy = {GlobalPoint(),
0079 GlobalPoint(),
0080 Measurement1D(-1.0, 1.0),
0081 Measurement1D(-1.0, 1.0),
0082 Measurement1D(-1.0, 1.0),
0083 Measurement1D(-1.0, 1.0),
0084 0.};
0085 return dummy;
0086 }
0087 }
0088
0089 const btag::TrackIPData &CombinedSVComputer::threshTrack(const TrackIPTagInfo &trackIPTagInfo,
0090 const btag::SortCriteria sort,
0091 const reco::Jet &jet,
0092 const GlobalPoint &pv) const {
0093 const edm::RefVector<TrackCollection> &tracks = trackIPTagInfo.selectedTracks();
0094 const std::vector<btag::TrackIPData> &ipData = trackIPTagInfo.impactParameterData();
0095 std::vector<std::size_t> indices = trackIPTagInfo.sortedIndexes(sort);
0096
0097 IterationRange range = flipIterate(indices.size(), false);
0098 TrackKinematics kin;
0099 range_for(i, range) {
0100 std::size_t idx = indices[i];
0101 const btag::TrackIPData &data = ipData[idx];
0102 const Track &track = *tracks[idx];
0103
0104 if (!trackNoDeltaRSelector(track, data, jet, pv))
0105 continue;
0106
0107 kin.add(track);
0108 if (kin.vectorSum().M() > charmCut)
0109 return data;
0110 }
0111
0112 if (trackFlip) {
0113 static const btag::TrackIPData dummy = {GlobalPoint(),
0114 GlobalPoint(),
0115 Measurement1D(1.0, 1.0),
0116 Measurement1D(1.0, 1.0),
0117 Measurement1D(1.0, 1.0),
0118 Measurement1D(1.0, 1.0),
0119 0.};
0120 return dummy;
0121 } else {
0122 static const btag::TrackIPData dummy = {GlobalPoint(),
0123 GlobalPoint(),
0124 Measurement1D(-1.0, 1.0),
0125 Measurement1D(-1.0, 1.0),
0126 Measurement1D(-1.0, 1.0),
0127 Measurement1D(-1.0, 1.0),
0128 0.};
0129 return dummy;
0130 }
0131 }
0132
0133 TaggingVariableList CombinedSVComputer::operator()(const TrackIPTagInfo &ipInfo,
0134 const SecondaryVertexTagInfo &svInfo) const {
0135 using namespace ROOT::Math;
0136
0137 edm::RefToBase<Jet> jet = ipInfo.jet();
0138 math::XYZVector jetDir = jet->momentum().Unit();
0139 TaggingVariableList vars;
0140
0141 TrackKinematics vertexKinematics;
0142
0143 double vtx_track_ptSum = 0.;
0144 double vtx_track_ESum = 0.;
0145
0146
0147 int vtx = -1;
0148 unsigned int numberofvertextracks = 0;
0149
0150 IterationRange range = flipIterate(svInfo.nVertices(), true);
0151 range_for(i, range) {
0152 numberofvertextracks = numberofvertextracks + (svInfo.secondaryVertex(i)).nTracks();
0153
0154 const Vertex &vertex = svInfo.secondaryVertex(i);
0155 bool hasRefittedTracks = vertex.hasRefittedTracks();
0156 for (reco::Vertex::trackRef_iterator track = vertex.tracks_begin(); track != vertex.tracks_end(); track++) {
0157 double w = vertex.trackWeight(*track);
0158 if (w < minTrackWeight)
0159 continue;
0160 if (hasRefittedTracks) {
0161 const Track actualTrack = vertex.refittedTrack(*track);
0162 vertexKinematics.add(actualTrack, w);
0163 vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir, actualTrack.momentum()), true);
0164 if (vtx < 0)
0165 {
0166 vtx_track_ptSum += std::sqrt(actualTrack.momentum().Perp2());
0167 vtx_track_ESum += std::sqrt(actualTrack.momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
0168 }
0169 } else {
0170 vertexKinematics.add(**track, w);
0171 vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir, (*track)->momentum()), true);
0172 if (vtx < 0)
0173 {
0174 vtx_track_ptSum += std::sqrt((*track)->momentum().Perp2());
0175 vtx_track_ESum += std::sqrt((*track)->momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
0176 }
0177 }
0178 }
0179
0180 if (vtx < 0)
0181 vtx = i;
0182 }
0183 if (vtx >= 0) {
0184 vars.insert(btau::vertexNTracks, numberofvertextracks, true);
0185 vars.insert(btau::vertexFitProb, (svInfo.secondaryVertex(vtx)).normalizedChi2(), true);
0186 }
0187
0188
0189 fillCommonVariables(vars, vertexKinematics, ipInfo, svInfo, vtx_track_ptSum, vtx_track_ESum);
0190
0191 vars.finalize();
0192 return vars;
0193 }
0194
0195 TaggingVariableList CombinedSVComputer::operator()(const CandIPTagInfo &ipInfo,
0196 const CandSecondaryVertexTagInfo &svInfo) const {
0197 using namespace ROOT::Math;
0198
0199 edm::RefToBase<Jet> jet = ipInfo.jet();
0200 math::XYZVector jetDir = jet->momentum().Unit();
0201 TaggingVariableList vars;
0202
0203 TrackKinematics vertexKinematics;
0204
0205 double vtx_track_ptSum = 0.;
0206 double vtx_track_ESum = 0.;
0207
0208
0209 int vtx = -1;
0210 unsigned int numberofvertextracks = 0;
0211
0212 IterationRange range = flipIterate(svInfo.nVertices(), true);
0213 range_for(i, range) {
0214 numberofvertextracks = numberofvertextracks + (svInfo.secondaryVertex(i)).numberOfSourceCandidatePtrs();
0215
0216 const reco::VertexCompositePtrCandidate &vertex = svInfo.secondaryVertex(i);
0217 const std::vector<CandidatePtr> &tracks = vertex.daughterPtrVector();
0218 for (std::vector<CandidatePtr>::const_iterator track = tracks.begin(); track != tracks.end(); ++track) {
0219 vertexKinematics.add(*track);
0220 vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir, (*track)->momentum()), true);
0221 if (vtx < 0)
0222 {
0223 vtx_track_ptSum += std::sqrt((*track)->momentum().Perp2());
0224 vtx_track_ESum += std::sqrt((*track)->momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
0225 }
0226 }
0227
0228 if (vtx < 0)
0229 vtx = i;
0230 }
0231 if (vtx >= 0) {
0232 vars.insert(btau::vertexNTracks, numberofvertextracks, true);
0233 vars.insert(btau::vertexFitProb, (svInfo.secondaryVertex(vtx)).vertexNormalizedChi2(), true);
0234 }
0235
0236
0237 fillCommonVariables(vars, vertexKinematics, ipInfo, svInfo, vtx_track_ptSum, vtx_track_ESum);
0238
0239 vars.finalize();
0240 return vars;
0241 }