CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MultiTrackSelector.cc
Go to the documentation of this file.
2 
8 
9 //#include "RecoTracker/DebugTools/interface/FixTrackHitPattern.h"
10 #include <Math/DistFunc.h>
11 #include "TMath.h"
12 #include "TFile.h"
13 
15 
16 MultiTrackSelector::MultiTrackSelector()
17 {
18  useForestFromDB_ = true;
19  forest_ = 0;
20  gbrVals_ = 0;
21 }
22 
24  src_( consumes<reco::TrackCollection>( cfg.getParameter<edm::InputTag>( "src" ) ) ),
25  beamspot_( consumes<reco::BeamSpot>( cfg.getParameter<edm::InputTag>( "beamspot" ) ) ),
26  useVertices_( cfg.getParameter<bool>( "useVertices" ) ),
27  useVtxError_( cfg.getParameter<bool>( "useVtxError" ) )
28  // now get the pset for each selector
29 {
30  if (useVertices_) vertices_ = consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>( "vertices" ));
31 
32  useAnyMVA_ = false;
33  forestLabel_ = "MVASelectorIter0";
34  std::string type = "BDTG";
35  useForestFromDB_ = true;
36  dbFileName_ = "";
37 
38  forest_ = 0;
39  gbrVals_ = 0;
40 
41  if(cfg.exists("useAnyMVA")) useAnyMVA_ = cfg.getParameter<bool>("useAnyMVA");
42 
43  if(useAnyMVA_){
44  gbrVals_ = new float[11];
45  if(cfg.exists("mvaType"))type = cfg.getParameter<std::string>("mvaType");
46  if(cfg.exists("GBRForestLabel"))forestLabel_ = cfg.getParameter<std::string>("GBRForestLabel");
47  if(cfg.exists("GBRForestFileName")){
48  dbFileName_ = cfg.getParameter<std::string>("GBRForestFileName");
49  useForestFromDB_ = false;
50  }
51 
52  mvaType_ = type;
53  }
54  std::vector<edm::ParameterSet> trkSelectors( cfg.getParameter<std::vector< edm::ParameterSet> >("trackSelectors") );
55  qualityToSet_.reserve(trkSelectors.size());
56  vtxNumber_.reserve(trkSelectors.size());
57  vertexCut_.reserve(trkSelectors.size());
58  res_par_.reserve(trkSelectors.size());
59  chi2n_par_.reserve(trkSelectors.size());
60  chi2n_no1Dmod_par_.reserve(trkSelectors.size());
61  d0_par1_.reserve(trkSelectors.size());
62  dz_par1_.reserve(trkSelectors.size());
63  d0_par2_.reserve(trkSelectors.size());
64  dz_par2_.reserve(trkSelectors.size());
65  applyAdaptedPVCuts_.reserve(trkSelectors.size());
66  max_d0_.reserve(trkSelectors.size());
67  max_z0_.reserve(trkSelectors.size());
68  nSigmaZ_.reserve(trkSelectors.size());
69  min_layers_.reserve(trkSelectors.size());
70  min_3Dlayers_.reserve(trkSelectors.size());
71  max_lostLayers_.reserve(trkSelectors.size());
72  min_hits_bypass_.reserve(trkSelectors.size());
73  applyAbsCutsIfNoPV_.reserve(trkSelectors.size());
74  max_d0NoPV_.reserve(trkSelectors.size());
75  max_z0NoPV_.reserve(trkSelectors.size());
76  preFilter_.reserve(trkSelectors.size());
77  max_relpterr_.reserve(trkSelectors.size());
78  min_nhits_.reserve(trkSelectors.size());
79  max_minMissHitOutOrIn_.reserve(trkSelectors.size());
80  max_lostHitFraction_.reserve(trkSelectors.size());
81  min_eta_.reserve(trkSelectors.size());
82  max_eta_.reserve(trkSelectors.size());
83  useMVA_.reserve(trkSelectors.size());
84  //mvaReaders_.reserve(trkSelectors.size());
85  min_MVA_.reserve(trkSelectors.size());
86  //mvaType_.reserve(trkSelectors.size());
87 
88  produces<edm::ValueMap<float> >("MVAVals");
89 
90  for ( unsigned int i=0; i<trkSelectors.size(); i++) {
91 
93  // parameters for vertex selection
94  vtxNumber_.push_back( useVertices_ ? trkSelectors[i].getParameter<int32_t>("vtxNumber") : 0 );
95  vertexCut_.push_back( useVertices_ ? trkSelectors[i].getParameter<std::string>("vertexCut") : 0);
96  // parameters for adapted optimal cuts on chi2 and primary vertex compatibility
97  res_par_.push_back(trkSelectors[i].getParameter< std::vector<double> >("res_par") );
98  chi2n_par_.push_back( trkSelectors[i].getParameter<double>("chi2n_par") );
99  chi2n_no1Dmod_par_.push_back( trkSelectors[i].getParameter<double>("chi2n_no1Dmod_par") );
100  d0_par1_.push_back(trkSelectors[i].getParameter< std::vector<double> >("d0_par1"));
101  dz_par1_.push_back(trkSelectors[i].getParameter< std::vector<double> >("dz_par1"));
102  d0_par2_.push_back(trkSelectors[i].getParameter< std::vector<double> >("d0_par2"));
103  dz_par2_.push_back(trkSelectors[i].getParameter< std::vector<double> >("dz_par2"));
104  // Boolean indicating if adapted primary vertex compatibility cuts are to be applied.
105  applyAdaptedPVCuts_.push_back(trkSelectors[i].getParameter<bool>("applyAdaptedPVCuts"));
106  // Impact parameter absolute cuts.
107  max_d0_.push_back(trkSelectors[i].getParameter<double>("max_d0"));
108  max_z0_.push_back(trkSelectors[i].getParameter<double>("max_z0"));
109  nSigmaZ_.push_back(trkSelectors[i].getParameter<double>("nSigmaZ"));
110  // Cuts on numbers of layers with hits/3D hits/lost hits.
111  min_layers_.push_back(trkSelectors[i].getParameter<uint32_t>("minNumberLayers") );
112  min_3Dlayers_.push_back(trkSelectors[i].getParameter<uint32_t>("minNumber3DLayers") );
113  max_lostLayers_.push_back(trkSelectors[i].getParameter<uint32_t>("maxNumberLostLayers"));
114  min_hits_bypass_.push_back(trkSelectors[i].getParameter<uint32_t>("minHitsToBypassChecks"));
115  // Flag to apply absolute cuts if no PV passes the selection
116  applyAbsCutsIfNoPV_.push_back(trkSelectors[i].getParameter<bool>("applyAbsCutsIfNoPV"));
117  keepAllTracks_.push_back( trkSelectors[i].getParameter<bool>("keepAllTracks"));
118  max_relpterr_.push_back(trkSelectors[i].getParameter<double>("max_relpterr"));
119  min_nhits_.push_back(trkSelectors[i].getParameter<uint32_t>("min_nhits"));
120  max_minMissHitOutOrIn_.push_back(
121  trkSelectors[i].existsAs<int32_t>("max_minMissHitOutOrIn") ?
122  trkSelectors[i].getParameter<int32_t>("max_minMissHitOutOrIn") : 99);
123  max_lostHitFraction_.push_back(
124  trkSelectors[i].existsAs<double>("max_lostHitFraction") ?
125  trkSelectors[i].getParameter<double>("max_lostHitFraction") : 1.0);
126  min_eta_.push_back(trkSelectors[i].existsAs<double>("min_eta") ?
127  trkSelectors[i].getParameter<double>("min_eta"):-9999);
128  max_eta_.push_back(trkSelectors[i].existsAs<double>("max_eta") ?
129  trkSelectors[i].getParameter<double>("max_eta"):9999);
130 
131  setQualityBit_.push_back( false );
132  std::string qualityStr = trkSelectors[i].getParameter<std::string>("qualityBit");
133  if (qualityStr != "") {
134  setQualityBit_[i] = true;
135  qualityToSet_[i] = TrackBase::qualityByName(trkSelectors[i].getParameter<std::string>("qualityBit"));
136  }
137 
138  if (setQualityBit_[i] && (qualityToSet_[i] == TrackBase::undefQuality)) throw cms::Exception("Configuration") <<
139  "You can't set the quality bit " << trkSelectors[i].getParameter<std::string>("qualityBit") << " as it is 'undefQuality' or unknown.\n";
140 
141  if (applyAbsCutsIfNoPV_[i]) {
142  max_d0NoPV_.push_back(trkSelectors[i].getParameter<double>("max_d0NoPV"));
143  max_z0NoPV_.push_back(trkSelectors[i].getParameter<double>("max_z0NoPV"));
144  }
145  else{//dummy values
146  max_d0NoPV_.push_back(0.);
147  max_z0NoPV_.push_back(0.);
148  }
149 
150  name_.push_back( trkSelectors[i].getParameter<std::string>("name") );
151 
152  preFilter_[i]=trkSelectors.size(); // no prefilter
153 
154  std::string pfName=trkSelectors[i].getParameter<std::string>("preFilterName");
155  if (pfName!="") {
156  bool foundPF=false;
157  for ( unsigned int j=0; j<i; j++)
158  if (name_[j]==pfName ) {
159  foundPF=true;
160  preFilter_[i]=j;
161  }
162  if ( !foundPF)
163  throw cms::Exception("Configuration") << "Invalid prefilter name in MultiTrackSelector "
164  << trkSelectors[i].getParameter<std::string>("preFilterName");
165 
166  }
167 
168  // produces<std::vector<int> >(name_[i]).setBranchAlias( name_[i] + "TrackQuals");
169  produces<edm::ValueMap<int> >(name_[i]).setBranchAlias( name_[i] + "TrackQuals");
170  if(useAnyMVA_){
171  bool thisMVA = false;
172  if(trkSelectors[i].exists("useMVA"))thisMVA = trkSelectors[i].getParameter<bool>("useMVA");
173  useMVA_.push_back(thisMVA);
174  if(thisMVA){
175  double minVal = -1;
176  if(trkSelectors[i].exists("minMVA"))minVal = trkSelectors[i].getParameter<double>("minMVA");
177  min_MVA_.push_back(minVal);
178 
179  }else{
180  min_MVA_.push_back(-9999.0);
181  }
182  }else{
183  min_MVA_.push_back(-9999.0);
184  }
185 
186  }
187 }
188 
190  if(gbrVals_)delete [] gbrVals_;
191  if(!useForestFromDB_ && forest_)delete forest_;
192 }
193 
195 {
196  using namespace std;
197  using namespace edm;
198  using namespace reco;
199 
200  // Get tracks
201  Handle<TrackCollection> hSrcTrack;
202  evt.getByToken( src_, hSrcTrack );
203  const TrackCollection& srcTracks(*hSrcTrack);
204 
205  // looking for the beam spot
207  evt.getByToken(beamspot_, hBsp);
208  const reco::BeamSpot& vertexBeamSpot(*hBsp);
209 
210 
211  // Select good primary vertices for use in subsequent track selection
213  if (useVertices_) evt.getByToken(vertices_, hVtx);
214 
215  unsigned int trkSize=srcTracks.size();
216  std::vector<int> selTracksSave( qualityToSet_.size()*trkSize,0);
217 
218  processMVA(evt,es);
219 
220  for (unsigned int i=0; i<qualityToSet_.size(); i++) {
221  std::vector<int> selTracks(trkSize,0);
222  auto_ptr<edm::ValueMap<int> > selTracksValueMap = auto_ptr<edm::ValueMap<int> >(new edm::ValueMap<int>);
223  edm::ValueMap<int>::Filler filler(*selTracksValueMap);
224 
225  std::vector<Point> points;
226  std::vector<float> vterr, vzerr;
227  if (useVertices_) selectVertices(i,*hVtx, points, vterr, vzerr);
228 
229  // Loop over tracks
230  size_t current = 0;
231  for (TrackCollection::const_iterator it = srcTracks.begin(), ed = srcTracks.end(); it != ed; ++it, ++current) {
232  const Track & trk = * it;
233  // Check if this track passes cuts
234 
235  LogTrace("TrackSelection") << "ready to check track with pt="<< trk.pt() ;
236 
237  //already removed
238  bool ok=true;
239  if (preFilter_[i]<i && selTracksSave[preFilter_[i]*trkSize+current] < 0) {
240  selTracks[current]=-1;
241  ok=false;
242  if ( !keepAllTracks_[i])
243  continue;
244  }
245  else {
246  double mvaVal = 0;
247  if(useAnyMVA_)mvaVal = mvaVals_[current];
248  ok = select(i,vertexBeamSpot, trk, points, vterr, vzerr,mvaVal);
249  if (!ok) {
250  LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " NOT selected";
251  if (!keepAllTracks_[i]) {
252  selTracks[current]=-1;
253  continue;
254  }
255  }
256  else
257  LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " selected";
258  }
259 
260  if (preFilter_[i]<i ) {
261  selTracks[current]=selTracksSave[preFilter_[i]*trkSize+current];
262  }
263  else {
264  selTracks[current]=trk.qualityMask();
265  }
266  if ( ok && setQualityBit_[i]) {
267  selTracks[current]= (selTracks[current] | (1<<qualityToSet_[i]));
269  selTracks[current]=(selTracks[current] | (1<<TrackBase::loose));
270  }
271  else if (qualityToSet_[i]==TrackBase::highPurity) {
272  selTracks[current]=(selTracks[current] | (1<<TrackBase::loose));
273  selTracks[current]=(selTracks[current] | (1<<TrackBase::tight));
274  }
275  if (!points.empty()) {
277  selTracks[current]=(selTracks[current] | (1<<TrackBase::looseSetWithPV));
278  }
279  else if (qualityToSet_[i]==TrackBase::highPurity) {
280  selTracks[current]=(selTracks[current] | (1<<TrackBase::looseSetWithPV));
281  selTracks[current]=(selTracks[current] | (1<<TrackBase::highPuritySetWithPV));
282  }
283  }
284  }
285  }
286  for ( unsigned int j=0; j< trkSize; j++ ) selTracksSave[j+i*trkSize]=selTracks[j];
287  filler.insert(hSrcTrack, selTracks.begin(),selTracks.end());
288  filler.fill();
289 
290  // evt.put(selTracks,name_[i]);
291  evt.put(selTracksValueMap,name_[i]);
292  }
293 }
294 
295 
296  bool MultiTrackSelector::select(unsigned int tsNum,
297  const reco::BeamSpot &vertexBeamSpot,
298  const reco::Track &tk,
299  const std::vector<Point> &points,
300  std::vector<float> &vterr,
301  std::vector<float> &vzerr,
302  double mvaVal) {
303  // Decide if the given track passes selection cuts.
304 
305  using namespace std;
306 
307  if(tk.found()>=min_hits_bypass_[tsNum]) return true;
308  if ( tk.ndof() < 1E-5 ) return false;
309 
310  // Cuts on numbers of layers with hits/3D hits/lost hits.
311  uint32_t nlayers = tk.hitPattern().trackerLayersWithMeasurement();
312  uint32_t nlayers3D = tk.hitPattern().pixelLayersWithMeasurement() +
314  uint32_t nlayersLost = tk.hitPattern().trackerLayersWithoutMeasurement();
315  LogDebug("TrackSelection") << "cuts on nlayers: " << nlayers << " " << nlayers3D << " " << nlayersLost << " vs "
316  << min_layers_[tsNum] << " " << min_3Dlayers_[tsNum] << " " << max_lostLayers_[tsNum];
317  if (nlayers < min_layers_[tsNum]) return false;
318  if (nlayers3D < min_3Dlayers_[tsNum]) return false;
319  if (nlayersLost > max_lostLayers_[tsNum]) return false;
320  LogTrace("TrackSelection") << "cuts on nlayers passed";
321 
322  float chi2n = tk.normalizedChi2();
323  float chi2n_no1Dmod = chi2n;
324 
325  int count1dhits = 0;
326  for (trackingRecHit_iterator ith = tk.recHitsBegin(), edh = tk.recHitsEnd(); ith != edh; ++ith) {
327  const TrackingRecHit * hit = ith->get();
328  if (hit->isValid()) {
329  if (typeid(*hit) == typeid(SiStripRecHit1D)) ++count1dhits;
330  }
331  }
332  if (count1dhits > 0) {
333  float chi2 = tk.chi2();
334  float ndof = tk.ndof();
335  chi2n = (chi2+count1dhits)/float(ndof+count1dhits);
336  }
337  // For each 1D rechit, the chi^2 and ndof is increased by one. This is a way of retaining approximately
338  // the same normalized chi^2 distribution as with 2D rechits.
339  if (chi2n > chi2n_par_[tsNum]*nlayers) return false;
340 
341  if (chi2n_no1Dmod > chi2n_no1Dmod_par_[tsNum]*nlayers) return false;
342 
343  // Get track parameters
344  float pt = std::max(float(tk.pt()),0.000001f);
345  float eta = tk.eta();
346  if (eta<min_eta_[tsNum] || eta>max_eta_[tsNum]) return false;
347 
348  //cuts on relative error on pt and number of valid hits
349  float relpterr = float(tk.ptError())/pt;
350  uint32_t nhits = tk.numberOfValidHits();
351  if(relpterr > max_relpterr_[tsNum]) return false;
352  if(nhits < min_nhits_[tsNum]) return false;
353 
356  int minLost = std::min(lostIn,lostOut);
357  if (minLost > max_minMissHitOutOrIn_[tsNum]) return false;
358  float lostMidFrac = tk.numberOfLostHits() / (tk.numberOfValidHits() + tk.numberOfLostHits());
359  if (lostMidFrac > max_lostHitFraction_[tsNum]) return false;
360 
361 
362 
364  //Adding the MVA selection before vertex cuts//
366 
367  if(useAnyMVA_ && useMVA_[tsNum]){
368  if(mvaVal < min_MVA_[tsNum])return false;
369  }
370 
372  //End of MVA selection section//
374 
375  //other track parameters
376  float d0 = -tk.dxy(vertexBeamSpot.position()), d0E = tk.d0Error(),
377  dz = tk.dz(vertexBeamSpot.position()), dzE = tk.dzError();
378 
379  // parametrized d0 resolution for the track pt
380  float nomd0E = sqrt(res_par_[tsNum][0]*res_par_[tsNum][0]+(res_par_[tsNum][1]/pt)*(res_par_[tsNum][1]/pt));
381  // parametrized z0 resolution for the track pt and eta
382  float nomdzE = nomd0E*(std::cosh(eta));
383 
384  float dzCut = min( pow(dz_par1_[tsNum][0]*nlayers,dz_par1_[tsNum][1])*nomdzE,
385  pow(dz_par2_[tsNum][0]*nlayers,dz_par2_[tsNum][1])*dzE );
386  float d0Cut = min( pow(d0_par1_[tsNum][0]*nlayers,d0_par1_[tsNum][1])*nomd0E,
387  pow(d0_par2_[tsNum][0]*nlayers,d0_par2_[tsNum][1])*d0E );
388 
389 
390  // ---- PrimaryVertex compatibility cut
391  bool primaryVertexZCompatibility(false);
392  bool primaryVertexD0Compatibility(false);
393 
394  if (points.empty()) { //If not primaryVertices are reconstructed, check just the compatibility with the BS
395  //z0 within (n sigma + dzCut) of the beam spot z, if no good vertex is found
396  if ( abs(dz) < hypot(vertexBeamSpot.sigmaZ()*nSigmaZ_[tsNum],dzCut) ) primaryVertexZCompatibility = true;
397  // d0 compatibility with beam line
398  if (abs(d0) < d0Cut) primaryVertexD0Compatibility = true;
399  }
400 
401  int iv=0;
402  for (std::vector<Point>::const_iterator point = points.begin(), end = points.end(); point != end; ++point) {
403  LogTrace("TrackSelection") << "Test track w.r.t. vertex with z position " << point->z();
404  if(primaryVertexZCompatibility && primaryVertexD0Compatibility) break;
405  float dzPV = tk.dz(*point); //re-evaluate the dz with respect to the vertex position
406  float d0PV = tk.dxy(*point); //re-evaluate the dxy with respect to the vertex position
407  if(useVtxError_){
408  float dzErrPV = sqrt(dzE*dzE+vzerr[iv]*vzerr[iv]); // include vertex error in z
409  float d0ErrPV = sqrt(d0E*d0E+vterr[iv]*vterr[iv]); // include vertex error in xy
410  iv++;
411  if (abs(dzPV) < dz_par1_[tsNum][0]*pow(nlayers,dz_par1_[tsNum][1])*nomdzE &&
412  abs(dzPV) < dz_par2_[tsNum][0]*pow(nlayers,dz_par2_[tsNum][1])*dzErrPV &&
413  abs(dzPV) < max_z0_[tsNum]) primaryVertexZCompatibility = true;
414  if (abs(d0PV) < d0_par1_[tsNum][0]*pow(nlayers,d0_par1_[tsNum][1])*nomd0E &&
415  abs(d0PV) < d0_par2_[tsNum][0]*pow(nlayers,d0_par2_[tsNum][1])*d0ErrPV &&
416  abs(d0PV) < max_d0_[tsNum]) primaryVertexD0Compatibility = true;
417  }else{
418  if (abs(dzPV) < dzCut) primaryVertexZCompatibility = true;
419  if (abs(d0PV) < d0Cut) primaryVertexD0Compatibility = true;
420  }
421  LogTrace("TrackSelection") << "distances " << dzPV << " " << d0PV << " vs " << dzCut << " " << d0Cut;
422  }
423 
424  if (points.empty() && applyAbsCutsIfNoPV_[tsNum]) {
425  if ( abs(dz) > max_z0NoPV_[tsNum] || abs(d0) > max_d0NoPV_[tsNum]) return false;
426  } else {
427  // Absolute cuts on all tracks impact parameters with respect to beam-spot.
428  // If BS is not compatible, verify if at least the reco-vertex is compatible (useful for incorrect BS settings)
429  if (abs(d0) > max_d0_[tsNum] && !primaryVertexD0Compatibility) return false;
430  LogTrace("TrackSelection") << "absolute cuts on d0 passed";
431  if (abs(dz) > max_z0_[tsNum] && !primaryVertexZCompatibility) return false;
432  LogTrace("TrackSelection") << "absolute cuts on dz passed";
433  }
434 
435  LogTrace("TrackSelection") << "cuts on PV: apply adapted PV cuts? " << applyAdaptedPVCuts_[tsNum]
436  << " d0 compatibility? " << primaryVertexD0Compatibility
437  << " z compatibility? " << primaryVertexZCompatibility ;
438 
439  if (applyAdaptedPVCuts_[tsNum]) {
440  return (primaryVertexD0Compatibility && primaryVertexZCompatibility);
441  } else {
442  return true;
443  }
444 
445 }
446 
447  void MultiTrackSelector::selectVertices(unsigned int tsNum,
448  const reco::VertexCollection &vtxs,
449  std::vector<Point> &points,
450  std::vector<float> &vterr,
451  std::vector<float> &vzerr) {
452  // Select good primary vertices
453  using namespace reco;
454  int32_t toTake = vtxNumber_[tsNum];
455  for (VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end(); it != ed; ++it) {
456 
457  LogDebug("SelectVertex") << " select vertex with z position " << it->z() << " "
458  << it->chi2() << " " << it->ndof() << " " << TMath::Prob(it->chi2(), static_cast<int32_t>(it->ndof()));
459  Vertex vtx = *it;
460  bool pass = vertexCut_[tsNum]( vtx );
461  if( pass ) {
462  points.push_back(it->position());
463  vterr.push_back(sqrt(it->yError()*it->xError()));
464  vzerr.push_back(it->zError());
465  LogTrace("SelectVertex") << " SELECTED vertex with z position " << it->z();
466  toTake--; if (toTake == 0) break;
467  }
468  }
469 }
470 
472 {
473 
474  using namespace std;
475  using namespace edm;
476  using namespace reco;
477 
478  // Get tracks
479  Handle<TrackCollection> hSrcTrack;
480  evt.getByToken( src_, hSrcTrack );
481  const TrackCollection& srcTracks(*hSrcTrack);
482 
483  auto_ptr<edm::ValueMap<float> >mvaValValueMap = auto_ptr<edm::ValueMap<float> >(new edm::ValueMap<float>);
484  edm::ValueMap<float>::Filler mvaFiller(*mvaValValueMap);
485 
486  mvaVals_.clear();
487 
488  if(!useAnyMVA_){
489  size_t current = 0;
490  for (TrackCollection::const_iterator it = srcTracks.begin(), ed = srcTracks.end(); it != ed; ++it, ++current) {
491  mvaVals_.push_back(-99.0);
492  }
493  mvaFiller.insert(hSrcTrack,mvaVals_.begin(),mvaVals_.end());
494  mvaFiller.fill();
495  evt.put(mvaValValueMap,"MVAVals");
496  return;
497  }
498 
499  if(!forest_){
500  if(useForestFromDB_){
501  edm::ESHandle<GBRForest> forestHandle;
502  es.get<GBRWrapperRcd>().get(forestLabel_,forestHandle);
503  forest_ = (GBRForest*)forestHandle.product();
504  }else{
505  TFile gbrfile(dbFileName_.c_str());
506  forest_ = (GBRForest*)gbrfile.Get(forestLabel_.c_str());
507  }
508  }
509 
510 
511 
512  size_t current = 0;
513  for (TrackCollection::const_iterator it = srcTracks.begin(), ed = srcTracks.end(); it != ed; ++it, ++current) {
514  const Track & trk = * it;
515  tmva_ndof_ = trk.ndof();
519  float chi2n = trk.normalizedChi2();
520  float chi2n_no1Dmod = chi2n;
521 
522  int count1dhits = 0;
523  for (trackingRecHit_iterator ith = trk.recHitsBegin(), edh = trk.recHitsEnd(); ith != edh; ++ith) {
524  const TrackingRecHit * hit = ith->get();
525  if (hit->isValid()) {
526  if (typeid(*hit) == typeid(SiStripRecHit1D)) ++count1dhits;
527  }
528  }
529  if (count1dhits > 0) {
530  float chi2 = trk.chi2();
531  float ndof = trk.ndof();
532  chi2n = (chi2+count1dhits)/float(ndof+count1dhits);
533  }
534  tmva_chi2n_ = chi2n;
535  tmva_chi2n_no1dmod_ = chi2n_no1Dmod;
536  tmva_eta_ = trk.eta();
537  tmva_relpterr_ = float(trk.ptError())/std::max(float(trk.pt()),0.000001f);
541  int minLost = std::min(lostIn,lostOut); tmva_minlost_ = minLost;
543 
545  gbrVals_[1] = tmva_minlost_;
546  gbrVals_[2] = tmva_nhits_;
548  gbrVals_[4] = tmva_eta_;
550  gbrVals_[6] = tmva_chi2n_;
553  gbrVals_[9] = tmva_nlayers_;
554  gbrVals_[10] = tmva_ndof_;
555 
556  double gbrVal = forest_->GetClassifier(gbrVals_);
557  mvaVals_.push_back((float)gbrVal);
558  }
559  mvaFiller.insert(hSrcTrack,mvaVals_.begin(),mvaVals_.end());
560  mvaFiller.fill();
561  evt.put(mvaValValueMap,"MVAVals");
562 
563 }
#define LogDebug(id)
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
std::vector< double > max_d0_
Impact parameter absolute cuts.
int i
Definition: DBlmapReader.cc:9
double d0Error() const
error on d0
Definition: TrackBase.h:209
std::vector< std::vector< double > > dz_par1_
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:109
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
std::vector< std::vector< double > > d0_par1_
int numberOfLostTrackerHits() const
Definition: HitPattern.h:643
std::vector< std::vector< double > > d0_par2_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
bool exists(std::string const &parameterName) const
checks if a parameter exists
int pixelLayersWithMeasurement() const
Definition: HitPattern.h:738
int trackerLayersWithoutMeasurement() const
Definition: HitPattern.h:751
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:232
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< uint32_t > min_3Dlayers_
T eta() const
std::vector< int32_t > vtxNumber_
vertex cuts
std::vector< std::string > name_
const HitPattern & trackerExpectedHitsOuter() const
Access the hit pattern counting (in the Tracker) the number of expected crossed layers after the last...
Definition: TrackBase.h:225
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:139
const T & max(const T &a, const T &b)
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:105
std::vector< std::vector< double > > res_par_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:107
int trackerLayersWithMeasurement() const
Definition: HitPattern.h:733
T sqrt(T t)
Definition: SSEVec.h:48
void processMVA(edm::Event &evt, const edm::EventSetup &es)
double pt() const
track transverse momentum
Definition: TrackBase.h:129
std::vector< uint32_t > min_nhits_
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:192
int qualityMask() const
Definition: TrackBase.h:284
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
std::vector< unsigned int > preFilter_
double f[11][100]
const HitPattern & trackerExpectedHitsInner() const
Access the hit pattern counting (in the Tracker) the number of expected crossed layers before the fir...
Definition: TrackBase.h:223
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:230
#define end
Definition: vmac.h:37
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:62
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:221
#define LogTrace(id)
std::vector< uint32_t > min_hits_bypass_
void selectVertices(unsigned int tsNum, const reco::VertexCollection &vtxs, std::vector< Point > &points, std::vector< float > &vterr, std::vector< float > &vzerr)
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:125
double dzError() const
error on dz
Definition: TrackBase.h:213
virtual void produce(edm::Event &evt, const edm::EventSetup &es) override
process one event
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:46
std::vector< double > chi2n_no1Dmod_par_
std::vector< StringCutObjectSelector< reco::Vertex > > vertexCut_
const T & get() const
Definition: EventSetup.h:55
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
edm::EDGetTokenT< reco::BeamSpot > beamspot_
edm::EDGetTokenT< reco::TrackCollection > src_
source collection label
T const * product() const
Definition: ESHandle.h:62
bool isValid() const
bool select(unsigned tsNum, const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk, const std::vector< Point > &points, std::vector< float > &vterr, std::vector< float > &vzerr, double mvaVal)
return class, or -1 if rejected
std::vector< TrackBase::TrackQuality > qualityToSet_
std::vector< bool > setQualityBit_
do I have to set a quality bit?
edm::EDGetTokenT< reco::VertexCollection > vertices_
std::vector< uint32_t > min_layers_
Cuts on numbers of layers with hits/3D hits/lost hits.
unsigned short found() const
Number of valid hits on track.
Definition: Track.h:99
const Point & position() const
position
Definition: BeamSpot.h:62
std::vector< uint32_t > max_lostLayers_
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:119
double GetClassifier(const float *vector) const
Definition: GBRForest.h:60
std::vector< std::vector< double > > dz_par2_
int numberOfValidStripLayersWithMonoAndStereo(uint32_t stripdet, uint32_t layer) const
Definition: HitPattern.cc:223
std::vector< int32_t > max_minMissHitOutOrIn_
std::vector< int32_t > max_lostHitFraction_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:64