CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Attributes
L1TrackQuality Class Reference

#include <L1TrackQuality.h>

Public Types

enum  QualityAlgorithm { QualityAlgorithm::Cut, QualityAlgorithm::GBDT, QualityAlgorithm::NN, QualityAlgorithm::None }
 

Public Member Functions

std::vector< float > featureTransform (TTTrack< Ref_Phase2TrackerDigi_ > &aTrack, std::vector< std::string > const &featureNames)
 
 L1TrackQuality ()
 
 L1TrackQuality (const edm::ParameterSet &qualityParams)
 
void setBonusFeatures (std::vector< float > bonusFeatures)
 
void setCutParameters (std::string const &AlgorithmString, float maxZ0, float maxEta, float chi2dofMax, float bendchi2Max, float minPt, int nStubmin)
 
void setL1TrackQuality (TTTrack< Ref_Phase2TrackerDigi_ > &aTrack)
 
void setONNXModel (std::string const &AlgorithmString, edm::FileInPath const &ONNXmodel, std::string const &ONNXInputName, std::vector< std::string > const &featureNames)
 
 ~L1TrackQuality ()=default
 

Private Attributes

float bendchi2Max_
 
std::vector< float > bonusFeatures_
 
float chi2dofMax_
 
std::vector< std::string > featureNames_
 
float maxEta_
 
float maxZ0_
 
float minPt_
 
int nStubsmin_
 
std::string ONNXInputName_
 
edm::FileInPath ONNXmodel_
 
QualityAlgorithm qualityAlgorithm_ = QualityAlgorithm::None
 
std::unique_ptr< cms::Ort::ONNXRuntimerunTime_
 
bool useHPH_
 

Detailed Description

Definition at line 26 of file L1TrackQuality.h.

Member Enumeration Documentation

◆ QualityAlgorithm

Enumerator
Cut 
GBDT 
NN 
None 

Definition at line 29 of file L1TrackQuality.h.

29 { Cut, GBDT, NN, None };
constexpr Matriplex::idx_t NN
Definition: Matrix.h:43

Constructor & Destructor Documentation

◆ L1TrackQuality() [1/2]

L1TrackQuality::L1TrackQuality ( )

Definition at line 11 of file L1TrackQuality.cc.

11 {}

◆ L1TrackQuality() [2/2]

L1TrackQuality::L1TrackQuality ( const edm::ParameterSet qualityParams)

Definition at line 13 of file L1TrackQuality.cc.

References edm::FileInPath::fullPath(), edm::ParameterSet::getParameter(), ONNXmodel_, runTime_, setCutParameters(), setONNXModel(), and AlCaHLTBitMon_QueryRunRegistry::string.

