CMS 3D CMS Logo

HitPatternHelper.h
Go to the documentation of this file.
1 // This is a helper function that can be used to decode hitpattern, which is a 7-bit integer produced by the Kalman filter (KF).
2 // Hitpattern is stored at TTTrack objects (DataFormats/L1TrackTrigger/interface/TTTrack.h)
3 // It can be accessed via TTTrack: hitPattern()
4 //
5 // There are two classes declared in HitPatternHelper (hph) namesapce:
6 // 1)Setup: This is used to produce a layermap and a collection of <tt::SensorModule> needed by HitPatternHelper.
7 // 2)HitPatternHelper: This function returns more specific information (e.g. module type, layer id,...) about each stub on the TTTrack objects.
8 // This function needs three variables from TTTrack: hitPattern(),tanL() and z0().
9 // It makes predictions in two different ways depending on which version of the KF is deployed:
10 //
11 // Old KF (L1Trigger/TrackFindingTMTT/plugins/TMTrackProducer.cc) is a CMSSW emulation of KF firmware that
12 // ignores truncation effect and is not bit/clock cycle accurate.
13 // With this version of the KF, HitPatternHelper relys on a hard-coded layermap to deterimne layer IDs of each stub.
14 //
15 // New KF (L1Trigger/TrackerTFP/plugins/ProducerKF.cc) is a new CMSSW emulation of KF firmware that
16 // considers truncaton effect and is bit/clock cycle accurate.
17 // With this version of the KF, HitPatternHelper makes predictions based on the spatial coordinates of tracks and detector modules.
18 //
19 // Created by J.Li on 1/23/21.
20 //
21 
22 #ifndef L1Trigger_TrackTrigger_interface_HitPatternHelper_h
23 #define L1Trigger_TrackTrigger_interface_HitPatternHelper_h
24 
35 
36 #include <bitset>
37 #include <iostream>
38 #include <vector>
39 #include <utility>
40 #include <map>
41 
42 namespace hph {
43 
44  //Class that stores configuration for HitPatternHelper
45  class Setup {
46  public:
47  Setup(const edm::ParameterSet& iConfig,
48  const tt::Setup& setupTT,
49  const trackerTFP::DataFormats& dataFormats,
51  ~Setup() {}
52 
53  bool hphDebug() const { return hphDebug_; }
54  bool useNewKF() const { return useNewKF_; }
55  double chosenRofZ() const { return chosenRofZ_; }
56  std::vector<double> etaRegions() const { return etaRegions_; }
57  std::map<int, std::map<int, std::vector<int>>> layermap() const { return layermap_; }
58  int nKalmanLayers() const { return nKalmanLayers_; }
59  int etaRegion(double z0, double cot, bool useNewKF) const;
60  int digiCot(double cot, int binEta) const;
61  int digiZT(double z0, double cot, int binEta) const;
62  const std::vector<int>& layerEncoding(int binEta, int binZT, int binCot) const {
63  return layerEncoding_.layerEncoding(binEta, binZT, binCot);
64  }
65  const std::map<int, const tt::SensorModule*>& layerEncodingMap(int binEta, int binZT, int binCot) const {
66  return layerEncoding_.layerEncodingMap(binEta, binZT, binCot);
67  }
68 
69  private:
72  const tt::Setup setupTT_; // Helper class to store TrackTrigger configuration
77  bool hphDebug_;
78  bool useNewKF_;
80  std::vector<double> etaRegionsNewKF_;
81  double chosenRofZ_;
82  std::vector<double> etaRegions_;
83  std::map<int, std::map<int, std::vector<int>>> layermap_; // Hard-coded layermap in Old KF
84  int nEtaRegions_; // # of eta regions
85  int nKalmanLayers_; // # of maximum KF layers allowed
86  }; // Only needed by Old KF
87 
88  //Class that returns decoded information from hitpattern
90  public:
91  HitPatternHelper(const Setup* setup, int hitpattern, double cot, double z0);
93 
94  int etaSector() { return etaSector_; } //Eta sectors defined in KF
95  int numExpLayer() { return numExpLayer_; } //The number of layers KF expects
96  int numMissingPS() {
97  return numMissingPS_;
98  } //The number of PS layers that are missing. It includes layers that are missing:
99  //1)before the innermost stub on the track,
100  //2)after the outermost stub on the track.
101  int numMissing2S() {
102  return numMissing2S_;
103  } //The number of 2S layers that are missing. It includes the two types of layers mentioned above.
104  int numPS() { return numPS_; } //The number of PS layers are found in hitpattern
105  int num2S() { return num2S_; } //The number of 2S layers are found in hitpattern
107  return numMissingInterior1_;
108  } //The number of missing interior layers (using only hitpattern)
110  return numMissingInterior2_;
111  } //The number of missing interior layers (using hitpattern, layermap from Old KF and sensor modules)
112  std::vector<int> binary() { return binary_; } //11-bit hitmask needed by TrackQuality.cc (0~5->L1~L6;6~10->D1~D5)
113  std::vector<float> bonusFeatures() { return bonusFeatures_; } //bonus features for track quality
114 
115  int reducedId(
116  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)
117  int findLayer(int layerId); //Search for a layer ID from sensor modules
118 
119  private:
120  const Setup* setup_;
121  bool hphDebug_;
122  bool useNewKF_;
123  std::vector<double> etaRegions_;
124  std::map<int, std::map<int, std::vector<int>>> layermap_;
126  int etaBin_;
127  int cotBin_;
128  int zTBin_;
129  std::vector<int> layerEncoding_;
130  std::map<int, const tt::SensorModule*> layerEncodingMap_;
137  int numPS_;
138  int num2S_;
141  std::vector<int> binary_;
142  std::vector<float> bonusFeatures_;
143  };
144 
145 } // namespace hph
146 
148 
149 #endif
std::vector< float > bonusFeatures()
edm::ParameterSet iConfig_
int nKalmanLayers() const
double chosenRofZ_
Class to process and provide run-time constants used by Track Trigger emulators.
Definition: Setup.h:44
bool useNewKF() const
std::vector< float > bonusFeatures_
const std::vector< int > & layerEncoding(int binEta, int binZT, int binCot) const
Definition: LayerEncoding.h:25
const trackerTFP::DataFormat dfcot_
double chosenRofZ() const
double chosenRofZNewKF_
int digiZT(double z0, double cot, int binEta) const
bool hphDebug() const
const std::map< int, const tt::SensorModule * > & layerEncodingMap(int binEta, int binZT, int binCot) const
const trackerTFP::LayerEncoding layerEncoding_
const trackerTFP::DataFormat dfzT_
edm::ParameterSet oldKFPSet_
std::vector< int > binary_
std::map< int, const tt::SensorModule * > layerEncodingMap_
#define EVENTSETUP_DATA_DEFAULT_RECORD(_data_, _record_)
std::map< int, std::map< int, std::vector< int > > > layermap() const
std::vector< double > etaRegions_
std::vector< double > etaRegionsNewKF_
const trackerTFP::DataFormats dataFormats_
std::vector< int > layerEncoding_
const std::map< int, const tt::SensorModule * > & layerEncodingMap(int binEta, int binZT, int binCot) const
Definition: LayerEncoding.h:28
int findLayer(int layerId)
std::vector< int > binary()
const tt::Setup setupTT_
std::vector< double > etaRegions() const
std::vector< double > etaRegions_
Class to encode layer ids for Kalman Filter Layers consitent with rough r-z track parameters are coun...
Definition: LayerEncoding.h:19
const std::vector< int > & layerEncoding(int binEta, int binZT, int binCot) const
Class to calculate and provide dataformats used by Track Trigger emulator.
Definition: DataFormats.h:216
Setup(const edm::ParameterSet &iConfig, const tt::Setup &setupTT, const trackerTFP::DataFormats &dataFormats, const trackerTFP::LayerEncoding &layerEncoding)
int reducedId(int layerId)
int digiCot(double cot, int binEta) const
std::map< int, std::map< int, std::vector< int > > > layermap_
HitPatternHelper(const Setup *setup, int hitpattern, double cot, double z0)
int etaRegion(double z0, double cot, bool useNewKF) const
std::map< int, std::map< int, std::vector< int > > > layermap_