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::helper::VertexingHelper Class Reference

Produces and/or checks pat::VertexAssociation's. More...

#include "PhysicsTools/PatAlgos/interface/VertexingHelper.h"

Public Member Functions

bool enabled () const
 returns true if this was given a non dummy configuration More...
 
void newEvent (const edm::Event &event)
 To be called for each new event, reads in the vertex collection. More...
 
void newEvent (const edm::Event &event, const edm::EventSetup &setup)
 
template<typename AnyCandRef >
pat::VertexAssociation operator() (const AnyCandRef &) const
 
 VertexingHelper ()
 
 VertexingHelper (const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
 

Private Member Functions

pat::VertexAssociation associate (const reco::Candidate &) const
 
reco::TrackBaseRef getTrack_ (const reco::Candidate &c) const
 Get out the track from the Candidate / RecoCandidate / PFCandidate. More...
 

Private Attributes

pat::VertexAssociationSelector assoSelector_
 selector of associations More...
 
bool enabled_
 true if it has non null configuration More...
 
bool playback_
 true if it's just reading the associations from the event More...
 
edm::ESHandle
< TransientTrackBuilder
ttBuilder_
 
bool useTracks_
 use tracks inside candidates More...
 
edm::EDGetTokenT
< edm::ValueMap
< pat::VertexAssociation > > 
vertexAssociationsToken_
 
edm::Handle< edm::ValueMap
< pat::VertexAssociation > > 
vertexAssoMap_
 
edm::Handle
< reco::VertexCollection
vertexHandle_
 
edm::EDGetTokenT
< reco::VertexCollection
verticesToken_
 

Detailed Description

Produces and/or checks pat::VertexAssociation's.

The VertexingHelper produces pat::VertexAssociation, or reads them from the event, and can use them to select if a candidate is good or not.

Author
Giovanni Petrucciani
Version
Id:
VertexingHelper.h,v 1.3 2008/06/06 14:13:40 gpetrucc Exp

Definition at line 46 of file VertexingHelper.h.

Constructor & Destructor Documentation

pat::helper::VertexingHelper::VertexingHelper ( )
inline

Definition at line 48 of file VertexingHelper.h.

48 : enabled_(false) {}
bool enabled_
true if it has non null configuration
pat::helper::VertexingHelper::VertexingHelper ( const edm::ParameterSet iConfig,
edm::ConsumesCollector &&  iC 
)

Definition at line 11 of file VertexingHelper.cc.

References assoSelector_, edm::ParameterSet::empty(), enabled_, edm::hlt::Exception, edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), playback_, useTracks_, vertexAssociationsToken_, and verticesToken_.

