CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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::VertexCollection
vertexCollectionToken_
 

Detailed Description

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

Definition at line 24 of file HIProtoTrackSelector.h.

Member Typedef Documentation

Definition at line 29 of file HIProtoTrackSelector.h.

typedef container::const_iterator HIProtoTrackSelector::const_iterator

Definition at line 35 of file HIProtoTrackSelector.h.

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

Definition at line 32 of file HIProtoTrackSelector.h.

Constructor & Destructor Documentation

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

Definition at line 38 of file HIProtoTrackSelector.h.

38  :
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"))
47  {};
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
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

const_iterator HIProtoTrackSelector::begin ( void  ) const
inline

Definition at line 125 of file HIProtoTrackSelector.h.

References selected_.

125 { return selected_.begin(); }
const_iterator HIProtoTrackSelector::end ( void  ) const
inline

Definition at line 128 of file HIProtoTrackSelector.h.

References selected_.

128 { return selected_.end(); }
void HIProtoTrackSelector::select ( edm::Handle< reco::TrackCollection > &  TCH,
const edm::Event iEvent,
const edm::EventSetup iSetup 
)
inline

Definition at line 50 of file HIProtoTrackSelector.h.

References SiPixelRawToDigiRegional_cfi::beamSpot, beamSpot_, beamSpotToken_, reco::BeamSpot::BeamWidthX(), reco::BeamSpot::BeamWidthY(), EnergyCorrector::c, edm::Event::getByToken(), edm::HandleBase::isValid(), bookConverter::max, maxD0Significance_, minZCut_, nSigmaZ_, reco::BeamSpot::position(), edm::Handle< T >::product(), ptMin_, selected_, mathSSE::sqrt(), and vertexCollectionToken_.

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

Definition at line 131 of file HIProtoTrackSelector.h.

References selected_.

131 { return selected_.size(); }

Member Data Documentation

edm::InputTag HIProtoTrackSelector::beamSpot_
private

Definition at line 138 of file HIProtoTrackSelector.h.

Referenced by select().

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

Definition at line 139 of file HIProtoTrackSelector.h.

Referenced by select().

double HIProtoTrackSelector::maxD0Significance_
private

Definition at line 143 of file HIProtoTrackSelector.h.

Referenced by select().

double HIProtoTrackSelector::minZCut_
private

Definition at line 142 of file HIProtoTrackSelector.h.

Referenced by select().

double HIProtoTrackSelector::nSigmaZ_
private

Definition at line 141 of file HIProtoTrackSelector.h.

Referenced by select().

double HIProtoTrackSelector::ptMin_
private

Definition at line 140 of file HIProtoTrackSelector.h.

Referenced by select().

container HIProtoTrackSelector::selected_
private

Definition at line 135 of file HIProtoTrackSelector.h.

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

edm::InputTag HIProtoTrackSelector::vertexCollection_
private

Definition at line 136 of file HIProtoTrackSelector.h.

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

Definition at line 137 of file HIProtoTrackSelector.h.

Referenced by select().