CMS 3D CMS Logo

TrackConverter.cc
Go to the documentation of this file.
2 
3 using namespace Phase2L1GMT;
4 
5 TrackConverter::TrackConverter(const edm::ParameterSet& iConfig) : verbose_(iConfig.getParameter<int>("verbose")) {}
6 
7 std::vector<ConvertedTTTrack> TrackConverter::convertTracks(
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 }
15 
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  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 
50  return convertedTrack;
51 }
const float LSBphi
Definition: Constants.h:93
std::vector< ConvertedTTTrack > convertTracks(const std::vector< edm::Ptr< l1t::TrackerMuon::L1TTTrackType > > &tracks)
TrackConverter(const edm::ParameterSet &iConfig)
const ap_uint< BITSPT > ptLUT[2251]
Definition: TPSLUTs.h:18
ConvertedTTTrack convert(const edm::Ptr< TTTrack< Ref_Phase2TrackerDigi_ > > &track)
const int BITSZ0
Definition: Constants.h:27
void setTrkPtr(const edm::Ptr< TTTrack< Ref_Phase2TrackerDigi_ > > &trkPtr)
T curvature(T InversePt, const MagneticField &field)
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
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
const int BITSTTTANL
Definition: Constants.h:13
void setOfflineQuantities(float pt, float eta, float phi)
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