12 {
13  if (!iConfig.empty()) {
14  enabled_ = true;
15  if ( iConfig.existsAs<edm::InputTag>("vertexAssociations") == iConfig.existsAs<edm::InputTag>("vertices")) {
16  throw cms::Exception("Configuration") <<
17  "VertexingHelper: you must configure either 'vertices' (to produce associations) or 'vertexAssociations' (to read them from disk), " <<
18  "you can't specify both, nor you can specify none!\n";
19  }
20 
21  if (iConfig.existsAs<edm::InputTag>("vertexAssociations")) {
22  playback_ = true;
24  }
25  if (iConfig.existsAs<edm::InputTag>("vertices")) { // vertex have been specified, so run on the fly
26  playback_ = false;
28  // ------ MODE ------------------
29  useTracks_ = iConfig.getParameter<bool>("useTracks");
30  // ------ CUTS (fully optional) ------------------
31  }
32  assoSelector_ = reco::modules::make<pat::VertexAssociationSelector>(iConfig);
33  } else {
34  enabled_ = false;
35  }
36 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
bool empty() const
Definition: ParameterSet.h:216
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:184
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool useTracks_
use tracks inside candidates
edm::EDGetTokenT< reco::VertexCollection > verticesToken_
edm::EDGetTokenT< edm::ValueMap< pat::VertexAssociation > > vertexAssociationsToken_
bool playback_
true if it&#39;s just reading the associations from the event
bool enabled_
true if it has non null configuration
pat::VertexAssociationSelector assoSelector_
selector of associations

Member Function Documentation

pat::VertexAssociation pat::helper::VertexingHelper::associate ( const reco::Candidate c) const
private

Try to associated this candidate to a vertex. If no association is found passing all cuts, return a null association

Definition at line 55 of file VertexingHelper.cc.

References end, edm::hlt::Exception, edm::RefToBase< T >::isNull(), pat::VertexAssociation::setDistances(), reco::TransientTrack::trajectoryStateClosestToPoint(), groupFilesInBlocks::tt, and reco::Candidate::vertex().

55  {
56  if (playback_) throw cms::Exception("Configuration") << "VertexingHelper: if this module was configured to read associations from the event," <<
57  " you must use 'operator()' passing a candidate ref, and not 'associate()' directly!\n";
58 
59  reco::VertexCollection::const_iterator vtx, end;
60  size_t ivtx;
63  if (useTracks_) {
64  if (!ttBuilder_.isValid()) throw cms::Exception("Configuration") << "VertexingHelper: If you use 'useTracks', you must call newEvent(iEvent,iSetup)!\n";
65  tk = getTrack_(c);
66  if (tk.isNull()) return pat::VertexAssociation();
67  tt = ttBuilder_->build(*tk);
68  }
69  for (vtx = vertexHandle_->begin(), end = vertexHandle_->end(), ivtx = 0; vtx != end; ++vtx, ++ivtx) {
71  if (useTracks_ == false) {
72  association.setDistances(c.vertex(), vtx->position(), vtx->error());
73  } else {
74  GlobalPoint vtxGP(vtx->x(), vtx->y(), vtx->z()); // need to convert XYZPoint to GlobalPoint
76  GlobalPoint trackPos = tscp.theState().position();
77  AlgebraicSymMatrix33 trackErr = tscp.theState().cartesianError().matrix().Sub<AlgebraicSymMatrix33>(0,0);
78  association.setDistances(trackPos, vtx->position(), trackErr + vtx->error());
79  }
80  if (assoSelector_(association)) return association;
81  }
82  return pat::VertexAssociation();
83 }
reco::TrackBaseRef getTrack_(const reco::Candidate &c) const
Get out the track from the Candidate / RecoCandidate / PFCandidate.
bool isNull() const
Checks for null.
Definition: RefToBase.h:270
edm::ESHandle< TransientTrackBuilder > ttBuilder_
bool useTracks_
use tracks inside candidates
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
virtual const Point & vertex() const =0
vertex position
#define end
Definition: vmac.h:37
bool playback_
true if it&#39;s just reading the associations from the event
edm::Handle< reco::VertexCollection > vertexHandle_
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
Analysis-level structure for vertex-related information.
Definition: Vertexing.h:26
bool isValid() const
Definition: ESHandle.h:37
pat::VertexAssociationSelector assoSelector_
selector of associations
bool pat::helper::VertexingHelper::enabled ( ) const
inline

returns true if this was given a non dummy configuration

Definition at line 52 of file VertexingHelper.h.

References enabled_.

Referenced by pat::PATGenericParticleProducer::produce().

52 { return enabled_; }
bool enabled_
true if it has non null configuration
reco::TrackBaseRef pat::helper::VertexingHelper::getTrack_ ( const reco::Candidate c) const
private

Get out the track from the Candidate / RecoCandidate / PFCandidate.

Definition at line 85 of file VertexingHelper.cc.

References reco::RecoCandidate::bestTrackRef(), trackerHits::c, and reco::PFCandidate::trackRef().

85  {
86  const reco::RecoCandidate *rc = dynamic_cast<const reco::RecoCandidate *>(&c);
87  if (rc != 0) { return rc->bestTrackRef(); }
88  const reco::PFCandidate *pfc = dynamic_cast<const reco::PFCandidate *>(&c);
89  if (pfc != 0) { return reco::TrackBaseRef(pfc->trackRef()); }
90 
91  return reco::TrackBaseRef();
92 }
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:429
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:22
virtual TrackBaseRef bestTrackRef() const
best track RefToBase
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:38
void pat::helper::VertexingHelper::newEvent ( const edm::Event event)

To be called for each new event, reads in the vertex collection.

Definition at line 39 of file VertexingHelper.cc.

References edm::Event::getByToken().

Referenced by pat::PATVertexAssociationProducer::produce(), and pat::PATGenericParticleProducer::produce().

39  {
40  if (playback_) {
42  } else {
43  iEvent.getByToken(verticesToken_, vertexHandle_);
44  }
45 }
edm::Handle< edm::ValueMap< pat::VertexAssociation > > vertexAssoMap_
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::VertexCollection > verticesToken_
edm::EDGetTokenT< edm::ValueMap< pat::VertexAssociation > > vertexAssociationsToken_
bool playback_
true if it&#39;s just reading the associations from the event
edm::Handle< reco::VertexCollection > vertexHandle_
void pat::helper::VertexingHelper::newEvent ( const edm::Event event,
const edm::EventSetup setup 
)

To be called for each new event, reads in the vertex collection and the tracking info You need this if 'useTrack' is true

Definition at line 48 of file VertexingHelper.cc.

References edm::EventSetup::get().

48  {
50  if (!playback_) iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", ttBuilder_);
51 }
edm::ESHandle< TransientTrackBuilder > ttBuilder_
int iEvent
Definition: GenABIO.cc:230
bool playback_
true if it&#39;s just reading the associations from the event
void newEvent(const edm::Event &event)
To be called for each new event, reads in the vertex collection.
template<typename AnyCandRef >
pat::VertexAssociation pat::helper::VertexingHelper::operator() ( const AnyCandRef &  cand) const

Return true if this candidate is associated to a valid vertex AnyCandRef should be a Ref<>, RefToBase<> or Ptr to a Candidate object

Definition at line 98 of file VertexingHelper.h.

99  {
100  if (playback_) {
101  const pat::VertexAssociation &assoc = (*vertexAssoMap_)[cand];
102  return assoSelector_(assoc) ? assoc : pat::VertexAssociation();
103  } else {
104  return associate( *cand );
105  }
106 
107  }
pat::VertexAssociation associate(const reco::Candidate &) const
bool playback_
true if it&#39;s just reading the associations from the event
Analysis-level structure for vertex-related information.
Definition: Vertexing.h:26
pat::VertexAssociationSelector assoSelector_
selector of associations

Member Data Documentation

pat::VertexAssociationSelector pat::helper::VertexingHelper::assoSelector_
private

selector of associations

Definition at line 74 of file VertexingHelper.h.

Referenced by VertexingHelper().

bool pat::helper::VertexingHelper::enabled_
private

true if it has non null configuration

Definition at line 68 of file VertexingHelper.h.

Referenced by enabled(), and VertexingHelper().

bool pat::helper::VertexingHelper::playback_
private

true if it's just reading the associations from the event

Definition at line 71 of file VertexingHelper.h.

Referenced by VertexingHelper().

edm::ESHandle<TransientTrackBuilder> pat::helper::VertexingHelper::ttBuilder_
private

Definition at line 81 of file VertexingHelper.h.

bool pat::helper::VertexingHelper::useTracks_
private

use tracks inside candidates

Definition at line 80 of file VertexingHelper.h.

Referenced by VertexingHelper().

edm::EDGetTokenT<edm::ValueMap<pat::VertexAssociation> > pat::helper::VertexingHelper::vertexAssociationsToken_
private

Definition at line 84 of file VertexingHelper.h.

Referenced by VertexingHelper().

edm::Handle<edm::ValueMap<pat::VertexAssociation> > pat::helper::VertexingHelper::vertexAssoMap_
private

Definition at line 85 of file VertexingHelper.h.

edm::Handle<reco::VertexCollection> pat::helper::VertexingHelper::vertexHandle_
private

Definition at line 78 of file VertexingHelper.h.

edm::EDGetTokenT<reco::VertexCollection> pat::helper::VertexingHelper::verticesToken_
private

Definition at line 77 of file VertexingHelper.h.

Referenced by VertexingHelper().