CMS 3D CMS Logo

TrackCutClassifier.cc
Go to the documentation of this file.
2 
3 
5 
7 
8 #include <cassert>
9 
10 #include "getBestVertex.h"
11 #include "powN.h"
12 
13 
14 namespace {
15 
16  void fillArrayF(float * x,const edm::ParameterSet & cfg, const char * name) {
17  auto v = cfg.getParameter< std::vector<double> >(name);
18  assert(v.size()==3);
20  }
21 
22  void fillArrayI(int * x,const edm::ParameterSet & cfg, const char * name) {
23  auto v = cfg.getParameter< std::vector<int> >(name);
24  assert(v.size()==3);
26  }
27 
28  // fake mva value to return for loose,tight,hp
29  constexpr float mvaVal[3] = {-.5,.5,1.};
30 
31  template<typename T,typename Comp>
32  inline float cut(T val, const T * cuts, Comp comp) {
33  for (int i=2; i>=0; --i) {
34  if ( comp(val,cuts[i]) ) return mvaVal[i];
35  }
36  return -1.f;
37  }
38 
39  inline float chi2n(reco::Track const & tk) { return tk.normalizedChi2();}
40 
41  inline float relPtErr(reco::Track const & tk) {
42  return (tk.pt() != 0. ? float(tk.ptError())/float(tk.pt()) : 9999999.);
43  }
44 
45  inline int lostLayers(reco::Track const & tk) {
47  }
48 
49  inline int n3DLayers(reco::Track const & tk, bool isHLT) {
50  uint32_t nlayers3D = tk.hitPattern().pixelLayersWithMeasurement();
51  if (!isHLT)
53  else {
54  size_t count3D = 0;
55  for ( auto it = tk.recHitsBegin(), et = tk.recHitsEnd(); it!=et; ++it) {
56  const TrackingRecHit* hit = (*it);
57  if (!trackerHitRTTI::isFromDetOrFast(*hit) ) continue;
58 
59  if ( hit->dimension()==2 ) {
60  auto const & thit = static_cast<BaseTrackerRecHit const&>(*hit);
61  if (thit.isMatched()) count3D++;
62  }
63  }
64  nlayers3D += count3D;
65  }
66  return nlayers3D;
67  }
68 
69  inline int nHits(reco::Track const & tk) {
70  return tk.numberOfValidHits();
71  }
72 
73  inline int nPixelHits(reco::Track const & tk) {
74  return tk.hitPattern().numberOfValidPixelHits();
75  }
76 
77  inline float dz(reco::Track const & trk, Point const & bestVertex) {
78  return std::abs(trk.dz(bestVertex));
79  }
80  inline float dr(reco::Track const & trk, Point const & bestVertex) {
81  return std::abs(trk.dxy(bestVertex));
82  }
83 
84  inline void dzCut_par1(reco::Track const & trk, int & nLayers, const float * par, const int * exp, float dzCut[]) {
85  float dzE = trk.dzError();
86  for (int i=2; i>=0; --i) {
87  dzCut[i] = powN(par[i]*nLayers,exp[i])*dzE;
88  }
89  }
90  inline void drCut_par1(reco::Track const & trk, int & nLayers, const float * par, const int * exp, float drCut[]) {
91  float drE = trk.d0Error();
92  for (int i=2; i>=0; --i) {
93  drCut[i] = powN(par[i]*nLayers,exp[i])*drE;
94  }
95  }
96 
97  inline void dzCut_par2(reco::Track const & trk, int & nLayers, const float * par, const int * exp, const float * d0err, const float * d0err_par, float dzCut[]) {
98  float pt = float(trk.pt());
99  float p = float(trk.p());
100 
101  for (int i=2; i>=0; --i) {
102  // parametrized d0 resolution for the track pt
103  float nomd0E = sqrt(d0err[i]*d0err[i]+(d0err_par[i]/pt)*(d0err_par[i]/pt));
104  // parametrized z0 resolution for the track pt and eta
105  float nomdzE = nomd0E*(p/pt); // cosh(eta):=abs(p)/pt
106 
107  dzCut[i] = powN(par[i]*nLayers,exp[i])*nomdzE;
108  }
109  }
110  inline void drCut_par2(reco::Track const & trk, int & nLayers, const float* par, const int * exp, const float * d0err, const float * d0err_par, float drCut[]) {
111  float pt = trk.pt();
112 
113  for (int i=2; i>=0; --i) {
114  // parametrized d0 resolution for the track pt
115  float nomd0E = sqrt(d0err[i]*d0err[i]+(d0err_par[i]/pt)*(d0err_par[i]/pt));
116 
117  drCut[i] = powN(par[i]*nLayers,exp[i])*nomd0E;
118  }
119  }
120 
121  inline void dzCut_wPVerror_par(reco::Track const & trk, int & nLayers, const float * par, const int * exp, Point const & bestVertexError, float dzCut[]) {
122  float dzE = trk.dzError();
123  float zPVerr = bestVertexError.z();
124 
125  float dzErrPV = std::sqrt(dzE*dzE+zPVerr*zPVerr);
126  for (int i=2; i>=0; --i) {
127  dzCut[i] = par[i]*dzErrPV;
128  if (exp[i] != 0)
129  dzCut[i] *= pow(nLayers,exp[i]);
130  }
131  }
132 
133  inline void drCut_wPVerror_par(reco::Track const & trk, int & nLayers, const float* par, const int * exp, Point const & bestVertexError, float drCut[]) {
134  float drE = trk.d0Error();
135  float rPVerr = sqrt(bestVertexError.x()*bestVertexError.y()); // shouldn't it be bestVertex.xError()*bestVertex.xError()+bestVertex.yError()*bestVertex.yError() ?!?!?
136 
137  float drErrPV = std::sqrt(drE*drE+rPVerr*rPVerr);
138  for (int i=2; i>=0; --i) {
139  drCut[i] = par[i]*drErrPV;
140  if (exp[i] != 0)
141  drCut[i] *= pow(nLayers,exp[i]);
142  }
143 
144  }
145 
146  struct Cuts {
147 
148  Cuts(const edm::ParameterSet & cfg) {
149  isHLT = cfg.getParameter<bool>("isHLT");
150  fillArrayF(minNdof, cfg,"minNdof");
151  fillArrayF(maxChi2, cfg,"maxChi2");
152  fillArrayF(maxChi2n, cfg,"maxChi2n");
153  fillArrayI(minHits4pass, cfg,"minHits4pass");
154  fillArrayI(minHits, cfg,"minHits");
155  fillArrayI(minPixelHits, cfg,"minPixelHits");
156  fillArrayI(min3DLayers, cfg,"min3DLayers");
157  fillArrayI(minLayers, cfg,"minLayers");
158  fillArrayI(maxLostLayers,cfg,"maxLostLayers");
159  fillArrayF(maxRelPtErr, cfg,"maxRelPtErr");
160  minNVtxTrk = cfg.getParameter<int>("minNVtxTrk");
161  fillArrayF(maxDz, cfg,"maxDz");
162  fillArrayF(maxDzWrtBS, cfg,"maxDzWrtBS");
163  fillArrayF(maxDr, cfg,"maxDr");
164  edm::ParameterSet dz_par = cfg.getParameter<edm::ParameterSet>("dz_par");
165  fillArrayI(dz_exp, dz_par,"dz_exp");
166  fillArrayF(dz_par1, dz_par,"dz_par1");
167  fillArrayF(dz_par2, dz_par,"dz_par2");
168  fillArrayF(dzWPVerr_par, dz_par,"dzWPVerr_par");
169  edm::ParameterSet dr_par = cfg.getParameter<edm::ParameterSet>("dr_par");
170  fillArrayI(dr_exp, dr_par,"dr_exp");
171  fillArrayF(dr_par1, dr_par,"dr_par1");
172  fillArrayF(dr_par2, dr_par,"dr_par2");
173  fillArrayF(d0err, dr_par,"d0err");
174  fillArrayF(d0err_par, dr_par,"d0err_par");
175  fillArrayF(drWPVerr_par, dr_par,"drWPVerr_par");
176  }
177 
178  void beginStream() {}
179  void initEvent(const edm::EventSetup&) {}
180 
181  float operator()(reco::Track const & trk,
182  reco::BeamSpot const & beamSpot,
183  reco::VertexCollection const & vertices) const {
184 
185  float ret = 1.f;
186  // minimum number of hits for by-passing the other checks
187  if ( minHits4pass[0] < std::numeric_limits<int>::max() ) {
188  ret = std::min(ret,cut(nHits(trk),minHits4pass,std::greater_equal<int>()));
189  if (ret==1.f) return ret;
190  }
191 
193  ret = std::min(ret,cut(relPtErr(trk),maxRelPtErr,std::less_equal<float>()) );
194  if (ret==-1.f) return ret;
195  }
196 
197  ret = std::min(ret,cut(float(trk.ndof()),minNdof,std::greater_equal<float>()) );
198  if (ret==-1.f) return ret;
199 
200  auto nLayers = trk.hitPattern().trackerLayersWithMeasurement();
201  ret = std::min(ret,cut(nLayers,minLayers,std::greater_equal<int>()));
202  if (ret==-1.f) return ret;
203 
204  ret = std::min(ret,cut(chi2n(trk)/float(nLayers),maxChi2n,std::less_equal<float>()));
205  if (ret==-1.f) return ret;
206 
207  ret = std::min(ret,cut(chi2n(trk),maxChi2,std::less_equal<float>()));
208  if (ret==-1.f) return ret;
209 
210  ret = std::min(ret,cut(n3DLayers(trk,isHLT),min3DLayers,std::greater_equal<int>()));
211  if (ret==-1.f) return ret;
212 
213  ret = std::min(ret,cut(nHits(trk),minHits,std::greater_equal<int>()));
214  if (ret==-1.f) return ret;
215 
216  ret = std::min(ret,cut(nPixelHits(trk),minPixelHits,std::greater_equal<int>()));
217  if (ret==-1.f) return ret;
218 
219  ret = std::min(ret,cut(lostLayers(trk),maxLostLayers,std::less_equal<int>()));
220  if (ret==-1.f) return ret;
221 
222  // original dz and dr cut
224 
225  // if not primaryVertices are reconstructed, check compatibility w.r.t. beam spot
226  Point bestVertex = getBestVertex(trk,vertices,minNVtxTrk); // min number of tracks [2 (=default) for offline, 3 for HLT]
227  float maxDzcut[3];
229  if (bestVertex.z() < -99998.) {
230  bestVertex = beamSpot.position();
231  std::copy(std::begin(maxDzWrtBS),std::end(maxDzWrtBS),std::begin(maxDzcut));
232  }
233  ret = std::min(ret,cut(dr(trk,bestVertex), maxDr,std::less<float>()));
234  if (ret==-1.f) return ret;
235 
236  ret = std::min(ret,cut(dz(trk,bestVertex), maxDzcut,std::less<float>()));
237  if (ret==-1.f) return ret;
238 
239  }
240 
241 
242  // parametrized dz and dr cut by using PV error
243  if (dzWPVerr_par[2]<std::numeric_limits<float>::max() || drWPVerr_par[2]<std::numeric_limits<float>::max()) {
244 
245  Point bestVertexError(-1.,-1.,-1.);
246  Point bestVertex = getBestVertex_withError(trk,vertices,bestVertexError,minNVtxTrk); // min number of tracks [2 (=default) for offline, 3 for HLT]
247 
248  float maxDz_par[3];
249  float maxDr_par[3];
250  dzCut_wPVerror_par(trk,nLayers,dzWPVerr_par,dz_exp,bestVertexError, maxDz_par);
251  drCut_wPVerror_par(trk,nLayers,drWPVerr_par,dr_exp,bestVertexError, maxDr_par);
252 
253  ret = std::min(ret,cut(dr(trk,bestVertex), maxDr_par,std::less<float>()));
254  if (ret==-1.f) return ret;
255 
256  ret = std::min(ret,cut(dz(trk,bestVertex), maxDr_par,std::less<float>()));
257  if (ret==-1.f) return ret;
258 
259  }
260 
261  // parametrized dz and dr cut by using their error
263 
264  float maxDz_par1[3];
265  float maxDr_par1[3];
266  dzCut_par1(trk,nLayers,dz_par1,dz_exp, maxDz_par1);
267  drCut_par1(trk,nLayers,dr_par1,dr_exp, maxDr_par1);
268 
269  float maxDz_par[3];
270  float maxDr_par[3];
271  std::copy(std::begin(maxDz_par1),std::end(maxDz_par1),std::begin(maxDz_par));
272  std::copy(std::begin(maxDr_par1),std::end(maxDr_par1),std::begin(maxDr_par));
273 
274  // parametrized dz and dr cut by using d0 and z0 resolution
276 
277  float maxDz_par2[3];
278  float maxDr_par2[3];
279  dzCut_par2(trk,nLayers,dz_par2,dz_exp,d0err,d0err_par, maxDz_par2);
280  drCut_par2(trk,nLayers,dr_par2,dr_exp,d0err,d0err_par, maxDr_par2);
281 
282 
283  for (int i=2; i>=0; --i) {
284  if (maxDr_par2[i]<maxDr_par[i]) maxDr_par[i] = maxDr_par2[i];
285  if (maxDz_par2[i]<maxDz_par[i]) maxDz_par[i] = maxDz_par2[i];
286  }
287  }
288 
289  Point bestVertex = getBestVertex(trk,vertices,minNVtxTrk); // min number of tracks 3 @HLT
290  if (bestVertex.z() < -99998.) {
291  bestVertex = beamSpot.position();
292  }
293 
294  ret = std::min(ret,cut(dz(trk,bestVertex), maxDz_par,std::less<float>()));
295  if (ret==-1.f) return ret;
296  ret = std::min(ret,cut(dr(trk,bestVertex), maxDr_par,std::less<float>()));
297  if (ret==-1.f) return ret;
298 
299 
300  }
301  if (ret==-1.f) return ret;
302 
303  return ret;
304 
305  }
306 
307 
308 
309  static const char * name() { return "TrackCutClassifier";}
310 
311  static void fillDescriptions(edm::ParameterSetDescription & desc) {
312  desc.add<bool>("isHLT",false);
313  desc.add<std::vector<int>>("minHits4pass", { std::numeric_limits<int>::max(), std::numeric_limits<int>::max(), std::numeric_limits<int>::max() } );
314  desc.add<std::vector<int>>("minHits", { 0, 0, 1});
315  desc.add<std::vector<int>>("minPixelHits", { 0, 0, 1});
316  desc.add<std::vector<int>>("minLayers", { 3, 4, 5});
317  desc.add<std::vector<int>>("min3DLayers", { 1, 2, 3});
318  desc.add<std::vector<int>>("maxLostLayers",{99, 3, 3});
319  desc.add<std::vector<double>>("maxRelPtErr", { std::numeric_limits<float>::max(),std::numeric_limits<float>::max(),std::numeric_limits<float>::max() } );
320  desc.add<std::vector<double>>("minNdof", {-1., -1., -1.});
321  desc.add<std::vector<double>>("maxChi2", {9999.,25., 16. });
322  desc.add<std::vector<double>>("maxChi2n", {9999., 1.0, 0.4});
323 
324  desc.add<int>("minNVtxTrk", 2);
325 
327  desc.add<std::vector<double>>("maxDzWrtBS",{std::numeric_limits<float>::max(),24.,15.});
329 
331  dz_par.add<std::vector<int>> ("dz_exp", {std::numeric_limits<int>::max(), std::numeric_limits<int>::max(), std::numeric_limits<int>::max()} ); // par = 4
332  dz_par.add<std::vector<double>>("dz_par1", {std::numeric_limits<float>::max(),std::numeric_limits<float>::max(),std::numeric_limits<float>::max()}); // par = 0.4
333  dz_par.add<std::vector<double>>("dz_par2", {std::numeric_limits<float>::max(),std::numeric_limits<float>::max(),std::numeric_limits<float>::max()}); // par = 0.35
334  dz_par.add<std::vector<double>>("dzWPVerr_par",{std::numeric_limits<float>::max(),std::numeric_limits<float>::max(),std::numeric_limits<float>::max()}); // par = 3.
335  desc.add<edm::ParameterSetDescription>("dz_par", dz_par);
336 
338  dr_par.add<std::vector<int>> ("dr_exp", {std::numeric_limits<int>::max(), std::numeric_limits<int>::max(), std::numeric_limits<int>::max()} ); // par = 4
339  dr_par.add<std::vector<double>>("dr_par1",{std::numeric_limits<float>::max(),std::numeric_limits<float>::max(),std::numeric_limits<float>::max()}); // par = 0.4
340  dr_par.add<std::vector<double>>("dr_par2",{std::numeric_limits<float>::max(),std::numeric_limits<float>::max(),std::numeric_limits<float>::max()}); // par = 0.3
341  dr_par.add<std::vector<double>>("d0err", {0.003, 0.003, 0.003});
342  dr_par.add<std::vector<double>>("d0err_par", {0.001, 0.001, 0.001});
343  dr_par.add<std::vector<double>>("drWPVerr_par",{std::numeric_limits<float>::max(),std::numeric_limits<float>::max(),std::numeric_limits<float>::max()}); // par = 3.
344  desc.add<edm::ParameterSetDescription>("dr_par", dr_par);
345 
346  }
347 
348  bool isHLT;
349  float maxRelPtErr[3];
350  float minNdof[3];
351  float maxChi2[3];
352  float maxChi2n[3];
353  int minLayers[3];
354  int min3DLayers[3];
355  int minHits4pass[3];
356  int minHits[3];
357  int minPixelHits[3];
358  int maxLostLayers[3];
359  int minNVtxTrk;
360  float maxDz[3];
361  float maxDzWrtBS[3];
362  float maxDr[3];
363  int dz_exp[3];
364  float dz_par1[3];
365  float dz_par2[3];
366  float dzWPVerr_par[3];
367  int dr_exp[3];
368  float dr_par1[3];
369  float dr_par2[3];
370  float d0err[3];
371  float d0err_par[3];
372  float drWPVerr_par[3];
373  };
374 
375 
376  using TrackCutClassifier = TrackMVAClassifier<Cuts>;
377 
378 }
379 
382 
383 DEFINE_FWK_MODULE(TrackCutClassifier);
384 
double p() const
momentum vector magnitude
Definition: TrackBase.h:654
T getParameter(std::string const &) const
double d0Error() const
error on d0
Definition: TrackBase.h:853
def copy(args, dbName)
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:600
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:538
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:557
Point getBestVertex_withError(reco::Track const &trk, reco::VertexCollection const &vertices, Point &error, const size_t minNtracks=2)
Definition: getBestVertex.h:30
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:379
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:594
T sqrt(T t)
Definition: SSEVec.h:18
virtual int dimension() const =0
double pt() const
track transverse momentum
Definition: TrackBase.h:660
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:814
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:901
#define end
Definition: vmac.h:39
T min(T a, T b)
Definition: MathUtil.h:58
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:109
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:648
double dzError() const
error on dz
Definition: TrackBase.h:865
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:485
et
define resolution functions of each parameter
int trackerLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:574
#define begin
Definition: vmac.h:32
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:55
Point getBestVertex(reco::Track const &trk, reco::VertexCollection const &vertices, const size_t minNtracks=2)
Definition: getBestVertex.h:9
int numberOfValidPixelHits() const
Definition: HitPattern.h:916
bool isFromDetOrFast(TrackingRecHit const &hit)
const Point & position() const
position
Definition: BeamSpot.h:62
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:630
long double T
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
#define constexpr
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:114