CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Attributes
AlignmentTrackFromVertexSelector Class Reference

#include <AlignmentTracksFromVertexSelector.h>

Public Types

typedef std::vector< const reco::Track * > Tracks
 

Public Member Functions

 AlignmentTrackFromVertexSelector (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 constructor More...
 
const reco::VertexfindClosestVertex (const reco::TrackCollection &leptonTracks, const reco::VertexCollection *vertices, const edm::EventSetup &setup) const
 
Tracks select (const edm::Handle< reco::TrackCollection > &tc, const edm::Event &evt, const edm::EventSetup &setup) const
 select tracks More...
 
 ~AlignmentTrackFromVertexSelector ()
 destructor More...
 

Private Attributes

edm::EDGetTokenT< reco::TrackCollectiondiLeptonToken_
 
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordttbESToken_
 
bool useClosestVertex_
 
unsigned int vertexIndex_
 
edm::EDGetTokenT< reco::VertexCollectionvertexToken_
 

Detailed Description

Definition at line 21 of file AlignmentTracksFromVertexSelector.h.

Member Typedef Documentation

◆ Tracks

Definition at line 23 of file AlignmentTracksFromVertexSelector.h.

Constructor & Destructor Documentation

◆ AlignmentTrackFromVertexSelector()

AlignmentTrackFromVertexSelector::AlignmentTrackFromVertexSelector ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

constructor

Definition at line 17 of file AlignmentTracksFromVertexSelector.cc.

19  : ttbESToken_(
20  iC.esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", "TransientTrackBuilder"))),
21  vertexToken_(iC.consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("vertices"))),
22  diLeptonToken_(iC.consumes<reco::TrackCollection>(cfg.getParameter<edm::InputTag>("leptonTracks"))),
23  useClosestVertex_(cfg.getParameter<bool>("useClosestVertexToDilepton")),
24  vertexIndex_(cfg.getParameter<unsigned int>("vertexIndex")) {}
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< reco::TrackCollection > diLeptonToken_
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > ttbESToken_
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9

◆ ~AlignmentTrackFromVertexSelector()

AlignmentTrackFromVertexSelector::~AlignmentTrackFromVertexSelector ( )

destructor

Definition at line 27 of file AlignmentTracksFromVertexSelector.cc.

27 {}

Member Function Documentation

◆ findClosestVertex()

const reco::Vertex * AlignmentTrackFromVertexSelector::findClosestVertex ( const reco::TrackCollection leptonTracks,
const reco::VertexCollection vertices,
const edm::EventSetup setup 
) const

Definition at line 30 of file AlignmentTracksFromVertexSelector.cc.

References TransientTrackBuilder::build(), counter, VertexDistance3D::distance(), TransientVertex::isValid(), singleTopDQM_cfi::setup, HLT_2024v12_cff::track, ttbESToken_, Measurement1D::value(), KalmanVertexFitter::vertex(), AlignmentTracksFromVertexSelector_cfi::vertices, and L1BJetProducer_cff::vtx.

Referenced by select().

32  {
33  reco::Vertex* defaultVtx = nullptr;
34 
35  // fill the transient track collection with the lepton tracks
36  const TransientTrackBuilder* theB = &setup.getData(ttbESToken_);
37  std::vector<reco::TransientTrack> tks;
38  for (const auto& track : leptonTracks) {
39  reco::TransientTrack trajectory = theB->build(track);
40  tks.push_back(trajectory);
41  }
42 
43  // compute the secondary vertex
44  TransientVertex aTransVtx;
45  KalmanVertexFitter kalman(true);
46  aTransVtx = kalman.vertex(tks);
47 
48  if (!aTransVtx.isValid())
49  return defaultVtx;
50 
51  // find the closest vertex to the secondary vertex in 3D
52  VertexDistance3D vertTool3D;
53  float minD = 9999.;
54  int closestVtxIndex = 0;
55  int counter = 0;
56  for (const auto& vtx : *vertices) {
57  double dist3D = vertTool3D.distance(aTransVtx, vtx).value();
58  if (dist3D < minD) {
59  minD = dist3D;
60  closestVtxIndex = counter;
61  }
62  counter++;
63  }
64  if ((*vertices).at(closestVtxIndex).isValid()) {
65  return &(vertices->at(closestVtxIndex));
66  } else {
67  return defaultVtx;
68  }
69 }
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > ttbESToken_
reco::TransientTrack build(const reco::Track *p) const
bool isValid() const
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
double value() const
Definition: Measurement1D.h:25
static std::atomic< unsigned int > counter

◆ select()

AlignmentTrackFromVertexSelector::Tracks AlignmentTrackFromVertexSelector::select ( const edm::Handle< reco::TrackCollection > &  tc,
const edm::Event evt,
const edm::EventSetup setup 
) const

select tracks

Definition at line 72 of file AlignmentTracksFromVertexSelector.cc.

References diLeptonToken_, spr::find(), findClosestVertex(), edm::Event::get(), edm::Event::getHandle(), edm::HandleBase::isValid(), submitPVResolutionJobs::key, edm::Ref< C, T, F >::key(), LogDebug, edm::Handle< T >::product(), mps_fire::result, singleTopDQM_cfi::setup, jetUpdater_cfi::sort, reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), useClosestVertex_, vertexIndex_, vertexToken_, and AlignmentTracksFromVertexSelector_cfi::vertices.

