CMS 3D CMS Logo

MatchCalculator.cc
Go to the documentation of this file.
1 // MatchCalculator
3 // This class loads pairs for tracklet/stubs and looks for the
4 // best possible match.
5 // Variables such as `best_ideltaphi_barrel` store the "global"
6 // best value for delta phi, r, z, and r*phi, for instances
7 // where the same tracklet has multiple stub pairs. This allows
8 // us to find the truly best match
10 
21 
25 
26 #include <filesystem>
27 #include <algorithm>
28 
29 using namespace std;
30 using namespace trklet;
31 
32 MatchCalculator::MatchCalculator(string name, Settings const& settings, Globals* global)
33  : ProcessBase(name, settings, global),
34  phimatchcuttable_(settings),
35  zmatchcuttable_(settings),
36  rphicutPStable_(settings),
37  rphicut2Stable_(settings),
38  rcutPStable_(settings),
39  rcut2Stable_(settings),
40  alphainner_(settings),
41  alphaouter_(settings),
42  rSSinner_(settings),
43  rSSouter_(settings) {
44  phiregion_ = name[8] - 'A';
46 
47  fullMatches_.resize(12, nullptr);
48 
49  //TODO - need to sort out constants here
50  icorrshift_ = 7;
51 
52  if (layerdisk_ < N_PSLAYER) {
54  } else {
56  }
57  phi0shift_ = 3;
58  fact_ = 1;
59  if (layerdisk_ >= N_PSLAYER && layerdisk_ < N_LAYER) {
60  fact_ = (1 << (settings_.nzbitsstub(0) - settings_.nzbitsstub(5)));
64  phi0shift_ = 0;
65  }
66 
67  unsigned int region = getName()[8] - 'A';
69 
70  if (layerdisk_ < N_LAYER) {
71  phimatchcuttable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::barrelphi, region);
72  zmatchcuttable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::barrelz, region);
73  } else {
74  rphicutPStable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::diskPSphi, region);
75  rphicut2Stable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::disk2Sphi, region);
76  rcutPStable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::diskPSr, region);
77  rcut2Stable_.initmatchcut(layerdisk_, TrackletLUT::MatchType::disk2Sr, region);
78  alphainner_.initmatchcut(layerdisk_, TrackletLUT::MatchType::alphainner, region);
79  alphaouter_.initmatchcut(layerdisk_, TrackletLUT::MatchType::alphaouter, region);
80  rSSinner_.initmatchcut(layerdisk_, TrackletLUT::MatchType::rSSinner, region);
81  rSSouter_.initmatchcut(layerdisk_, TrackletLUT::MatchType::rSSouter, region);
82  }
83 
84  for (unsigned int i = 0; i < N_DSS_MOD * 2; i++) {
86  (1 << (settings_.nbitsalpha() - 1)) / (settings_.rDSSinner(i) * settings_.rDSSinner(i)) /
87  settings_.kphi();
89  (1 << (settings_.nbitsalpha() - 1)) / (settings_.rDSSouter(i) * settings_.rDSSouter(i)) /
90  settings_.kphi();
91  }
92 }
93 
95  if (settings_.writetrace()) {
96  edm::LogVerbatim("Tracklet") << "In " << name_ << " adding output to " << memory->getName() << " to output "
97  << output;
98  }
99  if (output.substr(0, 8) == "matchout") {
100  auto* tmp = dynamic_cast<FullMatchMemory*>(memory);
101  assert(tmp != nullptr);
102  unsigned int iSeed = getISeed(memory->getName());
103  fullMatches_[iSeed] = tmp;
104  return;
105  }
106  throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " could not find output " << output;
107 }
108 
110  if (settings_.writetrace()) {
111  edm::LogVerbatim("Tracklet") << "In " << name_ << " adding input from " << memory->getName() << " to input "
112  << input;
113  }
114  if (input == "allstubin") {
115  auto* tmp = dynamic_cast<AllStubsMemory*>(memory);
116  assert(tmp != nullptr);
117  allstubs_ = tmp;
118  return;
119  }
120  if (input == "allprojin") {
121  auto* tmp = dynamic_cast<AllProjectionsMemory*>(memory);
122  assert(tmp != nullptr);
123  allprojs_ = tmp;
124  return;
125  }
126  if (input.substr(0, 5) == "match" && input.substr(input.size() - 2, 2) == "in") {
127  auto* tmp = dynamic_cast<CandidateMatchMemory*>(memory);
128  assert(tmp != nullptr);
129  matches_.push_back(tmp);
130  return;
131  }
132  throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " could not find input " << input;
133 }
134 
135 void MatchCalculator::execute(unsigned int iSector, double phioffset) {
136  unsigned int countall = 0;
137  unsigned int countsel = 0;
138 
139  //bool print = getName() == "MC_L4PHIC" && iSector == 3;
140 
141  Tracklet* oldTracklet = nullptr;
142 
143  // Get all tracklet/stub pairs
144  std::vector<std::pair<std::pair<Tracklet*, int>, const Stub*> > mergedMatches = mergeMatches(matches_);
145 
146  // Number of clock cycles the pipeline in HLS takes to process the projection merging to
147  // produce the first projection
148  unsigned int mergedepth = 3;
149 
150  unsigned int maxProc = std::min(settings_.maxStep("MC") - mergedepth, (unsigned int)mergedMatches.size());
151 
152  // Pick some initial large values
153  int best_ideltaphi_barrel = 0xFFFF;
154  int best_ideltaz_barrel = 0xFFFF;
155  int best_ideltaphi_disk = 0xFFFF;
156  int best_ideltar_disk = 0xFFFF;
157  unsigned int curr_projid = -1;
158  unsigned int next_projid = -1;
159 
160  for (unsigned int j = 0; j < maxProc; j++) {
161  if (settings_.debugTracklet() && j == 0) {
162  edm::LogVerbatim("Tracklet") << getName() << " has " << mergedMatches.size() << " candidate matches";
163  }
164 
165  countall++;
166 
167  const Stub* fpgastub = mergedMatches[j].second;
168  Tracklet* tracklet = mergedMatches[j].first.first;
169  const L1TStub* stub = fpgastub->l1tstub();
170 
171  //check that the matches are ordered correctly
172  //allow equal here since we can have more than one candidate match per tracklet projection
173  if (oldTracklet != nullptr) {
174  assert(oldTracklet->TCID() <= tracklet->TCID());
175  }
176  oldTracklet = tracklet;
177 
178  if (layerdisk_ < N_LAYER) {
179  //Integer calculation
180 
181  const Projection& proj = tracklet->proj(layerdisk_);
182 
183  int ir = fpgastub->r().value();
184  int iphi = proj.fpgaphiproj().value();
185  int icorr = (ir * proj.fpgaphiprojder().value()) >> icorrshift_;
186  iphi += icorr;
187 
188  int iz = proj.fpgarzproj().value();
189  int izcor = (ir * proj.fpgarzprojder().value() + (1 << (icorzshift_ - 1))) >> icorzshift_;
190  iz += izcor;
191 
192  int ideltaz = fpgastub->z().value() - iz;
193  int ideltaphi = (fpgastub->phi().value() << phi0shift_) - (iphi << (settings_.phi0bitshift() - 1 + phi0shift_));
194 
195  //Floating point calculations
196 
197  double phi = stub->phi() - phioffset;
198  double r = stub->r();
199  double z = stub->z();
200 
201  if (settings_.useapprox()) {
202  double dphi = reco::reduceRange(phi - fpgastub->phiapprox(0.0, 0.0));
203  assert(std::abs(dphi) < 0.001);
204  phi = fpgastub->phiapprox(0.0, 0.0);
205  z = fpgastub->zapprox();
206  r = fpgastub->rapprox();
207  }
208 
209  if (phi < 0)
210  phi += 2 * M_PI;
211 
212  double dr = r - settings_.rmean(layerdisk_);
213  assert(std::abs(dr) < settings_.drmax());
214 
215  double dphi = reco::reduceRange(phi - (proj.phiproj() + dr * proj.phiprojder()));
216 
217  double dz = z - (proj.rzproj() + dr * proj.rzprojder());
218 
219  double dphiapprox = reco::reduceRange(phi - (proj.phiprojapprox() + dr * proj.phiprojderapprox()));
220 
221  double dzapprox = z - (proj.rzprojapprox() + dr * proj.rzprojderapprox());
222 
223  int seedindex = tracklet->getISeed();
224  unsigned int projindex = mergedMatches[j].first.second; // Allproj index
225  curr_projid = next_projid;
226  next_projid = projindex;
227 
228  // Do we have a new tracklet?
229  bool newtracklet = (j == 0 || projindex != curr_projid);
230  if (j == 0)
231  best_ideltar_disk = (1 << (fpgastub->r().nbits() - 1)); // Set to the maximum possible
232  // If so, replace the "best" values with the cut tables
233  if (newtracklet) {
234  best_ideltaphi_barrel = (int)phimatchcuttable_.lookup(seedindex);
235  best_ideltaz_barrel = (int)zmatchcuttable_.lookup(seedindex);
236  }
237 
238  assert(phimatchcuttable_.lookup(seedindex) > 0);
239  assert(zmatchcuttable_.lookup(seedindex) > 0);
240 
241  if (settings_.bookHistos()) {
242  bool truthmatch = tracklet->stubtruthmatch(stub);
243 
245  hists->FillLayerResidual(layerdisk_ + 1,
246  seedindex,
247  dphiapprox * settings_.rmean(layerdisk_),
248  ideltaphi * settings_.kphi1() * settings_.rmean(layerdisk_),
249  ideltaz * fact_ * settings_.kz(),
250  dz,
251  truthmatch);
252  }
253 
254  //This would catch significant consistency problems in the configuration - helps to debug if there are problems.
255  if (std::abs(dphi) > 0.5 * settings_.dphisectorHG() || std::abs(dphiapprox) > 0.5 * settings_.dphisectorHG()) {
256  throw cms::Exception("LogicError")
257  << "WARNING dphi and/or dphiapprox too large : " << dphi << " " << dphiapprox << endl;
258  }
259 
260  if (settings_.writeMonitorData("Residuals")) {
261  double pt = 0.01 * settings_.c() * settings_.bfield() / std::abs(tracklet->rinv());
262 
263  globals_->ofstream("layerresiduals.txt")
264  << layerdisk_ + 1 << " " << seedindex << " " << pt << " "
265  << ideltaphi * settings_.kphi1() * settings_.rmean(layerdisk_) << " "
266  << dphiapprox * settings_.rmean(layerdisk_) << " "
267  << phimatchcuttable_.lookup(seedindex) * settings_.kphi1() * settings_.rmean(layerdisk_) << " "
268  << ideltaz * fact_ * settings_.kz() << " " << dz << " "
269  << zmatchcuttable_.lookup(seedindex) * settings_.kz() << endl;
270  }
271 
272  // integer match
273  bool imatch = (std::abs(ideltaphi) <= best_ideltaphi_barrel) && (ideltaz * fact_ < best_ideltaz_barrel) &&
274  (ideltaz * fact_ >= -best_ideltaz_barrel);
275  // Update the "best" values
276  if (imatch) {
277  best_ideltaphi_barrel = std::abs(ideltaphi);
278  best_ideltaz_barrel = std::abs(ideltaz * fact_);
279  }
280 
281  if (settings_.debugTracklet()) {
282  edm::LogVerbatim("Tracklet") << getName() << " imatch = " << imatch << " ideltaphi cut " << ideltaphi << " "
283  << phimatchcuttable_.lookup(seedindex) << " ideltaz*fact cut " << ideltaz * fact_
284  << " " << zmatchcuttable_.lookup(seedindex);
285  }
286 
287  if (imatch) {
288  countsel++;
289 
290  tracklet->addMatch(layerdisk_,
291  ideltaphi,
292  ideltaz,
293  dphi,
294  dz,
295  dphiapprox,
296  dzapprox,
297  (phiregion_ << 7) + fpgastub->stubindex().value(),
298  mergedMatches[j].second);
299 
300  if (settings_.debugTracklet()) {
301  edm::LogVerbatim("Tracklet") << "Accepted full match in layer " << getName() << " " << tracklet;
302  }
303 
304  fullMatches_[seedindex]->addMatch(tracklet, mergedMatches[j].second);
305  }
306  } else { //disk matches
307  //check that stubs and projections in same half of detector
308  assert(stub->z() * tracklet->t() > 0.0);
309 
310  int sign = (tracklet->t() > 0.0) ? 1 : -1;
311  int disk = sign * (layerdisk_ - (N_LAYER - 1));
312  assert(disk != 0);
313 
314  //Perform integer calculations here
315 
316  const Projection& proj = tracklet->proj(layerdisk_);
317 
318  int iz = fpgastub->z().value();
319  int iphi = proj.fpgaphiproj().value();
320 
321  //TODO - need to express in terms of constants
322  int shifttmp = 6;
323  int iphicorr = (iz * proj.fpgaphiprojder().value()) >> shifttmp;
324 
325  iphi += iphicorr;
326 
327  int ir = proj.fpgarzproj().value();
328 
329  //TODO - need to express in terms of constants
330  int shifttmp2 = 7;
331  int ircorr = (iz * proj.fpgarzprojder().value()) >> shifttmp2;
332 
333  ir += ircorr;
334 
335  int ideltaphi = fpgastub->phi().value() * settings_.kphi() / settings_.kphi() - iphi;
336 
337  int irstub = fpgastub->r().value();
338  int ialphafact = 0;
339  if (!stub->isPSmodule()) {
340  assert(irstub < (int)N_DSS_MOD * 2);
341  if (abs(disk) <= 2) {
342  ialphafact = ialphafactinner_[irstub];
343  irstub = settings_.rDSSinner(irstub) / settings_.kr();
344  } else {
345  ialphafact = ialphafactouter_[irstub];
346  irstub = settings_.rDSSouter(irstub) / settings_.kr();
347  }
348  }
349 
350  //TODO stub and projection r should not use different # bits...
351  int ideltar = (irstub >> 1) - ir;
352 
353  if (!stub->isPSmodule()) {
354  int ialpha = fpgastub->alpha().value();
355  int iphialphacor = ((ideltar * ialpha * ialphafact) >> settings_.alphashift());
356  ideltaphi += iphialphacor;
357  }
358 
359  //Perform floating point calculations here
360 
361  double phi = stub->phi() - phioffset;
362  double z = stub->z();
363  double r = stub->r();
364 
365  if (settings_.useapprox()) {
366  double dphi = reco::reduceRange(phi - fpgastub->phiapprox(0.0, 0.0));
367  assert(std::abs(dphi) < 0.001);
368  phi = fpgastub->phiapprox(0.0, 0.0);
369  z = fpgastub->zapprox();
370  r = fpgastub->rapprox();
371  }
372 
373  if (phi < 0)
374  phi += 2 * M_PI;
375 
376  double dz = z - sign * settings_.zmean(layerdisk_ - N_LAYER);
377 
378  if (std::abs(dz) > settings_.dzmax()) {
379  throw cms::Exception("LogicError") << __FILE__ << " " << __LINE__ << " " << name_ << " " << tracklet->getISeed()
380  << "\n stub " << stub->z() << " disk " << disk << " " << dz;
381  }
382 
383  double phiproj = proj.phiproj() + dz * proj.phiprojder();
384 
385  double rproj = proj.rzproj() + dz * proj.rzprojder();
386 
387  double deltar = r - rproj;
388 
389  double dr = stub->r() - rproj;
390 
391  double dphi = reco::reduceRange(phi - phiproj);
392 
393  double dphiapprox = reco::reduceRange(phi - (proj.phiprojapprox() + dz * proj.phiprojderapprox()));
394 
395  double drapprox = stub->r() - (proj.rzprojapprox() + dz * proj.rzprojderapprox());
396 
397  double drphi = dphi * stub->r();
398  double drphiapprox = dphiapprox * stub->r();
399 
400  if (!stub->isPSmodule()) {
401  double alphanorm = stub->alphanorm();
402  dphi += dr * alphanorm * settings_.half2SmoduleWidth() / stub->r2();
403  dphiapprox += drapprox * alphanorm * settings_.half2SmoduleWidth() / stub->r2();
404 
405  drphi += dr * alphanorm * settings_.half2SmoduleWidth() / stub->r();
406  drphiapprox += dr * alphanorm * settings_.half2SmoduleWidth() / stub->r();
407  }
408 
409  int seedindex = tracklet->getISeed();
410 
411  int idrphicut = rphicutPStable_.lookup(seedindex);
412  int idrcut = rcutPStable_.lookup(seedindex);
413  if (!stub->isPSmodule()) {
414  idrphicut = rphicut2Stable_.lookup(seedindex);
415  idrcut = rcut2Stable_.lookup(seedindex);
416  }
417 
418  unsigned int projindex = mergedMatches[j].first.second; // Allproj index
419  curr_projid = next_projid;
420  next_projid = projindex;
421  // Do we have a new tracklet?
422  bool newtracklet = (j == 0 || projindex != curr_projid);
423  // If so, replace the "best" values with the cut tables
424  if (newtracklet) {
425  best_ideltaphi_disk = idrphicut;
426  best_ideltar_disk = idrcut;
427  }
428 
429  double drphicut = idrphicut * settings_.kphi() * settings_.kr();
430  double drcut = idrcut * settings_.krprojshiftdisk();
431 
432  bool match, imatch;
433  if (std::abs(dphi) < third * settings_.dphisectorHG() &&
434  std::abs(dphiapprox) < third * settings_.dphisectorHG()) { //1/3 of sector size to catch errors
435  if (settings_.writeMonitorData("Residuals")) {
436  double pt = 0.01 * settings_.c() * settings_.bfield() / std::abs(tracklet->rinv());
437 
438  globals_->ofstream("diskresiduals.txt")
439  << disk << " " << stub->isPSmodule() << " " << tracklet->layer() << " " << abs(tracklet->disk()) << " "
440  << pt << " " << ideltaphi * settings_.kphi() * stub->r() << " " << drphiapprox << " " << drphicut << " "
441  << ideltar * settings_.krprojshiftdisk() << " " << deltar << " " << drcut << " " << endl;
442  }
443 
444  // floating point match
445  match = (std::abs(drphi) < drphicut) && (std::abs(deltar) < drcut);
446 
447  // integer match
448  imatch = (std::abs(ideltaphi) * irstub < best_ideltaphi_disk) && (std::abs(ideltar) < best_ideltar_disk);
449  // Update the "best" values
450  if (imatch) {
451  best_ideltaphi_disk = std::abs(ideltaphi) * irstub;
452  best_ideltar_disk = std::abs(ideltar);
453  }
454  } else {
455  edm::LogProblem("Tracklet") << "WARNING dphi and/or dphiapprox too large : " << dphi << " " << dphiapprox
456  << "dphi " << dphi << " Seed / ISeed " << tracklet->getISeed() << endl;
457  match = false;
458  imatch = false;
459  }
460  if (settings_.debugTracklet()) {
461  edm::LogVerbatim("Tracklet") << "imatch match disk: " << imatch << " " << match << " " << std::abs(ideltaphi)
462  << " " << drphicut / (settings_.kphi() * stub->r()) << " " << std::abs(ideltar)
463  << " " << drcut / settings_.krprojshiftdisk() << " r = " << stub->r();
464  }
465 
466  if (imatch) {
467  countsel++;
468 
469  if (settings_.debugTracklet()) {
470  edm::LogVerbatim("Tracklet") << "MatchCalculator found match in disk " << getName();
471  }
472 
473  tracklet->addMatch(layerdisk_,
474  ideltaphi,
475  ideltar,
476  drphi / stub->r(),
477  dr,
478  drphiapprox / stub->r(),
479  drapprox,
480  (phiregion_ << 7) + fpgastub->stubindex().value(),
481  fpgastub);
482 
483  if (settings_.debugTracklet()) {
484  edm::LogVerbatim("Tracklet") << "Accepted full match in disk " << getName() << " " << tracklet;
485  }
486 
487  fullMatches_[seedindex]->addMatch(tracklet, mergedMatches[j].second);
488  }
489  }
490  if (countall >= settings_.maxStep("MC"))
491  break;
492  }
493 
494  if (settings_.writeMonitorData("MC")) {
495  globals_->ofstream("matchcalculator.txt") << getName() << " " << countall << " " << countsel << endl;
496  }
497 }
498 
499 // Combines all tracklet/stub pairs into a vector
500 std::vector<std::pair<std::pair<Tracklet*, int>, const Stub*> > MatchCalculator::mergeMatches(
501  vector<CandidateMatchMemory*>& candmatch) {
502  std::vector<std::pair<std::pair<Tracklet*, int>, const Stub*> > tmp;
503 
504  std::vector<unsigned int> indexArray;
505  indexArray.reserve(candmatch.size());
506  for (unsigned int i = 0; i < candmatch.size(); i++) {
507  indexArray.push_back(0);
508  }
509 
510  int bestIndex = -1;
511  do {
512  int bestSector = 100;
513  int bestTCID = -1;
514  bestIndex = -1;
515  for (unsigned int i = 0; i < candmatch.size(); i++) {
516  if (indexArray[i] >= candmatch[i]->nMatches()) {
517  // skip as we were at the end
518  continue;
519  }
520  int TCID = candmatch[i]->getMatch(indexArray[i]).first.first->TCID();
521  int dSector = 0;
522  if (dSector > 2)
523  dSector -= N_SECTOR;
524  if (dSector < -2)
525  dSector += N_SECTOR;
526  assert(abs(dSector) < 2);
527  if (dSector == -1)
528  dSector = 2;
529  if (dSector < bestSector) {
530  bestSector = dSector;
531  bestTCID = TCID;
532  bestIndex = i;
533  }
534  if (dSector == bestSector) {
535  if (TCID < bestTCID || bestTCID < 0) {
536  bestTCID = TCID;
537  bestIndex = i;
538  }
539  }
540  }
541  if (bestIndex != -1) {
542  tmp.push_back(candmatch[bestIndex]->getMatch(indexArray[bestIndex]));
543  indexArray[bestIndex]++;
544  }
545  } while (bestIndex != -1);
546 
547  if (layerdisk_ < N_LAYER) {
548  int lastTCID = -1;
549  bool error = false;
550 
551  //Allow equal TCIDs since we can have multiple candidate matches
552  for (unsigned int i = 1; i < tmp.size(); i++) {
553  if (lastTCID > tmp[i].first.first->TCID()) {
554  edm::LogProblem("Tracklet") << "Wrong TCID ordering for projections in " << getName() << " last " << lastTCID
555  << " " << tmp[i].first.first->TCID();
556  error = true;
557  } else {
558  lastTCID = tmp[i].first.first->TCID();
559  }
560  }
561 
562  if (error) {
563  for (unsigned int i = 1; i < tmp.size(); i++) {
564  edm::LogProblem("Tracklet") << "Wrong order for in " << getName() << " " << i << " " << tmp[i].first.first
565  << " " << tmp[i].first.first->TCID();
566  }
567  }
568  }
569 
570  for (unsigned int i = 0; i < tmp.size(); i++) {
571  if (i > 0) {
572  //This allows for equal TCIDs. This means that we can e.g. have a track seeded
573  //in L1L2 that projects to both L3 and D4. The algorithm will pick up the first hit and
574  //drop the second
575 
576  assert(tmp[i - 1].first.first->TCID() <= tmp[i].first.first->TCID());
577  }
578  }
579 
580  return tmp;
581 }
Log< level::Info, true > LogVerbatim
double kz() const
Definition: Settings.h:342
double phi() const
Definition: L1TStub.h:65
double t() const
Definition: Tracklet.h:123
double zapprox() const
Definition: Stub.cc:166
unsigned int maxStep(std::string module) const
Definition: Settings.h:125
const FPGAWord & r() const
Definition: Stub.h:65
int PS_zderL_shift() const
Definition: Settings.h:392
bool bookHistos() const
Definition: Settings.h:219
unsigned int nrbitsstub(unsigned int layerdisk) const
Definition: Settings.h:93
std::string name_
Definition: ProcessBase.h:38
int disk() const
Definition: Tracklet.cc:801
double kphi1() const
Definition: Settings.h:339
void initmatchcut(unsigned int layerdisk, MatchType type, unsigned int region)
Definition: TrackletLUT.cc:190
double rDSSinner(unsigned int iBin) const
Definition: Settings.h:183
constexpr T reduceRange(T x)
Definition: deltaPhi.h:18
double phiapprox(double phimin, double) const
Definition: Stub.cc:185
const FPGAWord & z() const
Definition: Stub.h:66
double dphisectorHG() const
Definition: Settings.h:320
Projection & proj(int layerdisk)
Definition: Tracklet.h:87
Settings const & settings_
Definition: ProcessBase.h:40
double z() const
Definition: L1TStub.h:59
Globals * globals_
Definition: ProcessBase.h:41
int lookup(unsigned int index) const
bool writetrace() const
Definition: Settings.h:195
double dzmax() const
Definition: Settings.h:138
constexpr unsigned int N_DSS_MOD
Definition: Settings.h:31
assert(be >=bs)
int TCID() const
Definition: Tracklet.h:214
constexpr double third
Definition: Settings.h:46
static std::string const input
Definition: EdmProvDump.cc:50
unsigned int isPSmodule() const
Definition: L1TStub.h:103
void addMatch(unsigned int layerdisk, int ideltaphi, int ideltarz, double dphi, double drz, double dphiapprox, double drzapprox, int stubid, const trklet::Stub *stubptr)
Definition: Tracklet.cc:303
U second(std::pair< T, U > const &p)
void addOutput(MemoryBase *memory, std::string output) override
double half2SmoduleWidth() const
Definition: Settings.h:140
unsigned int nzbitsstub(unsigned int layerdisk) const
Definition: Settings.h:91
double rmean(unsigned int iLayer) const
Definition: Settings.h:173
AllStubsMemory * allstubs_
double rinv() const
Definition: Tracklet.h:120
std::vector< CandidateMatchMemory * > matches_
AllProjectionsMemory * allprojs_
double rDSSouter(unsigned int iBin) const
Definition: Settings.h:186
double bfield() const
Definition: Settings.h:277
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int value() const
Definition: FPGAWord.h:24
int alphashift() const
Definition: Settings.h:228
L1TStub * l1tstub()
Definition: Stub.h:83
bool writeMonitorData(std::string module) const
Definition: Settings.h:118
double rapprox() const
Definition: Stub.cc:152
double zmean(unsigned int iDisk) const
Definition: Settings.h:176
void initLayerDisk(unsigned int pos, int &layer, int &disk)
Definition: ProcessBase.cc:33
const FPGAWord & stubindex() const
Definition: Stub.h:72
#define M_PI
unsigned int nallstubs(unsigned int layerdisk) const
Definition: Settings.h:116
bool debugTracklet() const
Definition: Settings.h:194
constexpr unsigned int N_SECTOR
Definition: Settings.h:23
std::vector< std::pair< std::pair< Tracklet *, int >, const Stub * > > mergeMatches(std::vector< CandidateMatchMemory *> &candmatch)
constexpr unsigned int N_PSLAYER
Definition: Settings.h:27
double kr() const
Definition: Settings.h:344
void execute(unsigned int iSector, double phioffset)
double drmax() const
Definition: Settings.h:137
double alphanorm() const
Definition: L1TStub.cc:91
bool stubtruthmatch(const L1TStub *stub)
Definition: Tracklet.cc:138
int nbits() const
Definition: FPGAWord.h:25
Definition: deltar.py:1
int getISeed() const
Definition: Tracklet.cc:820
void addInput(MemoryBase *memory, std::string input) override
int ialphafactouter_[N_DSS_MOD *2]
double c() const
Definition: Settings.h:224
unsigned int getISeed(const std::string &name)
Definition: ProcessBase.cc:119
const FPGAWord & alpha() const
Definition: Stub.h:70
double r() const
Definition: L1TStub.h:60
std::ofstream & ofstream(std::string fname)
Definition: Globals.cc:44
int nbitsalpha() const
Definition: Settings.h:229
HistBase *& histograms()
Definition: Globals.h:38
Definition: output.py:1
double krprojshiftdisk() const
Definition: Settings.h:441
const FPGAWord & phi() const
Definition: Stub.h:68
int ialphafactinner_[N_DSS_MOD *2]
int layer() const
Definition: Tracklet.cc:792
double kphi() const
Definition: Settings.h:338
tmp
align.sh
Definition: createJobs.py:716
double r2() const
Definition: L1TStub.h:62
std::string const & getName() const
Definition: ProcessBase.h:22
int SS_zderL_shift() const
Definition: Settings.h:393
std::vector< FullMatchMemory * > fullMatches_
bool useapprox() const
Definition: Settings.h:247
Log< level::Error, true > LogProblem
int phi0bitshift() const
Definition: Settings.h:403
constexpr int N_LAYER
Definition: Settings.h:25