Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-03 00:12:15

0001 #include "L1Trigger/TrackerTFP/interface/CleanTrackBuilder.h"
0002 
0003 #include <numeric>
0004 #include <algorithm>
0005 #include <iterator>
0006 #include <deque>
0007 #include <vector>
0008 #include <cmath>
0009 
0010 namespace trackerTFP {
0011 
0012   CleanTrackBuilder::CleanTrackBuilder(const tt::Setup* setup,
0013                                        const DataFormats* dataFormats,
0014                                        const LayerEncoding* layerEncoding,
0015                                        const DataFormat& cot,
0016                                        std::vector<StubCTB>& stubs,
0017                                        std::vector<TrackCTB>& tracks)
0018       : setup_(setup),
0019         dataFormats_(dataFormats),
0020         layerEncoding_(layerEncoding),
0021         cot_(cot),
0022         stubsCTB_(stubs),
0023         tracksCTB_(tracks),
0024         r_(dataFormats_->format(Variable::r, Process::ctb)),
0025         phi_(dataFormats_->format(Variable::phi, Process::ctb)),
0026         z_(dataFormats_->format(Variable::z, Process::ctb)),
0027         phiT_(dataFormats_->format(Variable::phiT, Process::ctb)),
0028         zT_(dataFormats_->format(Variable::zT, Process::ctb)) {
0029     stubs_.reserve(stubs.capacity());
0030     tracks_.reserve(tracks.capacity());
0031     numChannelOut_ = dataFormats_->numChannel(Process::ctb);
0032     numChannel_ = dataFormats_->numChannel(Process::ht) / numChannelOut_;
0033     numLayers_ = setup_->numLayers();
0034     wlayer_ = dataFormats_->width(Variable::layer, Process::ctb);
0035     numBinsInv2R_ = setup_->ctbNumBinsInv2R();
0036     numBinsPhiT_ = setup_->ctbNumBinsPhiT();
0037     numBinsCot_ = setup_->ctbNumBinsCot();
0038     numBinsZT_ = setup_->ctbNumBinsZT();
0039     baseInv2R_ = dataFormats_->base(Variable::inv2R, Process::ctb) / numBinsInv2R_;
0040     basePhiT_ = phiT_.base() / numBinsPhiT_;
0041     baseCot_ = cot_.base();
0042     baseZT_ = zT_.base() / numBinsZT_;
0043   }
0044 
0045   // fill output products
0046   void CleanTrackBuilder::produce(const std::vector<std::vector<StubHT*>>& streamsIn,
0047                                   std::vector<std::deque<TrackCTB*>>& regionTracks,
0048                                   std::vector<std::vector<std::deque<StubCTB*>>>& regionStubs) {
0049     // loop over worker
0050     for (int channelOut = 0; channelOut < numChannelOut_; channelOut++) {
0051       // clean input tracks
0052       std::vector<std::deque<Track*>> streamsT(numChannel_);
0053       std::vector<std::deque<Stub*>> streamsS(numChannel_);
0054       for (int cin = 0; cin < numChannel_; cin++) {
0055         const int index = numChannel_ * cin + channelOut;
0056         cleanStream(streamsIn[index], streamsT[cin], streamsS[cin], index);
0057       }
0058       // route
0059       std::deque<Track*> tracks;
0060       std::vector<std::deque<Stub*>> stubs(numLayers_);
0061       route(streamsT, tracks);
0062       route(streamsS, stubs);
0063       // sort
0064       sort(tracks, stubs);
0065       // convert
0066       std::deque<TrackCTB*>& channelTracks = regionTracks[channelOut];
0067       std::vector<std::deque<StubCTB*>>& channelStubs = regionStubs[channelOut];
0068       convert(tracks, stubs, channelTracks, channelStubs);
0069     }
0070   }
0071 
0072   //
0073   void CleanTrackBuilder::cleanStream(const std::vector<StubHT*>& input,
0074                                       std::deque<Track*>& tracks,
0075                                       std::deque<Stub*>& stubs,
0076                                       int channelId) {
0077     const DataFormat& dfInv2R = dataFormats_->format(Variable::inv2R, Process::ht);
0078     const double inv2R = dfInv2R.floating(dfInv2R.toSigned(channelId));
0079     const int offset = channelId * setup_->ctbMaxTracks();
0080     int trackId = offset;
0081     // identify tracks in input container
0082     int id;
0083     auto toTrkId = [this](StubHT* stub) {
0084       const DataFormat& phiT = dataFormats_->format(Variable::phiT, Process::ht);
0085       const DataFormat& zT = dataFormats_->format(Variable::zT, Process::ht);
0086       return (phiT.ttBV(stub->phiT()) + zT.ttBV(stub->zT())).val();
0087     };
0088     auto different = [&id, toTrkId](StubHT* stub) { return id != toTrkId(stub); };
0089     int delta = -setup_->htMinLayers() + 1;
0090     int old = 0;
0091     for (auto it = input.begin(); it != input.end();) {
0092       id = toTrkId(*it);
0093       const auto start = it;
0094       const auto end = std::find_if(start, input.end(), different);
0095       const std::vector<StubHT*> track(start, end);
0096       // restore clock accurancy
0097       delta += track.size() - old;
0098       old = track.size();
0099       if (delta > 0) {
0100         stubs.insert(stubs.end(), delta, nullptr);
0101         tracks.insert(tracks.end(), delta, nullptr);
0102         delta = 0;
0103       }
0104       // run single track through r-phi and r-z hough transform
0105       cleanTrack(track, tracks, stubs, inv2R, (*start)->zT(), trackId++);
0106       if (trackId - offset == setup_->ctbMaxTracks())
0107         break;
0108       // set begin of next track
0109       it = end;
0110     }
0111   }
0112 
0113   // run single track through r-phi and r-z hough transform
0114   void CleanTrackBuilder::cleanTrack(const std::vector<StubHT*>& track,
0115                                      std::deque<Track*>& tracks,
0116                                      std::deque<Stub*>& stubs,
0117                                      double inv2R,
0118                                      int zT,
0119                                      int trackId) {
0120     const TTBV& maybePattern = layerEncoding_->maybePattern(zT);
0121     auto noTrack = [this, &maybePattern](const TTBV& pattern) {
0122       // not enough seeding layer
0123       if (pattern.count(0, setup_->kfMaxSeedingLayer()) < 2)
0124         return true;
0125       int nHits(0);
0126       int nGaps(0);
0127       bool doubleGap = false;
0128       for (int layer = 0; layer < numLayers_; layer++) {
0129         if (pattern.test(layer)) {
0130           doubleGap = false;
0131           if (++nHits == setup_->ctbMinLayers())
0132             return false;
0133         } else if (!maybePattern.test(layer)) {
0134           if (++nGaps == setup_->kfMaxGaps() || doubleGap)
0135             break;
0136           doubleGap = true;
0137         }
0138       }
0139       return true;
0140     };
0141     auto toLayerId = [this](StubHT* stub) { return stub->layer().val(wlayer_); };
0142     auto toDPhi = [this, inv2R](StubHT* stub) {
0143       const bool barrel = stub->layer()[5];
0144       const bool ps = stub->layer()[4];
0145       const bool tilt = stub->layer()[3];
0146       const double pitchRow = ps ? setup_->pitchRowPS() : setup_->pitchRow2S();
0147       const double pitchCol = ps ? setup_->pitchColPS() : setup_->pitchCol2S();
0148       const double pitchColR = barrel ? (tilt ? setup_->tiltUncertaintyR() : 0.0) : pitchCol;
0149       const double r = stub->r() + setup_->chosenRofPhi();
0150       const double dPhi = pitchRow / r + (setup_->scattering() + pitchColR) * std::abs(inv2R);
0151       return phi_.digi(dPhi / 2.);
0152     };
0153     auto toDZ = [this](StubHT* stub) {
0154       const double m = setup_->tiltApproxSlope();
0155       const double c = setup_->tiltApproxIntercept();
0156       const bool barrel = stub->layer()[5];
0157       const bool ps = stub->layer()[4];
0158       const bool tilt = stub->layer()[3];
0159       const double pitchCol = ps ? setup_->pitchColPS() : setup_->pitchCol2S();
0160       const double zT = zT_.floating(stub->zT());
0161       const double cot = std::abs(zT) / setup_->chosenRofZ();
0162       const double dZ = (barrel ? (tilt ? m * cot + c : 1.) : cot) * pitchCol;
0163       return z_.digi(dZ / 2.);
0164     };
0165     std::vector<Stub*> tStubs;
0166     tStubs.reserve(track.size());
0167     std::vector<TTBV> hitPatternPhi(numBinsInv2R_ * numBinsPhiT_, TTBV(0, numLayers_));
0168     std::vector<TTBV> hitPatternZ(numBinsCot_ * numBinsZT_, TTBV(0, numLayers_));
0169     TTBV tracksPhi(0, numBinsInv2R_ * numBinsPhiT_);
0170     TTBV tracksZ(0, numBinsCot_ * numBinsZT_);
0171     // identify finer tracks each stub is consistent with
0172     for (StubHT* stub : track) {
0173       const int layerId = toLayerId(stub);
0174       const double dPhi = toDPhi(stub);
0175       const double dZ = toDZ(stub);
0176       // r - phi HT
0177       auto phiT = [stub](double inv2R, double dPhi) { return inv2R * stub->r() + stub->phi() + dPhi; };
0178       TTBV hitsPhi(0, numBinsInv2R_ * numBinsPhiT_);
0179       for (int binInv2R = 0; binInv2R < numBinsInv2R_; binInv2R++) {
0180         const int offset = binInv2R * numBinsPhiT_;
0181         const double inv2RMin = (binInv2R - numBinsInv2R_ / 2.) * baseInv2R_;
0182         const double inv2RMax = inv2RMin + baseInv2R_;
0183         const auto phiTs = {phiT(inv2RMin, -dPhi), phiT(inv2RMax, -dPhi), phiT(inv2RMin, dPhi), phiT(inv2RMax, dPhi)};
0184         const int binPhiTMin =
0185             std::floor(*std::min_element(phiTs.begin(), phiTs.end()) / basePhiT_ + 1.e-11) + numBinsPhiT_ / 2;
0186         const int binPhiTMax =
0187             std::floor(*std::max_element(phiTs.begin(), phiTs.end()) / basePhiT_ + 1.e-11) + numBinsPhiT_ / 2;
0188         for (int binPhiT = 0; binPhiT < numBinsPhiT_; binPhiT++)
0189           if (binPhiT >= binPhiTMin && binPhiT <= binPhiTMax)
0190             hitsPhi.set(offset + binPhiT);
0191       }
0192       // check for tracks on the fly
0193       for (int phi : hitsPhi.ids()) {
0194         hitPatternPhi[phi].set(layerId);
0195         if (!noTrack(hitPatternPhi[phi]))
0196           tracksPhi.set(phi);
0197       }
0198       // r - z HT
0199       auto zT = [this, stub](double cot, double dZ) {
0200         const double r = r_.digi(stub->r() + r_.digi(setup_->chosenRofPhi() - setup_->chosenRofZ()));
0201         return cot * r + stub->z() + dZ;
0202       };
0203       TTBV hitsZ(0, numBinsCot_ * numBinsZT_);
0204       for (int binCot = 0; binCot < numBinsCot_; binCot++) {
0205         const int offset = binCot * numBinsZT_;
0206         const double cotMin = (binCot - numBinsCot_ / 2.) * baseCot_;
0207         const double cotMax = cotMin + baseCot_;
0208         const auto zTs = {zT(cotMin, -dZ), zT(cotMax, -dZ), zT(cotMin, dZ), zT(cotMax, dZ)};
0209         const int binZTMin = std::floor(*std::min_element(zTs.begin(), zTs.end()) / baseZT_ + 1.e-11) + numBinsZT_ / 2;
0210         const int binZTMax = std::floor(*std::max_element(zTs.begin(), zTs.end()) / baseZT_ + 1.e-11) + numBinsZT_ / 2;
0211         for (int binZT = 0; binZT < numBinsZT_; binZT++)
0212           if (binZT >= binZTMin && binZT <= binZTMax)
0213             hitsZ.set(offset + binZT);
0214       }
0215       // check for tracks on the fly
0216       for (int z : hitsZ.ids()) {
0217         hitPatternZ[z].set(layerId);
0218         if (!noTrack(hitPatternZ[z]))
0219           tracksZ.set(z);
0220       }
0221       // store stubs consistent finer tracks
0222       stubs_.emplace_back(stub, trackId, hitsPhi, hitsZ, layerId, dPhi, dZ);
0223       tStubs.push_back(&stubs_.back());
0224     }
0225     // clean
0226     tracks.insert(tracks.end(), tStubs.size() - 1, nullptr);
0227     tracks_.emplace_back(setup_, trackId, tracksPhi, tracksZ, tStubs, inv2R);
0228     tracks.push_back(&tracks_.back());
0229     stubs.insert(stubs.end(), tStubs.begin(), tStubs.end());
0230   }
0231 
0232   //
0233   void CleanTrackBuilder::Stub::update(const TTBV& phi, const TTBV& z, std::vector<int>& ids, int max) {
0234     auto consistent = [](const TTBV& stub, const TTBV& track) {
0235       for (int id : track.ids())
0236         if (stub[id])
0237           return true;
0238       return false;
0239     };
0240     if (consistent(hitsPhi_, phi) && consistent(hitsZ_, z) && ids[layerId_] < max)
0241       stubId_ = ids[layerId_]++;
0242     else
0243       valid_ = false;
0244   }
0245 
0246   // construct Track from Stubs
0247   CleanTrackBuilder::Track::Track(const tt::Setup* setup,
0248                                   int trackId,
0249                                   const TTBV& hitsPhi,
0250                                   const TTBV& hitsZ,
0251                                   const std::vector<Stub*>& stubs,
0252                                   double inv2R)
0253       : valid_(true), stubs_(stubs), trackId_(trackId), hitsPhi_(hitsPhi), hitsZ_(hitsZ), inv2R_(inv2R) {
0254     std::vector<int> stubIds(setup->numLayers(), 0);
0255     for (Stub* stub : stubs_)
0256       stub->update(hitsPhi_, hitsZ_, stubIds, setup->ctbMaxStubs());
0257     const int nLayer =
0258         std::accumulate(stubIds.begin(), stubIds.end(), 0, [](int sum, int i) { return sum + (i > 0 ? 1 : 0); });
0259     if (nLayer < setup->ctbMinLayers())
0260       valid_ = false;
0261     size_ = *std::max_element(stubIds.begin(), stubIds.end());
0262   }
0263 
0264   //
0265   void CleanTrackBuilder::route(std::vector<std::deque<Stub*>>& input, std::vector<std::deque<Stub*>>& outputs) const {
0266     for (int channelOut = 0; channelOut < (int)outputs.size(); channelOut++) {
0267       std::deque<Stub*>& output = outputs[channelOut];
0268       std::vector<std::deque<Stub*>> inputs(input);
0269       for (std::deque<Stub*>& stream : inputs) {
0270         for (Stub*& stub : stream)
0271           if (stub && (!stub->valid_ || stub->layerId_ != channelOut))
0272             stub = nullptr;
0273         for (auto it = stream.end(); it != stream.begin();)
0274           it = (*--it) ? stream.begin() : stream.erase(it);
0275       }
0276       std::vector<std::deque<Stub*>> stacks(input.size());
0277       // clock accurate firmware emulation, each while trip describes one clock tick, one stub in and one stub out per tick
0278       while (!std::all_of(inputs.begin(), inputs.end(), [](const std::deque<Stub*>& stubs) { return stubs.empty(); }) ||
0279              !std::all_of(stacks.begin(), stacks.end(), [](const std::deque<Stub*>& stubs) { return stubs.empty(); })) {
0280         // fill input fifos
0281         for (int channel = 0; channel < static_cast<int>(input.size()); channel++) {
0282           std::deque<Stub*>& stack = stacks[channel];
0283           Stub* stub = pop_front(inputs[channel]);
0284           if (stub) {
0285             if (setup_->enableTruncation() && (int)stack.size() == setup_->ctbDepthMemory() - 1)
0286               pop_front(stack);
0287             stack.push_back(stub);
0288           }
0289         }
0290         // merge input fifos to one stream, prioritizing lower input channel over higher channel
0291         bool nothingToRoute(true);
0292         for (int channel = 0; channel < static_cast<int>(input.size()); channel++) {
0293           Stub* stub = pop_front(stacks[channel]);
0294           if (stub) {
0295             nothingToRoute = false;
0296             output.push_back(stub);
0297             break;
0298           }
0299         }
0300         if (nothingToRoute)
0301           output.push_back(nullptr);
0302       }
0303     }
0304   }
0305 
0306   //
0307   void CleanTrackBuilder::route(std::vector<std::deque<Track*>>& inputs, std::deque<Track*>& output) const {
0308     std::vector<std::deque<Track*>> stacks(inputs.size());
0309     // clock accurate firmware emulation, each while trip describes one clock tick, one stub in and one stub out per tick
0310     while (!std::all_of(inputs.begin(), inputs.end(), [](const std::deque<Track*>& tracks) {
0311       return tracks.empty();
0312     }) || !std::all_of(stacks.begin(), stacks.end(), [](const std::deque<Track*>& tracks) { return tracks.empty(); })) {
0313       // fill input fifos
0314       for (int channel = 0; channel < static_cast<int>(inputs.size()); channel++) {
0315         std::deque<Track*>& stack = stacks[channel];
0316         Track* track = pop_front(inputs[channel]);
0317         if (track && track->valid_) {
0318           if (setup_->enableTruncation() && static_cast<int>(stack.size()) == setup_->ctbDepthMemory() - 1)
0319             pop_front(stack);
0320           stack.push_back(track);
0321         }
0322       }
0323       // merge input fifos to one stream, prioritizing lower input channel over higher channel
0324       bool nothingToRoute(true);
0325       for (int channel = 0; channel < static_cast<int>(inputs.size()); channel++) {
0326         Track* track = pop_front(stacks[channel]);
0327         if (track) {
0328           nothingToRoute = false;
0329           output.push_back(track);
0330           break;
0331         }
0332       }
0333       if (nothingToRoute)
0334         output.push_back(nullptr);
0335     }
0336   }
0337 
0338   // sort
0339   void CleanTrackBuilder::sort(std::deque<Track*>& tracks, std::vector<std::deque<Stub*>>& stubs) const {
0340     // aplly truncation
0341     if (setup_->enableTruncation()) {
0342       if (static_cast<int>(tracks.size()) > setup_->numFramesHigh())
0343         tracks.resize(setup_->numFramesHigh());
0344       for (std::deque<Stub*>& stream : stubs)
0345         if (static_cast<int>(stream.size()) > setup_->numFramesHigh())
0346           stream.resize(setup_->numFramesHigh());
0347     }
0348     // cycle event, remove all gaps
0349     tracks.erase(std::remove(tracks.begin(), tracks.end(), nullptr), tracks.end());
0350     for (std::deque<Stub*>& stream : stubs)
0351       stream.erase(std::remove(stream.begin(), stream.end(), nullptr), stream.end());
0352     // prepare sort according to track id arrival order
0353     std::vector<int> trackIds;
0354     trackIds.reserve(tracks.size());
0355     std::transform(
0356         tracks.begin(), tracks.end(), std::back_inserter(trackIds), [](Track* track) { return track->trackId_; });
0357     auto cleaned = [&trackIds](Stub* stub) {
0358       return std::find(trackIds.begin(), trackIds.end(), stub->trackId_) == trackIds.end();
0359     };
0360     auto order = [&trackIds](auto lhs, auto rhs) {
0361       const auto l = std::find(trackIds.begin(), trackIds.end(), lhs->trackId_);
0362       const auto r = std::find(trackIds.begin(), trackIds.end(), rhs->trackId_);
0363       return std::distance(r, l) < 0;
0364     };
0365     for (std::deque<Stub*>& stream : stubs) {
0366       // remove stubs from removed tracks
0367       stream.erase(std::remove_if(stream.begin(), stream.end(), cleaned), stream.end());
0368       // sort according to stub id on layer
0369       std::stable_sort(stream.begin(), stream.end(), [](Stub* lhs, Stub* rhs) { return lhs->stubId_ < rhs->stubId_; });
0370       // sort according to track id arrival order
0371       std::stable_sort(stream.begin(), stream.end(), order);
0372     }
0373     // add all gaps
0374     const int size =
0375         std::accumulate(tracks.begin(), tracks.end(), 0, [](int sum, Track* track) { return sum + track->size_; });
0376     for (int frame = 0; frame < size;) {
0377       const int trackId = tracks[frame]->trackId_;
0378       const int length = tracks[frame]->size_;
0379       tracks.insert(std::next(tracks.begin(), frame + 1), length - 1, nullptr);
0380       for (int layer = 0; layer < numLayers_; layer++) {
0381         std::deque<Stub*>& stream = stubs[layer];
0382         if (frame >= static_cast<int>(stream.size())) {
0383           stream.insert(stream.end(), length, nullptr);
0384           continue;
0385         }
0386         const auto begin = std::next(stream.begin(), frame);
0387         const auto end = std::find_if(begin, stream.end(), [trackId](Stub* stub) { return stub->trackId_ != trackId; });
0388         stream.insert(end, length - std::distance(begin, end), nullptr);
0389       }
0390       frame += length;
0391     }
0392   }
0393 
0394   //
0395   void CleanTrackBuilder::convert(const std::deque<Track*>& iTracks,
0396                                   const std::vector<std::deque<Stub*>>& iStubs,
0397                                   std::deque<TrackCTB*>& oTracks,
0398                                   std::vector<std::deque<StubCTB*>>& oStubs) {
0399     for (int iFrame = 0; iFrame < static_cast<int>(iTracks.size());) {
0400       Track* track = iTracks[iFrame];
0401       if (!track) {
0402         oTracks.push_back(nullptr);
0403         for (std::deque<StubCTB*>& stubs : oStubs)
0404           stubs.push_back(nullptr);
0405         iFrame++;
0406         continue;
0407       }
0408       StubHT* s = nullptr;
0409       for (int layer = 0; layer < numLayers_; layer++) {
0410         for (int n = 0; n < track->size_; n++) {
0411           Stub* stub = iStubs[layer][iFrame + n];
0412           if (!stub) {
0413             oStubs[layer].push_back(nullptr);
0414             continue;
0415           }
0416           s = stub->stubHT_;
0417           const double r = s->r();
0418           const double phi = s->phi();
0419           const double z = s->z();
0420           const double dPhi = stub->dPhi_;
0421           const double dZ = stub->dZ_;
0422           stubsCTB_.emplace_back(*s, r, phi, z, dPhi, dZ);
0423           oStubs[layer].push_back(&stubsCTB_.back());
0424         }
0425       }
0426       const double inv2R = track->inv2R_;
0427       const double phiT = dataFormats_->format(Variable::phiT, Process::ctb).floating(s->phiT());
0428       const double zT = dataFormats_->format(Variable::zT, Process::ctb).floating(s->zT());
0429       tracksCTB_.emplace_back(TTTrackRef(), dataFormats_, inv2R, phiT, zT);
0430       oTracks.push_back(&tracksCTB_.back());
0431       oTracks.insert(oTracks.end(), track->size_ - 1, nullptr);
0432       iFrame += track->size_;
0433     }
0434   }
0435 
0436   // remove and return first element of deque, returns nullptr if empty
0437   template <class T>
0438   T* CleanTrackBuilder::pop_front(std::deque<T*>& ts) const {
0439     T* t = nullptr;
0440     if (!ts.empty()) {
0441       t = ts.front();
0442       ts.pop_front();
0443     }
0444     return t;
0445   }
0446 
0447   void CleanTrackBuilder::put(TrackCTB* track,
0448                               const std::vector<std::vector<StubCTB*>>& stubs,
0449                               int region,
0450                               tt::TTTracks& ttTracks) const {
0451     const double dPhi = dataFormats_->format(Variable::phiT, Process::ctb).range();
0452     const double invR = -track->inv2R() * 2.;
0453     const double phi0 = tt::deltaPhi(track->phiT() - track->inv2R() * setup_->chosenRofPhi() + region * dPhi);
0454     const double zT = track->zT();
0455     const double cot = zT / setup_->chosenRofZ();
0456     TTBV hits(0, numLayers_);
0457     double chi2phi(0.);
0458     double chi2z(0.);
0459     const int nStubs = accumulate(
0460         stubs.begin(), stubs.end(), 0, [](int sum, const std::vector<StubCTB*>& layer) { return sum += layer.size(); });
0461     std::vector<TTStubRef> ttStubRefs;
0462     ttStubRefs.reserve(nStubs);
0463     for (int layer = 0; layer < numLayers_; layer++) {
0464       for (StubCTB* stub : stubs[layer]) {
0465         hits.set(layer);
0466         chi2phi += std::pow(stub->phi(), 2) / pow(stub->dPhi(), 2);
0467         chi2z += std::pow(stub->z(), 2) / pow(stub->dZ(), 2);
0468         ttStubRefs.push_back(stub->frame().first);
0469       }
0470     }
0471     static constexpr int nPar = 4;
0472     static constexpr double d0 = 0.;
0473     static constexpr double z0 = 0;
0474     static constexpr double trkMVA1 = 0.;
0475     static constexpr double trkMVA2 = 0.;
0476     static constexpr double trkMVA3 = 0.;
0477     const int hitPattern = hits.val();
0478     const double bField = setup_->bField();
0479     TTTrack<Ref_Phase2TrackerDigi_> ttTrack(
0480         invR, phi0, cot, z0, d0, chi2phi, chi2z, trkMVA1, trkMVA2, trkMVA3, hitPattern, nPar, bField);
0481     ttTrack.setStubRefs(ttStubRefs);
0482     ttTrack.setPhiSector(region);
0483     ttTracks.emplace_back(ttTrack);
0484   }
0485 
0486 }  // namespace trackerTFP