CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
Phase2L1GMT::TrackConverter Class Reference

#include <TrackConverter.h>

Public Member Functions

std::vector< ConvertedTTTrackconvertTracks (const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &tracks)
 
 TrackConverter (const edm::ParameterSet &iConfig)
 
 ~TrackConverter ()=default
 

Private Types

typedef ap_uint< 96 > wordtype
 

Private Member Functions

ConvertedTTTrack convert (const edm::Ptr< TTTrack< Ref_Phase2TrackerDigi_ > > &track)
 
uint etaLookup (uint absTanL)
 
uint generateQuality (const edm::Ptr< TTTrack< Ref_Phase2TrackerDigi_ > > &track)
 
uint ptLookup (uint absCurv)
 

Private Attributes

int verbose_
 

Detailed Description

Definition at line 11 of file TrackConverter.h.

Member Typedef Documentation

◆ wordtype

typedef ap_uint<96> Phase2L1GMT::TrackConverter::wordtype
private

Definition at line 20 of file TrackConverter.h.

Constructor & Destructor Documentation

◆ TrackConverter()

TrackConverter::TrackConverter ( const edm::ParameterSet iConfig)

Definition at line 5 of file TrackConverter.cc.

5 : verbose_(iConfig.getParameter<int>("verbose")) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:307

◆ ~TrackConverter()

Phase2L1GMT::TrackConverter::~TrackConverter ( )
default

Member Function Documentation

◆ convert()

ConvertedTTTrack TrackConverter::convert ( const edm::Ptr< TTTrack< Ref_Phase2TrackerDigi_ > > &  track)
private

Definition at line 16 of file TrackConverter.cc.

References Phase2L1GMT::BITSD0, Phase2L1GMT::BITSETA, Phase2L1GMT::BITSPHI, Phase2L1GMT::BITSTTCURV, Phase2L1GMT::BITSTTD0, Phase2L1GMT::BITSTTTANL, Phase2L1GMT::BITSTTZ0, Phase2L1GMT::BITSZ0, ALCARECOTkAlJpsiMuMu_cff::charge, PixelRecoUtilities::curvature(), d0, PVValHelper::eta, etaLookup(), Phase2L1GMT::etaLUT, generateQuality(), Phase2L1GMT::LSBeta, Phase2L1GMT::LSBphi, Phase2L1GMT::LSBpt, Phase2L1GMT::ConvertedTTTrack::print(), DiDispStaMuonMonitor_cfi::pt, ptLookup(), Phase2L1GMT::ptLUT, quality, Phase2L1GMT::ConvertedTTTrack::setOfflineQuantities(), Phase2L1GMT::ConvertedTTTrack::setTrkPtr(), HLT_2024v13_cff::track, parallelization::uint, and verbose_.

Referenced by convertTracks().

