CMS 3D CMS Logo

TripletEngine.cc
Go to the documentation of this file.
4 
7 
8 #include <algorithm>
9 
10 using namespace std;
11 using namespace trklet;
12 
13 TripletEngine::TripletEngine(string name, Settings const &settings, Globals *global)
14  : ProcessBase(name, settings, global) {
15  stubpairs_.clear();
16  thirdvmstubs_.clear();
17  layer1_ = 0;
18  layer2_ = 0;
19  layer3_ = 0;
20  disk1_ = 0;
21  disk2_ = 0;
22  disk3_ = 0;
23  dct1_ = 0;
24  dct2_ = 0;
25  dct3_ = 0;
26  phi1_ = 0;
27  phi2_ = 0;
28  phi3_ = 0;
29  z1_ = 0;
30  z2_ = 0;
31  z3_ = 0;
32  r1_ = 0;
33  r2_ = 0;
34  r3_ = 0;
35 
36  if (name_[4] == 'L')
37  layer1_ = name_[5] - '0';
38  if (name_[4] == 'D')
39  disk1_ = name_[5] - '0';
40  if (name_[7] == 'L')
41  layer2_ = name_[8] - '0';
42  if (name_[7] == 'D')
43  disk2_ = name_[8] - '0';
44 
45  if (layer1_ == 3 && layer2_ == 4) {
46  layer3_ = 2;
47  iSeed_ = 8;
48  } else if (layer1_ == 5 && layer2_ == 6) {
49  layer3_ = 4;
50  iSeed_ = 9;
51  } else if (layer1_ == 2 && layer2_ == 3) {
52  disk3_ = 1;
53  iSeed_ = 10;
54  } else if (disk1_ == 1 && disk2_ == 2) {
55  layer3_ = 2;
56  iSeed_ = 11;
57  } else
58  throw cms::Exception("LogicError") << __FILE__ << " " << __LINE__ << " Invalid seeding!";
59 
60  if ((layer2_ == 4 && layer3_ == 2) || (layer2_ == 6 && layer3_ == 4)) {
63  }
64  if ((layer2_ == 3 && disk3_ == 1) || (disk2_ == 2 && layer3_ == 2)) {
67  }
69  readTables();
70 }
71 
74  writeTables();
75 }
76 
78  if (settings_.writetrace()) {
79  edm::LogVerbatim("Tracklet") << "In " << name_ << " adding output to " << memory->getName() << " to output "
80  << output;
81  }
82  if (output == "stubtripout") {
83  auto *tmp = dynamic_cast<StubTripletsMemory *>(memory);
84  assert(tmp != nullptr);
86  return;
87  }
88  throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " Could not find output : " << output;
89 }
90 
92  if (settings_.writetrace()) {
93  edm::LogVerbatim("Tracklet") << "In " << name_ << " adding input from " << memory->getName() << " to input "
94  << input;
95  }
96  if (input == "thirdvmstubin") {
97  auto *tmp = dynamic_cast<VMStubsTEMemory *>(memory);
98  assert(tmp != nullptr);
99  thirdvmstubs_.push_back(tmp);
100  return;
101  }
102  if (input.substr(0, 8) == "stubpair") {
103  auto *tmp = dynamic_cast<StubPairsMemory *>(memory);
104  assert(tmp != nullptr);
105  stubpairs_.push_back(tmp);
106  return;
107  }
108  throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " Could not find input : " << input;
109 }
110 
112  unsigned int countall = 0;
113  unsigned int countpass = 0;
114  unsigned int nThirdStubs = 0;
115  count_ = 0;
116 
117  for (unsigned int iThirdMem = 0; iThirdMem < thirdvmstubs_.size();
118  nThirdStubs += thirdvmstubs_.at(iThirdMem)->nVMStubs(), iThirdMem++)
119  ;
120 
121  assert(!thirdvmstubs_.empty());
122  assert(!stubpairs_.empty());
123 
124  bool print = false && (getName().substr(0, 10) == "TRE_L2cL3c");
125 
126  print = print && nThirdStubs > 0;
127 
128  if (print) {
129  edm::LogVerbatim("Tracklet") << "In TripletEngine::execute : " << getName() << " " << nThirdStubs << ":";
130  for (unsigned int i = 0; i < thirdvmstubs_.size(); ++i) {
131  edm::LogVerbatim("Tracklet") << thirdvmstubs_.at(i)->getName() << " " << thirdvmstubs_.at(i)->nVMStubs();
132  }
133  std::string oss = "";
134  for (unsigned int i = 0; i < stubpairs_.size(); ++i) {
135  oss += std::to_string(stubpairs_.at(i)->nStubPairs());
136  oss += " ";
137  }
138  edm::LogVerbatim("Tracklet") << oss;
139  for (unsigned int i = 0; i < stubpairs_.size(); ++i) {
140  edm::LogVerbatim("Tracklet") << " " << stubpairs_.at(i)->getName();
141  }
142  }
143 
144  tmpSPTable_.clear();
145 
146  for (unsigned int i = 0; i < stubpairs_.size(); ++i) {
147  for (unsigned int j = 0; j < stubpairs_.at(i)->nStubPairs(); ++j) {
148  if (print)
149  edm::LogVerbatim("Tracklet") << " ***** " << stubpairs_.at(i)->getName() << " "
150  << stubpairs_.at(i)->nStubPairs();
151 
152  auto firstvmstub = stubpairs_.at(i)->getVMStub1(j);
153  auto secondvmstub = stubpairs_.at(i)->getVMStub2(j);
154 
155  if ((layer2_ == 4 && layer3_ == 2) || (layer2_ == 6 && layer3_ == 4)) {
156  constexpr unsigned int vmbitshift = 10;
157  int lookupbits = (int)((firstvmstub.vmbits().value() >> vmbitshift) & 1023); //1023=2^vmbitshift-1
158  int newbin = (lookupbits & 127);
159  int bin = newbin / 8;
160 
161  int start = (bin >> 1);
162  int last = start + (bin & 1);
163 
164  for (int ibin = start; ibin <= last; ibin++) {
165  for (unsigned int k = 0; k < thirdvmstubs_.size(); k++) {
166  string vmsteSuffix = thirdvmstubs_.at(k)->getLastPartOfName();
167  vmsteSuffix = vmsteSuffix.substr(0, vmsteSuffix.find_last_of('n'));
168  if (stubpairs_.at(i)->getLastPartOfName() != vmsteSuffix)
169  continue;
170  for (unsigned int l = 0; l < thirdvmstubs_.at(k)->nVMStubsBinned(ibin); l++) {
171  if (settings_.debugTracklet()) {
172  edm::LogVerbatim("Tracklet") << "In " << getName() << " have third stub";
173  }
174 
175  if (countall >= settings_.maxStep("TRE"))
176  break;
177  countall++;
178 
179  const VMStubTE &thirdvmstub = thirdvmstubs_.at(k)->getVMStubTEBinned(ibin, l);
180 
181  assert(secondphibits_ != -1);
182  assert(thirdphibits_ != -1);
183 
184  unsigned int nvmsecond = settings_.nallstubs(layer2_ - 1) * settings_.nvmte(1, iSeed_);
185  unsigned int nvmbitssecond = nbits(nvmsecond);
186 
187  FPGAWord iphisecondbin = secondvmstub.stub()->iphivmFineBins(nvmbitssecond, secondphibits_);
188 
189  //currently not using same number of bits as in the TED
190  //assert(iphisecondbin==(int)secondvmstub.finephi());
191  FPGAWord iphithirdbin = thirdvmstub.finephi();
192 
193  unsigned int index = (iphisecondbin.value() << thirdphibits_) + iphithirdbin.value();
194 
195  FPGAWord secondbend = secondvmstub.bend();
196  FPGAWord thirdbend = thirdvmstub.bend();
197 
198  index = (index << secondbend.nbits()) + secondbend.value();
199  index = (index << thirdbend.nbits()) + thirdbend.value();
200 
202  (index >= table_.size() || !table_[index])) {
203  if (settings_.debugTracklet()) {
204  edm::LogVerbatim("Tracklet")
205  << "Stub pair rejected because of stub pt cut bends : "
206  << settings_.benddecode(secondvmstub.bend().value(), layer2_ - 1, secondvmstub.isPSmodule())
207  << " " << settings_.benddecode(thirdvmstub.bend().value(), layer3_ - 1, thirdvmstub.isPSmodule());
208  }
209 
210  //FIXME temporarily commented out until bend table fixed
211  //if (!settings_.writeTripletTables())
212  // continue;
213  }
215  if (index >= table_.size())
216  table_.resize(index + 1, false);
217  table_[index] = true;
218 
219  const unsigned spIndex = stubpairs_.at(i)->getIndex(j);
220  const string &tedName = stubpairs_.at(i)->getTEDName(j);
221  if (!tmpSPTable_.count(tedName))
222  tmpSPTable_[tedName];
223  if (spIndex >= tmpSPTable_.at(tedName).size())
224  tmpSPTable_.at(tedName).resize(spIndex + 1);
225  tmpSPTable_.at(tedName).at(spIndex).push_back(stubpairs_.at(i)->getName());
226  }
227 
228  if (settings_.debugTracklet())
229  edm::LogVerbatim("Tracklet") << "Adding layer-layer pair in " << getName();
230  if (settings_.writeMonitorData("Seeds")) {
231  ofstream fout("seeds.txt", ofstream::app);
232  fout << __FILE__ << ":" << __LINE__ << " " << name_ << " " << iSeed_ << endl;
233  fout.close();
234  }
235  stubtriplets_->addStubs(thirdvmstub.stub(),
236  (stubpairs_.at(i))->getVMStub1(j).stub(),
237  (stubpairs_.at(i))->getVMStub2(j).stub());
238 
239  countpass++;
240  }
241  }
242  }
243 
244  }
245 
246  else if (disk2_ == 2 && layer3_ == 2) {
247  int lookupbits = (int)((firstvmstub.vmbits().value() >> 10) & 1023);
248  int newbin = (lookupbits & 127);
249  int bin = newbin / 8;
250 
251  int start = (bin >> 1);
252  int last = start + (bin & 1);
253 
254  if (firstvmstub.stub()->disk().value() < 0) { //TODO - negative disk should come from memory
255  start = settings_.NLONGVMBINS() - last - 1;
256  last = settings_.NLONGVMBINS() - start - 1;
257  }
258 
259  for (int ibin = start; ibin <= last; ibin++) {
260  for (unsigned int k = 0; k < thirdvmstubs_.size(); k++) {
261  string vmsteSuffix = thirdvmstubs_.at(k)->getLastPartOfName();
262  vmsteSuffix = vmsteSuffix.substr(0, vmsteSuffix.find_last_of('n'));
263  if (stubpairs_.at(i)->getLastPartOfName() != vmsteSuffix)
264  continue;
265 
266  for (unsigned int l = 0; l < thirdvmstubs_.at(k)->nVMStubsBinned(ibin); l++) {
267  if (countall >= settings_.maxStep("TRE"))
268  break;
269  countall++;
270 
271  const VMStubTE &thirdvmstub = thirdvmstubs_.at(k)->getVMStubTEBinned(ibin, l);
272 
273  assert(secondphibits_ != -1);
274  assert(thirdphibits_ != -1);
275 
276  FPGAWord iphisecondbin = secondvmstub.finephi();
277  FPGAWord iphithirdbin = thirdvmstub.finephi();
278 
279  unsigned int index = (iphisecondbin.value() << thirdphibits_) + iphithirdbin.value();
280 
281  FPGAWord secondbend = secondvmstub.bend();
282  FPGAWord thirdbend = thirdvmstub.bend();
283 
284  index = (index << secondbend.nbits()) + secondbend.value();
285  index = (index << thirdbend.nbits()) + thirdbend.value();
286 
288  (index >= table_.size() || !table_[index])) {
289  if (settings_.debugTracklet()) {
290  edm::LogVerbatim("Tracklet")
291  << "Stub triplet rejected because of stub pt cut bends : "
292  << settings_.benddecode(secondvmstub.bend().value(), disk2_ + 5, secondvmstub.isPSmodule()) << " "
293  << settings_.benddecode(thirdvmstub.bend().value(), layer3_ - 1, thirdvmstub.isPSmodule());
294  }
295  continue;
296  }
298  if (index >= table_.size())
299  table_.resize(index + 1, false);
300  table_[index] = true;
301 
302  const unsigned spIndex = stubpairs_.at(i)->getIndex(j);
303  const string &tedName = stubpairs_.at(i)->getTEDName(j);
304  if (!tmpSPTable_.count(tedName))
305  tmpSPTable_[tedName];
306  if (spIndex >= tmpSPTable_.at(tedName).size())
307  tmpSPTable_.at(tedName).resize(spIndex + 1);
308  tmpSPTable_.at(tedName).at(spIndex).push_back(stubpairs_.at(i)->getName());
309  }
310 
311  if (settings_.debugTracklet())
312  edm::LogVerbatim("Tracklet") << "Adding layer-disk pair in " << getName();
313  if (settings_.writeMonitorData("Seeds")) {
314  ofstream fout("seeds.txt", ofstream::app);
315  fout << __FILE__ << ":" << __LINE__ << " " << name_ << " " << iSeed_ << endl;
316  fout.close();
317  }
318  stubtriplets_->addStubs(thirdvmstub.stub(),
319  (stubpairs_.at(i))->getVMStub1(j).stub(),
320  (stubpairs_.at(i))->getVMStub2(j).stub());
321  countpass++;
322  }
323  }
324  }
325  }
326 
327  else if (layer2_ == 3 && disk3_ == 1) {
328  int lookupbits = (int)((firstvmstub.vmbits().value() >> 10) & 1023);
329 
330  int newbin = (lookupbits & 127);
331  int bin = newbin / 8;
332 
333  int start = (bin >> 1);
334  int last = start + (bin & 1);
335 
336  for (int ibin = start; ibin <= last; ibin++) {
337  for (unsigned int k = 0; k < thirdvmstubs_.size(); k++) {
338  string vmsteSuffix = thirdvmstubs_.at(k)->getLastPartOfName();
339  vmsteSuffix = vmsteSuffix.substr(0, vmsteSuffix.find_last_of('n'));
340  if (stubpairs_.at(i)->getLastPartOfName() != vmsteSuffix)
341  continue;
342  for (unsigned int l = 0; l < thirdvmstubs_.at(k)->nVMStubsBinned(ibin); l++) {
343  if (countall >= settings_.maxStep("TRE"))
344  break;
345  countall++;
346 
347  const VMStubTE &thirdvmstub = thirdvmstubs_.at(k)->getVMStubTEBinned(ibin, l);
348 
349  assert(secondphibits_ != -1);
350  assert(thirdphibits_ != -1);
351 
352  unsigned int nvmsecond;
353 
354  nvmsecond = settings_.nallstubs(layer2_ - 1) * settings_.nvmte(1, iSeed_);
355  unsigned int nvmbitssecond = nbits(nvmsecond);
356 
357  FPGAWord iphisecondbin = secondvmstub.stub()->iphivmFineBins(nvmbitssecond, secondphibits_);
358 
359  //currentlty not using same number of bits as in the TED
360  //assert(iphisecondbin==(int)secondvmstub.finephi());
361  FPGAWord iphithirdbin = thirdvmstub.finephi();
362 
363  unsigned int index = (iphisecondbin.value() << thirdphibits_) + iphithirdbin.value();
364 
365  FPGAWord secondbend = secondvmstub.bend();
366  FPGAWord thirdbend = thirdvmstub.bend();
367 
368  index = (index << secondbend.nbits()) + secondbend.value();
369  index = (index << thirdbend.nbits()) + thirdbend.value();
370 
372  (index >= table_.size() || !table_[index])) {
373  if (settings_.debugTracklet()) {
374  edm::LogVerbatim("Tracklet")
375  << "Stub pair rejected because of stub pt cut bends : "
376  << settings_.benddecode(secondvmstub.bend().value(), layer2_ - 1, secondvmstub.isPSmodule())
377  << " " << settings_.benddecode(thirdvmstub.bend().value(), disk3_ + 5, thirdvmstub.isPSmodule());
378  }
379  continue;
380  }
382  if (index >= table_.size())
383  table_.resize(index + 1, false);
384  table_[index] = true;
385 
386  const unsigned spIndex = stubpairs_.at(i)->getIndex(j);
387  const string &tedName = stubpairs_.at(i)->getTEDName(j);
388  if (!tmpSPTable_.count(tedName))
389  tmpSPTable_[tedName];
390  if (spIndex >= tmpSPTable_.at(tedName).size())
391  tmpSPTable_.at(tedName).resize(spIndex + 1);
392  tmpSPTable_.at(tedName).at(spIndex).push_back(stubpairs_.at(i)->getName());
393  }
394 
395  if (settings_.debugTracklet())
396  edm::LogVerbatim("Tracklet") << "Adding layer-disk pair in " << getName();
397  if (settings_.writeMonitorData("Seeds")) {
398  ofstream fout("seeds.txt", ofstream::app);
399  fout << __FILE__ << ":" << __LINE__ << " " << name_ << " " << iSeed_ << endl;
400  fout.close();
401  }
402  stubtriplets_->addStubs(thirdvmstub.stub(),
403  (stubpairs_.at(i))->getVMStub1(j).stub(),
404  (stubpairs_.at(i))->getVMStub2(j).stub());
405  countpass++;
406  }
407  }
408  }
409  }
410  }
411  }
412 
413  for (const auto &tedName : tmpSPTable_) {
414  for (unsigned spIndex = 0; spIndex < tedName.second.size(); spIndex++) {
415  if (tedName.second.at(spIndex).empty())
416  continue;
417  vector<string> entry(tedName.second.at(spIndex));
418  sort(entry.begin(), entry.end());
419  entry.erase(unique(entry.begin(), entry.end()), entry.end());
420  const string &spName = entry.at(0);
421 
422  if (!spTable_.count(tedName.first))
423  spTable_[tedName.first];
424  if (spIndex >= spTable_.at(tedName.first).size())
425  spTable_.at(tedName.first).resize(spIndex + 1);
426  if (!spTable_.at(tedName.first).at(spIndex).count(spName))
427  spTable_.at(tedName.first).at(spIndex)[spName] = 0;
428  spTable_.at(tedName.first).at(spIndex)[spName]++;
429  }
430  }
431 
432  if (settings_.writeMonitorData("TRE")) {
433  globals_->ofstream("tripletengine.txt") << getName() << " " << countall << " " << countpass << endl;
434  }
435 }
436 
438  ifstream fin;
439  string tableName, word;
440  unsigned num;
441 
442  string tablePath = settings_.tableTREFile();
443  unsigned int finddir = tablePath.find("table_TRE_");
444  tableName = tablePath.substr(0, finddir) + "table_" + name_ + ".txt";
445 
446  fin.open(tableName, ifstream::in);
447  if (!fin) {
448  throw cms::Exception("BadConfig") << "TripletEngine::readTables, file " << tableName << " not known";
449  }
450  while (!fin.eof()) {
451  fin >> word;
452  num = atoi(word.c_str());
453  table_.push_back(num > 0 ? true : false);
454  }
455  fin.close();
456 }
457 
459  ofstream fout;
460  stringstream tableName;
461 
462  tableName << "table/table_" << name_ << ".txt";
463 
464  fout.open(tableName.str(), ofstream::out);
465  for (const auto entry : table_)
466  fout << entry << endl;
467  fout.close();
468 
469  for (const auto &tedName : spTable_) {
470  tableName.str("");
471  tableName << "table/table_" << tedName.first << "_" << name_ << ".txt";
472 
473  fout.open(tableName.str(), ofstream::out);
474  for (const auto &entry : tedName.second) {
475  for (const auto &spName : entry)
476  fout << spName.first << ":" << spName.second << " ";
477  fout << endl;
478  }
479  fout.close();
480  }
481 }
Definition: start.py:1
Log< level::Info, true > LogVerbatim
std::map< std::string, std::vector< std::vector< std::string > > > tmpSPTable_
Definition: TripletEngine.h:59
unsigned int maxStep(std::string module) const
Definition: Settings.h:123
bool isPSmodule() const
Definition: VMStubTE.h:31
std::string name_
Definition: ProcessBase.h:38
double benddecode(int ibend, int layerdisk, bool isPSmodule) const
Definition: Settings.h:437
bool enableTripletTables() const
Definition: Settings.h:212
const FPGAWord & bend() const
Definition: VMStubTE.h:25
std::string const & tableTREFile() const
Definition: Settings.h:78
Settings const & settings_
Definition: ProcessBase.h:40
Globals * globals_
Definition: ProcessBase.h:41
void addInput(MemoryBase *memory, std::string input) override
bool writetrace() const
Definition: Settings.h:193
assert(be >=bs)
static std::string to_string(const XMLCh *ch)
static std::string const input
Definition: EdmProvDump.cc:50
unsigned int NLONGVMBINS() const
Definition: Settings.h:363
const Stub * stub() const
Definition: VMStubTE.h:29
uint64_t word
unsigned int nbits(unsigned int power)
Definition: ProcessBase.cc:17
std::vector< bool > table_
Definition: TripletEngine.h:61
void addOutput(MemoryBase *memory, std::string output) override
std::vector< VMStubsTEMemory * > thirdvmstubs_
Definition: TripletEngine.h:54
def unique(seq, keepstr=True)
Definition: tier0.py:24
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
Definition: Utilities.cc:47
int value() const
Definition: FPGAWord.h:24
void addStubs(const Stub *stub1, const Stub *stub2, const Stub *stub3)
StubTripletsMemory * stubtriplets_
Definition: TripletEngine.h:57
bool writeMonitorData(std::string module) const
Definition: Settings.h:116
unsigned int nvmte(unsigned int inner, unsigned int iSeed) const
Definition: Settings.h:108
unsigned int nallstubs(unsigned int layerdisk) const
Definition: Settings.h:114
const FPGAWord & finephi() const
Definition: VMStubTE.h:23
bool debugTracklet() const
Definition: Settings.h:192
std::map< std::string, std::vector< std::map< std::string, unsigned > > > spTable_
Definition: TripletEngine.h:60
int nbits() const
Definition: FPGAWord.h:25
int nfinephi(unsigned int inner, unsigned int iSeed) const
Definition: Settings.h:140
std::ofstream & ofstream(std::string fname)
Definition: Globals.cc:44
std::vector< StubPairsMemory * > stubpairs_
Definition: TripletEngine.h:55
Definition: output.py:1
tmp
align.sh
Definition: createJobs.py:716
std::string const & getName() const
Definition: ProcessBase.h:22
bool writeTripletTables() const
Definition: Settings.h:213