CMS 3D CMS Logo

L1TrackQuality.cc
Go to the documentation of this file.
1 /*
2 Track Quality Body file
3 
4 C.Brown & C.Savard 07/2020
5 */
6 
8 
9 //Constructors
10 
12 
13 L1TrackQuality::L1TrackQuality(const edm::ParameterSet& qualityParams) : setupHPH_(), useHPH_(false) {
14  std::string AlgorithmString = qualityParams.getParameter<std::string>("qualityAlgorithm");
15  // Unpacks EDM parameter set itself to save unecessary processing within TrackProducers
16  if (AlgorithmString == "Cut") {
17  setCutParameters(AlgorithmString,
18  (float)qualityParams.getParameter<double>("maxZ0"),
19  (float)qualityParams.getParameter<double>("maxEta"),
20  (float)qualityParams.getParameter<double>("chi2dofMax"),
21  (float)qualityParams.getParameter<double>("bendchi2Max"),
22  (float)qualityParams.getParameter<double>("minPt"),
23  qualityParams.getParameter<int>("nStubsmin"));
24  }
25 
26  else {
27  setONNXModel(AlgorithmString,
28  qualityParams.getParameter<edm::FileInPath>("ONNXmodel"),
29  qualityParams.getParameter<std::string>("ONNXInputName"),
30  qualityParams.getParameter<std::vector<std::string>>("featureNames"));
31  runTime_ = std::make_unique<cms::Ort::ONNXRuntime>(this->ONNXmodel_.fullPath());
32  }
33 }
34 
36  std::vector<std::string> const& featureNames) {
37  // List input features for MVA in proper order below, the current features options are
38  // {"phi", "eta", "z0", "bendchi2_bin", "nstub", "nlaymiss_interior", "chi2rphi_bin",
39  // "chi2rz_bin"}
40  //
41  // To use more features, they must be created here and added to feature_map below
42 
43  std::vector<float> transformedFeatures;
44 
45  // Define feature map, filled as features are generated
46  std::map<std::string, float> feature_map;
47 
48  // -------- calculate feature variables --------
49 
50  // calculate number of missed interior layers from hitpattern
51  int tmp_trk_hitpattern = aTrack.hitPattern();
52  int nbits = floor(log2(tmp_trk_hitpattern)) + 1;
53  int lay_i = 0;
54  int tmp_trk_nlaymiss_interior = 0;
55  bool seq = false;
56  for (int i = 0; i < nbits; i++) {
57  lay_i = ((1 << i) & tmp_trk_hitpattern) >> i; //0 or 1 in ith bit (right to left)
58 
59  if (lay_i && !seq)
60  seq = true; //sequence starts when first 1 found
61  if (!lay_i && seq)
62  tmp_trk_nlaymiss_interior++;
63  }
64 
65  // binned chi2 variables
66  int tmp_trk_bendchi2_bin = aTrack.getBendChi2Bits();
67  int tmp_trk_chi2rphi_bin = aTrack.getChi2RPhiBits();
68  int tmp_trk_chi2rz_bin = aTrack.getChi2RZBits();
69 
70  // Extra variables from HitPatternHelper don't improve performance, so commented out.
71  //if (useHPH_) {
72  // hph::HitPatternHelper hph(setupHPH_, tmp_trk_hitpattern, tmp_trk_tanL, tmp_trk_z0);
73  // hitpattern_expanded_binary = hph.binary();
74  // tmp_trk_nlaymiss_PS = hph.numMissingPS();
75  // tmp_trk_nlaymiss_2S = hph.numMissing2S();
76  //}
77 
78  // get the nstub
79  std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>> stubRefs =
80  aTrack.getStubRefs();
81  int tmp_trk_nstub = stubRefs.size();
82 
83  // get other variables directly from TTTrack
84  float tmp_trk_z0 = aTrack.z0();
85  float tmp_trk_phi = aTrack.phi();
86  float tmp_trk_eta = aTrack.eta();
87 
88  // -------- fill the feature map ---------
89 
90  feature_map["nstub"] = float(tmp_trk_nstub);
91  feature_map["z0"] = tmp_trk_z0;
92  feature_map["phi"] = tmp_trk_phi;
93  feature_map["eta"] = tmp_trk_eta;
94  feature_map["nlaymiss_interior"] = float(tmp_trk_nlaymiss_interior);
95  feature_map["bendchi2_bin"] = tmp_trk_bendchi2_bin;
96  feature_map["chi2rphi_bin"] = tmp_trk_chi2rphi_bin;
97  feature_map["chi2rz_bin"] = tmp_trk_chi2rz_bin;
98 
99  // fill tensor with track params
100  transformedFeatures.reserve(featureNames.size());
101  for (const std::string& feature : featureNames)
102  transformedFeatures.push_back(feature_map[feature]);
103 
104  return transformedFeatures;
105 }
106 
109  // Get Track parameters
110  float trk_pt = aTrack.momentum().perp();
111  float trk_bend_chi2 = aTrack.stubPtConsistency();
112  float trk_z0 = aTrack.z0();
113  float trk_eta = aTrack.momentum().eta();
114  float trk_chi2 = aTrack.chi2();
115  const auto& stubRefs = aTrack.getStubRefs();
116  int nStubs = stubRefs.size();
117 
118  float classification = 0.0; // Default classification is 0
119 
120  if (trk_pt >= this->minPt_ && abs(trk_z0) < this->maxZ0_ && abs(trk_eta) < this->maxEta_ &&
121  trk_chi2 < this->chi2dofMax_ && trk_bend_chi2 < this->bendchi2Max_ && nStubs >= this->nStubsmin_)
122  classification = 1.0;
123  // Classification updated to 1 if conditions are met
124 
125  aTrack.settrkMVA1(classification);
126  }
127 
129  // Setup ONNX input and output names and arrays
130  std::vector<std::string> ortinput_names;
131  std::vector<std::string> ortoutput_names;
132 
133  cms::Ort::FloatArrays ortinput;
134  cms::Ort::FloatArrays ortoutputs;
135 
136  std::vector<float> Transformed_features = featureTransform(aTrack, this->featureNames_);
137  // cms::Ort::ONNXRuntime runTime(this->ONNXmodel_.fullPath()); //Setup ONNX runtime
138 
139  ortinput_names.push_back(this->ONNXInputName_);
140  ortoutput_names = runTime_->getOutputNames();
141 
142  //ONNX runtime recieves a vector of vectors of floats so push back the input
143  // vector of float to create a 1,1,21 ortinput
144  ortinput.push_back(Transformed_features);
145 
146  // batch_size 1 as only one set of transformed features is being processed
147  int batch_size = 1;
148  // Run classification
149  ortoutputs = runTime_->run(ortinput_names, ortinput, {}, ortoutput_names, batch_size);
150 
152  aTrack.settrkMVA1(ortoutputs[0][0]);
153  }
154 
155  else if (this->qualityAlgorithm_ == QualityAlgorithm::GBDT) {
156  aTrack.settrkMVA1(ortoutputs[1][1]);
157  }
158  // Slight differences in the ONNX models of the GBDTs and NNs mean different
159  // indices of the ortoutput need to be accessed
160  }
161 
162  else {
163  aTrack.settrkMVA1(-999);
164  }
165 }
166 
167 void L1TrackQuality::setCutParameters(std::string const& AlgorithmString,
168  float maxZ0,
169  float maxEta,
170  float chi2dofMax,
171  float bendchi2Max,
172  float minPt,
173  int nStubmin) {
175  maxZ0_ = maxZ0;
176  maxEta_ = maxEta;
179  minPt_ = minPt;
180  nStubsmin_ = nStubmin;
181 }
182 
183 void L1TrackQuality::setONNXModel(std::string const& AlgorithmString,
184  edm::FileInPath const& ONNXmodel,
185  std::string const& ONNXInputName,
186  std::vector<std::string> const& featureNames) {
187  //Convert algorithm string to Enum class for track by track comparison
188  if (AlgorithmString == "NN") {
190  } else if (AlgorithmString == "GBDT") {
192  } else {
194  }
198 }
199 
200 void L1TrackQuality::beginRun(const hph::Setup* setupHPH) {
201  // Initialize helper for decoding hit patterns.
202  setupHPH_ = setupHPH;
203  useHPH_ = true;
204 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T perp() const
Definition: PV3DBase.h:69
double eta() const
Track eta.
Definition: TTTrack.h:310
std::string fullPath() const
Definition: FileInPath.cc:161
const hph::Setup * setupHPH_
unsigned int getChi2RZBits() const
void setL1TrackQuality(TTTrack< Ref_Phase2TrackerDigi_ > &aTrack)
void setCutParameters(std::string const &AlgorithmString, float maxZ0, float maxEta, float chi2dofMax, float bendchi2Max, float minPt, int nStubmin)
T eta() const
Definition: PV3DBase.h:73
unsigned int getChi2RPhiBits() const
double phi() const
Track phi.
Definition: TTTrack.h:315
minPt
Definition: PV_cfg.py:223
unsigned int hitPattern() const
Hit Pattern.
Definition: TTTrack.h:423
GlobalVector momentum() const
Track momentum.
Definition: TTTrack.h:295
std::vector< std::vector< float > > FloatArrays
Definition: ONNXRuntime.h:23
void setONNXModel(std::string const &AlgorithmString, edm::FileInPath const &ONNXmodel, std::string const &ONNXInputName, std::vector< std::string > const &featureNames)
void settrkMVA1(double atrkMVA1)
Definition: TTTrack.h:381
QualityAlgorithm qualityAlgorithm_
double chi2() const
Chi2.
Definition: TTTrack.h:341
void beginRun(const hph::Setup *setup)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
std::vector< edm::Ref< edmNew::DetSetVector< TTStub< T > >, TTStub< T > > > getStubRefs() const
Track components.
Definition: TTTrack.h:93
std::vector< float > featureTransform(TTTrack< Ref_Phase2TrackerDigi_ > &aTrack, std::vector< std::string > const &featureNames)
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
std::unique_ptr< cms::Ort::ONNXRuntime > runTime_
double stubPtConsistency() const
StubPtConsistency.
Definition: TTTrack.h:417
std::vector< std::string > featureNames_
unsigned int getBendChi2Bits() const
edm::FileInPath ONNXmodel_
double z0() const
Track z0.
Definition: TTTrack.h:330
std::string ONNXInputName_