16  {
17  uint charge = (track->rInv() < 0) ? 1 : 0;
18  ap_int<BITSTTCURV> curvature = ap_int<BITSTTCURV>(track->getRinvBits());
19  ap_int<BITSPHI> phisec = ap_int<BITSPHI>(ap_int<BITSTTPHI>(track->getPhiBits()) / 2);
20  ap_int<BITSTTTANL> tanLambda = ap_int<BITSTTTANL>(track->getTanlBits());
21  ap_int<BITSZ0> z0 = ap_int<BITSZ0>(ap_int<BITSTTZ0>(track->getZ0Bits()) / (1 << (BITSTTZ0 - BITSZ0)));
22  ap_int<BITSD0> d0 = ap_int<BITSD0>(ap_int<BITSTTD0>(track->getD0Bits()) / (1 << (BITSTTD0 - BITSD0)));
23  //calculate pt
24  ap_uint<BITSTTCURV - 1> absCurv =
25  curvature > 0 ? ap_uint<BITSTTCURV - 1>(curvature) : ap_uint<BITSTTCURV - 1>(-curvature);
26  ap_uint<BITSPT> pt = ptLUT[ptLookup(absCurv)];
27  ap_uint<1> quality = generateQuality(track);
28  ap_uint<BITSTTTANL - 1> absTanL =
29  tanLambda > 0 ? ap_uint<BITSTTTANL - 1>(tanLambda) : ap_uint<BITSTTTANL - 1>(-tanLambda);
30  ap_uint<BITSETA - 1> absEta = etaLUT[etaLookup(absTanL)];
31  ap_int<BITSETA> eta = tanLambda > 0 ? ap_int<BITSETA>(absEta) : ap_int<BITSETA>(-absEta);
32 
33  ap_int<BITSPHI> phi = ap_int<BITSPHI>(phisec + track->phiSector() * 910);
34 
35  wordtype word = 0;
36  int bstart = 0;
37  bstart = wordconcat<wordtype>(word, bstart, curvature, BITSTTCURV);
38  bstart = wordconcat<wordtype>(word, bstart, phi, BITSPHI); //was phiSec
39  bstart = wordconcat<wordtype>(word, bstart, tanLambda, BITSTTTANL);
40  bstart = wordconcat<wordtype>(word, bstart, z0, BITSZ0);
41  bstart = wordconcat<wordtype>(word, bstart, d0, BITSD0);
42  bstart = wordconcat<wordtype>(word, bstart, uint(track->chi2()), 4);
43 
44  ConvertedTTTrack convertedTrack(charge, curvature, absEta, pt, eta, phi, z0, d0, quality, word);
45  convertedTrack.setOfflineQuantities(LSBpt * pt, LSBeta * eta, LSBphi * phi);
46  if (verbose_)
47  convertedTrack.print();
48  convertedTrack.setTrkPtr(track);
49  return convertedTrack;
50 }
const float LSBphi
Definition: Constants.h:93
const ap_uint< BITSPT > ptLUT[2251]
Definition: TPSLUTs.h:18
const int BITSZ0
Definition: Constants.h:27
T curvature(T InversePt, const MagneticField &field)
ap_uint< 64 > wordtype
Definition: Constants.h:101
string quality
uint64_t word
const float LSBpt
Definition: Constants.h:92
const int BITSETA
Definition: Constants.h:26
uint etaLookup(uint absTanL)
const int BITSPHI
Definition: Constants.h:25
const int BITSD0
Definition: Constants.h:28
const int BITSTTCURV
Definition: Constants.h:10
uint ptLookup(uint absCurv)
const float LSBeta
Definition: Constants.h:94
static constexpr float d0
const int BITSTTTANL
Definition: Constants.h:13
uint generateQuality(const edm::Ptr< TTTrack< Ref_Phase2TrackerDigi_ > > &track)
const int BITSTTD0
Definition: Constants.h:16
const int BITSTTZ0
Definition: Constants.h:15
const uint etaLUT[4082]
Definition: TPSLUTs.h:142

◆ convertTracks()

std::vector< ConvertedTTTrack > TrackConverter::convertTracks ( const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &  tracks)

Definition at line 7 of file TrackConverter.cc.

References convert(), MillePedeFileConverter_cfg::out, submitPVValidationJobs::t, and DiMuonV_cfg::tracks.

8  {
9  std::vector<ConvertedTTTrack> out;
10  out.reserve(tracks.size());
11  for (const auto& t : tracks)
12  out.push_back(convert(t));
13  return out;
14 }
ConvertedTTTrack convert(const edm::Ptr< TTTrack< Ref_Phase2TrackerDigi_ > > &track)

◆ etaLookup()

uint Phase2L1GMT::TrackConverter::etaLookup ( uint  absTanL)
inlineprivate

Definition at line 36 of file TrackConverter.h.

References Phase2L1GMT::etaShifts, mps_fire::i, and parallelization::uint.

Referenced by convert().

36  {
37  for (auto i : etaShifts) {
38  if (absTanL >= uint(i[0]) && absTanL < uint(i[1])) {
39  if (i[2] < 0)
40  return i[4];
41  else
42  return (absTanL >> i[2]) + i[3];
43  }
44  }
45  return 0;
46  }
const int etaShifts[4][5]
Definition: TPSLUTs.h:139

◆ generateQuality()

uint Phase2L1GMT::TrackConverter::generateQuality ( const edm::Ptr< TTTrack< Ref_Phase2TrackerDigi_ > > &  track)
inlineprivate

Definition at line 22 of file TrackConverter.h.

Referenced by convert().

22 { return 1; }

◆ ptLookup()

uint Phase2L1GMT::TrackConverter::ptLookup ( uint  absCurv)
inlineprivate

Definition at line 24 of file TrackConverter.h.

References mps_fire::i, Phase2L1GMT::ptShifts, and parallelization::uint.

Referenced by convert().

24  {
25  for (auto i : ptShifts) {
26  if (absCurv >= uint(i[0]) && absCurv < uint(i[1])) {
27  if (i[2] < 0)
28  return i[4];
29  else
30  return (absCurv >> i[2]) + i[3];
31  }
32  }
33  return 0;
34  }
const int ptShifts[8][5]
Definition: TPSLUTs.h:9

Member Data Documentation

◆ verbose_

int Phase2L1GMT::TrackConverter::verbose_
private

Definition at line 19 of file TrackConverter.h.

Referenced by convert().