CMS 3D CMS Logo

KFbase.cc
Go to the documentation of this file.
1 
4 
14 #include "TMatrixD.h"
15 
16 #include <algorithm>
17 #include <functional>
18 #include <fstream>
19 #include <iomanip>
20 #include <atomic>
21 #include <sstream>
22 
23 using namespace std;
24 
25 namespace tmtt {
26 
27  /* Initialize cfg parameters */
28 
29  KFbase::KFbase(const Settings *settings, const uint nHelixPar, const string &fitterName, const uint nMeas)
30  : TrackFitGeneric(settings, fitterName) {
31  nHelixPar_ = nHelixPar;
32  nMeas_ = nMeas;
33  numEtaRegions_ = settings->numEtaRegions();
34  }
35 
36  /* Do track fit */
37 
38  L1fittedTrack KFbase::fit(const L1track3D &l1track3D) {
39  iPhiSec_ = l1track3D.iPhiSec();
40  iEtaReg_ = l1track3D.iEtaReg();
41  resetStates();
42  numUpdateCalls_ = 0;
43 
44  vector<Stub *> stubs = l1track3D.stubs();
45 
46  auto orderByLayer = [](const Stub *a, const Stub *b) { return bool(a->layerId() < b->layerId()); };
47  sort(stubs.begin(), stubs.end(), orderByLayer); // Makes debug printout pretty.
48 
49  //TP
50  const TP *tpa(nullptr);
51  if (l1track3D.matchedTP()) {
52  tpa = l1track3D.matchedTP();
53  }
54  tpa_ = tpa;
55 
56  //track information dump
57  if (settings_->kalmanDebugLevel() >= 1) {
58  PrintL1trk() << "===============================================================================";
59  std::stringstream text;
60  text << std::fixed << std::setprecision(4);
61  text << "Input track cand: [phiSec,etaReg]=[" << l1track3D.iPhiSec() << "," << l1track3D.iEtaReg() << "]";
62  text << " HT(m,c)=(" << l1track3D.cellLocationHT().first << "," << l1track3D.cellLocationHT().second
63  << ") q/pt=" << l1track3D.qOverPt() << " tanL=" << l1track3D.tanLambda() << " z0=" << l1track3D.z0()
64  << " phi0=" << l1track3D.phi0() << " nStubs=" << l1track3D.numStubs() << " d0=" << l1track3D.d0();
65  PrintL1trk() << text.str();
66  if (not settings_->hybrid())
67  printTP(tpa);
68  if (settings_->kalmanDebugLevel() >= 2) {
69  printStubLayers(stubs, l1track3D.iEtaReg());
70  printStubs(stubs);
71  }
72  }
73 
74  //Kalman Filter
75  const KalmanState *cand = doKF(l1track3D, stubs, tpa);
76 
77  //return L1fittedTrk for the selected state (if KF produced one it was happy with).
78  if (cand != nullptr) {
79  // Get track helix params.
80  TVectorD trackPars = trackParams(cand);
81  double d0 = (nHelixPar_ == 5) ? trackPars[D0] : 0.;
82 
83  L1fittedTrack fitTrk(settings_,
84  &l1track3D,
85  cand->stubs(),
86  cand->hitPattern(),
87  trackPars[QOVERPT],
88  d0,
89  trackPars[PHI0],
90  trackPars[Z0],
91  trackPars[T],
92  cand->chi2rphi(),
93  cand->chi2rz(),
94  nHelixPar_);
95 
96  // Store supplementary info, specific to KF fitter.
97  fitTrk.setInfoKF(cand->nSkippedLayers(), numUpdateCalls_);
98 
99  // If doing 5 parameter fit, optionally also calculate helix params & chi2 with beam-spot constraint applied,
100  // and store inside L1fittedTrack object.
102  if (nHelixPar_ == 5) {
103  double chi2rphi_bcon = 0.;
104  TVectorD trackPars_bcon = trackParams_BeamConstr(cand, chi2rphi_bcon);
105  fitTrk.setBeamConstr(trackPars_bcon[QOVERPT], trackPars_bcon[PHI0], chi2rphi_bcon);
106  }
107  }
108 
109  // Fitted track params must lie in same sector as HT originally found track in.
110  if (!settings_->hybrid()) { // consistentSector() function not yet working for Hybrid.
111 
112  // Bodge to take into account digitisation in sector consistency check.
113  if (settings_->enableDigitize())
114  fitTrk.digitizeTrack("KF4ParamsComb");
115 
116  if (!fitTrk.consistentSector()) {
117  if (settings_->kalmanDebugLevel() >= 1)
118  PrintL1trk() << "Track rejected by sector consistency test";
119  L1fittedTrack rejectedTrk;
120  return rejectedTrk;
121  }
122  }
123 
124  return fitTrk;
125 
126  } else { // Track rejected by fitter
127 
128  if (settings_->kalmanDebugLevel() >= 1) {
129  bool goodTrack = (tpa && tpa->useForAlgEff()); // Matches truth particle.
130  if (goodTrack) {
131  int tpin = tpa->index();
132  PrintL1trk() << "TRACK LOST: eta=" << l1track3D.iEtaReg() << " pt=" << l1track3D.pt() << " tp=" << tpin;
133 
134  for (auto stub : stubs) {
135  int kalmanLay =
136  this->kalmanLayer(l1track3D.iEtaReg(), stub->layerIdReduced(), stub->barrel(), stub->r(), stub->z());
137  std::stringstream text;
138  text << std::fixed << std::setprecision(4);
139  text << " Stub: lay_red=" << stub->layerIdReduced() << " KFlay=" << kalmanLay << " r=" << stub->r()
140  << " z=" << stub->z() << " assoc TPs =";
141  for (const TP *tp_i : stub->assocTPs())
142  text << " " << tp_i->index();
143  PrintL1trk() << text.str();
144  if (stub->assocTPs().empty())
145  PrintL1trk() << " none";
146  }
147  PrintL1trk() << "=====================";
148  }
149  }
150 
151  //dump on the missed TP for efficiency calculation.
152  if (settings_->kalmanDebugLevel() >= 3) {
153  if (tpa && tpa->useForAlgEff()) {
154  PrintL1trk() << "TP for eff. missed addr. index : " << tpa << " " << tpa->index();
155  printStubs(stubs);
156  }
157  }
158 
159  L1fittedTrack rejectedTrk;
160  return rejectedTrk;
161  }
162  }
163 
164  /* Do track fit (internal function) */
165 
166  const KalmanState *KFbase::doKF(const L1track3D &l1track3D, const vector<Stub *> &stubs, const TP *tpa) {
167  const KalmanState *finished_state = nullptr;
168 
169  map<unsigned int, const KalmanState *, std::greater<unsigned int>>
170  best_state_by_nstubs; // Best state (if any) for each viable no. of stubs on track value.
171 
172  // seed helix params & their covariance.
173  TVectorD x0 = seedX(l1track3D);
174  TMatrixD pxx0 = seedC(l1track3D);
175  TMatrixD K(nHelixPar_, 2);
176  TMatrixD dcov(2, 2);
177 
178  const KalmanState *state0 = mkState(l1track3D, 0, -1, nullptr, x0, pxx0, K, dcov, nullptr, 0, 0);
179 
180  // internal containers - i.e. the state FIFO. Contains estimate of helix params in last/next layer, with multiple entries if there were multiple stubs, yielding multiple states.
181  vector<const KalmanState *> new_states;
182  vector<const KalmanState *> prev_states;
183  prev_states.push_back(state0);
184 
185  // Get dead layers, if any.
186  bool remove2PSCut = settings_->kalmanRemove2PScut();
187  set<unsigned> kfDeadLayers = kalmanDeadLayers(remove2PSCut);
188 
189  // arrange stubs into Kalman layers according to eta region
190  int etaReg = l1track3D.iEtaReg();
191  map<int, vector<Stub *>> layerStubs;
192 
193  for (auto stub : stubs) {
194  // Get Kalman encoded layer ID for this stub.
195  int kalmanLay = this->kalmanLayer(etaReg, stub->layerIdReduced(), stub->barrel(), stub->r(), stub->z());
196 
197  constexpr unsigned int invalidKFlayer = 7;
198  if (kalmanLay != invalidKFlayer) {
199  if (layerStubs[kalmanLay].size() < settings_->kalmanMaxStubsPerLayer()) {
200  layerStubs[kalmanLay].push_back(stub);
201  } else {
202  // If too many stubs, FW keeps the last stub.
203  layerStubs[kalmanLay].back() = stub;
204  }
205  }
206  }
207 
208  // iterate using state->nextLayer() to determine next Kalman layer(s) to add stubs from
209  constexpr unsigned int nTypicalLayers = 6; // Number of tracker layers a typical track can pass through.
210  // If user asked to add up to 7 layers to track, increase number of iterations by 1.
211  const unsigned int maxIterations = std::max(nTypicalLayers, settings_->kalmanMaxNumStubs());
212  for (unsigned iteration = 0; iteration < maxIterations; iteration++) {
213  int combinations_per_iteration = 0;
214 
215  bool easy = (l1track3D.numStubs() < settings_->kalmanMaxStubsEasy());
216  unsigned int kalmanMaxSkipLayers =
218 
219  // update each state from previous iteration (or seed) using stubs in next Kalman layer
220  vector<const KalmanState *>::const_iterator i_state = prev_states.begin();
221  for (; i_state != prev_states.end(); i_state++) {
222  const KalmanState *the_state = *i_state;
223 
224  unsigned int layer = the_state->nextLayer(); // Get KF layer where stubs to be searched for next
225  unsigned nSkipped = the_state->nSkippedLayers();
226 
227  // If this layer is known to be dead, skip to the next layer (layer+1)
228  // The next_states_skipped will then look at layer+2
229  // However, if there are stubs in this layer, then don't skip (e.g. our phi/eta boundaries might not line up exactly with a dead region)
230  // Continue to skip until you reach a functioning layer (or a layer with stubs)
231  unsigned nSkippedDeadLayers = 0;
232  unsigned nSkippedAmbiguousLayers = 0;
233  while (kfDeadLayers.find(layer) != kfDeadLayers.end() && layerStubs[layer].empty()) {
234  layer += 1;
235  ++nSkippedDeadLayers;
236  }
237  while (this->kalmanAmbiguousLayer(etaReg, layer) && layerStubs[layer].empty()) {
238  layer += 1;
239  ++nSkippedAmbiguousLayers;
240  }
241 
242  // containers for updated state+stub combinations
243  vector<const KalmanState *> next_states;
244  vector<const KalmanState *> next_states_skipped;
245 
246  // find stubs for this layer
247  // (If layer > 6, this will return empty vector, so safe).
248  vector<Stub *> thislay_stubs = layerStubs[layer];
249 
250  // find stubs for next layer if we skip a layer, except when we are on the penultimate layer,
251  // or we have exceeded the max skipped layers
252  vector<Stub *> nextlay_stubs;
253 
254  // If the next layer (layer+1) is a dead layer, then proceed to the layer after next (layer+2), if possible
255  // Also note if we need to increase "skipped" by one more for these states
256  unsigned nSkippedDeadLayers_nextStubs = 0;
257  unsigned nSkippedAmbiguousLayers_nextStubs = 0;
258  if (nSkipped < kalmanMaxSkipLayers) {
259  if (kfDeadLayers.find(layer + 1) != kfDeadLayers.end() && layerStubs[layer + 1].empty()) {
260  nextlay_stubs = layerStubs[layer + 2];
261  nSkippedDeadLayers_nextStubs++;
262  } else if (this->kalmanAmbiguousLayer(etaReg, layer) && layerStubs[layer + 1].empty()) {
263  nextlay_stubs = layerStubs[layer + 2];
264  nSkippedAmbiguousLayers_nextStubs++;
265  } else {
266  nextlay_stubs = layerStubs[layer + 1];
267  }
268  }
269 
270  // If track was not rejected by isGoodState() is previous iteration, failure here usually means the tracker ran out of layers to explore.
271  // (Due to "kalmanLay" not having unique ID for each layer within a given eta sector).
272  if (settings_->kalmanDebugLevel() >= 2 && best_state_by_nstubs.empty() && thislay_stubs.empty() &&
273  nextlay_stubs.empty())
274  PrintL1trk() << "State is lost by start of iteration " << iteration
275  << " : #thislay_stubs=" << thislay_stubs.size() << " #nextlay_stubs=" << nextlay_stubs.size()
276  << " layer=" << layer << " eta=" << l1track3D.iEtaReg();
277 
278  // If we skipped over a dead layer, only increment "nSkipped" after the stubs in next+1 layer have been obtained
279  nSkipped += nSkippedDeadLayers;
280  nSkipped += nSkippedAmbiguousLayers;
281 
282  // check to guarantee no fewer than 2PS hits per state at iteration 1
283  // (iteration 0 will always include a PS hit, but iteration 1 could use 2S hits
284  // unless we include this)
285  if (iteration == 1 && !remove2PSCut) {
286  vector<Stub *> temp_thislaystubs;
287  vector<Stub *> temp_nextlaystubs;
288  for (auto stub : thislay_stubs) {
289  if (stub->psModule())
290  temp_thislaystubs.push_back(stub);
291  }
292  for (auto stub : nextlay_stubs) {
293  if (stub->psModule())
294  temp_nextlaystubs.push_back(stub);
295  }
296  thislay_stubs = temp_thislaystubs;
297  nextlay_stubs = temp_nextlaystubs;
298  }
299 
300  combinations_per_iteration += thislay_stubs.size() + nextlay_stubs.size();
301 
302  // loop over each stub in this layer and check for compatibility with this state
303  for (unsigned i = 0; i < thislay_stubs.size(); i++) {
304  Stub *stub = thislay_stubs[i];
305 
306  // Update helix params by adding this stub.
307  const KalmanState *new_state = kalmanUpdate(nSkipped, layer, stub, the_state, tpa);
308 
309  // Cut on track chi2, pt etc.
310  if (isGoodState(*new_state))
311  next_states.push_back(new_state);
312  }
313 
314  // loop over each stub in next layer if we skip, and check for compatibility with this state
315  for (unsigned i = 0; i < nextlay_stubs.size(); i++) {
316  Stub *stub = nextlay_stubs[i];
317 
318  const KalmanState *new_state =
319  kalmanUpdate(nSkipped + 1 + nSkippedDeadLayers_nextStubs + nSkippedAmbiguousLayers_nextStubs,
320  layer + 1 + nSkippedDeadLayers_nextStubs + nSkippedAmbiguousLayers_nextStubs,
321  stub,
322  the_state,
323  tpa);
324 
325  if (isGoodState(*new_state))
326  next_states_skipped.push_back(new_state);
327  }
328 
329  // post Kalman filter local sorting per state
330  auto orderByChi2 = [](const KalmanState *a, const KalmanState *b) {
331  return bool(a->chi2scaled() < b->chi2scaled());
332  };
333  sort(next_states.begin(), next_states.end(), orderByChi2);
334  sort(next_states_skipped.begin(), next_states_skipped.end(), orderByChi2);
335 
336  new_states.insert(new_states.end(), next_states.begin(), next_states.end());
337  new_states.insert(new_states.end(), next_states_skipped.begin(), next_states_skipped.end());
338  /*
339  i = 0;
340  for (auto state : next_states) {
341  new_states.push_back(state);
342  i++;
343  }
344 
345  i = 0;
346  for (auto state : next_states_skipped) {
347  new_states.push_back(state);
348  i++;
349  }
350 */
351  } //end of state loop
352 
353  // copy new_states into prev_states for next iteration or end if we are on
354  // last iteration by clearing all states and making final state selection
355 
356  auto orderByMinSkipChi2 = [](const KalmanState *a, const KalmanState *b) {
357  return bool((a->chi2scaled()) * (a->nSkippedLayers() + 1) < (b->chi2scaled()) * (b->nSkippedLayers() + 1));
358  };
359  sort(new_states.begin(), new_states.end(), orderByMinSkipChi2); // Sort by chi2*(skippedLayers+1)
360 
361  unsigned int nStubs = iteration + 1;
362  // Success. We have at least one state that passes all cuts. Save best state found with this number of stubs.
363  if (nStubs >= settings_->kalmanMinNumStubs() && not new_states.empty())
364  best_state_by_nstubs[nStubs] = new_states[0];
365 
366  if (nStubs == settings_->kalmanMaxNumStubs()) {
367  // We're done.
368  prev_states.clear();
369  new_states.clear();
370 
371  } else {
372  // Continue iterating.
373  prev_states = new_states;
374  new_states.clear();
375  }
376  }
377 
378  if (not best_state_by_nstubs.empty()) {
379  // Select state with largest number of stubs.
380  finished_state = best_state_by_nstubs.begin()->second; // First element has largest number of stubs.
381  if (settings_->kalmanDebugLevel() >= 1) {
382  std::stringstream text;
383  text << std::fixed << std::setprecision(4);
384  text << "Track found! final state selection: nLay=" << finished_state->nStubLayers()
385  << " hitPattern=" << std::hex << finished_state->hitPattern() << std::dec
386  << " phiSec=" << l1track3D.iPhiSec() << " etaReg=" << l1track3D.iEtaReg() << " HT(m,c)=("
387  << l1track3D.cellLocationHT().first << "," << l1track3D.cellLocationHT().second << ")";
388  TVectorD y = trackParams(finished_state);
389  text << " q/pt=" << y[QOVERPT] << " tanL=" << y[T] << " z0=" << y[Z0] << " phi0=" << y[PHI0];
390  if (nHelixPar_ == 5)
391  text << " d0=" << y[D0];
392  text << " chosen from states:";
393  for (const auto &p : best_state_by_nstubs)
394  text << " " << p.second->chi2() << "/" << p.second->nStubLayers();
395  PrintL1trk() << text.str();
396  }
397  } else {
398  if (settings_->kalmanDebugLevel() >= 1) {
399  PrintL1trk() << "Track lost";
400  }
401  }
402 
403  return finished_state;
404  }
405 
406  /*--- Update a helix state by adding a stub. */
407 
409  unsigned nSkipped, unsigned int layer, Stub *stub, const KalmanState *state, const TP *tpa) {
410  if (settings_->kalmanDebugLevel() >= 4) {
411  PrintL1trk() << "---------------";
412  PrintL1trk() << "kalmanUpdate";
413  PrintL1trk() << "---------------";
414  printStub(stub);
415  }
416 
417  numUpdateCalls_++; // For monitoring, count calls to updator per track.
418 
419  // Helix params & their covariance.
420  TVectorD vecX = state->vectorX();
421  TMatrixD matC = state->matrixC();
422  if (state->barrel() && !stub->barrel()) {
423  if (settings_->kalmanDebugLevel() >= 4) {
424  PrintL1trk() << "STATE BARREL TO ENDCAP BEFORE ";
425  PrintL1trk() << "state : " << vecX[0] << " " << vecX[1] << " " << vecX[2] << " " << vecX[3];
426  PrintL1trk() << "cov(x): ";
427  matC.Print();
428  }
429  if (settings_->kalmanDebugLevel() >= 4) {
430  PrintL1trk() << "STATE BARREL TO ENDCAP AFTER ";
431  PrintL1trk() << "state : " << vecX[0] << " " << vecX[1] << " " << vecX[2] << " " << vecX[3];
432  PrintL1trk() << "cov(x): ";
433  matC.Print();
434  }
435  }
436  // Matrix to propagate helix reference point from one layer to next.
437  TMatrixD matF = matrixF(stub, state);
438  TMatrixD matFtrans(TMatrixD::kTransposed, matF);
439  if (settings_->kalmanDebugLevel() >= 4) {
440  PrintL1trk() << "matF";
441  matF.Print();
442  }
443 
444  // Multiply matrices to get helix params relative to reference point at next layer.
445  TVectorD vecXref = matF * vecX;
446  if (settings_->kalmanDebugLevel() >= 4) {
447  PrintL1trk() << "vecFref = [";
448  for (unsigned i = 0; i < nHelixPar_; i++)
449  PrintL1trk() << vecXref[i] << ", ";
450  PrintL1trk() << "]";
451  }
452 
453  // Get stub residuals.
454  TVectorD delta = residual(stub, vecXref, state->candidate().qOverPt());
455  if (settings_->kalmanDebugLevel() >= 4) {
456  PrintL1trk() << "delta = " << delta[0] << ", " << delta[1];
457  }
458 
459  // Derivative of predicted (phi,z) intercept with layer w.r.t. helix params.
460  TMatrixD matH = matrixH(stub);
461  if (settings_->kalmanDebugLevel() >= 4) {
462  PrintL1trk() << "matH";
463  matH.Print();
464  }
465 
466  if (settings_->kalmanDebugLevel() >= 4) {
467  PrintL1trk() << "previous state covariance";
468  matC.Print();
469  }
470  // Get scattering contribution to helix parameter covariance (currently zero).
471  TMatrixD matScat(nHelixPar_, nHelixPar_);
472 
473  // Get covariance on helix parameters at new reference point including scattering..
474  TMatrixD matCref = matF * matC * matFtrans + matScat;
475  if (settings_->kalmanDebugLevel() >= 4) {
476  PrintL1trk() << "matCref";
477  matCref.Print();
478  }
479  // Get hit position covariance matrix.
480  TMatrixD matV = matrixV(stub, state);
481  if (settings_->kalmanDebugLevel() >= 4) {
482  PrintL1trk() << "matV";
483  matV.Print();
484  }
485 
486  TMatrixD matRinv = matrixRinv(matH, matCref, matV);
487  if (settings_->kalmanDebugLevel() >= 4) {
488  PrintL1trk() << "matRinv";
489  matRinv.Print();
490  }
491 
492  // Calculate Kalman Gain matrix.
493  TMatrixD matK = getKalmanGainMatrix(matH, matCref, matRinv);
494  if (settings_->kalmanDebugLevel() >= 4) {
495  PrintL1trk() << "matK";
496  matK.Print();
497  }
498 
499  // Update helix state & its covariance matrix with new stub.
500  TVectorD new_vecX(nHelixPar_);
501  TMatrixD new_matC(nHelixPar_, nHelixPar_);
502  adjustState(matK, matCref, vecXref, matH, delta, new_vecX, new_matC);
503 
504  // Update track fit chi2 with new stub.
505  double new_chi2rphi = 0, new_chi2rz = 0;
506  this->adjustChi2(state, matRinv, delta, new_chi2rphi, new_chi2rz);
507 
508  if (settings_->kalmanDebugLevel() >= 4) {
509  if (nHelixPar_ == 4)
510  PrintL1trk() << "adjusted x = " << new_vecX[0] << ", " << new_vecX[1] << ", " << new_vecX[2] << ", "
511  << new_vecX[3];
512  else if (nHelixPar_ == 5)
513  PrintL1trk() << "adjusted x = " << new_vecX[0] << ", " << new_vecX[1] << ", " << new_vecX[2] << ", "
514  << new_vecX[3] << ", " << new_vecX[4];
515  PrintL1trk() << "adjusted C ";
516  new_matC.Print();
517  PrintL1trk() << "adjust chi2rphi=" << new_chi2rphi << " chi2rz=" << new_chi2rz;
518  }
519 
520  const KalmanState *new_state = mkState(
521  state->candidate(), nSkipped, layer, state, new_vecX, new_matC, matK, matV, stub, new_chi2rphi, new_chi2rz);
522 
523  return new_state;
524  }
525 
526  /* Create a KalmanState, containing a helix state & next stub it is to be updated with. */
527 
528  const KalmanState *KFbase::mkState(const L1track3D &candidate,
529  unsigned nSkipped,
530  unsigned layer,
531  const KalmanState *last_state,
532  const TVectorD &vecX,
533  const TMatrixD &matC,
534  const TMatrixD &matK,
535  const TMatrixD &matV,
536  Stub *stub,
537  double chi2rphi,
538  double chi2rz) {
539  auto new_state = std::make_unique<const KalmanState>(
540  settings_, candidate, nSkipped, layer, last_state, vecX, matC, matK, matV, stub, chi2rphi, chi2rz);
541 
542  const KalmanState *p_new_state = new_state.get();
543  listAllStates_.push_back(std::move(new_state)); // Vector keeps ownership of all states.
544  return p_new_state;
545  }
546 
547  /* Product of H*C*H(transpose) (where C = helix covariance matrix) */
548 
549  TMatrixD KFbase::matrixHCHt(const TMatrixD &matH, const TMatrixD &matC) const {
550  TMatrixD matHtrans(TMatrixD::kTransposed, matH);
551  return matH * matC * matHtrans;
552  }
553 
554  /* Get inverted Kalman R matrix: inverse(V + HCHt) */
555 
556  TMatrixD KFbase::matrixRinv(const TMatrixD &matH, const TMatrixD &matCref, const TMatrixD &matV) const {
557  TMatrixD matHCHt = matrixHCHt(matH, matCref);
558  TMatrixD matR = matV + matHCHt;
559  TMatrixD matRinv(2, 2);
560  if (matR.Determinant() > 0) {
561  matRinv = TMatrixD(TMatrixD::kInverted, matR);
562  } else {
563  // Protection against rare maths instability.
564  const TMatrixD unitMatrix(TMatrixD::kUnit, TMatrixD(nHelixPar_, nHelixPar_));
565  const double big = 9.9e9;
566  matRinv = big * unitMatrix;
567  }
568  if (settings_->kalmanDebugLevel() >= 4) {
569  PrintL1trk() << "matHCHt";
570  matHCHt.Print();
571  PrintL1trk() << "matR";
572  matR.Print();
573  }
574  return matRinv;
575  }
576 
577  /* Determine Kalman gain matrix K */
578 
579  TMatrixD KFbase::getKalmanGainMatrix(const TMatrixD &matH, const TMatrixD &matCref, const TMatrixD &matRinv) const {
580  TMatrixD matHtrans(TMatrixD::kTransposed, matH);
581  TMatrixD matCrefht = matCref * matHtrans;
582  TMatrixD matK = matCrefht * matRinv;
583  return matK;
584  }
585 
586  /* Calculate stub residual w.r.t. helix */
587 
588  TVectorD KFbase::residual(const Stub *stub, const TVectorD &vecX, double candQoverPt) const {
589  TVectorD vd = vectorM(stub); // Get (phi relative to sector, z) of hit.
590  TMatrixD h = matrixH(stub);
591  TVectorD hx = h * vecX; // Get intercept of helix with layer (linear approx).
592  TVectorD delta = vd - hx;
593 
594  // Calculate higher order corrections to residuals.
595 
596  if (not settings_->kalmanHOfw()) {
597  TVectorD correction(2);
598 
599  float inv2R = (settings_->invPtToInvR()) * 0.5 * candQoverPt;
600  float tanL = vecX[T];
601  float z0 = vecX[Z0];
602 
603  float deltaS = 0.;
604  if (settings_->kalmanHOhelixExp()) {
605  // Higher order correction correction to circle expansion for improved accuracy at low Pt.
606  double corr = stub->r() * inv2R;
607 
608  // N.B. In endcap 2S, this correction to correction[0] is exactly cancelled by the deltaS-dependent correction to it below.
609  correction[0] += (1. / 6.) * pow(corr, 3);
610 
611  deltaS = (1. / 6.) * (stub->r()) * pow(corr, 2);
612  correction[1] -= deltaS * tanL;
613  }
614 
615  if ((not stub->barrel()) && not(stub->psModule())) {
616  // These corrections rely on inside --> outside tracking, so r-z track params in 2S modules known.
617  float rShift = (stub->z() - z0) / tanL - stub->r();
618 
620  rShift -= deltaS;
621 
622  if (settings_->kalmanHOprojZcorr() == 1) {
623  // Add correlation term related to conversion of stub residuals from (r,phi) to (z,phi).
624  correction[0] += inv2R * rShift;
625  }
626 
627  if (settings_->kalmanHOalpha() == 1) {
628  // Add alpha correction for non-radial 2S endcap strips..
629  correction[0] += stub->alpha() * rShift;
630  }
631  }
632 
633  // Apply correction to residuals.
634  delta += correction;
635  }
636 
637  delta[0] = reco::deltaPhi(delta[0], 0.);
638 
639  return delta;
640  }
641 
642  /* Update helix state & its covariance matrix with new stub */
643 
644  void KFbase::adjustState(const TMatrixD &matK,
645  const TMatrixD &matCref,
646  const TVectorD &vecXref,
647  const TMatrixD &matH,
648  const TVectorD &delta,
649  TVectorD &new_vecX,
650  TMatrixD &new_matC) const {
651  new_vecX = vecXref + matK * delta;
652  const TMatrixD unitMatrix(TMatrixD::kUnit, TMatrixD(nHelixPar_, nHelixPar_));
653  TMatrixD tmp = unitMatrix - matK * matH;
654  new_matC = tmp * matCref;
655  }
656 
657  /* Update track fit chi2 with new stub */
658 
660  const TMatrixD &matRinv,
661  const TVectorD &delta,
662  double &chi2rphi,
663  double &chi2rz) const {
664  // Change in chi2 (with r-phi/r-z correlation term included in r-phi component)
665  double delChi2rphi = delta[PHI] * delta[PHI] * matRinv[PHI][PHI] + 2 * delta[PHI] * delta[Z] * matRinv[PHI][Z];
666  double delChi2rz = delta[Z] * delta[Z] * matRinv[Z][Z];
667 
668  if (settings_->kalmanDebugLevel() >= 4) {
669  PrintL1trk() << "delta(chi2rphi)=" << delChi2rphi << " delta(chi2rz)= " << delChi2rz;
670  }
671  chi2rphi = state->chi2rphi() + delChi2rphi;
672  chi2rz = state->chi2rz() + delChi2rz;
673  return;
674  }
675 
676  /* Reset internal data ready for next track. */
677 
679 
680  /* Get Kalman layer mapping (i.e. layer order in which stubs should be processed) */
681 
682  unsigned int KFbase::kalmanLayer(
683  unsigned int iEtaReg, unsigned int layerIDreduced, bool barrel, float r, float z) const {
684  // index across is GP encoded layer ID (where barrel layers=1,2,7,5,4,3 & endcap wheels=3,4,5,6,7 & 0 never occurs)
685  // index down is eta reg
686  // element is kalman layer, where 7 is invalid
687 
688  // If stub with given GP encoded layer ID can have different KF layer ID depending on whether it
689  // is barrel or endcap, then in layerMap, the the barrel case is assumed.
690  // The endcap case is fixed by hand later in this function.
691 
692  const unsigned int nEta = 16;
693  const unsigned int nGPlayID = 7;
694 
695  if (nEta != numEtaRegions_)
696  throw cms::Exception("LogicError")
697  << "ERROR KFbase::getKalmanLayer hardwired value of nEta differs from NumEtaRegions cfg param";
698 
699  // In cases where identical GP encoded layer ID present in this sector from both barrel & endcap, this array filled considering barrel. The endcap is fixed by subsequent code.
700 
701  constexpr unsigned layerMap[nEta / 2][nGPlayID + 1] = {
702  {7, 0, 1, 5, 4, 3, 7, 2}, // B1 B2 B3 B4 B5 B6 -- current FW
703  {7, 0, 1, 5, 4, 3, 7, 2}, // B1 B2 B3 B4 B5 B6
704  {7, 0, 1, 5, 4, 3, 7, 2}, // B1 B2 B3 B4 B5 B6
705  {7, 0, 1, 5, 4, 3, 7, 2}, // B1 B2 B3 B4 B5 B6
706  {7, 0, 1, 5, 4, 3, 7, 2}, // B1 B2 B3 B4(/D3) B5(/D2) B6(/D1)
707  {7, 0, 1, 3, 4, 2, 6, 2}, // B1 B2 B3(/D5)+B4(/D3) D1 D2 X D4
708  {7, 0, 1, 1, 2, 3, 4, 5}, // B1 B2+D1 D2 D3 D5 D6
709  {7, 0, 7, 1, 2, 3, 4, 5}, // B1 D1 D2 D3 D4 D5
710  };
711 
712  unsigned int kfEtaReg; // KF VHDL eta sector def: small in barrel & large in endcap.
713  if (iEtaReg < numEtaRegions_ / 2) {
714  kfEtaReg = numEtaRegions_ / 2 - 1 - iEtaReg;
715  } else {
716  kfEtaReg = iEtaReg - numEtaRegions_ / 2;
717  }
718 
719  unsigned int kalmanLay = layerMap[kfEtaReg][layerIDreduced];
720 
721  // Fixes to layermap when "maybe layer" used
722  if (settings_->kfUseMaybeLayers()) {
723  switch (kfEtaReg) {
724  case 5: //case 5: B1 B2 (B3+B4)* D1 D2 D3+D4 D5+D6 -- B3 is combined with B4 and is flagged as "maybe layer"
725  if (layerIDreduced == 6) {
726  kalmanLay = 5;
727  }
728  break;
729  case 6: //case 6: B1* B2* D1 D2 D3 D4 D5 -- B1 and B2 are flagged as "maybe layer"
730  if (layerIDreduced > 2) {
731  kalmanLay++;
732  }
733  break;
734  default:
735  break;
736  }
737  }
738 
739  // Fixes to endcap stubs, for cases where identical GP encoded layer ID present in this sector from both barrel & endcap.
740 
741  if (not barrel) {
742  switch (kfEtaReg) {
743  case 4: // B1 B2 B3 B4 B5/D1 B6/D2 D3
744  if (layerIDreduced == 3) {
745  kalmanLay = 4;
746  } else if (layerIDreduced == 4) {
747  kalmanLay = 5;
748  } else if (layerIDreduced == 5) {
749  kalmanLay = 6;
750  }
751  break;
752  //case 5: // B1 B2 B3+B4 D1 D2 D3 D4/D5
753  case 5: // B1 B2 B3 D1+B4 D2 D3 D4/D5
754  if (layerIDreduced == 5) {
755  kalmanLay = 5;
756  } else if (layerIDreduced == 7) {
757  kalmanLay = 6;
758  }
759  break;
760  default:
761  break;
762  }
763  }
764 
765  /*
766  // Fix cases where a barrel layer only partially crosses the eta sector.
767  // (Logically should work, but actually reduces efficiency -- INVESTIGATE).
768 
769  const float barrelHalfLength = 120.;
770  const float barrel4Radius = 68.8;
771  const float barrel5Radius = 86.1;
772 
773  if ( not barrel) {
774  switch ( kfEtaReg ) {
775  case 4:
776  if (layerIDreduced==3) { // D1
777  float disk1_rCut = barrel5Radius*(std::abs(z)/barrelHalfLength);
778  if (r > disk1_rCut) kalmanLay++;
779  }
780  break;
781  case 5:
782  if (layerIDreduced==3) { // D1
783  float disk1_rCut = barrel4Radius*(std::abs(z)/barrelHalfLength);
784  if (r > disk1_rCut) kalmanLay++;
785  }
786  if (layerIDreduced==4) { // D2
787  float disk2_rCut = barrel4Radius*(std::abs(z)/barrelHalfLength);
788  if (r > disk2_rCut) kalmanLay++;
789  }
790  break;
791  default:
792  break;
793  }
794  }
795  */
796 
797  return kalmanLay;
798  }
799 
800  /*=== Check if particles in given eta sector are uncertain to go through the given KF layer. */
801  /*=== (If so, count layer for numbers of hit layers, but not for number of skipped layers). */
802 
803  bool KFbase::kalmanAmbiguousLayer(unsigned int iEtaReg, unsigned int kfLayer) {
804  // Only helps in extreme forward sector, and there not significantly.
805  // UNDERSTAND IF CAN BE USED ELSEWHERE.
806 
807  const unsigned int nEta = 16;
808  const unsigned int nKFlayer = 7;
809  constexpr bool ambiguityMap[nEta / 2][nKFlayer] = {
810  {false, false, false, false, false, false, false},
811  {false, false, false, false, false, false, false},
812  {false, false, false, false, false, false, false},
813  {false, false, false, false, false, false, false},
814  {false, false, false, false, false, false, false},
815  {false, false, true, false, false, false, false},
816  {true, true, false, false, false, false, false},
817  {true, false, false, false, false, false, false},
818  };
819 
820  unsigned int kfEtaReg; // KF VHDL eta sector def: small in barrel & large in endcap.
821  if (iEtaReg < numEtaRegions_ / 2) {
822  kfEtaReg = numEtaRegions_ / 2 - 1 - iEtaReg;
823  } else {
824  kfEtaReg = iEtaReg - numEtaRegions_ / 2;
825  }
826 
827  bool ambiguous = false;
829  ambiguous = ambiguityMap[kfEtaReg][kfLayer];
830 
831  return ambiguous;
832  }
833 
834  /* Adjust KF algorithm to allow for any dead tracker layers */
835 
836  set<unsigned> KFbase::kalmanDeadLayers(bool &remove2PSCut) const {
837  // Kill scenarios described StubKiller.cc
838 
839  // By which Stress Test scenario (if any) are dead modules being emulated?
840  const StubKiller::KillOptions killScenario = static_cast<StubKiller::KillOptions>(settings_->killScenario());
841  // Should TMTT tracking be modified to reduce efficiency loss due to dead modules?
842  const bool killRecover = settings_->killRecover();
843 
844  set<pair<unsigned, bool>> deadGPlayers; // GP layer ID & boolean indicating if in barrel.
845 
846  // Range of sectors chosen to cover dead regions from StubKiller.
847  if (killRecover) {
848  if (killScenario == StubKiller::KillOptions::layer5) { // barrel layer 5
849  if (iEtaReg_ >= 3 && iEtaReg_ <= 7 && iPhiSec_ >= 1 && iPhiSec_ <= 5) {
850  deadGPlayers.insert(pair<unsigned, bool>(4, true));
851  }
852  } else if (killScenario == StubKiller::KillOptions::layer1) { // barrel layer 1
853  if (iEtaReg_ <= 7 && iPhiSec_ >= 1 && iPhiSec_ <= 5) {
854  deadGPlayers.insert(pair<unsigned, bool>(1, true));
855  }
856  remove2PSCut = true;
857  } else if (killScenario == StubKiller::KillOptions::layer1layer2) { // barrel layers 1 & 2
858  if (iEtaReg_ <= 7 && iPhiSec_ >= 1 && iPhiSec_ <= 5) {
859  deadGPlayers.insert(pair<unsigned, bool>(1, true));
860  }
861  if (iEtaReg_ >= 1 && iEtaReg_ <= 7 && iPhiSec_ >= 1 && iPhiSec_ <= 5) {
862  deadGPlayers.insert(pair<unsigned, bool>(2, true));
863  }
864  remove2PSCut = true;
865  } else if (killScenario == StubKiller::KillOptions::layer1disk1) { // barrel layer 1 & disk 1
866  if (iEtaReg_ <= 7 && iPhiSec_ >= 1 && iPhiSec_ <= 5) {
867  deadGPlayers.insert(pair<unsigned, bool>(1, true));
868  }
869  if (iEtaReg_ <= 3 && iPhiSec_ >= 1 && iPhiSec_ <= 5) {
870  deadGPlayers.insert(pair<unsigned, bool>(3, false));
871  }
872  remove2PSCut = true;
873  }
874  }
875 
876  set<unsigned> kfDeadLayers;
877  for (const auto &p : deadGPlayers) {
878  unsigned int layer = p.first;
879  bool barrel = p.second;
880  float r = 0.; // This fails for r-dependent parts of kalmanLayer(). FIX
881  float z = 999.;
882  unsigned int kalmanLay = this->kalmanLayer(iEtaReg_, layer, barrel, r, z);
883  kfDeadLayers.insert(kalmanLay);
884  }
885 
886  return kfDeadLayers;
887  }
888 
889  //=== Function to calculate approximation for tilted barrel modules (aka B) copied from Stub class.
890 
891  float KFbase::approxB(float z, float r) const {
893  }
894 
895  /* Print truth particle */
896 
897  void KFbase::printTP(const TP *tp) const {
898  TVectorD tpParams(5);
899  bool useForAlgEff(false);
900  if (tp) {
901  useForAlgEff = tp->useForAlgEff();
902  tpParams[QOVERPT] = tp->qOverPt();
903  tpParams[PHI0] = tp->phi0();
904  tpParams[Z0] = tp->z0();
905  tpParams[T] = tp->tanLambda();
906  tpParams[D0] = tp->d0();
907  }
908  std::stringstream text;
909  text << std::fixed << std::setprecision(4);
910  if (tp) {
911  text << " TP index = " << tp->index() << " useForAlgEff = " << useForAlgEff << " ";
912  const string helixNames[5] = {"qOverPt", "phi0", "z0", "tanL", "d0"};
913  for (int i = 0; i < tpParams.GetNrows(); i++) {
914  text << helixNames[i] << ":" << tpParams[i] << ", ";
915  }
916  text << " inv2R = " << tp->qOverPt() * settings_->invPtToInvR() * 0.5;
917  } else {
918  text << " Fake";
919  }
920  PrintL1trk() << text.str();
921  }
922 
923  /* Print tracker layers with stubs */
924 
925  void KFbase::printStubLayers(const vector<Stub *> &stubs, unsigned int iEtaReg) const {
926  std::stringstream text;
927  text << std::fixed << std::setprecision(4);
928  if (stubs.empty())
929  text << "stub layers = []\n";
930  else {
931  text << "stub layers = [ ";
932  for (unsigned i = 0; i < stubs.size(); i++) {
933  text << stubs[i]->layerId();
934  if (i != stubs.size() - 1)
935  text << ", ";
936  }
937  text << " ] ";
938  text << "KF stub layers = [ ";
939  for (unsigned j = 0; j < stubs.size(); j++) {
940  unsigned int kalmanLay =
941  this->kalmanLayer(iEtaReg, stubs[j]->layerIdReduced(), stubs[j]->barrel(), stubs[j]->r(), stubs[j]->z());
942  text << kalmanLay;
943  if (j != stubs.size() - 1)
944  text << ", ";
945  }
946  text << " ]\n";
947  }
948  PrintL1trk() << text.str();
949  }
950 
951  /* Print a stub */
952 
953  void KFbase::printStub(const Stub *stub) const {
954  std::stringstream text;
955  text << std::fixed << std::setprecision(4);
956  text << "stub ";
957  text << "index=" << stub->index() << " ";
958  text << "layerId=" << stub->layerId() << " ";
959  text << "r=" << stub->r() << " ";
960  text << "phi=" << stub->phi() << " ";
961  text << "z=" << stub->z() << " ";
962  text << "sigmaX=" << stub->sigmaPerp() << " ";
963  text << "sigmaZ=" << stub->sigmaPar() << " ";
964  text << "TPids=";
965  std::set<const TP *> tps = stub->assocTPs();
966  for (auto tp : tps)
967  text << tp->index() << ",";
968  PrintL1trk() << text.str();
969  }
970 
971  /* Print all stubs */
972 
973  void KFbase::printStubs(const vector<Stub *> &stubs) const {
974  for (auto &stub : stubs) {
975  printStub(stub);
976  }
977  }
978 
979 } // namespace tmtt
tmtt::KFbase::matrixF
virtual TMatrixD matrixF(const Stub *stub=nullptr, const KalmanState *state=nullptr) const =0
tmtt::Settings::kalmanHOfw
bool kalmanHOfw() const
Definition: Settings.h:337
tmtt::Settings::numEtaRegions
unsigned int numEtaRegions() const
Definition: Settings.h:125
alignBH_cfg.fixed
fixed
Definition: alignBH_cfg.py:54
tmtt::L1track3D::qOverPt
float qOverPt() const override
Definition: L1track3D.h:149
tmtt::Stub::phi
float phi() const
Definition: Stub.h:107
electrons_cff.bool
bool
Definition: electrons_cff.py:366
mps_fire.i
i
Definition: mps_fire.py:428
tmtt::Settings::kalmanHOalpha
unsigned int kalmanHOalpha() const
Definition: Settings.h:333
tmtt::KFbase::adjustChi2
virtual void adjustChi2(const KalmanState *state, const TMatrixD &covRinv, const TVectorD &delta, double &chi2rphi, double &chi2rz) const
Definition: KFbase.cc:659
tmtt::L1track3D::z0
float z0() const
Definition: L1track3D.h:159
tmtt::KFbase::isGoodState
virtual bool isGoodState(const KalmanState &state) const =0
Reference_intrackfit_cff.barrel
list barrel
Definition: Reference_intrackfit_cff.py:37
tmtt::KFbase::printStubs
void printStubs(const std::vector< Stub * > &stubs) const
Definition: KFbase.cc:973
tmtt::KFbase::doKF
const KalmanState * doKF(const L1track3D &l1track3D, const std::vector< Stub * > &stubs, const TP *tpa)
Definition: KFbase.cc:166
tmtt::KFbase::iPhiSec_
unsigned int iPhiSec_
Definition: KFbase.h:157
tmtt::Stub::psModule
bool psModule() const
Definition: Stub.h:196
HLT_FULL_cff.nEta
nEta
Definition: HLT_FULL_cff.py:6586
tmtt::Settings::kalmanDebugLevel
unsigned kalmanDebugLevel() const
Definition: Settings.h:297
tmtt::L1track3D::stubs
const std::vector< Stub * > & stubs() const override
Definition: L1track3D.h:95
reco::deltaPhi
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
tmtt::StubKiller::KillOptions::layer1
deltaPhi.h
tmtt::Stub::assocTPs
const std::set< const TP * > & assocTPs() const
Definition: Stub.h:164
tmtt::KFbase::resetStates
void resetStates()
Definition: KFbase.cc:678
tmtt::KFbase::QOVERPT
Definition: KFbase.h:34
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
tmtt::L1track3D::iEtaReg
unsigned int iEtaReg() const override
Definition: L1track3D.h:175
tmtt::StubKiller::KillOptions::layer1layer2
tmtt::KFbase::numEtaRegions_
unsigned numEtaRegions_
Definition: KFbase.h:155
tmtt::KFbase::D0
Definition: KFbase.h:33
tmtt::Settings::enableDigitize
bool enableDigitize() const
Definition: Settings.h:80
tmtt::KFbase::residual
virtual TVectorD residual(const Stub *stub, const TVectorD &x, double candQoverPt) const
Definition: KFbase.cc:588
tmtt::Stub::r
float r() const
Definition: Stub.h:108
tmtt::KFbase::numUpdateCalls_
unsigned int numUpdateCalls_
Definition: KFbase.h:160
spr::goodTrack
bool goodTrack(const reco::Track *pTrack, math::XYZPoint leadPV, trackSelectionParameters parameters, bool debug=false)
Definition: TrackSelection.cc:10
tmtt::L1track3D::cellLocationHT
std::pair< unsigned int, unsigned int > cellLocationHT() const override
Definition: L1track3D.h:101
tmtt::KFbase::PHI0
Definition: KFbase.h:33
vd
std::vector< DeviationSensor2D * > vd
Definition: DeviationsFromFileSensor2D.h:21
tmtt::TP
Definition: TP.h:23
tmtt::L1track3D::tanLambda
float tanLambda() const
Definition: L1track3D.h:160
tmtt::Settings::kalmanMaxSkipLayersEasy
unsigned int kalmanMaxSkipLayersEasy() const
Definition: Settings.h:308
tmtt::KFbase::matrixHCHt
TMatrixD matrixHCHt(const TMatrixD &h, const TMatrixD &c) const
Definition: KFbase.cc:549
createJobs.tmp
tmp
align.sh
Definition: createJobs.py:716
tmtt::Settings::kalmanMaxStubsPerLayer
unsigned int kalmanMaxStubsPerLayer() const
Definition: Settings.h:322
parallelization.uint
uint
Definition: parallelization.py:124
tmtt::KFbase::printStub
void printStub(const Stub *stub) const
Definition: KFbase.cc:953
StubKiller.h
tmtt::Settings::kfUseMaybeLayers
bool kfUseMaybeLayers() const
Definition: Settings.h:312
tmtt::KFbase::Z0
Definition: KFbase.h:33
tmtt::L1track3D::pt
float pt() const
Definition: L1track3D.h:153
tmtt::KFbase::PHI
Definition: KFbase.h:35
Utility.h
tmtt::KFbase::trackParams
virtual TVectorD trackParams(const KalmanState *state) const =0
tmtt::TrackFitGeneric::settings_
const Settings * settings_
Definition: TrackFitGeneric.h:31
tmtt::L1track3D::iPhiSec
unsigned int iPhiSec() const override
Definition: L1track3D.h:174
tmtt::KFbase::matrixH
virtual TMatrixD matrixH(const Stub *stub) const =0
tmtt::L1track3D::matchedTP
const TP * matchedTP() const override
Definition: L1track3D.h:186
h
tmtt::KFbase::fit
L1fittedTrack fit(const L1track3D &l1track3D) override
Definition: KFbase.cc:38
alignCSCRings.corr
dictionary corr
Definition: alignCSCRings.py:124
Service.h
tmtt::KalmanState::nSkippedLayers
unsigned nSkippedLayers() const
Definition: KalmanState.h:43
tmtt::KFbase::mkState
const KalmanState * mkState(const L1track3D &candidate, unsigned nSkipped, unsigned layer, const KalmanState *last_state, const TVectorD &x, const TMatrixD &pxx, const TMatrixD &K, const TMatrixD &dcov, Stub *stub, double chi2rphi, double chi2rz)
Definition: KFbase.cc:528
tmtt::L1fittedTrack::setInfoKF
void setInfoKF(unsigned int nSkippedLayers, unsigned int numUpdateCalls)
Definition: L1fittedTrack.h:124
tmtt::KFbase::iEtaReg_
unsigned int iEtaReg_
Definition: KFbase.h:158
HLTMuonOfflineAnalyzer_cfi.z0
z0
Definition: HLTMuonOfflineAnalyzer_cfi.py:98
tmtt::Stub::index
unsigned int index() const
Definition: Stub.h:101
tmtt::PrintL1trk
Definition: PrintL1trk.h:12
tmtt::KFbase::T
Definition: KFbase.h:33
b
double b
Definition: hdecay.h:118
tmtt::KFbase::kalmanUpdate
virtual const KalmanState * kalmanUpdate(unsigned nSkipped, unsigned layer, Stub *stub, const KalmanState *state, const TP *tp)
Definition: KFbase.cc:408
cmsswSequenceInfo.tp
tp
Definition: cmsswSequenceInfo.py:17
tmtt::KFbase::nHelixPar_
unsigned nHelixPar_
Definition: KFbase.h:153
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
tmtt::Stub::z
float z() const
Definition: Stub.h:109
TFileService.h
tmtt::Settings::bApprox_intercept
double bApprox_intercept() const
Definition: Settings.h:105
TP.h
tmtt::L1fittedTrack
Definition: L1fittedTrack.h:30
KalmanState.h
tmtt::Stub::layerId
unsigned int layerId() const
Definition: Stub.h:198
a
double a
Definition: hdecay.h:119
tmtt::Settings::kalmanMaxStubsEasy
unsigned int kalmanMaxStubsEasy() const
Definition: Settings.h:310
HLT_FULL_cff.maxIterations
maxIterations
Definition: HLT_FULL_cff.py:13241
tmtt::KFbase::vectorM
virtual TVectorD vectorM(const Stub *stub) const =0
tmtt::Settings::killRecover
bool killRecover() const
Definition: Settings.h:345
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
tmtt::StubKiller::KillOptions
KillOptions
Definition: StubKiller.h:22
tmtt::Settings::kalmanMaxSkipLayersHard
unsigned int kalmanMaxSkipLayersHard() const
Definition: Settings.h:307
tmtt::Settings::kalmanRemove2PScut
bool kalmanRemove2PScut() const
Definition: Settings.h:305
pfMETCorrectionType0_cfi.correction
correction
Definition: pfMETCorrectionType0_cfi.py:39
tmtt::KalmanState
Definition: KalmanState.h:22
jetUpdater_cfi.sort
sort
Definition: jetUpdater_cfi.py:29
tmtt::KFbase::approxB
float approxB(float z, float r) const
Definition: KFbase.cc:891
dumpMFGeometry_cfg.delta
delta
Definition: dumpMFGeometry_cfg.py:25
tmtt::KalmanState::settings
const Settings * settings() const
Definition: KalmanState.h:37
tmtt::Settings::kalmanHOprojZcorr
unsigned int kalmanHOprojZcorr() const
Definition: Settings.h:335
cand
Definition: decayParser.h:32
tmtt::KFbase::getKalmanGainMatrix
TMatrixD getKalmanGainMatrix(const TMatrixD &h, const TMatrixD &pxcov, const TMatrixD &covRinv) const
Definition: KFbase.cc:579
tmtt::KFbase::printStubLayers
void printStubLayers(const std::vector< Stub * > &stubs, unsigned int iEtaReg) const
Definition: KFbase.cc:925
tmtt::KFbase::Z
Definition: KFbase.h:35
tmtt::KFbase::kalmanDeadLayers
std::set< unsigned > kalmanDeadLayers(bool &remove2PSCut) const
Definition: KFbase.cc:836
tmtt::KFbase::trackParams_BeamConstr
virtual TVectorD trackParams_BeamConstr(const KalmanState *state, double &chi2rphi_bcon) const =0
tmtt::Settings::kalmanMinNumStubs
unsigned int kalmanMinNumStubs() const
Definition: Settings.h:299
tmtt::Stub::sigmaPerp
float sigmaPerp() const
Definition: Stub.h:189
alignCSCRings.r
r
Definition: alignCSCRings.py:93
tmtt::KFbase::kalmanAmbiguousLayer
virtual bool kalmanAmbiguousLayer(unsigned int iEtaReg, unsigned int kfLayer)
Definition: KFbase.cc:803
tmtt::KFbase::matrixV
virtual TMatrixD matrixV(const Stub *stub, const KalmanState *state) const =0
tmtt::KalmanState::nextLayer
unsigned nextLayer() const
Definition: KalmanState.h:39
tmtt::Stub::barrel
bool barrel() const
Definition: Stub.h:201
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
tmtt::Settings::kalmanHOhelixExp
bool kalmanHOhelixExp() const
Definition: Settings.h:331
tmtt::Stub::alpha
float alpha() const
Definition: Stub.h:119
RunInfoPI::state
state
Definition: RunInfoPayloadInspectoHelper.h:16
tmtt::KFbase::seedX
virtual TVectorD seedX(const L1track3D &l1track3D) const =0
tmtt::KalmanState::hitPattern
unsigned int hitPattern() const
Definition: KalmanState.h:66
tmtt::StubKiller::KillOptions::layer1disk1
Exception
Definition: hltDiff.cc:245
tmtt::Settings::kalmanAddBeamConstr
bool kalmanAddBeamConstr() const
Definition: Settings.h:303
tmtt::L1track3D::d0
float d0() const
Definition: L1track3D.h:157
tmtt::KFbase::tpa_
const TP * tpa_
Definition: KFbase.h:165
tmtt::Stub::sigmaPar
float sigmaPar() const
Definition: Stub.h:191
KFbase.h
tmtt::KFbase::printTP
void printTP(const TP *tp) const
Definition: KFbase.cc:897
tmtt::KFbase::kalmanLayer
virtual unsigned int kalmanLayer(unsigned int iEtaReg, unsigned int layerIDreduced, bool barrel, float r, float z) const
Definition: KFbase.cc:682
tmtt::Settings::hybrid
bool hybrid() const
Definition: Settings.h:409
PrintL1trk.h
big
Definition: big.h:8
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
tmtt::KalmanState::nStubLayers
unsigned nStubLayers() const
Definition: KalmanState.h:65
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
runonSM.text
text
Definition: runonSM.py:43
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
tmtt::Settings::kalmanMaxNumStubs
unsigned int kalmanMaxNumStubs() const
Definition: Settings.h:301
tmtt::Stub
Definition: Stub.h:43
d0
static constexpr float d0
Definition: L1EGammaCrystalsEmulatorProducer.cc:85
align_cfg.iteration
iteration
Definition: align_cfg.py:5
tmtt::L1track3D
Definition: L1track3D.h:24
tmtt::Settings::bApprox_gradient
double bApprox_gradient() const
Definition: Settings.h:104
tmtt
=== This is the base class for the linearised chi-squared track fit algorithms.
Definition: Array2D.h:16
tmtt::KFbase::seedC
virtual TMatrixD seedC(const L1track3D &l1track3D) const =0
TauDecayModes.dec
dec
Definition: TauDecayModes.py:143
tmtt::StubKiller::KillOptions::layer5
tmtt::KFbase::adjustState
void adjustState(const TMatrixD &K, const TMatrixD &pxcov, const TVectorD &x, const TMatrixD &h, const TVectorD &delta, TVectorD &new_x, TMatrixD &new_xcov) const
Definition: KFbase.cc:644
tmtt::Settings::invPtToInvR
double invPtToInvR() const
Definition: Settings.h:395
tmtt::Settings::killScenario
unsigned int killScenario() const
Definition: Settings.h:343
tmtt::L1track3D::phi0
float phi0() const override
Definition: L1track3D.h:158
tmtt::KFbase::matrixRinv
TMatrixD matrixRinv(const TMatrixD &matH, const TMatrixD &matCref, const TMatrixD &matV) const
Definition: KFbase.cc:556
tmtt::L1track3D::numStubs
unsigned int numStubs() const override
Definition: L1track3D.h:97
tmtt::KFbase::listAllStates_
std::vector< std::unique_ptr< const KalmanState > > listAllStates_
Definition: KFbase.h:163
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443