CMS 3D CMS Logo

L1TrackQuality.cc
Go to the documentation of this file.
1 /*
2 Track Quality Body file
3 C.Brown & C.Savard 07/2020
4 */
5 
7 
8 //Constructors
9 
11 
12 L1TrackQuality::L1TrackQuality(const edm::ParameterSet& qualityParams) : useHPH_(false), bonusFeatures_() {
13  std::string AlgorithmString = qualityParams.getParameter<std::string>("qualityAlgorithm");
14  // Unpacks EDM parameter set itself to save unecessary processing within TrackProducers
15  if (AlgorithmString == "Cut") {
16  setCutParameters(AlgorithmString,
17  (float)qualityParams.getParameter<double>("maxZ0"),
18  (float)qualityParams.getParameter<double>("maxEta"),
19  (float)qualityParams.getParameter<double>("chi2dofMax"),
20  (float)qualityParams.getParameter<double>("bendchi2Max"),
21  (float)qualityParams.getParameter<double>("minPt"),
22  qualityParams.getParameter<int>("nStubsmin"));
23  }
24 
25  else {
26  setONNXModel(AlgorithmString,
27  qualityParams.getParameter<edm::FileInPath>("ONNXmodel"),
28  qualityParams.getParameter<std::string>("ONNXInputName"),
29  qualityParams.getParameter<std::vector<std::string>>("featureNames"));
30  if ((AlgorithmString == "GBDT") || (AlgorithmString == "NN"))
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  // get the nstub
71  std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>> stubRefs =
72  aTrack.getStubRefs();
73  int tmp_trk_nstub = stubRefs.size();
74 
75  // get other variables directly from TTTrack
76  float tmp_trk_z0 = aTrack.z0();
77  float tmp_trk_z0_scaled = tmp_trk_z0 / abs(aTrack.minZ0);
78  float tmp_trk_phi = aTrack.phi();
79  float tmp_trk_eta = aTrack.eta();
80  float tmp_trk_tanl = aTrack.tanL();
81 
82  // -------- fill the feature map ---------
83 
84  feature_map["nstub"] = float(tmp_trk_nstub);
85  feature_map["z0"] = tmp_trk_z0;
86  feature_map["z0_scaled"] = tmp_trk_z0_scaled;
87  feature_map["phi"] = tmp_trk_phi;
88  feature_map["eta"] = tmp_trk_eta;
89  feature_map["nlaymiss_interior"] = float(tmp_trk_nlaymiss_interior);
90  feature_map["bendchi2_bin"] = tmp_trk_bendchi2_bin;
91  feature_map["chi2rphi_bin"] = tmp_trk_chi2rphi_bin;
92  feature_map["chi2rz_bin"] = tmp_trk_chi2rz_bin;
93  feature_map["tanl"] = tmp_trk_tanl;
94 
95  // fill tensor with track params
96  transformedFeatures.reserve(featureNames.size());
97  for (const std::string& feature : featureNames)
98  transformedFeatures.push_back(feature_map[feature]);
99 
100  return transformedFeatures;
101 }
102 
105  // Get Track parameters
106  float trk_pt = aTrack.momentum().perp();
107  float trk_bend_chi2 = aTrack.stubPtConsistency();
108  float trk_z0 = aTrack.z0();
109  float trk_eta = aTrack.momentum().eta();
110  float trk_chi2 = aTrack.chi2();
111  const auto& stubRefs = aTrack.getStubRefs();
112  int nStubs = stubRefs.size();
113 
114  float classification = 0.0; // Default classification is 0
115 
116  if (trk_pt >= this->minPt_ && abs(trk_z0) < this->maxZ0_ && abs(trk_eta) < this->maxEta_ &&
117  trk_chi2 < this->chi2dofMax_ && trk_bend_chi2 < this->bendchi2Max_ && nStubs >= this->nStubsmin_)
118  classification = 1.0;
119  // Classification updated to 1 if conditions are met
120 
121  aTrack.settrkMVA1(classification);
122  }
123 
124  else if (this->qualityAlgorithm_ == QualityAlgorithm::GBDT_cpp) {
125  // load in bdt
126  conifer::BDT<float, float> bdt(this->ONNXmodel_.fullPath());
127 
128  // collect features and classify using bdt
129  std::vector<float> inputs = featureTransform(aTrack, this->featureNames_);
130  std::vector<float> output = bdt.decision_function(inputs);
131  aTrack.settrkMVA1(1. / (1. + exp(-output.at(0)))); // need logistic sigmoid fcn applied to xgb output
132  }
133 
135  // Setup ONNX input and output names and arrays
136  std::vector<std::string> ortinput_names;
137  std::vector<std::string> ortoutput_names;
138 
139  cms::Ort::FloatArrays ortinput;
140  cms::Ort::FloatArrays ortoutputs;
141 
142  std::vector<float> Transformed_features = featureTransform(aTrack, this->featureNames_);
143  // cms::Ort::ONNXRuntime runTime(this->ONNXmodel_.fullPath()); //Setup ONNX runtime
144 
145  ortinput_names.push_back(this->ONNXInputName_);
146  ortoutput_names = runTime_->getOutputNames();
147 
148  //ONNX runtime recieves a vector of vectors of floats so push back the input
149  // vector of float to create a 1,1,21 ortinput
150  ortinput.push_back(Transformed_features);
151 
152  // batch_size 1 as only one set of transformed features is being processed
153  int batch_size = 1;
154  // Run classification
155  ortoutputs = runTime_->run(ortinput_names, ortinput, {}, ortoutput_names, batch_size);
156 
158  aTrack.settrkMVA1(ortoutputs[0][0]);
159  }
160 
161  else if (this->qualityAlgorithm_ == QualityAlgorithm::GBDT) {
162  aTrack.settrkMVA1(ortoutputs[1][1]);
163  }
164  // Slight differences in the ONNX models of the GBDTs and NNs mean different
165  // indices of the ortoutput need to be accessed
166  }
167 
168  else {
169  aTrack.settrkMVA1(-999);
170  }
171 }
172 
173 float L1TrackQuality::runEmulatedTQ(std::vector<ap_fixed<10, 5>> inputFeatures) {
174  // load in bdt
175 
176  conifer::BDT<ap_fixed<10, 5>, ap_fixed<10, 5>> bdt(this->ONNXmodel_.fullPath());
177 
178  // collect features and classify using bdt
179  std::vector<ap_fixed<10, 5>> output = bdt.decision_function(inputFeatures);
180  return output.at(0).to_float(); // need logistic sigmoid fcn applied to xgb output
181 }
182 
183 void L1TrackQuality::setCutParameters(std::string const& AlgorithmString,
184  float maxZ0,
185  float maxEta,
186  float chi2dofMax,
187  float bendchi2Max,
188  float minPt,
189  int nStubmin) {
191  maxZ0_ = maxZ0;
192  maxEta_ = maxEta;
195  minPt_ = minPt;
196  nStubsmin_ = nStubmin;
197 }
198 
199 void L1TrackQuality::setONNXModel(std::string const& AlgorithmString,
200  edm::FileInPath const& ONNXmodel,
201  std::string const& ONNXInputName,
202  std::vector<std::string> const& featureNames) {
203  //Convert algorithm string to Enum class for track by track comparison
204  if (AlgorithmString == "NN") {
206  } else if (AlgorithmString == "GBDT") {
208  } else if (AlgorithmString == "GBDT_cpp") {
210  } else {
212  }
216 }
217 
218 void L1TrackQuality::setBonusFeatures(std::vector<float> bonusFeatures) {
219  bonusFeatures_ = bonusFeatures;
220  useHPH_ = true;
221 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T perp() const
Definition: PV3DBase.h:69
double eta() const
Track eta.
Definition: TTTrack.h:310
std::string fullPath() const
Definition: FileInPath.cc:161
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)
static constexpr double minZ0
T eta() const
Definition: PV3DBase.h:73
unsigned int getChi2RPhiBits() const
double phi() const
Track phi.
Definition: TTTrack.h:315
double tanL() const
Track tanL.
Definition: TTTrack.h:305
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float runEmulatedTQ(std::vector< ap_fixed< 10, 5 >> inputFeatures)
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::vector< float > bonusFeatures_
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_
Definition: output.py:1
void setBonusFeatures(std::vector< float > bonusFeatures)