CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions
pat::LeptonVertexSignificance Class Reference

Calculates a lepton's vertex association significance. More...

#include "PhysicsTools/PatUtils/interface/LeptonVertexSignificance.h"

Public Member Functions

float calculate (const Electron &anElectron, const reco::VertexCollection &vertices, const TransientTrackBuilder &builder)
 
float calculate (const Muon &aMuon, const reco::VertexCollection &vertices, const TransientTrackBuilder &builder)
 
 LeptonVertexSignificance ()=default
 
 ~LeptonVertexSignificance ()=default
 

Static Public Member Functions

static edm::InputTag vertexCollectionTag ()
 

Private Member Functions

float calculate (const reco::Track &track, const reco::VertexCollection &vertices, const TransientTrackBuilder &builder)
 

Detailed Description

Calculates a lepton's vertex association significance.

LeptonVertexSignificance calculates the significance of the association of the lepton to a given vertex, as defined in CMS Note 2006/024

Author
Steven Lowette

Definition at line 31 of file LeptonVertexSignificance.h.

Constructor & Destructor Documentation

◆ LeptonVertexSignificance()

pat::LeptonVertexSignificance::LeptonVertexSignificance ( )
default

◆ ~LeptonVertexSignificance()

pat::LeptonVertexSignificance::~LeptonVertexSignificance ( )
default

Member Function Documentation

◆ calculate() [1/3]

float LeptonVertexSignificance::calculate ( const Electron anElectron,
const reco::VertexCollection vertices,
const TransientTrackBuilder builder 
)

Definition at line 18 of file LeptonVertexSignificance.cc.

References pat::Electron::gsfTrack(), and AlignmentTracksFromVertexSelector_cfi::vertices.

Referenced by calculate().

20  {
21  return this->calculate(*theElectron.gsfTrack(), vertices, builder);
22 }
float calculate(const Electron &anElectron, const reco::VertexCollection &vertices, const TransientTrackBuilder &builder)

◆ calculate() [2/3]

float LeptonVertexSignificance::calculate ( const Muon aMuon,
const reco::VertexCollection vertices,
const TransientTrackBuilder builder 
)

Definition at line 24 of file LeptonVertexSignificance.cc.

References calculate(), pat::Muon::track(), and AlignmentTracksFromVertexSelector_cfi::vertices.

26  {
27  return this->calculate(*theMuon.track(), vertices, builder);
28 }
float calculate(const Electron &anElectron, const reco::VertexCollection &vertices, const TransientTrackBuilder &builder)

◆ calculate() [3/3]

float LeptonVertexSignificance::calculate ( const reco::Track track,
const reco::VertexCollection vertices,
const TransientTrackBuilder builder 
)
private

Definition at line 31 of file LeptonVertexSignificance.cc.

References TransientTrackBuilder::build(), reco::Vertex::position(), funct::pow(), mathSSE::sqrt(), TrajectoryStateClosestToPoint::theState(), reco::TransientTrack::trajectoryStateClosestToPoint(), AlignmentTracksFromVertexSelector_cfi::vertices, and reco::Vertex::zError().

33  {
34  if (vertices.empty())
35  return 0;
36  reco::Vertex const& theVertex = vertices.front();
37  // calculate the track-vertex association significance
38  reco::TransientTrack theTrTrack = builder.build(&theTrack);
39  GlobalPoint theVertexPoint(theVertex.position().x(), theVertex.position().y(), theVertex.position().z());
40  FreeTrajectoryState theLeptonNearVertex = theTrTrack.trajectoryStateClosestToPoint(theVertexPoint).theState();
41  return fabs(theVertex.position().z() - theLeptonNearVertex.position().z()) /
42  sqrt(std::pow(theVertex.zError(), 2) + theLeptonNearVertex.cartesianError().position().czz());
43 }
const Point & position() const
position
Definition: Vertex.h:128
double zError() const
error on z
Definition: Vertex.h:142
reco::TransientTrack build(const reco::Track *p) const
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
T sqrt(T t)
Definition: SSEVec.h:19
const FreeTrajectoryState & theState() const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ vertexCollectionTag()

edm::InputTag LeptonVertexSignificance::vertexCollectionTag ( )
static

Definition at line 13 of file LeptonVertexSignificance.cc.

References ProducerED_cfi::InputTag.

13  {
14  return edm::InputTag("offlinePrimaryVerticesFromCTFTracks");
15 }