Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-07-07 22:33:33

0001 #include "L1Trigger/TrackFindingTracklet/interface/MatchProcessor.h"
0002 #include "L1Trigger/TrackFindingTracklet/interface/Globals.h"
0003 #include "L1Trigger/TrackFindingTracklet/interface/Util.h"
0004 #include "L1Trigger/TrackFindingTracklet/interface/HistBase.h"
0005 
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "FWCore/Utilities/interface/Exception.h"
0008 #include "DataFormats/Math/interface/deltaPhi.h"
0009 #include "L1Trigger/TrackFindingTracklet/interface/IMATH_TrackletCalculator.h"
0010 
0011 #include <filesystem>
0012 
0013 using namespace std;
0014 using namespace trklet;
0015 
0016 MatchProcessor::MatchProcessor(string name, Settings const& settings, Globals* global)
0017     : ProcessBase(name, settings, global),
0018       phimatchcuttable_(settings),
0019       zmatchcuttable_(settings),
0020       rphicutPStable_(settings),
0021       rphicut2Stable_(settings),
0022       rcutPStable_(settings),
0023       rcut2Stable_(settings),
0024       fullmatches_(12),
0025       rinvbendlut_(settings),
0026       luttable_(settings),
0027       inputProjBuffer_(3) {
0028   phiregion_ = name[8] - 'A';
0029 
0030   layerdisk_ = initLayerDisk(3);
0031 
0032   barrel_ = layerdisk_ < N_LAYER;
0033 
0034   phishift_ = settings_.nphibitsstub(N_LAYER - 1) - settings_.nphibitsstub(layerdisk_);
0035   dzshift_ = settings_.nzbitsstub(0) - settings_.nzbitsstub(layerdisk_);
0036 
0037   if (barrel_) {
0038     icorrshift_ = ilog2(settings_.kphi(layerdisk_) / (settings_.krbarrel() * settings_.kphider()));
0039     icorzshift_ = ilog2(settings_.kz(layerdisk_) / (settings_.krbarrel() * settings_.kzder()));
0040   } else {
0041     icorrshift_ = ilog2(settings_.kphi(layerdisk_) / (settings_.kz() * settings_.kphiderdisk()));
0042     icorzshift_ = ilog2(settings_.krprojshiftdisk() / (settings_.kz() * settings_.krder()));
0043   }
0044 
0045   luttable_.initBendMatch(layerdisk_);
0046 
0047   nrbits_ = 5;
0048   nphiderbits_ = 6;
0049 
0050   if (!barrel_) {
0051     rinvbendlut_.initProjectionBend(
0052         global->ITC_L1L2()->der_phiD_final.K(), layerdisk_ - N_LAYER, nrbits_, nphiderbits_);
0053   }
0054 
0055   nrinv_ = NRINVBITS;
0056 
0057   unsigned int region = getName()[8] - 'A';
0058   assert(region < settings_.nallstubs(layerdisk_));
0059 
0060   if (barrel_) {
0061     phimatchcuttable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::barrelphi, region);
0062     zmatchcuttable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::barrelz, region);
0063   } else {
0064     rphicutPStable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::diskPSphi, region);
0065     rphicut2Stable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::disk2Sphi, region);
0066     rcutPStable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::diskPSr, region);
0067     rcut2Stable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::disk2Sr, region);
0068   }
0069 
0070   for (unsigned int i = 0; i < N_DSS_MOD * 2; i++) {
0071     ialphafactinner_[i] = (1 << settings_.alphashift()) * settings_.krprojshiftdisk() * settings_.half2SmoduleWidth() /
0072                           (1 << (settings_.nbitsalpha() - 1)) / (settings_.rDSSinner(i) * settings_.rDSSinner(i)) /
0073                           settings_.kphi();
0074     ialphafactouter_[i] = (1 << settings_.alphashift()) * settings_.krprojshiftdisk() * settings_.half2SmoduleWidth() /
0075                           (1 << (settings_.nbitsalpha() - 1)) / (settings_.rDSSouter(i) * settings_.rDSSouter(i)) /
0076                           settings_.kphi();
0077   }
0078 
0079   nvm_ = settings_.nvmme(layerdisk_) * settings_.nallstubs(layerdisk_);
0080   nvmbins_ = settings_.nvmme(layerdisk_);
0081   nvmbits_ = settings_.nbitsvmme(layerdisk_) + settings_.nbitsallstubs(layerdisk_);
0082 
0083   nMatchEngines_ = 4;
0084   for (unsigned int iME = 0; iME < nMatchEngines_; iME++) {
0085     MatchEngineUnit tmpME(barrel_, layerdisk_, luttable_);
0086     matchengines_.push_back(tmpME);
0087   }
0088 }
0089 
0090 void MatchProcessor::addOutput(MemoryBase* memory, string output) {
0091   if (settings_.writetrace()) {
0092     edm::LogVerbatim("Tracklet") << "In " << name_ << " adding output to " << memory->getName() << " to output "
0093                                  << output;
0094   }
0095   if (output.find("matchout") != std::string::npos) {
0096     auto* tmp = dynamic_cast<FullMatchMemory*>(memory);
0097     assert(tmp != nullptr);
0098     unsigned int iSeed = getISeed(tmp->getName());
0099     assert(iSeed < fullmatches_.size());
0100     assert(fullmatches_[iSeed] == nullptr);
0101     fullmatches_[iSeed] = tmp;
0102     return;
0103   }
0104   throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " could not find output: " << output;
0105 }
0106 
0107 void MatchProcessor::addInput(MemoryBase* memory, string input) {
0108   if (settings_.writetrace()) {
0109     edm::LogVerbatim("Tracklet") << "In " << name_ << " adding input from " << memory->getName() << " to input "
0110                                  << input;
0111   }
0112   if (input == "allstubin") {
0113     auto* tmp = dynamic_cast<AllStubsMemory*>(memory);
0114     assert(tmp != nullptr);
0115     allstubs_ = tmp;
0116     return;
0117   }
0118   if (input == "vmstubin") {
0119     auto* tmp = dynamic_cast<VMStubsMEMemory*>(memory);
0120     assert(tmp != nullptr);
0121     vmstubs_.push_back(tmp);  //to allow more than one stub in?  vmstubs_=tmp;
0122     return;
0123   }
0124   if (input == "projin") {
0125     auto* tmp = dynamic_cast<TrackletProjectionsMemory*>(memory);
0126     assert(tmp != nullptr);
0127     inputprojs_.push_back(tmp);
0128     return;
0129   }
0130   throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " could not find input: " << input;
0131 }
0132 
0133 void MatchProcessor::execute(unsigned int iSector, double phimin) {
0134   assert(vmstubs_.size() == 1);
0135 
0136   /*
0137     The code is organized in three 'steps' corresponding to the PR, ME, and MC functions. The output from
0138     the PR step is buffered in a 'circular' buffer, and similarly the ME output is put in a circular buffer. 
0139     
0140     The implementation is done in steps, emulating what can be done in firmware. One each step we do:
0141     
0142     1) A projection is read and if there is space it is insert into the inputProjBuffer_
0143     
0144     2) Process next match in the ME - if there is an idle ME the next projection is inserted
0145     
0146     3) Readout match from ME and send to match calculator
0147     
0148   */
0149 
0150   bool print = getName() == "MP_L3PHIC" && iSector == 3;
0151   print = false;
0152 
0153   phimin_ = phimin;
0154 
0155   Tracklet* oldTracklet = nullptr;
0156 
0157   unsigned int countme = 0;
0158   unsigned int countall = 0;
0159   unsigned int countsel = 0;
0160   unsigned int countinputproj = 0;
0161 
0162   unsigned int iprojmem = 0;
0163   while (iprojmem < inputprojs_.size() && inputprojs_[iprojmem]->nTracklets() == 0) {
0164     iprojmem++;
0165   }
0166 
0167   unsigned int iproj = 0;
0168 
0169   inputProjBuffer_.reset();
0170 
0171   for (const auto& inputproj : inputprojs_) {
0172     countinputproj += inputproj->nTracklets();
0173   }
0174 
0175   for (auto& matchengine : matchengines_) {
0176     matchengine.reset();
0177   }
0178 
0179   ProjectionTemp tmpProj_, tmpProj__;
0180   bool good_ = false;
0181   bool good__ = false;
0182 
0183   for (unsigned int istep = 0; istep < settings_.maxStep("MP"); istep++) {
0184     bool projBuffNearFull = inputProjBuffer_.nearfull();
0185 
0186     // This print statement is useful for detailed comparison with the HLS code
0187     // It prints out detailed status information for each clock step
0188     /* 
0189     if (print) {
0190       cout << "istep = "<<istep<<" projBuff: "<<inputProjBuffer_.rptr()<<" "<<inputProjBuffer_.wptr()<<" "<<projBuffNearFull;
0191       unsigned int iMEU = 0;
0192       for (auto& matchengine : matchengines_) {
0193     cout <<" MEU"<<iMEU<<": "<<matchengine.rptr()<<" "<<matchengine.wptr()
0194          <<" "<<matchengine.idle()<<" "<<matchengine.empty()
0195          <<" "<<matchengine.have_()<<matchengine.have__()
0196          <<" "<<matchengine.TCID();
0197     iMEU++;
0198       }
0199       cout << std::endl;
0200     }
0201     */
0202 
0203     bool projdone = false;
0204 
0205     //Step 3
0206     //Check if we have candidate match to process
0207 
0208     unsigned int iMEbest = 0;
0209     int bestTCID = matchengines_[0].TCID();
0210     bool meactive = matchengines_[0].active();
0211     for (unsigned int iME = 1; iME < nMatchEngines_; iME++) {
0212       meactive = meactive || matchengines_[iME].active();
0213       int tcid = matchengines_[iME].TCID();
0214       if (tcid < bestTCID) {
0215         bestTCID = tcid;
0216         iMEbest = iME;
0217       }
0218     }
0219 
0220     if (!matchengines_[iMEbest].empty()) {
0221       std::pair<Tracklet*, const Stub*> candmatch = matchengines_[iMEbest].read();
0222 
0223       const Stub* fpgastub = candmatch.second;
0224       Tracklet* tracklet = candmatch.first;
0225 
0226       if (oldTracklet != nullptr) {
0227         //allow equal here since we can have more than one cadidate match per tracklet projection
0228         //cout << "old new : "<<oldTracklet->TCID()<<" "<<tracklet->TCID()<<" "<<iMEbest<<endl;
0229         assert(oldTracklet->TCID() <= tracklet->TCID());
0230       }
0231       oldTracklet = tracklet;
0232 
0233       bool match = matchCalculator(tracklet, fpgastub, print, istep);
0234 
0235       if (settings_.debugTracklet() && match) {
0236         edm::LogVerbatim("Tracklet") << getName() << " have match";
0237       }
0238 
0239       countall++;
0240       if (match)
0241         countsel++;
0242     }
0243 
0244     //Step 2
0245     //Check if we have ME that can process projection
0246 
0247     bool addedProjection = false;
0248     for (unsigned int iME = 0; iME < nMatchEngines_; iME++) {
0249       if (!matchengines_[iME].idle())
0250         countme++;
0251       //if match engine empty and we have queued projections add to match engine
0252       if ((!addedProjection) && matchengines_[iME].idle() && (!inputProjBuffer_.empty())) {
0253         ProjectionTemp tmpProj = inputProjBuffer_.read();
0254         VMStubsMEMemory* stubmem = vmstubs_[0];
0255 
0256         if (settings_.debugTracklet()) {
0257           edm::LogVerbatim("Tracklet") << getName() << " adding projection to match engine";
0258         }
0259 
0260         int nbins = (1 << N_RZBITS);
0261         if (layerdisk_ >= N_LAYER) {
0262           nbins *= 2;  //twice as many bins in disks (since there are two disks)
0263         }
0264 
0265         matchengines_[iME].init(stubmem,
0266                                 nbins,
0267                                 tmpProj.slot(),
0268                                 tmpProj.iphi(),
0269                                 tmpProj.shift(),
0270                                 tmpProj.projrinv(),
0271                                 tmpProj.projfinerz(),
0272                                 tmpProj.projfinephi(),
0273                                 tmpProj.use(0, 0),
0274                                 tmpProj.use(0, 1),
0275                                 tmpProj.use(1, 0),
0276                                 tmpProj.use(1, 1),
0277                                 tmpProj.isPSseed(),
0278                                 tmpProj.proj(),
0279                                 print && (iME == 0));
0280         addedProjection = true;
0281       } else {
0282         matchengines_[iME].step(print && iME == 0);
0283       }
0284     }
0285 
0286     //Step 1
0287     //First step here checks if we have more input projections to put into
0288     //the input puffer for projections
0289 
0290     if (good__) {
0291       inputProjBuffer_.store(tmpProj__);
0292     }
0293 
0294     good__ = good_;
0295     tmpProj__ = tmpProj_;
0296 
0297     good_ = false;
0298 
0299     if (istep < settings_.maxStep("MP")) {
0300       if (iprojmem < inputprojs_.size()) {
0301         TrackletProjectionsMemory* projMem = inputprojs_[iprojmem];
0302         if (!projBuffNearFull) {
0303           if (settings_.debugTracklet()) {
0304             edm::LogVerbatim("Tracklet") << getName() << " have projection in memory : " << projMem->getName();
0305           }
0306 
0307           Tracklet* proj = projMem->getTracklet(iproj);
0308 
0309           FPGAWord fpgaphi = proj->proj(layerdisk_).fpgaphiproj();
0310 
0311           unsigned int iphi = (fpgaphi.value() >> (fpgaphi.nbits() - nvmbits_)) & (nvmbins_ - 1);
0312 
0313           int nextrabits = 2;
0314           int overlapbits = nvmbits_ + nextrabits;
0315 
0316           unsigned int extrabits = fpgaphi.bits(fpgaphi.nbits() - overlapbits, nextrabits);
0317 
0318           unsigned int ivmPlus = iphi;
0319 
0320           int shift = 0;
0321 
0322           if (extrabits == ((1U << nextrabits) - 1) && iphi != ((1U << settings_.nbitsvmme(layerdisk_)) - 1)) {
0323             shift = 1;
0324             ivmPlus++;
0325           }
0326           unsigned int ivmMinus = iphi;
0327           if (extrabits == 0 && iphi != 0) {
0328             shift = -1;
0329             ivmMinus--;
0330           }
0331 
0332           int projrinv = -1;
0333           if (barrel_) {
0334             FPGAWord phider = proj->proj(layerdisk_).fpgaphiprojder();
0335             projrinv = (1 << (nrinv_ - 1)) - 1 - (phider.value() >> (phider.nbits() - nrinv_));
0336           } else {
0337             //The next lines looks up the predicted bend based on:
0338             // 1 - r projections
0339             // 2 - phi derivative
0340             // 3 - the sign - i.e. if track is forward or backward
0341 
0342             int rindex = (proj->proj(layerdisk_).fpgarzproj().value() >>
0343                           (proj->proj(layerdisk_).fpgarzproj().nbits() - nrbits_)) &
0344                          ((1 << nrbits_) - 1);
0345 
0346             int phiderindex = (proj->proj(layerdisk_).fpgaphiprojder().value() >>
0347                                (proj->proj(layerdisk_).fpgaphiprojder().nbits() - nphiderbits_)) &
0348                               ((1 << nphiderbits_) - 1);
0349 
0350             int signindex = proj->proj(layerdisk_).fpgarzprojder().value() < 0;
0351 
0352             int bendindex = (signindex << (nphiderbits_ + nrbits_)) + (rindex << (nphiderbits_)) + phiderindex;
0353 
0354             projrinv = rinvbendlut_.lookup(bendindex);
0355 
0356             proj->proj(layerdisk_).setBendIndex(projrinv);
0357           }
0358           assert(projrinv >= 0);
0359 
0360           unsigned int slot = proj->proj(layerdisk_).fpgarzbin1projvm().value();
0361           bool second = proj->proj(layerdisk_).fpgarzbin2projvm().value();
0362 
0363           unsigned int projfinephi =
0364               (fpgaphi.value() >> (fpgaphi.nbits() - (nvmbits_ + NFINEPHIBITS))) & ((1 << NFINEPHIBITS) - 1);
0365           int projfinerz = proj->proj(layerdisk_).fpgafinerzvm().value();
0366 
0367           bool isPSseed = proj->PSseed();
0368 
0369           int nbins = (1 << N_RZBITS);
0370           if (layerdisk_ >= N_LAYER) {
0371             nbins *= 2;  //twice as many bins in disks (since there are two disks)
0372           }
0373 
0374           VMStubsMEMemory* stubmem = vmstubs_[0];
0375           bool usefirstMinus = stubmem->nStubsBin(ivmMinus * nbins + slot) != 0;
0376           bool usesecondMinus = (second && (stubmem->nStubsBin(ivmMinus * nbins + slot + 1) != 0));
0377           bool usefirstPlus = ivmPlus != ivmMinus && stubmem->nStubsBin(ivmPlus * nbins + slot) != 0;
0378           bool usesecondPlus = ivmPlus != ivmMinus && (second && (stubmem->nStubsBin(ivmPlus * nbins + slot + 1) != 0));
0379 
0380           good_ = usefirstPlus || usesecondPlus || usefirstMinus || usesecondMinus;
0381 
0382           if (good_) {
0383             ProjectionTemp tmpProj(proj,
0384                                    slot,
0385                                    projrinv,
0386                                    projfinerz,
0387                                    projfinephi,
0388                                    ivmMinus,
0389                                    shift,
0390                                    usefirstMinus,
0391                                    usefirstPlus,
0392                                    usesecondMinus,
0393                                    usesecondPlus,
0394                                    isPSseed);
0395             tmpProj_ = tmpProj;
0396           }
0397 
0398           iproj++;
0399           if (iproj == projMem->nTracklets()) {
0400             iproj = 0;
0401             do {
0402               iprojmem++;
0403             } while (iprojmem < inputprojs_.size() && inputprojs_[iprojmem]->nTracklets() == 0);
0404           }
0405         }
0406 
0407       } else {
0408         projdone = true && !good_ && !good__;
0409       }
0410     }
0411 
0412     //
0413     //  Check if done
0414     //
0415     //
0416     //
0417 
0418     if ((projdone && !meactive) || (istep == settings_.maxStep("MP") - 1)) {
0419       if (settings_.writeMonitorData("MP")) {
0420         globals_->ofstream("matchprocessor.txt") << getName() << " " << istep << " " << countall << " " << countsel
0421                                                  << " " << countme << " " << countinputproj << endl;
0422       }
0423       break;
0424     }
0425   }
0426 
0427   if (settings_.writeMonitorData("MC")) {
0428     globals_->ofstream("matchcalculator.txt") << getName() << " " << countall << " " << countsel << endl;
0429   }
0430 }
0431 
0432 bool MatchProcessor::matchCalculator(Tracklet* tracklet, const Stub* fpgastub, bool, unsigned int) {
0433   const L1TStub* stub = fpgastub->l1tstub();
0434 
0435   if (layerdisk_ < N_LAYER) {
0436     const Projection& proj = tracklet->proj(layerdisk_);
0437     int ir = fpgastub->r().value();
0438     int iphi = proj.fpgaphiproj().value();
0439     int icorr = (ir * proj.fpgaphiprojder().value()) >> icorrshift_;
0440     iphi += icorr;
0441 
0442     int iz = proj.fpgarzproj().value();
0443     int izcor = (ir * proj.fpgarzprojder().value() + (1 << (icorzshift_ - 1))) >> icorzshift_;
0444     iz += izcor;
0445 
0446     int ideltaz = fpgastub->z().value() - iz;
0447     int ideltaphi = (fpgastub->phi().value() - iphi) << phishift_;
0448 
0449     //Floating point calculations
0450 
0451     double phi = stub->phi();
0452     double r = stub->r();
0453     double z = stub->z();
0454 
0455     if (settings_.useapprox()) {
0456       double dphi = reco::reduceRange(phi - fpgastub->phiapprox(phimin_, 0.0));
0457       assert(std::abs(dphi) < 0.001);
0458       phi = fpgastub->phiapprox(phimin_, 0.0);
0459       z = fpgastub->zapprox();
0460       r = fpgastub->rapprox();
0461     }
0462 
0463     if (phi < 0)
0464       phi += 2 * M_PI;
0465     phi -= phimin_;
0466 
0467     double dr = r - settings_.rmean(layerdisk_);
0468     assert(std::abs(dr) < settings_.drmax());
0469 
0470     double dphi = reco::reduceRange(phi - (proj.phiproj() + dr * proj.phiprojder()));
0471 
0472     double dz = z - (proj.rzproj() + dr * proj.rzprojder());
0473 
0474     double dphiapprox = reco::reduceRange(phi - (proj.phiprojapprox() + dr * proj.phiprojderapprox()));
0475 
0476     double dzapprox = z - (proj.rzprojapprox() + dr * proj.rzprojderapprox());
0477 
0478     int seedindex = tracklet->getISeed();
0479 
0480     assert(phimatchcuttable_.lookup(seedindex) > 0);
0481     assert(zmatchcuttable_.lookup(seedindex) > 0);
0482 
0483     if (settings_.bookHistos()) {
0484       bool truthmatch = tracklet->stubtruthmatch(stub);
0485 
0486       HistBase* hists = globals_->histograms();
0487       hists->FillLayerResidual(layerdisk_ + 1,
0488                                seedindex,
0489                                dphiapprox * settings_.rmean(layerdisk_),
0490                                ideltaphi * settings_.kphi1() * settings_.rmean(layerdisk_),
0491                                (ideltaz << dzshift_) * settings_.kz(),
0492                                dz,
0493                                truthmatch);
0494     }
0495 
0496     if (settings_.writeMonitorData("Residuals")) {
0497       double pt = 0.01 * settings_.c() * settings_.bfield() / std::abs(tracklet->rinv());
0498 
0499       globals_->ofstream("layerresiduals.txt")
0500           << layerdisk_ + 1 << " " << seedindex << " " << pt << " "
0501           << ideltaphi * settings_.kphi1() * settings_.rmean(layerdisk_) << " "
0502           << dphiapprox * settings_.rmean(layerdisk_) << " "
0503           << phimatchcuttable_.lookup(seedindex) * settings_.kphi1() * settings_.rmean(layerdisk_) << "   "
0504           << (ideltaz << dzshift_) * settings_.kz() << " " << dz << " "
0505           << zmatchcuttable_.lookup(seedindex) * settings_.kz() << endl;
0506     }
0507 
0508     bool imatch = (std::abs(ideltaphi) <= phimatchcuttable_.lookup(seedindex)) &&
0509                   (std::abs(ideltaz << dzshift_) <= zmatchcuttable_.lookup(seedindex));
0510 
0511     if (settings_.debugTracklet()) {
0512       edm::LogVerbatim("Tracklet") << getName() << " imatch = " << imatch << " ideltaphi cut " << ideltaphi << " "
0513                                    << phimatchcuttable_.lookup(seedindex) << " ideltaz<<dzshift cut "
0514                                    << (ideltaz << dzshift_) << " " << zmatchcuttable_.lookup(seedindex);
0515     }
0516 
0517     //This would catch significant consistency problems in the configuration - helps to debug if there are problems.
0518     if (std::abs(dphi) > 0.5 * settings_.dphisectorHG() || std::abs(dphiapprox) > 0.5 * settings_.dphisectorHG()) {
0519       throw cms::Exception("LogicError") << "WARNING dphi and/or dphiapprox too large : " << dphi << " " << dphiapprox
0520                                          << endl;
0521     }
0522 
0523     if (imatch) {
0524       tracklet->addMatch(layerdisk_,
0525                          ideltaphi,
0526                          ideltaz,
0527                          dphi,
0528                          dz,
0529                          dphiapprox,
0530                          dzapprox,
0531                          (phiregion_ << N_BITSMEMADDRESS) + fpgastub->stubindex().value(),
0532                          fpgastub);
0533 
0534       if (settings_.debugTracklet()) {
0535         edm::LogVerbatim("Tracklet") << "Accepted full match in layer " << getName() << " " << tracklet;
0536       }
0537 
0538       int iSeed = tracklet->getISeed();
0539       assert(fullmatches_[iSeed] != nullptr);
0540       fullmatches_[iSeed]->addMatch(tracklet, fpgastub);
0541 
0542       return true;
0543     } else {
0544       return false;
0545     }
0546   } else {  //disk matches
0547 
0548     //check that stubs and projections in same half of detector
0549     assert(stub->z() * tracklet->t() > 0.0);
0550 
0551     int sign = (tracklet->t() > 0.0) ? 1 : -1;
0552     int disk = sign * (layerdisk_ - N_LAYER + 1);
0553     assert(disk != 0);
0554 
0555     //Perform integer calculations here
0556 
0557     int iz = fpgastub->z().value();
0558 
0559     const Projection& proj = tracklet->proj(layerdisk_);
0560 
0561     int iphi = proj.fpgaphiproj().value();
0562     int iphicorr = (iz * proj.fpgaphiprojder().value()) >> icorrshift_;
0563 
0564     iphi += iphicorr;
0565 
0566     int ir = proj.fpgarzproj().value();
0567     int ircorr = (iz * proj.fpgarzprojder().value()) >> icorzshift_;
0568     ir += ircorr;
0569 
0570     int ideltaphi = fpgastub->phi().value() - iphi;
0571 
0572     int irstub = fpgastub->r().value();
0573     int ialphafact = 0;
0574     if (!stub->isPSmodule()) {
0575       assert(irstub < (int)N_DSS_MOD * 2);
0576       if (layerdisk_ - N_LAYER <= 1) {
0577         ialphafact = ialphafactinner_[irstub];
0578         irstub = settings_.rDSSinner(irstub) / settings_.kr();
0579       } else {
0580         ialphafact = ialphafactouter_[irstub];
0581         irstub = settings_.rDSSouter(irstub) / settings_.kr();
0582       }
0583     }
0584 
0585     int ideltar = (irstub * settings_.kr()) / settings_.krprojshiftdisk() - ir;
0586 
0587     if (!stub->isPSmodule()) {
0588       int ialpha = fpgastub->alpha().value();
0589       int iphialphacor = ((ideltar * ialpha * ialphafact) >> settings_.alphashift());
0590       ideltaphi += iphialphacor;
0591     }
0592 
0593     //Perform floating point calculations here
0594 
0595     double phi = stub->phi();
0596     double z = stub->z();
0597     double r = stub->r();
0598 
0599     if (settings_.useapprox()) {
0600       double dphi = reco::reduceRange(phi - fpgastub->phiapprox(phimin_, 0.0));
0601       assert(std::abs(dphi) < 0.001);
0602       phi = fpgastub->phiapprox(phimin_, 0.0);
0603       z = fpgastub->zapprox();
0604       r = fpgastub->rapprox();
0605     }
0606 
0607     if (phi < 0)
0608       phi += 2 * M_PI;
0609     phi -= phimin_;
0610 
0611     double dz = z - sign * settings_.zmean(layerdisk_ - N_LAYER);
0612 
0613     if (std::abs(dz) > settings_.dzmax()) {
0614       edm::LogProblem("Tracklet") << __FILE__ << ":" << __LINE__ << " " << name_ << " " << tracklet->getISeed();
0615       edm::LogProblem("Tracklet") << "stub " << stub->z() << " disk " << disk << " " << dz;
0616       assert(std::abs(dz) < settings_.dzmax());
0617     }
0618 
0619     double phiproj = proj.phiproj() + dz * proj.phiprojder();
0620     double rproj = proj.rzproj() + dz * proj.rzprojder();
0621     double deltar = r - rproj;
0622 
0623     double dr = stub->r() - rproj;
0624     double drapprox = stub->r() - (proj.rzprojapprox() + dz * proj.rzprojderapprox());
0625 
0626     double dphi = reco::reduceRange(phi - phiproj);
0627 
0628     double dphiapprox = reco::reduceRange(phi - (proj.phiprojapprox() + dz * proj.phiprojderapprox()));
0629 
0630     double drphi = dphi * stub->r();
0631     double drphiapprox = dphiapprox * stub->r();
0632 
0633     if (!stub->isPSmodule()) {
0634       double alphanorm = stub->alphanorm();
0635       dphi += dr * alphanorm * settings_.half2SmoduleWidth() / stub->r2();
0636       ;
0637       dphiapprox += drapprox * alphanorm * settings_.half2SmoduleWidth() / stub->r2();
0638 
0639       drphi += dr * alphanorm * settings_.half2SmoduleWidth() / stub->r();
0640       drphiapprox += dr * alphanorm * settings_.half2SmoduleWidth() / stub->r();
0641     }
0642 
0643     int seedindex = tracklet->getISeed();
0644 
0645     int idrphicut = rphicutPStable_.lookup(seedindex);
0646     int idrcut = rcutPStable_.lookup(seedindex);
0647     if (!stub->isPSmodule()) {
0648       idrphicut = rphicut2Stable_.lookup(seedindex);
0649       idrcut = rcut2Stable_.lookup(seedindex);
0650     }
0651 
0652     double drphicut = idrphicut * settings_.kphi() * settings_.kr();
0653     double drcut = idrcut * settings_.krprojshiftdisk();
0654 
0655     if (settings_.writeMonitorData("Residuals")) {
0656       double pt = 0.01 * settings_.c() * settings_.bfield() / std::abs(tracklet->rinv());
0657 
0658       globals_->ofstream("diskresiduals.txt")
0659           << layerdisk_ - N_LAYER + 1 << " " << stub->isPSmodule() << " " << tracklet->layer() << " "
0660           << abs(tracklet->disk()) << " " << pt << " " << ideltaphi * settings_.kphi() * stub->r() << " " << drphiapprox
0661           << " " << drphicut << " " << ideltar * settings_.krprojshiftdisk() << " " << deltar << " " << drcut << " "
0662           << endl;
0663     }
0664 
0665     bool match = (std::abs(drphi) < drphicut) && (std::abs(deltar) < drcut);
0666     bool imatch = (std::abs(ideltaphi * irstub) < idrphicut) && (std::abs(ideltar) < idrcut);
0667 
0668     if (settings_.debugTracklet()) {
0669       edm::LogVerbatim("Tracklet") << "imatch match disk: " << imatch << " " << match << " " << std::abs(ideltaphi)
0670                                    << " " << drphicut / (settings_.kphi() * stub->r()) << " " << std::abs(ideltar)
0671                                    << " " << drcut / settings_.krprojshiftdisk() << " r = " << stub->r();
0672     }
0673 
0674     if (imatch) {
0675       if (settings_.debugTracklet()) {
0676         edm::LogVerbatim("Tracklet") << "MatchCalculator found match in disk " << getName();
0677       }
0678 
0679       if (std::abs(dphi) >= third * settings_.dphisectorHG()) {
0680         edm::LogPrint("Tracklet") << "dphi " << dphi << " ISeed " << tracklet->getISeed();
0681       }
0682       assert(std::abs(dphi) < third * settings_.dphisectorHG());
0683       assert(std::abs(dphiapprox) < third * settings_.dphisectorHG());
0684 
0685       tracklet->addMatch(layerdisk_,
0686                          ideltaphi,
0687                          ideltar,
0688                          drphi / stub->r(),
0689                          dr,
0690                          drphiapprox / stub->r(),
0691                          drapprox,
0692                          (phiregion_ << N_BITSMEMADDRESS) + fpgastub->stubindex().value(),
0693                          fpgastub);
0694 
0695       if (settings_.debugTracklet()) {
0696         edm::LogVerbatim("Tracklet") << "Accepted full match in disk " << getName() << " " << tracklet;
0697       }
0698 
0699       int iSeed = tracklet->getISeed();
0700       assert(fullmatches_[iSeed] != nullptr);
0701       fullmatches_[iSeed]->addMatch(tracklet, fpgastub);
0702 
0703       return true;
0704     } else {
0705       return false;
0706     }
0707   }
0708 }