00001 #include <iostream>
00002 #include <cstddef>
00003 #include <string>
00004 #include <cmath>
00005 #include <vector>
00006
00007 #include <Math/VectorUtil.h>
00008
00009 #include "FWCore/Utilities/interface/Exception.h"
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011
00012 #include "DataFormats/Math/interface/Vector3D.h"
00013 #include "DataFormats/Math/interface/LorentzVector.h"
00014 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
00015 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00016 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00017 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
00018 #include "DataFormats/TrackReco/interface/Track.h"
00019 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00020 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
00021 #include "DataFormats/BTauReco/interface/SecondaryVertexTagInfo.h"
00022 #include "DataFormats/BTauReco/interface/TaggingVariable.h"
00023 #include "DataFormats/BTauReco/interface/VertexTypes.h"
00024
00025 #include "RecoVertex/VertexPrimitives/interface/ConvertToFromReco.h"
00026
00027 #include "RecoBTag/SecondaryVertex/interface/ParticleMasses.h"
00028 #include "RecoBTag/SecondaryVertex/interface/TrackSorting.h"
00029 #include "RecoBTag/SecondaryVertex/interface/TrackSelector.h"
00030 #include "RecoBTag/SecondaryVertex/interface/TrackKinematics.h"
00031 #include "RecoBTag/SecondaryVertex/interface/V0Filter.h"
00032
00033 #include "RecoBTag/SecondaryVertex/interface/CombinedSVComputer.h"
00034
00035 using namespace reco;
00036
00037 struct CombinedSVComputer::IterationRange {
00038 int begin, end, increment;
00039 };
00040
00041 #define range_for(i, x) \
00042 for(int i = (x).begin; i != (x).end; i += (x).increment)
00043
00044 static edm::ParameterSet dropDeltaR(const edm::ParameterSet &pset)
00045 {
00046 edm::ParameterSet psetCopy(pset);
00047 psetCopy.addParameter<double>("jetDeltaRMax", 99999.0);
00048 return psetCopy;
00049 }
00050
00051 CombinedSVComputer::CombinedSVComputer(const edm::ParameterSet ¶ms) :
00052 trackFlip(params.getParameter<bool>("trackFlip")),
00053 vertexFlip(params.getParameter<bool>("vertexFlip")),
00054 charmCut(params.getParameter<double>("charmCut")),
00055 sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))),
00056 trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
00057 trackNoDeltaRSelector(dropDeltaR(params.getParameter<edm::ParameterSet>("trackSelection"))),
00058 trackPseudoSelector(params.getParameter<edm::ParameterSet>("trackPseudoSelection")),
00059 pseudoMultiplicityMin(params.getParameter<unsigned int>("pseudoMultiplicityMin")),
00060 trackMultiplicityMin(params.getParameter<unsigned int>("trackMultiplicityMin")),
00061 minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
00062 useTrackWeights(params.getParameter<bool>("useTrackWeights")),
00063 vertexMassCorrection(params.getParameter<bool>("correctVertexMass")),
00064 pseudoVertexV0Filter(params.getParameter<edm::ParameterSet>("pseudoVertexV0Filter")),
00065 trackPairV0Filter(params.getParameter<edm::ParameterSet>("trackPairV0Filter"))
00066 {
00067 }
00068
00069 inline double CombinedSVComputer::flipValue(double value, bool vertex) const
00070 {
00071 return (vertex ? vertexFlip : trackFlip) ? -value : value;
00072 }
00073
00074 inline CombinedSVComputer::IterationRange CombinedSVComputer::flipIterate(
00075 int size, bool vertex) const
00076 {
00077 IterationRange range;
00078 if (vertex ? vertexFlip : trackFlip) {
00079 range.begin = size - 1;
00080 range.end = -1;
00081 range.increment = -1;
00082 } else {
00083 range.begin = 0;
00084 range.end = size;
00085 range.increment = +1;
00086 }
00087
00088 return range;
00089 }
00090
00091 const TrackIPTagInfo::TrackIPData &
00092 CombinedSVComputer::threshTrack(const TrackIPTagInfo &trackIPTagInfo,
00093 const TrackIPTagInfo::SortCriteria sort,
00094 const reco::Jet &jet,
00095 const GlobalPoint &pv) const
00096 {
00097 const edm::RefVector<TrackCollection> &tracks =
00098 trackIPTagInfo.selectedTracks();
00099 const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
00100 trackIPTagInfo.impactParameterData();
00101 std::vector<std::size_t> indices = trackIPTagInfo.sortedIndexes(sort);
00102
00103 IterationRange range = flipIterate(indices.size(), false);
00104 TrackKinematics kin;
00105 range_for(i, range) {
00106 std::size_t idx = indices[i];
00107 const TrackIPTagInfo::TrackIPData &data = ipData[idx];
00108 const Track &track = *tracks[idx];
00109
00110 if (!trackNoDeltaRSelector(track, data, jet, pv))
00111 continue;
00112
00113 kin.add(track);
00114 if (kin.vectorSum().M() > charmCut)
00115 return data;
00116 }
00117
00118 static const TrackIPTagInfo::TrackIPData dummy = {
00119 Measurement1D(-1.0, 1.0),
00120 Measurement1D(-1.0, 1.0),
00121 };
00122 return dummy;
00123 }
00124
00125 static double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
00126 {
00127 double momPar = dir.Dot(track);
00128 double energy = std::sqrt(track.Mag2() +
00129 ROOT::Math::Square(ParticleMasses::piPlus));
00130
00131 return 0.5 * std::log((energy + momPar) / (energy - momPar));
00132 }
00133
00134 TaggingVariableList
00135 CombinedSVComputer::operator () (const TrackIPTagInfo &ipInfo,
00136 const SecondaryVertexTagInfo &svInfo) const
00137 {
00138 using namespace ROOT::Math;
00139
00140 edm::RefToBase<Jet> jet = ipInfo.jet();
00141 math::XYZVector jetDir = jet->momentum().Unit();
00142 bool havePv = ipInfo.primaryVertex().isNonnull();
00143 GlobalPoint pv;
00144 if (havePv)
00145 pv = GlobalPoint(ipInfo.primaryVertex()->x(),
00146 ipInfo.primaryVertex()->y(),
00147 ipInfo.primaryVertex()->z());
00148
00149 btag::Vertices::VertexType vtxType = btag::Vertices::NoVertex;
00150
00151 TaggingVariableList vars;
00152
00153 vars.insert(btau::jetPt, jet->pt(), true);
00154 vars.insert(btau::jetEta, jet->eta(), true);
00155
00156 if (ipInfo.tracks().size() < trackMultiplicityMin)
00157 return vars;
00158
00159 TrackKinematics allKinematics;
00160 TrackKinematics vertexKinematics;
00161
00162 int vtx = -1;
00163 IterationRange range = flipIterate(svInfo.nVertices(), true);
00164 range_for(i, range) {
00165 if (vtx < 0)
00166 vtx = i;
00167
00168 TrackRefVector tracks = svInfo.vertexTracks(i);
00169 for(TrackRefVector::const_iterator track = tracks.begin();
00170 track != tracks.end(); track++) {
00171 double w = svInfo.trackWeight(i, *track);
00172 if (w >= minTrackWeight) {
00173 vertexKinematics.add(**track, w);
00174 vars.insert(btau::trackEtaRel, etaRel(jetDir,
00175 (*track)->momentum()), true);
00176 }
00177 }
00178 }
00179
00180 if (vtx >= 0) {
00181 vtxType = btag::Vertices::RecoVertex;
00182
00183 vars.insert(btau::flightDistance2dVal,
00184 flipValue(
00185 svInfo.flightDistance(vtx, true).value(),
00186 true),
00187 true);
00188 vars.insert(btau::flightDistance2dSig,
00189 flipValue(
00190 svInfo.flightDistance(vtx, true).significance(),
00191 true),
00192 true);
00193 vars.insert(btau::flightDistance3dVal,
00194 flipValue(
00195 svInfo.flightDistance(vtx, false).value(),
00196 true),
00197 true);
00198 vars.insert(btau::flightDistance3dSig,
00199 flipValue(
00200 svInfo.flightDistance(vtx, false).significance(),
00201 true),
00202 true);
00203 vars.insert(btau::vertexJetDeltaR,
00204 Geom::deltaR(svInfo.flightDirection(vtx), jetDir),true);
00205 vars.insert(btau::jetNSecondaryVertices, svInfo.nVertices(), true);
00206 vars.insert(btau::vertexNTracks, svInfo.nVertexTracks(), true);
00207 }
00208
00209 std::vector<std::size_t> indices = ipInfo.sortedIndexes(sortCriterium);
00210 const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
00211 ipInfo.impactParameterData();
00212 const edm::RefVector<TrackCollection> &tracks =
00213 ipInfo.selectedTracks();
00214 std::vector<TrackRef> pseudoVertexTracks;
00215
00216 TrackRef trackPairV0Test[2];
00217 range = flipIterate(indices.size(), false);
00218 range_for(i, range) {
00219 std::size_t idx = indices[i];
00220 const TrackIPTagInfo::TrackIPData &data = ipData[idx];
00221 const TrackRef &trackRef = tracks[idx];
00222 const Track &track = *trackRef;
00223
00224
00225
00226 if (!trackSelector(track, data, *jet, pv))
00227 continue;
00228
00229
00230
00231 allKinematics.add(track);
00232
00233
00234
00235 if (vtxType == btag::Vertices::NoVertex &&
00236 trackPseudoSelector(track, data, *jet, pv)) {
00237 pseudoVertexTracks.push_back(trackRef);
00238 vertexKinematics.add(track);
00239 }
00240
00241
00242
00243 trackPairV0Test[0] = tracks[idx];
00244 bool ok = true;
00245 range_for(j, range) {
00246 if (i == j)
00247 continue;
00248
00249 std::size_t pairIdx = indices[j];
00250 const TrackIPTagInfo::TrackIPData &pairTrackData =
00251 ipData[pairIdx];
00252 const TrackRef &pairTrackRef = tracks[pairIdx];
00253 const Track &pairTrack = *pairTrackRef;
00254
00255 if (!trackSelector(pairTrack, pairTrackData, *jet, pv))
00256 continue;
00257
00258 trackPairV0Test[1] = pairTrackRef;
00259 if (!trackPairV0Filter(trackPairV0Test, 2)) {
00260 ok = false;
00261 break;
00262 }
00263 }
00264 if (!ok)
00265 continue;
00266
00267
00268
00269 math::XYZVector trackMom = track.momentum();
00270 double trackMag = std::sqrt(trackMom.Mag2());
00271
00272 vars.insert(btau::trackSip3dVal,
00273 flipValue(data.ip3d.value(), false), true);
00274 vars.insert(btau::trackSip3dSig,
00275 flipValue(data.ip3d.significance(), false), true);
00276 vars.insert(btau::trackSip2dVal,
00277 flipValue(data.ip2d.value(), false), true);
00278 vars.insert(btau::trackSip2dSig,
00279 flipValue(data.ip2d.significance(), false), true);
00280 vars.insert(btau::trackJetDist, data.distanceToJetAxis, true);
00281 vars.insert(btau::trackFirstTrackDist, data.distanceToFirstTrack, true);
00282 vars.insert(btau::trackDecayLenVal, havePv ? (data.closestToJetAxis - pv).mag() : -1.0, true);
00283
00284 vars.insert(btau::trackMomentum, trackMag, true);
00285 vars.insert(btau::trackEta, trackMom.Eta(), true);
00286
00287 vars.insert(btau::trackPtRel, VectorUtil::Perp(trackMom, jetDir), true);
00288 vars.insert(btau::trackPPar, jetDir.Dot(trackMom), true);
00289 vars.insert(btau::trackDeltaR, VectorUtil::DeltaR(trackMom, jetDir), true);
00290 vars.insert(btau::trackPtRatio, VectorUtil::Perp(trackMom, jetDir) / trackMag, true);
00291 vars.insert(btau::trackPParRatio, jetDir.Dot(trackMom) / trackMag, true);
00292 }
00293
00294 if (vtxType == btag::Vertices::NoVertex &&
00295 vertexKinematics.numberOfTracks() >= pseudoMultiplicityMin &&
00296 pseudoVertexV0Filter(pseudoVertexTracks)) {
00297 vtxType = btag::Vertices::PseudoVertex;
00298 for(std::vector<TrackRef>::const_iterator track =
00299 pseudoVertexTracks.begin();
00300 track != pseudoVertexTracks.end(); ++track)
00301 vars.insert(btau::trackEtaRel, etaRel(jetDir,
00302 (*track)->momentum()), true);
00303 }
00304
00305 vars.insert(btau::vertexCategory, vtxType, true);
00306 vars.insert(btau::trackSumJetDeltaR,
00307 VectorUtil::DeltaR(allKinematics.vectorSum(), jetDir), true);
00308 vars.insert(btau::trackSumJetEtRatio,
00309 allKinematics.vectorSum().Et() / ipInfo.jet()->et(), true);
00310 vars.insert(btau::trackSip3dSigAboveCharm,
00311 flipValue(
00312 threshTrack(ipInfo, TrackIPTagInfo::IP3DSig, *jet, pv)
00313 .ip3d.significance(),
00314 false),
00315 true);
00316 vars.insert(btau::trackSip2dSigAboveCharm,
00317 flipValue(
00318 threshTrack(ipInfo, TrackIPTagInfo::IP2DSig, *jet, pv)
00319 .ip2d.significance(),
00320 false),
00321 true);
00322
00323 if (vtxType != btag::Vertices::NoVertex) {
00324 math::XYZTLorentzVector allSum = useTrackWeights
00325 ? allKinematics.weightedVectorSum()
00326 : allKinematics.vectorSum();
00327 math::XYZTLorentzVector vertexSum = useTrackWeights
00328 ? vertexKinematics.weightedVectorSum()
00329 : vertexKinematics.vectorSum();
00330
00331 if (vtxType != btag::Vertices::RecoVertex) {
00332 vars.insert(btau::vertexNTracks,
00333 vertexKinematics.numberOfTracks(), true);
00334 vars.insert(btau::vertexJetDeltaR,
00335 VectorUtil::DeltaR(vertexSum, jetDir), true);
00336 }
00337
00338 double vertexMass = vertexSum.M();
00339 if (vtxType == btag::Vertices::RecoVertex &&
00340 vertexMassCorrection) {
00341 GlobalVector dir = svInfo.flightDirection(vtx);
00342 double vertexPt2 =
00343 math::XYZVector(dir.x(), dir.y(), dir.z()).
00344 Cross(vertexSum).Mag2() / dir.mag2();
00345 vertexMass = std::sqrt(vertexMass * vertexMass +
00346 vertexPt2) + std::sqrt(vertexPt2);
00347 }
00348 vars.insert(btau::vertexMass, vertexMass, true);
00349 vars.insert(btau::vertexEnergyRatio, vertexSum.E() / allSum.E(), true);
00350 }
00351
00352 vars.finalize();
00353
00354 return vars;
00355 }