Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:22:00

0001 #include "L1Trigger/TrackFindingTracklet/interface/MatchEngine.h"
0002 #include "L1Trigger/TrackFindingTracklet/interface/Globals.h"
0003 #include "L1Trigger/TrackFindingTracklet/interface/Util.h"
0004 #include "L1Trigger/TrackFindingTracklet/interface/MemoryBase.h"
0005 #include "L1Trigger/TrackFindingTracklet/interface/VMStubsMEMemory.h"
0006 #include "L1Trigger/TrackFindingTracklet/interface/VMProjectionsMemory.h"
0007 #include "L1Trigger/TrackFindingTracklet/interface/CandidateMatchMemory.h"
0008 #include "L1Trigger/TrackFindingTracklet/interface/Tracklet.h"
0009 
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/Utilities/interface/Exception.h"
0012 
0013 #include <filesystem>
0014 
0015 using namespace std;
0016 using namespace trklet;
0017 
0018 MatchEngine::MatchEngine(string name, Settings const& settings, Globals* global)
0019     : ProcessBase(name, settings, global), luttable_(settings) {
0020   layerdisk_ = initLayerDisk(3);
0021 
0022   barrel_ = layerdisk_ < N_LAYER;
0023 
0024   luttable_.initBendMatch(layerdisk_);
0025 
0026   nrinv_ = NRINVBITS;
0027 }
0028 
0029 void MatchEngine::addOutput(MemoryBase* memory, string output) {
0030   if (settings_.writetrace()) {
0031     edm::LogVerbatim("Tracklet") << "In " << name_ << " adding output to " << memory->getName() << " to output "
0032                                  << output;
0033   }
0034   if (output == "matchout") {
0035     auto* tmp = dynamic_cast<CandidateMatchMemory*>(memory);
0036     assert(tmp != nullptr);
0037     candmatches_ = tmp;
0038     return;
0039   }
0040   throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " could not find output: " << output;
0041 }
0042 
0043 void MatchEngine::addInput(MemoryBase* memory, string input) {
0044   if (settings_.writetrace()) {
0045     edm::LogVerbatim("Tracklet") << "In " << name_ << " adding input from " << memory->getName() << " to input "
0046                                  << input;
0047   }
0048   if (input == "vmstubin") {
0049     auto* tmp = dynamic_cast<VMStubsMEMemory*>(memory);
0050     assert(tmp != nullptr);
0051     vmstubs_ = tmp;
0052     return;
0053   }
0054   if (input == "vmprojin") {
0055     auto* tmp = dynamic_cast<VMProjectionsMemory*>(memory);
0056     assert(tmp != nullptr);
0057     vmprojs_ = tmp;
0058     return;
0059   }
0060   throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " could not find input: " << input;
0061 }
0062 
0063 void MatchEngine::execute(unsigned int iSector) {
0064   unsigned int countall = 0;
0065   unsigned int countpass = 0;
0066 
0067   //bool print = (iSector == 3 && getName() == "ME_L3PHIC20");
0068   bool print = false;
0069 
0070   constexpr unsigned int kNBitsBuffer = 3;
0071 
0072   int writeindex = 0;
0073   int readindex = 0;
0074   std::pair<int, int> projbuffer[1 << kNBitsBuffer];  //iproj zbin
0075 
0076   //The next projection to read, the number of projections and flag if we have more projections to read
0077   int iproj = 0;
0078   int nproj = vmprojs_->nTracklets();
0079   bool moreproj = iproj < nproj;
0080 
0081   //Projection that is read from the buffer and compared to stubs
0082   int rzbin = 0;
0083   int projfinerz = 0;
0084   int projfinerzadj = 0;
0085   unsigned int projfinephi = 0;
0086 
0087   int projindex = 0;
0088   int projrinv = 0;
0089   bool isPSseed = false;
0090 
0091   //Number of stubs for current zbin and the stub being processed on this clock
0092   int nstubs = 0;
0093   int istub = 0;
0094 
0095   //Main processing loops starts here
0096   for (unsigned int istep = 0; istep < settings_.maxStep("ME"); istep++) {
0097     countall++;
0098 
0099     int writeindexplus = (writeindex + 1) % (1 << kNBitsBuffer);
0100     int writeindexplusplus = (writeindex + 2) % (1 << kNBitsBuffer);
0101 
0102     //Determine if buffer is full - or near full as a projection
0103     //can point to two z bins we might fill two slots in the buffer
0104     bool bufferfull = (writeindexplus == readindex) || (writeindexplusplus == readindex);
0105 
0106     //Determin if buffer is empty
0107     bool buffernotempty = (writeindex != readindex);
0108 
0109     //If we have more projections and the buffer is not full we read
0110     //next projection and put in buffer if there are stubs in the
0111     //memory the projection points to
0112 
0113     if ((!moreproj) && (!buffernotempty))
0114       break;
0115 
0116     if (moreproj && (!bufferfull)) {
0117       Tracklet* proj = vmprojs_->getTracklet(iproj);
0118 
0119       int iprojtmp = iproj;
0120 
0121       iproj++;
0122       moreproj = iproj < nproj;
0123 
0124       unsigned int rzfirst = proj->proj(layerdisk_).fpgarzbin1projvm().value();
0125       unsigned int rzlast = rzfirst;
0126 
0127       bool second = proj->proj(layerdisk_).fpgarzbin2projvm().value();
0128 
0129       if (second)
0130         rzlast += 1;
0131 
0132       //Check if there are stubs in the memory
0133       int nstubfirst = vmstubs_->nStubsBin(rzfirst);
0134       int nstublast = vmstubs_->nStubsBin(rzlast);
0135       bool savefirst = nstubfirst != 0;
0136       bool savelast = second && (nstublast != 0);
0137 
0138       int writeindextmp = writeindex;
0139       int writeindextmpplus = (writeindex + 1) % (1 << kNBitsBuffer);
0140 
0141       if (savefirst && savelast) {
0142         writeindex = writeindexplusplus;
0143       } else if (savefirst || savelast) {
0144         writeindex = writeindexplus;
0145       }
0146 
0147       if (savefirst) {  //TODO for HLS (make code logic simpler)
0148         std::pair<int, int> tmp(iprojtmp, rzfirst);
0149         projbuffer[writeindextmp] = tmp;
0150       }
0151       if (savelast) {
0152         std::pair<int, int> tmp(iprojtmp, rzlast + 100);  //TODO for HLS (fix flagging that this is second bin)
0153         if (savefirst) {
0154           projbuffer[writeindextmpplus] = tmp;
0155         } else {
0156           projbuffer[writeindextmp] = tmp;
0157         }
0158       }
0159     }
0160 
0161     //If the buffer is not empty we have a projection that we need to process.
0162 
0163     if (buffernotempty) {
0164       int istubtmp = istub;
0165 
0166       //New projection
0167       if (istub == 0) {
0168         projindex = projbuffer[readindex].first;
0169         rzbin = projbuffer[readindex].second;
0170         bool second = false;
0171         if (rzbin >= 100) {
0172           rzbin -= 100;
0173           second = true;
0174         }
0175 
0176         Tracklet* proj = vmprojs_->getTracklet(projindex);
0177 
0178         FPGAWord fpgafinephi = proj->proj(layerdisk_).fpgafinephivm();
0179 
0180         projfinephi = fpgafinephi.value();
0181 
0182         nstubs = vmstubs_->nStubsBin(rzbin);
0183 
0184         projfinerz = proj->proj(layerdisk_).fpgafinerzvm().value();
0185 
0186         projrinv = barrel_ ? ((1 << (nrinv_ - 1)) + ((-2 * proj->proj(layerdisk_).fpgaphiprojder().value()) >>
0187                                                      (proj->proj(layerdisk_).fpgaphiprojder().nbits() - (nrinv_ - 1))))
0188                            : proj->proj(layerdisk_).getBendIndex().value();
0189         assert(projrinv >= 0);
0190         if (settings_.extended() && projrinv == (1 << nrinv_)) {
0191           if (settings_.debugTracklet()) {
0192             edm::LogVerbatim("Tracklet") << "Extended tracking, projrinv:" << projrinv;
0193           }
0194           projrinv = (1 << nrinv_) - 1;
0195         }
0196         assert(projrinv < (1 << nrinv_));
0197 
0198         isPSseed = proj->PSseed();
0199 
0200         //Calculate fine z position
0201         if (second) {
0202           projfinerzadj = projfinerz - (1 << NFINERZBITS);
0203         } else {
0204           projfinerzadj = projfinerz;
0205         }
0206         if (nstubs == 1) {
0207           istub = 0;
0208           readindex = (readindex + 1) % (1 << kNBitsBuffer);
0209         } else {
0210           istub++;
0211         }
0212       } else {
0213         //Check if last stub, if so, go to next buffer entry
0214         if (istub + 1 >= nstubs) {
0215           istub = 0;
0216           readindex = (readindex + 1) % (1 << kNBitsBuffer);
0217         } else {
0218           istub++;
0219         }
0220       }
0221 
0222       //Read vmstub memory and extract data fields
0223       const VMStubME& vmstub = vmstubs_->getVMStubMEBin(rzbin, istubtmp);
0224 
0225       int stubfinerz = vmstub.finerz().value();
0226 
0227       bool isPSmodule;
0228 
0229       if (barrel_) {
0230         isPSmodule = layerdisk_ < N_PSLAYER;
0231       } else {
0232         if (layerdisk_ < N_LAYER + 2) {
0233           isPSmodule = ((rzbin & 7) < 3) || ((rzbin & 7) == 3 && stubfinerz <= 3);
0234         } else {
0235           isPSmodule = ((rzbin & 7) < 3) || ((rzbin & 7) == 3 && stubfinerz <= 2);
0236         }
0237       }
0238 
0239       assert(isPSmodule == vmstub.isPSmodule());
0240 
0241       int nbits = isPSmodule ? N_BENDBITS_PS : N_BENDBITS_2S;
0242 
0243       int deltaphi = projfinephi - vmstub.finephi().value();
0244 
0245       constexpr int mindeltaphicut = 3;
0246       constexpr int maxdeltaphicut = 5;
0247       bool passphi = (std::abs(deltaphi) < mindeltaphicut) || (std::abs(deltaphi) > maxdeltaphicut);
0248 
0249       unsigned int index = (projrinv << nbits) + vmstub.bend().value();
0250       if (!barrel_ && isPSmodule) {
0251         index += (1 << (nrinv_ + N_BENDBITS_2S));
0252       }
0253 
0254       //Check if stub z position consistent
0255       int idrz = stubfinerz - projfinerzadj;
0256       bool passz;
0257 
0258       if (barrel_) {
0259         if (isPSseed) {
0260           constexpr int drzcut = 1;
0261           passz = std::abs(idrz) <= drzcut;
0262         } else {
0263           constexpr int drzcut = 5;
0264           passz = std::abs(idrz) <= drzcut;
0265         }
0266       } else {
0267         if (isPSmodule) {
0268           constexpr int drzcut = 1;
0269           passz = std::abs(idrz) <= drzcut;
0270         } else {
0271           constexpr int drzcut = 3;
0272           passz = std::abs(idrz) <= drzcut;
0273         }
0274       }
0275 
0276       if (print) {
0277         edm::LogVerbatim("Tracklet") << "istep index : " << istep << " " << index << " " << vmstub.bend().value()
0278                                      << " rzbin istubtmp : " << rzbin << " " << istubtmp << " dz " << stubfinerz << " "
0279                                      << projfinerzadj << "  dphi: " << deltaphi;
0280       }
0281 
0282       //Check if stub bend and proj rinv consistent
0283       if (passz && passphi) {
0284         if (luttable_.lookup(index)) {
0285           Tracklet* proj = vmprojs_->getTracklet(projindex);
0286           std::pair<Tracklet*, int> tmp(proj, vmprojs_->getAllProjIndex(projindex));
0287           if (settings_.writeMonitorData("Seeds")) {
0288             ofstream fout("seeds.txt", ofstream::app);
0289             fout << __FILE__ << ":" << __LINE__ << " " << name_ << " " << proj->getISeed() << endl;
0290             fout.close();
0291           }
0292           candmatches_->addMatch(tmp, vmstub.stub());
0293           countpass++;
0294         }
0295       }
0296     }
0297   }
0298 
0299   if (settings_.writeMonitorData("ME")) {
0300     globals_->ofstream("matchengine.txt") << getName() << " " << countall << " " << countpass << endl;
0301   }
0302 }