CMS 3D CMS Logo

FitTrack.cc
Go to the documentation of this file.
9 
12 
13 using namespace std;
14 using namespace trklet;
15 
16 FitTrack::FitTrack(string name, Settings const& settings, Globals* global)
17  : ProcessBase(name, settings, global), trackfit_(nullptr) {}
18 
20  if (settings_.writetrace()) {
21  edm::LogVerbatim("Tracklet") << "In " << name_ << " adding output to " << memory->getName() << " to output "
22  << output;
23  }
24  if (output == "trackout") {
25  TrackFitMemory* tmp = dynamic_cast<TrackFitMemory*>(memory);
26  assert(tmp != nullptr);
27  trackfit_ = tmp;
28  return;
29  }
30 
31  throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " addOutput, output = " << output << " not known";
32 }
33 
35  if (settings_.writetrace()) {
36  edm::LogVerbatim("Tracklet") << "In " << name_ << " adding input from " << memory->getName() << " to input "
37  << input;
38  }
39  if (input.substr(0, 4) == "tpar") {
40  auto* tmp = dynamic_cast<TrackletParametersMemory*>(memory);
41  assert(tmp != nullptr);
42  seedtracklet_.push_back(tmp);
43  return;
44  }
45  if (input.substr(0, 10) == "fullmatch1") {
46  auto* tmp = dynamic_cast<FullMatchMemory*>(memory);
47  assert(tmp != nullptr);
48  fullmatch1_.push_back(tmp);
49  return;
50  }
51  if (input.substr(0, 10) == "fullmatch2") {
52  auto* tmp = dynamic_cast<FullMatchMemory*>(memory);
53  assert(tmp != nullptr);
54  fullmatch2_.push_back(tmp);
55  return;
56  }
57  if (input.substr(0, 10) == "fullmatch3") {
58  auto* tmp = dynamic_cast<FullMatchMemory*>(memory);
59  assert(tmp != nullptr);
60  fullmatch3_.push_back(tmp);
61  return;
62  }
63  if (input.substr(0, 10) == "fullmatch4") {
64  auto* tmp = dynamic_cast<FullMatchMemory*>(memory);
65  assert(tmp != nullptr);
66  fullmatch4_.push_back(tmp);
67  return;
68  }
69 
70  throw cms::Exception("BadConfig") << __FILE__ << " " << __LINE__ << " input = " << input << " not found";
71 }
72 
73 #ifdef USEHYBRID
74 void FitTrack::trackFitKF(Tracklet* tracklet,
75  std::vector<const Stub*>& trackstublist,
76  std::vector<std::pair<int, int>>& stubidslist) {
77  if (settings_.doKF()) {
78  // From full match lists, collect all the stubs associated with the tracklet seed
79 
80  // Get seed stubs first
81  trackstublist.emplace_back(tracklet->innerFPGAStub());
82  if (tracklet->getISeed() >= (int)N_TRKLSEED + 1)
83  trackstublist.emplace_back(tracklet->middleFPGAStub());
84  trackstublist.emplace_back(tracklet->outerFPGAStub());
85 
86  // Now get ALL matches (can have multiple per layer)
87  for (const auto& i : fullmatch1_) {
88  for (unsigned int j = 0; j < i->nMatches(); j++) {
89  if (i->getTracklet(j)->TCID() == tracklet->TCID()) {
90  trackstublist.push_back(i->getMatch(j).second);
91  }
92  }
93  }
94 
95  for (const auto& i : fullmatch2_) {
96  for (unsigned int j = 0; j < i->nMatches(); j++) {
97  if (i->getTracklet(j)->TCID() == tracklet->TCID()) {
98  trackstublist.push_back(i->getMatch(j).second);
99  }
100  }
101  }
102 
103  for (const auto& i : fullmatch3_) {
104  for (unsigned int j = 0; j < i->nMatches(); j++) {
105  if (i->getTracklet(j)->TCID() == tracklet->TCID()) {
106  trackstublist.push_back(i->getMatch(j).second);
107  }
108  }
109  }
110 
111  for (const auto& i : fullmatch4_) {
112  for (unsigned int j = 0; j < i->nMatches(); j++) {
113  if (i->getTracklet(j)->TCID() == tracklet->TCID()) {
114  trackstublist.push_back(i->getMatch(j).second);
115  }
116  }
117  }
118 
119  // For merge removal, loop through the resulting list of stubs to calculate their stubids
120  if (settings_.removalType() == "merge") {
121  for (const auto& it : trackstublist) {
122  int layer = it->layer().value() + 1; // Assume layer (1-6) stub first
123  if (it->layer().value() < 0) { // if disk stub, though...
124  layer = it->disk().value() + 10 * it->disk().value() / abs(it->disk().value()); //disk = +/- 11-15
125  }
126  stubidslist.push_back(std::make_pair(layer, it->phiregionaddress()));
127  }
128 
129  // And that's all we need! The rest is just for track fit (done in PurgeDuplicate)
130 
131  } else {
132  // Track fit only called here if not running duplicate removal
133  // before fit. (e.g. If skipping duplicate removal).
134  HybridFit hybridFitter(iSector_, settings_, globals_);
135  hybridFitter.Fit(tracklet, trackstublist);
136  }
137  }
138 }
139 
140 #else // Code for pure Tracklet algo.
141 
142 void FitTrack::trackFitChisq(Tracklet* tracklet, std::vector<const Stub*>&, std::vector<std::pair<int, int>>&) {
143  if (globals_->trackDerTable() == nullptr) {
144  TrackDerTable* derTablePtr = new TrackDerTable(settings_);
145 
146  derTablePtr->readPatternFile(settings_.fitPatternFile());
147  derTablePtr->fillTable();
148  if (settings_.debugTracklet()) {
149  edm::LogVerbatim("Tracklet") << "Number of entries in derivative table: " << derTablePtr->getEntries();
150  }
151  assert(derTablePtr->getEntries() != 0);
152 
153  globals_->trackDerTable() = derTablePtr;
154  }
155 
156  const TrackDerTable& derTable = *globals_->trackDerTable();
157 
158  //First step is to build list of layers and disks.
159  int layers[N_LAYER];
160  double r[N_LAYER];
161  unsigned int nlayers = 0; // layers with found stub-projections
162  int disks[N_DISK];
163  double z[N_DISK];
164  unsigned int ndisks = 0; // disks with found stub-projections
165 
166  // residuals for each stub
167  double phiresid[N_FITSTUB];
168  double zresid[N_FITSTUB];
169  double phiresidexact[N_FITSTUB];
170  double zresidexact[N_FITSTUB];
171  int iphiresid[N_FITSTUB];
172  int izresid[N_FITSTUB];
173  double alpha[N_FITSTUB];
174 
175  for (unsigned int i = 0; i < N_FITSTUB; i++) {
176  iphiresid[i] = 0;
177  izresid[i] = 0;
178  alpha[i] = 0.0;
179 
180  phiresid[i] = 0.0;
181  zresid[i] = 0.0;
182  phiresidexact[i] = 0.0;
183  zresidexact[i] = 0.0;
184  }
185 
186  std::bitset<N_LAYER> lmatches; //layer matches
187  std::bitset<N_DISK * 2> dmatches; //disk matches (2 per disk to separate 2S from PS)
188 
189  int mult = 1;
190 
191  unsigned int layermask = 0;
192  unsigned int diskmask = 0;
193  unsigned int alphaindex = 0;
194  unsigned int power = 1;
195 
196  double t = tracklet->t();
197  double rinv = tracklet->rinv();
198 
199  if (tracklet->isBarrel()) {
200  for (unsigned int l = 1; l <= N_LAYER; l++) {
201  if (l == (unsigned int)tracklet->layer() || l == (unsigned int)tracklet->layer() + 1) {
202  lmatches.set(N_LAYER - l);
203  layermask |= (1 << (N_LAYER - l));
204  layers[nlayers++] = l;
205  continue;
206  }
207  if (tracklet->match(l - 1)) {
208  const Residual& resid = tracklet->resid(l - 1);
209  lmatches.set(N_LAYER - l);
210  layermask |= (1 << (N_LAYER - l));
211  phiresid[nlayers] = resid.phiresidapprox();
212  zresid[nlayers] = resid.rzresidapprox();
213  phiresidexact[nlayers] = resid.phiresid();
214  zresidexact[nlayers] = resid.rzresid();
215  iphiresid[nlayers] = resid.fpgaphiresid().value();
216  izresid[nlayers] = resid.fpgarzresid().value();
217 
218  layers[nlayers++] = l;
219  }
220  }
221 
222  for (unsigned int d = 1; d <= N_DISK; d++) {
223  if (layermask & (1 << (d - 1)))
224  continue;
225 
226  if (mult == 1 << (3 * settings_.alphaBitsTable()))
227  continue;
228 
229  if (ndisks + nlayers >= N_FITSTUB)
230  continue;
231  if (tracklet->match(N_LAYER + d - 1)) {
232  const Residual& resid = tracklet->resid(N_LAYER + d - 1);
233  double pitch = settings_.stripPitch(resid.stubptr()->l1tstub()->isPSmodule());
234  if (std::abs(resid.stubptr()->l1tstub()->alpha(pitch)) < 1e-20) {
235  dmatches.set(2 * d - 1);
236  diskmask |= (1 << (2 * (N_DISK - d) + 1));
237  } else {
238  int ialpha = resid.stubptr()->alpha().value();
239  int nalpha = resid.stubptr()->alpha().nbits();
240  nalpha = nalpha - settings_.alphaBitsTable();
241  ialpha = (1 << (settings_.alphaBitsTable() - 1)) + (ialpha >> nalpha);
242 
243  alphaindex += ialpha * power;
245  dmatches.set(2 * (N_DISK - d));
246  diskmask |= (1 << (2 * (N_DISK - d)));
248  }
249  alpha[ndisks] = resid.stubptr()->l1tstub()->alpha(pitch);
250  phiresid[nlayers + ndisks] = resid.phiresidapprox();
251  zresid[nlayers + ndisks] = resid.rzresidapprox();
252  phiresidexact[nlayers + ndisks] = resid.phiresid();
253  zresidexact[nlayers + ndisks] = resid.rzresid();
254  iphiresid[nlayers + ndisks] = resid.fpgaphiresid().value();
255  izresid[nlayers + ndisks] = resid.fpgarzresid().value();
256 
257  disks[ndisks++] = d;
258  }
259  }
260 
261  if (settings_.writeMonitorData("HitPattern")) {
262  if (mult <= 1 << (3 * settings_.alphaBitsTable())) {
263  globals_->ofstream("hitpattern.txt")
264  << lmatches.to_string() << " " << dmatches.to_string() << " " << mult << endl;
265  }
266  }
267  }
268 
269  if (tracklet->isDisk()) {
270  for (unsigned int l = 1; l <= 2; l++) {
271  if (tracklet->match(l - 1)) {
272  lmatches.set(N_LAYER - l);
273 
274  layermask |= (1 << (N_LAYER - l));
275  const Residual& resid = tracklet->resid(l - 1);
276  phiresid[nlayers] = resid.phiresidapprox();
277  zresid[nlayers] = resid.rzresidapprox();
278  phiresidexact[nlayers] = resid.phiresid();
279  zresidexact[nlayers] = resid.rzresid();
280  iphiresid[nlayers] = resid.fpgaphiresid().value();
281  izresid[nlayers] = resid.fpgarzresid().value();
282 
283  layers[nlayers++] = l;
284  }
285  }
286 
287  for (int d1 = 1; d1 <= N_DISK; d1++) {
288  int d = d1;
289 
290  // skip F/B5 if there's already a L2 match
291  if (d == 5 and layermask & (1 << 4))
292  continue;
293 
294  if (tracklet->fpgat().value() < 0.0)
295  d = -d1;
296  if (d1 == abs(tracklet->disk()) || d1 == abs(tracklet->disk()) + 1) {
297  dmatches.set(2 * d1 - 1);
298  diskmask |= (1 << (2 * (N_DISK - d1) + 1));
299  alpha[ndisks] = 0.0;
300  disks[ndisks++] = d;
301  continue;
302  }
303 
304  if (ndisks + nlayers >= N_FITSTUB)
305  continue;
306  if (tracklet->match(N_LAYER + abs(d) - 1)) {
307  const Residual& resid = tracklet->resid(N_LAYER + abs(d) - 1);
308  double pitch = settings_.stripPitch(resid.stubptr()->l1tstub()->isPSmodule());
309  if (std::abs(resid.stubptr()->l1tstub()->alpha(pitch)) < 1e-20) {
310  dmatches.set(2 * d1 - 1);
311  diskmask |= (1 << (2 * (N_DISK - d1) + 1));
312  } else {
313  int ialpha = resid.stubptr()->alpha().value();
314  int nalpha = resid.stubptr()->alpha().nbits();
315  nalpha = nalpha - settings_.alphaBitsTable();
316  ialpha = (1 << (settings_.alphaBitsTable() - 1)) + (ialpha >> nalpha);
317 
318  alphaindex += ialpha * power;
320  dmatches.set(2 * (N_DISK - d1));
321  diskmask |= (1 << (2 * (N_DISK - d1)));
323  }
324 
325  alpha[ndisks] = resid.stubptr()->l1tstub()->alpha(pitch);
326  assert(std::abs(resid.phiresidapprox()) < 0.2);
327  phiresid[nlayers + ndisks] = resid.phiresidapprox();
328  zresid[nlayers + ndisks] = resid.rzresidapprox();
329  assert(std::abs(resid.phiresid()) < 0.2);
330  phiresidexact[nlayers + ndisks] = resid.phiresid();
331  zresidexact[nlayers + ndisks] = resid.rzresid();
332  iphiresid[nlayers + ndisks] = resid.fpgaphiresid().value();
333  izresid[nlayers + ndisks] = resid.fpgarzresid().value();
334 
335  disks[ndisks++] = d;
336  }
337  }
338  }
339 
340  if (tracklet->isOverlap()) {
341  for (unsigned int l = 1; l <= 2; l++) {
342  if (l == (unsigned int)tracklet->layer()) {
343  lmatches.set(N_LAYER - l);
344  layermask |= (1 << (N_LAYER - l));
345  layers[nlayers++] = l;
346  continue;
347  }
348  if (tracklet->match(l - 1)) {
349  lmatches.set(N_LAYER - l);
350  layermask |= (1 << (N_LAYER - l));
351  const Residual& resid = tracklet->resid(l - 1);
352  assert(std::abs(resid.phiresidapprox()) < 0.2);
353  phiresid[nlayers] = resid.phiresidapprox();
354  zresid[nlayers] = resid.rzresidapprox();
355  assert(std::abs(resid.phiresid()) < 0.2);
356  phiresidexact[nlayers] = resid.phiresid();
357  zresidexact[nlayers] = resid.rzresid();
358  iphiresid[nlayers] = resid.fpgaphiresid().value();
359  izresid[nlayers] = resid.fpgarzresid().value();
360 
361  layers[nlayers++] = l;
362  }
363  }
364 
365  for (unsigned int d1 = 1; d1 <= N_DISK; d1++) {
366  if (mult == 1 << (3 * settings_.alphaBitsTable()))
367  continue;
368  int d = d1;
369  if (tracklet->fpgat().value() < 0.0)
370  d = -d1;
371  if (d == tracklet->disk()) { //All seeds in PS modules
372  disks[ndisks] = tracklet->disk();
373  dmatches.set(2 * d1 - 1);
374  diskmask |= (1 << (2 * (N_DISK - d1) + 1));
375  ndisks++;
376  continue;
377  }
378 
379  if (ndisks + nlayers >= N_FITSTUB)
380  continue;
381  if (tracklet->match(N_LAYER + abs(d) - 1)) {
382  const Residual& resid = tracklet->resid(N_LAYER + abs(d) - 1);
383  double pitch = settings_.stripPitch(resid.stubptr()->l1tstub()->isPSmodule());
384  if (std::abs(resid.stubptr()->l1tstub()->alpha(pitch)) < 1e-20) {
385  dmatches.set(2 * (N_DISK - d1));
386  diskmask |= (1 << (2 * (N_DISK - d1) + 1));
387  FPGAWord tmp;
388  tmp.set(diskmask, 10);
389  } else {
390  int ialpha = resid.stubptr()->alpha().value();
391  int nalpha = resid.stubptr()->alpha().nbits();
392  nalpha = nalpha - settings_.alphaBitsTable();
393  ialpha = (1 << (settings_.alphaBitsTable() - 1)) + (ialpha >> nalpha);
394 
395  alphaindex += ialpha * power;
397  dmatches.set(2 * (N_DISK - d1));
398  diskmask |= (1 << (2 * (N_DISK - d1)));
399  FPGAWord tmp;
400  tmp.set(diskmask, 10);
402  }
403 
404  alpha[ndisks] = resid.stubptr()->l1tstub()->alpha(pitch);
405  assert(std::abs(resid.phiresidapprox()) < 0.2);
406  phiresid[nlayers + ndisks] = resid.phiresidapprox();
407  zresid[nlayers + ndisks] = resid.rzresidapprox();
408  assert(std::abs(resid.phiresid()) < 0.2);
409  phiresidexact[nlayers + ndisks] = resid.phiresid();
410  zresidexact[nlayers + ndisks] = resid.rzresid();
411  iphiresid[nlayers + ndisks] = resid.fpgaphiresid().value();
412  izresid[nlayers + ndisks] = resid.fpgarzresid().value();
413 
414  disks[ndisks++] = d;
415  }
416  }
417  }
418 
419  int rinvindex =
420  (1 << (settings_.nrinvBitsTable() - 1)) * rinv / settings_.rinvmax() + (1 << (settings_.nrinvBitsTable() - 1));
421  if (rinvindex < 0)
422  rinvindex = 0;
423  if (rinvindex >= (1 << settings_.nrinvBitsTable()))
424  rinvindex = (1 << settings_.nrinvBitsTable()) - 1;
425 
426  const TrackDer* derivatives = derTable.getDerivatives(layermask, diskmask, alphaindex, rinvindex);
427 
428  if (derivatives == nullptr) {
429  if (settings_.warnNoDer()) {
430  FPGAWord tmpl, tmpd;
431  tmpl.set(layermask, 6);
432  tmpd.set(diskmask, 10);
433  edm::LogVerbatim("Tracklet") << "No derivative for layermask, diskmask : " << layermask << " " << tmpl.str()
434  << " " << diskmask << " " << tmpd.str() << " eta = " << asinh(t);
435  }
436  return;
437  }
438 
439  double ttabi = TrackDerTable::tpar(settings_, diskmask, layermask);
440  if (t < 0.0)
441  ttabi = -ttabi;
442  double ttab = ttabi;
443 
444  if (settings_.debugTracklet()) {
445  edm::LogVerbatim("Tracklet") << "Doing trackfit in " << getName();
446  }
447 
448  int sign = 1;
449  if (t < 0.0)
450  sign = -1;
451 
452  double rstub[6];
453 
454  double realrstub[3];
455  realrstub[0] = -1.0;
456  realrstub[1] = -1.0;
457  realrstub[2] = -1.0;
458 
459  for (unsigned i = 0; i < nlayers; i++) {
460  r[i] = settings_.rmean(layers[i] - 1);
461  if (layers[i] == tracklet->layer()) {
462  if (tracklet->isOverlap()) {
463  realrstub[i] = tracklet->outerFPGAStub()->l1tstub()->r();
464  } else {
465  realrstub[i] = tracklet->innerFPGAStub()->l1tstub()->r();
466  }
467  }
468  if (layers[i] == tracklet->layer() + 1) {
469  realrstub[i] = tracklet->outerFPGAStub()->l1tstub()->r();
470  }
471  if (tracklet->match(layers[i] - 1) && layers[i] < 4) {
472  const Stub* stubptr = tracklet->resid(layers[i] - 1).stubptr();
473  realrstub[i] = stubptr->l1tstub()->r();
474  assert(std::abs(realrstub[i] - r[i]) < 5.0);
475  }
476  rstub[i] = r[i];
477  }
478  for (unsigned i = 0; i < ndisks; i++) {
479  z[i] = sign * settings_.zmean(abs(disks[i]) - 1);
480  rstub[i + nlayers] = z[i] / ttabi;
481  }
482 
483  double D[N_FITPARAM][N_FITSTUB * 2];
484  double MinvDt[N_FITPARAM][N_FITSTUB * 2];
485  int iD[N_FITPARAM][N_FITSTUB * 2];
486  int iMinvDt[N_FITPARAM][N_FITSTUB * 2];
487  double sigma[N_FITSTUB * 2];
488  double kfactor[N_FITSTUB * 2];
489 
490  unsigned int n = nlayers + ndisks;
491 
492  if (settings_.exactderivatives()) {
494  settings_, nlayers, r, ndisks, z, alpha, t, rinv, D, iD, MinvDt, iMinvDt, sigma, kfactor);
495  ttabi = t;
496  ttab = t;
497  } else {
500  settings_, nlayers, r, ndisks, z, alpha, t, rinv, D, iD, MinvDt, iMinvDt, sigma, kfactor);
501 
502  double MinvDtDummy[N_FITPARAM][N_FITSTUB * 2];
503  derivatives->fill(tracklet->fpgat().value(), MinvDtDummy, iMinvDt);
504  ttab = t;
505  } else {
506  derivatives->fill(tracklet->fpgat().value(), MinvDt, iMinvDt);
507  }
508  }
509 
510  if (!settings_.exactderivatives()) {
511  for (unsigned int i = 0; i < nlayers; i++) {
512  if (r[i] > settings_.rPS2S())
513  continue;
514  for (unsigned int ii = 0; ii < nlayers; ii++) {
515  if (r[ii] > settings_.rPS2S())
516  continue;
517 
518  double tder = derivatives->tdzcorr(i, ii);
519  double zder = derivatives->z0dzcorr(i, ii);
520 
521  double dr = realrstub[i] - r[i];
522 
523  MinvDt[2][2 * ii + 1] += dr * tder;
524  MinvDt[3][2 * ii + 1] += dr * zder;
525 
526  int itder = derivatives->itdzcorr(i, ii);
527  int izder = derivatives->iz0dzcorr(i, ii);
528 
529  int idr = dr / settings_.kr();
530 
531  iMinvDt[2][2 * ii + 1] += ((idr * itder) >> settings_.rcorrbits());
532  iMinvDt[3][2 * ii + 1] += ((idr * izder) >> settings_.rcorrbits());
533  }
534  }
535  }
536 
537  double rinvseed = tracklet->rinvapprox();
538  double phi0seed = tracklet->phi0approx();
539  double tseed = tracklet->tapprox();
540  double z0seed = tracklet->z0approx();
541 
542  double rinvseedexact = tracklet->rinv();
543  double phi0seedexact = tracklet->phi0();
544  double tseedexact = tracklet->t();
545  double z0seedexact = tracklet->z0();
546 
547  double chisqseedexact = 0.0;
548 
549  double delta[2 * N_FITSTUB];
550  double deltaexact[2 * N_FITSTUB];
551  int idelta[2 * N_FITSTUB];
552 
553  for (unsigned int i = 0; i < 2 * N_FITSTUB; i++) {
554  delta[i] = 0.0;
555  deltaexact[i] = 0.0;
556  idelta[i] = 0;
557  }
558 
559  int j = 0;
560 
561  for (unsigned int i = 0; i < n; i++) {
562  if (i >= nlayers) {
563  iphiresid[i] *= (t / ttabi);
564  phiresid[i] *= (t / ttab);
565  phiresidexact[i] *= (t / ttab);
566  }
567 
568  idelta[j] = iphiresid[i];
569  delta[j] = phiresid[i];
570  if (std::abs(phiresid[i]) > 0.2) {
571  edm::LogWarning("Tracklet") << getName() << " WARNING too large phiresid: " << phiresid[i] << " "
572  << phiresidexact[i];
573  }
574  assert(std::abs(phiresid[i]) < 1.0);
575  assert(std::abs(phiresidexact[i]) < 1.0);
576  deltaexact[j++] = phiresidexact[i];
577 
578  idelta[j] = izresid[i];
579  delta[j] = zresid[i];
580  deltaexact[j++] = zresidexact[i];
581 
582  chisqseedexact += (deltaexact[j - 2] * deltaexact[j - 2] + deltaexact[j - 1] * deltaexact[j - 1]);
583  }
584  assert(j <= 12);
585 
586  double drinv = 0.0;
587  double dphi0 = 0.0;
588  double dt = 0.0;
589  double dz0 = 0.0;
590 
591  double drinvexact = 0.0;
592  double dphi0exact = 0.0;
593  double dtexact = 0.0;
594  double dz0exact = 0.0;
595 
596  int idrinv = 0;
597  int idphi0 = 0;
598  int idt = 0;
599  int idz0 = 0;
600 
601  double drinv_cov = 0.0;
602  double dphi0_cov = 0.0;
603  double dt_cov = 0.0;
604  double dz0_cov = 0.0;
605 
606  double drinv_covexact = 0.0;
607  double dphi0_covexact = 0.0;
608  double dt_covexact = 0.0;
609  double dz0_covexact = 0.0;
610 
611  for (unsigned int j = 0; j < 2 * n; j++) {
612  drinv -= MinvDt[0][j] * delta[j];
613  dphi0 -= MinvDt[1][j] * delta[j];
614  dt -= MinvDt[2][j] * delta[j];
615  dz0 -= MinvDt[3][j] * delta[j];
616 
617  drinv_cov += D[0][j] * delta[j];
618  dphi0_cov += D[1][j] * delta[j];
619  dt_cov += D[2][j] * delta[j];
620  dz0_cov += D[3][j] * delta[j];
621 
622  drinvexact -= MinvDt[0][j] * deltaexact[j];
623  dphi0exact -= MinvDt[1][j] * deltaexact[j];
624  dtexact -= MinvDt[2][j] * deltaexact[j];
625  dz0exact -= MinvDt[3][j] * deltaexact[j];
626 
627  drinv_covexact += D[0][j] * deltaexact[j];
628  dphi0_covexact += D[1][j] * deltaexact[j];
629  dt_covexact += D[2][j] * deltaexact[j];
630  dz0_covexact += D[3][j] * deltaexact[j];
631 
632  idrinv += ((iMinvDt[0][j] * idelta[j]));
633  idphi0 += ((iMinvDt[1][j] * idelta[j]));
634  idt += ((iMinvDt[2][j] * idelta[j]));
635  idz0 += ((iMinvDt[3][j] * idelta[j]));
636 
637  if (false && j % 2 == 0) {
638  edm::LogVerbatim("Tracklet") << "DEBUG CHI2FIT " << j << " " << rinvseed << " + " << MinvDt[0][j] * delta[j]
639  << " " << MinvDt[0][j] << " " << delta[j] * rstub[j / 2] * 10000 << " \n"
640  << j << " " << tracklet->fpgarinv().value() * settings_.krinvpars() << " + "
641  << ((iMinvDt[0][j] * idelta[j])) * settings_.krinvpars() / 1024.0 << " "
642  << iMinvDt[0][j] * settings_.krinvpars() / settings_.kphi() / 1024.0 << " "
643  << idelta[j] * settings_.kphi() * rstub[j / 2] * 10000 << " " << idelta[j];
644  }
645  }
646 
647  double deltaChisqexact =
648  drinvexact * drinv_covexact + dphi0exact * dphi0_covexact + dtexact * dt_covexact + dz0exact * dz0_covexact;
649 
650  int irinvseed = tracklet->fpgarinv().value();
651  int iphi0seed = tracklet->fpgaphi0().value();
652 
653  int itseed = tracklet->fpgat().value();
654  int iz0seed = tracklet->fpgaz0().value();
655 
656  int irinvfit = irinvseed + ((idrinv + (1 << settings_.fitrinvbitshift())) >> settings_.fitrinvbitshift());
657 
658  int iphi0fit = iphi0seed + (idphi0 >> settings_.fitphi0bitshift());
659  int itfit = itseed + (idt >> settings_.fittbitshift());
660 
661  int iz0fit = iz0seed + (idz0 >> settings_.fitz0bitshift());
662 
663  double rinvfit = rinvseed - drinv;
664  double phi0fit = phi0seed - dphi0;
665 
666  double tfit = tseed - dt;
667  double z0fit = z0seed - dz0;
668 
669  double rinvfitexact = rinvseedexact - drinvexact;
670  double phi0fitexact = phi0seedexact - dphi0exact;
671 
672  double tfitexact = tseedexact - dtexact;
673  double z0fitexact = z0seedexact - dz0exact;
674 
675  double chisqfitexact = chisqseedexact + deltaChisqexact;
676 
678  bool NewChisqDebug = false;
679  double chisqfit = 0.0;
680  uint ichisqfit = 0;
681 
682  double phifactor;
683  double rzfactor;
684  double iphifactor;
685  double irzfactor;
686  int k = 0; // column index of D matrix
687 
688  if (NewChisqDebug) {
689  edm::LogVerbatim("Tracklet") << "OG chisq: \n"
690  << "drinv/cov = " << drinv << "/" << drinv_cov << " \n"
691  << "dphi0/cov = " << drinv << "/" << dphi0_cov << " \n"
692  << "dt/cov = " << drinv << "/" << dt_cov << " \n"
693  << "dz0/cov = " << drinv << "/" << dz0_cov << "\n";
694  std::string myout = "D[0][k]= ";
695  for (unsigned int i = 0; i < 2 * n; i++) {
696  myout += std::to_string(D[0][i]);
697  myout += ", ";
698  }
699  edm::LogVerbatim("Tracklet") << myout;
700  }
701 
702  for (unsigned int i = 0; i < n; i++) { // loop over stubs
703 
704  phifactor = rstub[k / 2] * delta[k] / sigma[k] + D[0][k] * drinv + D[1][k] * dphi0 + D[2][k] * dt + D[3][k] * dz0;
705  iphifactor = kfactor[k] * rstub[k / 2] * idelta[k] * (1 << settings_.chisqphifactbits()) / sigma[k] -
706  iD[0][k] * idrinv - iD[1][k] * idphi0 - iD[2][k] * idt - iD[3][k] * idz0;
707 
708  if (NewChisqDebug) {
709  edm::LogVerbatim("Tracklet") << "delta[k]/sigma = " << delta[k] / sigma[k] << " delta[k] = " << delta[k] << "\n"
710  << "sum = " << phifactor - delta[k] / sigma[k] << " drinvterm = " << D[0][k] * drinv
711  << " dphi0term = " << D[1][k] * dphi0 << " dtterm = " << D[2][k] * dt
712  << " dz0term = " << D[3][k] * dz0 << "\n phifactor = " << phifactor;
713  }
714 
715  chisqfit += phifactor * phifactor;
716  ichisqfit += iphifactor * iphifactor / (1 << (2 * settings_.chisqphifactbits() - 4));
717 
718  k++;
719 
720  rzfactor = delta[k] / sigma[k] + D[0][k] * drinv + D[1][k] * dphi0 + D[2][k] * dt + D[3][k] * dz0;
721  irzfactor = kfactor[k] * idelta[k] * (1 << settings_.chisqzfactbits()) / sigma[k] - iD[0][k] * idrinv -
722  iD[1][k] * idphi0 - iD[2][k] * idt - iD[3][k] * idz0;
723 
724  if (NewChisqDebug) {
725  edm::LogVerbatim("Tracklet") << "delta[k]/sigma = " << delta[k] / sigma[k] << " delta[k] = " << delta[k] << "\n"
726  << "sum = " << rzfactor - delta[k] / sigma[k] << " drinvterm = " << D[0][k] * drinv
727  << " dphi0term = " << D[1][k] * dphi0 << " dtterm = " << D[2][k] * dt
728  << " dz0term = " << D[3][k] * dz0 << "\n rzfactor = " << rzfactor;
729  }
730 
731  chisqfit += rzfactor * rzfactor;
732  ichisqfit += irzfactor * irzfactor / (1 << (2 * settings_.chisqzfactbits() - 4));
733 
734  k++;
735  }
736 
737  if (settings_.writeMonitorData("ChiSq")) {
738  globals_->ofstream("chisq.txt") << asinh(itfit * settings_.ktpars()) << " " << chisqfit << " " << ichisqfit / 16.0
739  << endl;
740  }
741 
742  // Chisquare per DOF capped out at 11 bits, so 15 is an educated guess
743  if (ichisqfit >= (1 << 15)) {
744  if (NewChisqDebug) {
745  edm::LogVerbatim("Tracklet") << "CHISQUARE (" << ichisqfit << ") LARGER THAN 11 BITS!";
746  }
747  ichisqfit = (1 << 15) - 1;
748  }
749 
750  // Eliminate lower bits to fit in 8 bits
751  ichisqfit = ichisqfit >> 7;
752  // Probably redundant... enforce 8 bit cap
753  if (ichisqfit >= (1 << 8))
754  ichisqfit = (1 << 8) - 1;
755 
756  double phicrit = phi0fit - asin(0.5 * settings_.rcrit() * rinvfit);
757  bool keep = (phicrit > settings_.phicritmin()) && (phicrit < settings_.phicritmax());
758 
759  if (!keep) {
760  return;
761  }
762 
763  // NOTE: setFitPars in Tracklet.h now accepts chi2 r-phi and chi2 r-z values. This class only has access
764  // to the composite chi2. When setting fit parameters on a tracklet, this places all of the chi2 into the
765  // r-phi fit, and sets the r-z fit value to zero.
766  //
767  // This is also true for the call to setFitPars in trackFitFake.
768  tracklet->setFitPars(rinvfit,
769  phi0fit,
770  0.0,
771  tfit,
772  z0fit,
773  chisqfit,
774  0.0,
775  rinvfitexact,
776  phi0fitexact,
777  0.0,
778  tfitexact,
779  z0fitexact,
780  chisqfitexact,
781  0.0,
782  irinvfit,
783  iphi0fit,
784  0,
785  itfit,
786  iz0fit,
787  ichisqfit,
788  0,
789  0);
790 }
791 
792 #endif
793 
794 void FitTrack::trackFitFake(Tracklet* tracklet, std::vector<const Stub*>&, std::vector<std::pair<int, int>>&) {
795  tracklet->setFitPars(tracklet->rinvapprox(),
796  tracklet->phi0approx(),
797  tracklet->d0approx(),
798  tracklet->tapprox(),
799  tracklet->z0approx(),
800  0.0,
801  0.0,
802  tracklet->rinv(),
803  tracklet->phi0(),
804  tracklet->d0(),
805  tracklet->t(),
806  tracklet->z0(),
807  0.0,
808  0.0,
809  tracklet->fpgarinv().value(),
810  tracklet->fpgaphi0().value(),
811  tracklet->fpgad0().value(),
812  tracklet->fpgat().value(),
813  tracklet->fpgaz0().value(),
814  0,
815  0,
816  0);
817  return;
818 }
819 
820 std::vector<Tracklet*> FitTrack::orderedMatches(vector<FullMatchMemory*>& fullmatch) {
821  std::vector<Tracklet*> tmp;
822 
823  std::vector<unsigned int> indexArray;
824  for (auto& imatch : fullmatch) {
825  //check that we have correct order
826  if (imatch->nMatches() > 1) {
827  for (unsigned int j = 0; j < imatch->nMatches() - 1; j++) {
828  assert(imatch->getTracklet(j)->TCID() <= imatch->getTracklet(j + 1)->TCID());
829  }
830  }
831 
832  if (settings_.debugTracklet() && imatch->nMatches() != 0) {
833  edm::LogVerbatim("Tracklet") << "orderedMatches: " << imatch->getName() << " " << imatch->nMatches();
834  }
835 
836  indexArray.push_back(0);
837  }
838 
839  int bestIndex = -1;
840  do {
841  int bestTCID = -1;
842  bestIndex = -1;
843  for (unsigned int i = 0; i < fullmatch.size(); i++) {
844  if (indexArray[i] >= fullmatch[i]->nMatches()) {
845  //skip as we were at the end
846  continue;
847  }
848  int TCID = fullmatch[i]->getTracklet(indexArray[i])->TCID();
849  if (TCID < bestTCID || bestTCID < 0) {
850  bestTCID = TCID;
851  bestIndex = i;
852  }
853  }
854  if (bestIndex != -1) {
855  tmp.push_back(fullmatch[bestIndex]->getTracklet(indexArray[bestIndex]));
856  indexArray[bestIndex]++;
857  }
858  } while (bestIndex != -1);
859 
860  for (unsigned int i = 0; i < tmp.size(); i++) {
861  if (i > 0) {
862  //This allows for equal TCIDs. This means that we can e.g. have a track seeded in L1L2 that projects to both L3 and D4.
863  //The algorithm will pick up the first hit and drop the second.
864  if (tmp[i - 1]->TCID() > tmp[i]->TCID()) {
865  edm::LogVerbatim("Tracklet") << "Wrong TCID ordering in " << getName() << " : " << tmp[i - 1]->TCID() << " "
866  << tmp[i]->TCID();
867  }
868  }
869  }
870 
871  return tmp;
872 }
873 
874 // Adds the fitted track to the output memories to be used by pure Tracklet algo.
875 // (Also used by Hybrid algo with non-exact Old KF emulation)
876 // Also create output streams, that bypass these memories, (so can include gaps in time),
877 // to be used by Hybrid case with exact New KF emulation.
878 
879 void FitTrack::execute(deque<string>& streamTrackRaw,
880  vector<deque<StubStreamData>>& streamsStubRaw,
881  unsigned int iSector) {
882  // merge
883  const std::vector<Tracklet*>& matches1 = orderedMatches(fullmatch1_);
884  const std::vector<Tracklet*>& matches2 = orderedMatches(fullmatch2_);
885  const std::vector<Tracklet*>& matches3 = orderedMatches(fullmatch3_);
886  const std::vector<Tracklet*>& matches4 = orderedMatches(fullmatch4_);
887 
888  iSector_ = iSector;
889 
890  if (settings_.debugTracklet() && (matches1.size() + matches2.size() + matches3.size() + matches4.size()) > 0) {
891  for (auto& imatch : fullmatch1_) {
892  edm::LogVerbatim("Tracklet") << imatch->getName() << " " << imatch->nMatches();
893  }
894  edm::LogVerbatim("Tracklet") << getName() << " matches : " << matches1.size() << " " << matches2.size() << " "
895  << matches3.size() << " " << matches4.size();
896  }
897 
898  unsigned int indexArray[4];
899  for (unsigned int i = 0; i < 4; i++) {
900  indexArray[i] = 0;
901  }
902 
903  unsigned int countAll = 0;
904  unsigned int countFit = 0;
905 
906  Tracklet* bestTracklet = nullptr;
907  do {
908  countAll++;
909  bestTracklet = nullptr;
910 
911  if (indexArray[0] < matches1.size()) {
912  if (bestTracklet == nullptr) {
913  bestTracklet = matches1[indexArray[0]];
914  } else {
915  if (matches1[indexArray[0]]->TCID() < bestTracklet->TCID())
916  bestTracklet = matches1[indexArray[0]];
917  }
918  }
919 
920  if (indexArray[1] < matches2.size()) {
921  if (bestTracklet == nullptr) {
922  bestTracklet = matches2[indexArray[1]];
923  } else {
924  if (matches2[indexArray[1]]->TCID() < bestTracklet->TCID())
925  bestTracklet = matches2[indexArray[1]];
926  }
927  }
928 
929  if (indexArray[2] < matches3.size()) {
930  if (bestTracklet == nullptr) {
931  bestTracklet = matches3[indexArray[2]];
932  } else {
933  if (matches3[indexArray[2]]->TCID() < bestTracklet->TCID())
934  bestTracklet = matches3[indexArray[2]];
935  }
936  }
937 
938  if (indexArray[3] < matches4.size()) {
939  if (bestTracklet == nullptr) {
940  bestTracklet = matches4[indexArray[3]];
941  } else {
942  if (matches4[indexArray[3]]->TCID() < bestTracklet->TCID())
943  bestTracklet = matches4[indexArray[3]];
944  }
945  }
946 
947  if (bestTracklet == nullptr)
948  break;
949 
950  //Counts total number of matched hits
951  int nMatches = 0;
952 
953  //Counts unique hits in each layer
954  int nMatchesUniq = 0;
955  bool match = false;
956 
957  while (indexArray[0] < matches1.size() && matches1[indexArray[0]] == bestTracklet) {
958  indexArray[0]++;
959  nMatches++;
960  match = true;
961  }
962 
963  if (match)
964  nMatchesUniq++;
965  match = false;
966 
967  while (indexArray[1] < matches2.size() && matches2[indexArray[1]] == bestTracklet) {
968  indexArray[1]++;
969  nMatches++;
970  match = true;
971  }
972 
973  if (match)
974  nMatchesUniq++;
975  match = false;
976 
977  while (indexArray[2] < matches3.size() && matches3[indexArray[2]] == bestTracklet) {
978  indexArray[2]++;
979  nMatches++;
980  match = true;
981  }
982 
983  if (match)
984  nMatchesUniq++;
985  match = false;
986 
987  while (indexArray[3] < matches4.size() && matches4[indexArray[3]] == bestTracklet) {
988  indexArray[3]++;
989  nMatches++;
990  match = true;
991  }
992 
993  if (match)
994  nMatchesUniq++;
995 
996  if (settings_.debugTracklet()) {
997  edm::LogVerbatim("Tracklet") << getName() << " : nMatches = " << nMatches << " nMatchesUniq = " << nMatchesUniq
998  << " " << asinh(bestTracklet->t());
999  }
1000 
1001  std::vector<const Stub*> trackstublist;
1002  std::vector<std::pair<int, int>> stubidslist;
1003  // Track Builder cut of >= 4 layers with stubs.
1004  if ((bestTracklet->getISeed() >= (int)N_SEED_PROMPT && nMatchesUniq >= 1) ||
1005  nMatchesUniq >= 2) { //For seeds index >=8 (triplet seeds), there are three stubs associated from start.
1006  countFit++;
1007 
1008 #ifdef USEHYBRID
1009  if (settings_.fakefit()) {
1010  trackFitFake(bestTracklet, trackstublist, stubidslist);
1011  } else {
1012  trackFitKF(bestTracklet, trackstublist, stubidslist);
1013  }
1014 #else
1015  if (settings_.fakefit()) {
1016  trackFitFake(bestTracklet, trackstublist, stubidslist);
1017  } else {
1018  trackFitChisq(bestTracklet, trackstublist, stubidslist);
1019  }
1020 #endif
1021 
1022  if (settings_.removalType() == "merge") {
1023  trackfit_->addStubList(trackstublist);
1024  trackfit_->addStubidsList(stubidslist);
1025  bestTracklet->setTrackIndex(trackfit_->nTracks());
1026  trackfit_->addTrack(bestTracklet);
1027  } else if (bestTracklet->fit()) {
1028  assert(trackfit_ != nullptr);
1029  if (settings_.writeMonitorData("Seeds")) {
1030  ofstream fout("seeds.txt", ofstream::app);
1031  fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_"
1032  << " " << bestTracklet->getISeed() << endl;
1033  fout.close();
1034  }
1035  bestTracklet->setTrackIndex(trackfit_->nTracks());
1036  trackfit_->addTrack(bestTracklet);
1037  }
1038  }
1039 
1040  // store bit and clock accurate TB output
1041  if (settings_.storeTrackBuilderOutput() && bestTracklet) {
1042  // add gap if TrackBuilder rejected track (due to too few stub layers).
1043  if (!bestTracklet->fit()) {
1044  static const string invalid = "0";
1045  streamTrackRaw.emplace_back(invalid);
1046  for (auto& stream : streamsStubRaw)
1047  stream.emplace_back(StubStreamData());
1048  continue;
1049  }
1050  // convert Track word
1051  const string rinv = bestTracklet->fpgarinv().str();
1052  const string phi0 = bestTracklet->fpgaphi0().str();
1053  const string z0 = bestTracklet->fpgaz0().str();
1054  const string t = bestTracklet->fpgat().str();
1055  const int seedType = bestTracklet->getISeed();
1056  const string seed = TTBV(seedType, settings_.nbitsseed()).str();
1057  const string valid("1");
1058  streamTrackRaw.emplace_back(valid + seed + rinv + phi0 + z0 + t);
1059 
1060  // convert projected stubs
1061  unsigned int ihit(0);
1062  for (unsigned int ilayer = 0; ilayer < N_LAYER + N_DISK; ilayer++) {
1063  if (bestTracklet->match(ilayer)) {
1064  const Residual& resid = bestTracklet->resid(ilayer);
1065  // create bit accurate 64 bit word
1066  string r = resid.stubptr()->r().str();
1067  const string& phi = resid.fpgaphiresid().str();
1068  const string& rz = resid.fpgarzresid().str();
1069  const L1TStub* stub = resid.stubptr()->l1tstub();
1070  static constexpr int widthDisk2Sidentifier = 8;
1071  bool disk2S = (stub->disk() != 0) && (stub->isPSmodule() == 0);
1072  if (disk2S)
1073  r = string(widthDisk2Sidentifier, '0') + r;
1074  const string& stubId = resid.fpgastubid().str();
1075  // store seed, L1TStub, and bit accurate 64 bit word in clock accurate output
1076  streamsStubRaw[ihit++].emplace_back(seedType, *stub, valid + stubId + r + phi + rz);
1077  }
1078  }
1079  // convert seed stubs
1080  const string& stubId0 = bestTracklet->innerFPGAStub()->stubindex().str();
1081  const L1TStub* stub0 = bestTracklet->innerFPGAStub()->l1tstub();
1082  streamsStubRaw[ihit++].emplace_back(seedType, *stub0, valid + stubId0);
1083  const string& stubId1 = bestTracklet->outerFPGAStub()->stubindex().str();
1084  const L1TStub* stub1 = bestTracklet->outerFPGAStub()->l1tstub();
1085  streamsStubRaw[ihit++].emplace_back(seedType, *stub1, valid + stubId1);
1086  // fill all layers that have no stubs with gaps
1087  while (ihit < streamsStubRaw.size()) {
1088  streamsStubRaw[ihit++].emplace_back();
1089  }
1090  }
1091 
1092  } while (bestTracklet != nullptr && countAll < settings_.maxStep("TB"));
1093 
1094  if (settings_.writeMonitorData("FT")) {
1095  globals_->ofstream("fittrack.txt") << getName() << " " << countAll << " " << countFit << endl;
1096  }
1097 }
Log< level::Info, true > LogVerbatim
float dt
Definition: AMPTWrapper.h:136
const FPGAWord & fpgarinv() const
Definition: Tracklet.h:132
int chisqzfactbits() const
Definition: Settings.h:428
double t() const
Definition: Tracklet.h:123
double phiresid() const
Definition: Residual.h:47
double tdzcorr(int i, int j) const
Definition: TrackDer.h:66
int nrinvBitsTable() const
Definition: Settings.h:231
constexpr int N_DISK
Definition: Settings.h:26
unsigned int maxStep(std::string module) const
Definition: Settings.h:125
double z0dzcorr(int i, int j) const
Definition: TrackDer.h:67
const FPGAWord & r() const
Definition: Stub.h:65
double rPS2S() const
Definition: Settings.h:369
std::string name_
Definition: ProcessBase.h:42
const Residual & resid(unsigned int layerdisk)
Definition: Tracklet.h:110
static void calculateDerivatives(Settings const &settings, unsigned int nlayers, double r[N_LAYER], unsigned int ndisks, double z[N_DISK], double alpha[N_DISK], double t, double rinv, double D[N_FITPARAM][N_FITSTUB *2], int iD[N_FITPARAM][N_FITSTUB *2], double MinvDt[N_FITPARAM][N_FITSTUB *2], int iMinvDt[N_FITPARAM][N_FITSTUB *2], double sigma[N_FITSTUB *2], double kfactor[N_FITSTUB *2])
int disk() const
Definition: Tracklet.cc:801
bool isOverlap() const
Definition: Tracklet.h:203
int alphaBitsTable() const
Definition: Settings.h:230
tuple disks
Definition: alignBH_cfg.py:13
Bit vector used by Track Trigger emulators. Mainly used to convert integers into arbitrary (within ma...
Definition: TTBV.h:20
const FPGAWord & fpgastubid() const
Definition: Residual.h:42
void addTrack(Tracklet *tracklet)
std::vector< TrackletParametersMemory * > seedtracklet_
Definition: FitTrack.h:50
double rzresid() const
Definition: Residual.h:52
const FPGAWord & fpgaphiresid() const
Definition: Residual.h:32
int getEntries() const
Definition: TrackDerTable.h:37
const FPGAWord & fpgaz0() const
Definition: Tracklet.h:136
double rinvmax() const
Definition: Settings.h:226
unsigned int nbitsseed() const
Definition: Settings.h:320
const FPGAWord & fpgaphi0() const
Definition: Tracklet.h:133
Settings const & settings_
Definition: ProcessBase.h:44
double phi0() const
Definition: Tracklet.h:121
Globals * globals_
Definition: ProcessBase.h:45
bool exactderivatives() const
Definition: Settings.h:245
bool exactderivativesforfloating() const
Definition: Settings.h:246
double ktpars() const
Definition: Settings.h:437
bool writetrace() const
Definition: Settings.h:195
void trackFitKF(Tracklet *tracklet, std::vector< const Stub *> &trackstublist, std::vector< std::pair< int, int >> &stubidslist)
double krinvpars() const
Definition: Settings.h:432
void trackFitChisq(Tracklet *tracklet, std::vector< const Stub *> &, std::vector< std::pair< int, int >> &)
Definition: FitTrack.cc:142
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t stream
const FPGAWord & fpgat() const
Definition: Tracklet.h:135
void setTrackIndex(int index)
Definition: Tracklet.cc:832
assert(be >=bs)
void addOutput(MemoryBase *memory, std::string output) override
Definition: FitTrack.cc:19
int itdzcorr(int i, int j) const
Definition: TrackDer.h:78
int TCID() const
Definition: Tracklet.h:214
void addStubList(std::vector< const Stub *> stublist)
int fittbitshift() const
Definition: Settings.h:421
static std::string to_string(const XMLCh *ch)
int fitz0bitshift() const
Definition: Settings.h:422
static std::string const input
Definition: EdmProvDump.cc:50
unsigned int isPSmodule() const
Definition: L1TStub.h:103
TrackFitMemory * trackfit_
Definition: FitTrack.h:58
std::vector< Tracklet * > orderedMatches(std::vector< FullMatchMemory *> &fullmatch)
Definition: FitTrack.cc:820
bool fit() const
Definition: Tracklet.h:197
int rcorrbits() const
Definition: Settings.h:425
std::vector< FullMatchMemory * > fullmatch4_
Definition: FitTrack.h:54
const FPGAWord & fpgarzresid() const
Definition: Residual.h:37
double phi0approx() const
Definition: Tracklet.h:127
static double tpar(Settings const &settings, int diskmask, int layermask)
double rmean(unsigned int iLayer) const
Definition: Settings.h:173
double rinv() const
Definition: Tracklet.h:120
double phicritmin() const
Definition: Settings.h:338
const Stub * outerFPGAStub() const
Definition: Tracklet.h:65
constexpr unsigned int N_TRKLSEED
Definition: Settings.h:1100
double stripPitch(bool isPSmodule) const
Definition: Settings.h:291
void addInput(MemoryBase *memory, std::string input) override
Definition: FitTrack.cc:34
void fill(int t, double MinvDt[N_FITPARAM][N_FITSTUB *2], int iMinvDt[N_FITPARAM][N_FITSTUB *2]) const
Definition: TrackDer.cc:42
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int value() const
Definition: FPGAWord.h:24
unsigned int nTracks() const
const Stub * middleFPGAStub() const
Definition: Tracklet.h:63
std::string removalType() const
Definition: Settings.h:251
double tapprox() const
Definition: Tracklet.h:129
bool fakefit() const
Definition: Settings.h:255
unsigned int iSector_
Definition: FitTrack.h:56
L1TStub * l1tstub()
Definition: Stub.h:83
int disk() const
Definition: L1TStub.h:46
const TrackDer * getDerivatives(int index) const
Definition: TrackDerTable.h:24
std::string const & fitPatternFile() const
Definition: Settings.h:75
const FPGAWord & fpgad0() const
Definition: Tracklet.h:134
bool writeMonitorData(std::string module) const
Definition: Settings.h:118
bool isBarrel() const
Definition: Tracklet.h:202
double phiresidapprox() const
Definition: Residual.h:57
d
Definition: ztail.py:151
double zmean(unsigned int iDisk) const
Definition: Settings.h:176
const FPGAWord & stubindex() const
Definition: Stub.h:72
ii
Definition: cuy.py:589
std::vector< FullMatchMemory * > fullmatch1_
Definition: FitTrack.h:51
double alpha(double pitch) const
Definition: L1TStub.cc:79
bool debugTracklet() const
Definition: Settings.h:194
double rzresidapprox() const
Definition: Residual.h:62
double rinvapprox() const
Definition: Tracklet.h:126
double d0() const
Definition: Tracklet.h:122
double kr() const
Definition: Settings.h:351
constexpr unsigned int power(unsigned int base, unsigned int exponent)
std::vector< FullMatchMemory * > fullmatch2_
Definition: FitTrack.h:52
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
bool storeTrackBuilderOutput() const
Definition: Settings.h:257
double z0approx() const
Definition: Tracklet.h:130
int nbits() const
Definition: FPGAWord.h:25
void execute(std::deque< std::string > &streamTrackRaw, std::vector< std::deque< StubStreamData >> &stubStream, unsigned int iSector)
Definition: FitTrack.cc:879
double rinv(double phi1, double phi2, double r1, double r2)
Definition: Util.h:66
int iz0dzcorr(int i, int j) const
Definition: TrackDer.h:79
int getISeed() const
Definition: Tracklet.cc:820
void set(int value, int nbits, bool positive=true, int line=-1, const char *file=nullptr)
Definition: FPGAWord.cc:14
constexpr unsigned int N_FITPARAM
Definition: Settings.h:1104
void readPatternFile(std::string fileName)
const Stub * stubptr() const
Definition: Residual.h:67
double d0approx() const
Definition: Tracklet.h:128
double z0() const
Definition: Tracklet.h:124
const FPGAWord & alpha() const
Definition: Stub.h:70
double r() const
Definition: L1TStub.h:60
bool match(unsigned int layerdisk)
Definition: Tracklet.h:105
int isDisk() const
Definition: Tracklet.h:204
std::ofstream & ofstream(std::string fname)
Definition: Globals.cc:44
constexpr unsigned int N_FITSTUB
Definition: Settings.h:1105
int fitrinvbitshift() const
Definition: Settings.h:419
void setFitPars(double rinvfit, double phi0fit, double d0fit, double tfit, double z0fit, double chisqrphifit, double chisqrzfit, double rinvfitexact, double phi0fitexact, double d0fitexact, double tfitexact, double z0fitexact, double chisqrphifitexact, double chisqrzfitexact, int irinvfit, int iphi0fit, int id0fit, int itfit, int iz0fit, int ichisqrphifit, int ichisqrzfit, int hitpattern, const std::vector< const L1TStub *> &l1stubs=std::vector< const L1TStub *>())
Definition: Tracklet.cc:523
double phicritmax() const
Definition: Settings.h:339
int chisqphifactbits() const
Definition: Settings.h:427
TrackDerTable *& trackDerTable()
Definition: Globals.h:40
Definition: output.py:1
void addStubidsList(std::vector< std::pair< int, int >> stubidslist)
Log< level::Warning, false > LogWarning
const unsigned int kfactor
void trackFitFake(Tracklet *tracklet, std::vector< const Stub *> &, std::vector< std::pair< int, int >> &)
Definition: FitTrack.cc:794
int layer() const
Definition: Tracklet.cc:792
static constexpr float d1
double kphi() const
Definition: Settings.h:345
#define str(s)
tmp
align.sh
Definition: createJobs.py:716
int fitphi0bitshift() const
Definition: Settings.h:420
bool doKF() const
Definition: Settings.h:253
std::string const & getName() const
Definition: ProcessBase.h:22
std::vector< FullMatchMemory * > fullmatch3_
Definition: FitTrack.h:53
std::string str() const
Definition: FPGAWord.cc:54
bool warnNoDer() const
Definition: Settings.h:198
constexpr unsigned int N_SEED_PROMPT
Definition: Settings.h:29
constexpr int N_LAYER
Definition: Settings.h:25
double rcrit() const
Definition: Settings.h:334
const Stub * innerFPGAStub() const
Definition: Tracklet.h:61