CMS 3D CMS Logo

HIProtoTrackSelector.h
Go to the documentation of this file.
1 #ifndef HIProtoTrackSelection_h
2 #define HIProtoTrackSelection_h
3 
5 
8 
11 
14 
16 
17 #include <algorithm>
18 #include <iostream>
19 
25 {
26 
27  public:
28  // input collection type
30 
31  // output collection type
32  typedef std::vector<const reco::Track *> container;
33 
34  // iterator over result collection type.
35  typedef container::const_iterator const_iterator;
36 
37  // constructor from parameter set configurability
39  vertexCollection_(iConfig.getParameter<edm::InputTag>("VertexCollection")),
41  beamSpot_(iConfig.getParameter<edm::InputTag>("beamSpotLabel")),
42  beamSpotToken_(iC.consumes<reco::BeamSpot>(beamSpot_)),
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  {};
48 
49  // select object from a collection and possibly event content
51  {
52  selected_.clear();
53 
54  const collection & c = *(TCH.product());
55 
56  // Get fast reco vertex
60 
61  math::XYZPoint vtxPoint(0.0,0.0,0.0);
62  double vzErr =0.0;
63 
64  if(!vertices->empty()) {
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  }
123 
124  // iterators over selected objects: collection begin
125  const_iterator begin() const { return selected_.begin(); }
126 
127  // iterators over selected objects: collection end
128  const_iterator end() const { return selected_.end(); }
129 
130  // true if no object has been selected
131  size_t size() const { return selected_.size(); }
132 
133 
134  private:
135  container selected_;
140  double ptMin_;
141  double nSigmaZ_;
142  double minZCut_;
144 
145 };
146 
147 
148 #endif
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
const_iterator begin() const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< reco::VertexCollection > vertexCollectionToken_
int iEvent
Definition: GenABIO.cc:230
T sqrt(T t)
Definition: SSEVec.h:18
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:86
void select(edm::Handle< reco::TrackCollection > &TCH, const edm::Event &iEvent, const edm::EventSetup &iSetup)
bool isValid() const
Definition: HandleBase.h:74
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
std::vector< const reco::Track * > container
T const * product() const
Definition: Handle.h:81
HIProtoTrackSelector(const edm::ParameterSet &iConfig, edm::ConsumesCollector &&iC)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
reco::TrackCollection collection
const_iterator end() const
double BeamWidthY() const
beam width Y
Definition: BeamSpot.h:88
fixed size matrix
HLT enums.
const Point & position() const
position
Definition: BeamSpot.h:62
container::const_iterator const_iterator