Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:34

0001 #include <iostream>
0002 #include <cstddef>
0003 #include <string>
0004 #include <cmath>
0005 #include <vector>
0006 
0007 #include <Math/VectorUtil.h>
0008 
0009 #include "FWCore/Utilities/interface/Exception.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 
0012 #include "DataFormats/Math/interface/Vector3D.h"
0013 #include "DataFormats/Math/interface/LorentzVector.h"
0014 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
0015 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0016 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
0017 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
0018 #include "DataFormats/TrackReco/interface/Track.h"
0019 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0020 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
0021 #include "DataFormats/BTauReco/interface/SecondaryVertexTagInfo.h"
0022 #include "DataFormats/BTauReco/interface/TaggingVariable.h"
0023 #include "DataFormats/BTauReco/interface/VertexTypes.h"
0024 #include "DataFormats/BTauReco/interface/ParticleMasses.h"
0025 
0026 #include "RecoBTag/SecondaryVertex/interface/TrackSorting.h"
0027 #include "RecoBTag/SecondaryVertex/interface/TrackSelector.h"
0028 #include "RecoBTag/SecondaryVertex/interface/TrackKinematics.h"
0029 #include "RecoBTag/SecondaryVertex/interface/V0Filter.h"
0030 
0031 #include "RecoBTag/SecondaryVertex/interface/GhostTrackComputer.h"
0032 
0033 using namespace reco;
0034 
0035 static edm::ParameterSet dropDeltaR(const edm::ParameterSet &pset) {
0036   edm::ParameterSet psetCopy(pset);
0037   psetCopy.addParameter<double>("jetDeltaRMax", 99999.0);
0038   return psetCopy;
0039 }
0040 
0041 GhostTrackComputer::GhostTrackComputer(const edm::ParameterSet &params)
0042     : charmCut(params.getParameter<double>("charmCut")),
0043       sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))),
0044       trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
0045       trackNoDeltaRSelector(dropDeltaR(params.getParameter<edm::ParameterSet>("trackSelection"))),
0046       minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
0047       trackPairV0Filter(params.getParameter<edm::ParameterSet>("trackPairV0Filter")) {}
0048 
0049 const btag::TrackIPData &GhostTrackComputer::threshTrack(const TrackIPTagInfo &trackIPTagInfo,
0050                                                          const btag::SortCriteria sort,
0051                                                          const reco::Jet &jet,
0052                                                          const GlobalPoint &pv) const {
0053   const edm::RefVector<TrackCollection> &tracks = trackIPTagInfo.selectedTracks();
0054   const std::vector<btag::TrackIPData> &ipData = trackIPTagInfo.impactParameterData();
0055   std::vector<std::size_t> indices = trackIPTagInfo.sortedIndexes(sort);
0056 
0057   TrackKinematics kin;
0058   for (std::vector<std::size_t>::const_iterator iter = indices.begin(); iter != indices.end(); ++iter) {
0059     const btag::TrackIPData &data = ipData[*iter];
0060     const Track &track = *tracks[*iter];
0061 
0062     if (!trackNoDeltaRSelector(track, data, jet, pv))
0063       continue;
0064 
0065     kin.add(track);
0066     if (kin.vectorSum().M() > charmCut)
0067       return data;
0068   }
0069 
0070   static const btag::TrackIPData dummy = {GlobalPoint(),
0071                                           GlobalPoint(),
0072                                           Measurement1D(-1.0, 1.0),
0073                                           Measurement1D(-1.0, 1.0),
0074                                           Measurement1D(-1.0, 1.0),
0075                                           Measurement1D(-1.0, 1.0),
0076                                           0.};
0077   return dummy;
0078 }
0079 
0080 const btag::TrackIPData &GhostTrackComputer::threshTrack(const CandIPTagInfo &trackIPTagInfo,
0081                                                          const btag::SortCriteria sort,
0082                                                          const reco::Jet &jet,
0083                                                          const GlobalPoint &pv) const {
0084   const CandIPTagInfo::input_container &tracks = trackIPTagInfo.selectedTracks();
0085   const std::vector<btag::TrackIPData> &ipData = trackIPTagInfo.impactParameterData();
0086   std::vector<std::size_t> indices = trackIPTagInfo.sortedIndexes(sort);
0087 
0088   TrackKinematics kin;
0089   for (std::vector<std::size_t>::const_iterator iter = indices.begin(); iter != indices.end(); ++iter) {
0090     const btag::TrackIPData &data = ipData[*iter];
0091     const CandidatePtr &track = tracks[*iter];
0092 
0093     if (!trackNoDeltaRSelector(track, data, jet, pv))
0094       continue;
0095 
0096     kin.add(track);
0097     if (kin.vectorSum().M() > charmCut)
0098       return data;
0099   }
0100 
0101   static const btag::TrackIPData dummy = {GlobalPoint(),
0102                                           GlobalPoint(),
0103                                           Measurement1D(-1.0, 1.0),
0104                                           Measurement1D(-1.0, 1.0),
0105                                           Measurement1D(-1.0, 1.0),
0106                                           Measurement1D(-1.0, 1.0),
0107                                           0.};
0108   return dummy;
0109 }
0110 
0111 static void addMeas(std::pair<double, double> &sum, Measurement1D meas) {
0112   double weight = 1. / meas.error();
0113   weight *= weight;
0114   sum.first += weight * meas.value();
0115   sum.second += weight;
0116 }
0117 
0118 TaggingVariableList GhostTrackComputer::operator()(const TrackIPTagInfo &ipInfo,
0119                                                    const SecondaryVertexTagInfo &svInfo) const {
0120   using namespace ROOT::Math;
0121 
0122   edm::RefToBase<Jet> jet = ipInfo.jet();
0123   math::XYZVector jetDir = jet->momentum().Unit();
0124   bool havePv = ipInfo.primaryVertex().isNonnull();
0125   GlobalPoint pv;
0126   if (havePv)
0127     pv = GlobalPoint(ipInfo.primaryVertex()->x(), ipInfo.primaryVertex()->y(), ipInfo.primaryVertex()->z());
0128 
0129   btag::Vertices::VertexType vtxType = btag::Vertices::NoVertex;
0130 
0131   TaggingVariableList vars;
0132 
0133   vars.insert(btau::jetPt, jet->pt(), true);
0134   vars.insert(btau::jetEta, jet->eta(), true);
0135 
0136   TrackKinematics allKinematics;
0137   TrackKinematics vertexKinematics;
0138   TrackKinematics trackKinematics;
0139 
0140   std::pair<double, double> vertexDist2D, vertexDist3D;
0141   std::pair<double, double> tracksDist2D, tracksDist3D;
0142 
0143   unsigned int nVertices = 0;
0144   unsigned int nVertexTracks = 0;
0145   unsigned int nTracks = 0;
0146   for (unsigned int i = 0; i < svInfo.nVertices(); i++) {
0147     const Vertex &vertex = svInfo.secondaryVertex(i);
0148     bool hasRefittedTracks = vertex.hasRefittedTracks();
0149     TrackRefVector tracks = svInfo.vertexTracks(i);
0150     unsigned int n = 0;
0151     for (TrackRefVector::const_iterator track = tracks.begin(); track != tracks.end(); track++)
0152       if (svInfo.trackWeight(i, *track) >= minTrackWeight)
0153         n++;
0154 
0155     if (n < 1)
0156       continue;
0157     bool isTrackVertex = (n == 1);
0158     ++*(isTrackVertex ? &nTracks : &nVertices);
0159 
0160     addMeas(*(isTrackVertex ? &tracksDist2D : &vertexDist2D), svInfo.flightDistance(i, true));
0161     addMeas(*(isTrackVertex ? &tracksDist3D : &vertexDist3D), svInfo.flightDistance(i, false));
0162 
0163     TrackKinematics &kin = isTrackVertex ? trackKinematics : vertexKinematics;
0164     for (TrackRefVector::const_iterator track = tracks.begin(); track != tracks.end(); track++) {
0165       float w = svInfo.trackWeight(i, *track);
0166       if (w < minTrackWeight)
0167         continue;
0168       if (hasRefittedTracks) {
0169         Track actualTrack = vertex.refittedTrack(*track);
0170         kin.add(actualTrack, w);
0171         vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir, actualTrack.momentum()), true);
0172       } else {
0173         kin.add(**track, w);
0174         vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir, (*track)->momentum()), true);
0175       }
0176       if (!isTrackVertex)
0177         nVertexTracks++;
0178     }
0179   }
0180 
0181   Measurement1D dist2D, dist3D;
0182   if (nVertices) {
0183     vtxType = btag::Vertices::RecoVertex;
0184 
0185     if (nVertices == 1 && nTracks) {
0186       vertexDist2D.first += tracksDist2D.first;
0187       vertexDist2D.second += tracksDist2D.second;
0188       vertexDist3D.first += tracksDist3D.first;
0189       vertexDist3D.second += tracksDist3D.second;
0190       vertexKinematics += trackKinematics;
0191     }
0192 
0193     dist2D = Measurement1D(vertexDist2D.first / vertexDist2D.second, std::sqrt(1. / vertexDist2D.second));
0194     dist3D = Measurement1D(vertexDist3D.first / vertexDist3D.second, std::sqrt(1. / vertexDist3D.second));
0195 
0196     vars.insert(btau::jetNSecondaryVertices, nVertices, true);
0197     vars.insert(btau::vertexNTracks, nVertexTracks, true);
0198   } else if (nTracks) {
0199     vtxType = btag::Vertices::PseudoVertex;
0200     vertexKinematics = trackKinematics;
0201 
0202     dist2D = Measurement1D(tracksDist2D.first / tracksDist2D.second, std::sqrt(1. / tracksDist2D.second));
0203     dist3D = Measurement1D(tracksDist3D.first / tracksDist3D.second, std::sqrt(1. / tracksDist3D.second));
0204   }
0205 
0206   if (nVertices || nTracks) {
0207     vars.insert(btau::flightDistance2dVal, dist2D.value(), true);
0208     vars.insert(btau::flightDistance2dSig, dist2D.significance(), true);
0209     vars.insert(btau::flightDistance3dVal, dist3D.value(), true);
0210     vars.insert(btau::flightDistance3dSig, dist3D.significance(), true);
0211     vars.insert(btau::vertexJetDeltaR, Geom::deltaR(svInfo.flightDirection(0), jetDir), true);
0212     vars.insert(btau::jetNSingleTrackVertices, nTracks, true);
0213   }
0214 
0215   std::vector<std::size_t> indices = ipInfo.sortedIndexes(sortCriterium);
0216   const std::vector<btag::TrackIPData> &ipData = ipInfo.impactParameterData();
0217   const edm::RefVector<TrackCollection> &tracks = ipInfo.selectedTracks();
0218 
0219   TrackRef trackPairV0Test[2];
0220   for (unsigned int i = 0; i < indices.size(); i++) {
0221     std::size_t idx = indices[i];
0222     const btag::TrackIPData &data = ipData[idx];
0223     const TrackRef &trackRef = tracks[idx];
0224     const Track &track = *trackRef;
0225 
0226     // filter track
0227 
0228     if (!trackSelector(track, data, *jet, pv))
0229       continue;
0230 
0231     // add track to kinematics for all tracks in jet
0232 
0233     allKinematics.add(track);
0234 
0235     // check against all other tracks for V0 track pairs
0236 
0237     trackPairV0Test[0] = tracks[idx];
0238     bool ok = true;
0239     for (unsigned int j = 0; j < indices.size(); j++) {
0240       if (i == j)
0241         continue;
0242 
0243       std::size_t pairIdx = indices[j];
0244       const btag::TrackIPData &pairTrackData = ipData[pairIdx];
0245       const TrackRef &pairTrackRef = tracks[pairIdx];
0246       const Track &pairTrack = *pairTrackRef;
0247 
0248       if (!trackSelector(pairTrack, pairTrackData, *jet, pv))
0249         continue;
0250 
0251       trackPairV0Test[1] = pairTrackRef;
0252       if (!trackPairV0Filter(trackPairV0Test, 2)) {
0253         ok = false;
0254         break;
0255       }
0256     }
0257     if (!ok)
0258       continue;
0259 
0260     // add track variables
0261 
0262     const math::XYZVector &trackMom = track.momentum();
0263     double trackMag = std::sqrt(trackMom.Mag2());
0264 
0265     vars.insert(btau::trackSip3dVal, data.ip3d.value(), true);
0266     vars.insert(btau::trackSip3dSig, data.ip3d.significance(), true);
0267     vars.insert(btau::trackSip2dVal, data.ip2d.value(), true);
0268     vars.insert(btau::trackSip2dSig, data.ip2d.significance(), true);
0269     vars.insert(btau::trackJetDistVal, data.distanceToJetAxis.value(), true);
0270     vars.insert(btau::trackGhostTrackDistVal, data.distanceToGhostTrack.value(), true);
0271     vars.insert(btau::trackGhostTrackDistSig, data.distanceToGhostTrack.significance(), true);
0272     vars.insert(btau::trackGhostTrackWeight, data.ghostTrackWeight, true);
0273     vars.insert(btau::trackDecayLenVal, havePv ? (data.closestToGhostTrack - pv).mag() : -1.0, true);
0274 
0275     vars.insert(btau::trackMomentum, trackMag, true);
0276     vars.insert(btau::trackEta, trackMom.Eta(), true);
0277 
0278     vars.insert(btau::trackChi2, track.normalizedChi2(), true);
0279     vars.insert(btau::trackNPixelHits, track.hitPattern().pixelLayersWithMeasurement(), true);
0280     vars.insert(btau::trackNTotalHits, track.hitPattern().trackerLayersWithMeasurement(), true);
0281     vars.insert(btau::trackPtRel, VectorUtil::Perp(trackMom, jetDir), true);
0282     vars.insert(btau::trackDeltaR, VectorUtil::DeltaR(trackMom, jetDir), true);
0283   }
0284 
0285   vars.insert(btau::vertexCategory, vtxType, true);
0286   vars.insert(btau::trackSumJetDeltaR, VectorUtil::DeltaR(allKinematics.vectorSum(), jetDir), true);
0287   vars.insert(btau::trackSumJetEtRatio, allKinematics.vectorSum().Et() / ipInfo.jet()->et(), true);
0288   vars.insert(
0289       btau::trackSip3dSigAboveCharm, threshTrack(ipInfo, reco::btag::IP3DSig, *jet, pv).ip3d.significance(), true);
0290   vars.insert(
0291       btau::trackSip2dSigAboveCharm, threshTrack(ipInfo, reco::btag::IP2DSig, *jet, pv).ip2d.significance(), true);
0292 
0293   if (vtxType != btag::Vertices::NoVertex) {
0294     const math::XYZTLorentzVector &allSum = allKinematics.vectorSum();
0295     const math::XYZTLorentzVector &vertexSum = vertexKinematics.vectorSum();
0296 
0297     vars.insert(btau::vertexMass, vertexSum.M(), true);
0298     if (allKinematics.numberOfTracks())
0299       vars.insert(btau::vertexEnergyRatio, vertexSum.E() / allSum.E(), true);
0300     else
0301       vars.insert(btau::vertexEnergyRatio, 1, true);
0302   }
0303 
0304   vars.finalize();
0305 
0306   return vars;
0307 }
0308 
0309 TaggingVariableList GhostTrackComputer::operator()(const CandIPTagInfo &ipInfo,
0310                                                    const CandSecondaryVertexTagInfo &svInfo) const {
0311   using namespace ROOT::Math;
0312 
0313   edm::RefToBase<Jet> jet = ipInfo.jet();
0314   math::XYZVector jetDir = jet->momentum().Unit();
0315   bool havePv = ipInfo.primaryVertex().isNonnull();
0316   GlobalPoint pv;
0317   if (havePv)
0318     pv = GlobalPoint(ipInfo.primaryVertex()->x(), ipInfo.primaryVertex()->y(), ipInfo.primaryVertex()->z());
0319 
0320   btag::Vertices::VertexType vtxType = btag::Vertices::NoVertex;
0321 
0322   TaggingVariableList vars;
0323 
0324   vars.insert(btau::jetPt, jet->pt(), true);
0325   vars.insert(btau::jetEta, jet->eta(), true);
0326 
0327   TrackKinematics allKinematics;
0328   TrackKinematics vertexKinematics;
0329   TrackKinematics trackKinematics;
0330 
0331   std::pair<double, double> vertexDist2D, vertexDist3D;
0332   std::pair<double, double> tracksDist2D, tracksDist3D;
0333 
0334   unsigned int nVertices = 0;
0335   unsigned int nVertexTracks = 0;
0336   unsigned int nTracks = 0;
0337   for (unsigned int i = 0; i < svInfo.nVertices(); i++) {
0338     const reco::VertexCompositePtrCandidate &vertex = svInfo.secondaryVertex(i);
0339     const std::vector<CandidatePtr> &tracks = vertex.daughterPtrVector();
0340     unsigned int n = 0;
0341     for (std::vector<CandidatePtr>::const_iterator track = tracks.begin(); track != tracks.end(); ++track)
0342       n++;
0343 
0344     if (n < 1)
0345       continue;
0346     bool isTrackVertex = (n == 1);
0347     ++*(isTrackVertex ? &nTracks : &nVertices);
0348 
0349     addMeas(*(isTrackVertex ? &tracksDist2D : &vertexDist2D), svInfo.flightDistance(i, true));
0350     addMeas(*(isTrackVertex ? &tracksDist3D : &vertexDist3D), svInfo.flightDistance(i, false));
0351 
0352     TrackKinematics &kin = isTrackVertex ? trackKinematics : vertexKinematics;
0353 
0354     for (std::vector<CandidatePtr>::const_iterator track = tracks.begin(); track != tracks.end(); ++track) {
0355       kin.add(*track);
0356       vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir, (*track)->momentum()), true);
0357       if (!isTrackVertex)
0358         nVertexTracks++;
0359     }
0360   }
0361 
0362   Measurement1D dist2D, dist3D;
0363   if (nVertices) {
0364     vtxType = btag::Vertices::RecoVertex;
0365 
0366     if (nVertices == 1 && nTracks) {
0367       vertexDist2D.first += tracksDist2D.first;
0368       vertexDist2D.second += tracksDist2D.second;
0369       vertexDist3D.first += tracksDist3D.first;
0370       vertexDist3D.second += tracksDist3D.second;
0371       vertexKinematics += trackKinematics;
0372     }
0373 
0374     dist2D = Measurement1D(vertexDist2D.first / vertexDist2D.second, std::sqrt(1. / vertexDist2D.second));
0375     dist3D = Measurement1D(vertexDist3D.first / vertexDist3D.second, std::sqrt(1. / vertexDist3D.second));
0376 
0377     vars.insert(btau::jetNSecondaryVertices, nVertices, true);
0378     vars.insert(btau::vertexNTracks, nVertexTracks, true);
0379   } else if (nTracks) {
0380     vtxType = btag::Vertices::PseudoVertex;
0381     vertexKinematics = trackKinematics;
0382 
0383     dist2D = Measurement1D(tracksDist2D.first / tracksDist2D.second, std::sqrt(1. / tracksDist2D.second));
0384     dist3D = Measurement1D(tracksDist3D.first / tracksDist3D.second, std::sqrt(1. / tracksDist3D.second));
0385   }
0386 
0387   if (nVertices || nTracks) {
0388     vars.insert(btau::flightDistance2dVal, dist2D.value(), true);
0389     vars.insert(btau::flightDistance2dSig, dist2D.significance(), true);
0390     vars.insert(btau::flightDistance3dVal, dist3D.value(), true);
0391     vars.insert(btau::flightDistance3dSig, dist3D.significance(), true);
0392     vars.insert(btau::vertexJetDeltaR, Geom::deltaR(svInfo.flightDirection(0), jetDir), true);
0393     vars.insert(btau::jetNSingleTrackVertices, nTracks, true);
0394   }
0395 
0396   std::vector<std::size_t> indices = ipInfo.sortedIndexes(sortCriterium);
0397   const std::vector<btag::TrackIPData> &ipData = ipInfo.impactParameterData();
0398   const CandIPTagInfo::input_container &tracks = ipInfo.selectedTracks();
0399 
0400   const Track *trackPairV0Test[2];
0401   for (unsigned int i = 0; i < indices.size(); i++) {
0402     std::size_t idx = indices[i];
0403     const btag::TrackIPData &data = ipData[idx];
0404     const Track *trackPtr = reco::btag::toTrack(tracks[idx]);
0405     const Track &track = *trackPtr;
0406 
0407     // filter track
0408 
0409     if (!trackSelector(track, data, *jet, pv))
0410       continue;
0411 
0412     // add track to kinematics for all tracks in jet
0413 
0414     allKinematics.add(track);
0415 
0416     // check against all other tracks for V0 track pairs
0417 
0418     trackPairV0Test[0] = reco::btag::toTrack(tracks[idx]);
0419     bool ok = true;
0420     for (unsigned int j = 0; j < indices.size(); j++) {
0421       if (i == j)
0422         continue;
0423 
0424       std::size_t pairIdx = indices[j];
0425       const btag::TrackIPData &pairTrackData = ipData[pairIdx];
0426       const Track *pairTrackPtr = reco::btag::toTrack(tracks[pairIdx]);
0427       const Track &pairTrack = *pairTrackPtr;
0428 
0429       if (!trackSelector(pairTrack, pairTrackData, *jet, pv))
0430         continue;
0431 
0432       trackPairV0Test[1] = pairTrackPtr;
0433       if (!trackPairV0Filter(trackPairV0Test, 2)) {
0434         ok = false;
0435         break;
0436       }
0437     }
0438     if (!ok)
0439       continue;
0440 
0441     // add track variables
0442 
0443     const math::XYZVector &trackMom = track.momentum();
0444     double trackMag = std::sqrt(trackMom.Mag2());
0445 
0446     vars.insert(btau::trackSip3dVal, data.ip3d.value(), true);
0447     vars.insert(btau::trackSip3dSig, data.ip3d.significance(), true);
0448     vars.insert(btau::trackSip2dVal, data.ip2d.value(), true);
0449     vars.insert(btau::trackSip2dSig, data.ip2d.significance(), true);
0450     vars.insert(btau::trackJetDistVal, data.distanceToJetAxis.value(), true);
0451     vars.insert(btau::trackGhostTrackDistVal, data.distanceToGhostTrack.value(), true);
0452     vars.insert(btau::trackGhostTrackDistSig, data.distanceToGhostTrack.significance(), true);
0453     vars.insert(btau::trackGhostTrackWeight, data.ghostTrackWeight, true);
0454     vars.insert(btau::trackDecayLenVal, havePv ? (data.closestToGhostTrack - pv).mag() : -1.0, true);
0455 
0456     vars.insert(btau::trackMomentum, trackMag, true);
0457     vars.insert(btau::trackEta, trackMom.Eta(), true);
0458 
0459     vars.insert(btau::trackChi2, track.normalizedChi2(), true);
0460     vars.insert(btau::trackNPixelHits, track.hitPattern().pixelLayersWithMeasurement(), true);
0461     vars.insert(btau::trackNTotalHits, track.hitPattern().trackerLayersWithMeasurement(), true);
0462     vars.insert(btau::trackPtRel, VectorUtil::Perp(trackMom, jetDir), true);
0463     vars.insert(btau::trackDeltaR, VectorUtil::DeltaR(trackMom, jetDir), true);
0464   }
0465 
0466   vars.insert(btau::vertexCategory, vtxType, true);
0467   vars.insert(btau::trackSumJetDeltaR, VectorUtil::DeltaR(allKinematics.vectorSum(), jetDir), true);
0468   vars.insert(btau::trackSumJetEtRatio, allKinematics.vectorSum().Et() / ipInfo.jet()->et(), true);
0469   vars.insert(
0470       btau::trackSip3dSigAboveCharm, threshTrack(ipInfo, reco::btag::IP3DSig, *jet, pv).ip3d.significance(), true);
0471   vars.insert(
0472       btau::trackSip2dSigAboveCharm, threshTrack(ipInfo, reco::btag::IP2DSig, *jet, pv).ip2d.significance(), true);
0473 
0474   if (vtxType != btag::Vertices::NoVertex) {
0475     const math::XYZTLorentzVector &allSum = allKinematics.vectorSum();
0476     const math::XYZTLorentzVector &vertexSum = vertexKinematics.vectorSum();
0477 
0478     vars.insert(btau::vertexMass, vertexSum.M(), true);
0479     if (allKinematics.numberOfTracks())
0480       vars.insert(btau::vertexEnergyRatio, vertexSum.E() / allSum.E(), true);
0481     else
0482       vars.insert(btau::vertexEnergyRatio, 1, true);
0483   }
0484 
0485   vars.finalize();
0486 
0487   return vars;
0488 }