13  : useHPH_(false), bonusFeatures_() {
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 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::string fullPath() const
Definition: FileInPath.cc:161
void setCutParameters(std::string const &AlgorithmString, float maxZ0, float maxEta, float chi2dofMax, float bendchi2Max, float minPt, int nStubmin)
void setONNXModel(std::string const &AlgorithmString, edm::FileInPath const &ONNXmodel, std::string const &ONNXInputName, std::vector< std::string > const &featureNames)
std::vector< float > bonusFeatures_
std::unique_ptr< cms::Ort::ONNXRuntime > runTime_
edm::FileInPath ONNXmodel_

◆ ~L1TrackQuality()

L1TrackQuality::~L1TrackQuality ( )
default

Member Function Documentation

◆ featureTransform()

std::vector< float > L1TrackQuality::featureTransform ( TTTrack< Ref_Phase2TrackerDigi_ > &  aTrack,
std::vector< std::string > const &  featureNames 
)

Definition at line 35 of file L1TrackQuality.cc.

References TTTrack< T >::eta(), TrackQualityParams_cfi::featureNames, nano_mu_digi_cff::float, TTTrack_TrackWord::getBendChi2Bits(), TTTrack_TrackWord::getChi2RPhiBits(), TTTrack_TrackWord::getChi2RZBits(), TTTrack< T >::getStubRefs(), TTTrack< T >::hitPattern(), mps_fire::i, TTTrack< T >::phi(), cmsswSequenceInfo::seq, AlCaHLTBitMon_QueryRunRegistry::string, and TTTrack< T >::z0().

Referenced by setL1TrackQuality().

36  {
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_phi = aTrack.phi();
78  float tmp_trk_eta = aTrack.eta();
79 
80  // -------- fill the feature map ---------
81 
82  feature_map["nstub"] = float(tmp_trk_nstub);
83  feature_map["z0"] = tmp_trk_z0;
84  feature_map["phi"] = tmp_trk_phi;
85  feature_map["eta"] = tmp_trk_eta;
86  feature_map["nlaymiss_interior"] = float(tmp_trk_nlaymiss_interior);
87  feature_map["bendchi2_bin"] = tmp_trk_bendchi2_bin;
88  feature_map["chi2rphi_bin"] = tmp_trk_chi2rphi_bin;
89  feature_map["chi2rz_bin"] = tmp_trk_chi2rz_bin;
90 
91  // fill tensor with track params
92  transformedFeatures.reserve(featureNames.size());
93  for (const std::string& feature : featureNames)
94  transformedFeatures.push_back(feature_map[feature]);
95 
96  return transformedFeatures;
97 }
double eta() const
Track eta.
Definition: TTTrack.h:310
unsigned int getChi2RZBits() const
unsigned int getChi2RPhiBits() const
double phi() const
Track phi.
Definition: TTTrack.h:315
unsigned int hitPattern() const
Hit Pattern.
Definition: TTTrack.h:423
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
unsigned int getBendChi2Bits() const
double z0() const
Track z0.
Definition: TTTrack.h:330

◆ setBonusFeatures()

void L1TrackQuality::setBonusFeatures ( std::vector< float >  bonusFeatures)

Definition at line 192 of file L1TrackQuality.cc.

References bonusFeatures_, and useHPH_.

192  {
193  bonusFeatures_ = bonusFeatures;
194  useHPH_ = true;
195 }
std::vector< float > bonusFeatures_

◆ setCutParameters()

void L1TrackQuality::setCutParameters ( std::string const &  AlgorithmString,
float  maxZ0,
float  maxEta,
float  chi2dofMax,
float  bendchi2Max,
float  minPt,
int  nStubmin 
)

◆ setL1TrackQuality()

void L1TrackQuality::setL1TrackQuality ( TTTrack< Ref_Phase2TrackerDigi_ > &  aTrack)

Definition at line 99 of file L1TrackQuality.cc.

References funct::abs(), bendchi2Max_, TTTrack< T >::chi2(), chi2dofMax_, Cut, PV3DBase< T, PVType, FrameType >::eta(), featureNames_, featureTransform(), GBDT, TTTrack< T >::getStubRefs(), maxEta_, maxZ0_, minPt_, TTTrack< T >::momentum(), NN, nStubsmin_, ONNXInputName_, PV3DBase< T, PVType, FrameType >::perp(), qualityAlgorithm_, runTime_, TTTrack< T >::settrkMVA1(), TTTrack< T >::stubPtConsistency(), and TTTrack< T >::z0().

99  {
101  // Get Track parameters
102  float trk_pt = aTrack.momentum().perp();
103  float trk_bend_chi2 = aTrack.stubPtConsistency();
104  float trk_z0 = aTrack.z0();
105  float trk_eta = aTrack.momentum().eta();
106  float trk_chi2 = aTrack.chi2();
107  const auto& stubRefs = aTrack.getStubRefs();
108  int nStubs = stubRefs.size();
109 
110  float classification = 0.0; // Default classification is 0
111 
112  if (trk_pt >= this->minPt_ && abs(trk_z0) < this->maxZ0_ && abs(trk_eta) < this->maxEta_ &&
113  trk_chi2 < this->chi2dofMax_ && trk_bend_chi2 < this->bendchi2Max_ && nStubs >= this->nStubsmin_)
114  classification = 1.0;
115  // Classification updated to 1 if conditions are met
116 
117  aTrack.settrkMVA1(classification);
118  }
119 
121  // Setup ONNX input and output names and arrays
122  std::vector<std::string> ortinput_names;
123  std::vector<std::string> ortoutput_names;
124 
125  cms::Ort::FloatArrays ortinput;
126  cms::Ort::FloatArrays ortoutputs;
127 
128  std::vector<float> Transformed_features = featureTransform(aTrack, this->featureNames_);
129  // cms::Ort::ONNXRuntime runTime(this->ONNXmodel_.fullPath()); //Setup ONNX runtime
130 
131  ortinput_names.push_back(this->ONNXInputName_);
132  ortoutput_names = runTime_->getOutputNames();
133 
134  //ONNX runtime recieves a vector of vectors of floats so push back the input
135  // vector of float to create a 1,1,21 ortinput
136  ortinput.push_back(Transformed_features);
137 
138  // batch_size 1 as only one set of transformed features is being processed
139  int batch_size = 1;
140  // Run classification
141  ortoutputs = runTime_->run(ortinput_names, ortinput, {}, ortoutput_names, batch_size);
142 
144  aTrack.settrkMVA1(ortoutputs[0][0]);
145  }
146 
147  else if (this->qualityAlgorithm_ == QualityAlgorithm::GBDT) {
148  aTrack.settrkMVA1(ortoutputs[1][1]);
149  }
150  // Slight differences in the ONNX models of the GBDTs and NNs mean different
151  // indices of the ortoutput need to be accessed
152  }
153 
154  else {
155  aTrack.settrkMVA1(-999);
156  }
157 }
T perp() const
Definition: PV3DBase.h:69
T eta() const
Definition: PV3DBase.h:73
GlobalVector momentum() const
Track momentum.
Definition: TTTrack.h:295
std::vector< std::vector< float > > FloatArrays
Definition: ONNXRuntime.h:23
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
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)
std::unique_ptr< cms::Ort::ONNXRuntime > runTime_
double stubPtConsistency() const
StubPtConsistency.
Definition: TTTrack.h:417
std::vector< std::string > featureNames_
double z0() const
Track z0.
Definition: TTTrack.h:330
std::string ONNXInputName_

