CMS 3D CMS Logo

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

#include <HIProtoTrackSelector.h>

Public Types

typedef reco::TrackCollection collection
 
typedef container::const_iterator const_iterator
 
typedef std::vector< const reco::Track * > container
 

Public Member Functions

const_iterator begin () const
 
const_iterator end () const
 
 HIProtoTrackSelector (const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
 
void select (edm::Handle< reco::TrackCollection > &TCH, const edm::Event &iEvent, const edm::EventSetup &iSetup)
 
size_t size () const
 

Private Attributes

edm::InputTag beamSpot_
 
edm::EDGetTokenT< reco::BeamSpotbeamSpotToken_
 
double maxD0Significance_
 
double minZCut_
 
double nSigmaZ_
 
double ptMin_
 
container selected_
 
edm::InputTag vertexCollection_
 
edm::EDGetTokenT< reco::VertexCollectionvertexCollectionToken_
 

Detailed Description

Selector to select prototracks that pass certain kinematic cuts based on fast vertex

Definition at line 26 of file HIProtoTrackSelector.h.

Member Typedef Documentation

◆ collection

Definition at line 29 of file HIProtoTrackSelector.h.

◆ const_iterator

typedef container::const_iterator HIProtoTrackSelector::const_iterator

Definition at line 35 of file HIProtoTrackSelector.h.

◆ container

typedef std::vector<const reco::Track*> HIProtoTrackSelector::container

Definition at line 32 of file HIProtoTrackSelector.h.

Constructor & Destructor Documentation

◆ HIProtoTrackSelector()

HIProtoTrackSelector::HIProtoTrackSelector ( const edm::ParameterSet iConfig,
edm::ConsumesCollector &&  iC 
)
inline

Definition at line 38 of file HIProtoTrackSelector.h.

39  : vertexCollection_(iConfig.getParameter<edm::InputTag>("VertexCollection")),
41  beamSpot_(iConfig.getParameter<edm::InputTag>("beamSpotLabel")),
43  ptMin_(iConfig.getParameter<double>("ptMin")),
44  nSigmaZ_(iConfig.getParameter<double>("nSigmaZ")),
45  minZCut_(iConfig.getParameter<double>("minZCut")),
46  maxD0Significance_(iConfig.getParameter<double>("maxD0Significance")){};
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::VertexCollection > vertexCollectionToken_
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_

Member Function Documentation

◆ begin()

const_iterator HIProtoTrackSelector::begin ( void  ) const
inline

Definition at line 119 of file HIProtoTrackSelector.h.

References selected_.

119 { return selected_.begin(); }

◆ end()

const_iterator HIProtoTrackSelector::end ( void  ) const
inline

Definition at line 122 of file HIProtoTrackSelector.h.

References selected_.

122 { return selected_.end(); }

◆ select()

void HIProtoTrackSelector::select ( edm::Handle< reco::TrackCollection > &  TCH,
const edm::Event iEvent,
const edm::EventSetup iSetup 
)
inline

Definition at line 49 of file HIProtoTrackSelector.h.

References pwdgSkimBPark_cfi::beamSpot, beamSpot_, beamSpotToken_, HltBtagPostValidation_cff::c, d0, iEvent, edm::HandleBase::isValid(), SiStripPI::max, maxD0Significance_, minZCut_, nSigmaZ_, edm::Handle< T >::product(), ptMin_, selected_, mathSSE::sqrt(), vertexCollectionToken_, and AlignmentTracksFromVertexSelector_cfi::vertices.

49  {
50  selected_.clear();
51 
52  const collection& c = *(TCH.product());
53 
54  // Get fast reco vertex
56  iEvent.getByToken(vertexCollectionToken_, vc);
58 
59  math::XYZPoint vtxPoint(0.0, 0.0, 0.0);
60  double vzErr = 0.0;
61 
62  if (!vertices->empty()) {
63  vtxPoint = vertices->begin()->position();
64  vzErr = vertices->begin()->zError();
65  edm::LogInfo("HeavyIonVertexing") << "Select prototracks compatible with median vertex"
66  << "\n vz = " << vtxPoint.Z() << "\n " << nSigmaZ_
67  << " vz sigmas = " << vzErr * nSigmaZ_
68  << "\n cut at = " << std::max(vzErr * nSigmaZ_, minZCut_);
69  }
70  // Supress this warning, since events w/ no vertex are expected
71  //else {
72  //edm::LogError("HeavyIonVertexing") << "No vertex found in collection '" << vertexCollection_ << "'";
73  //}
74 
75  // Get beamspot
77  edm::Handle<reco::BeamSpot> beamSpotHandle;
78  iEvent.getByToken(beamSpotToken_, beamSpotHandle);
79 
80  math::XYZPoint bsPoint(0.0, 0.0, 0.0);
81  double bsWidth = 0.0;
82 
83  if (beamSpotHandle.isValid()) {
84  beamSpot = *beamSpotHandle;
85  bsPoint = beamSpot.position();
86  bsWidth = sqrt(beamSpot.BeamWidthX() * beamSpot.BeamWidthY());
87  edm::LogInfo("HeavyIonVertexing") << "Select prototracks compatible with beamspot"
88  << "\n (x,y,z) = (" << bsPoint.X() << "," << bsPoint.Y() << "," << bsPoint.Z()
89  << ")"
90  << "\n width = " << bsWidth
91  << "\n cut at d0/d0sigma = " << maxD0Significance_;
92  } else {
93  edm::LogError("HeavyIonVertexing") << "No beam spot available from '" << beamSpot_ << "\n";
94  }
95 
96  // Do selection
97  int nSelected = 0;
98  int nRejected = 0;
99  double d0 = 0.0;
100  double d0sigma = 0.0;
101  for (reco::TrackCollection::const_iterator trk = c.begin(); trk != c.end(); ++trk) {
102  d0 = -1. * trk->dxy(bsPoint);
103  d0sigma = sqrt(trk->d0Error() * trk->d0Error() + bsWidth * bsWidth);
104  if (trk->pt() > ptMin_ // keep only tracks above ptMin
105  && fabs(d0 / d0sigma) < maxD0Significance_ // keep only tracks with D0 significance less than cut
106  && fabs(trk->dz(vtxPoint)) < std::max(vzErr * nSigmaZ_, minZCut_) // within minZCut, nSigmaZ of fast vertex
107  ) {
108  nSelected++;
109  selected_.push_back(&*trk);
110  } else
111  nRejected++;
112  }
113 
114  edm::LogInfo("HeavyIonVertexing") << "selected " << nSelected << " prototracks out of " << nRejected + nSelected
115  << "\n";
116  }
T const * product() const
Definition: Handle.h:70
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
Log< level::Error, false > LogError
edm::EDGetTokenT< reco::VertexCollection > vertexCollectionToken_
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
Log< level::Info, false > LogInfo
static constexpr float d0
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
reco::TrackCollection collection
bool isValid() const
Definition: HandleBase.h:70

◆ size()

size_t HIProtoTrackSelector::size ( void  ) const
inline

Definition at line 125 of file HIProtoTrackSelector.h.

References selected_.

Referenced by ntupleDataFormat._Collection::__iter__(), and ntupleDataFormat._Collection::__len__().

125 { return selected_.size(); }

Member Data Documentation

◆ beamSpot_

edm::InputTag HIProtoTrackSelector::beamSpot_
private

Definition at line 131 of file HIProtoTrackSelector.h.

Referenced by select().

◆ beamSpotToken_

edm::EDGetTokenT<reco::BeamSpot> HIProtoTrackSelector::beamSpotToken_
private

Definition at line 132 of file HIProtoTrackSelector.h.

Referenced by select().

◆ maxD0Significance_

double HIProtoTrackSelector::maxD0Significance_
private

Definition at line 136 of file HIProtoTrackSelector.h.

Referenced by select().

◆ minZCut_

double HIProtoTrackSelector::minZCut_
private

Definition at line 135 of file HIProtoTrackSelector.h.

Referenced by select().

◆ nSigmaZ_

double HIProtoTrackSelector::nSigmaZ_
private

Definition at line 134 of file HIProtoTrackSelector.h.

Referenced by select().

◆ ptMin_

double HIProtoTrackSelector::ptMin_
private

Definition at line 133 of file HIProtoTrackSelector.h.

Referenced by select().

◆ selected_

container HIProtoTrackSelector::selected_
private

Definition at line 128 of file HIProtoTrackSelector.h.

Referenced by begin(), end(), select(), and size().

◆ vertexCollection_

edm::InputTag HIProtoTrackSelector::vertexCollection_
private

Definition at line 129 of file HIProtoTrackSelector.h.

◆ vertexCollectionToken_

edm::EDGetTokenT<reco::VertexCollection> HIProtoTrackSelector::vertexCollectionToken_
private

Definition at line 130 of file HIProtoTrackSelector.h.

Referenced by select().