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)
 
void select (edm::Handle< reco::TrackCollection > &TCH, const edm::Event &iEvent, const edm::EventSetup &iSetup)
 
size_t size () const
 

Private Attributes

edm::InputTag beamSpotLabel_
 
double maxD0Significance_
 
double minZCut_
 
double nSigmaZ_
 
double ptMin_
 
container selected_
 
edm::InputTag vertexCollection_
 

Detailed Description

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

Definition at line 22 of file HIProtoTrackSelector.h.

Member Typedef Documentation

Definition at line 27 of file HIProtoTrackSelector.h.

typedef container::const_iterator HIProtoTrackSelector::const_iterator

Definition at line 33 of file HIProtoTrackSelector.h.

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

Definition at line 30 of file HIProtoTrackSelector.h.

Constructor & Destructor Documentation

HIProtoTrackSelector::HIProtoTrackSelector ( const edm::ParameterSet iConfig)
inline

Definition at line 36 of file HIProtoTrackSelector.h.

36  :
37  vertexCollection_(iConfig.getParameter<edm::InputTag>("VertexCollection")),
38  beamSpotLabel_(iConfig.getParameter<edm::InputTag>("beamSpotLabel")),
39  ptMin_(iConfig.getParameter<double>("ptMin")),
40  nSigmaZ_(iConfig.getParameter<double>("nSigmaZ")),
41  minZCut_(iConfig.getParameter<double>("minZCut")),
42  maxD0Significance_(iConfig.getParameter<double>("maxD0Significance"))
43  {};
T getParameter(std::string const &) const

Member Function Documentation

const_iterator HIProtoTrackSelector::begin ( void  ) const
inline

Definition at line 119 of file HIProtoTrackSelector.h.

References selected_.

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

Definition at line 122 of file HIProtoTrackSelector.h.

References selected_.

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

Definition at line 46 of file HIProtoTrackSelector.h.

References beamSpotLabel_, reco::BeamSpot::BeamWidthX(), reco::BeamSpot::BeamWidthY(), trackerHits::c, debug_cff::d0, edm::Event::getByLabel(), edm::HandleBase::isValid(), max(), maxD0Significance_, minZCut_, nSigmaZ_, reco::BeamSpot::position(), edm::Handle< T >::product(), ptMin_, selected_, mathSSE::sqrt(), and vertexCollection_.

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

Definition at line 125 of file HIProtoTrackSelector.h.

References selected_.

125 { return selected_.size(); }

Member Data Documentation

edm::InputTag HIProtoTrackSelector::beamSpotLabel_
private

Definition at line 131 of file HIProtoTrackSelector.h.

Referenced by select().

double HIProtoTrackSelector::maxD0Significance_
private

Definition at line 135 of file HIProtoTrackSelector.h.

Referenced by select().

double HIProtoTrackSelector::minZCut_
private

Definition at line 134 of file HIProtoTrackSelector.h.

Referenced by select().

double HIProtoTrackSelector::nSigmaZ_
private

Definition at line 133 of file HIProtoTrackSelector.h.

Referenced by select().

double HIProtoTrackSelector::ptMin_
private

Definition at line 132 of file HIProtoTrackSelector.h.

Referenced by select().

container HIProtoTrackSelector::selected_
private

Definition at line 129 of file HIProtoTrackSelector.h.

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

edm::InputTag HIProtoTrackSelector::vertexCollection_
private

Definition at line 130 of file HIProtoTrackSelector.h.

Referenced by select().