Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-06-04 04:35:10

0001 /*
0002 Track Quality Body file
0003 C.Brown & C.Savard 07/2020
0004 */
0005 
0006 #include "L1Trigger/TrackTrigger/interface/L1TrackQuality.h"
0007 
0008 //Constructors
0009 
0010 L1TrackQuality::L1TrackQuality() {}
0011 
0012 L1TrackQuality::L1TrackQuality(const edm::ParameterSet& qualityParams) : useHPH_(false), bonusFeatures_() {
0013   // Unpacks EDM parameter set itself to save unecessary processing within TrackProducers
0014   setModel(qualityParams.getParameter<edm::FileInPath>("model"),
0015            qualityParams.getParameter<std::vector<std::string>>("featureNames"));
0016 }
0017 
0018 std::vector<float> L1TrackQuality::featureTransform(TTTrack<Ref_Phase2TrackerDigi_>& aTrack,
0019                                                     std::vector<std::string> const& featureNames) {
0020   // List input features for MVA in proper order below, the current features options are
0021   // {"phi", "eta", "z0", "bendchi2_bin", "nstub", "nlaymiss_interior", "chi2rphi_bin",
0022   // "chi2rz_bin"}
0023   //
0024   // To use more features, they must be created here and added to feature_map below
0025 
0026   std::vector<float> transformedFeatures;
0027 
0028   // Define feature map, filled as features are generated
0029   std::map<std::string, float> feature_map;
0030 
0031   // -------- calculate feature variables --------
0032 
0033   // calculate number of missed interior layers from hitpattern
0034   int tmp_trk_hitpattern = aTrack.hitPattern();
0035   int nbits = floor(log2(tmp_trk_hitpattern)) + 1;
0036   int lay_i = 0;
0037   int tmp_trk_nlaymiss_interior = 0;
0038   bool seq = false;
0039   for (int i = 0; i < nbits; i++) {
0040     lay_i = ((1 << i) & tmp_trk_hitpattern) >> i;  //0 or 1 in ith bit (right to left)
0041 
0042     if (lay_i && !seq)
0043       seq = true;  //sequence starts when first 1 found
0044     if (!lay_i && seq)
0045       tmp_trk_nlaymiss_interior++;
0046   }
0047 
0048   // binned chi2 variables
0049   int tmp_trk_bendchi2_bin = aTrack.getBendChi2Bits();
0050   int tmp_trk_chi2rphi_bin = aTrack.getChi2RPhiBits();
0051   int tmp_trk_chi2rz_bin = aTrack.getChi2RZBits();
0052 
0053   // get the nstub
0054   std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>> stubRefs =
0055       aTrack.getStubRefs();
0056   int tmp_trk_nstub = stubRefs.size();
0057 
0058   // get other variables directly from TTTrack
0059   float tmp_trk_z0 = aTrack.z0();
0060   float tmp_trk_z0_scaled = tmp_trk_z0 / abs(aTrack.minZ0);
0061   float tmp_trk_phi = aTrack.phi();
0062   float tmp_trk_eta = aTrack.eta();
0063   float tmp_trk_tanl = aTrack.tanL();
0064 
0065   // -------- fill the feature map ---------
0066 
0067   feature_map["nstub"] = float(tmp_trk_nstub);
0068   feature_map["z0"] = tmp_trk_z0;
0069   feature_map["z0_scaled"] = tmp_trk_z0_scaled;
0070   feature_map["phi"] = tmp_trk_phi;
0071   feature_map["eta"] = tmp_trk_eta;
0072   feature_map["nlaymiss_interior"] = float(tmp_trk_nlaymiss_interior);
0073   feature_map["bendchi2_bin"] = tmp_trk_bendchi2_bin;
0074   feature_map["chi2rphi_bin"] = tmp_trk_chi2rphi_bin;
0075   feature_map["chi2rz_bin"] = tmp_trk_chi2rz_bin;
0076   feature_map["tanl"] = tmp_trk_tanl;
0077 
0078   // fill tensor with track params
0079   transformedFeatures.reserve(featureNames.size());
0080   for (const std::string& feature : featureNames)
0081     transformedFeatures.push_back(feature_map[feature]);
0082 
0083   return transformedFeatures;
0084 }
0085 
0086 void L1TrackQuality::setL1TrackQuality(TTTrack<Ref_Phase2TrackerDigi_>& aTrack) {
0087   // load in bdt
0088   conifer::BDT<float, float> bdt(this->model_.fullPath());
0089 
0090   // collect features and classify using bdt
0091   std::vector<float> inputs = featureTransform(aTrack, this->featureNames_);
0092   std::vector<float> output = bdt.decision_function(inputs);
0093   aTrack.settrkMVA1(1. / (1. + exp(-output.at(0))));
0094 }
0095 
0096 float L1TrackQuality::runEmulatedTQ(std::vector<ap_fixed<10, 5>> inputFeatures) {
0097   // load in bdt
0098 
0099   conifer::BDT<ap_fixed<10, 5>, ap_fixed<10, 5>> bdt(this->model_.fullPath());
0100 
0101   // collect features and classify using bdt
0102   std::vector<ap_fixed<10, 5>> output = bdt.decision_function(inputFeatures);
0103   return output.at(0).to_float();  // need logistic sigmoid fcn applied to xgb output
0104 }
0105 
0106 void L1TrackQuality::setModel(edm::FileInPath const& model, std::vector<std::string> const& featureNames) {
0107   //Convert algorithm string to Enum class for track by track comparison
0108   model_ = model;
0109   featureNames_ = featureNames;
0110 }
0111 
0112 void L1TrackQuality::setBonusFeatures(std::vector<float> bonusFeatures) {
0113   bonusFeatures_ = bonusFeatures;
0114   useHPH_ = true;
0115 }