Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:21:45

0001 #include "L1Trigger/TrackerTFP/interface/DataFormats.h"
0002 #include "L1Trigger/TrackTrigger/interface/StubPtConsistency.h"
0003 
0004 #include <vector>
0005 #include <deque>
0006 #include <cmath>
0007 #include <tuple>
0008 #include <iterator>
0009 #include <algorithm>
0010 #include <string>
0011 #include <iostream>
0012 
0013 using namespace std;
0014 using namespace edm;
0015 using namespace tt;
0016 
0017 namespace trackerTFP {
0018 
0019   // default constructor, trying to needs space as proper constructed object
0020   DataFormats::DataFormats()
0021       : numDataFormats_(0),
0022         formats_(+Variable::end, std::vector<DataFormat*>(+Process::end, nullptr)),
0023         numUnusedBitsStubs_(+Process::end, TTBV::S_),
0024         numUnusedBitsTracks_(+Process::end, TTBV::S_),
0025         numChannel_(+Process::end, 0) {
0026     setup_ = nullptr;
0027     countFormats();
0028     dataFormats_.reserve(numDataFormats_);
0029     numStreams_.reserve(+Process::end);
0030     numStreamsStubs_.reserve(+Process::end);
0031     numStreamsTracks_.reserve(+Process::end);
0032   }
0033 
0034   // method to count number of unique data formats
0035   template <Variable v, Process p>
0036   void DataFormats::countFormats() {
0037     if constexpr (config_[+v][+p] == p)
0038       numDataFormats_++;
0039     if constexpr (++p != Process::end)
0040       countFormats<v, ++p>();
0041     else if constexpr (++v != Variable::end)
0042       countFormats<++v>();
0043   }
0044 
0045   // proper constructor
0046   DataFormats::DataFormats(const ParameterSet& iConfig, const Setup* setup) : DataFormats() {
0047     iConfig_ = iConfig;
0048     setup_ = setup;
0049     fillDataFormats();
0050     for (const Process p : Processes)
0051       for (const Variable v : stubs_[+p])
0052         numUnusedBitsStubs_[+p] -= formats_[+v][+p] ? formats_[+v][+p]->width() : 0;
0053     for (const Process p : Processes)
0054       for (const Variable v : tracks_[+p])
0055         numUnusedBitsTracks_[+p] -= formats_[+v][+p] ? formats_[+v][+p]->width() : 0;
0056     numChannel_[+Process::dtc] = setup_->numDTCsPerRegion();
0057     numChannel_[+Process::pp] = setup_->numDTCsPerTFP();
0058     numChannel_[+Process::gp] = setup_->numSectors();
0059     numChannel_[+Process::ht] = setup_->htNumBinsInv2R();
0060     numChannel_[+Process::mht] = setup_->htNumBinsInv2R();
0061     numChannel_[+Process::zht] = setup_->htNumBinsInv2R();
0062     numChannel_[+Process::kfin] = setup_->kfNumWorker() * setup_->numLayers();
0063     numChannel_[+Process::kf] = setup_->kfNumWorker();
0064     transform(numChannel_.begin(), numChannel_.end(), back_inserter(numStreams_), [this](int channel) {
0065       return channel * setup_->numRegions();
0066     });
0067     numStreamsStubs_ = numStreams_;
0068     numStreamsStubs_[+Process::kf] = numStreams_[+Process::kfin];
0069     numStreamsTracks_ = vector<int>(+Process::end, 0);
0070     numStreamsTracks_[+Process::kfin] = numStreams_[+Process::kf];
0071     numStreamsTracks_[+Process::kf] = numStreams_[+Process::kf];
0072     // Print digi data format of all variables of any specified algo step
0073     // (Look at DataFormat.h::tracks_ to see variable names).
0074     //for (const Variable v : tracks_[+Process::kf]) {
0075     //  const DataFormat& f = format(v, Process::kf);
0076     //  cout <<" KF "<< f.base() << " " << f.range() << " " << f.width() << endl;
0077     //}
0078   }
0079 
0080   // constructs data formats of all unique used variables and flavours
0081   template <Variable v, Process p>
0082   void DataFormats::fillDataFormats() {
0083     if constexpr (config_[+v][+p] == p) {
0084       dataFormats_.emplace_back(Format<v, p>(iConfig_, setup_));
0085       fillFormats<v, p>();
0086     }
0087     if constexpr (++p != Process::end)
0088       fillDataFormats<v, ++p>();
0089     else if constexpr (++v != Variable::end)
0090       fillDataFormats<++v>();
0091   }
0092 
0093   // helper (loop) data formats of all unique used variables and flavours
0094   template <Variable v, Process p, Process it>
0095   void DataFormats::fillFormats() {
0096     if (config_[+v][+it] == p) {
0097       formats_[+v][+it] = &dataFormats_.back();
0098     }
0099     if constexpr (++it != Process::end)
0100       fillFormats<v, p, ++it>();
0101   }
0102 
0103   // converts bits to ntuple of variables
0104   template <typename... Ts>
0105   void DataFormats::convertStub(Process p, const Frame& bv, tuple<Ts...>& data) const {
0106     TTBV ttBV(bv);
0107     extractStub(p, ttBV, data);
0108   }
0109 
0110   // helper (loop) to convert bits to ntuple of variables
0111   template <int it, typename... Ts>
0112   void DataFormats::extractStub(Process p, TTBV& ttBV, std::tuple<Ts...>& data) const {
0113     Variable v = *next(stubs_[+p].begin(), sizeof...(Ts) - 1 - it);
0114     formats_[+v][+p]->extract(ttBV, get<sizeof...(Ts) - 1 - it>(data));
0115     if constexpr (it + 1 != sizeof...(Ts))
0116       extractStub<it + 1>(p, ttBV, data);
0117   }
0118 
0119   // converts ntuple of variables to bits
0120   template <typename... Ts>
0121   void DataFormats::convertStub(Process p, const std::tuple<Ts...>& data, Frame& bv) const {
0122     TTBV ttBV(1, numUnusedBitsStubs_[+p]);
0123     attachStub(p, data, ttBV);
0124     bv = ttBV.bs();
0125   }
0126 
0127   // helper (loop) to convert ntuple of variables to bits
0128   template <int it, typename... Ts>
0129   void DataFormats::attachStub(Process p, const tuple<Ts...>& data, TTBV& ttBV) const {
0130     Variable v = *next(stubs_[+p].begin(), it);
0131     formats_[+v][+p]->attach(get<it>(data), ttBV);
0132     if constexpr (it + 1 != sizeof...(Ts))
0133       attachStub<it + 1>(p, data, ttBV);
0134   }
0135 
0136   // converts bits to ntuple of variables
0137   template <typename... Ts>
0138   void DataFormats::convertTrack(Process p, const Frame& bv, tuple<Ts...>& data) const {
0139     TTBV ttBV(bv);
0140     extractTrack(p, ttBV, data);
0141   }
0142 
0143   // helper (loop) to convert bits to ntuple of variables
0144   template <int it, typename... Ts>
0145   void DataFormats::extractTrack(Process p, TTBV& ttBV, std::tuple<Ts...>& data) const {
0146     Variable v = *next(tracks_[+p].begin(), sizeof...(Ts) - 1 - it);
0147     formats_[+v][+p]->extract(ttBV, get<sizeof...(Ts) - 1 - it>(data));
0148     if constexpr (it + 1 != sizeof...(Ts))
0149       extractTrack<it + 1>(p, ttBV, data);
0150   }
0151 
0152   // converts ntuple of variables to bits
0153   template <typename... Ts>
0154   void DataFormats::convertTrack(Process p, const std::tuple<Ts...>& data, Frame& bv) const {
0155     TTBV ttBV(1, numUnusedBitsTracks_[+p]);
0156     attachTrack(p, data, ttBV);
0157     bv = ttBV.bs();
0158   }
0159 
0160   // helper (loop) to convert ntuple of variables to bits
0161   template <int it, typename... Ts>
0162   void DataFormats::attachTrack(Process p, const tuple<Ts...>& data, TTBV& ttBV) const {
0163     Variable v = *next(tracks_[+p].begin(), it);
0164     formats_[+v][+p]->attach(get<it>(data), ttBV);
0165     if constexpr (it + 1 != sizeof...(Ts))
0166       attachTrack<it + 1>(p, data, ttBV);
0167   }
0168 
0169   // construct Stub from Frame
0170   template <typename... Ts>
0171   Stub<Ts...>::Stub(const FrameStub& frame, const DataFormats* dataFormats, Process p)
0172       : dataFormats_(dataFormats), p_(p), frame_(frame), trackId_(0) {
0173     dataFormats_->convertStub(p, frame.second, data_);
0174   }
0175 
0176   // construct Stub from other Stub
0177   template <typename... Ts>
0178   template <typename... Others>
0179   Stub<Ts...>::Stub(const Stub<Others...>& stub, Ts... data)
0180       : dataFormats_(stub.dataFormats()),
0181         p_(++stub.p()),
0182         frame_(stub.frame().first, Frame()),
0183         data_(data...),
0184         trackId_(0) {}
0185 
0186   // construct Stub from TTStubRef
0187   template <typename... Ts>
0188   Stub<Ts...>::Stub(const TTStubRef& ttStubRef, const DataFormats* dataFormats, Process p, Ts... data)
0189       : dataFormats_(dataFormats), p_(p), frame_(ttStubRef, Frame()), data_(data...), trackId_(0) {}
0190 
0191   // construct StubPP from Frame
0192   StubPP::StubPP(const FrameStub& frame, const DataFormats* formats) : Stub(frame, formats, Process::pp) {
0193     for (int sectorEta = sectorEtaMin(); sectorEta <= sectorEtaMax(); sectorEta++)
0194       for (int sectorPhi = 0; sectorPhi < width(Variable::sectorsPhi); sectorPhi++)
0195         sectors_[sectorEta * width(Variable::sectorsPhi) + sectorPhi] = sectorsPhi()[sectorPhi];
0196   }
0197 
0198   // construct StubGP from Frame
0199   StubGP::StubGP(const FrameStub& frame, const DataFormats* formats, int sectorPhi, int sectorEta)
0200       : Stub(frame, formats, Process::gp), sectorPhi_(sectorPhi), sectorEta_(sectorEta) {
0201     const Setup* setup = dataFormats_->setup();
0202     inv2RBins_ = TTBV(0, setup->htNumBinsInv2R());
0203     for (int inv2R = inv2RMin(); inv2R <= inv2RMax(); inv2R++)
0204       inv2RBins_.set(inv2R + inv2RBins_.size() / 2);
0205   }
0206 
0207   // construct StubGP from StubPP
0208   StubGP::StubGP(const StubPP& stub, int sectorPhi, int sectorEta)
0209       : Stub(stub, stub.r(), stub.phi(), stub.z(), stub.layer(), stub.inv2RMin(), stub.inv2RMax()),
0210         sectorPhi_(sectorPhi),
0211         sectorEta_(sectorEta) {
0212     const Setup* setup = dataFormats_->setup();
0213     get<1>(data_) -= (sectorPhi_ - .5) * setup->baseSector();
0214     get<2>(data_) -= (r() + dataFormats_->chosenRofPhi()) * setup->sectorCot(sectorEta_);
0215     dataFormats_->convertStub(p_, data_, frame_.second);
0216   }
0217 
0218   // construct StubHT from Frame
0219   StubHT::StubHT(const FrameStub& frame, const DataFormats* formats, int inv2R)
0220       : Stub(frame, formats, Process::ht), inv2R_(inv2R) {
0221     fillTrackId();
0222   }
0223 
0224   // construct StubHT from StubGP and HT cell assignment
0225   StubHT::StubHT(const StubGP& stub, int phiT, int inv2R)
0226       : Stub(stub, stub.r(), stub.phi(), stub.z(), stub.layer(), stub.sectorPhi(), stub.sectorEta(), phiT),
0227         inv2R_(inv2R) {
0228     get<1>(data_) -=
0229         format(Variable::inv2R).floating(this->inv2R()) * r() + format(Variable::phiT).floating(this->phiT());
0230     fillTrackId();
0231     dataFormats_->convertStub(p_, data_, frame_.second);
0232   }
0233 
0234   void StubHT::fillTrackId() {
0235     TTBV ttBV(bv());
0236     trackId_ = ttBV.extract(width(Variable::sectorPhi) + width(Variable::sectorEta) + width(Variable::phiT));
0237   }
0238 
0239   // construct StubMHT from Frame
0240   StubMHT::StubMHT(const FrameStub& frame, const DataFormats* formats) : Stub(frame, formats, Process::mht) {
0241     fillTrackId();
0242   }
0243 
0244   // construct StubMHT from StubHT and MHT cell assignment
0245   StubMHT::StubMHT(const StubHT& stub, int phiT, int inv2R)
0246       : Stub(stub,
0247              stub.r(),
0248              stub.phi(),
0249              stub.z(),
0250              stub.layer(),
0251              stub.sectorPhi(),
0252              stub.sectorEta(),
0253              stub.phiT(),
0254              stub.inv2R()) {
0255     const Setup* setup = dataFormats_->setup();
0256     // update track (phIT, inv2R), and phi residuals w.r.t. track, to reflect MHT cell assignment.
0257     get<6>(data_) = this->phiT() * setup->mhtNumBinsPhiT() + phiT;
0258     get<7>(data_) = this->inv2R() * setup->mhtNumBinsInv2R() + inv2R;
0259     get<1>(data_) -= base(Variable::inv2R) * (inv2R - .5) * r() + base(Variable::phiT) * (phiT - .5);
0260     dataFormats_->convertStub(p_, data_, frame_.second);
0261     fillTrackId();
0262   }
0263 
0264   // fills track id
0265   void StubMHT::fillTrackId() {
0266     TTBV ttBV(bv());
0267     trackId_ = ttBV.extract(width(Variable::sectorPhi) + width(Variable::sectorEta) + width(Variable::phiT) +
0268                             width(Variable::inv2R));
0269   }
0270 
0271   // construct StubZHT from Frame
0272   StubZHT::StubZHT(const FrameStub& frame, const DataFormats* formats) : Stub(frame, formats, Process::zht) {
0273     cot_ = 0.;
0274     zT_ = 0.;
0275     fillTrackId();
0276   }
0277 
0278   // construct StubZHT from StubMHT
0279   StubZHT::StubZHT(const StubMHT& stub)
0280       : Stub(stub,
0281              stub.r(),
0282              stub.phi(),
0283              stub.z(),
0284              stub.layer(),
0285              stub.sectorPhi(),
0286              stub.sectorEta(),
0287              stub.phiT(),
0288              stub.inv2R(),
0289              0,
0290              0) {
0291     cot_ = 0.;
0292     zT_ = 0.;
0293     r_ = format(Variable::r).digi(this->r() + dataFormats_->chosenRofPhi() - dataFormats_->setup()->chosenRofZ());
0294     chi_ = stub.z();
0295     trackId_ = stub.trackId();
0296   }
0297 
0298   //
0299   StubZHT::StubZHT(const StubZHT& stub, double zT, double cot, int id)
0300       : Stub(stub.frame().first,
0301              stub.dataFormats(),
0302              Process::zht,
0303              stub.r(),
0304              stub.phi(),
0305              stub.z(),
0306              stub.layer(),
0307              stub.sectorPhi(),
0308              stub.sectorEta(),
0309              stub.phiT(),
0310              stub.inv2R(),
0311              stub.zT(),
0312              stub.cot()) {
0313     // update track (zT, cot), and phi residuals w.r.t. track, to reflect ZHT cell assignment.
0314     r_ = stub.r_;
0315     cot_ = stub.cotf() + cot;
0316     zT_ = stub.ztf() + zT;
0317     chi_ = stub.z() - zT_ + r_ * cot_;
0318     get<8>(data_) = format(Variable::zT).integer(zT_);
0319     get<9>(data_) = format(Variable::cot).integer(cot_);
0320     dataFormats_->convertStub(p_, data_, frame_.second);
0321     trackId_ = stub.trackId() * 4 + id;
0322   }
0323 
0324   //
0325   StubZHT::StubZHT(const StubZHT& stub, int cot, int zT)
0326       : Stub(stub.frame().first,
0327              stub.dataFormats(),
0328              Process::zht,
0329              stub.r(),
0330              stub.phi(),
0331              0.,
0332              stub.layer(),
0333              stub.sectorPhi(),
0334              stub.sectorEta(),
0335              stub.phiT(),
0336              stub.inv2R(),
0337              zT,
0338              cot) {
0339     get<2>(data_) =
0340         format(Variable::z)
0341             .digi(stub.z() - (format(Variable::zT).floating(zT) - stub.r_ * format(Variable::cot).floating(cot)));
0342     dataFormats_->convertStub(p_, data_, frame_.second);
0343     fillTrackId();
0344   }
0345 
0346   // fills track id
0347   void StubZHT::fillTrackId() {
0348     TTBV ttBV(bv());
0349     trackId_ = ttBV.extract(width(Variable::sectorPhi) + width(Variable::sectorEta) + width(Variable::phiT) +
0350                             width(Variable::inv2R) + width(Variable::zT) + width(Variable::cot));
0351   }
0352 
0353   // construct StubKFin from Frame
0354   StubKFin::StubKFin(const FrameStub& frame, const DataFormats* formats, int layer)
0355       : Stub(frame, formats, Process::kfin), layer_(layer) {}
0356 
0357   // construct StubKFin from StubZHT
0358   StubKFin::StubKFin(const StubZHT& stub, double dPhi, double dZ, int layer)
0359       : Stub(stub, stub.r(), stub.phi(), stub.z(), dPhi, dZ), layer_(layer) {
0360     dataFormats_->convertStub(p_, data_, frame_.second);
0361   }
0362 
0363   // construct StubKFin from TTStubRef
0364   StubKFin::StubKFin(const TTStubRef& ttStubRef,
0365                      const DataFormats* dataFormats,
0366                      double r,
0367                      double phi,
0368                      double z,
0369                      double dPhi,
0370                      double dZ,
0371                      int layer)
0372       : Stub(ttStubRef, dataFormats, Process::kfin, r, phi, z, dPhi, dZ), layer_(layer) {
0373     dataFormats_->convertStub(p_, data_, frame_.second);
0374   }
0375 
0376   // construct StubKF from Frame
0377   StubKF::StubKF(const FrameStub& frame, const DataFormats* formats, int layer)
0378       : Stub(frame, formats, Process::kf), layer_(layer) {}
0379 
0380   // construct StubKF from StubKFin
0381   StubKF::StubKF(const StubKFin& stub, double inv2R, double phiT, double cot, double zT)
0382       : Stub(stub, stub.r(), 0., 0., stub.dPhi(), stub.dZ()), layer_(stub.layer()) {
0383     const Setup* setup = dataFormats_->setup();
0384     get<1>(data_) = format(Variable::phi).digi(stub.phi() - (phiT + this->r() * inv2R));
0385     const double d =
0386         (dataFormats_->hybrid() ? setup->hybridChosenRofPhi() : setup->chosenRofPhi()) - setup->chosenRofZ();
0387     const double rz = format(Variable::r).digi(this->r() + d);
0388     get<2>(data_) = format(Variable::z).digi(stub.z() - (zT + rz * cot));
0389     dataFormats_->convertStub(p_, data_, frame_.second);
0390   }
0391 
0392   // construct Track from Frame
0393   template <typename... Ts>
0394   Track<Ts...>::Track(const FrameTrack& frame, const DataFormats* dataFormats, Process p)
0395       : dataFormats_(dataFormats), p_(p), frame_(frame) {
0396     dataFormats_->convertTrack(p_, frame.second, data_);
0397   }
0398 
0399   // construct Track from other Track
0400   template <typename... Ts>
0401   template <typename... Others>
0402   Track<Ts...>::Track(const Track<Others...>& track, Ts... data)
0403       : dataFormats_(track.dataFormats()), p_(++track.p()), frame_(track.frame().first, Frame()), data_(data...) {}
0404 
0405   // construct Track from Stub
0406   template <typename... Ts>
0407   template <typename... Others>
0408   Track<Ts...>::Track(const Stub<Others...>& stub, const TTTrackRef& ttTrackRef, Ts... data)
0409       : dataFormats_(stub.dataFormats()), p_(++stub.p()), frame_(ttTrackRef, Frame()), data_(data...) {}
0410 
0411   // construct Track from TTTrackRef
0412   template <typename... Ts>
0413   Track<Ts...>::Track(const TTTrackRef& ttTrackRef, const DataFormats* dataFormats, Process p, Ts... data)
0414       : dataFormats_(dataFormats), p_(p), frame_(ttTrackRef, Frame()), data_(data...) {}
0415 
0416   // construct TrackKFin from Frame
0417   TrackKFin::TrackKFin(const FrameTrack& frame, const DataFormats* dataFormats, const vector<StubKFin*>& stubs)
0418       : Track(frame, dataFormats, Process::kfin), stubs_(setup()->numLayers()), hitPattern_(0, setup()->numLayers()) {
0419     vector<int> nStubs(stubs_.size(), 0);
0420     for (StubKFin* stub : stubs)
0421       nStubs[stub->layer()]++;
0422     for (int layer = 0; layer < dataFormats->setup()->numLayers(); layer++)
0423       stubs_[layer].reserve(nStubs[layer]);
0424     for (StubKFin* stub : stubs) {
0425       const int layer = stub->layer();
0426       stubs_[layer].push_back(stub);
0427       hitPattern_.set(layer);
0428     }
0429   }
0430 
0431   // construct TrackKFin from StubZHT
0432   TrackKFin::TrackKFin(const StubZHT& stub, const TTTrackRef& ttTrackRef, const TTBV& maybePattern)
0433       : Track(stub, ttTrackRef, maybePattern, stub.sectorPhi(), stub.sectorEta(), 0., 0., 0., 0.),
0434         stubs_(setup()->numLayers()),
0435         hitPattern_(0, setup()->numLayers()) {
0436     get<3>(data_) = format(Variable::phiT, Process::mht).floating(stub.phiT());
0437     get<4>(data_) = format(Variable::inv2R, Process::mht).floating(stub.inv2R());
0438     get<5>(data_) = format(Variable::zT, Process::zht).floating(stub.zT());
0439     get<6>(data_) = format(Variable::cot, Process::zht).floating(stub.cot());
0440     dataFormats_->convertTrack(p_, data_, frame_.second);
0441   }
0442 
0443   // construct TrackKFin from TTTrackRef
0444   TrackKFin::TrackKFin(const TTTrackRef& ttTrackRef,
0445                        const DataFormats* dataFormats,
0446                        const TTBV& maybePattern,
0447                        double phiT,
0448                        double inv2R,
0449                        double zT,
0450                        double cot,
0451                        int sectorPhi,
0452                        int sectorEta)
0453       : Track(ttTrackRef, dataFormats, Process::kfin, maybePattern, sectorPhi, sectorEta, phiT, inv2R, zT, cot),
0454         stubs_(setup()->numLayers()),
0455         hitPattern_(0, setup()->numLayers()) {
0456     dataFormats_->convertTrack(p_, data_, frame_.second);
0457   }
0458 
0459   vector<TTStubRef> TrackKFin::ttStubRefs(const TTBV& hitPattern, const vector<int>& layerMap) const {
0460     vector<TTStubRef> stubs;
0461     stubs.reserve(hitPattern.count());
0462     for (int layer = 0; layer < setup()->numLayers(); layer++)
0463       if (hitPattern[layer])
0464         stubs.push_back(stubs_[layer][layerMap[layer]]->ttStubRef());
0465     return stubs;
0466   }
0467 
0468   // construct TrackKF from Frame
0469   TrackKF::TrackKF(const FrameTrack& frame, const DataFormats* dataFormats) : Track(frame, dataFormats, Process::kf) {}
0470 
0471   // construct TrackKF from TrackKfin
0472   TrackKF::TrackKF(const TrackKFin& track, double phiT, double inv2R, double zT, double cot)
0473       : Track(
0474             track, false, track.sectorPhi(), track.sectorEta(), track.phiT(), track.inv2R(), track.cot(), track.zT()) {
0475     get<0>(data_) = abs(inv2R) < dataFormats_->format(Variable::inv2R, Process::zht).base() / 2. &&
0476                     abs(phiT) < dataFormats_->format(Variable::phiT, Process::zht).base() / 2.;
0477     get<3>(data_) += phiT;
0478     get<4>(data_) += inv2R;
0479     get<5>(data_) += cot;
0480     get<6>(data_) += zT;
0481     dataFormats_->convertTrack(p_, data_, frame_.second);
0482   }
0483 
0484   // conversion to TTTrack with given stubs
0485   TTTrack<Ref_Phase2TrackerDigi_> TrackKF::ttTrack(const vector<StubKF>& stubs) const {
0486     const double invR = -this->inv2R() * 2.;
0487     const double phi0 =
0488         deltaPhi(this->phiT() - this->inv2R() * dataFormats_->chosenRofPhi() +
0489                  setup()->baseSector() * (this->sectorPhi() - .5) + setup()->baseRegion() * frame_.first->phiSector());
0490     const double cot = this->cot() + setup()->sectorCot(this->sectorEta());
0491     const double z0 = this->zT() - this->cot() * setup()->chosenRofZ();
0492     TTBV hitVector(0, setup()->numLayers());
0493     double chi2phi(0.);
0494     double chi2z(0.);
0495     vector<TTStubRef> ttStubRefs;
0496     const int nLayer = stubs.size();
0497     ttStubRefs.reserve(nLayer);
0498     for (const StubKF& stub : stubs) {
0499       hitVector.set(stub.layer());
0500       const TTStubRef& ttStubRef = stub.ttStubRef();
0501       chi2phi += pow(stub.phi(), 2) / setup()->v0(ttStubRef, this->inv2R());
0502       chi2z += pow(stub.z(), 2) / setup()->v1(ttStubRef, cot);
0503       ttStubRefs.push_back(ttStubRef);
0504     }
0505     static constexpr int nPar = 4;
0506     static constexpr double d0 = 0.;
0507     static constexpr double trkMVA1 = 0.;
0508     static constexpr double trkMVA2 = 0.;
0509     static constexpr double trkMVA3 = 0.;
0510     const int hitPattern = hitVector.val();
0511     const double bField = setup()->bField();
0512     TTTrack<Ref_Phase2TrackerDigi_> ttTrack(
0513         invR, phi0, cot, z0, d0, chi2phi, chi2z, trkMVA1, trkMVA2, trkMVA3, hitPattern, nPar, bField);
0514     ttTrack.setStubRefs(ttStubRefs);
0515     ttTrack.setPhiSector(frame_.first->phiSector());
0516     ttTrack.setEtaSector(this->sectorEta());
0517     ttTrack.setTrackSeedType(frame_.first->trackSeedType());
0518     ttTrack.setStubPtConsistency(StubPtConsistency::getConsistency(
0519         ttTrack, setup()->trackerGeometry(), setup()->trackerTopology(), bField, nPar));
0520     return ttTrack;
0521   }
0522 
0523   // construct TrackDR from Frame
0524   TrackDR::TrackDR(const FrameTrack& frame, const DataFormats* dataFormats) : Track(frame, dataFormats, Process::dr) {}
0525 
0526   // construct TrackDR from TrackKF
0527   TrackDR::TrackDR(const TrackKF& track) : Track(track, 0., 0., 0., 0.) {
0528     get<0>(data_) = track.phiT() + track.inv2R() * dataFormats_->chosenRofPhi() +
0529                     dataFormats_->format(Variable::phi, Process::gp).range() * (track.sectorPhi() - .5);
0530     get<1>(data_) = track.inv2R();
0531     get<2>(data_) = track.zT() - track.cot() * setup()->chosenRofZ();
0532     get<3>(data_) = track.cot() + setup()->sectorCot(track.sectorEta());
0533     dataFormats_->convertTrack(p_, data_, frame_.second);
0534   }
0535 
0536   // conversion to TTTrack
0537   TTTrack<Ref_Phase2TrackerDigi_> TrackDR::ttTrack() const {
0538     const double inv2R = this->inv2R();
0539     const double phi0 = this->phi0();
0540     const double cot = this->cot();
0541     const double z0 = this->z0();
0542     static constexpr double d0 = 0.;
0543     static constexpr double chi2phi = 0.;
0544     static constexpr double chi2z = 0;
0545     static constexpr double trkMVA1 = 0.;
0546     static constexpr double trkMVA2 = 0.;
0547     static constexpr double trkMVA3 = 0.;
0548     static constexpr int hitPattern = 0.;
0549     static constexpr int nPar = 4;
0550     static constexpr double bField = 0.;
0551     const int sectorPhi = frame_.first->phiSector();
0552     const int sectorEta = frame_.first->etaSector();
0553     TTTrack<Ref_Phase2TrackerDigi_> ttTrack(
0554         inv2R, phi0, cot, z0, d0, chi2phi, chi2z, trkMVA1, trkMVA2, trkMVA3, hitPattern, nPar, bField);
0555     ttTrack.setPhiSector(sectorPhi);
0556     ttTrack.setEtaSector(sectorEta);
0557     return ttTrack;
0558   }
0559 
0560   template <>
0561   Format<Variable::phiT, Process::ht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0562     range_ = 2. * M_PI / (double)(setup->numRegions() * setup->numSectorsPhi());
0563     base_ = range_ / (double)setup->htNumBinsPhiT();
0564     width_ = ceil(log2(setup->htNumBinsPhiT()));
0565   }
0566 
0567   template <>
0568   Format<Variable::phiT, Process::mht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0569     const Format<Variable::phiT, Process::ht> ht(iConfig, setup);
0570     range_ = ht.range();
0571     base_ = ht.base() / setup->mhtNumBinsPhiT();
0572     width_ = ceil(log2(setup->htNumBinsPhiT() * setup->mhtNumBinsPhiT()));
0573   }
0574 
0575   template <>
0576   Format<Variable::inv2R, Process::ht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0577     const double mintPt = iConfig.getParameter<bool>("UseHybrid") ? setup->hybridMinPtCand() : setup->minPt();
0578     range_ = 2. * setup->invPtToDphi() / mintPt;
0579     base_ = range_ / (double)setup->htNumBinsInv2R();
0580     width_ = ceil(log2(setup->htNumBinsInv2R()));
0581   }
0582 
0583   template <>
0584   Format<Variable::inv2R, Process::mht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0585     const Format<Variable::inv2R, Process::ht> ht(iConfig, setup);
0586     range_ = ht.range();
0587     base_ = ht.base() / setup->mhtNumBinsInv2R();
0588     width_ = ceil(log2(setup->htNumBinsInv2R() * setup->mhtNumBinsInv2R()));
0589   }
0590 
0591   template <>
0592   Format<Variable::r, Process::ht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0593     const double chosenRofPhi =
0594         iConfig.getParameter<bool>("UseHybrid") ? setup->hybridChosenRofPhi() : setup->chosenRofPhi();
0595     width_ = setup->tmttWidthR();
0596     range_ = 2. * max(abs(setup->outerRadius() - chosenRofPhi), abs(setup->innerRadius() - chosenRofPhi));
0597     const Format<Variable::phiT, Process::ht> phiT(iConfig, setup);
0598     const Format<Variable::inv2R, Process::ht> inv2R(iConfig, setup);
0599     base_ = phiT.base() / inv2R.base();
0600     const int shift = ceil(log2(range_ / base_ / pow(2., width_)));
0601     base_ *= pow(2., shift);
0602   }
0603 
0604   template <>
0605   Format<Variable::phi, Process::gp>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0606     const Format<Variable::phiT, Process::ht> phiT(iConfig, setup);
0607     const Format<Variable::inv2R, Process::ht> inv2R(iConfig, setup);
0608     const Format<Variable::r, Process::ht> r(iConfig, setup);
0609     range_ = phiT.range() + inv2R.range() * r.base() * pow(2., r.width()) / 4.;
0610     const Format<Variable::phi, Process::dtc> dtc(iConfig, setup);
0611     base_ = dtc.base();
0612     width_ = ceil(log2(range_ / base_));
0613   }
0614 
0615   template <>
0616   Format<Variable::phi, Process::dtc>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0617     width_ = setup->tmttWidthPhi();
0618     const Format<Variable::phiT, Process::ht> phiT(iConfig, setup);
0619     const Format<Variable::inv2R, Process::ht> inv2R(iConfig, setup);
0620     const Format<Variable::r, Process::ht> r(iConfig, setup);
0621     range_ = 2. * M_PI / (double)setup->numRegions() + inv2R.range() * r.base() * pow(2., r.width()) / 4.;
0622     const int shift = ceil(log2(range_ / phiT.base() / pow(2., width_)));
0623     base_ = phiT.base() * pow(2., shift);
0624   }
0625 
0626   template <>
0627   Format<Variable::phi, Process::ht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0628     const Format<Variable::phiT, Process::ht> phiT(iConfig, setup);
0629     range_ = 2. * phiT.base();
0630     const Format<Variable::phi, Process::gp> gp(iConfig, setup);
0631     base_ = gp.base();
0632     width_ = ceil(log2(range_ / base_));
0633   }
0634 
0635   template <>
0636   Format<Variable::phi, Process::mht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0637     const Format<Variable::phiT, Process::mht> phiT(iConfig, setup);
0638     range_ = 2. * phiT.base();
0639     const Format<Variable::phi, Process::ht> ht(iConfig, setup);
0640     base_ = ht.base();
0641     width_ = ceil(log2(range_ / base_));
0642   }
0643 
0644   template <>
0645   Format<Variable::phi, Process::kf>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0646     const Format<Variable::phi, Process::zht> phi(iConfig, setup);
0647     const double rangeFactor = iConfig.getParameter<ParameterSet>("KalmanFilter").getParameter<double>("RangeFactor");
0648     range_ = rangeFactor * phi.range();
0649     base_ = phi.base();
0650     width_ = ceil(log2(range_ / base_));
0651   }
0652 
0653   template <>
0654   Format<Variable::z, Process::dtc>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0655     width_ = setup->tmttWidthZ();
0656     range_ = 2. * setup->halfLength();
0657     const Format<Variable::r, Process::ht> r(iConfig, setup);
0658     const int shift = ceil(log2(range_ / r.base())) - width_;
0659     base_ = r.base() * pow(2., shift);
0660   }
0661 
0662   template <>
0663   Format<Variable::z, Process::gp>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0664     range_ = setup->neededRangeChiZ();
0665     const Format<Variable::z, Process::dtc> dtc(iConfig, setup);
0666     base_ = dtc.base();
0667     width_ = ceil(log2(range_ / base_));
0668   }
0669 
0670   template <>
0671   Format<Variable::zT, Process::zht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0672     const int numBinsZT = iConfig.getParameter<ParameterSet>("ZHoughTransform").getParameter<int>("NumBinsZT");
0673     const int numStages = iConfig.getParameter<ParameterSet>("ZHoughTransform").getParameter<int>("NumStages");
0674     width_ = ceil(log2(pow(numBinsZT, numStages)));
0675     const Format<Variable::z, Process::dtc> z(iConfig, setup);
0676     range_ = -1.;
0677     for (int eta = 0; eta < setup->numSectorsEta(); eta++)
0678       range_ = max(range_, (sinh(setup->boundarieEta(eta + 1)) - sinh(setup->boundarieEta(eta))));
0679     range_ *= setup->chosenRofZ();
0680     const int shift = ceil(log2(range_ / z.base() / pow(2., width_)));
0681     base_ = z.base() * pow(2., shift);
0682   }
0683 
0684   template <>
0685   Format<Variable::cot, Process::zht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0686     const int numBinsCot = iConfig.getParameter<ParameterSet>("ZHoughTransform").getParameter<int>("NumBinsCot");
0687     const int numStages = iConfig.getParameter<ParameterSet>("ZHoughTransform").getParameter<int>("NumStages");
0688     width_ = ceil(log2(pow(numBinsCot, numStages)));
0689     const Format<Variable::zT, Process::zht> zT(iConfig, setup);
0690     range_ = (zT.range() + 2. * setup->beamWindowZ()) / setup->chosenRofZ();
0691     const int shift = ceil(log2(range_)) - width_;
0692     base_ = pow(2., shift);
0693   }
0694 
0695   template <>
0696   Format<Variable::z, Process::zht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0697     const Format<Variable::zT, Process::zht> zT(iConfig, setup);
0698     const Format<Variable::cot, Process::zht> cot(iConfig, setup);
0699     const double rangeR =
0700         2. * max(abs(setup->outerRadius() - setup->chosenRofZ()), abs(setup->innerRadius() - setup->chosenRofZ()));
0701     range_ = zT.base() + cot.base() * rangeR + setup->maxdZ();
0702     const Format<Variable::z, Process::dtc> dtc(iConfig, setup);
0703     base_ = dtc.base();
0704     width_ = ceil(log2(range_ / base_));
0705     /*const Format<Variable::z, Process::gp> z(iConfig, setup);
0706     width_ = z.width();
0707     range_ = z.range();
0708     base_ = z.base();*/
0709   }
0710 
0711   template <>
0712   Format<Variable::z, Process::kfin>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0713     const Format<Variable::z, Process::zht> zht(iConfig, setup);
0714     range_ = zht.range() * pow(2, setup->kfinShiftRangeZ());
0715     base_ = zht.base();
0716     width_ = ceil(log2(range_ / base_));
0717   }
0718 
0719   template <>
0720   Format<Variable::phi, Process::zht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0721     const Format<Variable::phiT, Process::mht> phiT(iConfig, setup);
0722     const Format<Variable::inv2R, Process::mht> inv2R(iConfig, setup);
0723     const double chosenRofPhi =
0724         iConfig.getParameter<bool>("UseHybrid") ? setup->hybridChosenRofPhi() : setup->chosenRofPhi();
0725     const double rangeR = 2. * max(abs(setup->outerRadius() - chosenRofPhi), abs(setup->innerRadius() - chosenRofPhi));
0726     range_ = phiT.base() + inv2R.base() * rangeR + setup->maxdPhi();
0727     const Format<Variable::phi, Process::dtc> dtc(iConfig, setup);
0728     base_ = dtc.base();
0729     width_ = ceil(log2(range_ / base_));
0730   }
0731 
0732   template <>
0733   Format<Variable::phi, Process::kfin>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0734     const Format<Variable::phi, Process::zht> zht(iConfig, setup);
0735     range_ = zht.range() * pow(2, setup->kfinShiftRangePhi());
0736     base_ = zht.base();
0737     width_ = ceil(log2(range_ / base_));
0738   }
0739 
0740   template <>
0741   Format<Variable::z, Process::kf>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0742     /*const Format<Variable::z, Process::zht> z(iConfig, setup);
0743     const double rangeFactor = iConfig.getParameter<ParameterSet>("KalmanFilter").getParameter<double>("RangeFactor");
0744     range_ = rangeFactor * z.range();
0745     base_ = z.base();
0746     width_ = ceil(log2(range_ / base_));*/
0747     const Format<Variable::zT, Process::zht> zT(iConfig, setup);
0748     const Format<Variable::cot, Process::zht> cot(iConfig, setup);
0749     const double rangeR =
0750         2. * max(abs(setup->outerRadius() - setup->chosenRofZ()), abs(setup->innerRadius() - setup->chosenRofZ()));
0751     range_ = zT.base() + cot.base() * rangeR + setup->maxdZ();
0752     const Format<Variable::z, Process::dtc> dtc(iConfig, setup);
0753     base_ = dtc.base();
0754     const double rangeFactor = iConfig.getParameter<ParameterSet>("KalmanFilter").getParameter<double>("RangeFactor");
0755     range_ *= rangeFactor;
0756     width_ = ceil(log2(range_ / base_));
0757   }
0758 
0759   template <>
0760   Format<Variable::layer, Process::ht>::Format(const ParameterSet& iConfig, const Setup* setup) : DataFormat(false) {
0761     range_ = setup->numLayers();
0762     width_ = ceil(log2(range_));
0763   }
0764 
0765   template <>
0766   Format<Variable::sectorEta, Process::gp>::Format(const ParameterSet& iConfig, const Setup* setup)
0767       : DataFormat(false) {
0768     range_ = setup->numSectorsEta();
0769     width_ = ceil(log2(range_));
0770   }
0771 
0772   template <>
0773   Format<Variable::sectorPhi, Process::gp>::Format(const ParameterSet& iConfig, const Setup* setup)
0774       : DataFormat(false) {
0775     range_ = setup->numSectorsPhi();
0776     width_ = ceil(log2(range_));
0777   }
0778 
0779   template <>
0780   Format<Variable::sectorsPhi, Process::dtc>::Format(const ParameterSet& iConfig, const Setup* setup)
0781       : DataFormat(false) {
0782     range_ = setup->numSectorsPhi();
0783     width_ = setup->numSectorsPhi();
0784   }
0785 
0786   template <>
0787   Format<Variable::match, Process::kf>::Format(const edm::ParameterSet& iConfig, const Setup* setup)
0788       : DataFormat(false) {
0789     width_ = 1;
0790     range_ = 1.;
0791   }
0792 
0793   template <>
0794   Format<Variable::hitPattern, Process::kfin>::Format(const edm::ParameterSet& iConfig, const Setup* setup)
0795       : DataFormat(false) {
0796     width_ = setup->numLayers();
0797   }
0798 
0799   template <>
0800   Format<Variable::phi0, Process::dr>::Format(const edm::ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0801     const Format<Variable::inv2R, Process::ht> inv2R(iConfig, setup);
0802     const Format<Variable::phiT, Process::ht> phiT(iConfig, setup);
0803     const double chosenRofPhi =
0804         iConfig.getParameter<bool>("UseHybrid") ? setup->hybridChosenRofPhi() : setup->chosenRofPhi();
0805     width_ = setup->tfpWidthPhi0();
0806     range_ = 2. * M_PI / (double)setup->numRegions() + inv2R.range() * chosenRofPhi;
0807     base_ = phiT.base();
0808     const int shift = ceil(log2(range_ / base_ / pow(2., width_)));
0809     base_ *= pow(2., shift);
0810   }
0811 
0812   template <>
0813   Format<Variable::inv2R, Process::dr>::Format(const edm::ParameterSet& iConfig, const Setup* setup)
0814       : DataFormat(true) {
0815     const Format<Variable::inv2R, Process::ht> inv2R(iConfig, setup);
0816     width_ = setup->tfpWidthInv2R();
0817     range_ = inv2R.range();
0818     base_ = inv2R.base();
0819     const int shift = ceil(log2(range_ / base_ / pow(2., width_)));
0820     base_ *= pow(2., shift);
0821   }
0822 
0823   template <>
0824   Format<Variable::z0, Process::dr>::Format(const edm::ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0825     const Format<Variable::zT, Process::zht> zT(iConfig, setup);
0826     width_ = setup->tfpWidthZ0();
0827     range_ = 2. * setup->beamWindowZ();
0828     base_ = zT.base();
0829     const int shift = ceil(log2(range_ / base_ / pow(2., width_)));
0830     base_ *= pow(2., shift);
0831   }
0832 
0833   template <>
0834   Format<Variable::cot, Process::dr>::Format(const edm::ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0835     const Format<Variable::cot, Process::zht> cot(iConfig, setup);
0836     width_ = setup->tfpWidthCot();
0837     range_ = 2. * setup->maxCot();
0838     base_ = cot.base();
0839     const int shift = ceil(log2(range_ / base_ / pow(2., width_)));
0840     base_ *= pow(2., shift);
0841   }
0842 
0843   template <>
0844   Format<Variable::phiT, Process::kf>::Format(const edm::ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0845     const Format<Variable::phi0, Process::dr> phi0(iConfig, setup);
0846     const Format<Variable::phiT, Process::ht> phiT(iConfig, setup);
0847     const double rangeFactor = iConfig.getParameter<ParameterSet>("KalmanFilter").getParameter<double>("RangeFactor");
0848     range_ = rangeFactor * phiT.range();
0849     base_ = phi0.base();
0850     width_ = ceil(log2(range_ / base_));
0851   }
0852 
0853   template <>
0854   Format<Variable::inv2R, Process::kf>::Format(const edm::ParameterSet& iConfig, const Setup* setup)
0855       : DataFormat(true) {
0856     const Format<Variable::inv2R, Process::dr> dr(iConfig, setup);
0857     const Format<Variable::inv2R, Process::mht> mht(iConfig, setup);
0858     const double rangeFactor = iConfig.getParameter<ParameterSet>("KalmanFilter").getParameter<double>("RangeFactor");
0859     range_ = mht.range() + rangeFactor * mht.base();
0860     base_ = dr.base();
0861     width_ = ceil(log2(range_ / base_));
0862   }
0863 
0864   template <>
0865   Format<Variable::zT, Process::kf>::Format(const edm::ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0866     const Format<Variable::z0, Process::dr> z0(iConfig, setup);
0867     const Format<Variable::zT, Process::zht> zT(iConfig, setup);
0868     const double rangeFactor = iConfig.getParameter<ParameterSet>("KalmanFilter").getParameter<double>("RangeFactor");
0869     range_ = zT.range() + rangeFactor * zT.base();
0870     base_ = z0.base();
0871     width_ = ceil(log2(range_ / base_));
0872   }
0873 
0874   template <>
0875   Format<Variable::cot, Process::kf>::Format(const edm::ParameterSet& iConfig, const Setup* setup) : DataFormat(true) {
0876     const Format<Variable::cot, Process::dr> dr(iConfig, setup);
0877     const Format<Variable::cot, Process::zht> zht(iConfig, setup);
0878     const double rangeFactor = iConfig.getParameter<ParameterSet>("KalmanFilter").getParameter<double>("RangeFactor");
0879     range_ = zht.range() + rangeFactor * zht.base();
0880     base_ = dr.base();
0881     width_ = ceil(log2(range_ / base_));
0882   }
0883 
0884   template <>
0885   Format<Variable::dPhi, Process::kfin>::Format(const edm::ParameterSet& iConfig, const Setup* setup)
0886       : DataFormat(false) {
0887     const Format<Variable::phi, Process::kfin> phi(iConfig, setup);
0888     range_ = setup->maxdPhi();
0889     base_ = phi.base();
0890     width_ = ceil(log2(range_ / base_));
0891   }
0892 
0893   template <>
0894   Format<Variable::dZ, Process::kfin>::Format(const edm::ParameterSet& iConfig, const Setup* setup)
0895       : DataFormat(false) {
0896     const Format<Variable::z, Process::kfin> z(iConfig, setup);
0897     range_ = setup->maxdZ();
0898     base_ = z.base();
0899     width_ = ceil(log2(range_ / base_));
0900   }
0901 
0902 }  // namespace trackerTFP