CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/RecoBTag/SecondaryVertex/src/GhostTrackComputer.cc

Go to the documentation of this file.
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 "RecoBTag/SecondaryVertex/interface/ParticleMasses.h"
00026 #include "RecoBTag/SecondaryVertex/interface/TrackSorting.h"
00027 #include "RecoBTag/SecondaryVertex/interface/TrackSelector.h"
00028 #include "RecoBTag/SecondaryVertex/interface/TrackKinematics.h"
00029 #include "RecoBTag/SecondaryVertex/interface/V0Filter.h"
00030 
00031 #include "RecoBTag/SecondaryVertex/interface/GhostTrackComputer.h"
00032 
00033 using namespace reco;
00034 
00035 static edm::ParameterSet dropDeltaR(const edm::ParameterSet &pset)
00036 {
00037         edm::ParameterSet psetCopy(pset);
00038         psetCopy.addParameter<double>("jetDeltaRMax", 99999.0);
00039         return psetCopy;
00040 }
00041 
00042 GhostTrackComputer::GhostTrackComputer(const edm::ParameterSet &params) :
00043         charmCut(params.getParameter<double>("charmCut")),
00044         sortCriterium(TrackSorting::getCriterium(params.getParameter<std::string>("trackSort"))),
00045         trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
00046         trackNoDeltaRSelector(dropDeltaR(params.getParameter<edm::ParameterSet>("trackSelection"))),
00047         minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
00048         trackPairV0Filter(params.getParameter<edm::ParameterSet>("trackPairV0Filter"))
00049 {
00050 }
00051 
00052 const TrackIPTagInfo::TrackIPData &
00053 GhostTrackComputer::threshTrack(const TrackIPTagInfo &trackIPTagInfo,
00054                                 const TrackIPTagInfo::SortCriteria sort,
00055                                 const reco::Jet &jet,
00056                                 const GlobalPoint &pv) const
00057 {
00058         const edm::RefVector<TrackCollection> &tracks =
00059                                         trackIPTagInfo.selectedTracks();
00060         const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
00061                                         trackIPTagInfo.impactParameterData();
00062         std::vector<std::size_t> indices = trackIPTagInfo.sortedIndexes(sort);
00063 
00064         TrackKinematics kin;
00065         for(std::vector<std::size_t>::const_iterator iter = indices.begin();
00066             iter != indices.end(); ++iter) {
00067                 const TrackIPTagInfo::TrackIPData &data = ipData[*iter];
00068                 const Track &track = *tracks[*iter];
00069 
00070                 if (!trackNoDeltaRSelector(track, data, jet, pv))
00071                         continue;
00072 
00073                 kin.add(track);
00074                 if (kin.vectorSum().M() > charmCut) 
00075                         return data;
00076         }
00077 
00078         static const TrackIPTagInfo::TrackIPData dummy = {
00079                 GlobalPoint(),
00080                 GlobalPoint(),
00081                 Measurement1D(-1.0, 1.0),
00082                 Measurement1D(-1.0, 1.0),
00083                 Measurement1D(-1.0, 1.0),
00084                 Measurement1D(-1.0, 1.0),
00085                 0.
00086         };
00087         return dummy;
00088 }
00089 
00090 static double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
00091 {
00092         double momPar = dir.Dot(track);
00093         double energy = std::sqrt(track.Mag2() +
00094                                   ROOT::Math::Square(ParticleMasses::piPlus));
00095 
00096         return 0.5 * std::log((energy + momPar) / (energy - momPar));
00097 }
00098 
00099 static void addMeas(std::pair<double, double> &sum, Measurement1D meas)
00100 {
00101         double weight = 1. / meas.error();
00102         weight *= weight;
00103         sum.first += weight * meas.value();
00104         sum.second += weight;
00105 }
00106 
00107 TaggingVariableList
00108 GhostTrackComputer::operator () (const TrackIPTagInfo &ipInfo,
00109                                  const SecondaryVertexTagInfo &svInfo) const
00110 {
00111         using namespace ROOT::Math;
00112 
00113         edm::RefToBase<Jet> jet = ipInfo.jet();
00114         math::XYZVector jetDir = jet->momentum().Unit();
00115         bool havePv = ipInfo.primaryVertex().isNonnull();
00116         GlobalPoint pv;
00117         if (havePv)
00118                 pv = GlobalPoint(ipInfo.primaryVertex()->x(),
00119                                  ipInfo.primaryVertex()->y(),
00120                                  ipInfo.primaryVertex()->z());
00121 
00122         btag::Vertices::VertexType vtxType = btag::Vertices::NoVertex;
00123 
00124         TaggingVariableList vars;
00125 
00126         vars.insert(btau::jetPt, jet->pt(), true);
00127         vars.insert(btau::jetEta, jet->eta(), true);
00128 
00129         TrackKinematics allKinematics;
00130         TrackKinematics vertexKinematics;
00131         TrackKinematics trackKinematics;
00132 
00133         std::pair<double, double> vertexDist2D, vertexDist3D;
00134         std::pair<double, double> tracksDist2D, tracksDist3D;
00135 
00136         unsigned int nVertices = 0;
00137         unsigned int nVertexTracks = 0;
00138         unsigned int nTracks = 0;
00139         for(unsigned int i = 0; i < svInfo.nVertices(); i++) {
00140                 const Vertex &vertex = svInfo.secondaryVertex(i);
00141                 bool hasRefittedTracks = vertex.hasRefittedTracks();
00142                 TrackRefVector tracks = svInfo.vertexTracks(i);
00143                 unsigned int n = 0;
00144                 for(TrackRefVector::const_iterator track = tracks.begin();
00145                     track != tracks.end(); track++)
00146                         if (svInfo.trackWeight(i, *track) >= minTrackWeight)
00147                                 n++;
00148 
00149                 if (n < 1)
00150                         continue;
00151                 bool isTrackVertex = (n == 1);
00152                 ++*(isTrackVertex ? &nTracks : &nVertices);
00153 
00154                 addMeas(*(isTrackVertex ? &tracksDist2D : &vertexDist2D),
00155                                         svInfo.flightDistance(i, true));
00156                 addMeas(*(isTrackVertex ? &tracksDist3D : &vertexDist3D),
00157                                         svInfo.flightDistance(i, false));
00158 
00159                 TrackKinematics &kin = isTrackVertex ? trackKinematics
00160                                                      : vertexKinematics;
00161                 for(TrackRefVector::const_iterator track = tracks.begin();
00162                     track != tracks.end(); track++) {
00163                         float w = svInfo.trackWeight(i, *track);
00164                         if (w < minTrackWeight)
00165                                 continue;
00166                         if (hasRefittedTracks) {
00167                                 Track actualTrack =
00168                                                 vertex.refittedTrack(*track);
00169                                 kin.add(actualTrack, w);
00170                                 vars.insert(btau::trackEtaRel, etaRel(jetDir,
00171                                                 actualTrack.momentum()), true);
00172                         } else {
00173                                 kin.add(**track, w);
00174                                 vars.insert(btau::trackEtaRel, etaRel(jetDir,
00175                                                 (*track)->momentum()), true);
00176                         }
00177                         if (!isTrackVertex)
00178                                 nVertexTracks++;
00179                 }
00180         }
00181 
00182         Measurement1D dist2D, dist3D;
00183         if (nVertices) {
00184                 vtxType = btag::Vertices::RecoVertex;
00185 
00186                 if (nVertices == 1 && nTracks) {
00187                         vertexDist2D.first += tracksDist2D.first;
00188                         vertexDist2D.second += tracksDist2D.second;
00189                         vertexDist3D.first += tracksDist3D.first;
00190                         vertexDist3D.second += tracksDist3D.second;
00191                         vertexKinematics += trackKinematics;
00192                 }
00193 
00194                 dist2D = Measurement1D(
00195                                 vertexDist2D.first / vertexDist2D.second,
00196                                 std::sqrt(1. / vertexDist2D.second));
00197                 dist3D = Measurement1D(
00198                                 vertexDist3D.first / vertexDist3D.second,
00199                                 std::sqrt(1. / vertexDist3D.second));
00200 
00201                 vars.insert(btau::jetNSecondaryVertices, nVertices, true);
00202                 vars.insert(btau::vertexNTracks, nVertexTracks, true);
00203         } else if (nTracks) {
00204                 vtxType = btag::Vertices::PseudoVertex;
00205                 vertexKinematics = trackKinematics;
00206 
00207                 dist2D = Measurement1D(
00208                                 tracksDist2D.first / tracksDist2D.second,
00209                                 std::sqrt(1. / tracksDist2D.second));
00210                 dist3D = Measurement1D(
00211                                 tracksDist3D.first / tracksDist3D.second,
00212                                 std::sqrt(1. / tracksDist3D.second));
00213         }
00214 
00215         if (nVertices || nTracks) {
00216                 vars.insert(btau::flightDistance2dVal, dist2D.value(), true);
00217                 vars.insert(btau::flightDistance2dSig, dist2D.significance(), true);
00218                 vars.insert(btau::flightDistance3dVal, dist3D.value(), true);
00219                 vars.insert(btau::flightDistance3dSig, dist3D.significance(), true);
00220                 vars.insert(btau::vertexJetDeltaR,
00221                             Geom::deltaR(svInfo.flightDirection(0), jetDir),true);
00222                 vars.insert(btau::jetNSingleTrackVertices, nTracks, true);
00223         }
00224 
00225         std::vector<std::size_t> indices = ipInfo.sortedIndexes(sortCriterium);
00226         const std::vector<TrackIPTagInfo::TrackIPData> &ipData =
00227                                                 ipInfo.impactParameterData();
00228         const edm::RefVector<TrackCollection> &tracks =
00229                                                 ipInfo.selectedTracks();
00230 
00231         TrackRef trackPairV0Test[2];
00232         for(unsigned int i = 0; i < indices.size(); i++) {
00233                 std::size_t idx = indices[i];
00234                 const TrackIPTagInfo::TrackIPData &data = ipData[idx];
00235                 const TrackRef &trackRef = tracks[idx];
00236                 const Track &track = *trackRef;
00237 
00238                 // filter track
00239 
00240                 if (!trackSelector(track, data, *jet, pv))
00241                         continue;
00242 
00243                 // add track to kinematics for all tracks in jet
00244 
00245                 allKinematics.add(track);
00246 
00247                 // check against all other tracks for V0 track pairs
00248 
00249                 trackPairV0Test[0] = tracks[idx];
00250                 bool ok = true;
00251                 for(unsigned int j = 0; j < indices.size(); j++) {
00252                         if (i == j)
00253                                 continue;
00254 
00255                         std::size_t pairIdx = indices[j];
00256                         const TrackIPTagInfo::TrackIPData &pairTrackData =
00257                                                         ipData[pairIdx];
00258                         const TrackRef &pairTrackRef = tracks[pairIdx];
00259                         const Track &pairTrack = *pairTrackRef;
00260 
00261                         if (!trackSelector(pairTrack, pairTrackData, *jet, pv))
00262                                 continue;
00263 
00264                         trackPairV0Test[1] = pairTrackRef;
00265                         if (!trackPairV0Filter(trackPairV0Test, 2)) {
00266                                 ok = false;
00267                                 break;
00268                         }
00269                 }
00270                 if (!ok)
00271                         continue;
00272 
00273                 // add track variables
00274 
00275                 math::XYZVector trackMom = track.momentum();
00276                 double trackMag = std::sqrt(trackMom.Mag2());
00277 
00278                 vars.insert(btau::trackSip3dVal, data.ip3d.value(), true);
00279                 vars.insert(btau::trackSip3dSig, data.ip3d.significance(), true);
00280                 vars.insert(btau::trackSip2dVal, data.ip2d.value(), true);
00281                 vars.insert(btau::trackSip2dSig, data.ip2d.significance(), true);
00282                 vars.insert(btau::trackJetDistVal, data.distanceToJetAxis.value(), true);
00283                 vars.insert(btau::trackGhostTrackDistVal, data.distanceToGhostTrack.value(), true);
00284                 vars.insert(btau::trackGhostTrackDistSig, data.distanceToGhostTrack.significance(), true);
00285                 vars.insert(btau::trackGhostTrackWeight, data.ghostTrackWeight, true);
00286                 vars.insert(btau::trackDecayLenVal, havePv ? (data.closestToGhostTrack - pv).mag() : -1.0, true);
00287 
00288                 vars.insert(btau::trackMomentum, trackMag, true);
00289                 vars.insert(btau::trackEta, trackMom.Eta(), true);
00290 
00291                 vars.insert(btau::trackChi2, track.normalizedChi2(), true);
00292                 vars.insert(btau::trackNPixelHits, track.hitPattern().pixelLayersWithMeasurement(), true);
00293                 vars.insert(btau::trackNTotalHits, track.hitPattern().trackerLayersWithMeasurement(), true);
00294                 vars.insert(btau::trackPtRel, VectorUtil::Perp(trackMom, jetDir), true);
00295                 vars.insert(btau::trackDeltaR, VectorUtil::DeltaR(trackMom, jetDir), true);
00296         } 
00297 
00298         vars.insert(btau::vertexCategory, vtxType, true);
00299         vars.insert(btau::trackSumJetDeltaR,
00300                     VectorUtil::DeltaR(allKinematics.vectorSum(), jetDir), true);
00301         vars.insert(btau::trackSumJetEtRatio,
00302                     allKinematics.vectorSum().Et() / ipInfo.jet()->et(), true);
00303         vars.insert(btau::trackSip3dSigAboveCharm,
00304                     threshTrack(ipInfo, TrackIPTagInfo::IP3DSig, *jet, pv)
00305                                                         .ip3d.significance(),
00306                     true);
00307         vars.insert(btau::trackSip2dSigAboveCharm,
00308                     threshTrack(ipInfo, TrackIPTagInfo::IP2DSig, *jet, pv)
00309                                                         .ip2d.significance(),
00310                     true);
00311 
00312         if (vtxType != btag::Vertices::NoVertex) {
00313                 math::XYZTLorentzVector allSum = allKinematics.vectorSum();
00314                 math::XYZTLorentzVector vertexSum = vertexKinematics.vectorSum();
00315 
00316                 vars.insert(btau::vertexMass, vertexSum.M(), true);
00317                 if (allKinematics.numberOfTracks())
00318                         vars.insert(btau::vertexEnergyRatio,
00319                                     vertexSum.E() / allSum.E(), true);
00320                 else
00321                         vars.insert(btau::vertexEnergyRatio, 1, true);
00322         }
00323 
00324         vars.finalize();
00325 
00326         return vars;
00327 }