Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:21:59

0001 //
0002 //  Created by J.Li on 1/23/21.
0003 //
0004 
0005 #include "L1Trigger/TrackFindingTracklet/interface/HitPatternHelper.h"
0006 #include "L1Trigger/TrackFindingTMTT/interface/KFbase.h"
0007 #include "L1Trigger/TrackFindingTMTT/interface/TrackerModule.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 
0010 #include <vector>
0011 #include <algorithm>
0012 #include <cmath>
0013 
0014 namespace hph {
0015 
0016   Setup::Setup(const edm::ParameterSet& iConfig,
0017                const tt::Setup& setupTT,
0018                const trackerTFP::DataFormats& dataFormats,
0019                const trackerTFP::LayerEncoding& layerEncoding)
0020       : iConfig_(iConfig),
0021         oldKFPSet_(iConfig.getParameter<edm::ParameterSet>("oldKFPSet")),
0022         setupTT_(setupTT),
0023         dataFormats_(dataFormats),
0024         dfcot_(dataFormats_.format(trackerTFP::Variable::cot, trackerTFP::Process::kfin)),
0025         dfzT_(dataFormats_.format(trackerTFP::Variable::zT, trackerTFP::Process::kfin)),
0026         layerEncoding_(layerEncoding),
0027         hphDebug_(iConfig.getParameter<bool>("hphDebug")),
0028         useNewKF_(iConfig.getParameter<bool>("useNewKF")),
0029         chosenRofZNewKF_(setupTT_.chosenRofZ()),
0030         etaRegionsNewKF_(setupTT_.boundarieEta()),
0031         layermap_(),
0032         nEtaRegions_(tmtt::KFbase::nEta_ / 2),
0033         nKalmanLayers_(tmtt::KFbase::nKFlayer_) {
0034     if (useNewKF_) {
0035       chosenRofZ_ = chosenRofZNewKF_;
0036       etaRegions_ = etaRegionsNewKF_;
0037     } else {
0038       chosenRofZ_ = oldKFPSet_.getParameter<double>("ChosenRofZ");
0039       etaRegions_ = oldKFPSet_.getParameter<std::vector<double>>("EtaRegions");
0040     }
0041     static constexpr auto layerIds = {1, 2, 3, 4, 5, 6, 11, 12, 13, 14, 15};  //layer ID 11~15 correspond to D1~D5
0042     // Converting tmtt::KFbase::layerMap_ to a format that is acceptatble by HitPatternHelper
0043     for (int i = 0; i < nEtaRegions_; i++) {
0044       for (int j : layerIds) {
0045         int layer;
0046         if (j < 7) {
0047           layer = tmtt::KFbase::layerMap_[i][tmtt::TrackerModule::calcLayerIdReduced(j)].first;
0048         } else {
0049           layer = tmtt::KFbase::layerMap_[i][tmtt::TrackerModule::calcLayerIdReduced(j)].second;
0050         }
0051         if (layer < nKalmanLayers_) {
0052           layermap_[i][layer].push_back(j);
0053         }
0054       }
0055     }
0056   }
0057 
0058   HitPatternHelper::HitPatternHelper(const Setup* setup, int hitpattern, double cot, double z0)
0059       : setup_(setup),
0060         hphDebug_(setup_->hphDebug()),
0061         useNewKF_(setup_->useNewKF()),
0062         etaRegions_(setup_->etaRegions()),
0063         layermap_(setup_->layermap()),
0064         nKalmanLayers_(setup_->nKalmanLayers()),
0065         etaBin_(setup_->etaRegion(z0, cot, true)),
0066         cotBin_(setup_->digiCot(cot, etaBin_)),
0067         zTBin_(setup_->digiZT(z0, cot, etaBin_)),
0068         layerEncoding_(setup->layerEncoding(etaBin_, zTBin_, cotBin_)),
0069         layerEncodingMap_(setup->layerEncodingMap(etaBin_, zTBin_, cotBin_)),
0070         numExpLayer_(layerEncoding_.size()),
0071         hitpattern_(hitpattern),
0072         etaSector_(setup_->etaRegion(z0, cot, useNewKF_)),
0073         numMissingLayer_(0),
0074         numMissingPS_(0),
0075         numMissing2S_(0),
0076         numPS_(0),
0077         num2S_(0),
0078         numMissingInterior1_(0),
0079         numMissingInterior2_(0),
0080         binary_(11, 0),  //there are 11 unique layer IDs, as defined in variable "layerIds"
0081         bonusFeatures_() {
0082     int kf_eta_reg = etaSector_;
0083     if (kf_eta_reg < ((int)etaRegions_.size() - 1) / 2) {
0084       kf_eta_reg = ((int)etaRegions_.size() - 1) / 2 - 1 - kf_eta_reg;
0085     } else {
0086       kf_eta_reg = kf_eta_reg - (int)(etaRegions_.size() - 1) / 2;
0087     }
0088 
0089     int nbits = floor(log2(hitpattern_)) + 1;
0090     int lay_i = 0;
0091     bool seq = false;
0092     for (int i = 0; i < nbits; i++) {
0093       lay_i = ((1 << i) & hitpattern_) >> i;  //0 or 1 in ith bit (right to left)
0094 
0095       if (lay_i && !seq)
0096         seq = true;  //sequence starts when first 1 found
0097       if (!lay_i && seq) {
0098         numMissingInterior1_++;  //This is the same as the "tmp_trk_nlaymiss_interior" calculated in Trackquality.cc
0099       }
0100       if (!lay_i) {
0101         bool realhit = false;
0102         if (layermap_[kf_eta_reg][i].empty())
0103           continue;
0104         for (int j : layermap_[kf_eta_reg][i]) {
0105           int k = findLayer(j);
0106           if (k > 0)
0107             realhit = true;
0108         }
0109         if (realhit)
0110           numMissingInterior2_++;  //This variable doesn't make sense for new KF because it uses the layermap from Old KF
0111       }
0112     }
0113 
0114     if (hphDebug_) {
0115       if (useNewKF_) {
0116         edm::LogVerbatim("TrackTriggerHPH") << "Running with New KF";
0117       } else {
0118         edm::LogVerbatim("TrackTriggerHPH") << "Running with Old KF";
0119       }
0120       edm::LogVerbatim("TrackTriggerHPH") << "======================================================";
0121       edm::LogVerbatim("TrackTriggerHPH")
0122           << "Looking at hitpattern " << std::bitset<7>(hitpattern_) << "; Looping over KF layers:";
0123     }
0124 
0125     if (useNewKF_) {
0126       //New KF uses sensor modules to determine the hitmask already
0127       for (int i = 0; i < numExpLayer_; i++) {
0128         if (hphDebug_) {
0129           edm::LogVerbatim("TrackTriggerHPH") << "--------------------------";
0130           edm::LogVerbatim("TrackTriggerHPH") << "Looking at KF layer " << i;
0131           if (layerEncoding_[i] < 10) {
0132             edm::LogVerbatim("TrackTriggerHPH") << "KF expects L" << layerEncoding_[i];
0133           } else {
0134             edm::LogVerbatim("TrackTriggerHPH") << "KF expects D" << layerEncoding_[i] - 10;
0135           }
0136         }
0137 
0138         if (((1 << i) & hitpattern_) >> i) {
0139           if (hphDebug_) {
0140             edm::LogVerbatim("TrackTriggerHPH") << "Layer found in hitpattern";
0141           }
0142 
0143           binary_[reducedId(layerEncoding_[i])] = 1;
0144           if (layerEncodingMap_[layerEncoding_[i]]->psModule()) {
0145             numPS_++;
0146           } else {
0147             num2S_++;
0148           }
0149         } else {
0150           if (hphDebug_) {
0151             edm::LogVerbatim("TrackTriggerHPH") << "Layer missing in hitpattern";
0152           }
0153 
0154           if (layerEncodingMap_[layerEncoding_[i]]->psModule()) {
0155             numMissingPS_++;
0156           } else {
0157             numMissing2S_++;
0158           }
0159         }
0160       }
0161 
0162     } else {
0163       //Old KF uses the hard coded layermap to determien hitmask
0164       for (int i = 0; i < nKalmanLayers_; i++) {  //Loop over each digit of hitpattern
0165 
0166         if (hphDebug_) {
0167           edm::LogVerbatim("TrackTriggerHPH") << "--------------------------";
0168           edm::LogVerbatim("TrackTriggerHPH") << "Looking at KF layer " << i;
0169         }
0170 
0171         if (layermap_[kf_eta_reg][i].empty()) {
0172           if (hphDebug_) {
0173             edm::LogVerbatim("TrackTriggerHPH") << "KF does not expect this layer";
0174           }
0175 
0176           continue;
0177         }
0178 
0179         for (int j :
0180              layermap_[kf_eta_reg][i]) {  //Find out which layer the Old KF is dealing with when hitpattern is encoded
0181 
0182           if (hphDebug_) {
0183             if (j < 10) {
0184               edm::LogVerbatim("TrackTriggerHPH") << "KF expects L" << j;
0185             } else {
0186               edm::LogVerbatim("TrackTriggerHPH") << "KF expects D" << j - 10;
0187             }
0188           }
0189 
0190           int k = findLayer(j);
0191           if (k < 0) {
0192             //k<0 means even though layer j is predicted by Old KF, this prediction is rejected because it contradicts
0193             if (hphDebug_) {  //a more accurate prediction made with the help of information from sensor modules
0194               edm::LogVerbatim("TrackTriggerHPH") << "Rejected by sensor modules";
0195             }
0196 
0197             continue;
0198           }
0199 
0200           if (hphDebug_) {
0201             edm::LogVerbatim("TrackTriggerHPH") << "Confirmed by sensor modules";
0202           }
0203           //prediction is accepted
0204           if (((1 << i) & hitpattern_) >> i) {
0205             if (hphDebug_) {
0206               edm::LogVerbatim("TrackTriggerHPH") << "Layer found in hitpattern";
0207             }
0208 
0209             binary_[reducedId(j)] = 1;
0210             if (layerEncodingMap_[layerEncoding_[k]]->psModule()) {
0211               numPS_++;
0212             } else {
0213               num2S_++;
0214             }
0215           } else {
0216             if (hphDebug_) {
0217               edm::LogVerbatim("TrackTriggerHPH") << "Layer missing in hitpattern";
0218             }
0219 
0220             if (layerEncodingMap_[layerEncoding_[k]]->psModule()) {
0221               numMissingPS_++;
0222             } else {
0223               numMissing2S_++;
0224             }
0225           }
0226         }
0227       }
0228     }
0229 
0230     if (hphDebug_) {
0231       edm::LogVerbatim("TrackTriggerHPH") << "------------------------------";
0232       edm::LogVerbatim("TrackTriggerHPH") << "numPS = " << numPS_ << ", num2S = " << num2S_
0233                                           << ", missingPS = " << numMissingPS_ << ", missing2S = " << numMissing2S_;
0234       edm::LogVerbatim("TrackTriggerHPH") << "======================================================";
0235     }
0236   }
0237 
0238   int Setup::etaRegion(double z0, double cot, bool useNewKF) const {
0239     //Calculating eta sector based on cot and z0
0240     double chosenRofZ;
0241     std::vector<double> etaRegions;
0242     if (useNewKF) {
0243       chosenRofZ = chosenRofZNewKF_;
0244       etaRegions = etaRegionsNewKF_;
0245     } else {
0246       chosenRofZ = chosenRofZ_;
0247       etaRegions = etaRegions_;
0248     }
0249     double kfzRef = z0 + chosenRofZ * cot;
0250     int kf_eta_reg = 0;
0251     for (int iEtaSec = 1; iEtaSec < ((int)etaRegions.size() - 1); iEtaSec++) {  // Doesn't apply eta < 2.4 cut.
0252       double etaMax = etaRegions[iEtaSec];
0253       double zRefMax = chosenRofZ / tan(2. * atan(exp(-etaMax)));
0254       if (kfzRef < zRefMax)
0255         break;
0256       kf_eta_reg = iEtaSec;
0257     }
0258     return kf_eta_reg;
0259   }
0260 
0261   int Setup::digiCot(double cot, int binEta) const {
0262     double cotLocal = dfcot_.digi(cot - setupTT_.sectorCot(binEta));
0263     return dfcot_.toUnsigned(dfcot_.integer(cotLocal));
0264   }
0265 
0266   int Setup::digiZT(double z0, double cot, int binEta) const {
0267     double zT = z0 + setupTT_.chosenRofZ() * cot;
0268     double zTLocal = dfzT_.digi(zT - setupTT_.sectorCot(binEta) * setupTT_.chosenRofZ());
0269     return dfzT_.toUnsigned(dfzT_.integer(zTLocal));
0270   }
0271 
0272   int HitPatternHelper::reducedId(int layerId) {
0273     if (hphDebug_ && (layerId > 15 || layerId < 1)) {
0274       edm::LogVerbatim("TrackTriggerHPH") << "Warning: invalid layer id !";
0275     }
0276     if (layerId <= 6) {
0277       layerId = layerId - 1;
0278       return layerId;
0279     } else {
0280       layerId = layerId - 5;
0281       return layerId;
0282     }
0283   }
0284 
0285   int HitPatternHelper::findLayer(int layerId) {
0286     for (int i = 0; i < (int)layerEncoding_.size(); i++) {
0287       if (layerId == (int)layerEncoding_[i]) {
0288         return i;
0289       }
0290     }
0291     return -1;
0292   }
0293 
0294 }  // namespace hph