File indexing completed on 2023-10-25 09:55:49
0001 #ifndef PHASE2GMT_TRACKCONVERTER
0002 #define PHASE2GMT_TRACKCONVERTER
0003
0004 #include "ConvertedTTTrack.h"
0005 #include "Constants.h"
0006 #include "DataFormats/L1TMuonPhase2/interface/TrackerMuon.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008
0009 namespace Phase2L1GMT {
0010
0011 class TrackConverter {
0012 public:
0013 TrackConverter(const edm::ParameterSet& iConfig) : verbose_(iConfig.getParameter<int>("verbose")) {}
0014 ~TrackConverter() {}
0015
0016 std::vector<ConvertedTTTrack> convertTracks(const std::vector<edm::Ptr<l1t::TrackerMuon::L1TTTrackType> >& tracks) {
0017 std::vector<ConvertedTTTrack> out;
0018 out.reserve(tracks.size());
0019 for (const auto& t : tracks)
0020 out.push_back(convert(t));
0021 return out;
0022 }
0023
0024 private:
0025 int verbose_;
0026 typedef ap_uint<96> wordtype;
0027
0028 uint generateQuality(const edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_> >& track) { return 1; }
0029
0030 uint ptLookup(uint absCurv) {
0031 for (auto i : ptShifts) {
0032 if (absCurv >= uint(i[0]) && absCurv < uint(i[1])) {
0033 if (i[2] < 0)
0034 return i[4];
0035 else
0036 return (absCurv >> i[2]) + i[3];
0037 }
0038 }
0039 return 0;
0040 }
0041
0042 uint etaLookup(uint absTanL) {
0043 for (auto i : etaShifts) {
0044 if (absTanL >= uint(i[0]) && absTanL < uint(i[1])) {
0045 if (i[2] < 0)
0046 return i[4];
0047 else
0048 return (absTanL >> i[2]) + i[3];
0049 }
0050 }
0051 return 0;
0052 }
0053
0054 ConvertedTTTrack convert(const edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_> >& track) {
0055 uint charge = (track->rInv() < 0) ? 1 : 0;
0056 int curvature = track->rInv() * (1 << (BITSTTCURV - 1)) / maxCurv_;
0057 int phi = track->phi() * (1 << (BITSPHI - 1)) / (M_PI);
0058 int tanLambda = track->tanL() * (1 << (BITSTTTANL - 1)) / maxTanl_;
0059 int z0 = track->z0() * (1 << (BITSZ0 - 1)) / maxZ0_;
0060 int d0 = track->d0() * (1 << (BITSD0 - 1)) / maxD0_;
0061
0062 uint absCurv = curvature > 0 ? (curvature) : (-curvature);
0063 uint pt = ptLUT[ptLookup(absCurv)];
0064 uint quality = generateQuality(track);
0065 uint absTanL = tanLambda > 0 ? (tanLambda) : (-tanLambda);
0066 uint absEta = etaLUT[etaLookup(absTanL)];
0067 int eta = tanLambda > 0 ? (absEta) : (-absEta);
0068
0069 ap_int<BITSPHI> phiSec = ap_int<BITSPHI>(phi) -
0070 ap_int<BITSPHI>((track->phiSector() * 40 * M_PI / 180.) * (1 << (BITSPHI - 1)) / (M_PI));
0071 ap_int<BITSPHI> phiCorrected = ap_int<BITSPHI>(phiSec + track->phiSector() * 910);
0072
0073 wordtype word = 0;
0074 int bstart = 0;
0075 bstart = wordconcat<wordtype>(word, bstart, curvature, BITSTTCURV);
0076 bstart = wordconcat<wordtype>(word, bstart, phiSec, BITSTTPHI);
0077 bstart = wordconcat<wordtype>(word, bstart, tanLambda, BITSTTTANL);
0078 bstart = wordconcat<wordtype>(word, bstart, z0, BITSZ0);
0079 bstart = wordconcat<wordtype>(word, bstart, d0, BITSD0);
0080 bstart = wordconcat<wordtype>(word, bstart, uint(track->chi2()), 4);
0081
0082 ConvertedTTTrack convertedTrack(charge, curvature, absEta, pt, eta, phiCorrected.to_int(), z0, d0, quality, word);
0083 convertedTrack.setOfflineQuantities(track->momentum().transverse(), track->eta(), track->phi());
0084 if (verbose_)
0085 convertedTrack.print();
0086 convertedTrack.setTrkPtr(track);
0087 return convertedTrack;
0088 }
0089 };
0090 }
0091
0092 #endif