CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TrackletLUT.cc
Go to the documentation of this file.
5 
6 #include <filesystem>
7 
8 using namespace std;
9 using namespace trklet;
10 
11 TrackletLUT::TrackletLUT(const Settings& settings) : settings_(settings) {}
12 
13 void TrackletLUT::initmatchcut(unsigned int layerdisk, MatchType type, unsigned int region) {
14  char cregion = 'A' + region;
15 
16  for (unsigned int iSeed = 0; iSeed < 12; iSeed++) {
17  if (type == barrelphi) {
18  table_.push_back(settings_.rphimatchcut(iSeed, layerdisk) / (settings_.kphi1() * settings_.rmean(layerdisk)));
19  }
20  if (type == barrelz) {
21  table_.push_back(settings_.zmatchcut(iSeed, layerdisk) / settings_.kz());
22  }
23  if (type == diskPSphi) {
24  table_.push_back(settings_.rphicutPS(iSeed, layerdisk - N_LAYER) / (settings_.kphi() * settings_.kr()));
25  }
26  if (type == disk2Sphi) {
27  table_.push_back(settings_.rphicut2S(iSeed, layerdisk - N_LAYER) / (settings_.kphi() * settings_.kr()));
28  }
29  if (type == disk2Sr) {
30  table_.push_back(settings_.rcut2S(iSeed, layerdisk - N_LAYER) / settings_.krprojshiftdisk());
31  }
32  if (type == diskPSr) {
33  table_.push_back(settings_.rcutPS(iSeed, layerdisk - N_LAYER) / settings_.krprojshiftdisk());
34  }
35  }
36 
37  name_ = settings_.combined() ? "MP_" : "MC_";
38 
39  if (type == barrelphi) {
40  name_ += TrackletConfigBuilder::LayerName(layerdisk) + "PHI" + cregion + "_phicut.tab";
41  }
42  if (type == barrelz) {
43  name_ += TrackletConfigBuilder::LayerName(layerdisk) + "PHI" + cregion + "_zcut.tab";
44  }
45  if (type == diskPSphi) {
46  name_ += TrackletConfigBuilder::LayerName(layerdisk) + "PHI" + cregion + "_PSphicut.tab";
47  }
48  if (type == disk2Sphi) {
49  name_ += TrackletConfigBuilder::LayerName(layerdisk) + "PHI" + cregion + "_2Sphicut.tab";
50  }
51  if (type == disk2Sr) {
52  name_ += TrackletConfigBuilder::LayerName(layerdisk) + "PHI" + cregion + "_2Srcut.tab";
53  }
54  if (type == diskPSr) {
55  name_ += TrackletConfigBuilder::LayerName(layerdisk) + "PHI" + cregion + "_PSrcut.tab";
56  }
57 
58  positive_ = false;
59 
60  writeTable();
61 }
62 
63 void TrackletLUT::initTPlut(bool fillInner,
64  unsigned int iSeed,
65  unsigned int layerdisk1,
66  unsigned int layerdisk2,
67  unsigned int nbitsfinephidiff,
68  unsigned int iTP) {
69  //number of fine phi bins in sector
70  int nfinephibins = settings_.nallstubs(layerdisk2) * settings_.nvmte(1, iSeed) * (1 << settings_.nfinephi(1, iSeed));
71  double dfinephi = settings_.dphisectorHG() / nfinephibins;
72 
73  int outerrbits = 3;
74 
75  if (iSeed == Seed::L1L2 || iSeed == Seed::L2L3 || iSeed == Seed::L3L4 || iSeed == Seed::L5L6) {
76  outerrbits = 0;
77  }
78 
79  int outerrbins = (1 << outerrbits);
80 
81  double dphi[2];
82  double router[2];
83 
84  unsigned int nbendbitsinner = 3;
85  unsigned int nbendbitsouter = 3;
86  if (iSeed == Seed::L3L4) {
87  nbendbitsouter = 4;
88  } else if (iSeed == Seed::L5L6) {
89  nbendbitsinner = 4;
90  nbendbitsouter = 4;
91  }
92 
93  int nbinsfinephidiff = (1 << nbitsfinephidiff);
94 
95  for (int iphibin = 0; iphibin < nbinsfinephidiff; iphibin++) {
96  int iphidiff = iphibin;
97  if (iphibin >= nbinsfinephidiff / 2) {
98  iphidiff = iphibin - nbinsfinephidiff;
99  }
100  //min and max dphi
101  dphi[0] = (iphidiff - 1.5) * dfinephi;
102  dphi[1] = (iphidiff + 1.5) * dfinephi;
103  for (int irouterbin = 0; irouterbin < outerrbins; irouterbin++) {
104  if (iSeed == Seed::D1D2 || iSeed == Seed::D3D4 || iSeed == Seed::L1D1 || iSeed == Seed::L2D1) {
105  router[0] =
106  settings_.rmindiskvm() + irouterbin * (settings_.rmaxdiskvm() - settings_.rmindiskvm()) / outerrbins;
107  router[1] =
108  settings_.rmindiskvm() + (irouterbin + 1) * (settings_.rmaxdiskvm() - settings_.rmindiskvm()) / outerrbins;
109  } else {
110  router[0] = settings_.rmean(layerdisk2);
111  router[1] = settings_.rmean(layerdisk2);
112  }
113 
114  double bendinnermin = 20.0;
115  double bendinnermax = -20.0;
116  double bendoutermin = 20.0;
117  double bendoutermax = -20.0;
118  double rinvmin = 1.0;
119  for (int i2 = 0; i2 < 2; i2++) {
120  for (int i3 = 0; i3 < 2; i3++) {
121  double rinner = 0.0;
122  if (iSeed == Seed::D1D2 || iSeed == Seed::D3D4) {
123  rinner = router[i3] * settings_.zmean(layerdisk1 - N_LAYER) / settings_.zmean(layerdisk2 - N_LAYER);
124  } else {
125  rinner = settings_.rmean(layerdisk1);
126  }
127  double rinv1 = rinv(0.0, -dphi[i2], rinner, router[i3]);
128  double pitchinner = (rinner < settings_.rcrit()) ? settings_.stripPitch(true) : settings_.stripPitch(false);
129  double pitchouter =
130  (router[i3] < settings_.rcrit()) ? settings_.stripPitch(true) : settings_.stripPitch(false);
131  double abendinner = bendstrip(rinner, rinv1, pitchinner);
132  double abendouter = bendstrip(router[i3], rinv1, pitchouter);
133  if (abendinner < bendinnermin)
134  bendinnermin = abendinner;
135  if (abendinner > bendinnermax)
136  bendinnermax = abendinner;
137  if (abendouter < bendoutermin)
138  bendoutermin = abendouter;
139  if (abendouter > bendoutermax)
140  bendoutermax = abendouter;
141  if (std::abs(rinv1) < rinvmin) {
142  rinvmin = std::abs(rinv1);
143  }
144  }
145  }
146 
147  bool passptcut = rinvmin < settings_.rinvcutte();
148 
149  if (fillInner) {
150  for (int ibend = 0; ibend < (1 << nbendbitsinner); ibend++) {
151  double bend = settings_.benddecode(ibend, layerdisk1, nbendbitsinner == 3);
152 
153  bool passinner = bend <= bendinnermax + settings_.bendcutte(ibend, layerdisk1, nbendbitsinner == 3) &&
154  bend >= bendinnermin - settings_.bendcutte(ibend, layerdisk1, nbendbitsinner == 3);
155  table_.push_back(passinner && passptcut);
156  }
157  } else {
158  for (int ibend = 0; ibend < (1 << nbendbitsouter); ibend++) {
159  double bend = settings_.benddecode(ibend, layerdisk2, nbendbitsouter == 3);
160 
161  bool passouter = bend <= bendoutermax + settings_.bendcutte(ibend, layerdisk2, nbendbitsouter == 3) &&
162  bend >= bendoutermin - settings_.bendcutte(ibend, layerdisk2, nbendbitsouter == 3);
163  table_.push_back(passouter && passptcut);
164  }
165  }
166  }
167  }
168 
169  positive_ = false;
170  char cTP = 'A' + iTP;
171 
172  name_ = "TP_" + TrackletConfigBuilder::LayerName(layerdisk1) + TrackletConfigBuilder::LayerName(layerdisk2) + cTP;
173 
174  if (fillInner) {
175  name_ += "_stubptinnercut.tab";
176  } else {
177  name_ += "_stubptoutercut.tab";
178  }
179 
180  writeTable();
181 }
182 
183 void TrackletLUT::initTPregionlut(unsigned int iSeed,
184  unsigned int layerdisk1,
185  unsigned int layerdisk2,
186  unsigned int iAllStub,
187  unsigned int nbitsfinephidiff,
188  unsigned int nbitsfinephi,
189  const TrackletLUT& tplutinner,
190  unsigned int iTP) {
191  int nirbits = 0;
192  if (iSeed == Seed::D1D2 || iSeed == Seed::D3D4 || iSeed == Seed::L1D1 || iSeed == Seed::L2D1) {
193  nirbits = 3;
194  }
195 
196  unsigned int nbendbitsinner = 3;
197 
198  if (iSeed == Seed::L5L6) {
199  nbendbitsinner = 4;
200  }
201 
202  for (int innerfinephi = 0; innerfinephi < (1 << nbitsfinephi); innerfinephi++) {
203  for (int innerbend = 0; innerbend < (1 << nbendbitsinner); innerbend++) {
204  for (int ir = 0; ir < (1 << nirbits); ir++) {
205  unsigned int usereg = 0;
206  for (unsigned int ireg = 0; ireg < settings_.nvmte(1, iSeed); ireg++) {
207  bool match = false;
208  for (int ifinephiouter = 0; ifinephiouter < (1 << settings_.nfinephi(1, iSeed)); ifinephiouter++) {
209  int outerfinephi = iAllStub * (1 << (nbitsfinephi - settings_.nbitsallstubs(layerdisk2))) +
210  ireg * (1 << settings_.nfinephi(1, iSeed)) + ifinephiouter;
211  int idphi = outerfinephi - innerfinephi;
212  bool inrange = (idphi < (1 << (nbitsfinephidiff - 1))) && (idphi >= -(1 << (nbitsfinephidiff - 1)));
213  if (idphi < 0)
214  idphi = idphi + (1 << nbitsfinephidiff);
215  int idphi1 = idphi;
216  if (iSeed >= 4)
217  idphi1 = (idphi << 3) + ir;
218  int ptinnerindexnew = (idphi1 << nbendbitsinner) + innerbend;
219  match = match || (inrange && tplutinner.lookup(ptinnerindexnew));
220  }
221  if (match) {
222  usereg = usereg | (1 << ireg);
223  }
224  }
225 
226  table_.push_back(usereg);
227  }
228  }
229  }
230 
231  positive_ = false;
232  char cTP = 'A' + iTP;
233 
234  name_ = "TP_" + TrackletConfigBuilder::LayerName(layerdisk1) + TrackletConfigBuilder::LayerName(layerdisk2) + cTP +
235  "_usereg.tab";
236 
237  writeTable();
238 }
239 
240 void TrackletLUT::initteptlut(bool fillInner,
241  bool fillTEMem,
242  unsigned int iSeed,
243  unsigned int layerdisk1,
244  unsigned int layerdisk2,
245  unsigned int innerphibits,
246  unsigned int outerphibits,
247  double innerphimin,
248  double innerphimax,
249  double outerphimin,
250  double outerphimax,
251  const std::string& innermem,
252  const std::string& outermem) {
253  int outerrbits = 0;
254  if (iSeed == Seed::D1D2 || iSeed == Seed::D3D4 || iSeed == Seed::L1D1 || iSeed == Seed::L2D1) {
255  outerrbits = 3;
256  }
257 
258  int outerrbins = (1 << outerrbits);
259  int innerphibins = (1 << innerphibits);
260  int outerphibins = (1 << outerphibits);
261 
262  double phiinner[2];
263  double phiouter[2];
264  double router[2];
265 
266  unsigned int nbendbitsinner = 3;
267  unsigned int nbendbitsouter = 3;
268  if (iSeed == Seed::L3L4) {
269  nbendbitsouter = 4;
270  }
271  if (iSeed == Seed::L5L6) {
272  nbendbitsinner = 4;
273  nbendbitsouter = 4;
274  }
275 
276  if (fillTEMem) {
277  if (fillInner) {
278  table_.resize((1 << nbendbitsinner), false);
279  } else {
280  table_.resize((1 << nbendbitsouter), false);
281  }
282  }
283 
284  for (int iphiinnerbin = 0; iphiinnerbin < innerphibins; iphiinnerbin++) {
285  phiinner[0] = innerphimin + iphiinnerbin * (innerphimax - innerphimin) / innerphibins;
286  phiinner[1] = innerphimin + (iphiinnerbin + 1) * (innerphimax - innerphimin) / innerphibins;
287  for (int iphiouterbin = 0; iphiouterbin < outerphibins; iphiouterbin++) {
288  phiouter[0] = outerphimin + iphiouterbin * (outerphimax - outerphimin) / outerphibins;
289  phiouter[1] = outerphimin + (iphiouterbin + 1) * (outerphimax - outerphimin) / outerphibins;
290  for (int irouterbin = 0; irouterbin < outerrbins; irouterbin++) {
291  if (iSeed == Seed::D1D2 || iSeed == Seed::D3D4 || iSeed == Seed::L1D1 || iSeed == Seed::L2D1) {
292  router[0] =
293  settings_.rmindiskvm() + irouterbin * (settings_.rmaxdiskvm() - settings_.rmindiskvm()) / outerrbins;
294  router[1] = settings_.rmindiskvm() +
295  (irouterbin + 1) * (settings_.rmaxdiskvm() - settings_.rmindiskvm()) / outerrbins;
296  } else {
297  router[0] = settings_.rmean(layerdisk2);
298  router[1] = settings_.rmean(layerdisk2);
299  }
300 
301  double bendinnermin = 20.0;
302  double bendinnermax = -20.0;
303  double bendoutermin = 20.0;
304  double bendoutermax = -20.0;
305  double rinvmin = 1.0;
306  for (int i1 = 0; i1 < 2; i1++) {
307  for (int i2 = 0; i2 < 2; i2++) {
308  for (int i3 = 0; i3 < 2; i3++) {
309  double rinner = 0.0;
310  if (iSeed == Seed::D1D2 || iSeed == Seed::D3D4) {
311  rinner = router[i3] * settings_.zmean(layerdisk1 - N_LAYER) / settings_.zmean(layerdisk2 - N_LAYER);
312  } else {
313  rinner = settings_.rmean(layerdisk1);
314  }
315  double rinv1 = -rinv(phiinner[i1], phiouter[i2], rinner, router[i3]);
316  double pitchinner =
317  (rinner < settings_.rcrit()) ? settings_.stripPitch(true) : settings_.stripPitch(false);
318  double pitchouter =
319  (router[i3] < settings_.rcrit()) ? settings_.stripPitch(true) : settings_.stripPitch(false);
320  double abendinner = bendstrip(rinner, rinv1, pitchinner);
321  double abendouter = bendstrip(router[i3], rinv1, pitchouter);
322  if (abendinner < bendinnermin)
323  bendinnermin = abendinner;
324  if (abendinner > bendinnermax)
325  bendinnermax = abendinner;
326  if (abendouter < bendoutermin)
327  bendoutermin = abendouter;
328  if (abendouter > bendoutermax)
329  bendoutermax = abendouter;
330  if (std::abs(rinv1) < rinvmin) {
331  rinvmin = std::abs(rinv1);
332  }
333  }
334  }
335  }
336 
337  bool passptcut = rinvmin < settings_.rinvcutte();
338 
339  if (fillInner) {
340  for (int ibend = 0; ibend < (1 << nbendbitsinner); ibend++) {
341  double bend = settings_.benddecode(ibend, layerdisk1, nbendbitsinner == 3);
342 
343  bool passinner = bend > bendinnermin - settings_.bendcutte(ibend, layerdisk1, nbendbitsinner == 3) &&
344  bend < bendinnermax + settings_.bendcutte(ibend, layerdisk1, nbendbitsinner == 3);
345 
346  if (fillTEMem) {
347  if (passinner) {
348  table_[ibend] = 1;
349  }
350  } else {
351  table_.push_back(passinner && passptcut);
352  }
353  }
354  } else {
355  for (int ibend = 0; ibend < (1 << nbendbitsouter); ibend++) {
356  double bend = settings_.benddecode(ibend, layerdisk2, nbendbitsouter == 3);
357 
358  bool passouter = bend > bendoutermin - settings_.bendcutte(ibend, layerdisk2, nbendbitsouter == 3) &&
359  bend < bendoutermax + settings_.bendcutte(ibend, layerdisk2, nbendbitsouter == 3);
360  if (fillTEMem) {
361  if (passouter) {
362  table_[ibend] = 1;
363  }
364  } else {
365  table_.push_back(passouter && passptcut);
366  }
367  }
368  }
369  }
370  }
371  }
372 
373  positive_ = false;
374 
375  if (fillTEMem) {
376  if (fillInner) {
377  name_ = "VMSTE_" + innermem + "_vmbendcut.tab";
378  } else {
379  name_ = "VMSTE_" + outermem + "_vmbendcut.tab";
380  }
381  } else {
382  name_ = "TE_" + innermem.substr(0, innermem.size() - 2) + "_" + outermem.substr(0, outermem.size() - 2);
383  if (fillInner) {
384  name_ += "_stubptinnercut.tab";
385  } else {
386  name_ += "_stubptoutercut.tab";
387  }
388  }
389 
390  writeTable();
391 }
392 
393 void TrackletLUT::initProjectionBend(double k_phider,
394  unsigned int idisk,
395  unsigned int nrbits,
396  unsigned int nphiderbits) {
397  unsigned int nsignbins = 2;
398  unsigned int nrbins = 1 << (nrbits);
399  unsigned int nphiderbins = 1 << (nphiderbits);
400 
401  for (unsigned int isignbin = 0; isignbin < nsignbins; isignbin++) {
402  for (unsigned int irbin = 0; irbin < nrbins; irbin++) {
403  int ir = irbin;
404  if (ir > (1 << (nrbits - 1)))
405  ir -= (1 << nrbits);
406  ir = ir << (settings_.nrbitsstub(N_LAYER) - nrbits);
407  for (unsigned int iphiderbin = 0; iphiderbin < nphiderbins; iphiderbin++) {
408  int iphider = iphiderbin;
409  if (iphider > (1 << (nphiderbits - 1)))
410  iphider -= (1 << nphiderbits);
411  iphider = iphider << (settings_.nbitsphiprojderL123() - nphiderbits);
412 
413  double rproj = ir * settings_.krprojshiftdisk();
414  double phider = iphider * k_phider;
415  double t = settings_.zmean(idisk) / rproj;
416 
417  if (isignbin)
418  t = -t;
419 
420  double rinv = -phider * (2.0 * t);
421 
422  double stripPitch = (rproj < settings_.rcrit()) ? settings_.stripPitch(true) : settings_.stripPitch(false);
423  double bendproj = bendstrip(rproj, rinv, stripPitch);
424 
425  int ibendproj = 2.0 * bendproj + 15.5;
426  if (ibendproj < 0)
427  ibendproj = 0;
428  if (ibendproj > 31)
429  ibendproj = 31;
430 
431  table_.push_back(ibendproj);
432  }
433  }
434  }
435 
436  positive_ = false;
437  name_ = settings_.combined() ? "MP_" : "PR_";
438  name_ += "ProjectionBend_" + TrackletConfigBuilder::LayerName(N_LAYER + idisk) + ".tab";
439 
440  writeTable();
441 }
442 
443 void TrackletLUT::initBendMatch(unsigned int layerdisk) {
444  unsigned int nrinv = NRINVBITS;
445  double rinvhalf = 0.5 * ((1 << nrinv) - 1);
446 
447  bool barrel = layerdisk < N_LAYER;
448  bool isPSmodule = layerdisk < N_PSLAYER;
449  double stripPitch = settings_.stripPitch(isPSmodule);
450 
451  if (barrel) {
452  unsigned int nbits = isPSmodule ? N_BENDBITS_PS : N_BENDBITS_2S;
453 
454  for (unsigned int irinv = 0; irinv < (1u << nrinv); irinv++) {
455  double rinv = (irinv - rinvhalf) * (1 << (settings_.nbitsrinv() - nrinv)) * settings_.krinvpars();
456 
457  double projbend = bendstrip(settings_.rmean(layerdisk), rinv, stripPitch);
458  for (unsigned int ibend = 0; ibend < (1u << nbits); ibend++) {
459  double stubbend = settings_.benddecode(ibend, layerdisk, isPSmodule);
460  bool pass = std::abs(stubbend - projbend) < settings_.bendcutme(ibend, layerdisk, isPSmodule);
461  table_.push_back(pass);
462  }
463  }
464  } else {
465  for (unsigned int iprojbend = 0; iprojbend < (1u << nrinv); iprojbend++) {
466  double projbend = 0.5 * (iprojbend - rinvhalf);
467  for (unsigned int ibend = 0; ibend < (1 << N_BENDBITS_2S); ibend++) {
468  double stubbend = settings_.benddecode(ibend, layerdisk, false);
469  bool pass = std::abs(stubbend - projbend) < settings_.bendcutme(ibend, layerdisk, false);
470  table_.push_back(pass);
471  }
472  }
473  for (unsigned int iprojbend = 0; iprojbend < (1u << nrinv); iprojbend++) {
474  double projbend = 0.5 * (iprojbend - rinvhalf);
475  for (unsigned int ibend = 0; ibend < (1 << N_BENDBITS_PS); ibend++) {
476  double stubbend = settings_.benddecode(ibend, layerdisk, true);
477  bool pass = std::abs(stubbend - projbend) < settings_.bendcutme(ibend, layerdisk, true);
478  table_.push_back(pass);
479  }
480  }
481  }
482 
483  positive_ = false;
484 
485  name_ = "METable_" + TrackletConfigBuilder::LayerName(layerdisk) + ".tab";
486 
487  writeTable();
488 }
489 
490 void TrackletLUT::initVMRTable(unsigned int layerdisk, VMRTableType type, int region) {
491  unsigned int zbits = settings_.vmrlutzbits(layerdisk);
492  unsigned int rbits = settings_.vmrlutrbits(layerdisk);
493 
494  unsigned int rbins = (1 << rbits);
495  unsigned int zbins = (1 << zbits);
496 
497  double zmin, zmax, rmin, rmax;
498 
499  if (layerdisk < N_LAYER) {
500  zmin = -settings_.zlength();
501  zmax = settings_.zlength();
502  rmin = settings_.rmean(layerdisk) - settings_.drmax();
503  rmax = settings_.rmean(layerdisk) + settings_.drmax();
504  } else {
505  rmin = 0;
506  rmax = settings_.rmaxdisk();
507  zmin = settings_.zmean(layerdisk - N_LAYER) - settings_.dzmax();
508  zmax = settings_.zmean(layerdisk - N_LAYER) + settings_.dzmax();
509  }
510 
511  double dr = (rmax - rmin) / rbins;
512  double dz = (zmax - zmin) / zbins;
513 
515 
516  for (unsigned int izbin = 0; izbin < zbins; izbin++) {
517  for (unsigned int irbin = 0; irbin < rbins; irbin++) {
518  double r = rmin + (irbin + 0.5) * dr;
519  double z = zmin + (izbin + 0.5) * dz;
520 
521  if (settings_.combined()) {
522  int iznew = izbin - (1 << (zbits - 1));
523  if (iznew < 0)
524  iznew += (1 << zbits);
525  assert(iznew >= 0);
526  assert(iznew < (1 << zbits));
527  z = zmin + (iznew + 0.5) * dz;
528  if (layerdisk < N_LAYER) {
529  int irnew = irbin - (1 << (rbits - 1));
530  if (irnew < 0)
531  irnew += (1 << rbits);
532  assert(irnew >= 0);
533  assert(irnew < (1 << rbits));
534  r = rmin + (irnew + 0.5) * dr;
535  }
536  }
537 
538  if (layerdisk >= N_LAYER && irbin < 10) //special case for the tabulated radii in 2S disks
539  r = (layerdisk < N_LAYER + 2) ? settings_.rDSSinner(irbin) : settings_.rDSSouter(irbin);
540 
541  int bin;
542  if (layerdisk < N_LAYER) {
543  double zproj = z * settings_.rmean(layerdisk) / r;
544  bin = NBINS * (zproj + settings_.zlength()) / (2 * settings_.zlength());
545  } else {
546  double rproj = r * settings_.zmean(layerdisk - N_LAYER) / z;
547  bin = NBINS * (rproj - settings_.rmindiskvm()) / (settings_.rmaxdisk() - settings_.rmindiskvm());
548  }
549  if (bin < 0)
550  bin = 0;
551  if (bin >= NBINS)
552  bin = NBINS - 1;
553 
554  if (type == VMRTableType::me) {
555  table_.push_back(bin);
556  }
557 
558  if (type == VMRTableType::disk) {
559  if (layerdisk >= N_LAYER) {
560  double rproj = r * settings_.zmean(layerdisk - N_LAYER) / z;
561  bin = 0.5 * NBINS * (rproj - settings_.rmindiskvm()) / (settings_.rmaxdiskvm() - settings_.rmindiskvm());
562  //bin value of zero indicates that stub is out of range
563  if (bin < 0)
564  bin = 0;
565  if (bin >= NBINS / 2)
566  bin = 0;
567  table_.push_back(bin);
568  }
569  }
570 
571  if (type == VMRTableType::inner) {
572  if (layerdisk == LayerDisk::L1 || layerdisk == LayerDisk::L3 || layerdisk == LayerDisk::L5 ||
573  layerdisk == LayerDisk::D1 || layerdisk == LayerDisk::D3) {
574  table_.push_back(getVMRLookup(layerdisk + 1, z, r, dz, dr));
575  }
576  if (layerdisk == LayerDisk::L2) {
577  table_.push_back(getVMRLookup(layerdisk + 1, z, r, dz, dr, Seed::L2L3));
578  }
579  }
580 
581  if (type == VMRTableType::inneroverlap) {
582  if (layerdisk == LayerDisk::L1 || layerdisk == LayerDisk::L2) {
583  table_.push_back(getVMRLookup(6, z, r, dz, dr, layerdisk + 6));
584  }
585  }
586 
587  if (type == VMRTableType::innerthird) {
588  if (layerdisk == LayerDisk::L2) { //projection from L2 to D1 for L2L3D1 seeding
589  table_.push_back(getVMRLookup(LayerDisk::D1, z, r, dz, dr, Seed::L2L3D1));
590  }
591 
592  if (layerdisk == LayerDisk::L5) { //projection from L5 to L4 for L5L6L4 seeding
593  table_.push_back(getVMRLookup(LayerDisk::L4, z, r, dz, dr));
594  }
595 
596  if (layerdisk == LayerDisk::L3) { //projection from L3 to L5 for L3L4L2 seeding
597  table_.push_back(getVMRLookup(LayerDisk::L2, z, r, dz, dr));
598  }
599 
600  if (layerdisk == LayerDisk::D1) { //projection from D1 to L2 for D1D2L2 seeding
601  table_.push_back(getVMRLookup(LayerDisk::L2, z, r, dz, dr));
602  }
603  }
604  }
605  }
606 
607  if (settings_.combined()) {
608  if (type == VMRTableType::me) {
609  positive_ = false;
610  name_ = "VMRME_" + TrackletConfigBuilder::LayerName(layerdisk) + ".tab";
611  }
612  if (type == VMRTableType::disk) {
613  positive_ = false;
614  name_ = "VMRTE_" + TrackletConfigBuilder::LayerName(layerdisk) + ".tab";
615  }
616  if (type == VMRTableType::inner) {
617  positive_ = true;
618  nbits_ = 10;
619  name_ = "TP_" + TrackletConfigBuilder::LayerName(layerdisk) + TrackletConfigBuilder::LayerName(layerdisk + 1) +
620  ".tab";
621  }
622 
623  if (type == VMRTableType::inneroverlap) {
624  positive_ = true;
625  nbits_ = 10;
627  }
628 
629  } else {
630  if (type == VMRTableType::me) {
631  //This if a hack where the same memory is used in both ME and TE modules
632  if (layerdisk == 1 || layerdisk == 2 || layerdisk == 3 || layerdisk == 4) {
633  positive_ = false;
634  name_ = "VMTableOuter" + TrackletConfigBuilder::LayerName(layerdisk) + ".tab";
635  writeTable();
636  }
637 
638  assert(region >= 0);
639  char cregion = 'A' + region;
640  name_ = "VMR_" + TrackletConfigBuilder::LayerName(layerdisk) + "PHI" + cregion + "_finebin.tab";
641  positive_ = false;
642  }
643 
644  if (type == VMRTableType::inner) {
645  positive_ = false;
646  name_ = "VMTableInner" + TrackletConfigBuilder::LayerName(layerdisk) +
647  TrackletConfigBuilder::LayerName(layerdisk + 1) + ".tab";
648  }
649 
650  if (type == VMRTableType::inneroverlap) {
651  positive_ = false;
653  ".tab";
654  }
655 
656  if (type == VMRTableType::disk) {
657  positive_ = false;
658  name_ = "VMTableOuter" + TrackletConfigBuilder::LayerName(layerdisk) + ".tab";
659  }
660  }
661 
662  writeTable();
663 }
664 
665 int TrackletLUT::getVMRLookup(unsigned int layerdisk, double z, double r, double dz, double dr, int iseed) const {
666  double z0cut = settings_.z0cut();
667 
668  if (layerdisk < N_LAYER) {
669  if (iseed == Seed::L2L3 && std::abs(z) < 52.0)
670  return -1;
671 
672  double rmean = settings_.rmean(layerdisk);
673 
674  double rratio1 = rmean / (r + 0.5 * dr);
675  double rratio2 = rmean / (r - 0.5 * dr);
676 
677  double z1 = (z - 0.5 * dz) * rratio1 + z0cut * (rratio1 - 1.0);
678  double z2 = (z + 0.5 * dz) * rratio1 + z0cut * (rratio1 - 1.0);
679  double z3 = (z - 0.5 * dz) * rratio2 + z0cut * (rratio2 - 1.0);
680  double z4 = (z + 0.5 * dz) * rratio2 + z0cut * (rratio2 - 1.0);
681  double z5 = (z - 0.5 * dz) * rratio1 - z0cut * (rratio1 - 1.0);
682  double z6 = (z + 0.5 * dz) * rratio1 - z0cut * (rratio1 - 1.0);
683  double z7 = (z - 0.5 * dz) * rratio2 - z0cut * (rratio2 - 1.0);
684  double z8 = (z + 0.5 * dz) * rratio2 - z0cut * (rratio2 - 1.0);
685 
686  double zmin = std::min({z1, z2, z3, z4, z5, z6, z7, z8});
687  double zmax = std::max({z1, z2, z3, z4, z5, z6, z7, z8});
688 
690 
691  int zbin1 = NBINS * (zmin + settings_.zlength()) / (2 * settings_.zlength());
692  int zbin2 = NBINS * (zmax + settings_.zlength()) / (2 * settings_.zlength());
693 
694  if (zbin1 >= NBINS)
695  return -1;
696  if (zbin2 < 0)
697  return -1;
698 
699  if (zbin2 >= NBINS)
700  zbin2 = NBINS - 1;
701  if (zbin1 < 0)
702  zbin1 = 0;
703 
704  // This is a 10 bit word:
705  // xxx|yyy|z|rrr
706  // xxx is the delta z window
707  // yyy is the z bin
708  // z is flag to look in next bin
709  // rrr first fine z bin
710  // NOTE : this encoding is not efficient z is one if xxx+rrr is greater than 8
711  // and xxx is only 1,2, or 3
712  // should also reject xxx=0 as this means projection is outside range
713 
714  int value = zbin1 / 8;
715  value *= 2;
716  if (zbin2 / 8 - zbin1 / 8 > 0)
717  value += 1;
718  value *= 8;
719  value += (zbin1 & 7);
720  assert(value / 8 < 15);
721  int deltaz = zbin2 - zbin1;
722  if (deltaz > 7) {
723  deltaz = 7;
724  }
725  assert(deltaz < 8);
726  value += (deltaz << 7);
727 
728  return value;
729 
730  } else {
731  if (std::abs(z) < 2.0 * z0cut)
732  return -1;
733 
734  double zmean = settings_.zmean(layerdisk - N_LAYER);
735  if (z < 0.0)
736  zmean = -zmean;
737 
738  double r1 = (r + 0.5 * dr) * (zmean + z0cut) / (z + 0.5 * dz + z0cut);
739  double r2 = (r - 0.5 * dr) * (zmean - z0cut) / (z + 0.5 * dz - z0cut);
740  double r3 = (r + 0.5 * dr) * (zmean + z0cut) / (z - 0.5 * dz + z0cut);
741  double r4 = (r - 0.5 * dr) * (zmean - z0cut) / (z - 0.5 * dz - z0cut);
742  double r5 = (r + 0.5 * dr) * (zmean - z0cut) / (z + 0.5 * dz - z0cut);
743  double r6 = (r - 0.5 * dr) * (zmean + z0cut) / (z + 0.5 * dz + z0cut);
744  double r7 = (r + 0.5 * dr) * (zmean - z0cut) / (z - 0.5 * dz - z0cut);
745  double r8 = (r - 0.5 * dr) * (zmean + z0cut) / (z - 0.5 * dz + z0cut);
746 
747  double rmin = std::min({r1, r2, r3, r4, r5, r6, r7, r8});
748  double rmax = std::max({r1, r2, r3, r4, r5, r6, r7, r8});
749 
751 
752  double rmindisk = settings_.rmindiskvm();
753  double rmaxdisk = settings_.rmaxdiskvm();
754 
755  if (iseed == Seed::L1D1)
756  rmaxdisk = settings_.rmaxdiskl1overlapvm();
757  if (iseed == Seed::L2D1)
758  rmindisk = settings_.rmindiskl2overlapvm();
759  if (iseed == Seed::L2L3D1)
760  rmaxdisk = settings_.rmaxdisk();
761 
762  if (rmin > rmaxdisk)
763  return -1;
764  if (rmax > rmaxdisk)
765  rmax = rmaxdisk;
766 
767  if (rmax < rmindisk)
768  return -1;
769  if (rmin < rmindisk)
770  rmin = rmindisk;
771 
772  int rbin1 = NBINS * (rmin - settings_.rmindiskvm()) / (settings_.rmaxdiskvm() - settings_.rmindiskvm());
773  int rbin2 = NBINS * (rmax - settings_.rmindiskvm()) / (settings_.rmaxdiskvm() - settings_.rmindiskvm());
774 
775  if (iseed == Seed::L2L3D1) {
776  constexpr double rminspec = 40.0;
777  rbin1 = NBINS * (rmin - rminspec) / (settings_.rmaxdisk() - rminspec);
778  rbin2 = NBINS * (rmax - rminspec) / (settings_.rmaxdisk() - rminspec);
779  }
780 
781  if (rbin2 >= NBINS)
782  rbin2 = NBINS - 1;
783  if (rbin1 < 0)
784  rbin1 = 0;
785 
786  // This is a 9 bit word:
787  // xxx|yy|z|rrr
788  // xxx is the delta r window
789  // yy is the r bin yy is three bits for overlaps
790  // z is flag to look in next bin
791  // rrr fine r bin
792  // NOTE : this encoding is not efficient z is one if xxx+rrr is greater than 8
793  // and xxx is only 1,2, or 3
794  // should also reject xxx=0 as this means projection is outside range
795 
796  bool overlap = iseed == Seed::L1D1 || iseed == Seed::L2D1 || iseed == Seed::L2L3D1;
797 
798  int value = rbin1 / 8;
799  if (overlap) {
800  if (z < 0.0)
801  value += 4;
802  }
803  value *= 2;
804  if (rbin2 / 8 - rbin1 / 8 > 0)
805  value += 1;
806  value *= 8;
807  value += (rbin1 & 7);
808  assert(value / 8 < 15);
809  int deltar = rbin2 - rbin1;
810  if (deltar > 7)
811  deltar = 7;
812  if (overlap) {
813  value += (deltar << 7);
814  } else {
815  value += (deltar << 6);
816  }
817 
818  return value;
819  }
820 }
821 
822 void TrackletLUT::initPhiCorrTable(unsigned int layerdisk, unsigned int rbits) {
823  bool psmodule = layerdisk < N_PSLAYER;
824 
825  unsigned int bendbits = psmodule ? N_BENDBITS_PS : N_BENDBITS_2S;
826 
827  unsigned int rbins = (1 << rbits);
828 
829  double rmean = settings_.rmean(layerdisk);
830  double drmax = settings_.drmax();
831 
832  double dr = 2.0 * drmax / rbins;
833 
834  unsigned int bendbins = (1 << bendbits);
835 
836  for (unsigned int ibend = 0; ibend < bendbins; ibend++) {
837  for (unsigned int irbin = 0; irbin < rbins; irbin++) {
838  int value = getphiCorrValue(layerdisk, ibend, irbin, rmean, dr, drmax);
839  table_.push_back(value);
840  }
841  }
842 
843  name_ = "VMPhiCorrL" + std::to_string(layerdisk + 1) + ".tab";
844  nbits_ = 14;
845  positive_ = false;
846 
847  writeTable();
848 }
849 
851  unsigned int layerdisk, unsigned int ibend, unsigned int irbin, double rmean, double dr, double drmax) const {
852  bool psmodule = layerdisk < N_PSLAYER;
853 
854  double bend = -settings_.benddecode(ibend, layerdisk, psmodule);
855 
856  //for the rbin - calculate the distance to the nominal layer radius
857  double Delta = (irbin + 0.5) * dr - drmax;
858 
859  //calculate the phi correction - this is a somewhat approximate formula
860  double dphi = (Delta / 0.18) * bend * settings_.stripPitch(psmodule) / rmean;
861 
862  double kphi = psmodule ? settings_.kphi() : settings_.kphi1();
863 
864  int idphi = dphi / kphi;
865 
866  return idphi;
867 }
868 
869 // Write LUT table.
871  if (!settings_.writeTable()) {
872  return;
873  }
874 
875  if (name_.empty()) {
876  return;
877  }
878 
879  ofstream out = openfile(settings_.tablePath(), name_, __FILE__, __LINE__);
880 
881  out << "{" << endl;
882  for (unsigned int i = 0; i < table_.size(); i++) {
883  if (i != 0) {
884  out << "," << endl;
885  }
886 
887  int itable = table_[i];
888  if (positive_) {
889  if (table_[i] < 0) {
890  itable = (1 << nbits_) - 1;
891  }
892  }
893 
894  out << itable;
895  }
896  out << endl << "};" << endl;
897  out.close();
898 }
899 
900 int TrackletLUT::lookup(unsigned int index) const {
901  assert(index < table_.size());
902  return table_[index];
903 }
void initBendMatch(unsigned int layerdisk)
Definition: TrackletLUT.cc:443
double dphisectorHG() const
Definition: Settings.h:281
double rmindiskl2overlapvm() const
Definition: Settings.h:319
double z0cut() const
Definition: Settings.h:324
unsigned int vmrlutrbits(unsigned int layerdisk) const
Definition: Settings.h:179
unsigned int vmrlutzbits(unsigned int layerdisk) const
Definition: Settings.h:178
void initmatchcut(unsigned int layerdisk, MatchType type, unsigned int region)
Definition: TrackletLUT.cc:13
constexpr unsigned int NRINVBITS
Definition: Settings.h:32
void initProjectionBend(double k_phider, unsigned int idisk, unsigned int nrbits, unsigned int nphiderbits)
Definition: TrackletLUT.cc:393
void initteptlut(bool fillInner, bool fillTEMem, unsigned int iSeed, unsigned int layerdisk1, unsigned int layerdisk2, unsigned int innerphibits, unsigned int outerphibits, double innerphimin, double innerphimax, double outerphimin, double outerphimax, const std::string &innermem, const std::string &outermem)
Definition: TrackletLUT.cc:240
double zmatchcut(unsigned int iSeed, unsigned int ilayer) const
Definition: Settings.h:158
int getVMRLookup(unsigned int layerdisk, double z, double r, double dz, double dr, int iseed=-1) const
Definition: TrackletLUT.cc:665
const Settings & settings_
Definition: TrackletLUT.h:85
double bendcutme(int ibend, int layerdisk, bool isPSmodule) const
Definition: Settings.h:428
double rmaxdisk() const
Definition: Settings.h:125
double zlength() const
Definition: Settings.h:124
double kphi() const
Definition: Settings.h:298
static std::string LayerName(unsigned int ilayer)
constexpr unsigned int N_BENDBITS_2S
Definition: Settings.h:30
double rcrit() const
Definition: Settings.h:288
assert(be >=bs)
double bendcutte(int ibend, int layerdisk, bool isPSmodule) const
Definition: Settings.h:426
constexpr unsigned int N_BENDBITS_PS
Definition: Settings.h:29
void initPhiCorrTable(unsigned int layerdisk, unsigned int rbits)
Definition: TrackletLUT.cc:822
int getphiCorrValue(unsigned int layerdisk, unsigned int ibend, unsigned int irbin, double rmean, double dr, double drmax) const
Definition: TrackletLUT.cc:850
bool overlap(const reco::Muon &muon1, const reco::Muon &muon2, double pullX=1.0, double pullY=1.0, bool checkAdjacentChambers=false)
double rmean(unsigned int iLayer) const
Definition: Settings.h:164
double krinvpars() const
Definition: Settings.h:384
double zmean(unsigned int iDisk) const
Definition: Settings.h:167
std::string name_
Definition: TrackletLUT.h:87
double rmaxdiskl1overlapvm() const
Definition: Settings.h:318
unsigned int NLONGVMBINS() const
Definition: Settings.h:329
Divides< A, C > D1
Definition: Factorize.h:136
unsigned int nvmte(unsigned int inner, unsigned int iSeed) const
Definition: Settings.h:101
unsigned int nrbitsstub(unsigned int layerdisk) const
Definition: Settings.h:84
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int nallstubs(unsigned int layerdisk) const
Definition: Settings.h:107
double rphicutPS(unsigned int iSeed, unsigned int idisk) const
Definition: Settings.h:159
T min(T a, T b)
Definition: MathUtil.h:58
double rDSSouter(unsigned int iBin) const
Definition: Settings.h:174
double kphi1() const
Definition: Settings.h:299
double dzmax() const
Definition: Settings.h:129
int lookup(unsigned int index) const
Definition: TrackletLUT.cc:900
constexpr unsigned int N_PSLAYER
Definition: Settings.h:23
unsigned int nbitsphiprojderL123() const
Definition: Settings.h:87
int iseed
Definition: AMPTWrapper.h:134
double rcut2S(unsigned int iSeed, unsigned int idisk) const
Definition: Settings.h:162
double rDSSinner(unsigned int iBin) const
Definition: Settings.h:171
double rinv(double phi1, double phi2, double r1, double r2)
Definition: Util.h:49
void initTPlut(bool fillInner, unsigned int iSeed, unsigned int layerdisk1, unsigned int layerdisk2, unsigned int nbitsfinephidiff, unsigned int iTP)
Definition: TrackletLUT.cc:63
unsigned int nbits_
Definition: TrackletLUT.h:91
unsigned int nbitsallstubs(unsigned int layerdisk) const
Definition: Settings.h:106
void initTPregionlut(unsigned int iSeed, unsigned int layerdisk1, unsigned int layerdisk2, unsigned int iAllStub, unsigned int nbitsfinephidiff, unsigned int nbitsfinephi, const TrackletLUT &tplutinner, unsigned int iTP)
Definition: TrackletLUT.cc:183
double rinvcutte() const
Definition: Settings.h:313
double rphimatchcut(unsigned int iSeed, unsigned int ilayer) const
Definition: Settings.h:157
double rmindiskvm() const
Definition: Settings.h:315
double kz() const
Definition: Settings.h:302
int nfinephi(unsigned int inner, unsigned int iSeed) const
Definition: Settings.h:133
std::string tablePath() const
Definition: Settings.h:193
double drmax() const
Definition: Settings.h:128
double rphicut2S(unsigned int iSeed, unsigned int idisk) const
Definition: Settings.h:161
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
double kr() const
Definition: Settings.h:304
std::vector< int > table_
Definition: TrackletLUT.h:89
double krprojshiftdisk() const
Definition: Settings.h:400
int nbitsrinv() const
Definition: Settings.h:334
double benddecode(int ibend, int layerdisk, bool isPSmodule) const
Definition: Settings.h:402
double stripPitch(bool isPSmodule) const
Definition: Settings.h:260
double bendstrip(double r, double rinv, double stripPitch)
Definition: Util.h:42
double rcutPS(unsigned int iSeed, unsigned int idisk) const
Definition: Settings.h:160
double rmaxdiskvm() const
Definition: Settings.h:316
bool writeTable() const
Definition: Settings.h:189
void initVMRTable(unsigned int layerdisk, VMRTableType type, int region=-1)
Definition: TrackletLUT.cc:490
const int NBINS
bool combined() const
Definition: Settings.h:250
std::ofstream openfile(const std::string &dir, const std::string &fname, const char *file, int line)
Definition: Util.h:139
constexpr int N_LAYER
Definition: Settings.h:21
void writeTable() const
Definition: TrackletLUT.cc:870