File indexing completed on 2024-04-06 12:23:58
0001 #include "PhysicsTools/PatAlgos/interface/VertexingHelper.h"
0002 #include <algorithm>
0003
0004 #include <iostream>
0005
0006 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0007 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0008 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0009 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0010
0011 pat::helper::VertexingHelper::VertexingHelper(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC) {
0012 if (!iConfig.empty()) {
0013 enabled_ = true;
0014 if (iConfig.existsAs<edm::InputTag>("vertexAssociations") == iConfig.existsAs<edm::InputTag>("vertices")) {
0015 throw cms::Exception("Configuration") << "VertexingHelper: you must configure either 'vertices' (to produce "
0016 "associations) or 'vertexAssociations' (to read them from disk), "
0017 << "you can't specify both, nor you can specify none!\n";
0018 }
0019
0020 if (iConfig.existsAs<edm::InputTag>("vertexAssociations")) {
0021 playback_ = true;
0022 vertexAssociationsToken_ = iC.consumes<edm::ValueMap<pat::VertexAssociation> >(
0023 iConfig.getParameter<edm::InputTag>("vertexAssociations"));
0024 }
0025 if (iConfig.existsAs<edm::InputTag>("vertices")) {
0026 playback_ = false;
0027 verticesToken_ = iC.consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"));
0028
0029 useTracks_ = iConfig.getParameter<bool>("useTracks");
0030
0031 }
0032 assoSelector_ = reco::modules::make<pat::VertexAssociationSelector>(iConfig);
0033 } else {
0034 enabled_ = false;
0035 }
0036 if (!playback_) {
0037 ttToken_ = iC.esConsumes(edm::ESInputTag("", "TransientTrackBuilder"));
0038 }
0039 }
0040
0041 void pat::helper::VertexingHelper::newEvent(const edm::Event &iEvent) {
0042 if (playback_) {
0043 iEvent.getByToken(vertexAssociationsToken_, vertexAssoMap_);
0044 } else {
0045 iEvent.getByToken(verticesToken_, vertexHandle_);
0046 }
0047 }
0048
0049 void pat::helper::VertexingHelper::newEvent(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0050 newEvent(iEvent);
0051 if (!playback_)
0052 ttBuilder_ = iSetup.getHandle(ttToken_);
0053 }
0054
0055 pat::VertexAssociation pat::helper::VertexingHelper::associate(const reco::Candidate &c) const {
0056 if (playback_)
0057 throw cms::Exception("Configuration")
0058 << "VertexingHelper: if this module was configured to read associations from the event,"
0059 << " you must use 'operator()' passing a candidate ref, and not 'associate()' directly!\n";
0060
0061 reco::VertexCollection::const_iterator vtx, end;
0062 size_t ivtx;
0063 reco::TrackBaseRef tk;
0064 reco::TransientTrack tt;
0065 if (useTracks_) {
0066 if (!ttBuilder_.isValid())
0067 throw cms::Exception("Configuration")
0068 << "VertexingHelper: If you use 'useTracks', you must call newEvent(iEvent,iSetup)!\n";
0069 tk = getTrack_(c);
0070 if (tk.isNull())
0071 return pat::VertexAssociation();
0072 tt = ttBuilder_->build(*tk);
0073 }
0074 for (vtx = vertexHandle_->begin(), end = vertexHandle_->end(), ivtx = 0; vtx != end; ++vtx, ++ivtx) {
0075 pat::VertexAssociation association(reco::VertexRef(vertexHandle_, ivtx), tk);
0076 if (useTracks_ == false) {
0077 association.setDistances(c.vertex(), vtx->position(), vtx->error());
0078 } else {
0079 GlobalPoint vtxGP(vtx->x(), vtx->y(), vtx->z());
0080 TrajectoryStateClosestToPoint tscp = tt.trajectoryStateClosestToPoint(vtxGP);
0081 GlobalPoint trackPos = tscp.theState().position();
0082 AlgebraicSymMatrix33 trackErr = tscp.theState().cartesianError().matrix().Sub<AlgebraicSymMatrix33>(0, 0);
0083 association.setDistances(trackPos, vtx->position(), trackErr + vtx->error());
0084 }
0085 if (assoSelector_(association))
0086 return association;
0087 }
0088 return pat::VertexAssociation();
0089 }
0090
0091 reco::TrackBaseRef pat::helper::VertexingHelper::getTrack_(const reco::Candidate &c) const {
0092 const reco::RecoCandidate *rc = dynamic_cast<const reco::RecoCandidate *>(&c);
0093 if (rc != nullptr) {
0094 return rc->bestTrackRef();
0095 }
0096 const reco::PFCandidate *pfc = dynamic_cast<const reco::PFCandidate *>(&c);
0097 if (pfc != nullptr) {
0098 return reco::TrackBaseRef(pfc->trackRef());
0099 }
0100
0101 return reco::TrackBaseRef();
0102 }