CMS 3D CMS Logo

HIProtoTrackSelector.h
Go to the documentation of this file.
1 #ifndef HIProtoTrackSelection_h
2 #define HIProtoTrackSelection_h
3 
7 
10 
13 
16 
18 
19 #include <algorithm>
20 #include <iostream>
21 
27 {
28 
29  public:
30  // input collection type
32 
33  // output collection type
34  typedef std::vector<const reco::Track *> container;
35 
36  // iterator over result collection type.
37  typedef container::const_iterator const_iterator;
38 
39  // constructor from parameter set configurability
41  vertexCollection_(iConfig.getParameter<edm::InputTag>("VertexCollection")),
43  beamSpot_(iConfig.getParameter<edm::InputTag>("beamSpotLabel")),
44  beamSpotToken_(iC.consumes<reco::BeamSpot>(beamSpot_)),
45  ptMin_(iConfig.getParameter<double>("ptMin")),
46  nSigmaZ_(iConfig.getParameter<double>("nSigmaZ")),
47  minZCut_(iConfig.getParameter<double>("minZCut")),
48  maxD0Significance_(iConfig.getParameter<double>("maxD0Significance"))
49  {};
50 
51  // select object from a collection and possibly event content
53  {
54  selected_.clear();
55 
56  const collection & c = *(TCH.product());
57 
58  // Get fast reco vertex
62 
63  math::XYZPoint vtxPoint(0.0,0.0,0.0);
64  double vzErr =0.0;
65 
66  if(!vertices->empty()) {
67  vtxPoint=vertices->begin()->position();
68  vzErr=vertices->begin()->zError();
69  edm::LogInfo("HeavyIonVertexing") << "Select prototracks compatible with median vertex"
70  << "\n vz = " << vtxPoint.Z()
71  << "\n " << nSigmaZ_ << " vz sigmas = " << vzErr*nSigmaZ_
72  << "\n cut at = " << std::max(vzErr*nSigmaZ_,minZCut_);
73  }
74  // Supress this warning, since events w/ no vertex are expected
75  //else {
76  //edm::LogError("HeavyIonVertexing") << "No vertex found in collection '" << vertexCollection_ << "'";
77  //}
78 
79  // Get beamspot
81  edm::Handle<reco::BeamSpot> beamSpotHandle;
82  iEvent.getByToken(beamSpotToken_, beamSpotHandle);
83 
84  math::XYZPoint bsPoint(0.0,0.0,0.0);
85  double bsWidth = 0.0;
86 
87  if ( beamSpotHandle.isValid() ) {
88  beamSpot = *beamSpotHandle;
89  bsPoint = beamSpot.position();
90  bsWidth = sqrt(beamSpot.BeamWidthX()*beamSpot.BeamWidthY());
91  edm::LogInfo("HeavyIonVertexing") << "Select prototracks compatible with beamspot"
92  << "\n (x,y,z) = (" << bsPoint.X() << "," << bsPoint.Y() << "," << bsPoint.Z() << ")"
93  << "\n width = " << bsWidth
94  << "\n cut at d0/d0sigma = " << maxD0Significance_;
95  } else {
96  edm::LogError("HeavyIonVertexing") << "No beam spot available from '" << beamSpot_ << "\n";
97  }
98 
99 
100  // Do selection
101  int nSelected=0;
102  int nRejected=0;
103  double d0=0.0;
104  double d0sigma=0.0;
105  for (reco::TrackCollection::const_iterator trk = c.begin(); trk != c.end(); ++ trk)
106  {
107 
108  d0 = -1.*trk->dxy(bsPoint);
109  d0sigma = sqrt(trk->d0Error()*trk->d0Error() + bsWidth*bsWidth);
110  if ( trk->pt() > ptMin_ // keep only tracks above ptMin
111  && fabs(d0/d0sigma) < maxD0Significance_ // keep only tracks with D0 significance less than cut
112  && fabs(trk->dz(vtxPoint)) < std::max(vzErr*nSigmaZ_,minZCut_) // within minZCut, nSigmaZ of fast vertex
113  )
114  {
115  nSelected++;
116  selected_.push_back( & * trk );
117  }
118  else
119  nRejected++;
120  }
121 
122  edm::LogInfo("HeavyIonVertexing") << "selected " << nSelected << " prototracks out of " << nRejected+nSelected << "\n";
123 
124  }
125 
126  // iterators over selected objects: collection begin
127  const_iterator begin() const { return selected_.begin(); }
128 
129  // iterators over selected objects: collection end
130  const_iterator end() const { return selected_.end(); }
131 
132  // true if no object has been selected
133  size_t size() const { return selected_.size(); }
134 
135 
136  private:
137  container selected_;
142  double ptMin_;
143  double nSigmaZ_;
144  double minZCut_;
146 
147 };
148 
149 
150 #endif
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:15
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:224
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:74
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