CMS 3D CMS Logo

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< TransientTrackBuilderttBuilder_
 
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordttToken_
 
bool useTracks_
 use tracks inside candidates More...
 
edm::EDGetTokenT< edm::ValueMap< pat::VertexAssociation > > vertexAssociationsToken_
 
edm::Handle< edm::ValueMap< pat::VertexAssociation > > vertexAssoMap_
 
edm::Handle< reco::VertexCollectionvertexHandle_
 
edm::EDGetTokenT< reco::VertexCollectionverticesToken_
 

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 51 of file VertexingHelper.h.

Constructor & Destructor Documentation

◆ VertexingHelper() [1/2]

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

Definition at line 53 of file VertexingHelper.h.

53 : enabled_(false) {}
bool enabled_
true if it has non null configuration

◆ VertexingHelper() [2/2]

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_, Exception, edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), playback_, ttToken_, useTracks_, vertexAssociationsToken_, and verticesToken_.

11  {
12  if (!iConfig.empty()) {
13  enabled_ = true;
14  if (iConfig.existsAs<edm::InputTag>("vertexAssociations") == iConfig.existsAs<edm::InputTag>("vertices")) {
15  throw cms::Exception("Configuration") << "VertexingHelper: you must configure either 'vertices' (to produce "
16  "associations) or 'vertexAssociations' (to read them from disk), "
17  << "you can't specify both, nor you can specify none!\n";
18  }
19 
20  if (iConfig.existsAs<edm::InputTag>("vertexAssociations")) {
21  playback_ = true;
23  iConfig.getParameter<edm::InputTag>("vertexAssociations"));
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  if (!playback_) {
37  ttToken_ = iC.esConsumes(edm::ESInputTag("", "TransientTrackBuilder"));
38  }
39 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > ttToken_
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
bool useTracks_
use tracks inside candidates
bool empty() const
Definition: ParameterSet.h:201
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

◆ associate()

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 HltBtagPostValidation_cff::c, Exception, edm::RefToBase< T >::isNull(), groupFilesInBlocks::tt, and L1BJetProducer_cff::vtx.

55  {
56  if (playback_)
57  throw cms::Exception("Configuration")
58  << "VertexingHelper: if this module was configured to read associations from the event,"
59  << " you must use 'operator()' passing a candidate ref, and not 'associate()' directly!\n";
60 
61  reco::VertexCollection::const_iterator vtx, end;
62  size_t ivtx;
65  if (useTracks_) {
66  if (!ttBuilder_.isValid())
67  throw cms::Exception("Configuration")
68  << "VertexingHelper: If you use 'useTracks', you must call newEvent(iEvent,iSetup)!\n";
69  tk = getTrack_(c);
70  if (tk.isNull())
71  return pat::VertexAssociation();
72  tt = ttBuilder_->build(*tk);
73  }
74  for (vtx = vertexHandle_->begin(), end = vertexHandle_->end(), ivtx = 0; vtx != end; ++vtx, ++ivtx) {
76  if (useTracks_ == false) {
77  association.setDistances(c.vertex(), vtx->position(), vtx->error());
78  } else {
79  GlobalPoint vtxGP(vtx->x(), vtx->y(), vtx->z()); // need to convert XYZPoint to GlobalPoint
80  TrajectoryStateClosestToPoint tscp = tt.trajectoryStateClosestToPoint(vtxGP);
81  GlobalPoint trackPos = tscp.theState().position();
82  AlgebraicSymMatrix33 trackErr = tscp.theState().cartesianError().matrix().Sub<AlgebraicSymMatrix33>(0, 0);
83  association.setDistances(trackPos, vtx->position(), trackErr + vtx->error());
84  }
86  return association;
87  }
88  return pat::VertexAssociation();
89 }
edm::ESHandle< TransientTrackBuilder > ttBuilder_
reco::TransientTrack build(const reco::Track *p) const
bool useTracks_
use tracks inside candidates
Definition: TTTypes.h:54
std::tuple< layerClusterToCaloParticle, caloParticleToLayerCluster > association
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:297
bool playback_
true if it&#39;s just reading the associations from the event
bool isValid() const
Definition: ESHandle.h:44
edm::Handle< reco::VertexCollection > vertexHandle_
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
Analysis-level structure for vertex-related information.
Definition: Vertexing.h:25
pat::VertexAssociationSelector assoSelector_
selector of associations

◆ enabled()

bool pat::helper::VertexingHelper::enabled ( ) const
inline

returns true if this was given a non dummy configuration

Definition at line 57 of file VertexingHelper.h.

References enabled_.

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

57 { return enabled_; }
bool enabled_
true if it has non null configuration

◆ getTrack_()

reco::TrackBaseRef pat::helper::VertexingHelper::getTrack_ ( const reco::Candidate c) const
private

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

Definition at line 91 of file VertexingHelper.cc.

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

91  {
92  const reco::RecoCandidate *rc = dynamic_cast<const reco::RecoCandidate *>(&c);
93  if (rc != nullptr) {
94  return rc->bestTrackRef();
95  }
96  const reco::PFCandidate *pfc = dynamic_cast<const reco::PFCandidate *>(&c);
97  if (pfc != nullptr) {
98  return reco::TrackBaseRef(pfc->trackRef());
99  }
100 
101  return reco::TrackBaseRef();
102 }
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:35
virtual TrackBaseRef bestTrackRef() const
best track RefToBase
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:437

◆ newEvent() [1/2]

void pat::helper::VertexingHelper::newEvent ( const edm::Event event)

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

Definition at line 41 of file VertexingHelper.cc.

References iEvent.

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

41  {
42  if (playback_) {
44  } else {
45  iEvent.getByToken(verticesToken_, vertexHandle_);
46  }
47 }
edm::Handle< edm::ValueMap< pat::VertexAssociation > > vertexAssoMap_
int iEvent
Definition: GenABIO.cc:224
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_

◆ newEvent() [2/2]

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 49 of file VertexingHelper.cc.

References edm::EventSetup::getHandle(), and iEvent.

49  {
51  if (!playback_)
52  ttBuilder_ = iSetup.getHandle(ttToken_);
53 }
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > ttToken_
edm::ESHandle< TransientTrackBuilder > ttBuilder_
int iEvent
Definition: GenABIO.cc:224
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.

◆ operator()()

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 103 of file VertexingHelper.h.

References trackingPlots::assoc.

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

Member Data Documentation

◆ assoSelector_

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

selector of associations

Definition at line 79 of file VertexingHelper.h.

Referenced by VertexingHelper().

◆ enabled_

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

true if it has non null configuration

Definition at line 73 of file VertexingHelper.h.

Referenced by enabled(), and VertexingHelper().

◆ playback_

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

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

Definition at line 76 of file VertexingHelper.h.

Referenced by VertexingHelper().

◆ ttBuilder_

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

Definition at line 87 of file VertexingHelper.h.

◆ ttToken_

edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> pat::helper::VertexingHelper::ttToken_
private

Definition at line 86 of file VertexingHelper.h.

Referenced by VertexingHelper().

◆ useTracks_

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

use tracks inside candidates

Definition at line 85 of file VertexingHelper.h.

Referenced by VertexingHelper().

◆ vertexAssociationsToken_

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

Definition at line 90 of file VertexingHelper.h.

Referenced by VertexingHelper().

◆ vertexAssoMap_

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

Definition at line 91 of file VertexingHelper.h.

◆ vertexHandle_

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

Definition at line 83 of file VertexingHelper.h.

◆ verticesToken_

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

Definition at line 82 of file VertexingHelper.h.

Referenced by VertexingHelper().