Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-03 00:12:20

0001 // This is a helper function that can be used to decode hitpattern, which is a 7-bit integer produced by the Kalman filter (KF).
0002 // Hitpattern is stored at TTTrack objects (DataFormats/L1TrackTrigger/interface/TTTrack.h)
0003 // It can be accessed via TTTrack: hitPattern()
0004 //
0005 // There are two classes declared in HitPatternHelper (hph) namesapce:
0006 // 1)Setup: This is used to produce a layermap and a collection of <tt::SensorModule> needed by HitPatternHelper.
0007 // 2)HitPatternHelper: This function returns more specific information (e.g. module type, layer id,...) about each stub on the TTTrack objects.
0008 // This function needs three variables from TTTrack: hitPattern(),tanL() and z0().
0009 // It makes predictions in two different ways depending on which version of the KF is deployed:
0010 //
0011 // Old KF (L1Trigger/TrackFindingTMTT/plugins/TMTrackProducer.cc) is a CMSSW emulation of KF firmware that
0012 // ignores truncation effect and is not bit/clock cycle accurate.
0013 // With this version of the KF, HitPatternHelper relys on a hard-coded layermap to deterimne layer IDs of each stub.
0014 //
0015 // New KF (L1Trigger/TrackerTFP/plugins/ProducerKF.cc) is a new CMSSW emulation of KF firmware that
0016 // considers truncaton effect and is bit/clock cycle accurate.
0017 // With this version of the KF, HitPatternHelper makes predictions based on the spatial coordinates of tracks and detector modules.
0018 //
0019 //  Created by J.Li on 1/23/21.
0020 //
0021 
0022 #ifndef L1Trigger_TrackTrigger_interface_HitPatternHelper_h
0023 #define L1Trigger_TrackTrigger_interface_HitPatternHelper_h
0024 
0025 #include "FWCore/Framework/interface/data_default_record_trait.h"
0026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0027 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0028 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0029 #include "Geometry/CommonTopologies/interface/PixelGeomDetUnit.h"
0030 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0031 #include "L1Trigger/TrackFindingTracklet/interface/HitPatternHelperRcd.h"
0032 #include "L1Trigger/TrackTrigger/interface/Setup.h"
0033 #include "L1Trigger/TrackerTFP/interface/DataFormats.h"
0034 #include "L1Trigger/TrackerTFP/interface/LayerEncoding.h"
0035 
0036 #include <bitset>
0037 #include <iostream>
0038 #include <vector>
0039 #include <utility>
0040 #include <map>
0041 
0042 namespace hph {
0043 
0044   //Class that stores configuration for HitPatternHelper
0045   class Setup {
0046   public:
0047     // Configuration
0048     struct Config {
0049       bool hphDebug_;
0050       bool useNewKF_;
0051       double chosenRofZ_;
0052       std::vector<double> etaRegions_;
0053     };
0054     Setup(const Config& iConfig,
0055           const tt::Setup& setupTT,
0056           const trackerTFP::DataFormats& dataFormats,
0057           const trackerTFP::LayerEncoding& layerEncoding);
0058     ~Setup() {}
0059 
0060     bool hphDebug() const { return hphDebug_; }
0061     bool useNewKF() const { return useNewKF_; }
0062     double chosenRofZ() const { return chosenRofZ_; }
0063     const std::vector<double>& etaRegions() const { return etaRegions_; }
0064     const std::map<int, std::map<int, std::vector<int>>>& layermap() const { return layermap_; }
0065     int nKalmanLayers() const { return nKalmanLayers_; }
0066     int etaRegion(double z0, double cot, bool useNewKF) const;
0067     const std::vector<int>& layerEncoding(double zT) const { return layerEncoding_->layerEncoding(zT); }
0068     // LayerEncoding call filling numPS, num2S, numMissingPS and numMissingPS for given hitPattern and trajectory
0069     void analyze(
0070         int hitpattern, double cot, double z0, int& numPS, int& num2S, int& numMissingPS, int& numMissing2S) const {
0071       layerEncoding_->analyze(hitpattern, cot, z0, numPS, num2S, numMissingPS, numMissing2S);
0072     }
0073 
0074   private:
0075     const tt::Setup* setupTT_;  // Helper class to store TrackTrigger configuration
0076     const trackerTFP::LayerEncoding* layerEncoding_;
0077     bool hphDebug_;
0078     bool useNewKF_;
0079     double chosenRofZNewKF_;
0080     std::vector<double> etaRegionsNewKF_;
0081     double chosenRofZ_;
0082     std::vector<double> etaRegions_;
0083     std::map<int, std::map<int, std::vector<int>>> layermap_;  // Hard-coded layermap in Old KF
0084     int nEtaRegions_;                                          // # of eta regions
0085     int nKalmanLayers_;                                        // # of maximum KF layers allowed
0086   };  // Only needed by Old KF
0087 
0088   //Class that returns decoded information from hitpattern
0089   class HitPatternHelper {
0090   public:
0091     HitPatternHelper(const Setup* setup, int hitpattern, double cot, double z0);
0092     ~HitPatternHelper() {}
0093 
0094     int etaSector() { return etaSector_; }      //Eta sectors defined in KF
0095     int numExpLayer() { return numExpLayer_; }  //The number of layers KF expects
0096     int numMissingPS() {
0097       return numMissingPS_;
0098     }  //The number of PS layers that are missing. It includes layers that are missing:
0099     //1)before the innermost stub on the track,
0100     //2)after the outermost stub on the track.
0101     int numMissing2S() {
0102       return numMissing2S_;
0103     }  //The number of 2S layers that are missing. It includes the two types of layers mentioned above.
0104     int numPS() { return numPS_; }  //The number of PS layers are found in hitpattern
0105     int num2S() { return num2S_; }  //The number of 2S layers are found in hitpattern
0106     int numMissingInterior1() {
0107       return numMissingInterior1_;
0108     }  //The number of missing interior layers (using only hitpattern)
0109     int numMissingInterior2() {
0110       return numMissingInterior2_;
0111     }  //The number of missing interior layers (using hitpattern, layermap from Old KF and sensor modules)
0112     const std::vector<int>& binary() {
0113       return binary_;
0114     }  //11-bit hitmask needed by TrackQuality.cc (0~5->L1~L6;6~10->D1~D5)
0115     const std::vector<float>& bonusFeatures() { return bonusFeatures_; }  //bonus features for track quality
0116 
0117     int reducedId(
0118         int layerId);  //Converts layer ID (1~6->L1~L6;11~15->D1~D5) to reduced layer ID (0~5->L1~L6;6~10->D1~D5)
0119     int findLayer(int layerId);  //Search for a layer ID from sensor modules
0120 
0121   private:
0122     const Setup* setup_;
0123     bool hphDebug_;
0124     bool useNewKF_;
0125     std::vector<double> etaRegions_;
0126     std::map<int, std::map<int, std::vector<int>>> layermap_;
0127     int nKalmanLayers_;
0128     double zT_;
0129     std::vector<int> layerEncoding_;
0130     int numExpLayer_;
0131     int hitpattern_;
0132     int etaSector_;
0133     int numMissingLayer_;
0134     int numMissingPS_;
0135     int numMissing2S_;
0136     int numPS_;
0137     int num2S_;
0138     int numMissingInterior1_;
0139     int numMissingInterior2_;
0140     std::vector<int> binary_;
0141     std::vector<float> bonusFeatures_;
0142   };
0143 
0144 }  // namespace hph
0145 
0146 EVENTSETUP_DATA_DEFAULT_RECORD(hph::Setup, hph::SetupRcd);
0147 
0148 #endif