Referenced by TrackFromVertexConfigSelector::select().

73  {
74  Tracks result;
75 
76  std::vector<unsigned int> thePVkeys;
77 
78  // get collection of reconstructed vertices from event
80 
81  // get collection of the di-lepton traxks
82  const auto& leptonTracks = evt.get(diLeptonToken_);
83 
84  // fill the vector of keys
85  if (vertexHandle.isValid()) {
86  const reco::VertexCollection* vertices = vertexHandle.product();
87  const reco::Vertex* theVtx = nullptr;
88 
89  if (useClosestVertex_) {
90  theVtx = findClosestVertex(leptonTracks, vertices, setup);
91  } else {
92  if ((*vertices).at(vertexIndex_).isValid()) {
93  theVtx = &(vertices->at(vertexIndex_));
94  }
95  }
96 
97  if (theVtx) {
98  for (auto tv = theVtx->tracks_begin(); tv != theVtx->tracks_end(); tv++) {
99  if (tv->isNonnull()) {
100  const reco::TrackRef trackRef = tv->castTo<reco::TrackRef>();
101  thePVkeys.push_back(trackRef.key());
102  }
103  }
104  }
105  }
106 
107  std::sort(thePVkeys.begin(), thePVkeys.end());
108 
109  LogDebug("AlignmentTrackFromVertexSelector")
110  << "after collecting the PV keys: thePVkeys.size()" << thePVkeys.size() << std::endl;
111  for (const auto& key : thePVkeys) {
112  LogDebug("AlignmentTrackFromVertexSelector") << key << ", ";
113  }
114  LogDebug("AlignmentTrackFromVertexSelector") << std::endl;
115 
116  if (tc.isValid()) {
117  int indx(0);
118  // put the track in the collection is it was used for the vertex
119  for (reco::TrackCollection::const_iterator tk = tc->begin(); tk != tc->end(); ++tk, ++indx) {
120  reco::TrackRef trackRef = reco::TrackRef(tc, indx);
121  if (std::find(thePVkeys.begin(), thePVkeys.end(), trackRef.key()) != thePVkeys.end()) {
122  LogDebug("AlignmentTrackFromVertexSelector") << "track index: " << indx << "filling result vector" << std::endl;
123  result.push_back(&(*tk));
124  } // if a valid key is found
125  } // end loop over tracks
126  } // if the handle is valid
127  return result;
128 }
std::vector< const reco::Track * > Tracks
edm::EDGetTokenT< reco::TrackCollection > diLeptonToken_
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:344
edm::EDGetTokenT< reco::VertexCollection > vertexToken_
T const * product() const
Definition: Handle.h:70
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const reco::Vertex * findClosestVertex(const reco::TrackCollection &leptonTracks, const reco::VertexCollection *vertices, const edm::EventSetup &setup) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.h:111
key_type key() const
Accessor for product key.
Definition: Ref.h:250
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.h:109
key
prepare the HTCondor submission files and eventually submit them
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
bool isValid() const
Definition: HandleBase.h:70
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:552
#define LogDebug(id)

Member Data Documentation

◆ diLeptonToken_

edm::EDGetTokenT<reco::TrackCollection> AlignmentTrackFromVertexSelector::diLeptonToken_
private

Definition at line 44 of file AlignmentTracksFromVertexSelector.h.

Referenced by select().

◆ ttbESToken_

const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> AlignmentTrackFromVertexSelector::ttbESToken_
private

Definition at line 42 of file AlignmentTracksFromVertexSelector.h.

Referenced by findClosestVertex().

◆ useClosestVertex_

bool AlignmentTrackFromVertexSelector::useClosestVertex_
private

Definition at line 45 of file AlignmentTracksFromVertexSelector.h.

Referenced by select().

◆ vertexIndex_

unsigned int AlignmentTrackFromVertexSelector::vertexIndex_
private

Definition at line 46 of file AlignmentTracksFromVertexSelector.h.

Referenced by select().

◆ vertexToken_

edm::EDGetTokenT<reco::VertexCollection> AlignmentTrackFromVertexSelector::vertexToken_
private

Definition at line 43 of file AlignmentTracksFromVertexSelector.h.

Referenced by select().