◆ setONNXModel()

void L1TrackQuality::setONNXModel ( std::string const &  AlgorithmString,
edm::FileInPath const &  ONNXmodel,
std::string const &  ONNXInputName,
std::vector< std::string > const &  featureNames 
)

Definition at line 175 of file L1TrackQuality.cc.

References TrackQualityParams_cfi::featureNames, featureNames_, GBDT, NN, None, TrackQualityParams_cfi::ONNXInputName, ONNXInputName_, TrackQualityParams_cfi::ONNXmodel, ONNXmodel_, and qualityAlgorithm_.

Referenced by L1TrackQuality().

178  {
179  //Convert algorithm string to Enum class for track by track comparison
180  if (AlgorithmString == "NN") {
182  } else if (AlgorithmString == "GBDT") {
184  } else {
186  }
190 }
QualityAlgorithm qualityAlgorithm_
std::vector< std::string > featureNames_
edm::FileInPath ONNXmodel_
std::string ONNXInputName_

Member Data Documentation

◆ bendchi2Max_

float L1TrackQuality::bendchi2Max_
private

Definition at line 71 of file L1TrackQuality.h.

Referenced by setCutParameters(), and setL1TrackQuality().

◆ bonusFeatures_

std::vector<float> L1TrackQuality::bonusFeatures_
private

Definition at line 75 of file L1TrackQuality.h.

Referenced by setBonusFeatures().

◆ chi2dofMax_

float L1TrackQuality::chi2dofMax_
private

Definition at line 70 of file L1TrackQuality.h.

Referenced by setCutParameters(), and setL1TrackQuality().

◆ featureNames_

std::vector<std::string> L1TrackQuality::featureNames_
private

Definition at line 67 of file L1TrackQuality.h.

Referenced by setL1TrackQuality(), and setONNXModel().

◆ maxEta_

float L1TrackQuality::maxEta_
private

Definition at line 69 of file L1TrackQuality.h.

Referenced by setCutParameters(), and setL1TrackQuality().

◆ maxZ0_

float L1TrackQuality::maxZ0_
private

Definition at line 68 of file L1TrackQuality.h.

Referenced by setCutParameters(), and setL1TrackQuality().

◆ minPt_

float L1TrackQuality::minPt_
private

Definition at line 72 of file L1TrackQuality.h.

Referenced by setCutParameters(), and setL1TrackQuality().

◆ nStubsmin_

int L1TrackQuality::nStubsmin_
private

Definition at line 73 of file L1TrackQuality.h.

Referenced by setCutParameters(), and setL1TrackQuality().

◆ ONNXInputName_

std::string L1TrackQuality::ONNXInputName_
private

Definition at line 66 of file L1TrackQuality.h.

Referenced by setL1TrackQuality(), and setONNXModel().

◆ ONNXmodel_

edm::FileInPath L1TrackQuality::ONNXmodel_
private

Definition at line 65 of file L1TrackQuality.h.

Referenced by L1TrackQuality(), and setONNXModel().

◆ qualityAlgorithm_

QualityAlgorithm L1TrackQuality::qualityAlgorithm_ = QualityAlgorithm::None
private

Definition at line 64 of file L1TrackQuality.h.

Referenced by setCutParameters(), setL1TrackQuality(), and setONNXModel().

◆ runTime_

std::unique_ptr<cms::Ort::ONNXRuntime> L1TrackQuality::runTime_
private

Definition at line 76 of file L1TrackQuality.h.

Referenced by L1TrackQuality(), and setL1TrackQuality().

◆ useHPH_

bool L1TrackQuality::useHPH_
private

Definition at line 74 of file L1TrackQuality.h.

Referenced by setBonusFeatures().