CMS 3D CMS Logo

TrackletCalculatorDisplaced.cc
Go to the documentation of this file.
8 
12 
13 using namespace std;
14 using namespace trklet;
15 
16 TrackletCalculatorDisplaced::TrackletCalculatorDisplaced(string name,
17  Settings const& settings,
18  Globals* global,
19  unsigned int iSector)
20  : ProcessBase(name, settings, global, iSector) {
21  for (unsigned int ilayer = 0; ilayer < N_LAYER; ilayer++) {
22  vector<TrackletProjectionsMemory*> tmp(settings.nallstubs(ilayer), nullptr);
23  trackletprojlayers_.push_back(tmp);
24  }
25 
26  for (unsigned int idisk = 0; idisk < N_DISK; idisk++) {
27  vector<TrackletProjectionsMemory*> tmp(settings.nallstubs(idisk + N_LAYER), nullptr);
28  trackletprojdisks_.push_back(tmp);
29  }
30 
31  layer_ = 0;
32  disk_ = 0;
33 
34  string name1 = name.substr(1); //this is to correct for "TCD" having one more letter then "TC"
35  if (name1[3] == 'L')
36  layer_ = name1[4] - '0';
37  if (name1[3] == 'D')
38  disk_ = name1[4] - '0';
39 
40  // set TC index
41  int iSeed = -1;
42 
43  int iTC = name1[9] - 'A';
44 
45  if (name1.substr(3, 6) == "L3L4L2")
46  iSeed = 8;
47  else if (name1.substr(3, 6) == "L5L6L4")
48  iSeed = 9;
49  else if (name1.substr(3, 6) == "L2L3D1")
50  iSeed = 10;
51  else if (name1.substr(3, 6) == "D1D2L2")
52  iSeed = 11;
53 
54  assert(iSeed != -1);
55 
56  TCIndex_ = (iSeed << 4) + iTC;
57  assert(TCIndex_ >= 128 && TCIndex_ < 191);
58 
59  assert((layer_ != 0) || (disk_ != 0));
60 
61  toR_.clear();
62  toZ_.clear();
63 
64  if (iSeed == 8 || iSeed == 9) {
65  if (layer_ == 3) {
66  rzmeanInv_[0] = 1.0 / settings_.rmean(2 - 1);
67  rzmeanInv_[1] = 1.0 / settings_.rmean(3 - 1);
68  rzmeanInv_[2] = 1.0 / settings_.rmean(4 - 1);
69 
70  rproj_[0] = settings_.rmean(0);
71  rproj_[1] = settings_.rmean(4);
72  rproj_[2] = settings_.rmean(5);
73  lproj_[0] = 1;
74  lproj_[1] = 5;
75  lproj_[2] = 6;
76 
77  dproj_[0] = 1;
78  dproj_[1] = 2;
79  dproj_[2] = 0;
80  toZ_.push_back(settings_.zmean(0));
81  toZ_.push_back(settings_.zmean(1));
82  }
83  if (layer_ == 5) {
84  rzmeanInv_[0] = 1.0 / settings_.rmean(4 - 1);
85  rzmeanInv_[1] = 1.0 / settings_.rmean(5 - 1);
86  rzmeanInv_[2] = 1.0 / settings_.rmean(6 - 1);
87 
88  rproj_[0] = settings_.rmean(0);
89  rproj_[1] = settings_.rmean(1);
90  rproj_[2] = settings_.rmean(2);
91  lproj_[0] = 1;
92  lproj_[1] = 2;
93  lproj_[2] = 3;
94 
95  dproj_[0] = 0;
96  dproj_[1] = 0;
97  dproj_[2] = 0;
98  }
99  for (unsigned int i = 0; i < N_LAYER - 3; ++i)
100  toR_.push_back(rproj_[i]);
101  }
102 
103  if (iSeed == 10 || iSeed == 11) {
104  if (layer_ == 2) {
105  rzmeanInv_[0] = 1.0 / settings_.rmean(2 - 1);
106  rzmeanInv_[1] = 1.0 / settings_.rmean(3 - 1);
107  rzmeanInv_[2] = 1.0 / settings_.zmean(1 - 1);
108 
109  rproj_[0] = settings_.rmean(0);
110  lproj_[0] = 1;
111  lproj_[1] = -1;
112  lproj_[2] = -1;
113 
114  zproj_[0] = settings_.zmean(1);
115  zproj_[1] = settings_.zmean(2);
116  zproj_[2] = settings_.zmean(3);
117  dproj_[0] = 2;
118  dproj_[1] = 3;
119  dproj_[2] = 4;
120  }
121  if (disk_ == 1) {
122  rzmeanInv_[0] = 1.0 / settings_.rmean(2 - 1);
123  rzmeanInv_[1] = 1.0 / settings_.zmean(1 - 1);
124  rzmeanInv_[2] = 1.0 / settings_.zmean(2 - 1);
125 
126  rproj_[0] = settings_.rmean(0);
127  lproj_[0] = 1;
128  lproj_[1] = -1;
129  lproj_[2] = -1;
130 
131  zproj_[0] = settings_.zmean(2);
132  zproj_[1] = settings_.zmean(3);
133  zproj_[2] = settings_.zmean(4);
134  dproj_[0] = 3;
135  dproj_[1] = 4;
136  dproj_[2] = 5;
137  }
138  toR_.push_back(settings_.rmean(0));
139  for (unsigned int i = 0; i < N_DISK - 2; ++i)
140  toZ_.push_back(zproj_[i]);
141  }
142 }
143 
145  outputProj = dynamic_cast<TrackletProjectionsMemory*>(memory);
146  assert(outputProj != nullptr);
147 }
148 
150  if (settings_.writetrace()) {
151  edm::LogVerbatim("Tracklet") << "In " << name_ << " adding output to " << memory->getName() << " to output "
152  << output;
153  }
154  if (output == "trackpar") {
155  auto* tmp = dynamic_cast<TrackletParametersMemory*>(memory);
156  assert(tmp != nullptr);
157  trackletpars_ = tmp;
158  return;
159  }
160 
161  if (output.substr(0, 7) == "projout") {
162  //output is on the form 'projoutL2PHIC' or 'projoutD3PHIB'
163  auto* tmp = dynamic_cast<TrackletProjectionsMemory*>(memory);
164  assert(tmp != nullptr);
165 
166  unsigned int layerdisk = output[8] - '1'; //layer or disk counting from 0
167  unsigned int phiregion = output[12] - 'A'; //phiregion counting from 0
168 
169  if (output[7] == 'L') {
170  assert(layerdisk < N_LAYER);
171  assert(phiregion < trackletprojlayers_[layerdisk].size());
172  //check that phiregion not already initialized
173  assert(trackletprojlayers_[layerdisk][phiregion] == nullptr);
174  trackletprojlayers_[layerdisk][phiregion] = tmp;
175  return;
176  }
177 
178  if (output[7] == 'D') {
179  assert(layerdisk < N_DISK);
180  assert(phiregion < trackletprojdisks_[layerdisk].size());
181  //check that phiregion not already initialized
182  assert(trackletprojdisks_[layerdisk][phiregion] == nullptr);
183  trackletprojdisks_[layerdisk][phiregion] = tmp;
184  return;
185  }
186  }
187 
188  throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " Could not find output : " << output;
189 }
190 
192  if (settings_.writetrace()) {
193  edm::LogVerbatim("Tracklet") << "In " << name_ << " adding input from " << memory->getName() << " to input "
194  << input;
195  }
196  if (input == "thirdallstubin") {
197  auto* tmp = dynamic_cast<AllStubsMemory*>(memory);
198  assert(tmp != nullptr);
199  innerallstubs_.push_back(tmp);
200  return;
201  }
202  if (input == "firstallstubin") {
203  auto* tmp = dynamic_cast<AllStubsMemory*>(memory);
204  assert(tmp != nullptr);
205  middleallstubs_.push_back(tmp);
206  return;
207  }
208  if (input == "secondallstubin") {
209  auto* tmp = dynamic_cast<AllStubsMemory*>(memory);
210  assert(tmp != nullptr);
211  outerallstubs_.push_back(tmp);
212  return;
213  }
214  if (input.find("stubtriplet") == 0) {
215  auto* tmp = dynamic_cast<StubTripletsMemory*>(memory);
216  assert(tmp != nullptr);
217  stubtriplets_.push_back(tmp);
218  return;
219  }
220  throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " Could not find input : " << input;
221 }
222 
224  unsigned int countall = 0;
225  unsigned int countsel = 0;
226 
227  for (auto& stubtriplet : stubtriplets_) {
229  edm::LogVerbatim("Tracklet") << "Will break on too many tracklets in " << getName();
230  break;
231  }
232  for (unsigned int i = 0; i < stubtriplet->nStubTriplets(); i++) {
233  countall++;
234 
235  const Stub* innerFPGAStub = stubtriplet->getFPGAStub1(i);
236  const L1TStub* innerStub = innerFPGAStub->l1tstub();
237 
238  const Stub* middleFPGAStub = stubtriplet->getFPGAStub2(i);
239  const L1TStub* middleStub = middleFPGAStub->l1tstub();
240 
241  const Stub* outerFPGAStub = stubtriplet->getFPGAStub3(i);
242  const L1TStub* outerStub = outerFPGAStub->l1tstub();
243 
244  if (settings_.debugTracklet())
245  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced execute " << getName() << "[" << iSector_ << "]";
246 
247  if (innerFPGAStub->isBarrel() && middleFPGAStub->isBarrel() && outerFPGAStub->isBarrel()) {
248  //barrel+barrel seeding
249  bool accept = LLLSeeding(innerFPGAStub, innerStub, middleFPGAStub, middleStub, outerFPGAStub, outerStub);
250  if (accept)
251  countsel++;
252  } else if (innerFPGAStub->isDisk() && middleFPGAStub->isDisk() && outerFPGAStub->isDisk()) {
253  throw cms::Exception("LogicError") << __FILE__ << " " << __LINE__ << " Invalid seeding!";
254  } else {
255  //layer+disk seeding
256  if (innerFPGAStub->isBarrel() && middleFPGAStub->isDisk() && outerFPGAStub->isDisk()) { //D1D2L2
257  bool accept = DDLSeeding(innerFPGAStub, innerStub, middleFPGAStub, middleStub, outerFPGAStub, outerStub);
258  if (accept)
259  countsel++;
260  } else if (innerFPGAStub->isDisk() && middleFPGAStub->isBarrel() && outerFPGAStub->isBarrel()) { //L2L3D1
261  bool accept = LLDSeeding(innerFPGAStub, innerStub, middleFPGAStub, middleStub, outerFPGAStub, outerStub);
262  if (accept)
263  countsel++;
264  } else {
265  throw cms::Exception("LogicError") << __FILE__ << " " << __LINE__ << " Invalid seeding!";
266  }
267  }
268 
270  edm::LogVerbatim("Tracklet") << "Will break on number of tracklets in " << getName();
271  break;
272  }
273 
274  if (countall >= settings_.maxStep("TC")) {
275  if (settings_.debugTracklet())
276  edm::LogVerbatim("Tracklet") << "Will break on MAXTC 1";
277  break;
278  }
279  if (settings_.debugTracklet())
280  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced execute done";
281  }
282  if (countall >= settings_.maxStep("TC")) {
283  if (settings_.debugTracklet())
284  edm::LogVerbatim("Tracklet") << "Will break on MAXTC 2";
285  break;
286  }
287  }
288 
289  if (settings_.writeMonitorData("TPD")) {
290  globals_->ofstream("trackletcalculatordisplaced.txt") << getName() << " " << countall << " " << countsel << endl;
291  }
292 }
293 
295  FPGAWord fpgar = tracklet->fpgarprojdisk(disk);
296 
297  if (fpgar.value() * settings_.krprojshiftdisk() < settings_.rmindiskvm())
298  return;
299  if (fpgar.value() * settings_.krprojshiftdisk() > settings_.rmaxdisk())
300  return;
301 
302  FPGAWord fpgaphi = tracklet->fpgaphiprojdisk(disk);
303 
304  int iphivmRaw = fpgaphi.value() >> (fpgaphi.nbits() - 5);
305  int iphi = iphivmRaw / (32 / settings_.nallstubs(abs(disk) + N_DISK));
306 
307  addProjectionDisk(disk, iphi, trackletprojdisks_[abs(disk) - 1][iphi], tracklet);
308 }
309 
311  assert(layer > 0);
312 
313  FPGAWord fpgaz = tracklet->fpgazproj(layer);
314  FPGAWord fpgaphi = tracklet->fpgaphiproj(layer);
315 
316  if (fpgaz.atExtreme())
317  return false;
318 
319  if (std::abs(fpgaz.value() * settings_.kz()) > settings_.zlength())
320  return false;
321 
322  int iphivmRaw = fpgaphi.value() >> (fpgaphi.nbits() - 5);
323  int iphi = iphivmRaw / (32 / settings_.nallstubs(layer - 1));
324 
326 
327  return true;
328 }
329 
331  int iphi,
332  TrackletProjectionsMemory* trackletprojs,
333  Tracklet* tracklet) {
334  if (trackletprojs == nullptr) {
335  if (settings_.warnNoMem()) {
336  edm::LogVerbatim("Tracklet") << "No projection memory exists in " << getName() << " for layer = " << layer
337  << " iphi = " << iphi + 1;
338  }
339  return;
340  }
341  assert(trackletprojs != nullptr);
342  trackletprojs->addProj(tracklet);
343 }
344 
346  int iphi,
347  TrackletProjectionsMemory* trackletprojs,
348  Tracklet* tracklet) {
349  if (trackletprojs == nullptr) {
350  if (layer_ == 3 && abs(disk) == 3)
351  return; //L3L4 projections to D3 are not used.
352  if (settings_.warnNoMem()) {
353  edm::LogVerbatim("Tracklet") << "No projection memory exists in " << getName() << " for disk = " << abs(disk)
354  << " iphi = " << iphi + 1;
355  }
356  return;
357  }
358  assert(trackletprojs != nullptr);
359  if (settings_.debugTracklet())
360  edm::LogVerbatim("Tracklet") << getName() << " adding projection to " << trackletprojs->getName();
361  trackletprojs->addProj(tracklet);
362 }
363 
365  const L1TStub* innerStub,
366  const Stub* middleFPGAStub,
367  const L1TStub* middleStub,
368  const Stub* outerFPGAStub,
369  const L1TStub* outerStub) {
370  if (settings_.debugTracklet())
371  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced " << getName() << " " << layer_
372  << " trying stub triplet in layer (L L L): " << innerFPGAStub->layer().value() << " "
373  << middleFPGAStub->layer().value() << " " << outerFPGAStub->layer().value();
374 
375  assert(outerFPGAStub->isBarrel());
376 
377  double r1 = innerStub->r();
378  double z1 = innerStub->z();
379  double phi1 = innerStub->phi();
380 
381  double r2 = middleStub->r();
382  double z2 = middleStub->z();
383  double phi2 = middleStub->phi();
384 
385  double r3 = outerStub->r();
386  double z3 = outerStub->z();
387  double phi3 = outerStub->phi();
388 
389  int take3 = 0;
390  if (layer_ == 5)
391  take3 = 1;
392  unsigned ndisks = 0;
393 
394  double rinv, phi0, d0, t, z0;
395 
396  LayerProjection layerprojs[N_LAYER - 2];
397  DiskProjection diskprojs[N_DISK];
398 
399  double phiproj[N_LAYER - 2], zproj[N_LAYER - 2], phider[N_LAYER - 2], zder[N_LAYER - 2];
400  double phiprojdisk[N_DISK], rprojdisk[N_DISK], phiderdisk[N_DISK], rderdisk[N_DISK];
401 
403  z1,
404  phi1,
405  r2,
406  z2,
407  phi2,
408  r3,
409  z3,
410  phi3,
411  take3,
412  rinv,
413  phi0,
414  d0,
415  t,
416  z0,
417  phiproj,
418  zproj,
419  phiprojdisk,
420  rprojdisk,
421  phider,
422  zder,
423  phiderdisk,
424  rderdisk);
425  if (settings_.debugTracklet())
426  edm::LogVerbatim("Tracklet") << __LINE__ << ":" << __FILE__ << "LLL Exact values " << innerFPGAStub->isBarrel()
427  << middleFPGAStub->isBarrel() << outerFPGAStub->isBarrel() << " " << phi1 << ", " << z1
428  << ", " << r1 << ", " << phi2 << ", " << z2 << ", " << r2 << ", " << phi3 << ", " << z3
429  << ", " << r3 << endl;
430 
431  if (settings_.useapprox()) {
432  phi1 = innerFPGAStub->phiapprox(phimin_, phimax_);
433  z1 = innerFPGAStub->zapprox();
434  r1 = innerFPGAStub->rapprox();
435 
436  phi2 = middleFPGAStub->phiapprox(phimin_, phimax_);
437  z2 = middleFPGAStub->zapprox();
438  r2 = middleFPGAStub->rapprox();
439 
440  phi3 = outerFPGAStub->phiapprox(phimin_, phimax_);
441  z3 = outerFPGAStub->zapprox();
442  r3 = outerFPGAStub->rapprox();
443  }
444 
445  if (settings_.debugTracklet())
446  edm::LogVerbatim("Tracklet") << __LINE__ << ":" << __FILE__ << "LLL Approx values " << innerFPGAStub->isBarrel()
447  << middleFPGAStub->isBarrel() << outerFPGAStub->isBarrel() << " " << phi1 << ", " << z1
448  << ", " << r1 << ", " << phi2 << ", " << z2 << ", " << r2 << ", " << phi3 << ", " << z3
449  << ", " << r3 << endl;
450 
451  double rinvapprox, phi0approx, d0approx, tapprox, z0approx;
452  double phiprojapprox[N_LAYER - 2], zprojapprox[N_LAYER - 2], phiderapprox[N_LAYER - 2], zderapprox[N_LAYER - 2];
453  double phiprojdiskapprox[N_DISK], rprojdiskapprox[N_DISK];
454  double phiderdiskapprox[N_DISK], rderdiskapprox[N_DISK];
455 
456  //TODO: implement the actual integer calculation
457  if (settings_.useapprox()) {
459  z1,
460  phi1,
461  r2,
462  z2,
463  phi2,
464  r3,
465  z3,
466  phi3,
467  take3,
468  ndisks,
469  rinvapprox,
470  phi0approx,
471  d0approx,
472  tapprox,
473  z0approx,
474  phiprojapprox,
475  zprojapprox,
476  phiderapprox,
477  zderapprox,
478  phiprojdiskapprox,
479  rprojdiskapprox,
480  phiderdiskapprox,
481  rderdiskapprox);
482  } else {
483  rinvapprox = rinv;
484  phi0approx = phi0;
485  d0approx = d0;
486  tapprox = t;
487  z0approx = z0;
488 
489  for (unsigned int i = 0; i < toR_.size(); ++i) {
490  phiprojapprox[i] = phiproj[i];
491  zprojapprox[i] = zproj[i];
492  phiderapprox[i] = phider[i];
493  zderapprox[i] = zder[i];
494  }
495 
496  for (unsigned int i = 0; i < toZ_.size(); ++i) {
497  phiprojdiskapprox[i] = phiprojdisk[i];
498  rprojdiskapprox[i] = rprojdisk[i];
499  phiderdiskapprox[i] = phiderdisk[i];
500  rderdiskapprox[i] = rderdisk[i];
501  }
502  }
503 
504  //store the approcximate results
505 
506  if (settings_.debugTracklet()) {
507  edm::LogVerbatim("Tracklet") << "rinvapprox: " << rinvapprox << " rinv: " << rinv << endl;
508  edm::LogVerbatim("Tracklet") << "phi0approx: " << phi0approx << " phi0: " << phi0 << endl;
509  edm::LogVerbatim("Tracklet") << "d0approx: " << d0approx << " d0: " << d0 << endl;
510  edm::LogVerbatim("Tracklet") << "tapprox: " << tapprox << " t: " << t << endl;
511  edm::LogVerbatim("Tracklet") << "z0approx: " << z0approx << " z0: " << z0 << endl;
512  }
513 
514  for (unsigned int i = 0; i < toR_.size(); ++i) {
515  if (settings_.debugTracklet()) {
516  edm::LogVerbatim("Tracklet") << "phiprojapprox[" << i << "]: " << phiprojapprox[i] << " phiproj[" << i
517  << "]: " << phiproj[i] << endl;
518  edm::LogVerbatim("Tracklet") << "zprojapprox[" << i << "]: " << zprojapprox[i] << " zproj[" << i
519  << "]: " << zproj[i] << endl;
520  edm::LogVerbatim("Tracklet") << "phiderapprox[" << i << "]: " << phiderapprox[i] << " phider[" << i
521  << "]: " << phider[i] << endl;
522  edm::LogVerbatim("Tracklet") << "zderapprox[" << i << "]: " << zderapprox[i] << " zder[" << i << "]: " << zder[i]
523  << endl;
524  }
525  }
526 
527  for (unsigned int i = 0; i < toZ_.size(); ++i) {
528  if (settings_.debugTracklet()) {
529  edm::LogVerbatim("Tracklet") << "phiprojdiskapprox[" << i << "]: " << phiprojdiskapprox[i] << " phiprojdisk[" << i
530  << "]: " << phiprojdisk[i] << endl;
531  edm::LogVerbatim("Tracklet") << "rprojdiskapprox[" << i << "]: " << rprojdiskapprox[i] << " rprojdisk[" << i
532  << "]: " << rprojdisk[i] << endl;
533  edm::LogVerbatim("Tracklet") << "phiderdiskapprox[" << i << "]: " << phiderdiskapprox[i] << " phiderdisk[" << i
534  << "]: " << phiderdisk[i] << endl;
535  edm::LogVerbatim("Tracklet") << "rderdiskapprox[" << i << "]: " << rderdiskapprox[i] << " rderdisk[" << i
536  << "]: " << rderdisk[i] << endl;
537  }
538  }
539 
540  //now binary
541  double krinv = settings_.kphi1() / settings_.kr() * pow(2, settings_.rinv_shift()),
542  kphi0 = settings_.kphi1() * pow(2, settings_.phi0_shift()),
543  kt = settings_.kz() / settings_.kr() * pow(2, settings_.t_shift()),
544  kz0 = settings_.kz() * pow(2, settings_.z0_shift()),
545  kphiproj = settings_.kphi1() * pow(2, settings_.SS_phiL_shift()),
546  kphider = settings_.kphi1() / settings_.kr() * pow(2, settings_.SS_phiderL_shift()),
547  kzproj = settings_.kz() * pow(2, settings_.PS_zL_shift()),
548  kzder = settings_.kz() / settings_.kr() * pow(2, settings_.PS_zderL_shift()),
549  kphiprojdisk = settings_.kphi1() * pow(2, settings_.SS_phiD_shift()),
550  kphiderdisk = settings_.kphi1() / settings_.kr() * pow(2, settings_.SS_phiderD_shift()),
551  krprojdisk = settings_.kr() * pow(2, settings_.PS_rD_shift()),
552  krderdisk = settings_.kr() / settings_.kz() * pow(2, settings_.PS_rderD_shift());
553 
554  int irinv, iphi0, id0, it, iz0;
555  int iphiproj[N_LAYER - 2], izproj[N_LAYER - 2], iphider[N_LAYER - 2], izder[N_LAYER - 2];
556  int iphiprojdisk[N_DISK], irprojdisk[N_DISK], iphiderdisk[N_DISK], irderdisk[N_DISK];
557 
558  //store the binary results
559  irinv = rinvapprox / krinv;
560  iphi0 = phi0approx / kphi0;
561  id0 = d0approx / settings_.kd0();
562  it = tapprox / kt;
563  iz0 = z0approx / kz0;
564 
565  bool success = true;
566  if (std::abs(rinvapprox) > settings_.rinvcut()) {
567  if (settings_.debugTracklet())
568  edm::LogVerbatim("Tracklet") << "TrackletCalculator::LLL Seeding irinv too large: " << rinvapprox << "(" << irinv
569  << ")";
570  success = false;
571  }
572  if (std::abs(z0approx) > settings_.disp_z0cut()) {
573  if (settings_.debugTracklet())
574  edm::LogVerbatim("Tracklet") << "Failed tracklet z0 cut " << z0approx << " in layer " << layer_;
575  success = false;
576  }
577  if (std::abs(d0approx) > settings_.maxd0()) {
578  if (settings_.debugTracklet())
579  edm::LogVerbatim("Tracklet") << "Failed tracklet d0 cut " << d0approx;
580  success = false;
581  }
582 
583  if (!success)
584  return false;
585 
586  double phicritapprox = phi0approx - asin(0.5 * settings_.rcrit() * rinvapprox);
587  int phicrit = iphi0 - 2 * irinv;
588 
589  int iphicritmincut = settings_.phicritminmc() / globals_->ITC_L1L2()->phi0_final.K();
590  int iphicritmaxcut = settings_.phicritmaxmc() / globals_->ITC_L1L2()->phi0_final.K();
591 
592  bool keepapprox = (phicritapprox > settings_.phicritminmc()) && (phicritapprox < settings_.phicritmaxmc()),
593  keep = (phicrit > iphicritmincut) && (phicrit < iphicritmaxcut);
594 
595  if (settings_.debugTracklet())
596  if (keep && !keepapprox)
597  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced::LLLSeeding tracklet kept with exact phicrit cut "
598  "but not approximate, phicritapprox: "
599  << phicritapprox;
600  if (settings_.usephicritapprox()) {
601  if (!keepapprox)
602  return false;
603  } else {
604  if (!keep)
605  return false;
606  }
607 
608  for (unsigned int i = 0; i < toR_.size(); ++i) {
609  iphiproj[i] = phiprojapprox[i] / kphiproj;
610  izproj[i] = zprojapprox[i] / kzproj;
611 
612  iphider[i] = phiderapprox[i] / kphider;
613  izder[i] = zderapprox[i] / kzder;
614 
615  //check that z projection is in range
616  if (izproj[i] < -(1 << (settings_.nzbitsstub(0) - 1)))
617  continue;
618  if (izproj[i] >= (1 << (settings_.nzbitsstub(0) - 1)))
619  continue;
620 
621  //check that phi projection is in range
622  if (iphiproj[i] >= (1 << settings_.nphibitsstub(N_LAYER - 1)) - 1)
623  continue;
624  if (iphiproj[i] <= 0)
625  continue;
626 
627  //adjust number of bits for phi and z projection
628  if (rproj_[i] < settings_.rPS2S()) {
629  iphiproj[i] >>= (settings_.nphibitsstub(N_LAYER - 1) - settings_.nphibitsstub(0));
630  if (iphiproj[i] >= (1 << settings_.nphibitsstub(0)) - 1)
631  iphiproj[i] = (1 << settings_.nphibitsstub(0)) - 2; //-2 not to hit atExtreme
632  } else {
633  izproj[i] >>= (settings_.nzbitsstub(0) - settings_.nzbitsstub(N_LAYER - 1));
634  }
635 
636  if (rproj_[i] < settings_.rPS2S()) {
637  if (iphider[i] < -(1 << (settings_.nbitsphiprojderL123() - 1))) {
638  iphider[i] = -(1 << (settings_.nbitsphiprojderL123() - 1));
639  }
640  if (iphider[i] >= (1 << (settings_.nbitsphiprojderL123() - 1))) {
641  iphider[i] = (1 << (settings_.nbitsphiprojderL123() - 1)) - 1;
642  }
643  } else {
644  if (iphider[i] < -(1 << (settings_.nbitsphiprojderL456() - 1))) {
645  iphider[i] = -(1 << (settings_.nbitsphiprojderL456() - 1));
646  }
647  if (iphider[i] >= (1 << (settings_.nbitsphiprojderL456() - 1))) {
648  iphider[i] = (1 << (settings_.nbitsphiprojderL456() - 1)) - 1;
649  }
650  }
651 
652  layerprojs[i].init(settings_,
653  lproj_[i],
654  rproj_[i],
655  iphiproj[i],
656  izproj[i],
657  iphider[i],
658  izder[i],
659  phiproj[i],
660  zproj[i],
661  phider[i],
662  zder[i],
663  phiprojapprox[i],
664  zprojapprox[i],
665  phiderapprox[i],
666  zderapprox[i]);
667  }
668 
669  if (std::abs(it * kt) > 1.0) {
670  for (unsigned int i = 0; i < toZ_.size(); ++i) {
671  iphiprojdisk[i] = phiprojdiskapprox[i] / kphiprojdisk;
672  irprojdisk[i] = rprojdiskapprox[i] / krprojdisk;
673 
674  iphiderdisk[i] = phiderdiskapprox[i] / kphiderdisk;
675  irderdisk[i] = rderdiskapprox[i] / krderdisk;
676 
677  //check phi projection in range
678  if (iphiprojdisk[i] <= 0)
679  continue;
680  if (iphiprojdisk[i] >= (1 << settings_.nphibitsstub(0)) - 1)
681  continue;
682 
683  //check r projection in range
684  if (rprojdiskapprox[i] < settings_.rmindisk() || rprojdiskapprox[i] > settings_.rmaxdisk())
685  continue;
686 
687  diskprojs[i].init(settings_,
688  i + 1,
689  rproj_[i],
690  iphiprojdisk[i],
691  irprojdisk[i],
692  iphiderdisk[i],
693  irderdisk[i],
694  phiprojdisk[i],
695  rprojdisk[i],
696  phiderdisk[i],
697  rderdisk[i],
698  phiprojdiskapprox[i],
699  rprojdiskapprox[i],
700  phiderdisk[i],
701  rderdisk[i]);
702  }
703  }
704 
705  if (settings_.writeMonitorData("TrackletPars")) {
706  globals_->ofstream("trackletpars.txt")
707  << layer_ << " , " << rinv << " , " << rinvapprox << " , " << phi0 << " , " << phi0approx << " , " << t << " , "
708  << tapprox << " , " << z0 << " , " << z0approx << " , " << d0 << " , " << d0approx << endl;
709  }
710 
711  Tracklet* tracklet = new Tracklet(settings_,
712  innerStub,
713  middleStub,
714  outerStub,
715  innerFPGAStub,
716  middleFPGAStub,
717  outerFPGAStub,
718  rinv,
719  phi0,
720  d0,
721  z0,
722  t,
723  rinvapprox,
724  phi0approx,
725  d0approx,
726  z0approx,
727  tapprox,
728  irinv,
729  iphi0,
730  id0,
731  iz0,
732  it,
733  layerprojs,
734  diskprojs,
735  false);
736 
737  if (settings_.debugTracklet())
738  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced " << getName()
739  << " Found LLL tracklet in sector = " << iSector_ << " phi0 = " << phi0;
740 
742  tracklet->setTCIndex(TCIndex_);
743 
744  if (settings_.writeMonitorData("Seeds")) {
745  ofstream fout("seeds.txt", ofstream::app);
746  fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_" << iSector_ << " " << tracklet->getISeed() << endl;
747  fout.close();
748  }
749  trackletpars_->addTracklet(tracklet);
750 
751  bool addL5 = false;
752  bool addL6 = false;
753  for (unsigned int j = 0; j < toR_.size(); j++) {
754  bool added = false;
755  if (settings_.debugTracklet())
756  edm::LogVerbatim("Tracklet") << "adding layer projection " << j << "/" << toR_.size() << " " << lproj_[j];
757  if (tracklet->validProj(lproj_[j])) {
758  added = addLayerProj(tracklet, lproj_[j]);
759  if (added && lproj_[j] == 5)
760  addL5 = true;
761  if (added && lproj_[j] == 6)
762  addL6 = true;
763  }
764  }
765 
766  for (unsigned int j = 0; j < toZ_.size(); j++) {
767  int disk = dproj_[j];
768  if (disk == 0)
769  continue;
770  if (disk == 2 && addL5)
771  continue;
772  if (disk == 1 && addL6)
773  continue;
774  if (it < 0)
775  disk = -disk;
776  if (settings_.debugTracklet())
777  edm::LogVerbatim("Tracklet") << "adding disk projection " << j << "/" << toZ_.size() << " " << disk;
778  if (tracklet->validProjDisk(abs(disk))) {
779  addDiskProj(tracklet, disk);
780  }
781  }
782 
783  return true;
784 }
785 
787  const L1TStub* innerStub,
788  const Stub* middleFPGAStub,
789  const L1TStub* middleStub,
790  const Stub* outerFPGAStub,
791  const L1TStub* outerStub) {
792  if (settings_.debugTracklet())
793  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced " << getName() << " " << layer_
794  << " trying stub triplet in (L2 D1 D2): " << innerFPGAStub->layer().value() << " "
795  << middleFPGAStub->disk().value() << " " << outerFPGAStub->disk().value();
796 
797  int take3 = 1; //D1D2L2
798  unsigned ndisks = 2;
799 
800  double r1 = innerStub->r();
801  double z1 = innerStub->z();
802  double phi1 = innerStub->phi();
803 
804  double r2 = middleStub->r();
805  double z2 = middleStub->z();
806  double phi2 = middleStub->phi();
807 
808  double r3 = outerStub->r();
809  double z3 = outerStub->z();
810  double phi3 = outerStub->phi();
811 
812  double rinv, phi0, d0, t, z0;
813 
814  double phiproj[N_LAYER - 2], zproj[N_LAYER - 2], phider[N_LAYER - 2], zder[N_LAYER - 2];
815  double phiprojdisk[N_DISK], rprojdisk[N_DISK], phiderdisk[N_DISK], rderdisk[N_DISK];
816 
818  z1,
819  phi1,
820  r2,
821  z2,
822  phi2,
823  r3,
824  z3,
825  phi3,
826  take3,
827  rinv,
828  phi0,
829  d0,
830  t,
831  z0,
832  phiproj,
833  zproj,
834  phiprojdisk,
835  rprojdisk,
836  phider,
837  zder,
838  phiderdisk,
839  rderdisk);
840 
841  if (settings_.debugTracklet())
842  edm::LogVerbatim("Tracklet") << __LINE__ << ":" << __FILE__ << " DLL Exact values " << innerFPGAStub->isBarrel()
843  << middleFPGAStub->isBarrel() << outerFPGAStub->isBarrel() << " " << phi1 << ", " << z1
844  << ", " << r1 << ", " << phi2 << ", " << z2 << ", " << r2 << ", " << phi3 << ", " << z3
845  << ", " << r3 << endl;
846 
847  if (settings_.useapprox()) {
848  phi1 = innerFPGAStub->phiapprox(phimin_, phimax_);
849  z1 = innerFPGAStub->zapprox();
850  r1 = innerFPGAStub->rapprox();
851 
852  phi2 = middleFPGAStub->phiapprox(phimin_, phimax_);
853  z2 = middleFPGAStub->zapprox();
854  r2 = middleFPGAStub->rapprox();
855 
856  phi3 = outerFPGAStub->phiapprox(phimin_, phimax_);
857  z3 = outerFPGAStub->zapprox();
858  r3 = outerFPGAStub->rapprox();
859  }
860 
861  if (settings_.debugTracklet())
862  edm::LogVerbatim("Tracklet") << __LINE__ << ":" << __FILE__ << "DLL Approx values " << innerFPGAStub->isBarrel()
863  << middleFPGAStub->isBarrel() << outerFPGAStub->isBarrel() << " " << phi1 << ", " << z1
864  << ", " << r1 << ", " << phi2 << ", " << z2 << ", " << r2 << ", " << phi3 << ", " << z3
865  << ", " << r3 << endl;
866 
867  double rinvapprox, phi0approx, d0approx, tapprox, z0approx;
868  double phiprojapprox[N_LAYER - 2], zprojapprox[N_LAYER - 2], phiderapprox[N_LAYER - 2], zderapprox[N_LAYER - 2];
869  double phiprojdiskapprox[N_DISK], rprojdiskapprox[N_DISK];
870  double phiderdiskapprox[N_DISK], rderdiskapprox[N_DISK];
871 
872  //TODO: implement the actual integer calculation
873  if (settings_.useapprox()) {
875  z1,
876  phi1,
877  r2,
878  z2,
879  phi2,
880  r3,
881  z3,
882  phi3,
883  take3,
884  ndisks,
885  rinvapprox,
886  phi0approx,
887  d0approx,
888  tapprox,
889  z0approx,
890  phiprojapprox,
891  zprojapprox,
892  phiderapprox,
893  zderapprox,
894  phiprojdiskapprox,
895  rprojdiskapprox,
896  phiderdiskapprox,
897  rderdiskapprox);
898  } else {
899  rinvapprox = rinv;
900  phi0approx = phi0;
901  d0approx = d0;
902  tapprox = t;
903  z0approx = z0;
904 
905  for (unsigned int i = 0; i < toR_.size(); ++i) {
906  phiprojapprox[i] = phiproj[i];
907  zprojapprox[i] = zproj[i];
908  phiderapprox[i] = phider[i];
909  zderapprox[i] = zder[i];
910  }
911 
912  for (unsigned int i = 0; i < toZ_.size(); ++i) {
913  phiprojdiskapprox[i] = phiprojdisk[i];
914  rprojdiskapprox[i] = rprojdisk[i];
915  phiderdiskapprox[i] = phiderdisk[i];
916  rderdiskapprox[i] = rderdisk[i];
917  }
918  }
919 
920  //store the approcximate results
921  if (settings_.debugTracklet()) {
922  edm::LogVerbatim("Tracklet") << "rinvapprox: " << rinvapprox << " rinv: " << rinv << endl;
923  edm::LogVerbatim("Tracklet") << "phi0approx: " << phi0approx << " phi0: " << phi0 << endl;
924  edm::LogVerbatim("Tracklet") << "d0approx: " << d0approx << " d0: " << d0 << endl;
925  edm::LogVerbatim("Tracklet") << "tapprox: " << tapprox << " t: " << t << endl;
926  edm::LogVerbatim("Tracklet") << "z0approx: " << z0approx << " z0: " << z0 << endl;
927  }
928 
929  for (unsigned int i = 0; i < toR_.size(); ++i) {
930  if (settings_.debugTracklet()) {
931  edm::LogVerbatim("Tracklet") << "phiprojapprox[" << i << "]: " << phiprojapprox[i] << " phiproj[" << i
932  << "]: " << phiproj[i] << endl;
933  edm::LogVerbatim("Tracklet") << "zprojapprox[" << i << "]: " << zprojapprox[i] << " zproj[" << i
934  << "]: " << zproj[i] << endl;
935  edm::LogVerbatim("Tracklet") << "phiderapprox[" << i << "]: " << phiderapprox[i] << " phider[" << i
936  << "]: " << phider[i] << endl;
937  edm::LogVerbatim("Tracklet") << "zderapprox[" << i << "]: " << zderapprox[i] << " zder[" << i << "]: " << zder[i]
938  << endl;
939  }
940  }
941 
942  for (unsigned int i = 0; i < toZ_.size(); ++i) {
943  if (settings_.debugTracklet()) {
944  edm::LogVerbatim("Tracklet") << "phiprojdiskapprox[" << i << "]: " << phiprojdiskapprox[i] << " phiprojdisk[" << i
945  << "]: " << phiprojdisk[i] << endl;
946  edm::LogVerbatim("Tracklet") << "rprojdiskapprox[" << i << "]: " << rprojdiskapprox[i] << " rprojdisk[" << i
947  << "]: " << rprojdisk[i] << endl;
948  edm::LogVerbatim("Tracklet") << "phiderdiskapprox[" << i << "]: " << phiderdiskapprox[i] << " phiderdisk[" << i
949  << "]: " << phiderdisk[i] << endl;
950  edm::LogVerbatim("Tracklet") << "rderdiskapprox[" << i << "]: " << rderdiskapprox[i] << " rderdisk[" << i
951  << "]: " << rderdisk[i] << endl;
952  }
953  }
954 
955  //now binary
956  double krinv = settings_.kphi1() / settings_.kr() * pow(2, settings_.rinv_shift()),
957  kphi0 = settings_.kphi1() * pow(2, settings_.phi0_shift()),
958  kt = settings_.kz() / settings_.kr() * pow(2, settings_.t_shift()),
959  kz0 = settings_.kz() * pow(2, settings_.z0_shift()),
960  kphiproj = settings_.kphi1() * pow(2, settings_.SS_phiL_shift()),
961  kphider = settings_.kphi1() / settings_.kr() * pow(2, settings_.SS_phiderL_shift()),
962  kzproj = settings_.kz() * pow(2, settings_.PS_zL_shift()),
963  kzder = settings_.kz() / settings_.kr() * pow(2, settings_.PS_zderL_shift()),
964  kphiprojdisk = settings_.kphi1() * pow(2, settings_.SS_phiD_shift()),
965  kphiderdisk = settings_.kphi1() / settings_.kr() * pow(2, settings_.SS_phiderD_shift()),
966  krprojdisk = settings_.kr() * pow(2, settings_.PS_rD_shift()),
967  krderdisk = settings_.kr() / settings_.kz() * pow(2, settings_.PS_rderD_shift());
968 
969  int irinv, iphi0, id0, it, iz0;
970  int iphiproj[N_LAYER - 2], izproj[N_LAYER - 2], iphider[N_LAYER - 2], izder[N_LAYER - 2];
971  int iphiprojdisk[N_DISK], irprojdisk[N_DISK], iphiderdisk[N_DISK], irderdisk[N_DISK];
972 
973  //store the binary results
974  irinv = rinvapprox / krinv;
975  iphi0 = phi0approx / kphi0;
976  id0 = d0approx / settings_.kd0();
977  it = tapprox / kt;
978  iz0 = z0approx / kz0;
979 
980  bool success = true;
981  if (std::abs(rinvapprox) > settings_.rinvcut()) {
982  if (settings_.debugTracklet())
983  edm::LogVerbatim("Tracklet") << "TrackletCalculator::DDL Seeding irinv too large: " << rinvapprox << "(" << irinv
984  << ")";
985  success = false;
986  }
987  if (std::abs(z0approx) > settings_.disp_z0cut()) {
988  if (settings_.debugTracklet())
989  edm::LogVerbatim("Tracklet") << "Failed tracklet z0 cut " << z0approx;
990  success = false;
991  }
992  if (std::abs(d0approx) > settings_.maxd0()) {
993  if (settings_.debugTracklet())
994  edm::LogVerbatim("Tracklet") << "Failed tracklet d0 cut " << d0approx;
995  success = false;
996  }
997 
998  if (!success)
999  return false;
1000 
1001  double phicritapprox = phi0approx - asin(0.5 * settings_.rcrit() * rinvapprox);
1002  int phicrit = iphi0 - 2 * irinv;
1003 
1004  int iphicritmincut = settings_.phicritminmc() / globals_->ITC_L1L2()->phi0_final.K();
1005  int iphicritmaxcut = settings_.phicritmaxmc() / globals_->ITC_L1L2()->phi0_final.K();
1006 
1007  bool keepapprox = (phicritapprox > settings_.phicritminmc()) && (phicritapprox < settings_.phicritmaxmc()),
1008  keep = (phicrit > iphicritmincut) && (phicrit < iphicritmaxcut);
1009 
1010  if (settings_.debugTracklet())
1011  if (keep && !keepapprox)
1012  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced::DDLSeeding tracklet kept with exact phicrit cut "
1013  "but not approximate, phicritapprox: "
1014  << phicritapprox;
1015  if (settings_.usephicritapprox()) {
1016  if (!keepapprox)
1017  return false;
1018  } else {
1019  if (!keep)
1020  return false;
1021  }
1022 
1023  LayerProjection layerprojs[N_LAYER - 2];
1024  DiskProjection diskprojs[N_DISK];
1025 
1026  for (unsigned int i = 0; i < toR_.size(); ++i) {
1027  iphiproj[i] = phiprojapprox[i] / kphiproj;
1028  izproj[i] = zprojapprox[i] / kzproj;
1029 
1030  iphider[i] = phiderapprox[i] / kphider;
1031  izder[i] = zderapprox[i] / kzder;
1032 
1033  //check that z projection in range
1034  if (izproj[i] < -(1 << (settings_.nzbitsstub(0) - 1)))
1035  continue;
1036  if (izproj[i] >= (1 << (settings_.nzbitsstub(0) - 1)))
1037  continue;
1038 
1039  //check that phi projection in range
1040  if (iphiproj[i] >= (1 << settings_.nphibitsstub(5)) - 1)
1041  continue;
1042  if (iphiproj[i] <= 0)
1043  continue;
1044 
1045  if (rproj_[i] < settings_.rPS2S()) {
1046  iphiproj[i] >>= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
1047  } else {
1048  izproj[i] >>= (settings_.nzbitsstub(0) - settings_.nzbitsstub(5));
1049  }
1050 
1051  if (rproj_[i] < settings_.rPS2S()) {
1052  if (iphider[i] < -(1 << (settings_.nbitsphiprojderL123() - 1)))
1053  iphider[i] = -(1 << (settings_.nbitsphiprojderL123() - 1));
1054  if (iphider[i] >= (1 << (settings_.nbitsphiprojderL123() - 1)))
1055  iphider[i] = (1 << (settings_.nbitsphiprojderL123() - 1)) - 1;
1056  } else {
1057  if (iphider[i] < -(1 << (settings_.nbitsphiprojderL456() - 1)))
1058  iphider[i] = -(1 << (settings_.nbitsphiprojderL456() - 1));
1059  if (iphider[i] >= (1 << (settings_.nbitsphiprojderL456() - 1)))
1060  iphider[i] = (1 << (settings_.nbitsphiprojderL456() - 1)) - 1;
1061  }
1062 
1063  layerprojs[i].init(settings_,
1064  lproj_[i],
1065  rproj_[i],
1066  iphiproj[i],
1067  izproj[i],
1068  iphider[i],
1069  izder[i],
1070  phiproj[i],
1071  zproj[i],
1072  phider[i],
1073  zder[i],
1074  phiprojapprox[i],
1075  zprojapprox[i],
1076  phiderapprox[i],
1077  zderapprox[i]);
1078  }
1079 
1080  if (std::abs(it * kt) > 1.0) {
1081  for (unsigned int i = 0; i < toZ_.size(); ++i) {
1082  iphiprojdisk[i] = phiprojdiskapprox[i] / kphiprojdisk;
1083  irprojdisk[i] = rprojdiskapprox[i] / krprojdisk;
1084 
1085  iphiderdisk[i] = phiderdiskapprox[i] / kphiderdisk;
1086  irderdisk[i] = rderdiskapprox[i] / krderdisk;
1087 
1088  if (iphiprojdisk[i] <= 0)
1089  continue;
1090  if (iphiprojdisk[i] >= (1 << settings_.nphibitsstub(0)) - 1)
1091  continue;
1092 
1093  if (irprojdisk[i] < settings_.rmindisk() / krprojdisk || irprojdisk[i] > settings_.rmaxdisk() / krprojdisk)
1094  continue;
1095 
1096  diskprojs[i].init(settings_,
1097  i + 1,
1098  rproj_[i],
1099  iphiprojdisk[i],
1100  irprojdisk[i],
1101  iphiderdisk[i],
1102  irderdisk[i],
1103  phiprojdisk[i],
1104  rprojdisk[i],
1105  phiderdisk[i],
1106  rderdisk[i],
1107  phiprojdiskapprox[i],
1108  rprojdiskapprox[i],
1109  phiderdisk[i],
1110  rderdisk[i]);
1111  }
1112  }
1113 
1114  if (settings_.writeMonitorData("TrackletPars")) {
1115  globals_->ofstream("trackletpars.txt")
1116  << layer_ << " , " << rinv << " , " << rinvapprox << " , " << phi0 << " , " << phi0approx << " , " << t << " , "
1117  << tapprox << " , " << z0 << " , " << z0approx << " , " << d0 << " , " << d0approx << endl;
1118  }
1119 
1120  Tracklet* tracklet = new Tracklet(settings_,
1121  innerStub,
1122  middleStub,
1123  outerStub,
1124  innerFPGAStub,
1125  middleFPGAStub,
1126  outerFPGAStub,
1127  rinv,
1128  phi0,
1129  d0,
1130  z0,
1131  t,
1132  rinvapprox,
1133  phi0approx,
1134  d0approx,
1135  z0approx,
1136  tapprox,
1137  irinv,
1138  iphi0,
1139  id0,
1140  iz0,
1141  it,
1142  layerprojs,
1143  diskprojs,
1144  true);
1145 
1146  if (settings_.debugTracklet())
1147  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced " << getName()
1148  << " Found DDL tracklet in sector = " << iSector_ << " phi0 = " << phi0;
1149 
1151  tracklet->setTCIndex(TCIndex_);
1152 
1153  if (settings_.writeMonitorData("Seeds")) {
1154  ofstream fout("seeds.txt", ofstream::app);
1155  fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_" << iSector_ << " " << tracklet->getISeed() << endl;
1156  fout.close();
1157  }
1158  trackletpars_->addTracklet(tracklet);
1159 
1160  for (unsigned int j = 0; j < toR_.size(); j++) {
1161  if (settings_.debugTracklet())
1162  edm::LogVerbatim("Tracklet") << "adding layer projection " << j << "/" << toR_.size() << " " << lproj_[j] << " "
1163  << tracklet->validProj(lproj_[j]);
1164  if (tracklet->validProj(lproj_[j])) {
1165  addLayerProj(tracklet, lproj_[j]);
1166  }
1167  }
1168 
1169  for (unsigned int j = 0; j < toZ_.size(); j++) {
1170  int disk = dproj_[j];
1171  if (disk == 0)
1172  continue;
1173  if (it < 0)
1174  disk = -disk;
1175  if (settings_.debugTracklet())
1176  edm::LogVerbatim("Tracklet") << "adding disk projection " << j << "/" << toZ_.size() << " " << disk << " "
1177  << tracklet->validProjDisk(abs(disk));
1178  if (tracklet->validProjDisk(abs(disk))) {
1179  addDiskProj(tracklet, disk);
1180  }
1181  }
1182 
1183  return true;
1184 }
1185 
1187  const L1TStub* innerStub,
1188  const Stub* middleFPGAStub,
1189  const L1TStub* middleStub,
1190  const Stub* outerFPGAStub,
1191  const L1TStub* outerStub) {
1192  if (settings_.debugTracklet())
1193  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced " << getName() << " " << layer_
1194  << " trying stub triplet in (L2L3D1): " << middleFPGAStub->layer().value() << " "
1195  << outerFPGAStub->layer().value() << " " << innerFPGAStub->disk().value();
1196 
1197  int take3 = 0; //L2L3D1
1198  unsigned ndisks = 1;
1199 
1200  double r3 = innerStub->r();
1201  double z3 = innerStub->z();
1202  double phi3 = innerStub->phi();
1203 
1204  double r1 = middleStub->r();
1205  double z1 = middleStub->z();
1206  double phi1 = middleStub->phi();
1207 
1208  double r2 = outerStub->r();
1209  double z2 = outerStub->z();
1210  double phi2 = outerStub->phi();
1211 
1212  double rinv, phi0, d0, t, z0;
1213 
1214  double phiproj[N_LAYER - 2], zproj[N_LAYER - 2], phider[N_LAYER - 2], zder[N_LAYER - 2];
1215  double phiprojdisk[N_DISK], rprojdisk[N_DISK], phiderdisk[N_DISK], rderdisk[N_DISK];
1216 
1217  exacttracklet(r1,
1218  z1,
1219  phi1,
1220  r2,
1221  z2,
1222  phi2,
1223  r3,
1224  z3,
1225  phi3,
1226  take3,
1227  rinv,
1228  phi0,
1229  d0,
1230  t,
1231  z0,
1232  phiproj,
1233  zproj,
1234  phiprojdisk,
1235  rprojdisk,
1236  phider,
1237  zder,
1238  phiderdisk,
1239  rderdisk);
1240 
1241  if (settings_.debugTracklet())
1242  edm::LogVerbatim("Tracklet") << __LINE__ << ":" << __FILE__ << "LLD Exact values " << innerFPGAStub->isBarrel()
1243  << middleFPGAStub->isBarrel() << outerFPGAStub->isBarrel() << " " << phi3 << ", " << z3
1244  << ", " << r3 << ", " << phi1 << ", " << z1 << ", " << r1 << ", " << phi2 << ", " << z2
1245  << ", " << r2 << endl;
1246 
1247  if (settings_.useapprox()) {
1248  phi3 = innerFPGAStub->phiapprox(phimin_, phimax_);
1249  z3 = innerFPGAStub->zapprox();
1250  r3 = innerFPGAStub->rapprox();
1251 
1252  phi1 = middleFPGAStub->phiapprox(phimin_, phimax_);
1253  z1 = middleFPGAStub->zapprox();
1254  r1 = middleFPGAStub->rapprox();
1255 
1256  phi2 = outerFPGAStub->phiapprox(phimin_, phimax_);
1257  z2 = outerFPGAStub->zapprox();
1258  r2 = outerFPGAStub->rapprox();
1259  }
1260 
1261  if (settings_.debugTracklet())
1262  edm::LogVerbatim("Tracklet") << __LINE__ << ":" << __FILE__ << "LLD approx values " << innerFPGAStub->isBarrel()
1263  << middleFPGAStub->isBarrel() << outerFPGAStub->isBarrel() << " " << phi3 << ", " << z3
1264  << ", " << r3 << ", " << phi1 << ", " << z1 << ", " << r1 << ", " << phi2 << ", " << z2
1265  << ", " << r2 << endl;
1266 
1267  double rinvapprox, phi0approx, d0approx, tapprox, z0approx;
1268  double phiprojapprox[N_LAYER - 2], zprojapprox[N_LAYER - 2], phiderapprox[N_LAYER - 2], zderapprox[N_LAYER - 2];
1269  double phiprojdiskapprox[N_DISK], rprojdiskapprox[N_DISK];
1270  double phiderdiskapprox[N_DISK], rderdiskapprox[N_DISK];
1271 
1272  //TODO: implement the actual integer calculation
1273  if (settings_.useapprox()) {
1275  z1,
1276  phi1,
1277  r2,
1278  z2,
1279  phi2,
1280  r3,
1281  z3,
1282  phi3,
1283  take3,
1284  ndisks,
1285  rinvapprox,
1286  phi0approx,
1287  d0approx,
1288  tapprox,
1289  z0approx,
1290  phiprojapprox,
1291  zprojapprox,
1292  phiderapprox,
1293  zderapprox,
1294  phiprojdiskapprox,
1295  rprojdiskapprox,
1296  phiderdiskapprox,
1297  rderdiskapprox);
1298  } else {
1299  rinvapprox = rinv;
1300  phi0approx = phi0;
1301  d0approx = d0;
1302  tapprox = t;
1303  z0approx = z0;
1304 
1305  for (unsigned int i = 0; i < toR_.size(); ++i) {
1306  phiprojapprox[i] = phiproj[i];
1307  zprojapprox[i] = zproj[i];
1308  phiderapprox[i] = phider[i];
1309  zderapprox[i] = zder[i];
1310  }
1311 
1312  for (unsigned int i = 0; i < toZ_.size(); ++i) {
1313  phiprojdiskapprox[i] = phiprojdisk[i];
1314  rprojdiskapprox[i] = rprojdisk[i];
1315  phiderdiskapprox[i] = phiderdisk[i];
1316  rderdiskapprox[i] = rderdisk[i];
1317  }
1318  }
1319 
1320  //store the approcximate results
1321  if (settings_.debugTracklet()) {
1322  edm::LogVerbatim("Tracklet") << "rinvapprox: " << rinvapprox << " rinv: " << rinv << endl;
1323  edm::LogVerbatim("Tracklet") << "phi0approx: " << phi0approx << " phi0: " << phi0 << endl;
1324  edm::LogVerbatim("Tracklet") << "d0approx: " << d0approx << " d0: " << d0 << endl;
1325  edm::LogVerbatim("Tracklet") << "tapprox: " << tapprox << " t: " << t << endl;
1326  edm::LogVerbatim("Tracklet") << "z0approx: " << z0approx << " z0: " << z0 << endl;
1327  }
1328 
1329  for (unsigned int i = 0; i < toR_.size(); ++i) {
1330  if (settings_.debugTracklet()) {
1331  edm::LogVerbatim("Tracklet") << "phiprojapprox[" << i << "]: " << phiprojapprox[i] << " phiproj[" << i
1332  << "]: " << phiproj[i] << endl;
1333  edm::LogVerbatim("Tracklet") << "zprojapprox[" << i << "]: " << zprojapprox[i] << " zproj[" << i
1334  << "]: " << zproj[i] << endl;
1335  edm::LogVerbatim("Tracklet") << "phiderapprox[" << i << "]: " << phiderapprox[i] << " phider[" << i
1336  << "]: " << phider[i] << endl;
1337  edm::LogVerbatim("Tracklet") << "zderapprox[" << i << "]: " << zderapprox[i] << " zder[" << i << "]: " << zder[i]
1338  << endl;
1339  }
1340  }
1341 
1342  for (unsigned int i = 0; i < toZ_.size(); ++i) {
1343  if (settings_.debugTracklet()) {
1344  edm::LogVerbatim("Tracklet") << "phiprojdiskapprox[" << i << "]: " << phiprojdiskapprox[i] << " phiprojdisk[" << i
1345  << "]: " << phiprojdisk[i] << endl;
1346  edm::LogVerbatim("Tracklet") << "rprojdiskapprox[" << i << "]: " << rprojdiskapprox[i] << " rprojdisk[" << i
1347  << "]: " << rprojdisk[i] << endl;
1348  edm::LogVerbatim("Tracklet") << "phiderdiskapprox[" << i << "]: " << phiderdiskapprox[i] << " phiderdisk[" << i
1349  << "]: " << phiderdisk[i] << endl;
1350  edm::LogVerbatim("Tracklet") << "rderdiskapprox[" << i << "]: " << rderdiskapprox[i] << " rderdisk[" << i
1351  << "]: " << rderdisk[i] << endl;
1352  }
1353  }
1354 
1355  //now binary
1356  double krinv = settings_.kphi1() / settings_.kr() * pow(2, settings_.rinv_shift()),
1357  kphi0 = settings_.kphi1() * pow(2, settings_.phi0_shift()),
1358  kt = settings_.kz() / settings_.kr() * pow(2, settings_.t_shift()),
1359  kz0 = settings_.kz() * pow(2, settings_.z0_shift()),
1360  kphiproj = settings_.kphi1() * pow(2, settings_.SS_phiL_shift()),
1361  kphider = settings_.kphi1() / settings_.kr() * pow(2, settings_.SS_phiderL_shift()),
1362  kzproj = settings_.kz() * pow(2, settings_.PS_zL_shift()),
1363  kzder = settings_.kz() / settings_.kr() * pow(2, settings_.PS_zderL_shift()),
1364  kphiprojdisk = settings_.kphi1() * pow(2, settings_.SS_phiD_shift()),
1365  kphiderdisk = settings_.kphi1() / settings_.kr() * pow(2, settings_.SS_phiderD_shift()),
1366  krprojdisk = settings_.kr() * pow(2, settings_.PS_rD_shift()),
1367  krderdisk = settings_.kr() / settings_.kz() * pow(2, settings_.PS_rderD_shift());
1368 
1369  int irinv, iphi0, id0, it, iz0;
1370  int iphiproj[N_LAYER - 2], izproj[N_LAYER - 2], iphider[N_LAYER - 2], izder[N_LAYER - 2];
1371  int iphiprojdisk[N_DISK], irprojdisk[N_DISK], iphiderdisk[N_DISK], irderdisk[N_DISK];
1372 
1373  //store the binary results
1374  irinv = rinvapprox / krinv;
1375  iphi0 = phi0approx / kphi0;
1376  id0 = d0approx / settings_.kd0();
1377  it = tapprox / kt;
1378  iz0 = z0approx / kz0;
1379 
1380  bool success = true;
1381  if (std::abs(rinvapprox) > settings_.rinvcut()) {
1382  if (settings_.debugTracklet())
1383  edm::LogVerbatim("Tracklet") << "TrackletCalculator:: LLD Seeding irinv too large: " << rinvapprox << "(" << irinv
1384  << ")";
1385  success = false;
1386  }
1387  if (std::abs(z0approx) > settings_.disp_z0cut()) {
1388  if (settings_.debugTracklet())
1389  edm::LogVerbatim("Tracklet") << "Failed tracklet z0 cut " << z0approx;
1390  success = false;
1391  }
1392  if (std::abs(d0approx) > settings_.maxd0()) {
1393  if (settings_.debugTracklet())
1394  edm::LogVerbatim("Tracklet") << "Failed tracklet d0 cut " << d0approx;
1395  success = false;
1396  }
1397 
1398  if (!success)
1399  return false;
1400 
1401  double phicritapprox = phi0approx - asin(0.5 * settings_.rcrit() * rinvapprox);
1402  int phicrit = iphi0 - 2 * irinv;
1403 
1404  int iphicritmincut = settings_.phicritminmc() / globals_->ITC_L1L2()->phi0_final.K();
1405  int iphicritmaxcut = settings_.phicritmaxmc() / globals_->ITC_L1L2()->phi0_final.K();
1406 
1407  bool keepapprox = (phicritapprox > settings_.phicritminmc()) && (phicritapprox < settings_.phicritmaxmc()),
1408  keep = (phicrit > iphicritmincut) && (phicrit < iphicritmaxcut);
1409 
1410  if (settings_.debugTracklet())
1411  if (keep && !keepapprox)
1412  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced::LLDSeeding tracklet kept with exact phicrit cut "
1413  "but not approximate, phicritapprox: "
1414  << phicritapprox;
1415  if (settings_.usephicritapprox()) {
1416  if (!keepapprox)
1417  return false;
1418  } else {
1419  if (!keep)
1420  return false;
1421  }
1422 
1423  LayerProjection layerprojs[N_LAYER - 2];
1424  DiskProjection diskprojs[N_DISK];
1425 
1426  for (unsigned int i = 0; i < toR_.size(); ++i) {
1427  iphiproj[i] = phiprojapprox[i] / kphiproj;
1428  izproj[i] = zprojapprox[i] / kzproj;
1429 
1430  iphider[i] = phiderapprox[i] / kphider;
1431  izder[i] = zderapprox[i] / kzder;
1432 
1433  if (izproj[i] < -(1 << (settings_.nzbitsstub(0) - 1)))
1434  continue;
1435  if (izproj[i] >= (1 << (settings_.nzbitsstub(0) - 1)))
1436  continue;
1437 
1438  //this is left from the original....
1439  if (iphiproj[i] >= (1 << settings_.nphibitsstub(5)) - 1)
1440  continue;
1441  if (iphiproj[i] <= 0)
1442  continue;
1443 
1444  if (rproj_[i] < settings_.rPS2S()) {
1445  iphiproj[i] >>= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
1446  } else {
1447  izproj[i] >>= (settings_.nzbitsstub(0) - settings_.nzbitsstub(5));
1448  }
1449 
1450  if (rproj_[i] < settings_.rPS2S()) {
1451  if (iphider[i] < -(1 << (settings_.nbitsphiprojderL123() - 1)))
1452  iphider[i] = -(1 << (settings_.nbitsphiprojderL123() - 1));
1453  if (iphider[i] >= (1 << (settings_.nbitsphiprojderL123() - 1)))
1454  iphider[i] = (1 << (settings_.nbitsphiprojderL123() - 1)) - 1;
1455  } else {
1456  if (iphider[i] < -(1 << (settings_.nbitsphiprojderL456() - 1)))
1457  iphider[i] = -(1 << (settings_.nbitsphiprojderL456() - 1));
1458  if (iphider[i] >= (1 << (settings_.nbitsphiprojderL456() - 1)))
1459  iphider[i] = (1 << (settings_.nbitsphiprojderL456() - 1)) - 1;
1460  }
1461 
1462  layerprojs[i].init(settings_,
1463  lproj_[i],
1464  rproj_[i],
1465  iphiproj[i],
1466  izproj[i],
1467  iphider[i],
1468  izder[i],
1469  phiproj[i],
1470  zproj[i],
1471  phider[i],
1472  zder[i],
1473  phiprojapprox[i],
1474  zprojapprox[i],
1475  phiderapprox[i],
1476  zderapprox[i]);
1477  }
1478 
1479  if (std::abs(it * kt) > 1.0) {
1480  for (unsigned int i = 0; i < toZ_.size(); ++i) {
1481  iphiprojdisk[i] = phiprojdiskapprox[i] / kphiprojdisk;
1482  irprojdisk[i] = rprojdiskapprox[i] / krprojdisk;
1483 
1484  iphiderdisk[i] = phiderdiskapprox[i] / kphiderdisk;
1485  irderdisk[i] = rderdiskapprox[i] / krderdisk;
1486 
1487  //Check phi range of projection
1488  if (iphiprojdisk[i] <= 0)
1489  continue;
1490  if (iphiprojdisk[i] >= (1 << settings_.nphibitsstub(0)) - 1)
1491  continue;
1492 
1493  //Check r range of projection
1494  if (irprojdisk[i] < settings_.rmindisk() / krprojdisk || irprojdisk[i] > settings_.rmaxdisk() / krprojdisk)
1495  continue;
1496 
1497  diskprojs[i].init(settings_,
1498  i + 1,
1499  rproj_[i],
1500  iphiprojdisk[i],
1501  irprojdisk[i],
1502  iphiderdisk[i],
1503  irderdisk[i],
1504  phiprojdisk[i],
1505  rprojdisk[i],
1506  phiderdisk[i],
1507  rderdisk[i],
1508  phiprojdiskapprox[i],
1509  rprojdiskapprox[i],
1510  phiderdisk[i],
1511  rderdisk[i]);
1512  }
1513  }
1514 
1515  if (settings_.writeMonitorData("TrackletPars")) {
1516  globals_->ofstream("trackletpars.txt")
1517  << layer_ << " , " << rinv << " , " << rinvapprox << " , " << phi0 << " , " << phi0approx << " , " << t << " , "
1518  << tapprox << " , " << z0 << " , " << z0approx << " , " << d0 << " , " << d0approx << endl;
1519  }
1520 
1521  Tracklet* tracklet = new Tracklet(settings_,
1522  innerStub,
1523  middleStub,
1524  outerStub,
1525  innerFPGAStub,
1526  middleFPGAStub,
1527  outerFPGAStub,
1528  rinv,
1529  phi0,
1530  d0,
1531  z0,
1532  t,
1533  rinvapprox,
1534  phi0approx,
1535  d0approx,
1536  z0approx,
1537  tapprox,
1538  irinv,
1539  iphi0,
1540  id0,
1541  iz0,
1542  it,
1543  layerprojs,
1544  diskprojs,
1545  false);
1546 
1547  if (settings_.debugTracklet())
1548  edm::LogVerbatim("Tracklet") << "TrackletCalculatorDisplaced " << getName()
1549  << " Found LLD tracklet in sector = " << iSector_ << " phi0 = " << phi0;
1550 
1552  tracklet->setTCIndex(TCIndex_);
1553 
1554  if (settings_.writeMonitorData("Seeds")) {
1555  ofstream fout("seeds.txt", ofstream::app);
1556  fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_" << iSector_ << " " << tracklet->getISeed() << endl;
1557  fout.close();
1558  }
1559  trackletpars_->addTracklet(tracklet);
1560 
1561  for (unsigned int j = 0; j < toR_.size(); j++) {
1562  if (settings_.debugTracklet())
1563  edm::LogVerbatim("Tracklet") << "adding layer projection " << j << "/" << toR_.size() << " " << lproj_[j];
1564  if (tracklet->validProj(lproj_[j])) {
1565  addLayerProj(tracklet, lproj_[j]);
1566  }
1567  }
1568 
1569  for (unsigned int j = 0; j < toZ_.size(); j++) {
1570  int disk = dproj_[j];
1571  if (disk == 0)
1572  continue;
1573  if (it < 0)
1574  disk = -disk;
1575  if (settings_.debugTracklet())
1576  edm::LogVerbatim("Tracklet") << "adding disk projection " << j << "/" << toZ_.size() << " " << disk;
1577  if (tracklet->validProjDisk(abs(disk))) {
1578  addDiskProj(tracklet, disk);
1579  }
1580  }
1581 
1582  return true;
1583 }
1584 
1586  double rinv,
1587  double phi0,
1588  double d0,
1589  double t,
1590  double z0,
1591  double r0,
1592  double& phiproj,
1593  double& zproj,
1594  double& phider,
1595  double& zder) {
1596  double rho = 1 / rinv;
1597  if (rho < 0) {
1598  r0 = -r0;
1599  }
1600  phiproj = phi0 - asin((rproj * rproj + r0 * r0 - rho * rho) / (2 * rproj * r0));
1601  double beta = acos((rho * rho + r0 * r0 - rproj * rproj) / (2 * r0 * rho));
1602  zproj = z0 + t * std::abs(rho * beta);
1603 
1604  //not exact, but close
1605  phider = -0.5 * rinv / sqrt(1 - pow(0.5 * rproj * rinv, 2)) + d0 / (rproj * rproj);
1606  zder = t / sqrt(1 - pow(0.5 * rproj * rinv, 2));
1607 
1608  if (settings_.debugTracklet())
1609  edm::LogVerbatim("Tracklet") << "exact proj layer at " << rproj << " : " << phiproj << " " << zproj;
1610 }
1611 
1613  double rinv,
1614  double,
1615  double, //phi0 and d0 are not used.
1616  double t,
1617  double z0,
1618  double x0,
1619  double y0,
1620  double& phiproj,
1621  double& rproj,
1622  double& phider,
1623  double& rder) {
1624  //protect against t=0
1625  if (std::abs(t) < 0.1)
1626  t = 0.1;
1627  if (t < 0)
1628  zproj = -zproj;
1629  double rho = std::abs(1 / rinv);
1630  double beta = (zproj - z0) / (t * rho);
1631  double phiV = atan2(-y0, -x0);
1632  double c = rinv > 0 ? -1 : 1;
1633 
1634  double x = x0 + rho * cos(phiV + c * beta);
1635  double y = y0 + rho * sin(phiV + c * beta);
1636 
1637  phiproj = atan2(y, x);
1638 
1639  phiproj = reco::reduceRange(phiproj - phimin_);
1640 
1641  rproj = sqrt(x * x + y * y);
1642 
1643  phider = c / t / (x * x + y * y) * (rho + x0 * cos(phiV + c * beta) + y0 * sin(phiV + c * beta));
1644  rder = c / t / rproj * (y0 * cos(phiV + c * beta) - x0 * sin(phiV + c * beta));
1645 
1646  if (settings_.debugTracklet())
1647  edm::LogVerbatim("Tracklet") << "exact proj disk at" << zproj << " : " << phiproj << " " << rproj;
1648 }
1649 
1651  double z1,
1652  double phi1,
1653  double r2,
1654  double z2,
1655  double phi2,
1656  double r3,
1657  double z3,
1658  double phi3,
1659  int take3,
1660  double& rinv,
1661  double& phi0,
1662  double& d0,
1663  double& t,
1664  double& z0,
1665  double phiproj[N_DISK],
1666  double zproj[N_DISK],
1667  double phiprojdisk[N_DISK],
1668  double rprojdisk[N_DISK],
1669  double phider[N_DISK],
1670  double zder[N_DISK],
1671  double phiderdisk[N_DISK],
1672  double rderdisk[N_DISK]) {
1673  //two lines perpendicular to the 1->2 and 2->3
1674  double x1 = r1 * cos(phi1);
1675  double x2 = r2 * cos(phi2);
1676  double x3 = r3 * cos(phi3);
1677 
1678  double y1 = r1 * sin(phi1);
1679  double y2 = r2 * sin(phi2);
1680  double y3 = r3 * sin(phi3);
1681 
1682  double k1 = -(x2 - x1) / (y2 - y1);
1683  double k2 = -(x3 - x2) / (y3 - y2);
1684  double b1 = 0.5 * (y2 + y1) - 0.5 * (x1 + x2) * k1;
1685  double b2 = 0.5 * (y3 + y2) - 0.5 * (x2 + x3) * k2;
1686  //their intersection gives the center of the circle
1687  double y0 = (b1 * k2 - b2 * k1) / (k2 - k1);
1688  double x0 = (b1 - b2) / (k2 - k1);
1689  //get the radius three ways:
1690  double R1 = sqrt(pow(x1 - x0, 2) + pow(y1 - y0, 2));
1691  double R2 = sqrt(pow(x2 - x0, 2) + pow(y2 - y0, 2));
1692  double R3 = sqrt(pow(x3 - x0, 2) + pow(y3 - y0, 2));
1693  //check if the same
1694  double eps1 = std::abs(R1 / R2 - 1);
1695  double eps2 = std::abs(R3 / R2 - 1);
1696  if (eps1 > 1e-10 || eps2 > 1e-10)
1697  edm::LogVerbatim("Tracklet") << "&&&&&&&&&&&& bad circle! " << R1 << "\t" << R2 << "\t" << R3;
1698 
1699  if (settings_.debugTracklet())
1700  edm::LogVerbatim("Tracklet") << "phimin_: " << phimin_ << " phimax_: " << phimax_;
1701  //results
1702  rinv = 1. / R1;
1703  phi0 = 0.5 * M_PI + atan2(y0, x0);
1704 
1705  phi0 -= phimin_;
1706 
1707  d0 = -R1 + sqrt(x0 * x0 + y0 * y0);
1708  //sign of rinv:
1709  double dphi = reco::reduceRange(phi3 - atan2(y0, x0));
1710  if (dphi < 0) {
1711  rinv = -rinv;
1712  d0 = -d0;
1713  phi0 = phi0 + M_PI;
1714  }
1715  phi0 = angle0to2pi::make0To2pi(phi0);
1716 
1717  //now in RZ:
1718  //turning angle
1719  double beta1 = reco::reduceRange(atan2(y1 - y0, x1 - x0) - atan2(-y0, -x0));
1720  double beta2 = reco::reduceRange(atan2(y2 - y0, x2 - x0) - atan2(-y0, -x0));
1721  double beta3 = reco::reduceRange(atan2(y3 - y0, x3 - x0) - atan2(-y0, -x0));
1722 
1723  double t12 = (z2 - z1) / std::abs(beta2 - beta1) / R1;
1724  double z12 = (z1 * beta2 - z2 * beta1) / (beta2 - beta1);
1725  double t13 = (z3 - z1) / std::abs(beta3 - beta1) / R1;
1726  double z13 = (z1 * beta3 - z3 * beta1) / (beta3 - beta1);
1727 
1728  if (take3 > 0) {
1729  //take 13 (large lever arm)
1730  t = t13;
1731  z0 = z13;
1732  } else {
1733  //take 12 (pixel layers)
1734  t = t12;
1735  z0 = z12;
1736  }
1737 
1738  for (unsigned int i = 0; i < toR_.size(); i++) {
1739  exactproj(toR_[i], rinv, phi0, d0, t, z0, sqrt(x0 * x0 + y0 * y0), phiproj[i], zproj[i], phider[i], zder[i]);
1740  }
1741 
1742  for (unsigned int i = 0; i < toZ_.size(); i++) {
1743  exactprojdisk(toZ_[i], rinv, phi0, d0, t, z0, x0, y0, phiprojdisk[i], rprojdisk[i], phiderdisk[i], rderdisk[i]);
1744  }
1745 
1746  if (settings_.debugTracklet())
1747  edm::LogVerbatim("Tracklet") << "exact tracklet: " << rinv << " " << phi0 << " " << t << " " << z0 << " " << d0;
1748 }
1749 
1751  double phi0,
1752  double d0,
1753  double t,
1754  double z0,
1755  double halfRinv_0,
1756  double d0_0, // zeroeth order result for higher order terms calculation
1757  double rmean,
1758  double& phiproj,
1759  double& phiprojder,
1760  double& zproj,
1761  double& zprojder) {
1762  if (std::abs(2.0 * halfRinv) > settings_.rinvcut() || std::abs(z0) > settings_.disp_z0cut() ||
1763  std::abs(d0) > settings_.maxd0()) {
1764  phiproj = 0.0;
1765  return;
1766  }
1767  double rmeanInv = 1.0 / rmean;
1768 
1769  phiproj = phi0 + rmean * (-halfRinv + 2.0 * d0_0 * halfRinv_0 * halfRinv_0) +
1770  rmeanInv * (-d0 + halfRinv_0 * d0_0 * d0_0) + sixth * pow(-rmean * halfRinv_0 - rmeanInv * d0_0, 3);
1771  phiprojder = -halfRinv + d0 * rmeanInv * rmeanInv; //removed all high terms
1772 
1773  zproj = z0 + t * rmean - 0.5 * rmeanInv * t * d0_0 * d0_0 - t * rmean * halfRinv * d0 +
1774  sixth * pow(rmean, 3) * t * halfRinv_0 * halfRinv_0;
1775  zprojder = t; // removed all high terms
1776 
1777  phiproj = angle0to2pi::make0To2pi(phiproj);
1778 
1779  if (settings_.debugTracklet())
1780  edm::LogVerbatim("Tracklet") << "approx proj layer at " << rmean << " : " << phiproj << " " << zproj << endl;
1781 }
1782 
1784  double phi0,
1785  double d0,
1786  double t,
1787  double z0,
1788  double halfRinv_0,
1789  double d0_0, // zeroeth order result for higher order terms calculation
1790  double zmean,
1791  double& phiproj,
1792  double& phiprojder,
1793  double& rproj,
1794  double& rprojder) {
1795  if (std::abs(2.0 * halfRinv) > settings_.rinvcut() || std::abs(z0) > settings_.disp_z0cut() ||
1796  std::abs(d0) > settings_.maxd0()) {
1797  phiproj = 0.0;
1798  return;
1799  }
1800 
1801  if (t < 0)
1802  zmean = -zmean;
1803 
1804  double zmeanInv = 1.0 / zmean, rstar = (zmean - z0) / t,
1805  epsilon = 0.5 * zmeanInv * zmeanInv * d0_0 * d0_0 * t * t + halfRinv * d0 -
1806  sixth * rstar * rstar * halfRinv_0 * halfRinv_0;
1807 
1808  rproj = rstar * (1 + epsilon);
1809  rprojder = 1 / t;
1810 
1811  double A = rproj * halfRinv;
1812  double B = -d0 * t * zmeanInv * (1 + z0 * zmeanInv) * (1 - epsilon);
1813  double C = -d0 * halfRinv;
1814  double A_0 = rproj * halfRinv_0;
1815  double B_0 = -d0_0 * t * zmeanInv * (1 + z0 * zmeanInv) * (1 - epsilon);
1816  // double C_0 = -d0_0 * halfRinv_0;
1817 
1818  phiproj = phi0 - A + B * (1 + C - 2 * A_0 * A_0) + sixth * pow(-A_0 + B_0, 3);
1819  phiprojder = -halfRinv / t - d0 * t * t * zmeanInv * zmeanInv;
1820 
1821  phiproj = angle0to2pi::make0To2pi(phiproj);
1822 
1823  if (settings_.debugTracklet())
1824  edm::LogVerbatim("Tracklet") << "approx proj disk at" << zmean << " : " << phiproj << " " << rproj << endl;
1825 }
1826 
1828  double z1,
1829  double phi1,
1830  double r2,
1831  double z2,
1832  double phi2,
1833  double r3,
1834  double z3,
1835  double phi3,
1836  bool take3,
1837  unsigned ndisks,
1838  double& rinv,
1839  double& phi0,
1840  double& d0,
1841  double& t,
1842  double& z0,
1843  double phiproj[4],
1844  double zproj[4],
1845  double phider[4],
1846  double zder[4],
1847  double phiprojdisk[5],
1848  double rprojdisk[5],
1849  double phiderdisk[5],
1850  double rderdisk[5]) {
1851  double a = 1.0 / ((r1 - r2) * (r1 - r3));
1852  double b = 1.0 / ((r1 - r2) * (r2 - r3));
1853  double c = 1.0 / ((r1 - r3) * (r2 - r3));
1854 
1855  // first iteration in r-phi plane
1856  double halfRinv_0 = -phi1 * r1 * a + phi2 * r2 * b - phi3 * r3 * c;
1857  double d0_0 = r1 * r2 * r3 * (-phi1 * a + phi2 * b - phi3 * c);
1858 
1859  // corrections to phi1, phi2, and phi3
1860  double r = r2, z = z2;
1861  if (take3)
1862  r = r3, z = z3;
1863 
1864  double d0OverR1 = d0_0 * rzmeanInv_[0] * (ndisks > 2 ? std::abs((z - z1) / (r - r1)) : 1.0);
1865  double d0OverR2 = d0_0 * rzmeanInv_[1] * (ndisks > 1 ? std::abs((z - z1) / (r - r1)) : 1.0);
1866  double d0OverR3 = d0_0 * rzmeanInv_[2] * (ndisks > 0 ? std::abs((z - z1) / (r - r1)) : 1.0);
1867 
1868  double d0OverR = d0OverR2;
1869  if (take3)
1870  d0OverR = d0OverR3;
1871 
1872  double c1 = d0_0 * halfRinv_0 * d0OverR1 + 2.0 * d0_0 * halfRinv_0 * r1 * halfRinv_0 +
1873  sixth * pow(-r1 * halfRinv_0 - d0OverR1, 3);
1874  double c2 = d0_0 * halfRinv_0 * d0OverR2 + 2.0 * d0_0 * halfRinv_0 * r2 * halfRinv_0 +
1875  sixth * pow(-r2 * halfRinv_0 - d0OverR2, 3);
1876  double c3 = d0_0 * halfRinv_0 * d0OverR3 + 2.0 * d0_0 * halfRinv_0 * r3 * halfRinv_0 +
1877  sixth * pow(-r3 * halfRinv_0 - d0OverR3, 3);
1878 
1879  double phi1c = phi1 - c1;
1880  double phi2c = phi2 - c2;
1881  double phi3c = phi3 - c3;
1882 
1883  // second iteration in r-phi plane
1884  double halfRinv = -phi1c * r1 * a + phi2c * r2 * b - phi3c * r3 * c;
1885  phi0 = -phi1c * r1 * (r2 + r3) * a + phi2c * r2 * (r1 + r3) * b - phi3c * r3 * (r1 + r2) * c;
1886  d0 = r1 * r2 * r3 * (-phi1c * a + phi2c * b - phi3c * c);
1887 
1888  t = ((z - z1) / (r - r1)) *
1889  (1. + d0 * halfRinv - 0.5 * d0OverR1 * d0OverR - sixth * (r1 * r1 + r2 * r2 + r1 * r2) * halfRinv_0 * halfRinv_0);
1890  z0 = z1 - t * r1 * (1.0 - d0_0 * halfRinv_0 - 0.5 * d0OverR1 * d0OverR1 + sixth * r1 * r1 * halfRinv_0 * halfRinv_0);
1891 
1892  rinv = 2.0 * halfRinv;
1893  phi0 += -phimin_;
1894 
1895  phi0 = angle0to2pi::make0To2pi(phi0);
1896 
1897  for (unsigned int i = 0; i < toR_.size(); i++) {
1898  approxproj(halfRinv,
1899  phi0,
1900  d0,
1901  t,
1902  z0,
1903  halfRinv_0,
1904  d0_0, // added _0 version for high term calculations
1905  toR_.at(i),
1906  phiproj[i],
1907  phider[i],
1908  zproj[i],
1909  zder[i]);
1910  }
1911 
1912  for (unsigned int i = 0; i < toZ_.size(); i++) {
1913  approxprojdisk(halfRinv,
1914  phi0,
1915  d0,
1916  t,
1917  z0,
1918  halfRinv_0,
1919  d0_0, // added _0 version for high term calculations
1920  toZ_.at(i),
1921  phiprojdisk[i],
1922  phiderdisk[i],
1923  rprojdisk[i],
1924  rderdisk[i]);
1925  }
1926 
1928  std::abs(d0) > settings_.maxd0()) {
1929  phi0 = 0.0;
1930  return;
1931  }
1932 
1933  if (settings_.debugTracklet())
1934  edm::LogVerbatim("Tracklet") << "TCD approx tracklet: " << rinv << " " << phi0 << " " << t << " " << z0 << " " << d0
1935  << endl;
1936 }
trklet::TrackletCalculatorDisplaced::toZ_
std::vector< double > toZ_
Definition: TrackletCalculatorDisplaced.h:169
Settings.h
trklet::TrackletCalculatorDisplaced::outerallstubs_
std::vector< AllStubsMemory * > outerallstubs_
Definition: TrackletCalculatorDisplaced.h:173
trklet::Settings::phicritmaxmc
double phicritmaxmc() const
Definition: Settings.h:265
L1TStub.h
trklet::Settings::rmindisk
double rmindisk() const
Definition: Settings.h:113
trklet::Settings::maxd0
double maxd0() const
Definition: Settings.h:274
trklet::Settings::rmindiskvm
double rmindiskvm() const
Definition: Settings.h:281
mps_fire.i
i
Definition: mps_fire.py:428
trklet::Settings::ntrackletmax
unsigned int ntrackletmax() const
Definition: Settings.h:297
trklet::Settings::writetrace
bool writetrace() const
Definition: Settings.h:162
input
static const std::string input
Definition: EdmProvDump.cc:48
MessageLogger.h
trklet::TrackletCalculatorDisplaced::addProjection
void addProjection(int layer, int iphi, TrackletProjectionsMemory *trackletprojs, Tracklet *tracklet)
Definition: TrackletCalculatorDisplaced.cc:330
trklet::Settings::krprojshiftdisk
double krprojshiftdisk() const
Definition: Settings.h:364
trklet::TrackletCalculatorDisplaced::innerallstubs_
std::vector< AllStubsMemory * > innerallstubs_
Definition: TrackletCalculatorDisplaced.h:171
trklet::Tracklet::fpgaphiprojdisk
const FPGAWord & fpgaphiprojdisk(int disk) const
Definition: Tracklet.h:266
trklet::TrackletParametersMemory::addTracklet
void addTracklet(Tracklet *tracklet)
Definition: TrackletParametersMemory.h:23
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:72
trklet::TrackletCalculatorDisplaced::stubtriplets_
std::vector< StubTripletsMemory * > stubtriplets_
Definition: TrackletCalculatorDisplaced.h:174
deltaPhi.h
HLT_FULL_cff.beta
beta
Definition: HLT_FULL_cff.py:8686
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, bool isPSseed=false)
Definition: LayerProjection.cc:9
trklet::Settings::SS_phiL_shift
int SS_phiL_shift() const
Definition: Settings.h:314
testProducerWithPsetDescEmpty_cfi.x2
x2
Definition: testProducerWithPsetDescEmpty_cfi.py:28
memory
Definition: HGCRecHitSoA.h:20
trklet::Settings::phicritminmc
double phicritminmc() const
Definition: Settings.h:264
trklet::TrackletCalculatorDisplaced::disk_
int disk_
Definition: TrackletCalculatorDisplaced.h:161
trklet::TrackletProjectionsMemory
Definition: TrackletProjectionsMemory.h:15
trklet::Settings
Definition: Settings.h:31
trklet::L1TStub
Definition: L1TStub.h:12
trklet::TrackletCalculatorDisplaced::addLayerProj
bool addLayerProj(Tracklet *tracklet, int layer)
Definition: TrackletCalculatorDisplaced.cc:310
trklet::L1TStub::z
double z() const
Definition: L1TStub.h:56
trklet::Tracklet::getISeed
int getISeed() const
Definition: Tracklet.cc:849
trklet::TrackletCalculatorDisplaced::middleallstubs_
std::vector< AllStubsMemory * > middleallstubs_
Definition: TrackletCalculatorDisplaced.h:172
cms::cuda::assert
assert(be >=bs)
trklet::Settings::rmean
double rmean(unsigned int iLayer) const
Definition: Settings.h:143
trklet::N_DISK
constexpr int N_DISK
Definition: Settings.h:20
trklet::ProcessBase::settings_
Settings const & settings_
Definition: ProcessBase.h:44
trklet::Stub::phiapprox
double phiapprox(double phimin, double) const
Definition: Stub.cc:242
trklet::Settings::kphi1
double kphi1() const
Definition: Settings.h:268
trklet::TrackletCalculatorDisplaced::lproj_
int lproj_[N_LAYER - 2]
Definition: TrackletCalculatorDisplaced.h:163
trklet::FPGAWord::nbits
int nbits() const
Definition: FPGAWord.h:25
b2
static constexpr float b2
Definition: L1EGammaCrystalsEmulatorProducer.cc:83
createJobs.tmp
tmp
align.sh
Definition: createJobs.py:716
trklet::Globals
Definition: Globals.h:32
trklet::Tracklet::validProj
bool validProj(int layer) const
Definition: Tracklet.h:89
trklet::Settings::PS_rderD_shift
int PS_rderD_shift() const
Definition: Settings.h:326
trklet::Stub::isBarrel
bool isBarrel() const
Definition: Stub.h:60
trklet::LayerProjection
Definition: LayerProjection.h:10
trklet::Tracklet
Definition: Tracklet.h:28
geometryDiff.epsilon
int epsilon
Definition: geometryDiff.py:26
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::sixth
constexpr double sixth
Definition: Settings.h:29
trklet::TrackletCalculatorDisplaced::approxprojdisk
void approxprojdisk(double halfRinv, double phi0, double d0, double t, double z0, double halfRinv_0, double d0_0, double zmean, double &phiproj, double &phiprojder, double &rproj, double &rprojder)
Definition: TrackletCalculatorDisplaced.cc:1783
b1
static constexpr float b1
Definition: L1EGammaCrystalsEmulatorProducer.cc:83
accept
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
testProducerWithPsetDescEmpty_cfi.x1
x1
Definition: testProducerWithPsetDescEmpty_cfi.py:33
TrackletCalculatorDisplaced.h
testProducerWithPsetDescEmpty_cfi.y1
y1
Definition: testProducerWithPsetDescEmpty_cfi.py:29
trklet::Tracklet::fpgazproj
const FPGAWord & fpgazproj(int layer) const
Definition: Tracklet.h:99
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::MemoryBase::getName
std::string const & getName() const
Definition: MemoryBase.h:19
trklet::Settings::z0_shift
int z0_shift() const
Definition: Settings.h:309
trklet::Stub
Definition: Stub.h:16
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
trklet::DiskProjection
Definition: DiskProjection.h:12
trklet::TrackletCalculatorDisplaced::layer_
int layer_
Definition: TrackletCalculatorDisplaced.h:160
trklet::TrackletCalculatorDisplaced::trackletprojlayers_
std::vector< std::vector< TrackletProjectionsMemory * > > trackletprojlayers_
Definition: TrackletCalculatorDisplaced.h:179
trklet::N_LAYER
constexpr int N_LAYER
Definition: Settings.h:19
summarizeEdmComparisonLogfiles.success
success
Definition: summarizeEdmComparisonLogfiles.py:115
trklet::TrackletCalculatorDisplaced::rzmeanInv_
double rzmeanInv_[N_DISK - 2]
Definition: TrackletCalculatorDisplaced.h:166
trklet::Settings::zmean
double zmean(unsigned int iDisk) const
Definition: Settings.h:146
trklet::TrackletCalculatorDisplaced::exacttracklet
void exacttracklet(double r1, double z1, double phi1, double r2, double z2, double phi2, double r3, double z3, double phi3, int take3, double &rinv, double &phi0, double &d0, double &t, double &z0, double phiproj[N_LAYER - 2], double zproj[N_LAYER - 2], double phiprojdisk[N_DISK], double rprojdisk[N_DISK], double phider[N_LAYER - 2], double zder[N_LAYER - 2], double phiderdisk[N_DISK], double rderdisk[N_DISK])
Definition: TrackletCalculatorDisplaced.cc:1650
trklet::Stub::isDisk
bool isDisk() const
Definition: Stub.h:61
HLTMuonOfflineAnalyzer_cfi.z0
z0
Definition: HLTMuonOfflineAnalyzer_cfi.py:98
b
double b
Definition: hdecay.h:118
trklet::Settings::phi0_shift
int phi0_shift() const
Definition: Settings.h:307
trklet::Settings::SS_phiderD_shift
int SS_phiderD_shift() const
Definition: Settings.h:325
trklet::FPGAWord
Definition: FPGAWord.h:9
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
trklet::Settings::kr
double kr() const
Definition: Settings.h:271
trklet::TrackletCalculatorDisplaced::exactprojdisk
void exactprojdisk(double zproj, double rinv, double, double, double t, double z0, double x0, double y0, double &phiproj, double &rproj, double &phider, double &rder)
Definition: TrackletCalculatorDisplaced.cc:1612
trklet::TrackletCalculatorDisplaced::LLDSeeding
bool LLDSeeding(const Stub *innerFPGAStub, const L1TStub *innerStub, const Stub *middleFPGAStub, const L1TStub *middleStub, const Stub *outerFPGAStub, const L1TStub *outerStub)
Definition: TrackletCalculatorDisplaced.cc:1186
trklet::Settings::warnNoMem
bool warnNoMem() const
Definition: Settings.h:164
testProducerWithPsetDescEmpty_cfi.y2
y2
Definition: testProducerWithPsetDescEmpty_cfi.py:30
trklet::Settings::rmaxdisk
double rmaxdisk() const
Definition: Settings.h:112
trklet::MemoryBase
Definition: MemoryBase.h:13
trklet::Tracklet::setTrackletIndex
void setTrackletIndex(unsigned int index)
Definition: Tracklet.cc:844
a
double a
Definition: hdecay.h:119
trklet::Settings::nbitsphiprojderL456
unsigned int nbitsphiprojderL456() const
Definition: Settings.h:75
trklet::TrackletCalculatorDisplaced::addDiskProj
void addDiskProj(Tracklet *tracklet, int disk)
Definition: TrackletCalculatorDisplaced.cc:294
trklet::FPGAWord::atExtreme
bool atExtreme() const
Definition: FPGAWord.cc:79
trklet::rinv
double rinv(double phi1, double phi2, double r1, double r2)
Definition: Util.h:173
trklet::TrackletCalculatorDisplaced::LLLSeeding
bool LLLSeeding(const Stub *innerFPGAStub, const L1TStub *innerStub, const Stub *middleFPGAStub, const L1TStub *middleStub, const Stub *outerFPGAStub, const L1TStub *outerStub)
Definition: TrackletCalculatorDisplaced.cc:364
trklet::VarBase::K
double K() const
Definition: imath.h:246
trklet::Globals::ITC_L1L2
IMATH_TrackletCalculator * ITC_L1L2()
Definition: Globals.h:52
trklet::Settings::nzbitsstub
unsigned int nzbitsstub(unsigned int layerdisk) const
Definition: Settings.h:69
trklet::IMATH_TrackletCalculator::phi0_final
VarAdjustK phi0_final
Definition: IMATH_TrackletCalculator.h:218
trklet::Settings::SS_phiderL_shift
int SS_phiderL_shift() const
Definition: Settings.h:317
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::TrackletCalculatorDisplaced::execute
void execute()
Definition: TrackletCalculatorDisplaced.cc:223
alignmentValidation.c1
c1
do drawing
Definition: alignmentValidation.py:1025
Globals.h
trklet::Settings::nphibitsstub
unsigned int nphibitsstub(unsigned int layerdisk) const
Definition: Settings.h:70
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
trklet::Settings::t_shift
int t_shift() const
Definition: Settings.h:308
trklet::TrackletCalculatorDisplaced::zproj_
double zproj_[N_DISK - 2]
Definition: TrackletCalculatorDisplaced.h:164
trklet::TrackletCalculatorDisplaced::TCIndex_
int TCIndex_
Definition: TrackletCalculatorDisplaced.h:159
trklet::TrackletCalculatorDisplaced::rproj_
double rproj_[N_LAYER - 2]
Definition: TrackletCalculatorDisplaced.h:162
trklet
Definition: AllProjectionsMemory.h:9
trklet::FPGAWord::value
int value() const
Definition: FPGAWord.h:24
IMATH_TrackletCalculator.h
alignCSCRings.r
r
Definition: alignCSCRings.py:93
angle0to2pi::make0To2pi
constexpr valType make0To2pi(valType angle)
Definition: deltaPhi.h:67
trklet::TrackletCalculatorDisplaced::trackletprojdisks_
std::vector< std::vector< TrackletProjectionsMemory * > > trackletprojdisks_
Definition: TrackletCalculatorDisplaced.h:180
trklet::TrackletCalculatorDisplaced::approxtracklet
void approxtracklet(double r1, double z1, double phi1, double r2, double z2, double phi2, double r3, double z3, double phi3, bool take3, unsigned ndisks, double &rinv, double &phi0, double &d0, double &t, double &z0, double phiproj[4], double zproj[4], double phider[4], double zder[4], double phiprojdisk[5], double rprojdisk[5], double phiderdisk[5], double rderdisk[5])
Definition: TrackletCalculatorDisplaced.cc:1827
trklet::TrackletCalculatorDisplaced::approxproj
void approxproj(double halfRinv, double phi0, double d0, double t, double z0, double halfRinv_0, double d0_0, double rmean, double &phiproj, double &phiprojder, double &zproj, double &zprojder)
Definition: TrackletCalculatorDisplaced.cc:1750
trklet::TrackletCalculatorDisplaced::dproj_
int dproj_[N_DISK - 2]
Definition: TrackletCalculatorDisplaced.h:165
trklet::Settings::writeMonitorData
bool writeMonitorData(std::string module) const
Definition: Settings.h:96
TtFullHadDaughter::B
static const std::string B
Definition: TtFullHadronicEvent.h:9
std
Definition: JetResolutionObject.h:76
trklet::TrackletCalculatorDisplaced::toR_
std::vector< double > toR_
Definition: TrackletCalculatorDisplaced.h:168
trklet::Settings::usephicritapprox
bool usephicritapprox() const
Definition: Settings.h:212
trklet::TrackletParametersMemory::nTracklets
unsigned int nTracklets() const
Definition: TrackletParametersMemory.h:25
trklet::Tracklet::setTCIndex
void setTCIndex(int index)
Definition: Tracklet.h:496
trklet::ProcessBase
Definition: ProcessBase.h:12
trklet::Settings::rinvcut
double rinvcut() const
Definition: Settings.h:187
gen::C
C
Definition: PomwigHadronizer.cc:78
trklet::Globals::ofstream
std::ofstream & ofstream(std::string fname)
Definition: Globals.cc:44
trklet::Settings::SS_phiD_shift
int SS_phiD_shift() const
Definition: Settings.h:322
trklet::TrackletProjectionsMemory::addProj
void addProj(Tracklet *tracklet)
Definition: TrackletProjectionsMemory.cc:17
trklet::TrackletCalculatorDisplaced::exactproj
void exactproj(double rproj, double rinv, double phi0, double d0, double t, double z0, double r0, double &phiproj, double &zproj, double &phider, double &zder)
Definition: TrackletCalculatorDisplaced.cc:1585
trklet::Stub::l1tstub
L1TStub * l1tstub()
Definition: Stub.h:69
diffTwoXMLs.r1
r1
Definition: diffTwoXMLs.py:53
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
trklet::Settings::useapprox
bool useapprox() const
Definition: Settings.h:211
trklet::TrackletCalculatorDisplaced::DDLSeeding
bool DDLSeeding(const Stub *innerFPGAStub, const L1TStub *innerStub, const Stub *middleFPGAStub, const L1TStub *middleStub, const Stub *outerFPGAStub, const L1TStub *outerStub)
Definition: TrackletCalculatorDisplaced.cc:786
Exception
Definition: hltDiff.cc:245
trklet::TrackletCalculatorDisplaced::addOutput
void addOutput(MemoryBase *memory, std::string output) override
Definition: TrackletCalculatorDisplaced.cc:149
trklet::TrackletCalculatorDisplaced::trackletpars_
TrackletParametersMemory * trackletpars_
Definition: TrackletCalculatorDisplaced.h:176
trklet::ProcessBase::getName
std::string const & getName() const
Definition: ProcessBase.h:22
trklet::Tracklet::fpgaphiproj
const FPGAWord & fpgaphiproj(int layer) const
Definition: Tracklet.h:104
trklet::Settings::kz
double kz() const
Definition: Settings.h:270
trklet::ProcessBase::name_
std::string name_
Definition: ProcessBase.h:38
trklet::Settings::nallstubs
unsigned int nallstubs(unsigned int layerdisk) const
Definition: Settings.h:94
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:161
trklet::Settings::rcrit
double rcrit() const
Definition: Settings.h:257
Exception.h
trklet::TrackletCalculatorDisplaced::addInput
void addInput(MemoryBase *memory, std::string input) override
Definition: TrackletCalculatorDisplaced.cc:191
trklet::Tracklet::fpgarprojdisk
const FPGAWord & fpgarprojdisk(int disk) const
Definition: Tracklet.h:276
trklet::Settings::maxStep
unsigned int maxStep(std::string module) const
Definition: Settings.h:103
trklet::Settings::kd0
double kd0() const
Definition: Settings.h:277
trklet::L1TStub::r
double r() const
Definition: L1TStub.h:57
trklet::Settings::PS_zderL_shift
int PS_zderL_shift() const
Definition: Settings.h:318
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
trklet::Settings::PS_rD_shift
int PS_rD_shift() const
Definition: Settings.h:323
trklet::ProcessBase::phimin_
double phimin_
Definition: ProcessBase.h:41
trklet::Settings::PS_zL_shift
int PS_zL_shift() const
Definition: Settings.h:315
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
c
auto & c
Definition: CAHitNtupletGeneratorKernelsImpl.h:46
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
trklet::Settings::nbitsphiprojderL123
unsigned int nbitsphiprojderL123() const
Definition: Settings.h:74
trklet::Stub::layer
const FPGAWord & layer() const
Definition: Stub.h:56
trklet::Tracklet::validProjDisk
bool validProjDisk(int disk) const
Definition: Tracklet.h:206
trklet::ProcessBase::globals_
Globals * globals_
Definition: ProcessBase.h:45
trklet::Settings::rPS2S
double rPS2S() const
Definition: Settings.h:288
submitPVValidationJobs.t
string t
Definition: submitPVValidationJobs.py:644
d0
static constexpr float d0
Definition: L1EGammaCrystalsEmulatorProducer.cc:85
edm::Log
Definition: MessageLogger.h:70
trklet::ProcessBase::phimax_
double phimax_
Definition: ProcessBase.h:42
A
trklet::TrackletCalculatorDisplaced::addOutputProjection
void addOutputProjection(TrackletProjectionsMemory *&outputProj, MemoryBase *memory)
Definition: TrackletCalculatorDisplaced.cc:144
keep
const int keep
Definition: GenParticlePruner.cc:48
Stub.h
Tracklet.h
trklet::Stub::zapprox
double zapprox() const
Definition: Stub.cc:223
trklet::Settings::disp_z0cut
double disp_z0cut() const
Definition: Settings.h:292
trklet::L1TStub::phi
double phi() const
Definition: L1TStub.h:62
trklet::Settings::rinv_shift
int rinv_shift() const
Definition: Settings.h:306
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
trklet::TrackletCalculatorDisplaced::addProjectionDisk
void addProjectionDisk(int disk, int iphi, TrackletProjectionsMemory *trackletprojs, Tracklet *tracklet)
Definition: TrackletCalculatorDisplaced.cc:345
reco::reduceRange
constexpr T reduceRange(T x)
Definition: deltaPhi.h:18
trklet::Settings::zlength
double zlength() const
Definition: Settings.h:111
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