forked from cms-sw/cmssw
-
Notifications
You must be signed in to change notification settings - Fork 1
/
TrackCutClassifier.cc
384 lines (300 loc) · 14.8 KB
/
TrackCutClassifier.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
#include "RecoTracker/FinalTrackSelectors/interface/TrackMVAClassifier.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
#include <cassert>
#include "getBestVertex.h"
#include "powN.h"
namespace {
void fillArrayF(float * x,const edm::ParameterSet & cfg, const char * name) {
auto v = cfg.getParameter< std::vector<double> >(name);
assert(v.size()==3);
std::copy(std::begin(v),std::end(v),x);
}
void fillArrayI(int * x,const edm::ParameterSet & cfg, const char * name) {
auto v = cfg.getParameter< std::vector<int> >(name);
assert(v.size()==3);
std::copy(std::begin(v),std::end(v),x);
}
// fake mva value to return for loose,tight,hp
constexpr float mvaVal[3] = {-.5,.5,1.};
template<typename T,typename Comp>
inline float cut(T val, const T * cuts, Comp comp) {
for (int i=2; i>=0; --i) {
if ( comp(val,cuts[i]) ) return mvaVal[i];
}
return -1.f;
}
inline float chi2n(reco::Track const & tk) { return tk.normalizedChi2();}
inline float relPtErr(reco::Track const & tk) {
return (tk.pt() != 0. ? float(tk.ptError())/float(tk.pt()) : 9999999.);
}
inline int lostLayers(reco::Track const & tk) {
return tk.hitPattern().trackerLayersWithoutMeasurement(reco::HitPattern::TRACK_HITS);
}
inline int n3DLayers(reco::Track const & tk, bool isHLT) {
uint32_t nlayers3D = tk.hitPattern().pixelLayersWithMeasurement();
if (!isHLT)
nlayers3D += tk.hitPattern().numberOfValidStripLayersWithMonoAndStereo();
else {
size_t count3D = 0;
for ( auto it = tk.recHitsBegin(), et = tk.recHitsEnd(); it!=et; ++it) {
const TrackingRecHit* hit = (*it);
if ( trackerHitRTTI::isUndef(*hit) ) continue;
if ( hit->dimension()==2 ) {
auto const & thit = static_cast<BaseTrackerRecHit const&>(*hit);
if (thit.isMatched()) count3D++;
}
}
nlayers3D += count3D;
}
return nlayers3D;
}
inline int nHits(reco::Track const & tk) {
return tk.numberOfValidHits();
}
inline int nPixelHits(reco::Track const & tk) {
return tk.hitPattern().numberOfValidPixelHits();
}
inline float dz(reco::Track const & trk, Point const & bestVertex) {
return std::abs(trk.dz(bestVertex));
}
inline float dr(reco::Track const & trk, Point const & bestVertex) {
return std::abs(trk.dxy(bestVertex));
}
inline void dzCut_par1(reco::Track const & trk, int & nLayers, const float * par, const int * exp, float dzCut[]) {
float dzE = trk.dzError();
for (int i=2; i>=0; --i) {
dzCut[i] = powN(par[i]*nLayers,exp[i])*dzE;
}
}
inline void drCut_par1(reco::Track const & trk, int & nLayers, const float * par, const int * exp, float drCut[]) {
float drE = trk.d0Error();
for (int i=2; i>=0; --i) {
drCut[i] = powN(par[i]*nLayers,exp[i])*drE;
}
}
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[]) {
float pt = float(trk.pt());
float p = float(trk.p());
for (int i=2; i>=0; --i) {
// parametrized d0 resolution for the track pt
float nomd0E = sqrt(d0err[i]*d0err[i]+(d0err_par[i]/pt)*(d0err_par[i]/pt));
// parametrized z0 resolution for the track pt and eta
float nomdzE = nomd0E*(p/pt); // cosh(eta):=abs(p)/pt
dzCut[i] = powN(par[i]*nLayers,exp[i])*nomdzE;
}
}
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[]) {
float pt = trk.pt();
for (int i=2; i>=0; --i) {
// parametrized d0 resolution for the track pt
float nomd0E = sqrt(d0err[i]*d0err[i]+(d0err_par[i]/pt)*(d0err_par[i]/pt));
drCut[i] = powN(par[i]*nLayers,exp[i])*nomd0E;
}
}
inline void dzCut_wPVerror_par(reco::Track const & trk, int & nLayers, const float * par, const int * exp, Point const & bestVertexError, float dzCut[]) {
float dzE = trk.dzError();
float zPVerr = bestVertexError.z();
float dzErrPV = std::sqrt(dzE*dzE+zPVerr*zPVerr);
for (int i=2; i>=0; --i) {
dzCut[i] = par[i]*dzErrPV;
if (exp[i] != 0)
dzCut[i] *= pow(nLayers,exp[i]);
}
}
inline void drCut_wPVerror_par(reco::Track const & trk, int & nLayers, const float* par, const int * exp, Point const & bestVertexError, float drCut[]) {
float drE = trk.d0Error();
float rPVerr = sqrt(bestVertexError.x()*bestVertexError.y()); // shouldn't it be bestVertex.xError()*bestVertex.xError()+bestVertex.yError()*bestVertex.yError() ?!?!?
float drErrPV = std::sqrt(drE*drE+rPVerr*rPVerr);
for (int i=2; i>=0; --i) {
drCut[i] = par[i]*drErrPV;
if (exp[i] != 0)
drCut[i] *= pow(nLayers,exp[i]);
}
}
struct Cuts {
Cuts(const edm::ParameterSet & cfg) {
isHLT = cfg.getParameter<bool>("isHLT");
fillArrayF(minNdof, cfg,"minNdof");
fillArrayF(maxChi2, cfg,"maxChi2");
fillArrayF(maxChi2n, cfg,"maxChi2n");
fillArrayI(minHits4pass, cfg,"minHits4pass");
fillArrayI(minHits, cfg,"minHits");
fillArrayI(minPixelHits, cfg,"minPixelHits");
fillArrayI(min3DLayers, cfg,"min3DLayers");
fillArrayI(minLayers, cfg,"minLayers");
fillArrayI(maxLostLayers,cfg,"maxLostLayers");
fillArrayF(maxRelPtErr, cfg,"maxRelPtErr");
minNVtxTrk = cfg.getParameter<int>("minNVtxTrk");
fillArrayF(maxDz, cfg,"maxDz");
fillArrayF(maxDzWrtBS, cfg,"maxDzWrtBS");
fillArrayF(maxDr, cfg,"maxDr");
edm::ParameterSet dz_par = cfg.getParameter<edm::ParameterSet>("dz_par");
fillArrayI(dz_exp, dz_par,"dz_exp");
fillArrayF(dz_par1, dz_par,"dz_par1");
fillArrayF(dz_par2, dz_par,"dz_par2");
fillArrayF(dzWPVerr_par, dz_par,"dzWPVerr_par");
edm::ParameterSet dr_par = cfg.getParameter<edm::ParameterSet>("dr_par");
fillArrayI(dr_exp, dr_par,"dr_exp");
fillArrayF(dr_par1, dr_par,"dr_par1");
fillArrayF(dr_par2, dr_par,"dr_par2");
fillArrayF(d0err, dr_par,"d0err");
fillArrayF(d0err_par, dr_par,"d0err_par");
fillArrayF(drWPVerr_par, dr_par,"drWPVerr_par");
}
void beginStream() {}
void initEvent(const edm::EventSetup&) {}
float operator()(reco::Track const & trk,
reco::BeamSpot const & beamSpot,
reco::VertexCollection const & vertices) const {
float ret = 1.f;
// minimum number of hits for by-passing the other checks
if ( minHits4pass[0] < std::numeric_limits<int>::max() ) {
ret = std::min(ret,cut(nHits(trk),minHits4pass,std::greater_equal<int>()));
if (ret==1.f) return ret;
}
if ( maxRelPtErr[2] < std::numeric_limits<float>::max() ) {
ret = std::min(ret,cut(relPtErr(trk),maxRelPtErr,std::less_equal<float>()) );
if (ret==-1.f) return ret;
}
ret = std::min(ret,cut(float(trk.ndof()),minNdof,std::greater_equal<float>()) );
if (ret==-1.f) return ret;
auto nLayers = trk.hitPattern().trackerLayersWithMeasurement();
ret = std::min(ret,cut(nLayers,minLayers,std::greater_equal<int>()));
if (ret==-1.f) return ret;
ret = std::min(ret,cut(chi2n(trk)/float(nLayers),maxChi2n,std::less_equal<float>()));
if (ret==-1.f) return ret;
ret = std::min(ret,cut(chi2n(trk),maxChi2,std::less_equal<float>()));
if (ret==-1.f) return ret;
ret = std::min(ret,cut(n3DLayers(trk,isHLT),min3DLayers,std::greater_equal<int>()));
if (ret==-1.f) return ret;
ret = std::min(ret,cut(nHits(trk),minHits,std::greater_equal<int>()));
if (ret==-1.f) return ret;
ret = std::min(ret,cut(nPixelHits(trk),minPixelHits,std::greater_equal<int>()));
if (ret==-1.f) return ret;
ret = std::min(ret,cut(lostLayers(trk),maxLostLayers,std::less_equal<int>()));
if (ret==-1.f) return ret;
// original dz and dr cut
if (maxDz[2]<std::numeric_limits<float>::max() || maxDr[2]<std::numeric_limits<float>::max()) {
// if not primaryVertices are reconstructed, check compatibility w.r.t. beam spot
Point bestVertex = getBestVertex(trk,vertices,minNVtxTrk); // min number of tracks [2 (=default) for offline, 3 for HLT]
float maxDzcut[3];
std::copy(std::begin(maxDz),std::end(maxDz),std::begin(maxDzcut));
if (bestVertex.z() < -99998.) {
bestVertex = beamSpot.position();
std::copy(std::begin(maxDzWrtBS),std::end(maxDzWrtBS),std::begin(maxDzcut));
}
ret = std::min(ret,cut(dr(trk,bestVertex), maxDr,std::less<float>()));
if (ret==-1.f) return ret;
ret = std::min(ret,cut(dz(trk,bestVertex), maxDzcut,std::less<float>()));
if (ret==-1.f) return ret;
}
// parametrized dz and dr cut by using PV error
if (dzWPVerr_par[2]<std::numeric_limits<float>::max() || drWPVerr_par[2]<std::numeric_limits<float>::max()) {
Point bestVertexError(-1.,-1.,-1.);
Point bestVertex = getBestVertex_withError(trk,vertices,bestVertexError,minNVtxTrk); // min number of tracks [2 (=default) for offline, 3 for HLT]
float maxDz_par[3];
float maxDr_par[3];
dzCut_wPVerror_par(trk,nLayers,dzWPVerr_par,dz_exp,bestVertexError, maxDz_par);
drCut_wPVerror_par(trk,nLayers,drWPVerr_par,dr_exp,bestVertexError, maxDr_par);
ret = std::min(ret,cut(dr(trk,bestVertex), maxDr_par,std::less<float>()));
if (ret==-1.f) return ret;
ret = std::min(ret,cut(dz(trk,bestVertex), maxDr_par,std::less<float>()));
if (ret==-1.f) return ret;
}
// parametrized dz and dr cut by using their error
if (dz_par1[2]<std::numeric_limits<float>::max() || dr_par1[2]<std::numeric_limits<float>::max()) {
float maxDz_par1[3];
float maxDr_par1[3];
dzCut_par1(trk,nLayers,dz_par1,dz_exp, maxDz_par1);
drCut_par1(trk,nLayers,dr_par1,dr_exp, maxDr_par1);
float maxDz_par[3];
float maxDr_par[3];
std::copy(std::begin(maxDz_par1),std::end(maxDz_par1),std::begin(maxDz_par));
std::copy(std::begin(maxDr_par1),std::end(maxDr_par1),std::begin(maxDr_par));
// parametrized dz and dr cut by using d0 and z0 resolution
if (dz_par2[2]<std::numeric_limits<float>::max() || dr_par2[2]<std::numeric_limits<float>::max()) {
float maxDz_par2[3];
float maxDr_par2[3];
dzCut_par2(trk,nLayers,dz_par2,dz_exp,d0err,d0err_par, maxDz_par2);
drCut_par2(trk,nLayers,dr_par2,dr_exp,d0err,d0err_par, maxDr_par2);
for (int i=2; i>=0; --i) {
if (maxDr_par2[i]<maxDr_par[i]) maxDr_par[i] = maxDr_par2[i];
if (maxDz_par2[i]<maxDz_par[i]) maxDz_par[i] = maxDz_par2[i];
}
}
Point bestVertex = getBestVertex(trk,vertices,minNVtxTrk); // min number of tracks 3 @HLT
if (bestVertex.z() < -99998.) {
bestVertex = beamSpot.position();
}
ret = std::min(ret,cut(dz(trk,bestVertex), maxDz_par,std::less<float>()));
if (ret==-1.f) return ret;
ret = std::min(ret,cut(dr(trk,bestVertex), maxDr_par,std::less<float>()));
if (ret==-1.f) return ret;
}
if (ret==-1.f) return ret;
return ret;
}
static const char * name() { return "TrackCutClassifier";}
static void fillDescriptions(edm::ParameterSetDescription & desc) {
desc.add<bool>("isHLT",false);
desc.add<std::vector<int>>("minHits4pass", { std::numeric_limits<int>::max(), std::numeric_limits<int>::max(), std::numeric_limits<int>::max() } );
desc.add<std::vector<int>>("minHits", { 0, 0, 1});
desc.add<std::vector<int>>("minPixelHits", { 0, 0, 1});
desc.add<std::vector<int>>("minLayers", { 3, 4, 5});
desc.add<std::vector<int>>("min3DLayers", { 1, 2, 3});
desc.add<std::vector<int>>("maxLostLayers",{99, 3, 3});
desc.add<std::vector<double>>("maxRelPtErr", { std::numeric_limits<float>::max(),std::numeric_limits<float>::max(),std::numeric_limits<float>::max() } );
desc.add<std::vector<double>>("minNdof", {-1., -1., -1.});
desc.add<std::vector<double>>("maxChi2", {9999.,25., 16. });
desc.add<std::vector<double>>("maxChi2n", {9999., 1.0, 0.4});
desc.add<int>("minNVtxTrk", 2);
desc.add<std::vector<double>>("maxDz",{std::numeric_limits<float>::max(),std::numeric_limits<float>::max(),std::numeric_limits<float>::max()});
desc.add<std::vector<double>>("maxDzWrtBS",{std::numeric_limits<float>::max(),24.,15.});
desc.add<std::vector<double>>("maxDr",{std::numeric_limits<float>::max(),std::numeric_limits<float>::max(),std::numeric_limits<float>::max()});
edm::ParameterSetDescription dz_par;
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
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
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
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.
desc.add<edm::ParameterSetDescription>("dz_par", dz_par);
edm::ParameterSetDescription dr_par;
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
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
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
dr_par.add<std::vector<double>>("d0err", {0.003, 0.003, 0.003});
dr_par.add<std::vector<double>>("d0err_par", {0.001, 0.001, 0.001});
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.
desc.add<edm::ParameterSetDescription>("dr_par", dr_par);
}
bool isHLT;
float maxRelPtErr[3];
float minNdof[3];
float maxChi2[3];
float maxChi2n[3];
int minLayers[3];
int min3DLayers[3];
int minHits4pass[3];
int minHits[3];
int minPixelHits[3];
int maxLostLayers[3];
int minNVtxTrk;
float maxDz[3];
float maxDzWrtBS[3];
float maxDr[3];
int dz_exp[3];
float dz_par1[3];
float dz_par2[3];
float dzWPVerr_par[3];
int dr_exp[3];
float dr_par1[3];
float dr_par2[3];
float d0err[3];
float d0err_par[3];
float drWPVerr_par[3];
};
using TrackCutClassifier = TrackMVAClassifier<Cuts>;
}
#include "FWCore/PluginManager/interface/ModuleDef.h"
#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(TrackCutClassifier);