CMS 3D CMS Logo

codecs_tracks.cc
Go to the documentation of this file.
3 
4 namespace l1t::demo::codecs {
5 
6  // Return true if a track is contained within a collection
10  auto it = std::find_if(
11  trackRefCollection->begin(),
12  trackRefCollection->end(),
13  [&trackRef](edm::Ref<std::vector<TTTrack<Ref_Phase2TrackerDigi_>>> const& obj) { return obj == trackRef; });
14  if (it != trackRefCollection->end())
15  return true;
16  else
17  return false;
18  }
19 
20  // Encodes a single track into a 96-bit track word
21  ap_uint<96> encodeTrack(const TTTrack_TrackWord& t) { return t.getTrackWord(); }
22 
23  // Return the 96-bit track words from a given track collection and place them on the appropriate 18 'logical' links
24  std::array<std::vector<ap_uint<96>>, 18> getTrackWords(const edm::View<TTTrack<Ref_Phase2TrackerDigi_>>& tracks) {
25  std::array<std::vector<ap_uint<96>>, 18> trackWords;
26  for (const auto& track : tracks) {
27  trackWords.at(gttLinkID(track)).push_back(encodeTrack(track));
28  }
29  return trackWords;
30  }
31 
32  // Return the 96-bit track words from a given track collection and place them on the appropriate 18 'logical' links
33  std::array<std::vector<ap_uint<96>>, 18> getTrackWords(
36  std::array<std::vector<ap_uint<96>>, 18> trackWords;
37  for (unsigned int itrack = 0; itrack < referenceTracks->size(); itrack++) {
38  const auto& referenceTrack = referenceTracks->at(itrack);
39  edm::Ref<std::vector<TTTrack<Ref_Phase2TrackerDigi_>>> referenceTrackRef(referenceTracks, itrack);
40 
41  if (trackInCollection(referenceTrackRef, tracks)) {
42  trackWords.at(gttLinkID(referenceTrack)).push_back(encodeTrack(referenceTrack));
43  } else {
44  trackWords.at(gttLinkID(referenceTrack)).push_back(ap_uint<96>(0));
45  }
46  }
47  return trackWords;
48  }
49 
50  // Encodes a set of tracks onto a set of links
51  size_t encodeLinks(std::array<std::vector<ap_uint<96>>, 18>& trackWords,
52  std::array<std::vector<ap_uint<64>>, 18>& linkData) {
53  size_t counter = 0;
54  for (size_t i = 0; i < linkData.size(); i++) {
55  // Pad track vectors -> full packet length (156 frames = 104 tracks)
56  trackWords.at(i).resize(104, 0);
57  linkData.at(i).resize(156, {0});
58 
59  for (size_t j = 0; (j < trackWords.at(i).size()); j += 2) {
60  linkData.at(i).at(3 * j / 2) = trackWords.at(i).at(j)(63, 0);
61  linkData.at(i).at(3 * j / 2 + 1) =
62  (ap_uint<32>(trackWords.at(i).at(j + 1)(31, 0)), ap_uint<32>(trackWords.at(i).at(j)(95, 64)));
63  linkData.at(i).at(3 * j / 2 + 2) = trackWords.at(i).at(j + 1)(95, 32);
64  counter += trackWords.at(i).at(j)(95, 95) + trackWords.at(i).at(j + 1)(95, 95);
65  }
66  }
67  return counter;
68  }
69 
70  // Encodes track collection onto 18 output links (2x9 eta-phi sectors; , -/+ eta pairs)
71  std::array<std::vector<ap_uint<64>>, 18> encodeTracks(const edm::View<TTTrack<Ref_Phase2TrackerDigi_>>& tracks,
72  int debug) {
73  if (debug > 0) {
74  edm::LogInfo("l1t::demo::codecs") << "encodeTrack::Encoding " << tracks.size() << " tracks";
75  }
76 
77  std::array<std::vector<ap_uint<96>>, 18> trackWords = getTrackWords(tracks);
78  std::array<std::vector<ap_uint<64>>, 18> linkData;
79  size_t counter = encodeLinks(trackWords, linkData);
80 
81  if (debug > 0) {
82  edm::LogInfo("l1t::demo::codecs") << "encodeTrack::Encoded " << counter << " tracks";
83  }
84 
85  return linkData;
86  }
87 
88  // Encodes a track collection based off the ordering of another track collection
89  // Requirement: The second collection must be a subset of the first
90  std::array<std::vector<ap_uint<64>>, 18> encodeTracks(
93  int debug) {
94  if (debug > 0) {
95  edm::LogInfo("l1t::demo::codecs") << "encodeTrack::Encoding " << tracks->size() << " tracks";
96  }
97 
98  std::array<std::vector<ap_uint<96>>, 18> trackWords = getTrackWords(referenceTracks, tracks);
99  std::array<std::vector<ap_uint<64>>, 18> linkData;
100  size_t counter = encodeLinks(trackWords, linkData);
101 
102  if (debug > 0) {
103  edm::LogInfo("l1t::demo::codecs") << "encodeTrack::Encoded " << counter << " tracks";
104  }
105 
106  return linkData;
107  }
108 
109  std::vector<TTTrack_TrackWord> decodeTracks(const std::vector<ap_uint<64>>& frames) {
110  std::vector<TTTrack_TrackWord> tracks;
111 
112  if ((frames.size() % 3) != 0) {
113  std::stringstream message;
114  message << "The number of track frames (" << frames.size() << ") is not evenly divisible by 3";
115  throw std::runtime_error(message.str());
116  }
117 
118  for (size_t i = 0; i < frames.size(); i += 3) {
119  TTTrack_TrackWord::tkword_t combination1 = (ap_uint<32>(frames.at(i + 1)(31, 0)), frames.at(i)(63, 0));
120  TTTrack_TrackWord::tkword_t combination2 = (frames.at(i + 2)(63, 0), ap_uint<32>(frames.at(i + 1)(63, 32)));
121  TTTrack_TrackWord track1, track2;
122  track1.setTrackWord(
137  track2.setTrackWord(
152  tracks.push_back(track1);
153  tracks.push_back(track2);
154  }
155 
156  return tracks;
157  }
158 
159  // Decodes the tracks from 18 'logical' output links (2x9 eta-phi sectors; , -/+ eta pairs)
160  std::array<std::vector<TTTrack_TrackWord>, 18> decodeTracks(const std::array<std::vector<ap_uint<64>>, 18>& frames) {
161  std::array<std::vector<TTTrack_TrackWord>, 18> tracks;
162 
163  for (size_t i = 0; i < tracks.size(); i++) {
164  tracks.at(i) = decodeTracks(frames.at(i));
165  }
166 
167  return tracks;
168  }
169 
170 } // namespace l1t::demo::codecs
ap_uint< 96 > encodeTrack(const TTTrack_TrackWord &t)
ap_uint< TrackBitWidths::kMVAQualitySize > qualityMVA_t
size_t encodeLinks(std::array< std::vector< ap_uint< 96 >>, 18 > &trackWords, std::array< std::vector< ap_uint< 64 >>, 18 > &linkData)
ap_uint< kTrackWordSize > tkword_t
std::array< std::vector< ap_uint< 64 > >, 18 > encodeTracks(const edm::View< TTTrack< Ref_Phase2TrackerDigi_ >> &, int debug=0)
ap_uint< TrackBitWidths::kZ0Size > z0_t
void setTrackWord(unsigned int valid, const GlobalVector &momentum, const GlobalPoint &POCA, double rInv, double chi2RPhi, double chi2RZ, double bendChi2, unsigned int hitPattern, double mvaQuality, double mvaOther, unsigned int sector)
ap_uint< TrackBitWidths::kMVAOtherSize > otherMVA_t
unsigned int gttLinkID(T track)
Definition: tracks.h:20
ap_uint< TrackBitWidths::kValidSize > valid_t
ap_uint< TrackBitWidths::kRinvSize > rinv_t
ap_uint< TrackBitWidths::kHitPatternSize > hit_t
ap_uint< TrackBitWidths::kChi2RZSize > chi2rz_t
Log< level::Info, false > LogInfo
#define debug
Definition: HDRShower.cc:19
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
std::array< std::vector< ap_uint< 96 > >, 18 > getTrackWords(const edm::View< TTTrack< Ref_Phase2TrackerDigi_ >> &)
bool trackInCollection(const edm::Ref< std::vector< TTTrack< Ref_Phase2TrackerDigi_ >>> &, const edm::Handle< edm::RefVector< std::vector< TTTrack< Ref_Phase2TrackerDigi_ >>>> &)
Definition: codecs_tracks.cc:7
ap_uint< TrackBitWidths::kTanlSize > tanl_t
static std::atomic< unsigned int > counter
ap_uint< TrackBitWidths::kD0Size > d0_t
ap_uint< TrackBitWidths::kChi2RPhiSize > chi2rphi_t
std::vector< TTTrack_TrackWord > decodeTracks(const std::vector< ap_uint< 64 >> &)
ap_uint< TrackBitWidths::kBendChi2Size > bendChi2_t
ap_uint< TrackBitWidths::kPhiSize > phi_t