File indexing completed on 2024-04-06 12:22:05
0001 #include "L1Trigger/TrackFindingTracklet/interface/TripletEngine.h"
0002 #include "L1Trigger/TrackFindingTracklet/interface/Settings.h"
0003 #include "L1Trigger/TrackFindingTracklet/interface/Globals.h"
0004
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include "FWCore/Utilities/interface/Exception.h"
0007
0008 #include <algorithm>
0009
0010 using namespace std;
0011 using namespace trklet;
0012
0013 TripletEngine::TripletEngine(string name, Settings const &settings, Globals *global)
0014 : ProcessBase(name, settings, global) {
0015 stubpairs_.clear();
0016 thirdvmstubs_.clear();
0017 layer1_ = 0;
0018 layer2_ = 0;
0019 layer3_ = 0;
0020 disk1_ = 0;
0021 disk2_ = 0;
0022 disk3_ = 0;
0023 dct1_ = 0;
0024 dct2_ = 0;
0025 dct3_ = 0;
0026 phi1_ = 0;
0027 phi2_ = 0;
0028 phi3_ = 0;
0029 z1_ = 0;
0030 z2_ = 0;
0031 z3_ = 0;
0032 r1_ = 0;
0033 r2_ = 0;
0034 r3_ = 0;
0035
0036 if (name_[4] == 'L')
0037 layer1_ = name_[5] - '0';
0038 if (name_[4] == 'D')
0039 disk1_ = name_[5] - '0';
0040 if (name_[7] == 'L')
0041 layer2_ = name_[8] - '0';
0042 if (name_[7] == 'D')
0043 disk2_ = name_[8] - '0';
0044
0045 if (layer1_ == 3 && layer2_ == 4) {
0046 layer3_ = 2;
0047 iSeed_ = 8;
0048 } else if (layer1_ == 5 && layer2_ == 6) {
0049 layer3_ = 4;
0050 iSeed_ = 9;
0051 } else if (layer1_ == 2 && layer2_ == 3) {
0052 disk3_ = 1;
0053 iSeed_ = 10;
0054 } else if (disk1_ == 1 && disk2_ == 2) {
0055 layer3_ = 2;
0056 iSeed_ = 11;
0057 } else
0058 throw cms::Exception("LogicError") << __FILE__ << " " << __LINE__ << " Invalid seeding!";
0059
0060 if ((layer2_ == 4 && layer3_ == 2) || (layer2_ == 6 && layer3_ == 4)) {
0061 secondphibits_ = settings_.nfinephi(1, iSeed_);
0062 thirdphibits_ = settings_.nfinephi(2, iSeed_);
0063 }
0064 if ((layer2_ == 3 && disk3_ == 1) || (disk2_ == 2 && layer3_ == 2)) {
0065 secondphibits_ = settings_.nfinephi(1, iSeed_);
0066 thirdphibits_ = settings_.nfinephi(2, iSeed_);
0067 }
0068 if (settings_.enableTripletTables() && !settings_.writeTripletTables())
0069 readTables();
0070 }
0071
0072 TripletEngine::~TripletEngine() {
0073 if (settings_.writeTripletTables())
0074 writeTables();
0075 }
0076
0077 void TripletEngine::addOutput(MemoryBase *memory, string output) {
0078 if (settings_.writetrace()) {
0079 edm::LogVerbatim("Tracklet") << "In " << name_ << " adding output to " << memory->getName() << " to output "
0080 << output;
0081 }
0082 if (output == "stubtripout") {
0083 auto *tmp = dynamic_cast<StubTripletsMemory *>(memory);
0084 assert(tmp != nullptr);
0085 stubtriplets_ = tmp;
0086 return;
0087 }
0088 throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " Could not find output : " << output;
0089 }
0090
0091 void TripletEngine::addInput(MemoryBase *memory, string input) {
0092 if (settings_.writetrace()) {
0093 edm::LogVerbatim("Tracklet") << "In " << name_ << " adding input from " << memory->getName() << " to input "
0094 << input;
0095 }
0096 if (input == "thirdvmstubin") {
0097 auto *tmp = dynamic_cast<VMStubsTEMemory *>(memory);
0098 assert(tmp != nullptr);
0099 thirdvmstubs_.push_back(tmp);
0100 return;
0101 }
0102 if (input.substr(0, 8) == "stubpair") {
0103 auto *tmp = dynamic_cast<StubPairsMemory *>(memory);
0104 assert(tmp != nullptr);
0105 stubpairs_.push_back(tmp);
0106 return;
0107 }
0108 throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " Could not find input : " << input;
0109 }
0110
0111 void TripletEngine::execute() {
0112 unsigned int countall = 0;
0113 unsigned int countpass = 0;
0114 unsigned int nThirdStubs = 0;
0115 count_ = 0;
0116
0117 for (unsigned int iThirdMem = 0; iThirdMem < thirdvmstubs_.size();
0118 nThirdStubs += thirdvmstubs_.at(iThirdMem)->nVMStubs(), iThirdMem++)
0119 ;
0120
0121 assert(!thirdvmstubs_.empty());
0122 assert(!stubpairs_.empty());
0123
0124 bool print = false && (getName().substr(0, 10) == "TRE_L2cL3c");
0125
0126 print = print && nThirdStubs > 0;
0127
0128 if (print) {
0129 edm::LogVerbatim("Tracklet") << "In TripletEngine::execute : " << getName() << " " << nThirdStubs << ":";
0130 for (unsigned int i = 0; i < thirdvmstubs_.size(); ++i) {
0131 edm::LogVerbatim("Tracklet") << thirdvmstubs_.at(i)->getName() << " " << thirdvmstubs_.at(i)->nVMStubs();
0132 }
0133 std::string oss = "";
0134 for (unsigned int i = 0; i < stubpairs_.size(); ++i) {
0135 oss += std::to_string(stubpairs_.at(i)->nStubPairs());
0136 oss += " ";
0137 }
0138 edm::LogVerbatim("Tracklet") << oss;
0139 for (unsigned int i = 0; i < stubpairs_.size(); ++i) {
0140 edm::LogVerbatim("Tracklet") << " " << stubpairs_.at(i)->getName();
0141 }
0142 }
0143
0144 tmpSPTable_.clear();
0145
0146 for (unsigned int i = 0; i < stubpairs_.size(); ++i) {
0147 for (unsigned int j = 0; j < stubpairs_.at(i)->nStubPairs(); ++j) {
0148 if (print)
0149 edm::LogVerbatim("Tracklet") << " ***** " << stubpairs_.at(i)->getName() << " "
0150 << stubpairs_.at(i)->nStubPairs();
0151
0152 auto firstvmstub = stubpairs_.at(i)->getVMStub1(j);
0153 auto secondvmstub = stubpairs_.at(i)->getVMStub2(j);
0154
0155 if ((layer2_ == 4 && layer3_ == 2) || (layer2_ == 6 && layer3_ == 4)) {
0156 constexpr unsigned int vmbitshift = 10;
0157 int lookupbits = (int)((firstvmstub.vmbits().value() >> vmbitshift) & 1023);
0158 int newbin = (lookupbits & 127);
0159 int bin = newbin / 8;
0160
0161 int start = (bin >> 1);
0162 int last = start + (bin & 1);
0163
0164 for (int ibin = start; ibin <= last; ibin++) {
0165 for (unsigned int k = 0; k < thirdvmstubs_.size(); k++) {
0166 string vmsteSuffix = thirdvmstubs_.at(k)->getLastPartOfName();
0167 vmsteSuffix = vmsteSuffix.substr(0, vmsteSuffix.find_last_of('n'));
0168 if (stubpairs_.at(i)->getLastPartOfName() != vmsteSuffix)
0169 continue;
0170 for (unsigned int l = 0; l < thirdvmstubs_.at(k)->nVMStubsBinned(ibin); l++) {
0171 if (settings_.debugTracklet()) {
0172 edm::LogVerbatim("Tracklet") << "In " << getName() << " have third stub";
0173 }
0174
0175 if (countall >= settings_.maxStep("TRE"))
0176 break;
0177 countall++;
0178
0179 const VMStubTE &thirdvmstub = thirdvmstubs_.at(k)->getVMStubTEBinned(ibin, l);
0180
0181 assert(secondphibits_ != -1);
0182 assert(thirdphibits_ != -1);
0183
0184 unsigned int nvmsecond = settings_.nallstubs(layer2_ - 1) * settings_.nvmte(1, iSeed_);
0185 unsigned int nvmbitssecond = nbits(nvmsecond);
0186
0187 FPGAWord iphisecondbin = secondvmstub.stub()->iphivmFineBins(nvmbitssecond, secondphibits_);
0188
0189
0190
0191 FPGAWord iphithirdbin = thirdvmstub.finephi();
0192
0193 unsigned int index = (iphisecondbin.value() << thirdphibits_) + iphithirdbin.value();
0194
0195 FPGAWord secondbend = secondvmstub.bend();
0196 FPGAWord thirdbend = thirdvmstub.bend();
0197
0198 index = (index << secondbend.nbits()) + secondbend.value();
0199 index = (index << thirdbend.nbits()) + thirdbend.value();
0200
0201 if ((settings_.enableTripletTables() && !settings_.writeTripletTables()) &&
0202 (index >= table_.size() || !table_[index])) {
0203 if (settings_.debugTracklet()) {
0204 edm::LogVerbatim("Tracklet")
0205 << "Stub pair rejected because of stub pt cut bends : "
0206 << settings_.benddecode(secondvmstub.bend().value(), layer2_ - 1, secondvmstub.isPSmodule())
0207 << " " << settings_.benddecode(thirdvmstub.bend().value(), layer3_ - 1, thirdvmstub.isPSmodule());
0208 }
0209
0210
0211
0212
0213 }
0214 if (settings_.writeTripletTables()) {
0215 if (index >= table_.size())
0216 table_.resize(index + 1, false);
0217 table_[index] = true;
0218
0219 const unsigned spIndex = stubpairs_.at(i)->getIndex(j);
0220 const string &tedName = stubpairs_.at(i)->getTEDName(j);
0221 if (!tmpSPTable_.count(tedName))
0222 tmpSPTable_[tedName];
0223 if (spIndex >= tmpSPTable_.at(tedName).size())
0224 tmpSPTable_.at(tedName).resize(spIndex + 1);
0225 tmpSPTable_.at(tedName).at(spIndex).push_back(stubpairs_.at(i)->getName());
0226 }
0227
0228 if (settings_.debugTracklet())
0229 edm::LogVerbatim("Tracklet") << "Adding layer-layer pair in " << getName();
0230 if (settings_.writeMonitorData("Seeds")) {
0231 ofstream fout("seeds.txt", ofstream::app);
0232 fout << __FILE__ << ":" << __LINE__ << " " << name_ << " " << iSeed_ << endl;
0233 fout.close();
0234 }
0235 stubtriplets_->addStubs(thirdvmstub.stub(),
0236 (stubpairs_.at(i))->getVMStub1(j).stub(),
0237 (stubpairs_.at(i))->getVMStub2(j).stub());
0238
0239 countpass++;
0240 }
0241 }
0242 }
0243
0244 }
0245
0246 else if (disk2_ == 2 && layer3_ == 2) {
0247 int lookupbits = (int)((firstvmstub.vmbits().value() >> 10) & 1023);
0248 int newbin = (lookupbits & 127);
0249 int bin = newbin / 8;
0250
0251 int start = (bin >> 1);
0252 int last = start + (bin & 1);
0253
0254 if (firstvmstub.stub()->disk().value() < 0) {
0255 start = settings_.NLONGVMBINS() - last - 1;
0256 last = settings_.NLONGVMBINS() - start - 1;
0257 }
0258
0259 for (int ibin = start; ibin <= last; ibin++) {
0260 for (unsigned int k = 0; k < thirdvmstubs_.size(); k++) {
0261 string vmsteSuffix = thirdvmstubs_.at(k)->getLastPartOfName();
0262 vmsteSuffix = vmsteSuffix.substr(0, vmsteSuffix.find_last_of('n'));
0263 if (stubpairs_.at(i)->getLastPartOfName() != vmsteSuffix)
0264 continue;
0265
0266 for (unsigned int l = 0; l < thirdvmstubs_.at(k)->nVMStubsBinned(ibin); l++) {
0267 if (countall >= settings_.maxStep("TRE"))
0268 break;
0269 countall++;
0270
0271 const VMStubTE &thirdvmstub = thirdvmstubs_.at(k)->getVMStubTEBinned(ibin, l);
0272
0273 assert(secondphibits_ != -1);
0274 assert(thirdphibits_ != -1);
0275
0276 FPGAWord iphisecondbin = secondvmstub.finephi();
0277 FPGAWord iphithirdbin = thirdvmstub.finephi();
0278
0279 unsigned int index = (iphisecondbin.value() << thirdphibits_) + iphithirdbin.value();
0280
0281 FPGAWord secondbend = secondvmstub.bend();
0282 FPGAWord thirdbend = thirdvmstub.bend();
0283
0284 index = (index << secondbend.nbits()) + secondbend.value();
0285 index = (index << thirdbend.nbits()) + thirdbend.value();
0286
0287 if ((settings_.enableTripletTables() && !settings_.writeTripletTables()) &&
0288 (index >= table_.size() || !table_[index])) {
0289 if (settings_.debugTracklet()) {
0290 edm::LogVerbatim("Tracklet")
0291 << "Stub triplet rejected because of stub pt cut bends : "
0292 << settings_.benddecode(secondvmstub.bend().value(), disk2_ + 5, secondvmstub.isPSmodule()) << " "
0293 << settings_.benddecode(thirdvmstub.bend().value(), layer3_ - 1, thirdvmstub.isPSmodule());
0294 }
0295 continue;
0296 }
0297 if (settings_.writeTripletTables()) {
0298 if (index >= table_.size())
0299 table_.resize(index + 1, false);
0300 table_[index] = true;
0301
0302 const unsigned spIndex = stubpairs_.at(i)->getIndex(j);
0303 const string &tedName = stubpairs_.at(i)->getTEDName(j);
0304 if (!tmpSPTable_.count(tedName))
0305 tmpSPTable_[tedName];
0306 if (spIndex >= tmpSPTable_.at(tedName).size())
0307 tmpSPTable_.at(tedName).resize(spIndex + 1);
0308 tmpSPTable_.at(tedName).at(spIndex).push_back(stubpairs_.at(i)->getName());
0309 }
0310
0311 if (settings_.debugTracklet())
0312 edm::LogVerbatim("Tracklet") << "Adding layer-disk pair in " << getName();
0313 if (settings_.writeMonitorData("Seeds")) {
0314 ofstream fout("seeds.txt", ofstream::app);
0315 fout << __FILE__ << ":" << __LINE__ << " " << name_ << " " << iSeed_ << endl;
0316 fout.close();
0317 }
0318 stubtriplets_->addStubs(thirdvmstub.stub(),
0319 (stubpairs_.at(i))->getVMStub1(j).stub(),
0320 (stubpairs_.at(i))->getVMStub2(j).stub());
0321 countpass++;
0322 }
0323 }
0324 }
0325 }
0326
0327 else if (layer2_ == 3 && disk3_ == 1) {
0328 int lookupbits = (int)((firstvmstub.vmbits().value() >> 10) & 1023);
0329
0330 int newbin = (lookupbits & 127);
0331 int bin = newbin / 8;
0332
0333 int start = (bin >> 1);
0334 int last = start + (bin & 1);
0335
0336 for (int ibin = start; ibin <= last; ibin++) {
0337 for (unsigned int k = 0; k < thirdvmstubs_.size(); k++) {
0338 string vmsteSuffix = thirdvmstubs_.at(k)->getLastPartOfName();
0339 vmsteSuffix = vmsteSuffix.substr(0, vmsteSuffix.find_last_of('n'));
0340 if (stubpairs_.at(i)->getLastPartOfName() != vmsteSuffix)
0341 continue;
0342 for (unsigned int l = 0; l < thirdvmstubs_.at(k)->nVMStubsBinned(ibin); l++) {
0343 if (countall >= settings_.maxStep("TRE"))
0344 break;
0345 countall++;
0346
0347 const VMStubTE &thirdvmstub = thirdvmstubs_.at(k)->getVMStubTEBinned(ibin, l);
0348
0349 assert(secondphibits_ != -1);
0350 assert(thirdphibits_ != -1);
0351
0352 unsigned int nvmsecond;
0353
0354 nvmsecond = settings_.nallstubs(layer2_ - 1) * settings_.nvmte(1, iSeed_);
0355 unsigned int nvmbitssecond = nbits(nvmsecond);
0356
0357 FPGAWord iphisecondbin = secondvmstub.stub()->iphivmFineBins(nvmbitssecond, secondphibits_);
0358
0359
0360
0361 FPGAWord iphithirdbin = thirdvmstub.finephi();
0362
0363 unsigned int index = (iphisecondbin.value() << thirdphibits_) + iphithirdbin.value();
0364
0365 FPGAWord secondbend = secondvmstub.bend();
0366 FPGAWord thirdbend = thirdvmstub.bend();
0367
0368 index = (index << secondbend.nbits()) + secondbend.value();
0369 index = (index << thirdbend.nbits()) + thirdbend.value();
0370
0371 if ((settings_.enableTripletTables() && !settings_.writeTripletTables()) &&
0372 (index >= table_.size() || !table_[index])) {
0373 if (settings_.debugTracklet()) {
0374 edm::LogVerbatim("Tracklet")
0375 << "Stub pair rejected because of stub pt cut bends : "
0376 << settings_.benddecode(secondvmstub.bend().value(), layer2_ - 1, secondvmstub.isPSmodule())
0377 << " " << settings_.benddecode(thirdvmstub.bend().value(), disk3_ + 5, thirdvmstub.isPSmodule());
0378 }
0379 continue;
0380 }
0381 if (settings_.writeTripletTables()) {
0382 if (index >= table_.size())
0383 table_.resize(index + 1, false);
0384 table_[index] = true;
0385
0386 const unsigned spIndex = stubpairs_.at(i)->getIndex(j);
0387 const string &tedName = stubpairs_.at(i)->getTEDName(j);
0388 if (!tmpSPTable_.count(tedName))
0389 tmpSPTable_[tedName];
0390 if (spIndex >= tmpSPTable_.at(tedName).size())
0391 tmpSPTable_.at(tedName).resize(spIndex + 1);
0392 tmpSPTable_.at(tedName).at(spIndex).push_back(stubpairs_.at(i)->getName());
0393 }
0394
0395 if (settings_.debugTracklet())
0396 edm::LogVerbatim("Tracklet") << "Adding layer-disk pair in " << getName();
0397 if (settings_.writeMonitorData("Seeds")) {
0398 ofstream fout("seeds.txt", ofstream::app);
0399 fout << __FILE__ << ":" << __LINE__ << " " << name_ << " " << iSeed_ << endl;
0400 fout.close();
0401 }
0402 stubtriplets_->addStubs(thirdvmstub.stub(),
0403 (stubpairs_.at(i))->getVMStub1(j).stub(),
0404 (stubpairs_.at(i))->getVMStub2(j).stub());
0405 countpass++;
0406 }
0407 }
0408 }
0409 }
0410 }
0411 }
0412
0413 for (const auto &tedName : tmpSPTable_) {
0414 for (unsigned spIndex = 0; spIndex < tedName.second.size(); spIndex++) {
0415 if (tedName.second.at(spIndex).empty())
0416 continue;
0417 vector<string> entry(tedName.second.at(spIndex));
0418 sort(entry.begin(), entry.end());
0419 entry.erase(unique(entry.begin(), entry.end()), entry.end());
0420 const string &spName = entry.at(0);
0421
0422 if (!spTable_.count(tedName.first))
0423 spTable_[tedName.first];
0424 if (spIndex >= spTable_.at(tedName.first).size())
0425 spTable_.at(tedName.first).resize(spIndex + 1);
0426 if (!spTable_.at(tedName.first).at(spIndex).count(spName))
0427 spTable_.at(tedName.first).at(spIndex)[spName] = 0;
0428 spTable_.at(tedName.first).at(spIndex)[spName]++;
0429 }
0430 }
0431
0432 if (settings_.writeMonitorData("TRE")) {
0433 globals_->ofstream("tripletengine.txt") << getName() << " " << countall << " " << countpass << endl;
0434 }
0435 }
0436
0437 void TripletEngine::readTables() {
0438 ifstream fin;
0439 string tableName, word;
0440 unsigned num;
0441
0442 string tablePath = settings_.tableTREFile();
0443 unsigned int finddir = tablePath.find("table_TRE_");
0444 tableName = tablePath.substr(0, finddir) + "table_" + name_ + ".txt";
0445
0446 fin.open(tableName, ifstream::in);
0447 if (!fin) {
0448 throw cms::Exception("BadConfig") << "TripletEngine::readTables, file " << tableName << " not known";
0449 }
0450 while (!fin.eof()) {
0451 fin >> word;
0452 num = atoi(word.c_str());
0453 table_.push_back(num > 0 ? true : false);
0454 }
0455 fin.close();
0456 }
0457
0458 void TripletEngine::writeTables() {
0459 ofstream fout;
0460 stringstream tableName;
0461
0462 tableName << "table/table_" << name_ << ".txt";
0463
0464 fout.open(tableName.str(), ofstream::out);
0465 for (const auto entry : table_)
0466 fout << entry << endl;
0467 fout.close();
0468
0469 for (const auto &tedName : spTable_) {
0470 tableName.str("");
0471 tableName << "table/table_" << tedName.first << "_" << name_ << ".txt";
0472
0473 fout.open(tableName.str(), ofstream::out);
0474 for (const auto &entry : tedName.second) {
0475 for (const auto &spName : entry)
0476 fout << spName.first << ":" << spName.second << " ";
0477 fout << endl;
0478 }
0479 fout.close();
0480 }
0481 }