Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "L1Trigger/TrackFindingTracklet/interface/TrackletProcessorDisplaced.h"
0002 #include "L1Trigger/TrackFindingTracklet/interface/Settings.h"
0003 #include "L1Trigger/TrackFindingTracklet/interface/Globals.h"
0004 #include "L1Trigger/TrackFindingTracklet/interface/AllStubsMemory.h"
0005 #include "L1Trigger/TrackFindingTracklet/interface/AllInnerStubsMemory.h"
0006 #include "L1Trigger/TrackFindingTracklet/interface/StubPairsMemory.h"
0007 #include "L1Trigger/TrackFindingTracklet/interface/Tracklet.h"
0008 #include "L1Trigger/TrackFindingTracklet/interface/Util.h"
0009 #include "L1Trigger/TrackFindingTracklet/interface/IMATH_TrackletCalculator.h"
0010 
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "FWCore/Utilities/interface/Exception.h"
0013 #include "DataFormats/Math/interface/deltaPhi.h"
0014 
0015 #include <utility>
0016 #include <tuple>
0017 
0018 using namespace std;
0019 using namespace trklet;
0020 
0021 // TrackletProcessorDisplaced
0022 //
0023 // This module takes in collections of stubs within a phi region and a
0024 // displaced seed name and tries to create that displaced seed out of the stubs
0025 //
0026 // Update: Claire Savard, Nov. 2024
0027 
0028 TrackletProcessorDisplaced::TrackletProcessorDisplaced(string name, Settings const& settings, Globals* globals)
0029     : TrackletCalculatorDisplaced(name, settings, globals),
0030       trpbuffer_(CircularBuffer<TrpEData>(3), 0, 0, 0, 0),
0031       innerTable_(settings),
0032       innerThirdTable_(settings) {
0033   innerallstubs_.clear();
0034   middleallstubs_.clear();
0035   outerallstubs_.clear();
0036   innervmstubs_.clear();
0037   outervmstubs_.clear();
0038 
0039   // set layer/disk types based on input seed name
0040   initLayerDisksandISeedDisp(layerdisk1_, layerdisk2_, layerdisk3_, iSeed_);
0041 
0042   // get projection tables
0043   unsigned int region = name.back() - 'A';
0044   innerTable_.initVMRTable(
0045       layerdisk1_, TrackletLUT::VMRTableType::inner, region, false);  //projection to next layer/disk
0046   innerThirdTable_.initVMRTable(
0047       layerdisk1_, TrackletLUT::VMRTableType::innerthird, region, false);  //projection to third layer/disk
0048 
0049   nbitszfinebintable_ = settings_.vmrlutzbits(layerdisk1_);
0050   nbitsrfinebintable_ = settings_.vmrlutrbits(layerdisk1_);
0051 
0052   for (unsigned int ilayer = 0; ilayer < N_LAYER; ilayer++) {
0053     vector<TrackletProjectionsMemory*> tmp(settings_.nallstubs(ilayer), nullptr);
0054     trackletprojlayers_.push_back(tmp);
0055   }
0056 
0057   for (unsigned int idisk = 0; idisk < N_DISK; idisk++) {
0058     vector<TrackletProjectionsMemory*> tmp(settings_.nallstubs(idisk + N_LAYER), nullptr);
0059     trackletprojdisks_.push_back(tmp);
0060   }
0061 
0062   // set TC index
0063   iTC_ = region;
0064   TCIndex_ = (iSeed_ << settings.nbitsseed()) + iTC_;
0065 
0066   maxStep_ = settings_.maxStep("TPD");
0067 }
0068 
0069 void TrackletProcessorDisplaced::addOutputProjection(TrackletProjectionsMemory*& outputProj, MemoryBase* memory) {
0070   outputProj = dynamic_cast<TrackletProjectionsMemory*>(memory);
0071   assert(outputProj != nullptr);
0072 }
0073 
0074 void TrackletProcessorDisplaced::addOutput(MemoryBase* memory, string output) {
0075   if (settings_.writetrace()) {
0076     edm::LogVerbatim("Tracklet") << "In " << name_ << " adding output to " << memory->getName() << " to output "
0077                                  << output;
0078   }
0079 
0080   if (output == "trackpar") {
0081     auto* tmp = dynamic_cast<TrackletParametersMemory*>(memory);
0082     assert(tmp != nullptr);
0083     trackletpars_ = tmp;
0084     return;
0085   }
0086 
0087   if (output.substr(0, 7) == "projout") {
0088     //output is on the form 'projoutL2PHIC' or 'projoutD3PHIB'
0089     auto* tmp = dynamic_cast<TrackletProjectionsMemory*>(memory);
0090     assert(tmp != nullptr);
0091 
0092     constexpr unsigned layerdiskPosInprojout = 8;
0093     constexpr unsigned phiPosInprojout = 12;
0094 
0095     unsigned int layerdisk = output[layerdiskPosInprojout] - '1';  //layer or disk counting from 0
0096     unsigned int phiregion = output[phiPosInprojout] - 'A';        //phiregion counting from 0
0097 
0098     if (output[7] == 'L') {
0099       assert(layerdisk < N_LAYER);
0100       assert(phiregion < trackletprojlayers_[layerdisk].size());
0101       //check that phiregion not already initialized
0102       assert(trackletprojlayers_[layerdisk][phiregion] == nullptr);
0103       trackletprojlayers_[layerdisk][phiregion] = tmp;
0104       return;
0105     }
0106 
0107     if (output[7] == 'D') {
0108       assert(layerdisk < N_DISK);
0109       assert(phiregion < trackletprojdisks_[layerdisk].size());
0110       //check that phiregion not already initialized
0111       assert(trackletprojdisks_[layerdisk][phiregion] == nullptr);
0112       trackletprojdisks_[layerdisk][phiregion] = tmp;
0113       return;
0114     }
0115   }
0116 
0117   throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " Could not find output : " << output;
0118 }
0119 
0120 void TrackletProcessorDisplaced::addInput(MemoryBase* memory, string input) {
0121   if (settings_.writetrace()) {
0122     edm::LogVerbatim("Tracklet") << "In " << name_ << " adding input from " << memory->getName() << " to input "
0123                                  << input;
0124   }
0125 
0126   if (input == "thirdallstubin") {
0127     auto* tmp = dynamic_cast<AllStubsMemory*>(memory);
0128     assert(tmp != nullptr);
0129     innerallstubs_.push_back(tmp);
0130     return;
0131   }
0132   if (input == "firstallstubin") {
0133     auto* tmp = dynamic_cast<AllStubsMemory*>(memory);
0134     assert(tmp != nullptr);
0135     middleallstubs_.push_back(tmp);
0136     return;
0137   }
0138   if (input == "secondallstubin") {
0139     auto* tmp = dynamic_cast<AllStubsMemory*>(memory);
0140     assert(tmp != nullptr);
0141     outerallstubs_.push_back(tmp);
0142     return;
0143   }
0144   if (input == "thirdvmstubin") {
0145     auto* tmp = dynamic_cast<VMStubsTEMemory*>(memory);
0146     assert(tmp != nullptr);
0147     innervmstubs_.push_back(tmp);
0148     return;
0149   }
0150   if (input == "secondvmstubin") {
0151     auto* tmp = dynamic_cast<VMStubsTEMemory*>(memory);
0152     assert(tmp != nullptr);
0153     outervmstubs_.push_back(tmp);
0154     return;
0155   }
0156 
0157   throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " Could not find input : " << input;
0158 }
0159 
0160 void TrackletProcessorDisplaced::execute(unsigned int iSector, double phimin, double phimax) {
0161   phimin_ = phimin;
0162   phimax_ = phimax;
0163   iSector_ = iSector;
0164 
0165   unsigned int countall = 0;
0166   unsigned int countsel = 0;
0167   int donecount = 0;
0168 
0169   // set the triplet engine units and buffer
0170   TripletEngineUnit trpunit(&settings_, layerdisk1_, layerdisk2_, layerdisk3_, iSeed_, innervmstubs_, outervmstubs_);
0171   trpunits_.resize(settings_.trpunits(iSeed_), trpunit);
0172   trpbuffer_ = tuple<CircularBuffer<TrpEData>, unsigned int, unsigned int, unsigned int, unsigned int>(
0173       CircularBuffer<TrpEData>(3), 0, 0, 0, middleallstubs_.size());
0174 
0175   // reset the trpunits
0176   for (auto& trpunit : trpunits_) {
0177     trpunit.reset();
0178   }
0179 
0180   // reset the tebuffer
0181   std::get<0>(trpbuffer_).reset();
0182   std::get<1>(trpbuffer_) = 0;
0183   std::get<2>(trpbuffer_) = std::get<3>(trpbuffer_);
0184 
0185   TrpEData trpdata;
0186   TrpEData trpdata__;
0187   TrpEData trpdata___;
0188   bool goodtrpdata = false;
0189   bool goodtrpdata__ = false;
0190   bool goodtrpdata___ = false;
0191 
0192   bool trpbuffernearfull;
0193   for (unsigned int istep = 0; istep < maxStep_; istep++) {
0194     CircularBuffer<TrpEData>& trpdatabuffer = std::get<0>(trpbuffer_);
0195     trpbuffernearfull = trpdatabuffer.nearfull();
0196 
0197     //
0198     // First block here checks if there is a trpunit with data that should be used
0199     // to calculate the tracklet parameters
0200     //
0201 
0202     // set pointer to the last filled trpunit
0203     TripletEngineUnit* trpunitptr = nullptr;
0204     for (auto& trpunit : trpunits_) {
0205       trpunit.setNearFull();
0206       if (!trpunit.empty()) {
0207         trpunitptr = &trpunit;
0208       }
0209     }
0210 
0211     if (trpunitptr != nullptr) {
0212       auto stubtriplet = trpunitptr->read();
0213 
0214       countall++;
0215 
0216       const Stub* innerFPGAStub = std::get<0>(stubtriplet);
0217       const Stub* middleFPGAStub = std::get<1>(stubtriplet);
0218       const Stub* outerFPGAStub = std::get<2>(stubtriplet);
0219 
0220       const L1TStub* innerStub = innerFPGAStub->l1tstub();
0221       const L1TStub* middleStub = middleFPGAStub->l1tstub();
0222       const L1TStub* outerStub = outerFPGAStub->l1tstub();
0223 
0224       if (settings_.debugTracklet()) {
0225         edm::LogVerbatim("Tracklet") << "TrackletProcessorDisplaced execute " << getName() << "[" << iSector_ << "]";
0226       }
0227 
0228       // check if the seed made from the 3 stubs is valid
0229       bool accept = false;
0230       if (iSeed_ == Seed::L2L3L4 || iSeed_ == Seed::L4L5L6)
0231         accept = LLLSeeding(innerFPGAStub, innerStub, middleFPGAStub, middleStub, outerFPGAStub, outerStub);
0232       else if (iSeed_ == Seed::L2L3D1)
0233         accept = LLDSeeding(innerFPGAStub, innerStub, middleFPGAStub, middleStub, outerFPGAStub, outerStub);
0234       else if (iSeed_ == Seed::D1D2L2)
0235         accept = DDLSeeding(innerFPGAStub, innerStub, middleFPGAStub, middleStub, outerFPGAStub, outerStub);
0236 
0237       if (accept)
0238         countsel++;
0239 
0240       if (trackletpars_->nTracklets() >= settings_.ntrackletmax()) {
0241         edm::LogVerbatim("Tracklet") << "Will break on number of tracklets in " << getName();
0242         assert(0);
0243         break;
0244       }
0245 
0246       if (settings_.debugTracklet()) {
0247         edm::LogVerbatim("Tracklet") << "TrackletProcessor execute done";
0248       }
0249     }
0250 
0251     //
0252     // The second block fills the trpunit if data in buffer and process TripletEngineUnit step
0253     //
0254     //
0255 
0256     bool notemptytrpbuffer = !trpdatabuffer.empty();
0257     for (auto& trpunit : trpunits_) {
0258       if (trpunit.idle() && notemptytrpbuffer) {  // only fill one idle unit every step
0259         trpunit.init(std::get<0>(trpbuffer_).read());
0260         notemptytrpbuffer = false;  //prevent initializing another triplet engine unit
0261       }
0262       trpunit.step();
0263     }
0264 
0265     //
0266     // The third block here checks if we have input stubs to process
0267     //
0268     //
0269 
0270     if (goodtrpdata___)
0271       trpdatabuffer.store(trpdata___);
0272     goodtrpdata = false;
0273 
0274     unsigned int& istub = std::get<1>(trpbuffer_);
0275     unsigned int& midmem = std::get<2>(trpbuffer_);
0276     unsigned int midmemend = std::get<4>(trpbuffer_);
0277 
0278     if ((!trpbuffernearfull) && midmem < midmemend && istub < middleallstubs_[midmem]->nStubs()) {
0279       const Stub* stub = middleallstubs_[midmem]->getStub(istub);
0280 
0281       if (settings_.debugTracklet()) {
0282         edm::LogVerbatim("Tracklet") << "In " << getName() << " have middle stub";
0283       }
0284 
0285       // get r/z index of the middle stub
0286       int indexz = (((1 << (stub->z().nbits() - 1)) + stub->z().value()) >> (stub->z().nbits() - nbitszfinebintable_));
0287       int indexr = -1;
0288       bool negdisk = (stub->disk().value() < 0);  // check if disk in negative z region
0289       if (layerdisk1_ >= LayerDisk::D1) {         // if a disk
0290         if (negdisk)
0291           indexz = (1 << nbitszfinebintable_) - indexz;
0292         indexr = stub->r().value();
0293         if (stub->isPSmodule()) {
0294           indexr = stub->r().value() >> (stub->r().nbits() - nbitsrfinebintable_);
0295         }
0296       } else {  // else a layer
0297         indexr = (((1 << (stub->r().nbits() - 1)) + stub->r().value()) >> (stub->r().nbits() - nbitsrfinebintable_));
0298       }
0299 
0300       // create lookupbits that define projections from middle stub
0301       int lutval = -1;
0302       const auto& lutshift = innerTable_.nbits();
0303       lutval = innerTable_.lookup((indexz << nbitsrfinebintable_) + indexr);
0304       int lutval2 = innerThirdTable_.lookup((indexz << nbitsrfinebintable_) + indexr);
0305       if (lutval != -1 && lutval2 != -1)
0306         lutval += (lutval2 << lutshift);
0307 
0308       if (lutval != -1) {
0309         unsigned int lutwidth = settings_.lutwidthtabextended(0, iSeed_);
0310         FPGAWord lookupbits(lutval, lutwidth, true, __LINE__, __FILE__);
0311 
0312         // get r/z bins for projection into outer layer/disk
0313         int nbitsrzbin_out = N_RZBITS;
0314         if (iSeed_ == Seed::D1D2L2)
0315           nbitsrzbin_out--;
0316         int rzbinfirst_out = lookupbits.bits(0, NFINERZBITS);
0317         int rzdiffmax_out = lookupbits.bits(NFINERZBITS + 1 + nbitsrzbin_out, NFINERZBITS);
0318         int start_out = lookupbits.bits(NFINERZBITS + 1, nbitsrzbin_out);  // first rz bin projection
0319         int next_out = lookupbits.bits(NFINERZBITS, 1);
0320         if (iSeed_ == Seed::D1D2L2 && negdisk)  // if projecting into disk
0321           start_out += (1 << nbitsrzbin_out);
0322         int last_out = start_out + next_out;  // last rz bin projection
0323 
0324         // get r/z bins for projection into third (inner) layer/disk
0325         int nbitsrzbin_in = N_RZBITS;
0326         int start_in = lookupbits.bits(lutshift + NFINERZBITS + 1, nbitsrzbin_in);  // first rz bin projection
0327         int next_in = lookupbits.bits(lutshift + NFINERZBITS, 1);
0328         if (iSeed_ == Seed::D1D2L2 && negdisk)  // if projecting from disk into layer
0329           start_in = settings_.NLONGVMBINS() - 1 - start_in - next_in;
0330         int last_in = start_in + next_in;  // last rz bin projection
0331 
0332         // fill trpdata with projection info of middle stub
0333         trpdata.stub_ = stub;
0334         trpdata.rzbinfirst_out_ = rzbinfirst_out;
0335         trpdata.rzdiffmax_out_ = rzdiffmax_out;
0336         trpdata.start_out_ = start_out;
0337         trpdata.start_in_ = start_in;
0338 
0339         // fill projection bins info for single engine unit
0340         trpdata.projbin_out_.clear();
0341         trpdata.projbin_in_.clear();
0342         for (int ibin_out = start_out; ibin_out <= last_out; ibin_out++) {
0343           for (unsigned int outmem = 0; outmem < outervmstubs_.size(); outmem++) {
0344             int nstubs_out = outervmstubs_[outmem]->nVMStubsBinned(ibin_out);
0345             if (nstubs_out > 0)
0346               trpdata.projbin_out_.emplace_back(tuple<int, int, int>(ibin_out - start_out, outmem, nstubs_out));
0347           }
0348         }
0349         for (int ibin_in = start_in; ibin_in <= last_in; ibin_in++) {
0350           for (unsigned int inmem = 0; inmem < innervmstubs_.size(); inmem++) {
0351             int nstubs_in = innervmstubs_[inmem]->nVMStubsBinned(ibin_in);
0352             if (nstubs_in > 0)
0353               trpdata.projbin_in_.emplace_back(tuple<int, int, int>(ibin_in - start_in, inmem, nstubs_in));
0354           }
0355         }
0356 
0357         if (!trpdata.projbin_in_.empty() && !trpdata.projbin_out_.empty()) {
0358           goodtrpdata = true;
0359         }
0360       }
0361 
0362       istub++;
0363       if (istub >= middleallstubs_[midmem]->nStubs()) {
0364         istub = 0;
0365         midmem++;
0366       }
0367 
0368     } else if ((!trpbuffernearfull) && midmem < midmemend && istub == 0)
0369       midmem++;
0370 
0371     goodtrpdata___ = goodtrpdata__;
0372     goodtrpdata__ = goodtrpdata;
0373 
0374     trpdata___ = trpdata__;
0375     trpdata__ = trpdata;
0376 
0377     //
0378     // stop looping over istep if done
0379     //
0380 
0381     bool done = true;
0382 
0383     if (midmem < midmemend || (!trpdatabuffer.empty())) {
0384       done = false;
0385     }
0386 
0387     for (auto& trpunit : trpunits_) {
0388       if (!(trpunit.idle() && trpunit.empty()))
0389         done = false;
0390     }
0391 
0392     if (done) {
0393       donecount++;
0394     }
0395 
0396     //FIXME This should be done cleaner... Not too hard, but need to check fully the TEBuffer and TEUnit buffer.
0397     if (donecount > 4) {
0398       break;
0399     }
0400   }
0401 
0402   if (settings_.writeMonitorData("TPD")) {
0403     globals_->ofstream("trackletprocessordisplaced.txt")
0404         << getName() << " " << countall << " " << countsel << std::endl;
0405   }
0406 }