Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:58

0001 #include "L1Trigger/L1TMuonEndCap/interface/PtAssignmentEngineDxy.h"
0002 
0003 #include <cassert>
0004 #include <iostream>
0005 #include <sstream>
0006 
0007 #include "helper.h"  // assert_no_abort
0008 
0009 PtAssignmentEngineDxy::PtAssignmentEngineDxy() : graphDefDxy_(nullptr), sessionDxy_(nullptr) {}
0010 
0011 PtAssignmentEngineDxy::~PtAssignmentEngineDxy() {
0012   if (sessionDxy_ != nullptr) {
0013     tensorflow::closeSession(sessionDxy_);
0014   }
0015   delete graphDefDxy_;
0016 }
0017 
0018 void PtAssignmentEngineDxy::configure(int verbose, const std::string pbFileNameDxy) {
0019   verbose_ = verbose;
0020 
0021   pbFileNameDxy_ = pbFileNameDxy;
0022   std::string pbFilePathDxy_ = "L1Trigger/L1TMuon/data/emtf_luts/" + pbFileNameDxy_;
0023 
0024   inputNameDxy_ = "input1";
0025   outputNamesDxy_ = {"Identity"};
0026 
0027   if (graphDefDxy_ == nullptr) {
0028     graphDefDxy_ = tensorflow::loadGraphDef(edm::FileInPath(pbFilePathDxy_).fullPath());
0029   }
0030   emtf_assert(graphDefDxy_ != nullptr);
0031 
0032   if (sessionDxy_ == nullptr) {
0033     sessionDxy_ = tensorflow::createSession(graphDefDxy_);
0034   }
0035 
0036   emtf_assert(sessionDxy_ != nullptr);
0037 }
0038 
0039 const PtAssignmentEngineAux2017& PtAssignmentEngineDxy::aux() const {
0040   static const PtAssignmentEngineAux2017 instance;
0041   return instance;
0042 }
0043 
0044 void PtAssignmentEngineDxy::calculate_pt_dxy(const EMTFTrack& track,
0045                                              emtf::Feature& feature,
0046                                              emtf::Prediction& prediction) const {
0047   // This is called for each track instead of for entire track collection as was done in Phase-2 implementation
0048   preprocessing_dxy(track, feature);
0049   call_tensorflow_dxy(feature, prediction);
0050   return;
0051 }
0052 
0053 void PtAssignmentEngineDxy::preprocessing_dxy(const EMTFTrack& track, emtf::Feature& feature) const {
0054   // Mimic Phase-1 EMTF input calculations
0055   // 6 delta Phis: S1-S2, S1-S3, S1-S4, S2-S3, S2-S4, S3-S4
0056   // 6 delta Thetas: S1-S2, S1-S3, S1-S4, S2-S3, S2-S4, S3-S4
0057   // 6 delta Phi signs: S1-S2, S1-S3, S1-S4, S2-S3, S2-S4, S3-S4
0058   // 6 delta Theta signs: S1-S2, S1-S3, S1-S4, S2-S3, S2-S4, S3-S4
0059   // 1 track Theta taken from stub coordinate in ME2, ME3, ME4 (in this priority)
0060   // 4 CSC pattern values (Run 2 convention): S1, S2, S3, S4
0061   // Total: 29 variables
0062   std::array<float, 6> x_dphi;
0063   std::array<float, 6> x_dphi_sign;
0064   std::array<float, 6> x_dtheta;
0065   std::array<float, 6> x_dtheta_sign;
0066   std::array<float, 1> x_trk_theta;
0067   std::array<float, 4> x_csc_pattern;
0068 
0069   // Initialize x_csc_pattern to zeros
0070   x_csc_pattern.fill(0);
0071 
0072   EMTFPtLUT data = track.PtLUT();
0073 
0074   const int invalid_dtheta = 127;
0075   const int invalid_dphi = 8191;
0076 
0077   // // Which stations have hits
0078   bool st1 = (track.Mode() >= 8);
0079   bool st2 = ((track.Mode() % 8) >= 4);
0080   bool st3 = ((track.Mode() % 4) >= 2);
0081   bool st4 = ((track.Mode() % 2) == 1);
0082 
0083   // Get valid pattern values
0084   if (st1)
0085     x_csc_pattern[0] = data.cpattern[0];
0086   if (st2)
0087     x_csc_pattern[1] = data.cpattern[1];
0088   if (st3)
0089     x_csc_pattern[2] = data.cpattern[2];
0090   if (st4)
0091     x_csc_pattern[3] = data.cpattern[3];
0092 
0093   for (int i = 0; i < 6; ++i) {  // There are 6 deltas between 4 stations.
0094     // Calculate delta phi
0095     x_dphi[i] = (data.delta_ph[i] != invalid_dphi) ? data.delta_ph[i] : 0;
0096 
0097     // Calculate delta theta
0098     x_dtheta[i] = (data.delta_th[i] != invalid_dtheta) ? data.delta_th[i] : 0;
0099 
0100     // Get delta phi and theta signs
0101     x_dphi_sign[i] = data.sign_ph[i];
0102     x_dtheta_sign[i] = data.sign_th[i];
0103   }
0104 
0105   // Set dPhi and dTheta values to 0 if there was no hit in the station
0106   if (!st1) {
0107     x_dphi[0] = 0;
0108     x_dphi[1] = 0;
0109     x_dphi[2] = 0;
0110 
0111     x_dtheta[0] = 0;
0112     x_dtheta[1] = 0;
0113     x_dtheta[2] = 0;
0114   }
0115   if (!st2) {
0116     x_dphi[0] = 0;
0117     x_dphi[3] = 0;
0118     x_dphi[4] = 0;
0119 
0120     x_dtheta[0] = 0;
0121     x_dtheta[3] = 0;
0122     x_dtheta[4] = 0;
0123   }
0124   if (!st3) {
0125     x_dphi[1] = 0;
0126     x_dphi[3] = 0;
0127     x_dphi[5] = 0;
0128 
0129     x_dtheta[1] = 0;
0130     x_dtheta[3] = 0;
0131     x_dtheta[5] = 0;
0132   }
0133   if (!st4) {
0134     x_dphi[2] = 0;
0135     x_dphi[4] = 0;
0136     x_dphi[5] = 0;
0137 
0138     x_dtheta[2] = 0;
0139     x_dtheta[4] = 0;
0140     x_dtheta[5] = 0;
0141   }
0142 
0143   x_trk_theta[0] = track.Theta_fp();
0144 
0145   // Set NN inputs
0146   feature = {{x_dphi[0],        x_dphi[1],        x_dphi[2],        x_dphi[3],        x_dphi[4],
0147               x_dphi[5],        x_dphi_sign[0],   x_dphi_sign[1],   x_dphi_sign[2],   x_dphi_sign[3],
0148               x_dphi_sign[4],   x_dphi_sign[5],   x_dtheta[0],      x_dtheta[1],      x_dtheta[2],
0149               x_dtheta[3],      x_dtheta[4],      x_dtheta[5],      x_dtheta_sign[0], x_dtheta_sign[1],
0150               x_dtheta_sign[2], x_dtheta_sign[3], x_dtheta_sign[4], x_dtheta_sign[5], x_csc_pattern[0],
0151               x_csc_pattern[1], x_csc_pattern[2], x_csc_pattern[3], x_trk_theta[0]}};
0152 
0153   return;
0154 }
0155 
0156 void PtAssignmentEngineDxy::call_tensorflow_dxy(const emtf::Feature& feature, emtf::Prediction& prediction) const {
0157   tensorflow::Tensor input(tensorflow::DT_FLOAT, {1, emtf::NUM_FEATURES});
0158   std::vector<tensorflow::Tensor> outputs;
0159   emtf_assert(feature.size() == emtf::NUM_FEATURES);
0160 
0161   float* d = input.flat<float>().data();
0162   std::copy(feature.begin(), feature.end(), d);
0163   tensorflow::run(sessionDxy_, {{inputNameDxy_, input}}, outputNamesDxy_, &outputs);
0164   emtf_assert(outputs.size() == 1);
0165   emtf_assert(prediction.size() == emtf::NUM_PREDICTIONS);
0166 
0167   prediction.at(0) = outputs[0].matrix<float>()(0, 0);
0168   prediction.at(1) = outputs[0].matrix<float>()(0, 1);
0169 
0170   return;
0171 }