CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HIProtoTrackSelector.h
Go to the documentation of this file.
1 #ifndef HIProtoTrackSelection_h
2 #define HIProtoTrackSelection_h
3 
6 
9 
12 
14 
15 #include <algorithm>
16 #include <iostream>
17 
23 {
24 
25  public:
26  // input collection type
28 
29  // output collection type
30  typedef std::vector<const reco::Track *> container;
31 
32  // iterator over result collection type.
33  typedef container::const_iterator const_iterator;
34 
35  // constructor from parameter set configurability
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  {};
44 
45  // select object from a collection and possibly event content
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  }
117 
118  // iterators over selected objects: collection begin
119  const_iterator begin() const { return selected_.begin(); }
120 
121  // iterators over selected objects: collection end
122  const_iterator end() const { return selected_.end(); }
123 
124  // true if no object has been selected
125  size_t size() const { return selected_.size(); }
126 
127 
128  private:
132  double ptMin_;
133  double nSigmaZ_;
134  double minZCut_;
136 
137 };
138 
139 
140 #endif
HIProtoTrackSelector(const edm::ParameterSet &iConfig)
tuple d0
Definition: debug_cff.py:3
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
const_iterator begin() const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
int iEvent
Definition: GenABIO.cc:243
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
void select(edm::Handle< reco::TrackCollection > &TCH, const edm::Event &iEvent, const edm::EventSetup &iSetup)
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
std::vector< const reco::Track * > container
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
reco::TrackCollection collection
const_iterator end() const
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
container::const_iterator const_iterator