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 
33 
34 #include <bitset>
35 #include <iostream>
36 #include <vector>
37 #include <utility>
38 #include <map>
39 
40 namespace hph {
41 
42  //Class that stores configuration for HitPatternHelper
43  class Setup {
44  public:
45  Setup() {}
46  Setup(const edm::ParameterSet& iConfig, const tt::Setup& setupTT);
47  ~Setup() {}
48 
49  bool hphDebug() const { return iConfig_.getParameter<bool>("hphDebug"); }
50  bool useNewKF() const { return iConfig_.getParameter<bool>("useNewKF"); }
51  double deltaTanL() const { return iConfig_.getParameter<double>("deltaTanL"); }
52  double chosenRofZ() const { return setupTT_.chosenRofZ(); }
53  std::vector<double> etaRegions() const { return setupTT_.boundarieEta(); }
54  std::vector<tt::SensorModule> sensorModules() const { return setupTT_.sensorModules(); }
55  std::map<int, std::map<int, std::vector<int>>> layermap() const { return layermap_; }
56  int nKalmanLayers() const { return nKalmanLayers_; }
57  static auto smallerID(std::pair<int, bool> lhs, std::pair<int, bool> rhs) { return lhs.first < rhs.first; }
58  static auto equalID(std::pair<int, bool> lhs, std::pair<int, bool> rhs) { return lhs.first == rhs.first; }
59 
60  private:
62  const tt::Setup setupTT_; // Helper class to store TrackTrigger configuration
63  std::vector<std::pair<int, bool>>
64  layerIds_; // layer IDs (1~6->L1~L6;11~15->D1~D5) and whether or not they are from tracker barrel
65  // Only needed by Old KF
66  std::map<int, std::map<int, std::vector<int>>> layermap_; // Hard-coded layermap in Old KF
67  int nEtaRegions_; // # of eta regions
68  int nKalmanLayers_; // # of maximum KF layers allowed
69  }; // Only needed by Old KF
70 
71  //Class that returns decoded information from hitpattern
73  public:
75  HitPatternHelper(const Setup* setup, int hitpattern, double cot, double z0);
77 
78  int etaSector() { return etaSector_; } //Eta sectors defined in KF
79  int numExpLayer() { return numExpLayer_; } //The number of layers KF expects
80  int numMissingPS() {
81  return numMissingPS_;
82  } //The number of PS layers that are missing. It includes layers that are missing:
83  //1)before the innermost stub on the track,
84  //2)after the outermost stub on the track.
85  int numMissing2S() {
86  return numMissing2S_;
87  } //The number of 2S layers that are missing. It includes the two types of layers mentioned above.
88  int numPS() { return numPS_; } //The number of PS layers are found in hitpattern
89  int num2S() { return num2S_; } //The number of 2S layers are found in hitpattern
91  return numMissingInterior1_;
92  } //The number of missing interior layers (using only hitpattern)
94  return numMissingInterior2_;
95  } //The number of missing interior layers (using hitpattern, layermap from Old KF and sensor modules)
96  std::vector<int> binary() { return binary_; } //11-bit hitmask needed by TrackQuality.cc (0~5->L1~L6;6~10->D1~D5)
97  static auto smallerID(tt::SensorModule lhs, tt::SensorModule rhs) { return lhs.layerId() < rhs.layerId(); }
98  static auto equalID(tt::SensorModule lhs, tt::SensorModule rhs) { return lhs.layerId() == rhs.layerId(); }
99 
100  int ReducedId(
101  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)
102  int findLayer(int layerId); //Search for a layer ID from sensor modules
103 
104  private:
111  int numPS_;
112  int num2S_;
115  double cot_;
116  double z0_;
117  const Setup* setup_;
118  std::vector<tt::SensorModule> layers_; //Sensor modules that particles are expected to hit
119  std::vector<int> binary_;
120  bool hphDebug_;
121  bool useNewKF_;
122  float chosenRofZ_;
123  float deltaTanL_; // Uncertainty added to tanL (cot) when layermap in new KF is determined
124  std::vector<double> etaRegions_;
126  std::map<int, std::map<int, std::vector<int>>> layermap_;
127  };
128 
129 } // namespace hph
130 
132 
133 #endif
edm::ParameterSet iConfig_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
int nKalmanLayers() const
double chosenRofZ() const
Definition: Setup.h:413
Class to process and provide run-time constants used by Track Trigger emulators.
Definition: Setup.h:44
static auto equalID(tt::SensorModule lhs, tt::SensorModule rhs)
bool useNewKF() const
double boundarieEta(int eta) const
Definition: Setup.h:417
double deltaTanL() const
double chosenRofZ() const
bool hphDebug() const
std::vector< int > binary_
#define EVENTSETUP_DATA_DEFAULT_RECORD(_data_, _record_)
static auto smallerID(tt::SensorModule lhs, tt::SensorModule rhs)
std::map< int, std::map< int, std::vector< int > > > layermap() const
const std::vector< SensorModule > & sensorModules() const
Definition: Setup.h:126
std::vector< double > etaRegions_
std::vector< tt::SensorModule > layers_
static auto smallerID(std::pair< int, bool > lhs, std::pair< int, bool > rhs)
std::vector< std::pair< int, bool > > layerIds_
int findLayer(int layerId)
std::vector< int > binary()
std::vector< tt::SensorModule > sensorModules() const
const tt::Setup setupTT_
std::vector< double > etaRegions() const
int ReducedId(int layerId)
int layerId() const
Definition: SensorModule.h:43
std::map< int, std::map< int, std::vector< int > > > layermap_
static auto equalID(std::pair< int, bool > lhs, std::pair< int, bool > rhs)
std::map< int, std::map< int, std::vector< int > > > layermap_