CMS 3D CMS Logo

VertexingHelper.cc
Go to the documentation of this file.
2 #include <algorithm>
3 
4 #include <iostream>
5 
10 
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;
23  vertexAssociationsToken_ = iC.consumes<edm::ValueMap<pat::VertexAssociation> >(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;
27  verticesToken_ = iC.consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"));
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 }
37 
38 void
40  if (playback_) {
42  } else {
44  }
45 }
46 
47 void
49  newEvent(iEvent);
50  if (!playback_) iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", ttBuilder_);
51 }
52 
53 
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 }
84 
86  const reco::RecoCandidate *rc = dynamic_cast<const reco::RecoCandidate *>(&c);
87  if (rc != nullptr) { return rc->bestTrackRef(); }
88  const reco::PFCandidate *pfc = dynamic_cast<const reco::PFCandidate *>(&c);
89  if (pfc != nullptr) { return reco::TrackBaseRef(pfc->trackRef()); }
90 
91  return reco::TrackBaseRef();
92 }
void setDistances(const AlgebraicVector3 &dist, const AlgebraicSymMatrix33 &err)
Set dz and dr given the distance and the 3x3 total covariance matrix of the distance.
Definition: Vertexing.cc:5
T getParameter(std::string const &) const
bool empty() const
Definition: ParameterSet.h:218
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
edm::Handle< edm::ValueMap< pat::VertexAssociation > > vertexAssoMap_
reco::TrackBaseRef getTrack_(const reco::Candidate &c) const
Get out the track from the Candidate / RecoCandidate / PFCandidate.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
reco::TransientTrack build(const reco::Track *p) const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::ESHandle< TransientTrackBuilder > ttBuilder_
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:442
bool useTracks_
use tracks inside candidates
int iEvent
Definition: GenABIO.cc:230
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
edm::RefToBase< reco::Track > TrackBaseRef
persistent reference to a Track, using views
Definition: TrackFwd.h:32
virtual TrackBaseRef bestTrackRef() const
best track RefToBase
pat::VertexAssociation associate(const reco::Candidate &) const
edm::EDGetTokenT< reco::VertexCollection > verticesToken_
edm::EDGetTokenT< edm::ValueMap< pat::VertexAssociation > > vertexAssociationsToken_
#define end
Definition: vmac.h:39
bool playback_
true if it&#39;s just reading the associations from the event
edm::Handle< reco::VertexCollection > vertexHandle_
void newEvent(const edm::Event &event)
To be called for each new event, reads in the vertex collection.
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
bool isNull() const
Checks for null.
Definition: RefToBase.h:328
const T & get() const
Definition: EventSetup.h:59
bool enabled_
true if it has non null configuration
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
Analysis-level structure for vertex-related information.
Definition: Vertexing.h:26
virtual const Point & vertex() const =0
vertex position
bool isValid() const
Definition: ESHandle.h:47
pat::VertexAssociationSelector assoSelector_
selector of associations