CMS 3D CMS Logo

DefaultTrackMVAClassifier.cc
Go to the documentation of this file.
2 
3 
6 #include <limits>
7 
8 #include "getBestVertex.h"
9 
10 namespace {
11 
12 template<bool PROMPT>
13 struct mva {
14  mva(const edm::ParameterSet &){}
15  float operator()(reco::Track const & trk,
16  reco::BeamSpot const & beamSpot,
18  GBRForest const * forestP) const {
19 
20  auto const & forest = *forestP;
21  auto tmva_pt_ = trk.pt();
22  auto tmva_ndof_ = trk.ndof();
23  auto tmva_nlayers_ = trk.hitPattern().trackerLayersWithMeasurement();
24  auto tmva_nlayers3D_ = trk.hitPattern().pixelLayersWithMeasurement()
27  float chi2n = trk.normalizedChi2();
28  float chi2n_no1Dmod = chi2n;
29 
30  int count1dhits = 0;
31  for (auto ith =trk.recHitsBegin(); ith!=trk.recHitsEnd(); ++ith) {
32  const auto & hit = *(*ith);
33  if (hit.dimension()==1) ++count1dhits;
34  }
35 
36  if (count1dhits > 0) {
37  float chi2 = trk.chi2();
38  float ndof = trk.ndof();
39  chi2n = (chi2+count1dhits)/float(ndof+count1dhits);
40  }
41  auto tmva_chi2n_ = chi2n;
42  auto tmva_chi2n_no1dmod_ = chi2n_no1Dmod;
43  auto tmva_eta_ = trk.eta();
44  auto tmva_relpterr_ = float(trk.ptError())/std::max(float(trk.pt()),0.000001f);
45  auto tmva_nhits_ = trk.numberOfValidHits();
48  auto tmva_minlost_ = std::min(lostIn,lostOut);
49  auto tmva_lostmidfrac_ = static_cast<float>(trk.numberOfLostHits()) / static_cast<float>(trk.numberOfValidHits() + trk.numberOfLostHits());
50 
51  float gbrVals_[PROMPT ? 16 : 12];
52  gbrVals_[0] = tmva_pt_;
53  gbrVals_[1] = tmva_lostmidfrac_;
54  gbrVals_[2] = tmva_minlost_;
55  gbrVals_[3] = tmva_nhits_;
56  gbrVals_[4] = tmva_relpterr_;
57  gbrVals_[5] = tmva_eta_;
58  gbrVals_[6] = tmva_chi2n_no1dmod_;
59  gbrVals_[7] = tmva_chi2n_;
60  gbrVals_[8] = tmva_nlayerslost_;
61  gbrVals_[9] = tmva_nlayers3D_;
62  gbrVals_[10] = tmva_nlayers_;
63  gbrVals_[11] = tmva_ndof_;
64 
65  if (PROMPT) {
66  auto tmva_absd0_ = std::abs(trk.dxy(beamSpot.position()));
67  auto tmva_absdz_ = std::abs(trk.dz(beamSpot.position()));
68  Point bestVertex = getBestVertex(trk,vertices);
69  auto tmva_absd0PV_ = std::abs(trk.dxy(bestVertex));
70  auto tmva_absdzPV_ = std::abs(trk.dz(bestVertex));
71 
72  gbrVals_[12] = tmva_absd0PV_;
73  gbrVals_[13] = tmva_absdzPV_;
74  gbrVals_[14] = tmva_absdz_;
75  gbrVals_[15] = tmva_absd0_;
76  }
77 
78 
79 
80 
81  return forest.GetClassifier(gbrVals_);
82 
83  }
84 
85  static const char * name();
86 
87  static void fillDescriptions(edm::ParameterSetDescription & desc) {
88  }
89 
90 
91 };
92 
93  using TrackMVAClassifierDetached = TrackMVAClassifier<mva<false>>;
94  using TrackMVAClassifierPrompt = TrackMVAClassifier<mva<true>>;
95  template<>
96  const char * mva<false>::name() { return "TrackMVAClassifierDetached";}
97  template<>
98  const char * mva<true>::name() { return "TrackMVAClassifierPrompt";}
99 
100 }
101 
104 
105 DEFINE_FWK_MODULE(TrackMVAClassifierDetached);
106 DEFINE_FWK_MODULE(TrackMVAClassifierPrompt);
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:556
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:821
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:493
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:512
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:646
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:335
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:544
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:550
double pt() const
track transverse momentum
Definition: TrackBase.h:616
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:758
math::XYZPoint Point
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:815
T min(T a, T b)
Definition: MathUtil.h:58
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
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:604
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
int trackerLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:529
Point getBestVertex(reco::Track const &trk, reco::VertexCollection const &vertices, const size_t minNtracks=2)
Definition: getBestVertex.h:8
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:586
int numberOfHits(HitCategory category) const
Definition: HitPattern.h:807
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109