CMS 3D CMS Logo

TrackletCalculatorBase.cc
Go to the documentation of this file.
9 
13 
14 using namespace std;
15 using namespace trklet;
16 
17 TrackletCalculatorBase::TrackletCalculatorBase(string name,
18  Settings const& settings,
19  Globals* global,
20  unsigned int iSector)
21  : ProcessBase(name, settings, global, iSector) {}
22 
24  double z1,
25  double phi1,
26  double r2,
27  double z2,
28  double phi2,
29  double,
30  double& rinv,
31  double& phi0,
32  double& t,
33  double& z0,
34  double phiproj[N_LAYER - 2],
35  double zproj[N_LAYER - 2],
36  double phider[N_LAYER - 2],
37  double zder[N_LAYER - 2],
38  double phiprojdisk[N_DISK],
39  double rprojdisk[N_DISK],
40  double phiderdisk[N_DISK],
41  double rderdisk[N_DISK]) {
42  double deltaphi = reco::reduceRange(phi1 - phi2);
43 
44  double dist = sqrt(r2 * r2 + r1 * r1 - 2 * r1 * r2 * cos(deltaphi));
45 
46  rinv = 2 * sin(deltaphi) / dist;
47 
48  double phi1tmp = phi1 - phimin_;
49 
50  phi0 = reco::reduceRange(phi1tmp + asin(0.5 * r1 * rinv));
51 
52  double rhopsi1 = 2 * asin(0.5 * r1 * rinv) / rinv;
53  double rhopsi2 = 2 * asin(0.5 * r2 * rinv) / rinv;
54 
55  t = (z1 - z2) / (rhopsi1 - rhopsi2);
56 
57  z0 = z1 - t * rhopsi1;
58 
59  for (unsigned int i = 0; i < N_LAYER - 2; i++) {
61  rinv,
62  phi0,
63  t,
64  z0,
65  phiproj[i],
66  zproj[i],
67  phider[i],
68  zder[i]);
69  }
70 
71  for (unsigned int i = 0; i < N_DISK; i++) {
72  exactprojdisk(settings_.zmean(i), rinv, phi0, t, z0, phiprojdisk[i], rprojdisk[i], phiderdisk[i], rderdisk[i]);
73  }
74 }
75 
77  double z1,
78  double phi1,
79  double r2,
80  double z2,
81  double phi2,
82  double,
83  double& rinv,
84  double& phi0,
85  double& t,
86  double& z0,
87  double phiprojLayer[N_PSLAYER], //=3 (project to PS barrel layers only)
88  double zprojLayer[N_PSLAYER],
89  double phiderLayer[N_PSLAYER],
90  double zderLayer[N_PSLAYER],
91  double phiproj[N_DISK - 2], //=3 (max project to 3 other disks)
92  double rproj[N_DISK - 2],
93  double phider[N_DISK - 2],
94  double rder[N_DISK - 2]) {
95  double deltaphi = reco::reduceRange(phi1 - phi2);
96 
97  double dist = sqrt(r2 * r2 + r1 * r1 - 2 * r1 * r2 * cos(deltaphi));
98 
99  rinv = 2 * sin(deltaphi) / dist;
100 
101  double phi1tmp = phi1 - phimin_;
102 
103  phi0 = reco::reduceRange(phi1tmp + asin(0.5 * r1 * rinv));
104 
105  double rhopsi1 = 2 * asin(0.5 * r1 * rinv) / rinv;
106  double rhopsi2 = 2 * asin(0.5 * r2 * rinv) / rinv;
107 
108  t = (z1 - z2) / (rhopsi1 - rhopsi2);
109 
110  z0 = z1 - t * rhopsi1;
111 
112  for (unsigned int i = 0; i < N_DISK - 2; i++) {
114  rinv,
115  phi0,
116  t,
117  z0,
118  phiproj[i],
119  rproj[i],
120  phider[i],
121  rder[i]);
122  }
123 
124  for (unsigned int i = 0; i < N_DISK - 2; i++) {
125  exactproj(settings_.rmean(i), rinv, phi0, t, z0, phiprojLayer[i], zprojLayer[i], phiderLayer[i], zderLayer[i]);
126  }
127 }
128 
130  double z1,
131  double phi1,
132  double r2,
133  double z2,
134  double phi2,
135  double,
136  double& rinv,
137  double& phi0,
138  double& t,
139  double& z0,
140  double phiprojLayer[N_PSLAYER],
141  double zprojLayer[N_PSLAYER],
142  double phiderLayer[N_PSLAYER],
143  double zderLayer[N_PSLAYER],
144  double phiproj[N_DISK - 2],
145  double rproj[N_DISK - 2],
146  double phider[N_DISK - 2],
147  double rder[N_DISK - 2]) {
148  double deltaphi = reco::reduceRange(phi1 - phi2);
149 
150  double dist = sqrt(r2 * r2 + r1 * r1 - 2 * r1 * r2 * cos(deltaphi));
151 
152  rinv = 2 * sin(deltaphi) / dist;
153 
154  if (r1 > r2)
155  rinv = -rinv;
156 
157  double phi1tmp = phi1 - phimin_;
158 
159  phi0 = reco::reduceRange(phi1tmp + asin(0.5 * r1 * rinv));
160 
161  double rhopsi1 = 2 * asin(0.5 * r1 * rinv) / rinv;
162  double rhopsi2 = 2 * asin(0.5 * r2 * rinv) / rinv;
163 
164  t = (z1 - z2) / (rhopsi1 - rhopsi2);
165 
166  z0 = z1 - t * rhopsi1;
167 
168  for (int i = 0; i < 4; i++) {
169  exactprojdisk(settings_.zmean(i + 1), rinv, phi0, t, z0, phiproj[i], rproj[i], phider[i], rder[i]);
170  }
171 
172  for (int i = 0; i < 1; i++) {
173  exactproj(settings_.rmean(i), rinv, phi0, t, z0, phiprojLayer[i], zprojLayer[i], phiderLayer[i], zderLayer[i]);
174  }
175 }
176 
178  double rinv,
179  double phi0,
180  double t,
181  double z0,
182  double& phiproj,
183  double& zproj,
184  double& phider,
185  double& zder) {
186  phiproj = phi0 - asin(0.5 * rproj * rinv);
187  zproj = z0 + (2 * t / rinv) * asin(0.5 * rproj * rinv);
188 
189  phider = -0.5 * rinv / sqrt(1 - pow(0.5 * rproj * rinv, 2));
190  zder = t / sqrt(1 - pow(0.5 * rproj * rinv, 2));
191 }
192 
194  double rinv,
195  double phi0,
196  double t,
197  double z0,
198  double& phiproj,
199  double& rproj,
200  double& phider,
201  double& rder) {
202  if (t < 0)
203  zproj = -zproj;
204 
205  double tmp = rinv * (zproj - z0) / (2.0 * t);
206  rproj = (2.0 / rinv) * sin(tmp);
207  phiproj = phi0 - tmp;
208 
209  phider = -rinv / (2 * t);
210  rder = cos(tmp) / t;
211 }
212 
214  FPGAWord fpgar = tracklet->fpgarprojdisk(disk);
215 
216  if (fpgar.value() * settings_.krprojshiftdisk() < settings_.rmindiskvm())
217  return;
218  if (fpgar.value() * settings_.krprojshiftdisk() > settings_.rmaxdisk())
219  return;
220 
221  FPGAWord fpgaphi = tracklet->fpgaphiprojdisk(disk);
222 
223  int iphivmRaw = fpgaphi.value() >> (fpgaphi.nbits() - 5);
224 
225  int iphi = iphivmRaw / (32 / settings_.nallstubs(abs(disk) + N_DISK));
226 
227  addProjectionDisk(disk, iphi, trackletprojdisks_[abs(disk) - 1][iphi], tracklet);
228 }
229 
230 bool TrackletCalculatorBase::addLayerProj(Tracklet* tracklet, int layer) {
231  assert(layer > 0);
232 
233  FPGAWord fpgaz = tracklet->fpgazproj(layer);
234  FPGAWord fpgaphi = tracklet->fpgaphiproj(layer);
235 
236  if (fpgaphi.atExtreme())
237  edm::LogProblem("Tracklet") << "at extreme! " << fpgaphi.value();
238 
239  assert(!fpgaphi.atExtreme());
240 
241  if (fpgaz.atExtreme())
242  return false;
243 
244  if (std::abs(fpgaz.value() * settings_.kz()) > settings_.zlength())
245  return false;
246 
247  int iphivmRaw = fpgaphi.value() >> (fpgaphi.nbits() - 5);
248  int iphi = iphivmRaw / (32 / settings_.nallstubs(layer - 1));
249 
250  addProjection(layer, iphi, trackletprojlayers_[layer - 1][iphi], tracklet);
251 
252  return true;
253 }
254 
256  int iphi,
257  TrackletProjectionsMemory* trackletprojs,
258  Tracklet* tracklet) {
259  if (trackletprojs == nullptr) {
260  if (settings_.warnNoMem()) {
261  edm::LogVerbatim("Tracklet") << "No projection memory exists in " << getName() << " for layer = " << layer
262  << " iphi = " << iphi + 1;
263  }
264  return;
265  }
266  assert(trackletprojs != nullptr);
267  trackletprojs->addProj(tracklet);
268 }
269 
271  int iphi,
272  TrackletProjectionsMemory* trackletprojs,
273  Tracklet* tracklet) {
274  if (iSeed_ == 2 && abs(disk) == 4)
275  return; //L3L4 projections to D3 are not used. Should be in configuration
276  if (trackletprojs == nullptr) {
277  if (iSeed_ == 2 && abs(disk) == 3)
278  return; //L3L4 projections to D3 are not used.
279  if (settings_.warnNoMem()) {
280  edm::LogVerbatim("Tracklet") << "No projection memory exists in " << getName() << " for disk = " << abs(disk)
281  << " iphi = " << iphi + 1;
282  }
283  return;
284  }
285  assert(trackletprojs != nullptr);
286  trackletprojs->addProj(tracklet);
287 }
288 
289 bool TrackletCalculatorBase::goodTrackPars(bool goodrinv, bool goodz0) {
290  bool success = true;
291  if (!goodrinv) {
292  if (settings_.debugTracklet()) {
293  edm::LogVerbatim("Tracklet") << getName() << " TrackletCalculatorBase irinv too large";
294  }
295  success = false;
296  }
297  if (!goodz0) {
298  if (settings_.debugTracklet()) {
299  edm::LogVerbatim("Tracklet") << getName() << " TrackletCalculatorBase z0 cut to large";
300  }
301  success = false;
302  }
303  return success;
304 }
305 
306 bool TrackletCalculatorBase::inSector(int iphi0, int irinv, double phi0approx, double rinvapprox) {
307  double phicritapprox = phi0approx - asin(0.5 * settings_.rcrit() * rinvapprox);
308 
309  int ifactor = 0.5 * settings_.rcrit() * settings_.krinvpars() / settings_.kphi0pars() * (1 << 8);
310  int iphicrit = iphi0 - (irinv >> 8) * ifactor;
311 
312  int iphicritmincut = settings_.phicritminmc() / globals_->ITC_L1L2()->phi0_final.K();
313  int iphicritmaxcut = settings_.phicritmaxmc() / globals_->ITC_L1L2()->phi0_final.K();
314 
315  bool keepapprox = (phicritapprox > settings_.phicritminmc()) && (phicritapprox < settings_.phicritmaxmc()),
316  keep = (iphicrit > iphicritmincut) && (iphicrit < iphicritmaxcut);
317  if (settings_.debugTracklet())
318  if (keepapprox && !keep)
319  edm::LogVerbatim("Tracklet") << getName()
320  << " Tracklet kept with exact phicrit cut but not approximate, phicritapprox: "
321  << phicritapprox;
322  if (settings_.usephicritapprox()) {
323  return keepapprox;
324  } else {
325  return keep;
326  }
327 
328  return true;
329 }
330 
331 bool TrackletCalculatorBase::barrelSeeding(const Stub* innerFPGAStub,
332  const L1TStub* innerStub,
333  const Stub* outerFPGAStub,
334  const L1TStub* outerStub) {
335  if (settings_.debugTracklet()) {
336  edm::LogVerbatim("Tracklet") << "TrackletCalculator " << getName()
337  << " trying stub pair in layer (inner outer): " << innerFPGAStub->layer().value()
338  << " " << outerFPGAStub->layer().value();
339  }
340 
341  assert(outerFPGAStub->isBarrel());
342  assert(layerdisk1_ == (unsigned int)innerFPGAStub->layer().value());
344 
345  double r1 = innerStub->r();
346  double z1 = innerStub->z();
347  double phi1 = innerStub->phi();
348 
349  double r2 = outerStub->r();
350  double z2 = outerStub->z();
351  double phi2 = outerStub->phi();
352 
353  double rinv, phi0, t, z0;
354 
355  double phiproj[N_LAYER - 2], zproj[N_LAYER - 2], phider[N_LAYER - 2], zder[N_LAYER - 2];
356  double phiprojdisk[N_DISK], rprojdisk[N_DISK], phiderdisk[N_DISK], rderdisk[N_DISK];
357 
359  z1,
360  phi1,
361  r2,
362  z2,
363  phi2,
364  outerStub->sigmaz(),
365  rinv,
366  phi0,
367  t,
368  z0,
369  phiproj,
370  zproj,
371  phider,
372  zder,
373  phiprojdisk,
374  rprojdisk,
375  phiderdisk,
376  rderdisk);
377 
378  if (settings_.useapprox()) {
379  phi1 = innerFPGAStub->phiapprox(phimin_, phimax_);
380  z1 = innerFPGAStub->zapprox();
381  r1 = innerFPGAStub->rapprox();
382 
383  phi2 = outerFPGAStub->phiapprox(phimin_, phimax_);
384  z2 = outerFPGAStub->zapprox();
385  r2 = outerFPGAStub->rapprox();
386  }
387 
388  double rinvapprox, phi0approx, tapprox, z0approx;
389  double phiprojapprox[N_LAYER - 2], zprojapprox[N_LAYER - 2];
390  double phiprojdiskapprox[N_DISK], rprojdiskapprox[N_DISK];
391 
393  if (iSeed_ == 0)
394  ITC = globals_->ITC_L1L2();
395  else if (iSeed_ == 1)
396  ITC = globals_->ITC_L2L3();
397  else if (iSeed_ == 2)
398  ITC = globals_->ITC_L3L4();
399  else
400  ITC = globals_->ITC_L5L6();
401 
404  ITC->z1.set_fval(z1);
405  ITC->z2.set_fval(z2);
406  double sphi1 = angle0to2pi::make0To2pi(phi1 - phioffset_);
407  double sphi2 = angle0to2pi::make0To2pi(phi2 - phioffset_);
408 
409  ITC->phi1.set_fval(sphi1);
410  ITC->phi2.set_fval(sphi2);
411 
416 
417  ITC->zproj0.set_fval(t > 0 ? settings_.zmean(0) : -settings_.zmean(0));
418  ITC->zproj1.set_fval(t > 0 ? settings_.zmean(1) : -settings_.zmean(1));
419  ITC->zproj2.set_fval(t > 0 ? settings_.zmean(2) : -settings_.zmean(2));
420  ITC->zproj3.set_fval(t > 0 ? settings_.zmean(3) : -settings_.zmean(3));
421  ITC->zproj4.set_fval(t > 0 ? settings_.zmean(4) : -settings_.zmean(4));
422 
423  ITC->rinv_final.calculate();
424  ITC->phi0_final.calculate();
425  ITC->t_final.calculate();
426  ITC->z0_final.calculate();
427 
428  ITC->phiL_0_final.calculate();
429  ITC->phiL_1_final.calculate();
430  ITC->phiL_2_final.calculate();
431  ITC->phiL_3_final.calculate();
432 
433  ITC->zL_0_final.calculate();
434  ITC->zL_1_final.calculate();
435  ITC->zL_2_final.calculate();
436  ITC->zL_3_final.calculate();
437 
438  ITC->phiD_0_final.calculate();
439  ITC->phiD_1_final.calculate();
440  ITC->phiD_2_final.calculate();
441  ITC->phiD_3_final.calculate();
442  ITC->phiD_4_final.calculate();
443 
444  ITC->rD_0_final.calculate();
445  ITC->rD_1_final.calculate();
446  ITC->rD_2_final.calculate();
447  ITC->rD_3_final.calculate();
448  ITC->rD_4_final.calculate();
449 
450  ITC->der_phiL_final.calculate();
451  ITC->der_zL_final.calculate();
452  ITC->der_phiD_final.calculate();
453  ITC->der_rD_final.calculate();
454 
455  //store the approximate results
456  rinvapprox = ITC->rinv_final.fval();
457  phi0approx = ITC->phi0_final.fval();
458  tapprox = ITC->t_final.fval();
459  z0approx = ITC->z0_final.fval();
460 
461  phiprojapprox[0] = ITC->phiL_0_final.fval();
462  phiprojapprox[1] = ITC->phiL_1_final.fval();
463  phiprojapprox[2] = ITC->phiL_2_final.fval();
464  phiprojapprox[3] = ITC->phiL_3_final.fval();
465 
466  zprojapprox[0] = ITC->zL_0_final.fval();
467  zprojapprox[1] = ITC->zL_1_final.fval();
468  zprojapprox[2] = ITC->zL_2_final.fval();
469  zprojapprox[3] = ITC->zL_3_final.fval();
470 
471  phiprojdiskapprox[0] = ITC->phiD_0_final.fval();
472  phiprojdiskapprox[1] = ITC->phiD_1_final.fval();
473  phiprojdiskapprox[2] = ITC->phiD_2_final.fval();
474  phiprojdiskapprox[3] = ITC->phiD_3_final.fval();
475  phiprojdiskapprox[4] = ITC->phiD_4_final.fval();
476 
477  rprojdiskapprox[0] = ITC->rD_0_final.fval();
478  rprojdiskapprox[1] = ITC->rD_1_final.fval();
479  rprojdiskapprox[2] = ITC->rD_2_final.fval();
480  rprojdiskapprox[3] = ITC->rD_3_final.fval();
481  rprojdiskapprox[4] = ITC->rD_4_final.fval();
482 
483  //now binary
484 
485  int irinv, iphi0, it, iz0;
486  LayerProjection layerprojs[N_LAYER - 2];
487  DiskProjection diskprojs[N_DISK];
488 
489  int iphiproj[N_LAYER - 2], izproj[N_LAYER - 2];
490  int iphiprojdisk[N_DISK], irprojdisk[N_DISK];
491 
492  int ir1 = innerFPGAStub->r().value();
493  int iphi1 = innerFPGAStub->phi().value();
494  int iz1 = innerFPGAStub->z().value();
495 
496  int ir2 = outerFPGAStub->r().value();
497  int iphi2 = outerFPGAStub->phi().value();
498  int iz2 = outerFPGAStub->z().value();
499 
502  ir1 <<= (8 - settings_.nrbitsstub(layerdisk1_));
503  ir2 <<= (8 - settings_.nrbitsstub(layerdisk2_));
504 
507 
508  ITC->r1.set_ival(ir1);
509  ITC->r2.set_ival(ir2);
510  ITC->z1.set_ival(iz1);
511  ITC->z2.set_ival(iz2);
512  ITC->phi1.set_ival(iphi1);
513  ITC->phi2.set_ival(iphi2);
514 
515  ITC->rinv_final.calculate();
516  ITC->phi0_final.calculate();
517  ITC->t_final.calculate();
518  ITC->z0_final.calculate();
519 
520  ITC->phiL_0_final.calculate();
521  ITC->phiL_1_final.calculate();
522  ITC->phiL_2_final.calculate();
523  ITC->phiL_3_final.calculate();
524 
525  ITC->zL_0_final.calculate();
526  ITC->zL_1_final.calculate();
527  ITC->zL_2_final.calculate();
528  ITC->zL_3_final.calculate();
529 
530  ITC->phiD_0_final.calculate();
531  ITC->phiD_1_final.calculate();
532  ITC->phiD_2_final.calculate();
533  ITC->phiD_3_final.calculate();
534  ITC->phiD_4_final.calculate();
535 
536  ITC->rD_0_final.calculate();
537  ITC->rD_1_final.calculate();
538  ITC->rD_2_final.calculate();
539  ITC->rD_3_final.calculate();
540  ITC->rD_4_final.calculate();
541 
542  ITC->der_phiL_final.calculate();
543  ITC->der_zL_final.calculate();
544  ITC->der_phiD_final.calculate();
545  ITC->der_rD_final.calculate();
546 
547  //store the binary results
548  irinv = ITC->rinv_final.ival();
549  iphi0 = ITC->phi0_final.ival();
550  it = ITC->t_final.ival();
551  iz0 = ITC->z0_final.ival();
552 
553  iphiproj[0] = ITC->phiL_0_final.ival();
554  iphiproj[1] = ITC->phiL_1_final.ival();
555  iphiproj[2] = ITC->phiL_2_final.ival();
556  iphiproj[3] = ITC->phiL_3_final.ival();
557 
558  izproj[0] = ITC->zL_0_final.ival();
559  izproj[1] = ITC->zL_1_final.ival();
560  izproj[2] = ITC->zL_2_final.ival();
561  izproj[3] = ITC->zL_3_final.ival();
562 
564  return false;
565 
566  if (!inSector(iphi0, irinv, phi0approx, rinvapprox))
567  return false;
568 
569  for (unsigned int i = 0; i < N_LAYER - 2; ++i) {
570  //reject projection if z is out of range
571  if (izproj[i] < -(1 << (settings_.nzbitsstub(0) - 1)))
572  continue;
573  if (izproj[i] >= (1 << (settings_.nzbitsstub(0) - 1)))
574  continue;
575 
576  //reject projection if phi is out of range
577  if (iphiproj[i] >= (1 << settings_.nphibitsstub(5)) - 1)
578  continue;
579  if (iphiproj[i] <= 0)
580  continue;
581 
582  //Adjust bits for r and z projection depending on layer
583  if (settings_.projlayers(iSeed_, i) <= 3) { //TODO clean up logic
585  } else {
586  izproj[i] >>= (settings_.nzbitsstub(0) - settings_.nzbitsstub(5));
587  }
588 
589  layerprojs[i].init(settings_,
592  iphiproj[i],
593  izproj[i],
594  ITC->der_phiL_final.ival(),
595  ITC->der_zL_final.ival(),
596  phiproj[i],
597  zproj[i],
598  phider[i],
599  zder[i],
600  phiprojapprox[i],
601  zprojapprox[i],
602  ITC->der_phiL_final.fval(),
603  ITC->der_zL_final.fval());
604  }
605 
606  iphiprojdisk[0] = ITC->phiD_0_final.ival();
607  iphiprojdisk[1] = ITC->phiD_1_final.ival();
608  iphiprojdisk[2] = ITC->phiD_2_final.ival();
609  iphiprojdisk[3] = ITC->phiD_3_final.ival();
610  iphiprojdisk[4] = ITC->phiD_4_final.ival();
611 
612  irprojdisk[0] = ITC->rD_0_final.ival();
613  irprojdisk[1] = ITC->rD_1_final.ival();
614  irprojdisk[2] = ITC->rD_2_final.ival();
615  irprojdisk[3] = ITC->rD_3_final.ival();
616  irprojdisk[4] = ITC->rD_4_final.ival();
617 
618  if (std::abs(it * ITC->t_final.K()) > 1.0) {
619  for (unsigned int i = 0; i < N_DISK; ++i) {
620  if (iphiprojdisk[i] <= 0)
621  continue;
622  if (iphiprojdisk[i] >= (1 << settings_.nphibitsstub(0)) - 1)
623  continue;
624 
625  if (irprojdisk[i] < settings_.rmindisk() / ITC->rD_0_final.K() ||
626  irprojdisk[i] > settings_.rmaxdisk() / ITC->rD_0_final.K())
627  continue;
628 
629  diskprojs[i].init(settings_,
630  i + 1,
631  settings_.zmean(i),
632  iphiprojdisk[i],
633  irprojdisk[i],
634  ITC->der_phiD_final.ival(),
635  ITC->der_rD_final.ival(),
636  phiprojdisk[i],
637  rprojdisk[i],
638  phiderdisk[i],
639  rderdisk[i],
640  phiprojdiskapprox[i],
641  rprojdiskapprox[i],
642  ITC->der_phiD_final.fval(),
643  ITC->der_rD_final.fval());
644  }
645  }
646 
647  if (settings_.writeMonitorData("TPars")) {
648  globals_->ofstream("trackletpars.txt")
649  << "Trackpars " << layerdisk1_ + 1 << " " << rinv << " " << rinvapprox << " " << ITC->rinv_final.fval()
650  << " " << phi0 << " " << phi0approx << " " << ITC->phi0_final.fval() << " " << t << " " << tapprox << " "
651  << ITC->t_final.fval() << " " << z0 << " " << z0approx << " " << ITC->z0_final.fval() << endl;
652  }
653 
654  Tracklet* tracklet = new Tracklet(settings_,
655  innerStub,
656  nullptr,
657  outerStub,
658  innerFPGAStub,
659  nullptr,
660  outerFPGAStub,
661  rinv,
662  phi0,
663  0.0,
664  z0,
665  t,
666  rinvapprox,
667  phi0approx,
668  0.0,
669  z0approx,
670  tapprox,
671  irinv,
672  iphi0,
673  0,
674  iz0,
675  it,
676  layerprojs,
677  diskprojs,
678  false);
679 
680  if (settings_.debugTracklet()) {
681  edm::LogVerbatim("Tracklet") << "TrackletCalculator " << getName() << " Found tracklet for seed = " << iSeed_ << " "
682  << iSector_ << " phi0 = " << phi0;
683  }
684 
686  tracklet->setTCIndex(TCIndex_);
687 
688  if (settings_.writeMonitorData("Seeds")) {
689  ofstream fout("seeds.txt", ofstream::app);
690  fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_" << iSector_ << " " << tracklet->getISeed() << endl;
691  fout.close();
692  }
693  trackletpars_->addTracklet(tracklet);
694 
695  if (settings_.bookHistos()) {
697  int tp = tracklet->tpseed();
698  hists->fillTrackletParams(settings_,
699  globals_,
700  iSeed_,
701  iSector_,
702  rinvapprox,
703  irinv * ITC->rinv_final.K(),
704  phi0approx,
705  iphi0 * ITC->phi0_final.K(),
706  asinh(tapprox),
707  asinh(it * ITC->t_final.K()),
708  z0approx,
709  iz0 * ITC->z0_final.K(),
710  tp);
711  }
712 
713  bool addL3 = false;
714  bool addL4 = false;
715  bool addL5 = false;
716  bool addL6 = false;
717  for (unsigned int j = 0; j < N_LAYER - 2; j++) {
718  int lproj = settings_.projlayers(iSeed_, j);
719  bool added = false;
720  if (tracklet->validProj(lproj)) {
721  added = addLayerProj(tracklet, lproj);
722  if (added && lproj == 3)
723  addL3 = true;
724  if (added && lproj == 4)
725  addL4 = true;
726  if (added && lproj == 5)
727  addL5 = true;
728  if (added && lproj == 6)
729  addL6 = true;
730  }
731  }
732 
733  for (unsigned int j = 0; j < N_DISK - 1; j++) { //no projections to 5th disk!!
734  int disk = j + 1;
735  if (disk == 4 && addL3)
736  continue;
737  if (disk == 3 && addL4)
738  continue;
739  if (disk == 2 && addL5)
740  continue;
741  if (disk == 1 && addL6)
742  continue;
743  if (it < 0)
744  disk = -disk;
745  if (tracklet->validProjDisk(abs(disk))) {
746  addDiskProj(tracklet, disk);
747  }
748  }
749 
750  return true;
751 }
752 
753 bool TrackletCalculatorBase::diskSeeding(const Stub* innerFPGAStub,
754  const L1TStub* innerStub,
755  const Stub* outerFPGAStub,
756  const L1TStub* outerStub) {
757  if (settings_.debugTracklet()) {
758  edm::LogVerbatim("Tracklet") << "TrackletCalculator::execute calculate disk seeds";
759  }
760 
761  int sign = 1;
762  if (innerFPGAStub->disk().value() < 0)
763  sign = -1;
764 
765  int disk = innerFPGAStub->disk().value();
766  assert(abs(disk) == 1 || abs(disk) == 3);
767 
768  assert(innerStub->isPSmodule());
769  assert(outerStub->isPSmodule());
770 
771  double r1 = innerStub->r();
772  double z1 = innerStub->z();
773  double phi1 = innerStub->phi();
774 
775  double r2 = outerStub->r();
776  double z2 = outerStub->z();
777  double phi2 = outerStub->phi();
778 
779  if (r2 < r1 + 2.0) {
780  return false; //Protection... Should be handled cleaner to avoid problem with floating point calculation
781  }
782 
783  double rinv, phi0, t, z0;
784 
785  double phiproj[N_PSLAYER], zproj[N_PSLAYER], phider[N_PSLAYER], zder[N_PSLAYER];
786  double phiprojdisk[N_DISK - 2], rprojdisk[N_DISK - 2], phiderdisk[N_DISK - 2], rderdisk[N_DISK - 2];
787 
789  z1,
790  phi1,
791  r2,
792  z2,
793  phi2,
794  outerStub->sigmaz(),
795  rinv,
796  phi0,
797  t,
798  z0,
799  phiproj,
800  zproj,
801  phider,
802  zder,
803  phiprojdisk,
804  rprojdisk,
805  phiderdisk,
806  rderdisk);
807 
808  //Truncates floating point positions to integer representation precision
809  if (settings_.useapprox()) {
810  phi1 = innerFPGAStub->phiapprox(phimin_, phimax_);
811  z1 = innerFPGAStub->zapprox();
812  r1 = innerFPGAStub->rapprox();
813 
814  phi2 = outerFPGAStub->phiapprox(phimin_, phimax_);
815  z2 = outerFPGAStub->zapprox();
816  r2 = outerFPGAStub->rapprox();
817  }
818 
819  double rinvapprox, phi0approx, tapprox, z0approx;
820  double phiprojapprox[N_PSLAYER], zprojapprox[N_PSLAYER];
821  double phiprojdiskapprox[N_DISK - 2], rprojdiskapprox[N_DISK - 2];
822 
824  if (disk == 1)
825  ITC = globals_->ITC_F1F2();
826  else if (disk == 3)
827  ITC = globals_->ITC_F3F4();
828  else if (disk == -1)
829  ITC = globals_->ITC_B1B2();
830  else
831  ITC = globals_->ITC_B3B4();
832 
833  ITC->r1.set_fval(r1);
834  ITC->r2.set_fval(r2);
835  int signt = t > 0 ? 1 : -1;
836  ITC->z1.set_fval(z1 - signt * settings_.zmean(layerdisk1_ - N_LAYER));
837  ITC->z2.set_fval(z2 - signt * settings_.zmean(layerdisk2_ - N_LAYER));
838  double sphi1 = angle0to2pi::make0To2pi(phi1 - phioffset_);
839  double sphi2 = angle0to2pi::make0To2pi(phi2 - phioffset_);
840  ITC->phi1.set_fval(sphi1);
841  ITC->phi2.set_fval(sphi2);
842 
843  ITC->rproj0.set_fval(settings_.rmean(0));
844  ITC->rproj1.set_fval(settings_.rmean(1));
845  ITC->rproj2.set_fval(settings_.rmean(2));
846 
847  ITC->zproj0.set_fval(signt * settings_.zmean(settings_.projdisks(iSeed_, 0) - 1));
848  ITC->zproj1.set_fval(signt * settings_.zmean(settings_.projdisks(iSeed_, 1) - 1));
849  ITC->zproj2.set_fval(signt * settings_.zmean(settings_.projdisks(iSeed_, 2) - 1));
850 
851  ITC->rinv_final.calculate();
852  ITC->phi0_final.calculate();
853  ITC->t_final.calculate();
854  ITC->z0_final.calculate();
855 
856  ITC->phiL_0_final.calculate();
857  ITC->phiL_1_final.calculate();
858  ITC->phiL_2_final.calculate();
859 
860  ITC->zL_0_final.calculate();
861  ITC->zL_1_final.calculate();
862  ITC->zL_2_final.calculate();
863 
864  ITC->phiD_0_final.calculate();
865  ITC->phiD_1_final.calculate();
866  ITC->phiD_2_final.calculate();
867 
868  ITC->rD_0_final.calculate();
869  ITC->rD_1_final.calculate();
870  ITC->rD_2_final.calculate();
871 
872  ITC->der_phiL_final.calculate();
873  ITC->der_zL_final.calculate();
874  ITC->der_phiD_final.calculate();
875  ITC->der_rD_final.calculate();
876 
877  //store the approximate results
878  rinvapprox = ITC->rinv_final.fval();
879  phi0approx = ITC->phi0_final.fval();
880  tapprox = ITC->t_final.fval();
881  z0approx = ITC->z0_final.fval();
882 
883  phiprojapprox[0] = ITC->phiL_0_final.fval();
884  phiprojapprox[1] = ITC->phiL_1_final.fval();
885  phiprojapprox[2] = ITC->phiL_2_final.fval();
886 
887  zprojapprox[0] = ITC->zL_0_final.fval();
888  zprojapprox[1] = ITC->zL_1_final.fval();
889  zprojapprox[2] = ITC->zL_2_final.fval();
890 
891  phiprojdiskapprox[0] = ITC->phiD_0_final.fval();
892  phiprojdiskapprox[1] = ITC->phiD_1_final.fval();
893  phiprojdiskapprox[2] = ITC->phiD_2_final.fval();
894 
895  rprojdiskapprox[0] = ITC->rD_0_final.fval();
896  rprojdiskapprox[1] = ITC->rD_1_final.fval();
897  rprojdiskapprox[2] = ITC->rD_2_final.fval();
898 
899  //now binary
900 
901  int irinv, iphi0, it, iz0;
902  int iphiproj[N_PSLAYER], izproj[N_PSLAYER];
903 
904  int iphiprojdisk[N_DISK - 2], irprojdisk[N_DISK - 2];
905 
906  int ir1 = innerFPGAStub->r().value();
907  int iphi1 = innerFPGAStub->phi().value();
908  int iz1 = innerFPGAStub->z().value();
909 
910  int ir2 = outerFPGAStub->r().value();
911  int iphi2 = outerFPGAStub->phi().value();
912  int iz2 = outerFPGAStub->z().value();
913 
914  //To get same precission as for layers.
915  iphi1 <<= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
916  iphi2 <<= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
917 
918  ITC->r1.set_ival(ir1);
919  ITC->r2.set_ival(ir2);
920  ITC->z1.set_ival(iz1);
921  ITC->z2.set_ival(iz2);
922  ITC->phi1.set_ival(iphi1);
923  ITC->phi2.set_ival(iphi2);
924 
925  ITC->rinv_final.calculate();
926  ITC->phi0_final.calculate();
927  ITC->t_final.calculate();
928  ITC->z0_final.calculate();
929 
930  ITC->phiL_0_final.calculate();
931  ITC->phiL_1_final.calculate();
932  ITC->phiL_2_final.calculate();
933 
934  ITC->zL_0_final.calculate();
935  ITC->zL_1_final.calculate();
936  ITC->zL_2_final.calculate();
937 
938  ITC->phiD_0_final.calculate();
939  ITC->phiD_1_final.calculate();
940  ITC->phiD_2_final.calculate();
941 
942  ITC->rD_0_final.calculate();
943  ITC->rD_1_final.calculate();
944  ITC->rD_2_final.calculate();
945 
946  ITC->der_phiL_final.calculate();
947  ITC->der_zL_final.calculate();
948  ITC->der_phiD_final.calculate();
949  ITC->der_rD_final.calculate();
950 
951  //store the binary results
952  irinv = ITC->rinv_final.ival();
953  iphi0 = ITC->phi0_final.ival();
954  it = ITC->t_final.ival();
955  iz0 = ITC->z0_final.ival();
956 
957  iphiproj[0] = ITC->phiL_0_final.ival();
958  iphiproj[1] = ITC->phiL_1_final.ival();
959  iphiproj[2] = ITC->phiL_2_final.ival();
960 
961  izproj[0] = ITC->zL_0_final.ival();
962  izproj[1] = ITC->zL_1_final.ival();
963  izproj[2] = ITC->zL_2_final.ival();
964 
966  return false;
967 
968  if (!inSector(iphi0, irinv, phi0approx, rinvapprox))
969  return false;
970 
971  LayerProjection layerprojs[N_LAYER - 2];
972  DiskProjection diskprojs[N_DISK - 2];
973 
974  for (unsigned int i = 0; i < N_DISK - 2; ++i) {
975  //Check is outside z range
976  if (izproj[i] < -(1 << (settings_.nzbitsstub(0) - 1)))
977  continue;
978  if (izproj[i] >= (1 << (settings_.nzbitsstub(0) - 1)))
979  continue;
980 
981  //Check if outside phi range
982  if (iphiproj[i] >= (1 << settings_.nphibitsstub(5)) - 1)
983  continue;
984  if (iphiproj[i] <= 0)
985  continue;
986 
987  //shift bits - allways in PS modules for disk seeding
988  iphiproj[i] >>= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
989 
990  layerprojs[i].init(settings_,
991  i + 1,
992  settings_.rmean(i),
993  iphiproj[i],
994  izproj[i],
995  ITC->der_phiL_final.ival(),
996  ITC->der_zL_final.ival(),
997  phiproj[i],
998  zproj[i],
999  phider[i],
1000  zder[i],
1001  phiprojapprox[i],
1002  zprojapprox[i],
1003  ITC->der_phiL_final.fval(),
1004  ITC->der_zL_final.fval());
1005  }
1006 
1007  iphiprojdisk[0] = ITC->phiD_0_final.ival();
1008  iphiprojdisk[1] = ITC->phiD_1_final.ival();
1009  iphiprojdisk[2] = ITC->phiD_2_final.ival();
1010 
1011  irprojdisk[0] = ITC->rD_0_final.ival();
1012  irprojdisk[1] = ITC->rD_1_final.ival();
1013  irprojdisk[2] = ITC->rD_2_final.ival();
1014 
1015  for (unsigned int i = 0; i < N_DISK - 2; ++i) {
1016  //check that phi projection in range
1017  if (iphiprojdisk[i] <= 0)
1018  continue;
1019  if (iphiprojdisk[i] >= (1 << settings_.nphibitsstub(0)) - 1)
1020  continue;
1021 
1022  //check that r projection in range
1023  if (irprojdisk[i] <= 0 || irprojdisk[i] > settings_.rmaxdisk() / ITC->rD_0_final.K())
1024  continue;
1025 
1026  diskprojs[i].init(settings_,
1027  i + 1,
1029  iphiprojdisk[i],
1030  irprojdisk[i],
1031  ITC->der_phiD_final.ival(),
1032  ITC->der_rD_final.ival(),
1033  phiprojdisk[i],
1034  rprojdisk[i],
1035  phiderdisk[i],
1036  rderdisk[i],
1037  phiprojdiskapprox[i],
1038  rprojdiskapprox[i],
1039  ITC->der_phiD_final.fval(),
1040  ITC->der_rD_final.fval());
1041  }
1042 
1043  if (settings_.writeMonitorData("TPars")) {
1044  globals_->ofstream("trackletparsdisk.txt")
1045  << "Trackpars " << layerdisk1_ - 5 << " " << rinv << " " << rinvapprox << " "
1046  << ITC->rinv_final.fval() << " " << phi0 << " " << phi0approx << " " << ITC->phi0_final.fval() << " " << t
1047  << " " << tapprox << " " << ITC->t_final.fval() << " " << z0 << " " << z0approx << " " << ITC->z0_final.fval()
1048  << endl;
1049  }
1050 
1051  Tracklet* tracklet = new Tracklet(settings_,
1052  innerStub,
1053  nullptr,
1054  outerStub,
1055  innerFPGAStub,
1056  nullptr,
1057  outerFPGAStub,
1058  rinv,
1059  phi0,
1060  0.0,
1061  z0,
1062  t,
1063  rinvapprox,
1064  phi0approx,
1065  0.0,
1066  z0approx,
1067  tapprox,
1068  irinv,
1069  iphi0,
1070  0,
1071  iz0,
1072  it,
1073  layerprojs,
1074  diskprojs,
1075  true);
1076 
1077  if (settings_.debugTracklet()) {
1078  edm::LogVerbatim("Tracklet") << "Found tracklet for disk seed = " << iSeed_ << " " << tracklet << " " << iSector_;
1079  }
1080 
1082  tracklet->setTCIndex(TCIndex_);
1083 
1084  if (settings_.writeMonitorData("Seeds")) {
1085  ofstream fout("seeds.txt", ofstream::app);
1086  fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_" << iSector_ << " " << tracklet->getISeed() << endl;
1087  fout.close();
1088  }
1089  trackletpars_->addTracklet(tracklet);
1090 
1091  if (tracklet->validProj(1)) {
1092  addLayerProj(tracklet, 1);
1093  }
1094 
1095  if (tracklet->validProj(2)) {
1096  addLayerProj(tracklet, 2);
1097  }
1098 
1099  for (unsigned int j = 0; j < N_DISK - 2; j++) {
1100  if (tracklet->validProjDisk(sign * settings_.projdisks(iSeed_, j))) {
1101  addDiskProj(tracklet, sign * settings_.projdisks(iSeed_, j));
1102  }
1103  }
1104 
1105  return true;
1106 }
1107 
1109  const L1TStub* innerStub,
1110  const Stub* outerFPGAStub,
1111  const L1TStub* outerStub) {
1112  //Deal with overlap stubs here
1113  assert(outerFPGAStub->isBarrel());
1114 
1115  assert(innerFPGAStub->isDisk());
1116 
1117  int disk = innerFPGAStub->disk().value();
1118 
1119  if (settings_.debugTracklet()) {
1120  edm::LogVerbatim("Tracklet") << "trying to make overlap tracklet for seed = " << iSeed_ << " " << getName();
1121  }
1122 
1123  double r1 = innerStub->r();
1124  double z1 = innerStub->z();
1125  double phi1 = innerStub->phi();
1126 
1127  double r2 = outerStub->r();
1128  double z2 = outerStub->z();
1129  double phi2 = outerStub->phi();
1130 
1131  //Protection for wrong radii. Could be handled cleaner to avoid problem with floating point calculation and with overflows in the integer calculation.
1132  if (r1 < r2 + 1.5) {
1133  return false;
1134  }
1135 
1136  double rinv, phi0, t, z0;
1137 
1138  double phiproj[N_PSLAYER], zproj[N_PSLAYER], phider[N_PSLAYER], zder[N_PSLAYER];
1139  double phiprojdisk[N_DISK - 1], rprojdisk[N_DISK - 1], phiderdisk[N_DISK - 1], rderdisk[N_DISK - 1];
1140 
1142  z1,
1143  phi1,
1144  r2,
1145  z2,
1146  phi2,
1147  outerStub->sigmaz(),
1148  rinv,
1149  phi0,
1150  t,
1151  z0,
1152  phiproj,
1153  zproj,
1154  phider,
1155  zder,
1156  phiprojdisk,
1157  rprojdisk,
1158  phiderdisk,
1159  rderdisk);
1160 
1161  //Truncates floating point positions to integer representation precision
1162  if (settings_.useapprox()) {
1163  phi1 = innerFPGAStub->phiapprox(phimin_, phimax_);
1164  z1 = innerFPGAStub->zapprox();
1165  r1 = innerFPGAStub->rapprox();
1166 
1167  phi2 = outerFPGAStub->phiapprox(phimin_, phimax_);
1168  z2 = outerFPGAStub->zapprox();
1169  r2 = outerFPGAStub->rapprox();
1170  }
1171 
1172  double rinvapprox, phi0approx, tapprox, z0approx;
1173  double phiprojapprox[N_PSLAYER], zprojapprox[N_PSLAYER];
1174  double phiprojdiskapprox[N_DISK - 1], rprojdiskapprox[N_DISK - 1];
1175 
1177  int ll = outerFPGAStub->layer().value() + 1;
1178  if (ll == 1 && disk == 1)
1179  ITC = globals_->ITC_L1F1();
1180  else if (ll == 2 && disk == 1)
1181  ITC = globals_->ITC_L2F1();
1182  else if (ll == 1 && disk == -1)
1183  ITC = globals_->ITC_L1B1();
1184  else if (ll == 2 && disk == -1)
1185  ITC = globals_->ITC_L2B1();
1186  else
1187  throw cms::Exception("LogicError") << __FILE__ << " " << __LINE__ << " Invalid seeding!";
1188 
1189  ITC->r1.set_fval(r2 - settings_.rmean(ll - 1));
1190  ITC->r2.set_fval(r1);
1191  int signt = t > 0 ? 1 : -1;
1192  ITC->z1.set_fval(z2);
1193  ITC->z2.set_fval(z1 - signt * settings_.zmean(layerdisk2_ - N_LAYER));
1194  double sphi1 = angle0to2pi::make0To2pi(phi1 - phioffset_);
1195  double sphi2 = angle0to2pi::make0To2pi(phi2 - phioffset_);
1196  ITC->phi1.set_fval(sphi2);
1197  ITC->phi2.set_fval(sphi1);
1198 
1199  ITC->rproj0.set_fval(settings_.rmean(0));
1200  ITC->rproj1.set_fval(settings_.rmean(1));
1201  ITC->rproj2.set_fval(settings_.rmean(2));
1202 
1203  ITC->zproj0.set_fval(signt * settings_.zmean(1));
1204  ITC->zproj1.set_fval(signt * settings_.zmean(2));
1205  ITC->zproj2.set_fval(signt * settings_.zmean(3));
1206  ITC->zproj3.set_fval(signt * settings_.zmean(4));
1207 
1208  ITC->rinv_final.calculate();
1209  ITC->phi0_final.calculate();
1210  ITC->t_final.calculate();
1211  ITC->z0_final.calculate();
1212 
1213  ITC->phiL_0_final.calculate();
1214  ITC->phiL_1_final.calculate();
1215  ITC->phiL_2_final.calculate();
1216 
1217  ITC->zL_0_final.calculate();
1218  ITC->zL_1_final.calculate();
1219  ITC->zL_2_final.calculate();
1220 
1221  ITC->phiD_0_final.calculate();
1222  ITC->phiD_1_final.calculate();
1223  ITC->phiD_2_final.calculate();
1224  ITC->phiD_3_final.calculate();
1225 
1226  ITC->rD_0_final.calculate();
1227  ITC->rD_1_final.calculate();
1228  ITC->rD_2_final.calculate();
1229  ITC->rD_3_final.calculate();
1230 
1231  ITC->der_phiL_final.calculate();
1232  ITC->der_zL_final.calculate();
1233  ITC->der_phiD_final.calculate();
1234  ITC->der_rD_final.calculate();
1235 
1236  //store the approximate results
1237  rinvapprox = ITC->rinv_final.fval();
1238  phi0approx = ITC->phi0_final.fval();
1239  tapprox = ITC->t_final.fval();
1240  z0approx = ITC->z0_final.fval();
1241 
1242  phiprojapprox[0] = ITC->phiL_0_final.fval();
1243  phiprojapprox[1] = ITC->phiL_1_final.fval();
1244  phiprojapprox[2] = ITC->phiL_2_final.fval();
1245 
1246  zprojapprox[0] = ITC->zL_0_final.fval();
1247  zprojapprox[1] = ITC->zL_1_final.fval();
1248  zprojapprox[2] = ITC->zL_2_final.fval();
1249 
1250  phiprojdiskapprox[0] = ITC->phiD_0_final.fval();
1251  phiprojdiskapprox[1] = ITC->phiD_1_final.fval();
1252  phiprojdiskapprox[2] = ITC->phiD_2_final.fval();
1253  phiprojdiskapprox[3] = ITC->phiD_3_final.fval();
1254 
1255  rprojdiskapprox[0] = ITC->rD_0_final.fval();
1256  rprojdiskapprox[1] = ITC->rD_1_final.fval();
1257  rprojdiskapprox[2] = ITC->rD_2_final.fval();
1258  rprojdiskapprox[3] = ITC->rD_3_final.fval();
1259 
1260  //now binary
1261 
1262  int irinv, iphi0, it, iz0;
1263  int iphiproj[N_LAYER - 2], izproj[N_LAYER - 2];
1264  int iphiprojdisk[N_DISK], irprojdisk[N_DISK];
1265 
1266  int ir2 = innerFPGAStub->r().value();
1267  int iphi2 = innerFPGAStub->phi().value();
1268  int iz2 = innerFPGAStub->z().value();
1269 
1270  int ir1 = outerFPGAStub->r().value();
1271  int iphi1 = outerFPGAStub->phi().value();
1272  int iz1 = outerFPGAStub->z().value();
1273 
1274  //To get global precission
1275  ir1 <<= (8 - settings_.nrbitsstub(ll - 1));
1276  iphi1 <<= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
1277  iphi2 <<= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
1278 
1279  ITC->r1.set_ival(ir1);
1280  ITC->r2.set_ival(ir2);
1281  ITC->z1.set_ival(iz1);
1282  ITC->z2.set_ival(iz2);
1283  ITC->phi1.set_ival(iphi1);
1284  ITC->phi2.set_ival(iphi2);
1285 
1286  ITC->rinv_final.calculate();
1287  ITC->phi0_final.calculate();
1288  ITC->t_final.calculate();
1289  ITC->z0_final.calculate();
1290 
1291  ITC->phiL_0_final.calculate();
1292  ITC->phiL_1_final.calculate();
1293  ITC->phiL_2_final.calculate();
1294 
1295  ITC->zL_0_final.calculate();
1296  ITC->zL_1_final.calculate();
1297  ITC->zL_2_final.calculate();
1298 
1299  ITC->phiD_0_final.calculate();
1300  ITC->phiD_1_final.calculate();
1301  ITC->phiD_2_final.calculate();
1302  ITC->phiD_3_final.calculate();
1303 
1304  ITC->rD_0_final.calculate();
1305  ITC->rD_1_final.calculate();
1306  ITC->rD_2_final.calculate();
1307  ITC->rD_3_final.calculate();
1308 
1309  ITC->der_phiL_final.calculate();
1310  ITC->der_zL_final.calculate();
1311  ITC->der_phiD_final.calculate();
1312  ITC->der_rD_final.calculate();
1313 
1314  //store the binary results
1315  irinv = ITC->rinv_final.ival();
1316  iphi0 = ITC->phi0_final.ival();
1317  it = ITC->t_final.ival();
1318  iz0 = ITC->z0_final.ival();
1319 
1320  iphiproj[0] = ITC->phiL_0_final.ival();
1321  iphiproj[1] = ITC->phiL_1_final.ival();
1322  iphiproj[2] = ITC->phiL_2_final.ival();
1323 
1324  izproj[0] = ITC->zL_0_final.ival();
1325  izproj[1] = ITC->zL_1_final.ival();
1326  izproj[2] = ITC->zL_2_final.ival();
1327 
1328  iphiprojdisk[0] = ITC->phiD_0_final.ival();
1329  iphiprojdisk[1] = ITC->phiD_1_final.ival();
1330  iphiprojdisk[2] = ITC->phiD_2_final.ival();
1331  iphiprojdisk[3] = ITC->phiD_3_final.ival();
1332 
1333  irprojdisk[0] = ITC->rD_0_final.ival();
1334  irprojdisk[1] = ITC->rD_1_final.ival();
1335  irprojdisk[2] = ITC->rD_2_final.ival();
1336  irprojdisk[3] = ITC->rD_3_final.ival();
1337 
1339  return false;
1340 
1341  if (!inSector(iphi0, irinv, phi0approx, rinvapprox))
1342  return false;
1343 
1344  LayerProjection layerprojs[N_LAYER - 2];
1345  DiskProjection diskprojs[N_DISK];
1346 
1347  for (unsigned int i = 0; i < N_DISK - 2; ++i) {
1348  //check that zproj is in range
1349  if (izproj[i] < -(1 << (settings_.nzbitsstub(0) - 1)))
1350  continue;
1351  if (izproj[i] >= (1 << (settings_.nzbitsstub(0) - 1)))
1352  continue;
1353 
1354  //check that phiproj is in range
1355  if (iphiproj[i] >= (1 << settings_.nphibitsstub(5)) - 1)
1356  continue;
1357  if (iphiproj[i] <= 0)
1358  continue;
1359 
1360  //adjust bits for PS modules (no 2S modules in overlap seeds)
1361  iphiproj[i] >>= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
1362 
1363  layerprojs[i].init(settings_,
1364  i + 1,
1365  settings_.rmean(i),
1366  iphiproj[i],
1367  izproj[i],
1368  ITC->der_phiL_final.ival(),
1369  ITC->der_zL_final.ival(),
1370  phiproj[i],
1371  zproj[i],
1372  phider[i],
1373  zder[i],
1374  phiprojapprox[i],
1375  zprojapprox[i],
1376  ITC->der_phiL_final.fval(),
1377  ITC->der_zL_final.fval());
1378  }
1379 
1380  for (int i = 0; i < 4; ++i) {
1381  //check that phi projection in range
1382  if (iphiprojdisk[i] <= 0)
1383  continue;
1384  if (iphiprojdisk[i] >= (1 << settings_.nphibitsstub(0)) - 1)
1385  continue;
1386 
1387  //check that r projection in range
1388  if (irprojdisk[i] <= 0 || irprojdisk[i] > settings_.rmaxdisk() / ITC->rD_0_final.K())
1389  continue;
1390 
1391  diskprojs[i].init(settings_,
1392  i + 1,
1393  settings_.zmean(i),
1394  iphiprojdisk[i],
1395  irprojdisk[i],
1396  ITC->der_phiD_final.ival(),
1397  ITC->der_rD_final.ival(),
1398  phiprojdisk[i],
1399  rprojdisk[i],
1400  phiderdisk[i],
1401  rderdisk[i],
1402  phiprojdiskapprox[i],
1403  rprojdiskapprox[i],
1404  ITC->der_phiD_final.fval(),
1405  ITC->der_rD_final.fval());
1406  }
1407 
1408  if (settings_.writeMonitorData("TPars")) {
1409  globals_->ofstream("trackletparsoverlap.txt")
1410  << "Trackpars " << layerdisk1_ - 5 << " " << rinv << " " << irinv << " " << ITC->rinv_final.fval() << " "
1411  << phi0 << " " << iphi0 << " " << ITC->phi0_final.fval() << " " << t << " " << it << " "
1412  << ITC->t_final.fval() << " " << z0 << " " << iz0 << " " << ITC->z0_final.fval() << endl;
1413  }
1414 
1415  Tracklet* tracklet = new Tracklet(settings_,
1416  innerStub,
1417  nullptr,
1418  outerStub,
1419  innerFPGAStub,
1420  nullptr,
1421  outerFPGAStub,
1422  rinv,
1423  phi0,
1424  0.0,
1425  z0,
1426  t,
1427  rinvapprox,
1428  phi0approx,
1429  0.0,
1430  z0approx,
1431  tapprox,
1432  irinv,
1433  iphi0,
1434  0,
1435  iz0,
1436  it,
1437  layerprojs,
1438  diskprojs,
1439  false,
1440  true);
1441 
1442  if (settings_.debugTracklet()) {
1443  edm::LogVerbatim("Tracklet") << "Found tracklet in overlap seed = " << iSeed_ << " " << tracklet << " " << iSector_;
1444  }
1445 
1447  tracklet->setTCIndex(TCIndex_);
1448 
1449  if (settings_.writeMonitorData("Seeds")) {
1450  ofstream fout("seeds.txt", ofstream::app);
1451  fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_" << iSector_ << " " << tracklet->getISeed() << endl;
1452  fout.close();
1453  }
1454  trackletpars_->addTracklet(tracklet);
1455 
1456  int layer = outerFPGAStub->layer().value() + 1;
1457 
1458  if (layer == 2) {
1459  if (tracklet->validProj(1)) {
1460  addLayerProj(tracklet, 1);
1461  }
1462  }
1463 
1464  for (unsigned int disk = 2; disk < 6; disk++) {
1465  if (layer == 2 && disk == 5)
1466  continue;
1467  if (tracklet->validProjDisk(disk)) {
1468  addDiskProj(tracklet, disk);
1469  }
1470  }
1471 
1472  return true;
1473 }
trklet::VarBase::fval
double fval() const
Definition: imath.h:212
trklet::N_PSLAYER
constexpr unsigned int N_PSLAYER
Definition: Settings.h:21
trklet::IMATH_TrackletCalculatorDisk::zproj0
VarDef zproj0
Definition: IMATH_TrackletCalculatorDisk.h:155
trklet::IMATH_TrackletCalculatorDisk::der_phiD_final
VarAdjustK der_phiD_final
Definition: IMATH_TrackletCalculatorDisk.h:299
trklet::Settings::phicritmaxmc
double phicritmaxmc() const
Definition: Settings.h:244
trklet::IMATH_TrackletCalculator::rproj1
VarDef rproj1
Definition: IMATH_TrackletCalculator.h:165
trklet::IMATH_TrackletCalculator::r1
VarDef r1
Definition: IMATH_TrackletCalculator.h:155
trklet::TrackletCalculatorBase::iSeed_
unsigned int iSeed_
Definition: TrackletCalculatorBase.h:128
trklet::Settings::rmindisk
double rmindisk() const
Definition: Settings.h:103
trklet::IMATH_TrackletCalculator::rD_0_final
VarAdjustK rD_0_final
Definition: IMATH_TrackletCalculator.h:374
trklet::Settings::rmindiskvm
double rmindiskvm() const
Definition: Settings.h:260
mps_fire.i
i
Definition: mps_fire.py:428
trklet::Globals::ITC_L2L3
IMATH_TrackletCalculator * ITC_L2L3()
Definition: Globals.h:53
trklet::IMATH_TrackletCalculator::rD_3_final
VarAdjustK rD_3_final
Definition: IMATH_TrackletCalculator.h:377
trklet::TrackletCalculatorBase::addProjectionDisk
void addProjectionDisk(int disk, int iphi, TrackletProjectionsMemory *trackletprojs, Tracklet *tracklet)
Definition: TrackletCalculatorBase.cc:270
trklet::TrackletCalculatorBase::phioffset_
double phioffset_
Definition: TrackletCalculatorBase.h:134
trklet::IMATH_TrackletCalculatorOverlap::zL_0_final
VarAdjustK zL_0_final
Definition: IMATH_TrackletCalculatorOverlap.h:271
MessageLogger.h
trklet::IMATH_TrackletCalculatorOverlap::rD_0_final
VarAdjustK rD_0_final
Definition: IMATH_TrackletCalculatorOverlap.h:340
trklet::TrackletCalculatorBase::trackletpars_
TrackletParametersMemory * trackletpars_
Definition: TrackletCalculatorBase.h:139
trklet::IMATH_TrackletCalculator::z1
VarDef z1
Definition: IMATH_TrackletCalculator.h:157
trklet::Settings::krprojshiftdisk
double krprojshiftdisk() const
Definition: Settings.h:341
trklet::IMATH_TrackletCalculatorOverlap::zproj1
VarDef zproj1
Definition: IMATH_TrackletCalculatorOverlap.h:161
trklet::Tracklet::fpgaphiprojdisk
const FPGAWord & fpgaphiprojdisk(int disk) const
Definition: Tracklet.h:266
trklet::TrackletParametersMemory::addTracklet
void addTracklet(Tracklet *tracklet)
Definition: TrackletParametersMemory.h:23
trklet::Stub::phi
const FPGAWord & phi() const
Definition: Stub.h:51
trklet::Settings::krinvpars
double krinvpars() const
Definition: Settings.h:327
trklet::IMATH_TrackletCalculatorDisk::zL_0_final
VarAdjustK zL_0_final
Definition: IMATH_TrackletCalculatorDisk.h:265
trklet::IMATH_TrackletCalculatorOverlap::phiD_3_final
VarAdjustK phiD_3_final
Definition: IMATH_TrackletCalculatorOverlap.h:306
deltaPhi.h
trklet::LayerProjection::init
void init(Settings const &settings, int projlayer, double rproj, int iphiproj, int izproj, int iphider, int izder, double phiproj, double zproj, double phiprojder, double zprojder, double phiprojapprox, double zprojapprox, double phiprojderapprox, double zprojderapprox)
Definition: LayerProjection.cc:9
trklet::IMATH_TrackletCalculatorDisk::t_final
VarAdjustK t_final
Definition: IMATH_TrackletCalculatorDisk.h:203
trklet::IMATH_TrackletCalculatorDisk::phi2
VarDef phi2
Definition: IMATH_TrackletCalculatorDisk.h:149
trklet::Settings::phicritminmc
double phicritminmc() const
Definition: Settings.h:243
trklet::IMATH_TrackletCalculatorOverlap::z1
VarDef z1
Definition: IMATH_TrackletCalculatorOverlap.h:150
trklet::TrackletProjectionsMemory
Definition: TrackletProjectionsMemory.h:15
trklet::VarBase::calculate
bool calculate(int debug_level)
Definition: imath_calculate.cc:6
trklet::Settings
Definition: Settings.h:26
trklet::IMATH_TrackletCalculator::zL_2_final
VarAdjustKR zL_2_final
Definition: IMATH_TrackletCalculator.h:295
trklet::L1TStub
Definition: L1TStub.h:12
trklet::L1TStub::z
double z() const
Definition: L1TStub.h:56
trklet::IMATH_TrackletCalculatorOverlap::rD_3_final
VarAdjustK rD_3_final
Definition: IMATH_TrackletCalculatorOverlap.h:343
trklet::IMATH_TrackletCalculatorOverlap::rD_1_final
VarAdjustK rD_1_final
Definition: IMATH_TrackletCalculatorOverlap.h:341
trklet::Tracklet::getISeed
int getISeed() const
Definition: Tracklet.cc:849
trklet::IMATH_TrackletCalculatorOverlap::der_zL_final
VarAdjustK der_zL_final
Definition: IMATH_TrackletCalculatorOverlap.h:275
trklet::IMATH_TrackletCalculatorDisk::phi1
VarDef phi1
Definition: IMATH_TrackletCalculatorDisk.h:148
trklet::IMATH_TrackletCalculator::zL_1_final
VarAdjustKR zL_1_final
Definition: IMATH_TrackletCalculator.h:294
cms::cuda::assert
assert(be >=bs)
trklet::TrackletCalculatorBase::addProjection
void addProjection(int layer, int iphi, TrackletProjectionsMemory *trackletprojs, Tracklet *tracklet)
Definition: TrackletCalculatorBase.cc:255
trklet::IMATH_TrackletCalculator::phiD_0_final
VarAdjustK phiD_0_final
Definition: IMATH_TrackletCalculator.h:331
Validation_hcalonly_cfi.sign
sign
Definition: Validation_hcalonly_cfi.py:32
trklet::Settings::rmean
double rmean(unsigned int iLayer) const
Definition: Settings.h:128
trklet::IMATH_TrackletCalculatorDisk::rD_2_final
VarAdjustK rD_2_final
Definition: IMATH_TrackletCalculatorDisk.h:326
trklet::TrackletCalculatorBase::barrelSeeding
bool barrelSeeding(const Stub *innerFPGAStub, const L1TStub *innerStub, const Stub *outerFPGAStub, const L1TStub *outerStub)
Definition: TrackletCalculatorBase.cc:331
trklet::IMATH_TrackletCalculatorDisk
Definition: IMATH_TrackletCalculatorDisk.h:19
trklet::N_DISK
constexpr int N_DISK
Definition: Settings.h:20
trklet::IMATH_TrackletCalculatorDisk::rD_0_final
VarAdjustK rD_0_final
Definition: IMATH_TrackletCalculatorDisk.h:324
trklet::ProcessBase::settings_
Settings const & settings_
Definition: ProcessBase.h:44
compare.hists
hists
Definition: compare.py:319
trklet::Stub::phiapprox
double phiapprox(double phimin, double) const
Definition: Stub.cc:236
trklet::IMATH_TrackletCalculatorDisk::z2
VarDef z2
Definition: IMATH_TrackletCalculatorDisk.h:146
trklet::FPGAWord::nbits
int nbits() const
Definition: FPGAWord.h:25
createJobs.tmp
tmp
align.sh
Definition: createJobs.py:716
trklet::Globals::ITC_B1B2
IMATH_TrackletCalculatorDisk * ITC_B1B2()
Definition: Globals.h:59
trklet::Globals
Definition: Globals.h:32
trklet::IMATH_TrackletCalculatorDisk::phi0_final
VarAdjustK phi0_final
Definition: IMATH_TrackletCalculatorDisk.h:202
trklet::TrackletCalculatorBase::trackletprojlayers_
std::vector< std::vector< TrackletProjectionsMemory * > > trackletprojlayers_
Definition: TrackletCalculatorBase.h:142
trklet::Tracklet::validProj
bool validProj(int layer) const
Definition: Tracklet.h:89
trklet::Stub::isBarrel
bool isBarrel() const
Definition: Stub.h:60
trklet::Stub::r
const FPGAWord & r() const
Definition: Stub.h:49
trklet::LayerProjection
Definition: LayerProjection.h:10
trklet::Settings::bookHistos
bool bookHistos() const
Definition: Settings.h:166
trklet::Tracklet::tpseed
int tpseed()
Definition: Tracklet.cc:120
trklet::Settings::nrbitsstub
unsigned int nrbitsstub(unsigned int layerdisk) const
Definition: Settings.h:66
trklet::HistBase
Definition: HistBase.h:16
trklet::Tracklet
Definition: Tracklet.h:28
trklet::IMATH_TrackletCalculatorOverlap::rinv_final
VarAdjustK rinv_final
Definition: IMATH_TrackletCalculatorOverlap.h:206
trklet::IMATH_TrackletCalculator::rproj0
VarDef rproj0
Definition: IMATH_TrackletCalculator.h:164
trklet::IMATH_TrackletCalculatorOverlap::r1
VarDef r1
Definition: IMATH_TrackletCalculatorOverlap.h:148
trklet::IMATH_TrackletCalculatorDisk::r2
VarDef r2
Definition: IMATH_TrackletCalculatorDisk.h:144
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
testProducerWithPsetDescEmpty_cfi.z2
z2
Definition: testProducerWithPsetDescEmpty_cfi.py:41
LEDCalibrationChannels.iphi
iphi
Definition: LEDCalibrationChannels.py:64
trklet::Globals::ITC_L1F1
IMATH_TrackletCalculatorOverlap * ITC_L1F1()
Definition: Globals.h:62
trklet::IMATH_TrackletCalculator::phiL_0_final
VarAdjustK phiL_0_final
Definition: IMATH_TrackletCalculator.h:268
trklet::IMATH_TrackletCalculatorOverlap::phiD_1_final
VarAdjustK phiD_1_final
Definition: IMATH_TrackletCalculatorOverlap.h:304
trklet::Globals::ITC_L2B1
IMATH_TrackletCalculatorOverlap * ITC_L2B1()
Definition: Globals.h:65
trklet::IMATH_TrackletCalculatorDisk::der_rD_final
VarAdjustK der_rD_final
Definition: IMATH_TrackletCalculatorDisk.h:328
trklet::IMATH_TrackletCalculatorOverlap::phiD_0_final
VarAdjustK phiD_0_final
Definition: IMATH_TrackletCalculatorOverlap.h:303
trklet::IMATH_TrackletCalculatorOverlap::zL_2_final
VarAdjustK zL_2_final
Definition: IMATH_TrackletCalculatorOverlap.h:273
trklet::IMATH_TrackletCalculatorDisk::zL_1_final
VarAdjustK zL_1_final
Definition: IMATH_TrackletCalculatorDisk.h:266
trklet::VarDef::set_ival
void set_ival(int ival)
Definition: imath.h:501
trklet::Tracklet::fpgazproj
const FPGAWord & fpgazproj(int layer) const
Definition: Tracklet.h:99
trklet::Settings::projdisks
unsigned int projdisks(unsigned int iSeed, unsigned int i) const
Definition: Settings.h:120
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
trklet::Stub::disk
const FPGAWord & disk() const
Definition: Stub.h:57
trklet::Stub::z
const FPGAWord & z() const
Definition: Stub.h:50
trklet::Stub
Definition: Stub.h:16
trklet::L1TStub::isPSmodule
unsigned int isPSmodule() const
Definition: L1TStub.h:93
trklet::IMATH_TrackletCalculatorOverlap
Definition: IMATH_TrackletCalculatorOverlap.h:19
trklet::IMATH_TrackletCalculatorDisk::zproj1
VarDef zproj1
Definition: IMATH_TrackletCalculatorDisk.h:156
trklet::IMATH_TrackletCalculator::t_final
VarAdjustKR t_final
Definition: IMATH_TrackletCalculator.h:219
trklet::IMATH_TrackletCalculatorOverlap::phiL_1_final
VarAdjustK phiL_1_final
Definition: IMATH_TrackletCalculatorOverlap.h:251
trklet::IMATH_TrackletCalculatorOverlap::rproj1
VarDef rproj1
Definition: IMATH_TrackletCalculatorOverlap.h:157
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
trklet::IMATH_TrackletCalculator::zL_3_final
VarAdjustKR zL_3_final
Definition: IMATH_TrackletCalculator.h:296
trklet::Globals::ITC_F3F4
IMATH_TrackletCalculatorDisk * ITC_F3F4()
Definition: Globals.h:58
trklet::IMATH_TrackletCalculatorDisk::zproj2
VarDef zproj2
Definition: IMATH_TrackletCalculatorDisk.h:157
trklet::DiskProjection
Definition: DiskProjection.h:12
trklet::TrackletCalculatorBase::exactproj
void exactproj(double rproj, double rinv, double phi0, double t, double z0, double &phiproj, double &zproj, double &phider, double &zder)
Definition: TrackletCalculatorBase.cc:177
trklet::N_LAYER
constexpr int N_LAYER
Definition: Settings.h:19
summarizeEdmComparisonLogfiles.success
success
Definition: summarizeEdmComparisonLogfiles.py:115
trklet::IMATH_TrackletCalculatorOverlap::der_phiL_final
VarAdjustK der_phiL_final
Definition: IMATH_TrackletCalculatorOverlap.h:254
trklet::IMATH_TrackletCalculator::z2
VarDef z2
Definition: IMATH_TrackletCalculator.h:158
trklet::Settings::zmean
double zmean(unsigned int iDisk) const
Definition: Settings.h:131
trklet::IMATH_TrackletCalculator::phiL_3_final
VarAdjustK phiL_3_final
Definition: IMATH_TrackletCalculator.h:271
trklet::IMATH_TrackletCalculatorOverlap::zL_1_final
VarAdjustK zL_1_final
Definition: IMATH_TrackletCalculatorOverlap.h:272
trklet::TrackletCalculatorBase::addDiskProj
void addDiskProj(Tracklet *tracklet, int disk)
Definition: TrackletCalculatorBase.cc:213
trklet::IMATH_TrackletCalculator::rD_1_final
VarAdjustK rD_1_final
Definition: IMATH_TrackletCalculator.h:375
trklet::IMATH_TrackletCalculator::der_phiD_final
VarAdjustK der_phiD_final
Definition: IMATH_TrackletCalculator.h:339
trklet::IMATH_TrackletCalculatorDisk::zL_2_final
VarAdjustK zL_2_final
Definition: IMATH_TrackletCalculatorDisk.h:267
trklet::IMATH_TrackletCalculatorOverlap::z0_final
VarAdjustK z0_final
Definition: IMATH_TrackletCalculatorOverlap.h:210
trklet::IMATH_TrackletCalculatorDisk::phiD_2_final
VarAdjustK phiD_2_final
Definition: IMATH_TrackletCalculatorDisk.h:295
trklet::TrackletCalculatorBase::diskSeeding
bool diskSeeding(const Stub *innerFPGAStub, const L1TStub *innerStub, const Stub *outerFPGAStub, const L1TStub *outerStub)
Definition: TrackletCalculatorBase.cc:753
trklet::TrackletCalculatorBase::addLayerProj
bool addLayerProj(Tracklet *tracklet, int layer)
Definition: TrackletCalculatorBase.cc:230
trklet::Stub::isDisk
bool isDisk() const
Definition: Stub.h:61
HLTMuonOfflineAnalyzer_cfi.z0
z0
Definition: HLTMuonOfflineAnalyzer_cfi.py:98
trklet::IMATH_TrackletCalculatorDisk::phiD_0_final
VarAdjustK phiD_0_final
Definition: IMATH_TrackletCalculatorDisk.h:293
trklet::IMATH_TrackletCalculatorDisk::der_zL_final
VarAdjustK der_zL_final
Definition: IMATH_TrackletCalculatorDisk.h:269
TrackletCalculatorBase.h
trklet::IMATH_TrackletCalculatorOverlap::zproj0
VarDef zproj0
Definition: IMATH_TrackletCalculatorOverlap.h:160
trklet::IMATH_TrackletCalculator::phiD_4_final
VarAdjustK phiD_4_final
Definition: IMATH_TrackletCalculator.h:335
trklet::VarBase::local_passes
bool local_passes() const
Definition: imath.cc:315
cmsswSequenceInfo.tp
tp
Definition: cmsswSequenceInfo.py:17
trklet::FPGAWord
Definition: FPGAWord.h:9
trklet::Settings::warnNoMem
bool warnNoMem() const
Definition: Settings.h:149
trklet::IMATH_TrackletCalculatorOverlap::phiD_2_final
VarAdjustK phiD_2_final
Definition: IMATH_TrackletCalculatorOverlap.h:305
trklet::Settings::rmaxdisk
double rmaxdisk() const
Definition: Settings.h:102
trklet::IMATH_TrackletCalculatorOverlap::zproj2
VarDef zproj2
Definition: IMATH_TrackletCalculatorOverlap.h:162
trklet::TrackletCalculatorBase::exacttracklet
void exacttracklet(double r1, double z1, double phi1, double r2, double z2, double phi2, double, double &rinv, double &phi0, double &t, double &z0, double phiproj[N_LAYER - 2], double zproj[N_LAYER - 2], double phider[N_LAYER - 2], double zder[N_LAYER - 2], double phiprojdisk[N_DISK], double rprojdisk[N_DISK], double phiderdisk[N_DISK], double rderdisk[N_DISK])
Definition: TrackletCalculatorBase.cc:23
trklet::FPGAWord::atExtreme
bool atExtreme() const
Definition: FPGAWord.cc:79
trklet::Settings::kphi0pars
double kphi0pars() const
Definition: Settings.h:331
trklet::IMATH_TrackletCalculator::zL_0_final
VarAdjustKR zL_0_final
Definition: IMATH_TrackletCalculator.h:293
trklet::IMATH_TrackletCalculatorOverlap::z2
VarDef z2
Definition: IMATH_TrackletCalculatorOverlap.h:151
trklet::TrackletCalculatorBase::layerdisk1_
unsigned int layerdisk1_
Definition: TrackletCalculatorBase.h:129
trklet::rinv
double rinv(double phi1, double phi2, double r1, double r2)
Definition: Util.h:167
trklet::IMATH_TrackletCalculator::z0_final
VarAdjustKR z0_final
Definition: IMATH_TrackletCalculator.h:220
trklet::TrackletCalculatorBase::layerdisk2_
unsigned int layerdisk2_
Definition: TrackletCalculatorBase.h:130
trklet::IMATH_TrackletCalculator::zproj0
VarDef zproj0
Definition: IMATH_TrackletCalculator.h:169
trklet::VarBase::K
double K() const
Definition: imath.h:246
trklet::IMATH_TrackletCalculatorDisk::rproj0
VarDef rproj0
Definition: IMATH_TrackletCalculatorDisk.h:151
trklet::IMATH_TrackletCalculator::phi1
VarDef phi1
Definition: IMATH_TrackletCalculator.h:161
trklet::IMATH_TrackletCalculatorDisk::z0_final
VarAdjustK z0_final
Definition: IMATH_TrackletCalculatorDisk.h:204
trklet::IMATH_TrackletCalculatorOverlap::rproj2
VarDef rproj2
Definition: IMATH_TrackletCalculatorOverlap.h:158
HistBase.h
trklet::IMATH_TrackletCalculatorDisk::phiD_1_final
VarAdjustK phiD_1_final
Definition: IMATH_TrackletCalculatorDisk.h:294
trklet::Globals::ITC_L1L2
IMATH_TrackletCalculator * ITC_L1L2()
Definition: Globals.h:52
trklet::IMATH_TrackletCalculatorOverlap::zproj3
VarDef zproj3
Definition: IMATH_TrackletCalculatorOverlap.h:163
trklet::Settings::projlayers
unsigned int projlayers(unsigned int iSeed, unsigned int i) const
Definition: Settings.h:119
trklet::IMATH_TrackletCalculatorDisk::rD_1_final
VarAdjustK rD_1_final
Definition: IMATH_TrackletCalculatorDisk.h:325
trklet::Settings::nzbitsstub
unsigned int nzbitsstub(unsigned int layerdisk) const
Definition: Settings.h:64
trklet::IMATH_TrackletCalculator::phiL_1_final
VarAdjustK phiL_1_final
Definition: IMATH_TrackletCalculator.h:269
trklet::TrackletCalculatorBase::overlapSeeding
bool overlapSeeding(const Stub *innerFPGAStub, const L1TStub *innerStub, const Stub *outerFPGAStub, const L1TStub *outerStub)
Definition: TrackletCalculatorBase.cc:1108
trklet::IMATH_TrackletCalculator::phi0_final
VarAdjustK phi0_final
Definition: IMATH_TrackletCalculator.h:218
diffTwoXMLs.r2
r2
Definition: diffTwoXMLs.py:73
groupFilesInBlocks.fout
fout
Definition: groupFilesInBlocks.py:162
trklet::Stub::rapprox
double rapprox() const
Definition: Stub.cc:209
trklet::TrackletCalculatorBase::trackletprojdisks_
std::vector< std::vector< TrackletProjectionsMemory * > > trackletprojdisks_
Definition: TrackletCalculatorBase.h:143
trklet::IMATH_TrackletCalculator::phi2
VarDef phi2
Definition: IMATH_TrackletCalculator.h:162
Globals.h
trklet::Settings::nphibitsstub
unsigned int nphibitsstub(unsigned int layerdisk) const
Definition: Settings.h:65
trklet::IMATH_TrackletCalculatorDisk::rinv_final
VarAdjustK rinv_final
Definition: IMATH_TrackletCalculatorDisk.h:200
trklet::IMATH_TrackletCalculator::phiD_3_final
VarAdjustK phiD_3_final
Definition: IMATH_TrackletCalculator.h:334
trklet::IMATH_TrackletCalculatorDisk::r1
VarDef r1
Definition: IMATH_TrackletCalculatorDisk.h:143
trklet::IMATH_TrackletCalculator::rD_4_final
VarAdjustK rD_4_final
Definition: IMATH_TrackletCalculator.h:378
trklet::IMATH_TrackletCalculator::rD_2_final
VarAdjustK rD_2_final
Definition: IMATH_TrackletCalculator.h:376
trklet::Globals::ITC_F1F2
IMATH_TrackletCalculatorDisk * ITC_F1F2()
Definition: Globals.h:57
trklet::IMATH_TrackletCalculatorOverlap::rD_2_final
VarAdjustK rD_2_final
Definition: IMATH_TrackletCalculatorOverlap.h:342
trklet::IMATH_TrackletCalculatorOverlap::phiL_2_final
VarAdjustK phiL_2_final
Definition: IMATH_TrackletCalculatorOverlap.h:252
trklet::Tracklet::setTrackletIndex
void setTrackletIndex(int index)
Definition: Tracklet.cc:844
trklet::IMATH_TrackletCalculatorDisk::der_phiL_final
VarAdjustK der_phiL_final
Definition: IMATH_TrackletCalculatorDisk.h:248
trklet::VarBase::ival
long int ival() const
Definition: imath.h:213
trklet::IMATH_TrackletCalculatorDisk::rproj2
VarDef rproj2
Definition: IMATH_TrackletCalculatorDisk.h:153
trklet::IMATH_TrackletCalculatorOverlap::r2
VarDef r2
Definition: IMATH_TrackletCalculatorOverlap.h:149
trklet
Definition: AllProjectionsMemory.h:9
trklet::TrackletCalculatorBase::exacttrackletOverlap
void exacttrackletOverlap(double r1, double z1, double phi1, double r2, double z2, double phi2, double, double &rinv, double &phi0, double &t, double &z0, double phiprojLayer[N_PSLAYER], double zprojLayer[N_PSLAYER], double phiderLayer[N_PSLAYER], double zderLayer[N_PSLAYER], double phiproj[N_DISK - 2], double rproj[N_DISK - 2], double phider[N_DISK - 2], double rder[N_DISK - 2])
Definition: TrackletCalculatorBase.cc:129
trklet::FPGAWord::value
int value() const
Definition: FPGAWord.h:24
trklet::IMATH_TrackletCalculator::phiD_1_final
VarAdjustK phiD_1_final
Definition: IMATH_TrackletCalculator.h:332
trklet::Globals::ITC_L5L6
IMATH_TrackletCalculator * ITC_L5L6()
Definition: Globals.h:55
IMATH_TrackletCalculator.h
trklet::Globals::ITC_B3B4
IMATH_TrackletCalculatorDisk * ITC_B3B4()
Definition: Globals.h:60
trklet::IMATH_TrackletCalculator::zproj1
VarDef zproj1
Definition: IMATH_TrackletCalculator.h:170
angle0to2pi::make0To2pi
constexpr valType make0To2pi(valType angle)
Definition: deltaPhi.h:67
trklet::IMATH_TrackletCalculatorDisk::z1
VarDef z1
Definition: IMATH_TrackletCalculatorDisk.h:145
trklet::VarDef::set_fval
void set_fval(double fval)
Definition: imath.h:493
trklet::IMATH_TrackletCalculator::der_zL_final
VarAdjustK der_zL_final
Definition: IMATH_TrackletCalculator.h:298
trklet::IMATH_TrackletCalculator::rinv_final
VarAdjustK rinv_final
Definition: IMATH_TrackletCalculator.h:216
trklet::Settings::writeMonitorData
bool writeMonitorData(std::string module) const
Definition: Settings.h:86
trklet::TrackletCalculatorBase::TCIndex_
int TCIndex_
Definition: TrackletCalculatorBase.h:132
std
Definition: JetResolutionObject.h:76
trklet::IMATH_TrackletCalculatorOverlap::rproj0
VarDef rproj0
Definition: IMATH_TrackletCalculatorOverlap.h:156
trklet::TrackletParametersMemory::nTracklets
unsigned int nTracklets() const
Definition: TrackletParametersMemory.h:25
trklet::Settings::usephicritapprox
bool usephicritapprox() const
Definition: Settings.h:194
trklet::Tracklet::setTCIndex
void setTCIndex(int index)
Definition: Tracklet.h:496
trklet::ProcessBase
Definition: ProcessBase.h:12
trklet::IMATH_TrackletCalculator::zproj3
VarDef zproj3
Definition: IMATH_TrackletCalculator.h:172
trklet::Globals::ITC_L1B1
IMATH_TrackletCalculatorOverlap * ITC_L1B1()
Definition: Globals.h:63
trklet::IMATH_TrackletCalculatorDisk::phiL_0_final
VarAdjustK phiL_0_final
Definition: IMATH_TrackletCalculatorDisk.h:244
trklet::Globals::ofstream
std::ofstream & ofstream(std::string fname)
Definition: Globals.cc:44
trklet::TrackletProjectionsMemory::addProj
void addProj(Tracklet *tracklet)
Definition: TrackletProjectionsMemory.cc:16
trklet::IMATH_TrackletCalculator::zproj2
VarDef zproj2
Definition: IMATH_TrackletCalculator.h:171
trklet::IMATH_TrackletCalculator::phiL_2_final
VarAdjustK phiL_2_final
Definition: IMATH_TrackletCalculator.h:270
diffTwoXMLs.r1
r1
Definition: diffTwoXMLs.py:53
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
trklet::IMATH_TrackletCalculator::zproj4
VarDef zproj4
Definition: IMATH_TrackletCalculator.h:173
trklet::Settings::useapprox
bool useapprox() const
Definition: Settings.h:193
Exception
Definition: hltDiff.cc:246
trklet::TrackletCalculatorBase::inSector
bool inSector(int iphi0, int irinv, double phi0approx, double rinvapprox)
Definition: TrackletCalculatorBase.cc:306
trklet::ProcessBase::getName
std::string const & getName() const
Definition: ProcessBase.h:22
trklet::IMATH_TrackletCalculator::r2
VarDef r2
Definition: IMATH_TrackletCalculator.h:156
trklet::Tracklet::fpgaphiproj
const FPGAWord & fpgaphiproj(int layer) const
Definition: Tracklet.h:104
trklet::Settings::kz
double kz() const
Definition: Settings.h:249
trklet::ProcessBase::name_
std::string name_
Definition: ProcessBase.h:38
trklet::IMATH_TrackletCalculatorOverlap::t_final
VarAdjustK t_final
Definition: IMATH_TrackletCalculatorOverlap.h:209
trklet::Settings::nallstubs
unsigned int nallstubs(unsigned int layerdisk) const
Definition: Settings.h:84
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
trklet::ProcessBase::iSector_
unsigned int iSector_
Definition: ProcessBase.h:39
trklet::Settings::debugTracklet
bool debugTracklet() const
Definition: Settings.h:146
trklet::IMATH_TrackletCalculator::rproj3
VarDef rproj3
Definition: IMATH_TrackletCalculator.h:167
trklet::Settings::rcrit
double rcrit() const
Definition: Settings.h:236
trklet::Globals::ITC_L3L4
IMATH_TrackletCalculator * ITC_L3L4()
Definition: Globals.h:54
trklet::IMATH_TrackletCalculator
Definition: IMATH_TrackletCalculator.h:20
Exception.h
trklet::IMATH_TrackletCalculatorOverlap::phi0_final
VarAdjustK phi0_final
Definition: IMATH_TrackletCalculatorOverlap.h:208
trklet::IMATH_TrackletCalculatorOverlap::phi1
VarDef phi1
Definition: IMATH_TrackletCalculatorOverlap.h:153
trklet::TrackletCalculatorBase::exactprojdisk
void exactprojdisk(double zproj, double rinv, double phi0, double t, double z0, double &phiproj, double &rproj, double &phider, double &rder)
Definition: TrackletCalculatorBase.cc:193
trklet::Tracklet::fpgarprojdisk
const FPGAWord & fpgarprojdisk(int disk) const
Definition: Tracklet.h:276
trklet::L1TStub::sigmaz
double sigmaz() const
Definition: L1TStub.h:71
trklet::L1TStub::r
double r() const
Definition: L1TStub.h:57
trklet::IMATH_TrackletCalculatorOverlap::phiL_0_final
VarAdjustK phiL_0_final
Definition: IMATH_TrackletCalculatorOverlap.h:250
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
trklet::Globals::ITC_L2F1
IMATH_TrackletCalculatorOverlap * ITC_L2F1()
Definition: Globals.h:64
trklet::ProcessBase::phimin_
double phimin_
Definition: ProcessBase.h:41
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
trklet::IMATH_TrackletCalculator::rproj2
VarDef rproj2
Definition: IMATH_TrackletCalculator.h:166
trklet::IMATH_TrackletCalculator::phiD_2_final
VarAdjustK phiD_2_final
Definition: IMATH_TrackletCalculator.h:333
trklet::IMATH_TrackletCalculator::der_phiL_final
VarAdjustK der_phiL_final
Definition: IMATH_TrackletCalculator.h:273
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
trklet::Stub::layer
const FPGAWord & layer() const
Definition: Stub.h:56
trklet::Tracklet::validProjDisk
bool validProjDisk(int disk) const
Definition: Tracklet.h:206
trklet::TrackletCalculatorBase::exacttrackletdisk
void exacttrackletdisk(double r1, double z1, double phi1, double r2, double z2, double phi2, double, double &rinv, double &phi0, double &t, double &z0, double phiprojLayer[N_PSLAYER], double zprojLayer[N_PSLAYER], double phiderLayer[N_PSLAYER], double zderLayer[N_PSLAYER], double phiproj[N_DISK - 2], double rproj[N_DISK - 2], double phider[N_DISK - 2], double rder[N_DISK - 2])
Definition: TrackletCalculatorBase.cc:76
trklet::Globals::histograms
HistBase *& histograms()
Definition: Globals.h:40
trklet::IMATH_TrackletCalculatorDisk::rproj1
VarDef rproj1
Definition: IMATH_TrackletCalculatorDisk.h:152
trklet::IMATH_TrackletCalculatorOverlap::der_rD_final
VarAdjustK der_rD_final
Definition: IMATH_TrackletCalculatorOverlap.h:345
trklet::ProcessBase::globals_
Globals * globals_
Definition: ProcessBase.h:45
trklet::IMATH_TrackletCalculatorDisk::phiL_1_final
VarAdjustK phiL_1_final
Definition: IMATH_TrackletCalculatorDisk.h:245
submitPVValidationJobs.t
string t
Definition: submitPVValidationJobs.py:644
IMATH_TrackletCalculatorOverlap.h
edm::Log
Definition: MessageLogger.h:70
trklet::IMATH_TrackletCalculator::der_rD_final
VarAdjustK der_rD_final
Definition: IMATH_TrackletCalculator.h:380
trklet::ProcessBase::phimax_
double phimax_
Definition: ProcessBase.h:42
keep
const int keep
Definition: GenParticlePruner.cc:48
Stub.h
Tracklet.h
trklet::Stub::zapprox
double zapprox() const
Definition: Stub.cc:217
trklet::IMATH_TrackletCalculatorOverlap::phi2
VarDef phi2
Definition: IMATH_TrackletCalculatorOverlap.h:154
trklet::L1TStub::phi
double phi() const
Definition: L1TStub.h:62
trklet::IMATH_TrackletCalculatorOverlap::der_phiD_final
VarAdjustK der_phiD_final
Definition: IMATH_TrackletCalculatorOverlap.h:310
IMATH_TrackletCalculatorDisk.h
trklet::TrackletCalculatorBase::goodTrackPars
bool goodTrackPars(bool goodrinv, bool goodz0)
Definition: TrackletCalculatorBase.cc:289
trklet::IMATH_TrackletCalculatorDisk::phiL_2_final
VarAdjustK phiL_2_final
Definition: IMATH_TrackletCalculatorDisk.h:246
reco::reduceRange
constexpr T reduceRange(T x)
Definition: deltaPhi.h:18
trklet::Settings::zlength
double zlength() const
Definition: Settings.h:101
trklet::DiskProjection::init
void init(Settings const &settings, int projdisk, double zproj, int iphiproj, int irproj, int iphider, int irder, double phiproj, double rproj, double phiprojder, double rprojder, double phiprojapprox, double rprojapprox, double phiprojderapprox, double rprojderapprox)
Definition: DiskProjection.cc:12