Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "PhysicsTools/PatUtils/interface/LeptonVertexSignificance.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "FWCore/Utilities/interface/Exception.h"
0004 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0005 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0006 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
0007 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0008 #include "DataFormats/PatCandidates/interface/Electron.h"
0009 #include "DataFormats/PatCandidates/interface/Muon.h"
0010 
0011 using namespace pat;
0012 
0013 edm::InputTag LeptonVertexSignificance::vertexCollectionTag() {
0014   return edm::InputTag("offlinePrimaryVerticesFromCTFTracks");
0015 }
0016 
0017 // calculate the TrackIsoPt for the lepton object
0018 float LeptonVertexSignificance::calculate(const Electron& theElectron,
0019                                           const reco::VertexCollection& vertices,
0020                                           const TransientTrackBuilder& builder) {
0021   return this->calculate(*theElectron.gsfTrack(), vertices, builder);
0022 }
0023 
0024 float LeptonVertexSignificance::calculate(const Muon& theMuon,
0025                                           const reco::VertexCollection& vertices,
0026                                           const TransientTrackBuilder& builder) {
0027   return this->calculate(*theMuon.track(), vertices, builder);
0028 }
0029 
0030 // calculate the TrackIsoPt for the lepton's track
0031 float LeptonVertexSignificance::calculate(const reco::Track& theTrack,
0032                                           const reco::VertexCollection& vertices,
0033                                           const TransientTrackBuilder& builder) {
0034   if (vertices.empty())
0035     return 0;
0036   reco::Vertex const& theVertex = vertices.front();
0037   // calculate the track-vertex association significance
0038   reco::TransientTrack theTrTrack = builder.build(&theTrack);
0039   GlobalPoint theVertexPoint(theVertex.position().x(), theVertex.position().y(), theVertex.position().z());
0040   FreeTrajectoryState theLeptonNearVertex = theTrTrack.trajectoryStateClosestToPoint(theVertexPoint).theState();
0041   return fabs(theVertex.position().z() - theLeptonNearVertex.position().z()) /
0042          sqrt(std::pow(theVertex.zError(), 2) + theLeptonNearVertex.cartesianError().position().czz());
0043 }