CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
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 edm::Event &iEvent)
 
float calculate (const Muon &aMuon, const edm::Event &iEvent)
 
 LeptonVertexSignificance ()
 
 LeptonVertexSignificance (const edm::EventSetup &iSetup, edm::ConsumesCollector &&iC)
 
 ~LeptonVertexSignificance ()
 

Private Member Functions

float calculate (const reco::Track &track, const edm::Event &iEvent)
 

Private Attributes

TransientTrackBuildertheTrackBuilder_
 
edm::EDGetTokenT
< reco::VertexCollection
vertexToken_
 

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 37 of file LeptonVertexSignificance.h.

Constructor & Destructor Documentation

pat::LeptonVertexSignificance::LeptonVertexSignificance ( )
LeptonVertexSignificance::LeptonVertexSignificance ( const edm::EventSetup iSetup,
edm::ConsumesCollector &&  iC 
)

Definition at line 19 of file LeptonVertexSignificance.cc.

References edm::EventSetup::get(), edm::ESHandle< class >::product(), and theTrackBuilder_.

20 : vertexToken_( iC.consumes< reco::VertexCollection >( edm::InputTag( "offlinePrimaryVerticesFromCTFTracks" ) ) )
21 {
22  // instantiate the transient-track builder
24  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", builder);
26 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
TransientTrackBuilder * theTrackBuilder_
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
LeptonVertexSignificance::~LeptonVertexSignificance ( )

Definition at line 29 of file LeptonVertexSignificance.cc.

References theTrackBuilder_.

29  {
30  delete theTrackBuilder_;
31 }
TransientTrackBuilder * theTrackBuilder_

Member Function Documentation

float LeptonVertexSignificance::calculate ( const Electron anElectron,
const edm::Event iEvent 
)

Definition at line 34 of file LeptonVertexSignificance.cc.

References pat::Electron::gsfTrack(), and iEvent.

Referenced by calculate().

34  {
35  return this->calculate(*theElectron.gsfTrack(), iEvent);
36 }
int iEvent
Definition: GenABIO.cc:243
float calculate(const Electron &anElectron, const edm::Event &iEvent)
float LeptonVertexSignificance::calculate ( const Muon aMuon,
const edm::Event iEvent 
)

Definition at line 38 of file LeptonVertexSignificance.cc.

References calculate(), iEvent, and pat::Muon::track().

38  {
39  return this->calculate(*theMuon.track(), iEvent);
40 }
int iEvent
Definition: GenABIO.cc:243
float calculate(const Electron &anElectron, const edm::Event &iEvent)
float LeptonVertexSignificance::calculate ( const reco::Track track,
const edm::Event iEvent 
)
private

Definition at line 43 of file LeptonVertexSignificance.cc.

References TransientTrackBuilder::build(), edm::Event::getByToken(), reco::Vertex::position(), funct::pow(), edm::Handle< T >::product(), mathSSE::sqrt(), TrajectoryStateClosestToPoint::theState(), theTrackBuilder_, reco::TransientTrack::trajectoryStateClosestToPoint(), vertexToken_, and reco::Vertex::zError().

43  {
44  // FIXME: think more about how to handle events without vertices
45  // lepton LR calculation should have nothing to do with event selection
47  iEvent.getByToken(vertexToken_, vertexHandle);
48  if (vertexHandle.product()->size() == 0) return 0;
49  reco::Vertex theVertex = vertexHandle.product()->front();
50  // calculate the track-vertex association significance
51  reco::TransientTrack theTrTrack = theTrackBuilder_->build(&theTrack);
52  GlobalPoint theVertexPoint(theVertex.position().x(), theVertex.position().y(), theVertex.position().z());
53  FreeTrajectoryState theLeptonNearVertex = theTrTrack.trajectoryStateClosestToPoint(theVertexPoint).theState();
54  return fabs(theVertex.position().z() - theLeptonNearVertex.position().z())
55  / sqrt(std::pow(theVertex.zError(), 2) + theLeptonNearVertex.cartesianError().position().czz());
56 }
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
double zError() const
error on z
Definition: Vertex.h:104
const FreeTrajectoryState & theState() const
reco::TransientTrack build(const reco::Track *p) const
const Point & position() const
position
Definition: Vertex.h:92
TransientTrackBuilder * theTrackBuilder_
T sqrt(T t)
Definition: SSEVec.h:48
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
T const * product() const
Definition: Handle.h:81
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40

Member Data Documentation

TransientTrackBuilder* pat::LeptonVertexSignificance::theTrackBuilder_
private
edm::EDGetTokenT<reco::VertexCollection> pat::LeptonVertexSignificance::vertexToken_
private

Definition at line 49 of file LeptonVertexSignificance.h.

Referenced by calculate().