File indexing completed on 2024-04-06 12:22:01
0001 #include "L1Trigger/TrackFindingTracklet/interface/Sector.h"
0002 #include "L1Trigger/TrackFindingTracklet/interface/MemoryBase.h"
0003 #include "L1Trigger/TrackFindingTracklet/interface/Settings.h"
0004 #include "L1Trigger/TrackFindingTracklet/interface/Globals.h"
0005
0006 #include "L1Trigger/TrackFindingTracklet/interface/DTCLinkMemory.h"
0007 #include "L1Trigger/TrackFindingTracklet/interface/InputLinkMemory.h"
0008 #include "L1Trigger/TrackFindingTracklet/interface/AllStubsMemory.h"
0009 #include "L1Trigger/TrackFindingTracklet/interface/AllInnerStubsMemory.h"
0010 #include "L1Trigger/TrackFindingTracklet/interface/VMStubsTEMemory.h"
0011 #include "L1Trigger/TrackFindingTracklet/interface/VMStubsMEMemory.h"
0012 #include "L1Trigger/TrackFindingTracklet/interface/StubPairsMemory.h"
0013 #include "L1Trigger/TrackFindingTracklet/interface/StubTripletsMemory.h"
0014 #include "L1Trigger/TrackFindingTracklet/interface/TrackletParametersMemory.h"
0015 #include "L1Trigger/TrackFindingTracklet/interface/TrackletProjectionsMemory.h"
0016 #include "L1Trigger/TrackFindingTracklet/interface/AllProjectionsMemory.h"
0017 #include "L1Trigger/TrackFindingTracklet/interface/VMProjectionsMemory.h"
0018 #include "L1Trigger/TrackFindingTracklet/interface/CandidateMatchMemory.h"
0019 #include "L1Trigger/TrackFindingTracklet/interface/FullMatchMemory.h"
0020 #include "L1Trigger/TrackFindingTracklet/interface/TrackFitMemory.h"
0021 #include "L1Trigger/TrackFindingTracklet/interface/CleanTrackMemory.h"
0022
0023 #include "L1Trigger/TrackFindingTracklet/interface/InputRouter.h"
0024 #include "L1Trigger/TrackFindingTracklet/interface/VMRouterCM.h"
0025 #include "L1Trigger/TrackFindingTracklet/interface/VMRouter.h"
0026 #include "L1Trigger/TrackFindingTracklet/interface/TrackletEngine.h"
0027 #include "L1Trigger/TrackFindingTracklet/interface/TrackletEngineDisplaced.h"
0028 #include "L1Trigger/TrackFindingTracklet/interface/TripletEngine.h"
0029 #include "L1Trigger/TrackFindingTracklet/interface/TrackletCalculator.h"
0030 #include "L1Trigger/TrackFindingTracklet/interface/TrackletProcessor.h"
0031 #include "L1Trigger/TrackFindingTracklet/interface/TrackletProcessorDisplaced.h"
0032 #include "L1Trigger/TrackFindingTracklet/interface/TrackletCalculatorDisplaced.h"
0033 #include "L1Trigger/TrackFindingTracklet/interface/ProjectionRouter.h"
0034 #include "L1Trigger/TrackFindingTracklet/interface/MatchEngine.h"
0035 #include "L1Trigger/TrackFindingTracklet/interface/MatchCalculator.h"
0036 #include "L1Trigger/TrackFindingTracklet/interface/MatchProcessor.h"
0037 #include "L1Trigger/TrackFindingTracklet/interface/FitTrack.h"
0038 #include "L1Trigger/TrackFindingTracklet/interface/PurgeDuplicate.h"
0039 #include "L1Trigger/TrackFindingTracklet/interface/StubStreamData.h"
0040 #include "L1Trigger/TrackFindingTracklet/interface/Util.h"
0041
0042 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0043 #include "FWCore/Utilities/interface/Exception.h"
0044 #include "DataFormats/Math/interface/deltaPhi.h"
0045 #include "L1Trigger/TrackFindingTracklet/interface/TrackletLUT.h"
0046
0047 #include <deque>
0048
0049 using namespace std;
0050 using namespace trklet;
0051
0052 Sector::Sector(Settings const& settings, Globals* globals) : isector_(-1), settings_(settings), globals_(globals) {}
0053
0054 Sector::~Sector() = default;
0055
0056 void Sector::setSector(unsigned int isector) {
0057 assert(isector < N_SECTOR);
0058 isector_ = isector;
0059 double dphi = 2 * M_PI / N_SECTOR;
0060 double dphiHG = 0.5 * settings_.dphisectorHG() - M_PI / N_SECTOR;
0061 phimin_ = isector_ * dphi - dphiHG;
0062 phimax_ = phimin_ + dphi + 2 * dphiHG;
0063 phimin_ -= M_PI / N_SECTOR;
0064 phimax_ -= M_PI / N_SECTOR;
0065 phimin_ = reco::reduceRange(phimin_);
0066 phimax_ = reco::reduceRange(phimax_);
0067 if (phimin_ > phimax_) {
0068 phimin_ -= 2 * M_PI;
0069 }
0070 }
0071
0072 bool Sector::addStub(L1TStub stub, string dtc) {
0073 unsigned int layerdisk = stub.layerdisk();
0074 int nrbits = 3;
0075
0076 if (layerdisk < N_LAYER && globals_->phiCorr(layerdisk) == nullptr) {
0077 globals_->phiCorr(layerdisk) = new TrackletLUT(settings_);
0078 globals_->phiCorr(layerdisk)->initPhiCorrTable(layerdisk, nrbits);
0079 }
0080
0081 Stub fpgastub(stub, settings_, *globals_);
0082
0083 if (layerdisk < N_LAYER) {
0084 FPGAWord r = fpgastub.r();
0085 int bendbin = fpgastub.bend().value();
0086 int rbin = (r.value() + (1 << (r.nbits() - 1))) >> (r.nbits() - nrbits);
0087 const TrackletLUT& phiCorrTable = *globals_->phiCorr(layerdisk);
0088 int iphicorr = phiCorrTable.lookup(bendbin * (1 << nrbits) + rbin);
0089 fpgastub.setPhiCorr(iphicorr);
0090 }
0091
0092 int nadd = 0;
0093 for (unsigned int i = 0; i < DL_.size(); i++) {
0094 const string& name = DL_[i]->getName();
0095 if (name.find("_" + dtc) == string::npos)
0096 continue;
0097 DL_[i]->addStub(stub, fpgastub);
0098 nadd++;
0099 }
0100
0101 if (!(settings_.reduced()))
0102 assert(nadd == 1);
0103
0104 return true;
0105 }
0106
0107 void Sector::addMem(string memType, string memName) {
0108 if (memType == "DTCLink:") {
0109 addMemToVec(DL_, memName, settings_, phimin_, phimax_);
0110 } else if (memType == "InputLink:") {
0111 addMemToVec(IL_, memName, settings_, phimin_, phimax_);
0112 } else if (memType == "AllStubs:") {
0113 addMemToVec(AS_, memName, settings_);
0114 } else if (memType == "AllInnerStubs:") {
0115 addMemToVec(AIS_, memName, settings_);
0116 } else if (memType == "VMStubsTE:") {
0117 addMemToVec(VMSTE_, memName, settings_);
0118 } else if (memType == "VMStubsME:") {
0119 addMemToVec(VMSME_, memName, settings_);
0120 } else if (memType == "StubPairs:" || memType == "StubPairsDisplaced:") {
0121 addMemToVec(SP_, memName, settings_);
0122 } else if (memType == "StubTriplets:") {
0123 addMemToVec(ST_, memName, settings_);
0124 } else if (memType == "TrackletParameters:") {
0125 addMemToVec(TPAR_, memName, settings_);
0126 } else if (memType == "TrackletProjections:") {
0127 addMemToVec(TPROJ_, memName, settings_);
0128 } else if (memType == "AllProj:") {
0129 addMemToVec(AP_, memName, settings_);
0130 } else if (memType == "VMProjections:") {
0131 addMemToVec(VMPROJ_, memName, settings_);
0132 } else if (memType == "CandidateMatch:") {
0133 addMemToVec(CM_, memName, settings_);
0134 } else if (memType == "FullMatch:") {
0135 addMemToVec(FM_, memName, settings_);
0136 } else if (memType == "TrackFit:") {
0137 addMemToVec(TF_, memName, settings_, phimin_, phimax_);
0138 } else if (memType == "CleanTrack:") {
0139 addMemToVec(CT_, memName, settings_, phimin_, phimax_);
0140 } else {
0141 edm::LogPrint("Tracklet") << "Don't know of memory type: " << memType;
0142 exit(0);
0143 }
0144 }
0145
0146 void Sector::addProc(string procType, string procName) {
0147 if (procType == "InputRouter:") {
0148 addProcToVec(IR_, procName, settings_, globals_);
0149 } else if (procType == "VMRouter:") {
0150 addProcToVec(VMR_, procName, settings_, globals_);
0151 } else if (procType == "VMRouterCM:") {
0152 addProcToVec(VMRCM_, procName, settings_, globals_);
0153 } else if (procType == "TrackletEngine:") {
0154 addProcToVec(TE_, procName, settings_, globals_);
0155 } else if (procType == "TrackletEngineDisplaced:") {
0156 addProcToVec(TED_, procName, settings_, globals_);
0157 } else if (procType == "TripletEngine:") {
0158 addProcToVec(TRE_, procName, settings_, globals_);
0159 } else if (procType == "TrackletCalculator:") {
0160 addProcToVec(TC_, procName, settings_, globals_);
0161 } else if (procType == "TrackletProcessor:") {
0162 addProcToVec(TP_, procName, settings_, globals_);
0163 } else if (procType == "TrackletProcessorDisplaced:") {
0164 addProcToVec(TPD_, procName, settings_, globals_);
0165 } else if (procType == "TrackletCalculatorDisplaced:") {
0166 addProcToVec(TCD_, procName, settings_, globals_);
0167 } else if (procType == "ProjectionRouter:") {
0168 addProcToVec(PR_, procName, settings_, globals_);
0169 } else if (procType == "MatchEngine:") {
0170 addProcToVec(ME_, procName, settings_, globals_);
0171 } else if (procType == "MatchCalculator:" ||
0172 procType == "DiskMatchCalculator:") {
0173 addProcToVec(MC_, procName, settings_, globals_);
0174 } else if (procType == "MatchProcessor:") {
0175 addProcToVec(MP_, procName, settings_, globals_);
0176 } else if (procType == "FitTrack:") {
0177 addProcToVec(FT_, procName, settings_, globals_);
0178 } else if (procType == "PurgeDuplicate:") {
0179 addProcToVec(PD_, procName, settings_, globals_);
0180 } else {
0181 edm::LogPrint("Tracklet") << "Don't know of processing type: " << procType;
0182 exit(0);
0183 }
0184 }
0185
0186 void Sector::addWire(string mem, string procinfull, string procoutfull) {
0187 stringstream ss1(procinfull);
0188 string procin, output;
0189 getline(ss1, procin, '.');
0190 getline(ss1, output);
0191
0192 stringstream ss2(procoutfull);
0193 string procout, input;
0194 getline(ss2, procout, '.');
0195 getline(ss2, input);
0196
0197 MemoryBase* memory = getMem(mem);
0198
0199 if (!procin.empty()) {
0200 ProcessBase* inProc = getProc(procin);
0201 inProc->addOutput(memory, output);
0202 }
0203
0204 if (!procout.empty()) {
0205 ProcessBase* outProc = getProc(procout);
0206 outProc->addInput(memory, input);
0207 }
0208 }
0209
0210 ProcessBase* Sector::getProc(string procName) {
0211 auto it = Processes_.find(procName);
0212
0213 if (it != Processes_.end()) {
0214 return it->second;
0215 }
0216 throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " Could not find process : " << procName << endl;
0217 return nullptr;
0218 }
0219
0220 MemoryBase* Sector::getMem(string memName) {
0221 auto it = Memories_.find(memName);
0222
0223 if (it != Memories_.end()) {
0224 return it->second;
0225 }
0226 throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " Could not find memory : " << memName;
0227 return nullptr;
0228 }
0229
0230 void Sector::writeDTCStubs(bool first) {
0231 for (auto& i : DL_) {
0232 i->writeStubs(first, isector_);
0233 }
0234 }
0235
0236 void Sector::writeIRStubs(bool first) {
0237 for (auto& i : IL_) {
0238 i->writeStubs(first, isector_);
0239 }
0240 }
0241
0242 void Sector::writeVMSTE(bool first) {
0243 for (auto& i : VMSTE_) {
0244 i->writeStubs(first, isector_);
0245 }
0246 }
0247
0248 void Sector::writeVMSME(bool first) {
0249 for (auto& i : VMSME_) {
0250 i->writeStubs(first, isector_);
0251 }
0252 }
0253
0254 void Sector::writeAS(bool first) {
0255 for (auto& i : AS_) {
0256 i->writeStubs(first, isector_);
0257 }
0258 }
0259
0260 void Sector::writeAIS(bool first) {
0261 for (auto& i : AIS_) {
0262 i->writeStubs(first, isector_);
0263 }
0264 }
0265
0266 void Sector::writeSP(bool first) {
0267 for (auto& i : SP_) {
0268 i->writeSP(first, isector_);
0269 }
0270 }
0271
0272 void Sector::writeST(bool first) {
0273 for (auto& i : ST_) {
0274 i->writeST(first, isector_);
0275 }
0276 }
0277
0278 void Sector::writeTPAR(bool first) {
0279 for (auto& i : TPAR_) {
0280 i->writeTPAR(first, isector_);
0281 }
0282 }
0283
0284 void Sector::writeTPROJ(bool first) {
0285 for (auto& i : TPROJ_) {
0286 i->writeTPROJ(first, isector_);
0287 }
0288 }
0289
0290 void Sector::writeAP(bool first) {
0291 for (auto& i : AP_) {
0292 i->writeAP(first, isector_);
0293 }
0294 }
0295
0296 void Sector::writeVMPROJ(bool first) {
0297 for (auto& i : VMPROJ_) {
0298 i->writeVMPROJ(first, isector_);
0299 }
0300 }
0301
0302 void Sector::writeCM(bool first) {
0303 for (auto& i : CM_) {
0304 i->writeCM(first, isector_);
0305 }
0306 }
0307
0308 void Sector::writeMC(bool first) {
0309 for (auto& i : FM_) {
0310 i->writeMC(first, isector_);
0311 }
0312 }
0313
0314 void Sector::writeTF(bool first) {
0315 for (auto& i : TF_) {
0316 i->writeTF(first, isector_);
0317 }
0318 }
0319
0320 void Sector::writeCT(bool first) {
0321 for (auto& i : CT_) {
0322 i->writeCT(first, isector_);
0323 }
0324 }
0325
0326 void Sector::clean() {
0327 for (auto& mem : MemoriesV_) {
0328 mem->clean();
0329 }
0330 }
0331
0332 void Sector::executeIR() {
0333 for (auto& i : IR_) {
0334 i->execute();
0335 }
0336 }
0337
0338 void Sector::executeVMR() {
0339 if (settings_.writeMonitorData("IL")) {
0340 ofstream& out = globals_->ofstream("inputlink.txt");
0341 for (auto& i : IL_) {
0342 out << i->getName() << " " << i->nStubs() << endl;
0343 }
0344 }
0345 for (auto& i : VMR_) {
0346 i->execute();
0347 }
0348 for (auto& i : VMRCM_) {
0349 i->execute(isector_);
0350 }
0351 }
0352
0353 void Sector::executeTE() {
0354 for (auto& i : TE_) {
0355 i->execute();
0356 }
0357 }
0358
0359 void Sector::executeTED() {
0360 for (auto& i : TED_) {
0361 i->execute();
0362 }
0363 }
0364
0365 void Sector::executeTRE() {
0366 for (auto& i : TRE_) {
0367 i->execute();
0368 }
0369 }
0370
0371 void Sector::executeTP() {
0372 for (auto& i : TP_) {
0373 i->execute(isector_, phimin_, phimax_);
0374 }
0375 }
0376
0377 void Sector::executeTPD() {
0378 for (auto& i : TPD_) {
0379 i->execute(isector_, phimin_, phimax_);
0380 }
0381 }
0382
0383 void Sector::executeTC() {
0384 for (auto& i : TC_) {
0385 i->execute(isector_, phimin_, phimax_);
0386 }
0387
0388 if (settings_.writeMonitorData("TrackProjOcc")) {
0389 ofstream& out = globals_->ofstream("trackprojocc.txt");
0390 for (auto& i : TPROJ_) {
0391 out << i->getName() << " " << i->nTracklets() << endl;
0392 }
0393 }
0394 }
0395
0396 void Sector::executeTCD() {
0397 for (auto& i : TCD_) {
0398 i->execute(isector_, phimin_, phimax_);
0399 }
0400 }
0401
0402 void Sector::executePR() {
0403 for (auto& i : PR_) {
0404 i->execute();
0405 }
0406 }
0407
0408 void Sector::executeME() {
0409 for (auto& i : ME_) {
0410 i->execute(isector_);
0411 }
0412 }
0413
0414 void Sector::executeMC() {
0415 for (auto& i : MC_) {
0416 i->execute(isector_, phimin_);
0417 }
0418 }
0419
0420 void Sector::executeMP() {
0421 for (auto& i : MP_) {
0422 i->execute(isector_, phimin_);
0423 }
0424 }
0425
0426
0427
0428
0429
0430 void Sector::executeFT(vector<vector<string>>& streamsTrackRaw, vector<vector<StubStreamData>>& streamsStubRaw) {
0431 const int numChannels = streamsTrackRaw.size() / N_SECTOR;
0432 const int maxNumProjectionLayers = streamsStubRaw.size() / streamsTrackRaw.size();
0433 const int offsetTrack = isector_ * numChannels;
0434 int channelTrack(0);
0435
0436 for (auto& i : FT_) {
0437
0438 deque<string> streamTrackTmp;
0439 vector<deque<StubStreamData>> streamsStubTmp(maxNumProjectionLayers);
0440 i->execute(streamTrackTmp, streamsStubTmp, isector_);
0441 if (!settings_.storeTrackBuilderOutput())
0442 continue;
0443 const int offsetStub = (offsetTrack + channelTrack) * maxNumProjectionLayers;
0444 streamsTrackRaw[offsetTrack + channelTrack] = vector<string>(streamTrackTmp.begin(), streamTrackTmp.end());
0445 channelTrack++;
0446 int channelStub(0);
0447 for (auto& stream : streamsStubTmp)
0448 streamsStubRaw[offsetStub + channelStub++] = vector<StubStreamData>(stream.begin(), stream.end());
0449 }
0450 }
0451
0452
0453 void Sector::executePD(std::vector<Track>& tracks) {
0454 for (auto& i : PD_) {
0455 i->execute(tracks, isector_);
0456 }
0457 }
0458
0459 std::vector<Tracklet*> Sector::getAllTracklets() const {
0460 std::vector<Tracklet*> tmp;
0461 for (auto& tpar : TPAR_) {
0462 for (unsigned int j = 0; j < tpar->nTracklets(); j++) {
0463 tmp.push_back(tpar->getTracklet(j));
0464 }
0465 }
0466 return tmp;
0467 }
0468
0469 std::vector<const Stub*> Sector::getStubs() const {
0470 std::vector<const Stub*> tmp;
0471
0472 for (auto& imem : IL_) {
0473 for (unsigned int istub = 0; istub < imem->nStubs(); istub++) {
0474 tmp.push_back(imem->getStub(istub));
0475 }
0476 }
0477
0478 return tmp;
0479 }
0480
0481 std::unordered_set<int> Sector::seedMatch(int itp) const {
0482 std::unordered_set<int> tmpSeeds;
0483 for (auto& i : TPAR_) {
0484 unsigned int nTracklet = i->nTracklets();
0485 for (unsigned int j = 0; j < nTracklet; j++) {
0486 if (i->getTracklet(j)->tpseed() == itp) {
0487 tmpSeeds.insert(i->getTracklet(j)->getISeed());
0488 }
0489 }
0490 }
0491 return tmpSeeds;
0492 }