KalmanFilterAlg.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////
2 ///
3 /// \file KalmanFilterAlg.cxx
4 ///
5 /// \brief Kalman Filter.
6 ///
7 /// \author H. Greenlee
8 ///
9 ////////////////////////////////////////////////////////////////////////
10 
13 #include "boost/numeric/ublas/matrix_proxy.hpp"
14 #include "boost/numeric/ublas/vector_proxy.hpp"
15 #include "cetlib_except/exception.h"
22 
23 #include "Rtypes.h"
24 #include "TCanvas.h"
25 #include "TGaxis.h"
26 #include "TLegend.h"
27 #include "TLegendEntry.h"
28 #include "TMarker.h"
29 #include "TObject.h"
30 #include "TPad.h"
31 #include "TPaveText.h"
32 #include "TText.h"
33 #include "TVirtualPad.h"
34 
35 // Local functions.
36 
37 namespace {
38 
39  void
40  hit_position(const trkf::KHitBase& hit, double& z, double& x)
41  {
42 
43  // Map hit -> (z,x) for plotting.
44 
45  // Default result.
46 
47  z = 0.;
48  x = 0.;
49 
50  // Cast this hit to KHit<1>. Only know how to handle 1D hits.
51 
52  const trkf::KHit<1>* phit1 = dynamic_cast<const trkf::KHit<1>*>(&hit);
53  if (phit1) {
54  const std::shared_ptr<const trkf::Surface>& psurf = hit.getMeasSurface();
55 
56  // Handle SurfYZPlane.
57 
58  if (const trkf::SurfYZPlane* pyz = dynamic_cast<const trkf::SurfYZPlane*>(&*psurf)) {
59 
60  // Now finished doing casts.
61  // We have a kind of hit and measurement surface that we know how to handle.
62 
63  // Get x coordinate from hit and surface.
64 
65  x = pyz->x0() + phit1->getMeasVector()(0);
66 
67  // Get z position from surface parameters.
68  // The "z" position is actually calculated as the perpendicular distance
69  // from a corner, which is a proxy for wire number.
70 
71  double z0 = pyz->z0();
72  double y0 = pyz->y0();
73  double phi = pyz->phi();
75  double ymax = geom->DetHalfWidth();
76  if (phi > 0.)
77  z = z0 * std::cos(phi) + (ymax - y0) * std::sin(phi);
78  else
79  z = z0 * std::cos(phi) - (ymax + y0) * std::sin(phi);
80  }
81  else if (const trkf::SurfYZLine* pyz = dynamic_cast<const trkf::SurfYZLine*>(&*psurf)) {
82 
83  // Now finished doing casts.
84  // We have a kind of hit and measurement surface that we know how to handle.
85 
86  // Get x coordinate from surface.
87 
88  x = pyz->x0();
89 
90  // Get z position from surface parameters.
91  // The "z" position is actually calculated as the perpendicular distance
92  // from a corner, which is a proxy for wire number.
93 
94  double z0 = pyz->z0();
95  double y0 = pyz->y0();
96  double phi = pyz->phi();
98  double ymax = geom->DetHalfWidth();
99  if (phi > 0.)
100  z = z0 * std::cos(phi) + (ymax - y0) * std::sin(phi);
101  else
102  z = z0 * std::cos(phi) - (ymax + y0) * std::sin(phi);
103  }
104  }
105  }
106 
107  void
108  update_momentum(const trkf::KVector<2>::type& defl,
109  const trkf::KSymMatrix<2>::type& errc,
110  const trkf::KSymMatrix<2>::type& errn,
111  double mass,
112  double& invp,
113  double& var_invp)
114  // Momentum updater.
115  //
116  // Arguments: defl - Deflection (2D slope residual between two
117  // sampled points on track).
118  // errc - Error matrix of deflection residual
119  // exclusive of multiple scattering.
120  // errn - Part of deflection noise error matrix
121  // proportional to 1/(beta*momentum).
122  // mass - Mass of particle.
123  // invp - Modified 1/momentum track parameter.
124  // var_invp - Modified 1/momentum variance.
125  //
126  // Returns: True if success.
127  //
128  // The total deflection residual error matrix is
129  //
130  // err = errc + [1/(beta*momentum)] * errn.
131  //
132  // The inverse momentum and error are updated using using log-likelihood
133  //
134  // log(L) = -0.5 * log(det(err)) - 0.5 * defl^T * err^(-1) * defl.
135  // log(L) = -0.5 * tr(log(err)) - 0.5 * defl^T * err^(-1) * defl.
136  //
137  {
138  // Calculate original k = 1./(beta*p), and original variance of k.
139 
140  double invp2 = invp * invp;
141  double invp3 = invp * invp2;
142  double invp4 = invp2 * invp2;
143  double mass2 = mass * mass;
144  double k = std::sqrt(invp2 + mass2 * invp4);
145  double dkdinvp = (invp + 2. * mass2 * invp3) / k;
146  double vark = var_invp * dkdinvp * dkdinvp;
147 
148  // First, find current inverse error matrix using momentum hypothesis.
149 
150  trkf::KSymMatrix<2>::type inverr = errc + k * errn;
151  trkf::syminvert(inverr);
152 
153  // Find the first and second derivatives of the log likelihood
154  // with respact to k.
155 
156  trkf::KMatrix<2, 2>::type temp1 = prod(inverr, errn);
157  trkf::KMatrix<2, 2>::type temp2 = prod(temp1, temp1);
158 
159  // KJK-FIXME: Ideally, the explicit type specification below
160  // ('trk::KVector<T>::type') should be replaced with the 'auto'
161  // keyword. However, Boost 1.66.0a runs into a static-assert
162  // failure with the 'complexity' of vtemp2 not being 0. I am not
163  // sure what the complexity means.
164  //
165  // If 'trkf::KVector<2>::type vtemp2 = prod(temp1, vtemp1)' is
166  // used below, then the GCC 7.3 compiler thinks there may be an
167  // uninitialized variable in vtemp2. A workaround is to use the
168  // Boost interface that fills a default-constructed vector. Note
169  // that this workaround may hide an actual problem. It is the
170  // responsibility of the maintainer of this code to determine if
171  // an actual problem exists.
172 
173  trkf::KVector<2>::type vtemp1 = prod(inverr, defl);
174  trkf::KVector<2>::type vtemp2;
175  prod(temp1, vtemp1, vtemp2);
176  trkf::KVector<2>::type vtemp3 = prod(temp1, vtemp2);
177  double derivk1 = -0.5 * trkf::trace(temp1) + 0.5 * inner_prod(defl, vtemp2);
178  double derivk2 = 0.5 * trkf::trace(temp2) - inner_prod(defl, vtemp3);
179 
180  // We expect the log-likelihood to be most nearly Gaussian
181  // with respect to variable q = k^(-1/2) = std::sqrt(beta*p).
182  // Therefore, transform the original variables and log-likelihood
183  // derivatives to q-space.
184 
185  double q = 1. / std::sqrt(k);
186  double varq = vark / (4. * k * k * k);
187  double derivq1 = (-2. * k / q) * derivk1;
188  double derivq2 = 6. * k * k * derivk1 + 4. * k * k * k * derivk2;
189 
190  if (derivq2 < 0.) {
191 
192  // Estimate the measurement in q-space.
193 
194  double q1 = q - derivq1 / derivq2;
195  double varq1 = -1. / derivq2;
196 
197  // Get updated estimated q, combining original estimate and
198  // update from the current measurement.
199 
200  double newvarq = 1. / (1. / varq + 1. / varq1);
201  double newq = newvarq * (q / varq + q1 / varq1);
202  q = newq;
203  varq = newvarq;
204 
205  // Calculate updated c = 1./p and variance.
206 
207  double q2 = q * q;
208  double q4 = q2 * q2;
209  double c2 = 2. / (q4 + q2 * std::sqrt(q4 + 4. * mass2));
210  double c = std::sqrt(c2);
211  double dcdq = -2. * (c / q) * (1. + mass2 * c2) / (1. + 2. * mass2 * c2);
212  double varc = varq * dcdq * dcdq;
213 
214  // Update result.
215 
216  invp = c;
217  var_invp = varc;
218  }
219  }
220 }
221 
222 /// Constructor.
223 
225  : fTrace(false)
226  , fMaxPErr(0.)
227  , fGoodPErr(0.)
228  , fMaxIncChisq(0.)
229  , fMaxSeedIncChisq(0.)
230  , fMaxSmoothIncChisq(0.)
231  , fMaxEndChisq(0.)
232  , fMinLHits(0)
233  , fMaxLDist(0.)
234  , fMaxPredDist(0.)
235  , fMaxSeedPredDist(0.)
236  , fMaxPropDist(0.)
237  , fMinSortDist(0.)
238  , fMaxSortDist(0.)
239  , fMaxSamePlane(0)
240  , fGapDist(0.)
241  , fMaxNoiseHits(0)
242  , fMinSampleDist(0.)
243  , fFitMomRange(false)
244  , fFitMomMS(false)
245  , fGTrace(false)
246  , fGTraceWW(0)
247  , fGTraceWH(0)
248  , fPlane(-1)
249  , fInfoPad(0)
250  , fMessages(0)
251 {
252  mf::LogInfo("KalmanFilterAlg") << "KalmanFilterAlg instantiated.";
253 
254  fTrace = pset.get<bool>("Trace");
255  fMaxPErr = pset.get<double>("MaxPErr");
256  fGoodPErr = pset.get<double>("GoodPErr");
257  fMaxIncChisq = pset.get<double>("MaxIncChisq");
258  fMaxSeedIncChisq = pset.get<double>("MaxSeedIncChisq");
259  fMaxSmoothIncChisq = pset.get<double>("MaxSmoothIncChisq");
260  fMaxEndChisq = pset.get<double>("MaxEndChisq");
261  fMinLHits = pset.get<int>("MinLHits");
262  fMaxLDist = pset.get<double>("MaxLDist");
263  fMaxPredDist = pset.get<double>("MaxPredDist");
264  fMaxSeedPredDist = pset.get<double>("MaxSeedPredDist");
265  fMaxPropDist = pset.get<double>("MaxPropDist");
266  fMinSortDist = pset.get<double>("MinSortDist");
267  fMaxSortDist = pset.get<double>("MaxSortDist");
268  fMaxSamePlane = pset.get<int>("MaxSamePlane");
269  fGapDist = pset.get<double>("GapDist");
270  fMaxNoiseHits = pset.get<double>("MaxNoiseHits");
271  fMinSampleDist = pset.get<double>("MinSampleDist");
272  fFitMomRange = pset.get<bool>("FitMomRange");
273  fFitMomMS = pset.get<bool>("FitMomMS");
274  fGTrace = pset.get<bool>("GTrace");
275  if (fGTrace) {
276  fGTraceWW = pset.get<double>("GTraceWW");
277  fGTraceWH = pset.get<double>("GTraceWH");
278  fGTraceXMin = pset.get<double>("GTraceXMin");
279  fGTraceXMax = pset.get<double>("GTraceXMax");
280  fGTraceZMin = pset.get<std::vector<double>>("GTraceZMin");
281  fGTraceZMax = pset.get<std::vector<double>>("GTraceZMax");
282  }
283 }
284 
285 /// Add hits to track.
286 ///
287 /// Arguments:
288 ///
289 /// trk - Starting track.
290 /// trg - Result global track.
291 /// prop - Propagator.
292 /// dir - Direction.
293 /// hits - Candidate hits.
294 ///
295 /// Returns: True if success.
296 ///
297 /// This method makes a unidirectional Kalman fit in the specified
298 /// direction, visiting each surface of the passed KHitContainer and
299 /// updating the track. In case of multiple measurements on the same
300 /// surface, keep (at most) the one with the smallest incremental
301 /// chisquare. Any measurements that fail the incremental chisquare
302 /// cut are rejected. Resolved hits (accepted or rejected) are moved
303 /// to the unused list in KHitContainer. The container is sorted at
304 /// the start of the method, and may be resorted during the progress
305 /// of the fit.
306 ///
307 bool
309  KGTrack& trg,
310  const Propagator& prop,
312  KHitContainer& hits,
313  bool linear) const
314 {
315  // Direction must be forward or backward (unknown is not allowed).
316 
317  if (dir != Propagator::FORWARD && dir != Propagator::BACKWARD)
318  throw cet::exception("KalmanFilterAlg") << "No direction for Kalman fit.\n";
319 
320  // Set up canvas for graphical trace.
321 
322  if (fGTrace) {
323 
324  // Make a new canvas with a unique name.
325 
326  static int cnum = 0;
327  ++cnum;
328  std::ostringstream ostr;
329  ostr << "khit" << cnum;
330  fCanvases.emplace_back(
331  new TCanvas(ostr.str().c_str(), ostr.str().c_str(), fGTraceWW, fGTraceWH));
332  fPads.clear();
333  fMarkerMap.clear();
334  int nview = 3;
335  fCanvases.back()->Divide(2, nview / 2 + 1);
336 
337  // Make subpad for each view.
338 
339  for (int iview = 0; iview < nview; ++iview) {
340 
341  std::ostringstream ostr;
342  ostr << "Plane " << iview;
343 
344  fCanvases.back()->cd(iview + 1);
345  double zmin = 0.06;
346  double zmax = 0.94;
347  double xmin = 0.04;
348  double xmax = 0.95;
349  TPad* p = new TPad(ostr.str().c_str(), ostr.str().c_str(), zmin, xmin, zmax, xmax);
350  p->SetBit(kCanDelete); // Give away ownership.
351  p->Range(fGTraceZMin[iview], fGTraceXMin, fGTraceZMax[iview], fGTraceXMax);
352  p->SetFillStyle(4000); // Transparent.
353  p->Draw();
354  fPads.push_back(p);
355 
356  // Draw label.
357 
358  TText* t = new TText(zmax - 0.02, xmax - 0.03, ostr.str().c_str());
359  t->SetBit(kCanDelete); // Give away ownership.
360  t->SetTextAlign(33);
361  t->SetTextFont(42);
362  t->SetTextSize(0.04);
363  t->Draw();
364 
365  // Draw axes.
366 
367  TGaxis* pz1 =
368  new TGaxis(zmin, xmin, zmax, xmin, fGTraceZMin[iview], fGTraceZMax[iview], 510, "");
369  pz1->SetBit(kCanDelete); // Give away ownership.
370  pz1->Draw();
371 
372  TGaxis* px1 = new TGaxis(zmin, xmin, zmin, xmax, fGTraceXMin, fGTraceXMax, 510, "");
373  px1->SetBit(kCanDelete); // Give away ownership.
374  px1->Draw();
375 
376  TGaxis* pz2 =
377  new TGaxis(zmin, xmax, zmax, xmax, fGTraceZMin[iview], fGTraceZMax[iview], 510, "-");
378  pz2->SetBit(kCanDelete); // Give away ownership.
379  pz2->Draw();
380 
381  TGaxis* px2 = new TGaxis(zmax, xmin, zmax, xmax, fGTraceXMin, fGTraceXMax, 510, "L+");
382  px2->SetBit(kCanDelete); // Give away ownership.
383  px2->Draw();
384  }
385 
386  // Draw legend.
387 
388  fCanvases.back()->cd(nview + 1);
389  fInfoPad = gPad;
390  TLegend* leg = new TLegend(0.6, 0.5, 0.99, 0.99);
391  leg->SetBit(kCanDelete); // Give away ownership.
392 
393  TLegendEntry* entry = 0;
394  TMarker* m = 0;
395 
396  m = new TMarker(0., 0., 20);
397  m->SetBit(kCanDelete);
398  m->SetMarkerSize(1.2);
399  m->SetMarkerColor(kRed);
400  entry = leg->AddEntry(m, "Hits on Track", "P");
401  entry->SetBit(kCanDelete);
402 
403  m = new TMarker(0., 0., 20);
404  m->SetBit(kCanDelete);
405  m->SetMarkerSize(1.2);
406  m->SetMarkerColor(kOrange);
407  entry = leg->AddEntry(m, "Smoothed Hits on Track", "P");
408  entry->SetBit(kCanDelete);
409 
410  m = new TMarker(0., 0., 20);
411  m->SetBit(kCanDelete);
412  m->SetMarkerSize(1.2);
413  m->SetMarkerColor(kBlack);
414  entry = leg->AddEntry(m, "Available Hits", "P");
415  entry->SetBit(kCanDelete);
416 
417  m = new TMarker(0., 0., 20);
418  m->SetBit(kCanDelete);
419  m->SetMarkerSize(1.2);
420  m->SetMarkerColor(kBlue);
421  entry = leg->AddEntry(m, "Rejected Hits", "P");
422  entry->SetBit(kCanDelete);
423 
424  m = new TMarker(0., 0., 20);
425  m->SetBit(kCanDelete);
426  m->SetMarkerSize(1.2);
427  m->SetMarkerColor(kGreen);
428  entry = leg->AddEntry(m, "Removed Hits", "P");
429  entry->SetBit(kCanDelete);
430 
431  m = new TMarker(0., 0., 20);
432  m->SetBit(kCanDelete);
433  m->SetMarkerSize(1.2);
434  m->SetMarkerColor(kMagenta);
435  entry = leg->AddEntry(m, "Seed Hits", "P");
436  entry->SetBit(kCanDelete);
437 
438  m = new TMarker(0., 0., 20);
439  m->SetBit(kCanDelete);
440  m->SetMarkerSize(1.2);
441  m->SetMarkerColor(kCyan);
442  entry = leg->AddEntry(m, "Unreachable Seed Hits", "P");
443  entry->SetBit(kCanDelete);
444 
445  leg->Draw();
446 
447  // Draw message box.
448 
449  fMessages = new TPaveText(0.01, 0.01, 0.55, 0.99, "");
450  fMessages->SetBit(kCanDelete);
451  fMessages->SetTextFont(42);
452  fMessages->SetTextSize(0.04);
453  fMessages->SetTextAlign(12);
454  fMessages->Draw();
455  TText* t = fMessages->AddText("Enter buildTrack");
456  t->SetBit(kCanDelete);
457 
458  fCanvases.back()->Update();
459  }
460 
461  // Sort container using this seed track.
462 
463  hits.sort(trk, true, prop, Propagator::UNKNOWN);
464 
465  // Draw hits and populate hit->marker map.
466 
467  if (fGTrace && fCanvases.size() > 0) {
468 
469  // Loop over sorted KHitGroups.
470  // Paint sorted seed hits magenta.
471 
472  const std::list<KHitGroup>& groups = hits.getSorted();
473  for (auto const& gr : groups) {
474 
475  // Loop over hits in this group.
476 
477  const std::vector<std::shared_ptr<const KHitBase>>& phits = gr.getHits();
478  for (auto const& phit : phits) {
479  const KHitBase& hit = *phit;
480  int pl = hit.getMeasPlane();
481  if (pl >= 0 && pl < int(fPads.size())) {
482  double z = 0.;
483  double x = 0.;
484  hit_position(hit, z, x);
485  TMarker* marker = new TMarker(z, x, 20);
486  fMarkerMap[hit.getID()] = marker;
487  fPads[pl]->cd();
488  marker->SetBit(kCanDelete); // Give away ownership.
489  marker->SetMarkerSize(0.5);
490  marker->SetMarkerColor(kMagenta);
491  marker->Draw();
492  }
493  }
494  }
495 
496  // Loop over unsorted KHitGroups.
497  // Paint unsorted seed hits cyan.
498  // There should be few, if any, unsorted seed hits.
499 
500  const std::list<KHitGroup>& ugroups = hits.getUnsorted();
501  for (auto const& gr : ugroups) {
502 
503  // Loop over hits in this group.
504 
505  const std::vector<std::shared_ptr<const KHitBase>>& phits = gr.getHits();
506  for (auto const& phit : phits) {
507  const KHitBase& hit = *phit;
508  int pl = hit.getMeasPlane();
509  if (pl >= 0 && pl < int(fPads.size())) {
510  double z = 0.;
511  double x = 0.;
512  hit_position(hit, z, x);
513  TMarker* marker = new TMarker(z, x, 20);
514  fMarkerMap[hit.getID()] = marker;
515  fPads[pl]->cd();
516  marker->SetBit(kCanDelete); // Give away ownership.
517  marker->SetMarkerSize(0.5);
518  marker->SetMarkerColor(kCyan);
519  marker->Draw();
520  }
521  }
522  }
523  fCanvases.back()->Update();
524  }
525 
526  // Loop over measurements (KHitGroup) from sorted list.
527 
528  double tchisq = 0.; // Cumulative chisquare.
529  double path = 0.; // Cumulative path distance.
530  int step = 0; // Step count.
531  int nsame = 0; // Number of consecutive measurements in same plane.
532  int last_plane = -1; // Last plane.
533  bool has_pref_plane = false; // Set to true when first preferred plane hit is added to track.
534 
535  // Make a copy of the starting track, in the form of a KFitTrack,
536  // which we will update as we go.
537 
538  TrackError err;
539  trk.getSurface()->getStartingError(err);
540  KETrack tre(trk, err);
541  KFitTrack trf(tre, path, tchisq);
542 
543  // Also use the starting track as the reference track for linearized
544  // propagation, until the track is established with reasonably small
545  // errors.
546 
547  KTrack ref(trk);
548  KTrack* pref = &ref;
549 
550  mf::LogInfo log("KalmanFilterAlg");
551 
552  // Loop over measurement groups (KHitGroups).
553 
554  while (hits.getSorted().size() > 0) {
555  ++step;
556  if (fTrace) {
557  log << "Build Step " << step << "\n";
558  log << "KGTrack has " << trg.numHits() << " hits.\n";
559  log << trf;
560  }
561 
562  // Get an iterator for the next KHitGroup.
563 
565  if (dir == Propagator::FORWARD)
566  it = hits.getSorted().begin();
567  else {
568  it = hits.getSorted().end();
569  --it;
570  }
571  const KHitGroup& gr = *it;
572 
573  if (fTrace) {
574  double path_est = gr.getPath();
575  log << "Next surface: " << *(gr.getSurface()) << "\n";
576  log << " Estimated path length of hit group = " << path_est << "\n";
577  }
578 
579  // Get the next prediction surface. If this KHitGroup is on the
580  // preferred plane, use that as the prediction surface.
581  // Otherwise, use the current track surface as the prediction
582  // surface.
583 
584  std::shared_ptr<const Surface> psurf = trf.getSurface();
585  if (gr.getPlane() < 0) throw cet::exception("KalmanFilterAlg") << "negative plane?\n";
586  if (fPlane < 0 || gr.getPlane() < 0 || fPlane == gr.getPlane()) psurf = gr.getSurface();
587 
588  // Propagate track to the prediction surface.
589 
590  std::optional<double> dist = prop.noise_prop(trf, psurf, Propagator::UNKNOWN, true, pref);
591  if (dist && std::abs(*dist) > fMaxPropDist) dist = std::nullopt;
592  double ds = 0.;
593 
594  if (dist) {
595 
596  // Propagation succeeded.
597  // Update cumulative path distance and track status.
598 
599  ds = *dist;
600  path += ds;
601  trf.setPath(path);
602  if (dir == Propagator::FORWARD)
604  else {
606  }
607  if (fTrace) {
608  log << "After propagation\n";
609  log << " Incremental propagation distance = " << ds << "\n";
610  log << " Path length of prediction surface = " << path << "\n";
611  log << "KGTrack has " << trg.numHits() << " hits.\n";
612  log << trf;
613  }
614 
615  // Loop over measurements in this group.
616 
617  const std::vector<std::shared_ptr<const KHitBase>>& hits = gr.getHits();
618  double best_chisq = 0.;
619  std::shared_ptr<const KHitBase> best_hit;
620  for (std::vector<std::shared_ptr<const KHitBase>>::const_iterator ihit = hits.begin();
621  ihit != hits.end();
622  ++ihit) {
623  const KHitBase& hit = **ihit;
624 
625  // Turn this hit blue.
626 
627  if (fGTrace && fCanvases.size() > 0) {
628  auto marker_it = fMarkerMap.find(hit.getID());
629  if (marker_it != fMarkerMap.end()) {
630  TMarker* marker = marker_it->second;
631  marker->SetMarkerColor(kBlue);
632  }
633  //fCanvases.back()->Update();
634  }
635 
636  // Update predction using current track hypothesis and get
637  // incremental chisquare.
638 
639  if (hit.predict(trf, prop)) {
640  double chisq = hit.getChisq();
641  double preddist = hit.getPredDistance();
642  if (fTrace) {
643  log << "Trying Hit.\n"
644  << hit << "\nchisq = " << chisq << "\n"
645  << "prediction distance = " << preddist << "\n";
646  }
647  if ((!has_pref_plane || abs(preddist) < fMaxSeedPredDist) &&
648  (best_hit.get() == 0 || chisq < best_chisq)) {
649  best_chisq = chisq;
650  if (chisq < fMaxSeedIncChisq) best_hit = *ihit;
651  }
652  }
653  }
654  if (fTrace) {
655  log << "Best hit incremental chisquare = " << best_chisq << "\n";
656  if (best_hit.get() != 0) {
657  log << "Hit after prediction\n";
658  log << *best_hit;
659  }
660  else
661  log << "No hit passed chisquare cut.\n";
662  }
663  if (fGTrace && fCanvases.size() > 0) fCanvases.back()->Update();
664 
665  // If we found a best measurement, and if the incremental
666  // chisquare passes the cut, add it to the track and update
667  // fit information.
668 
669  bool update_ok = false;
670  if (best_hit.get() != 0) {
671  KFitTrack trf0(trf);
672  best_hit->update(trf);
673  update_ok = trf.isValid();
674  if (!update_ok) trf = trf0;
675  }
676  if (update_ok) {
677  ds += best_hit->getPredDistance();
678  tchisq += best_chisq;
679  trf.setChisq(tchisq);
680  if (dir == Propagator::FORWARD)
682  else {
684  }
685 
686  // If the pointing error got too large, and there is no
687  // reference track, quit.
688 
689  if (pref == 0 && trf.PointingError() > fMaxPErr) {
690  if (fTrace) log << "Quitting because pointing error got too large.\n";
691  break;
692  }
693 
694  // Test number of consecutive measurements in the same plane.
695 
696  if (gr.getPlane() >= 0) {
697  if (gr.getPlane() == last_plane)
698  ++nsame;
699  else {
700  nsame = 1;
701  last_plane = gr.getPlane();
702  }
703  }
704  else {
705  nsame = 0;
706  last_plane = -1;
707  }
708  if (nsame <= fMaxSamePlane) {
709 
710  // Turn best hit red.
711 
712  if (fGTrace && fCanvases.size() > 0) {
713  int pl = best_hit->getMeasPlane();
714  if (pl >= 0 && pl < int(fPads.size())) {
715  auto marker_it = fMarkerMap.find(best_hit->getID());
716  if (marker_it != fMarkerMap.end()) {
717  TMarker* marker = marker_it->second;
718  marker->SetMarkerColor(kRed);
719 
720  // Redraw marker so that it will be on top.
721 
722  fPads[pl]->cd();
723  marker->Draw();
724  }
725  }
726  fCanvases.back()->Update();
727  }
728 
729  // Make a KHitTrack and add it to the KGTrack.
730 
731  KHitTrack trh(trf, best_hit);
732  trg.addTrack(trh);
733  if (fPlane == gr.getPlane()) has_pref_plane = true;
734 
735  // Decide if we want to kill the reference track.
736 
737  if (!linear && pref != 0 && int(trg.numHits()) >= fMinLHits &&
738  (trf.PointingError() < fGoodPErr || path > fMaxLDist)) {
739  pref = 0;
740  if (fTrace) log << "Killing reference track.\n";
741  }
742 
743  if (fTrace) {
744  log << "After update\n";
745  log << "KGTrack has " << trg.numHits() << " hits.\n";
746  log << trf;
747  }
748  }
749  }
750  }
751 
752  // The current KHitGroup is now resolved.
753  // Move it to unused list.
754 
755  hits.getUnused().splice(hits.getUnused().end(), hits.getSorted(), it);
756 
757  // If the propagation distance was the wrong direction, resort the measurements.
758 
759  if (pref == 0 && dist &&
760  ((dir == Propagator::FORWARD && (ds < fMinSortDist || ds > fMaxSortDist)) ||
761  (dir == Propagator::BACKWARD && (-ds < fMinSortDist || -ds > fMaxSortDist)))) {
762  if (fTrace) log << "Resorting measurements.\n";
763  hits.sort(trf, true, prop, dir);
764  }
765  }
766 
767  // Clean track.
768 
769  cleanTrack(trg);
770 
771  // Set the fit status of the last added KHitTrack to optimal and get
772  // the final chisquare and path length.
773 
774  double fchisq = 0.;
775  path = 0.;
776  if (trg.isValid()) {
777  path = trg.endTrack().getPath() - trg.startTrack().getPath();
778  if (dir == Propagator::FORWARD) {
780  fchisq = trg.endTrack().getChisq();
781  }
782  else {
784  fchisq = trg.startTrack().getChisq();
785  }
786  }
787 
788  // Summary.
789 
790  log << "KalmanFilterAlg build track summary.\n"
791  << "Build direction = " << (dir == Propagator::FORWARD ? "FORWARD" : "BACKWARD") << "\n"
792  << "Track has " << trg.numHits() << " hits.\n"
793  << "Track length = " << path << "\n"
794  << "Track chisquare = " << fchisq << "\n";
795 
796  // Done. Return success if we added at least one measurement.
797 
798  bool ok = trg.numHits() > 0;
799 
800  if (fGTrace && fCanvases.size() > 0) {
801  TText* t = 0;
802  if (ok)
803  t = fMessages->AddText("Exit buildTrack, status success");
804  else
805  t = fMessages->AddText("Exit buildTrack, status fail");
806  t->SetBit(kCanDelete);
807  fInfoPad->Modified();
808  fCanvases.back()->Update();
809  }
810 
811  return ok;
812 }
813 
814 /// Smooth track.
815 ///
816 /// Arguments:
817 ///
818 /// trg - Track to be smoothed.
819 /// trg1 - Track to receive result of unidirectional fit.
820 /// prop - Propagator.
821 ///
822 /// Returns: True if success.
823 ///
824 /// The starting track should be a global track that has been fit in
825 /// one direction. Fit status should be optimal at (at least) one
826 /// end. It is an error if the fit status is not optimal at either
827 /// end. If the fit status is optimal at both ends, do nothing, but
828 /// return success.
829 ///
830 /// If the second argument is non-null, save the result of the
831 /// unidirectional track fit produced as a byproduct of the smoothing
832 /// operation. This track can be smoothed in order to iterate the
833 /// Kalman fit, etc.
834 ///
835 /// The Kalman smoothing algorithm starts at the optimal end and fits
836 /// the track in the reverse direction, calculating optimal track
837 /// parameters at each measurement surface.
838 ///
839 /// All measurements are included in the reverse fit. No incremental
840 /// chisquare cut is applied.
841 ///
842 /// If any measurement surface can not be reached because of a
843 /// measurement error, the entire smoothing operation is considered a
844 /// failure. In that case, false is returned and the track is left in
845 /// an undefined state.
846 ///
847 bool
849 {
850  if (not trg.isValid()) {
851  // It is an error if the KGTrack is not valid.
852  return false;
853  }
854 
855  bool result = false;
856 
857  // Examine the track endpoints and figure out which end of the track
858  // to start from. The fit always starts at the optimal end. It is
859  // an error if neither end point is optimal. Do nothing and return
860  // success if both end points are optimal.
861 
862  const KHitTrack& trh0 = trg.startTrack();
863  const KHitTrack& trh1 = trg.endTrack();
864  KFitTrack::FitStatus stat0 = trh0.getStat();
865  KFitTrack::FitStatus stat1 = trh1.getStat();
866  bool dofit = false;
867 
868  // Remember starting direction, track, and distance.
869 
871  const KTrack* trk = 0;
872  double path = 0.;
873 
874  if (stat0 == KFitTrack::OPTIMAL) {
875  if (stat1 == KFitTrack::OPTIMAL) {
876 
877  // Both ends optimal (do nothing, return success).
878 
879  dofit = false;
880  result = true;
881  }
882  else {
883 
884  // Start optimal.
885 
886  dofit = true;
887  dir = Propagator::FORWARD;
888  trk = &trh0;
889  path = 0.;
890  }
891  }
892  else {
893  if (stat1 == KFitTrack::OPTIMAL) {
894 
895  // End optimal.
896 
897  dofit = true;
898  dir = Propagator::BACKWARD;
899  trk = &trh1;
900  path = trh1.getPath();
901  }
902  else {
903 
904  // Neither end optimal (do nothing, return failure).
905 
906  dofit = false;
907  result = false;
908  }
909  }
910  if (dofit) {
911  if (!trk)
912  throw cet::exception("KalmanFilterAlg") << "KalmanFilterAlg::smoothTrack(): no track!\n";
913 
914  // Cumulative chisquare.
915 
916  double tchisq = 0.;
917 
918  // Construct starting KFitTrack with track information and distance
919  // taken from the optimal end, but with "infinite" errors.
920 
921  TrackError err;
922  trk->getSurface()->getStartingError(err);
923  KETrack tre(*trk, err);
924  KFitTrack trf(tre, path, tchisq);
925 
926  // Make initial reference track to be same as initial fit track.
927 
928  KTrack ref(trf);
929 
930  // Loop over KHitTracks contained in KGTrack.
931 
934  switch (dir) {
935  case Propagator::FORWARD:
936  it = trg.getTrackMap().begin();
937  itend = trg.getTrackMap().end();
938  break;
940  it = trg.getTrackMap().end();
941  itend = trg.getTrackMap().begin();
942  break;
943  default:
944  throw cet::exception("KalmanFilterAlg")
945  << "KalmanFilterAlg::smoothTrack(): invalid direction\n";
946  } // switch
947 
948  mf::LogInfo log("KalmanFilterAlg");
949 
950  // Loop starts here.
951 
952  result = true; // Result success unless we find an error.
953  int step = 0; // Step count.
954  while (dofit && it != itend) {
955  ++step;
956  if (fTrace) {
957  log << "Smooth Step " << step << "\n";
958  log << "Reverse fit track:\n";
959  log << trf;
960  }
961 
962  // For backward fit, decrement iterator at start of loop.
963 
964  if (dir == Propagator::BACKWARD) --it;
965 
966  KHitTrack& trh = (*it).second;
967  if (fTrace) {
968  log << "Forward track:\n";
969  log << trh;
970  }
971 
972  // Extract measurement.
973 
974  const KHitBase& hit = *(trh.getHit());
975 
976  // Propagate KFitTrack to the next track surface.
977 
978  std::shared_ptr<const Surface> psurf = trh.getSurface();
979  std::optional<double> dist = prop.noise_prop(trf, psurf, Propagator::UNKNOWN, true, &ref);
980 
981  // Check if propagation succeeded. If propagation fails, this
982  // measurement will be dropped from the unidirectional fit
983  // track. This measurement will still be in the original
984  // track, but with a status other than optimal.
985 
986  if (dist) {
987 
988  // Propagation succeeded.
989  // Update cumulative path distance and track status.
990 
991  double ds = *dist;
992  path += ds;
993  trf.setPath(path);
994  if (dir == Propagator::FORWARD)
996  else {
998  }
999  if (fTrace) {
1000  log << "Reverse fit track after propagation:\n";
1001  log << " Propagation distance = " << ds << "\n";
1002  log << trf;
1003  }
1004 
1005  // See if we have the proper information to calculate an optimal track
1006  // at this surface (should normally be possible).
1007 
1009  KFitTrack::FitStatus newstat = trf.getStat();
1010 
1011  if ((newstat == KFitTrack::FORWARD_PREDICTED && stat == KFitTrack::BACKWARD) ||
1012  (newstat == KFitTrack::BACKWARD_PREDICTED && stat == KFitTrack::FORWARD)) {
1013 
1014  // Update stored KHitTrack to be optimal.
1015 
1016  bool ok = trh.combineFit(trf);
1017 
1018  // Update the stored path distance to be from the currently fitting track.
1019 
1020  trh.setPath(trf.getPath());
1021 
1022  // Update reference track.
1023 
1024  ref = trh;
1025 
1026  // If combination failed, abandon the fit and return failure.
1027 
1028  if (!ok) {
1029  dofit = false;
1030  result = false;
1031  break;
1032  }
1033  if (fTrace) {
1034  log << "Combined track:\n";
1035  log << trh;
1036  }
1037  }
1038 
1039  // Update measurement predction using current track hypothesis.
1040 
1041  if (not hit.predict(trf, prop, &ref)) {
1042 
1043  // Abandon the fit and return failure.
1044 
1045  dofit = false;
1046  result = false;
1047  break;
1048  }
1049 
1050  // Prediction succeeded. Get incremental chisquare. If this
1051  // hit fails incremental chisquare cut, this hit will be
1052  // dropped from the unidirecitonal Kalman fit track, but may
1053  // still be in the smoothed track.
1054 
1055  double chisq = hit.getChisq();
1056  if (chisq < fMaxSmoothIncChisq) {
1057 
1058  // Update the reverse fitting track using the current measurement
1059  // (both track parameters and status).
1060 
1061  KFitTrack trf0(trf);
1062  hit.update(trf);
1063  bool update_ok = trf.isValid();
1064  if (!update_ok)
1065  trf = trf0;
1066  else {
1067  tchisq += chisq;
1068  trf.setChisq(tchisq);
1069 
1070  if (dir == Propagator::FORWARD)
1072  else {
1074  }
1075  if (fTrace) {
1076  log << "Reverse fit track after update:\n";
1077  log << trf;
1078  }
1079  if (fGTrace && fCanvases.size() > 0) {
1080  auto marker_it = fMarkerMap.find(hit.getID());
1081  if (marker_it != fMarkerMap.end()) {
1082  TMarker* marker = marker_it->second;
1083  marker->SetMarkerColor(kOrange);
1084  }
1085  fCanvases.back()->Update();
1086  }
1087 
1088  // If unidirectional track pointer is not null, make a
1089  // KHitTrack and save it in the unidirectional track.
1090 
1091  if (trg1 != 0) {
1092  KHitTrack trh1(trf, trh.getHit());
1093  trg1->addTrack(trh1);
1094  }
1095  }
1096  }
1097  }
1098 
1099  // For forward fit, increment iterator at end of loop.
1100 
1101  if (dir == Propagator::FORWARD) ++it;
1102 
1103  } // Loop over KHitTracks.
1104 
1105  // If fit was successful and the unidirectional track pointer
1106  // is not null and the track is valid, set the fit status of
1107  // the last added KHitTrack to optimal.
1108 
1109  if (result && trg1 != 0 && trg1->isValid()) {
1110  if (dir == Propagator::FORWARD)
1112  else {
1114  }
1115  }
1116 
1117  // Recalibrate track map.
1118 
1119  trg.recalibrate();
1120 
1121  } // Do fit.
1122 
1123  // Get the final chisquare.
1124 
1125  double fchisq = 0.5 * (trg.startTrack().getChisq() + trg.endTrack().getChisq());
1126 
1127  // Summary.
1128 
1129  mf::LogInfo log("KalmanFilterAlg");
1130  log << "KalmanFilterAlg smooth track summary.\n"
1131  << "Smooth direction = " << (dir == Propagator::FORWARD ? "FORWARD" : "BACKWARD") << "\n"
1132  << "Track has " << trg.numHits() << " hits.\n"
1133  << "Track length = " << trg.endTrack().getPath() - trg.startTrack().getPath() << "\n"
1134  << "Track chisquare = " << fchisq << "\n";
1135 
1136  return result;
1137 }
1138 
1139 /// Add hits to existing track.
1140 ///
1141 /// Arguments:
1142 ///
1143 /// trg - Track to extend.
1144 /// prop - Propagator.
1145 /// hits - Hit collection to choose hits from.
1146 ///
1147 /// This method extends a KGTrack by adding hits. The KGTrack must
1148 /// have previously been produced by a unidirectional Kalman fit (it
1149 /// should be optimal at one end). This method finds the optimal end
1150 /// and extends the track in that direction. If any hits are added,
1151 /// the originally optimal end has its status reset to forward or
1152 /// backward, and the new endpoint is optimal. In any case, the final
1153 /// result is unidirectionally fit KGTrack.
1154 ///
1155 bool
1157 {
1158  // Default result failure.
1159 
1160  bool result = false;
1161 
1162  if (fGTrace && fCanvases.size() > 0) {
1163  TText* t = fMessages->AddText("Enter extendTrack");
1164  t->SetBit(kCanDelete);
1165  fInfoPad->Modified();
1166  fCanvases.back()->Update();
1167  }
1168 
1169  // Remember the original number of measurement.
1170 
1171  unsigned int nhits0 = trg.numHits();
1172 
1173  // It is an error if the KGTrack is not valid.
1174 
1175  if (trg.isValid()) {
1176  mf::LogInfo log("KalmanFilterAlg");
1177 
1178  // Examine the track endpoints and figure out which end of the
1179  // track to extend. The track is always extended from the optimal
1180  // end. It is an error if neither end point is optimal, or both
1181  // endoints are optimal. Reset the status of the optimal, and
1182  // make a copy of the starting track fit. Also get starting path
1183  // and chisquare.
1184 
1185  KHitTrack& trh0 = trg.startTrack();
1186  KHitTrack& trh1 = trg.endTrack();
1187  KFitTrack::FitStatus stat0 = trh0.getStat();
1188  KFitTrack::FitStatus stat1 = trh1.getStat();
1189  bool dofit = false;
1191  KFitTrack trf;
1192  double path = 0.;
1193  double tchisq = 0.;
1194 
1195  if (stat0 == KFitTrack::OPTIMAL) {
1196  if (stat1 == KFitTrack::OPTIMAL) {
1197 
1198  // Both ends optimal (do nothing, return failure).
1199 
1200  dofit = false;
1201  result = false;
1202  return result;
1203  }
1204  else {
1205 
1206  // Start optimal. Extend backward.
1207 
1208  dofit = true;
1209  dir = Propagator::BACKWARD;
1211  trf = trh0;
1212  path = trh0.getPath();
1213  tchisq = trh0.getChisq();
1214  }
1215  }
1216  else {
1217  if (stat1 == KFitTrack::OPTIMAL) {
1218 
1219  // End optimal. Extend forward.
1220 
1221  dofit = true;
1222  dir = Propagator::FORWARD;
1224  trf = trh1;
1225  path = trh1.getPath();
1226  tchisq = trh1.getChisq();
1227 
1228  // Make sure forward extend track momentum is over some
1229  // minimum value.
1230 
1231  if (trf.getVector()(4) > 5.) {
1232  trf.getVector()(4) = 5.;
1233  trf.getError()(4, 4) = 5.;
1234  }
1235  }
1236  else {
1237 
1238  // Neither end optimal (do nothing, return failure).
1239 
1240  dofit = false;
1241  result = false;
1242  return result;
1243  }
1244  }
1245  if (dofit) {
1246 
1247  // Sort hit container using starting track.
1248 
1249  hits.sort(trf, true, prop, dir);
1250 
1251  // Draw and add hits in hit->marker map that are not already there.
1252 
1253  if (fGTrace && fCanvases.size() > 0) {
1254 
1255  // Loop over sorted KHitGroups.
1256  // Paint sorted hits black.
1257 
1258  const std::list<KHitGroup>& groups = hits.getSorted();
1259  for (auto const& gr : groups) {
1260 
1261  // Loop over hits in this group.
1262 
1263  const std::vector<std::shared_ptr<const KHitBase>>& phits = gr.getHits();
1264  for (auto const& phit : phits) {
1265  const KHitBase& hit = *phit;
1266  int pl = hit.getMeasPlane();
1267  if (pl >= 0 && pl < int(fPads.size())) {
1268 
1269  // Is this hit already in map?
1270 
1271  TMarker* marker = 0;
1272  auto marker_it = fMarkerMap.find(hit.getID());
1273  if (marker_it == fMarkerMap.end()) {
1274  double z = 0.;
1275  double x = 0.;
1276  hit_position(hit, z, x);
1277  marker = new TMarker(z, x, 20);
1278  fMarkerMap[hit.getID()] = marker;
1279  fPads[pl]->cd();
1280  marker->SetBit(kCanDelete); // Give away ownership.
1281  marker->SetMarkerSize(0.5);
1282  marker->Draw();
1283  }
1284  else
1285  marker = marker_it->second;
1286  marker->SetMarkerColor(kBlack);
1287  }
1288  }
1289  }
1290 
1291  // Loop over unsorted KHitGroups.
1292  // Paint unsorted hits blue.
1293 
1294  const std::list<KHitGroup>& ugroups = hits.getUnsorted();
1295  for (auto const& gr : ugroups) {
1296 
1297  // Loop over hits in this group.
1298 
1299  const std::vector<std::shared_ptr<const KHitBase>>& phits = gr.getHits();
1300  for (auto const& phit : phits) {
1301  const KHitBase& hit = *phit;
1302  int pl = hit.getMeasPlane();
1303  if (pl >= 0 && pl < int(fPads.size())) {
1304 
1305  // Is this hit already in map?
1306 
1307  TMarker* marker = 0;
1308  auto marker_it = fMarkerMap.find(hit.getID());
1309  if (marker_it == fMarkerMap.end()) {
1310  double z = 0.;
1311  double x = 0.;
1312  hit_position(hit, z, x);
1313  marker = new TMarker(z, x, 20);
1314  fMarkerMap[hit.getID()] = marker;
1315  fPads[pl]->cd();
1316  marker->SetBit(kCanDelete); // Give away ownership.
1317  marker->SetMarkerSize(0.5);
1318  marker->Draw();
1319  }
1320  else
1321  marker = marker_it->second;
1322  marker->SetMarkerColor(kBlue);
1323  }
1324  }
1325  }
1326  fCanvases.back()->Update();
1327  }
1328 
1329  //std::cout << "extendTrack: marker map has " << fMarkerMap.size() << " entries." << std::endl;
1330 
1331  // Extend loop starts here.
1332 
1333  int step = 0;
1334  int nsame = 0;
1335  int last_plane = -1;
1336  while (hits.getSorted().size() > 0) {
1337  ++step;
1338  if (fTrace) {
1339  log << "Extend Step " << step << "\n";
1340  log << "KGTrack has " << trg.numHits() << " hits.\n";
1341  log << trf;
1342  }
1343 
1344  // Get an iterator for the next KHitGroup.
1345 
1347  switch (dir) {
1348  case Propagator::FORWARD: it = hits.getSorted().begin(); break;
1349  case Propagator::BACKWARD:
1350  it = hits.getSorted().end();
1351  --it;
1352  break;
1353  default:
1354  throw cet::exception("KalmanFilterAlg")
1355  << "KalmanFilterAlg::extendTrack(): invalid direction\n";
1356  } // switch
1357  const KHitGroup& gr = *it;
1358 
1359  if (fTrace) {
1360  double path_est = gr.getPath();
1361  log << "Next surface: " << *(gr.getSurface()) << "\n";
1362  log << " Estimated path length of hit group = " << path_est << "\n";
1363  }
1364 
1365  // Get the next prediction surface. If this KHitGroup is on the
1366  // preferred plane, use that as the prediction surface.
1367  // Otherwise, use the current track surface as the prediction
1368  // surface.
1369 
1370  std::shared_ptr<const Surface> psurf = trf.getSurface();
1371  if (gr.getPlane() < 0)
1372  throw cet::exception("KalmanFilterAlg")
1373  << "KalmanFilterAlg::extendTrack(): negative plane?\n";
1374  if (fPlane < 0 || gr.getPlane() < 0 || fPlane == gr.getPlane()) psurf = gr.getSurface();
1375 
1376  // Propagate track to the prediction surface.
1377 
1378  std::optional<double> dist = prop.noise_prop(trf, psurf, Propagator::UNKNOWN, true);
1379  if (dist && std::abs(*dist) > fMaxPropDist) dist = std::nullopt;
1380  double ds = 0.;
1381 
1382  if (dist) {
1383 
1384  // Propagation succeeded.
1385  // Update cumulative path distance and track status.
1386 
1387  ds = *dist;
1388  path += ds;
1389  trf.setPath(path);
1390  if (dir == Propagator::FORWARD)
1392  else {
1394  }
1395  if (fTrace) {
1396  log << "After propagation\n";
1397  log << " Incremental distance = " << ds << "\n";
1398  log << " Track path length = " << path << "\n";
1399  log << "KGTrack has " << trg.numHits() << " hits.\n";
1400  log << trf;
1401  }
1402 
1403  // Loop over measurements in this group.
1404 
1405  const std::vector<std::shared_ptr<const KHitBase>>& hits = gr.getHits();
1406  double best_chisq = 0.;
1407  std::shared_ptr<const KHitBase> best_hit;
1408  for (std::vector<std::shared_ptr<const KHitBase>>::const_iterator ihit = hits.begin();
1409  ihit != hits.end();
1410  ++ihit) {
1411  const KHitBase& hit = **ihit;
1412 
1413  // Turn this hit blue.
1414 
1415  if (fGTrace && fCanvases.size() > 0) {
1416  auto marker_it = fMarkerMap.find(hit.getID());
1417  if (marker_it != fMarkerMap.end()) {
1418  TMarker* marker = marker_it->second;
1419  marker->SetMarkerColor(kBlue);
1420  }
1421  //fCanvases.back()->Update();
1422  }
1423 
1424  // Update predction using current track hypothesis and get
1425  // incremental chisquare.
1426 
1427  bool ok = hit.predict(trf, prop);
1428  if (ok) {
1429  double chisq = hit.getChisq();
1430  double preddist = hit.getPredDistance();
1431  if (abs(preddist) < fMaxPredDist && (best_hit.get() == 0 || chisq < best_chisq)) {
1432  best_chisq = chisq;
1433  if (chisq < fMaxIncChisq) best_hit = *ihit;
1434  }
1435  }
1436  }
1437  if (fTrace) {
1438  log << "Best hit incremental chisquare = " << best_chisq << "\n";
1439  if (best_hit.get() != 0) {
1440  log << "Hit after prediction\n";
1441  log << *best_hit;
1442  }
1443  else
1444  log << "No hit passed chisquare cut.\n";
1445  }
1446  if (fGTrace && fCanvases.size() > 0) fCanvases.back()->Update();
1447 
1448  // If we found a best measurement, and if the incremental
1449  // chisquare passes the cut, add it to the track and update
1450  // fit information.
1451 
1452  bool update_ok = false;
1453  if (best_hit.get() != 0) {
1454  KFitTrack trf0(trf);
1455  best_hit->update(trf);
1456  update_ok = trf.isValid();
1457  if (!update_ok) trf = trf0;
1458  }
1459  if (update_ok) {
1460  ds += best_hit->getPredDistance();
1461  tchisq += best_chisq;
1462  trf.setChisq(tchisq);
1463  if (dir == Propagator::FORWARD)
1465  else {
1467  }
1468 
1469  // If the pointing error got too large, quit.
1470 
1471  if (trf.PointingError() > fMaxPErr) {
1472  if (fTrace) log << "Quitting because pointing error got too large.\n";
1473  break;
1474  }
1475 
1476  // Test number of consecutive measurements in the same plane.
1477 
1478  if (gr.getPlane() >= 0) {
1479  if (gr.getPlane() == last_plane)
1480  ++nsame;
1481  else {
1482  nsame = 1;
1483  last_plane = gr.getPlane();
1484  }
1485  }
1486  else {
1487  nsame = 0;
1488  last_plane = -1;
1489  }
1490  if (nsame <= fMaxSamePlane) {
1491 
1492  // Turn best hit red.
1493 
1494  if (fGTrace && fCanvases.size() > 0) {
1495  int pl = best_hit->getMeasPlane();
1496  if (pl >= 0 && pl < int(fPads.size())) {
1497  auto marker_it = fMarkerMap.find(best_hit->getID());
1498  if (marker_it != fMarkerMap.end()) {
1499  TMarker* marker = marker_it->second;
1500  marker->SetMarkerColor(kRed);
1501 
1502  // Redraw marker so that it will be on top.
1503 
1504  fPads[pl]->cd();
1505  marker->Draw();
1506  }
1507  }
1508  fCanvases.back()->Update();
1509  }
1510 
1511  // Make a KHitTrack and add it to the KGTrack.
1512 
1513  KHitTrack trh(trf, best_hit);
1514  trg.addTrack(trh);
1515 
1516  if (fTrace) {
1517  log << "After update\n";
1518  log << "KGTrack has " << trg.numHits() << " hits.\n";
1519  log << trf;
1520  }
1521  }
1522  }
1523  }
1524 
1525  // The current KHitGroup is now resolved.
1526  // Move it to unused list.
1527 
1528  hits.getUnused().splice(hits.getUnused().end(), hits.getSorted(), it);
1529 
1530  // If the propagation distance was the wrong direction, resort the measurements.
1531 
1532  if (dist && ((dir == Propagator::FORWARD && (ds < fMinSortDist || ds > fMaxSortDist)) ||
1533  (dir == Propagator::BACKWARD && (-ds < fMinSortDist || -ds > fMaxSortDist)))) {
1534  if (fTrace) log << "Resorting measurements.\n";
1535  hits.sort(trf, true, prop, dir);
1536  }
1537  }
1538  }
1539 
1540  // Clean track.
1541 
1542  cleanTrack(trg);
1543 
1544  // Set the fit status of the last added KHitTrack to optimal and
1545  // get the final chisquare and path length.
1546 
1547  double fchisq = 0.;
1548  path = 0.;
1549  if (trg.isValid()) {
1550  path = trg.endTrack().getPath() - trg.startTrack().getPath();
1551  switch (dir) {
1552  case Propagator::FORWARD:
1554  fchisq = trg.endTrack().getChisq();
1555  break;
1556  case Propagator::BACKWARD:
1558  fchisq = trg.startTrack().getChisq();
1559  break;
1560  default:
1561  throw cet::exception("KalmanFilterAlg")
1562  << "KalmanFilterAlg::extendTrack(): invalid direction [II]\n";
1563  } // switch
1564  }
1565 
1566  // Summary.
1567 
1568  log << "KalmanFilterAlg extend track summary.\n"
1569  << "Extend direction = " << (dir == Propagator::FORWARD ? "FORWARD" : "BACKWARD") << "\n"
1570  << "Track has " << trg.numHits() << " hits.\n"
1571  << "Track length = " << path << "\n"
1572  << "Track chisquare = " << fchisq << "\n";
1573  // log << trg << "\n";
1574  }
1575 
1576  // Done.
1577 
1578  result = (trg.numHits() > nhits0);
1579 
1580  if (fGTrace && fCanvases.size() > 0) {
1581  TText* t = 0;
1582  if (result)
1583  t = fMessages->AddText("Exit extendTrack, status success");
1584  else
1585  t = fMessages->AddText("Exit extendTrack, status fail");
1586  t->SetBit(kCanDelete);
1587  fInfoPad->Modified();
1588  fCanvases.back()->Update();
1589  }
1590 
1591  return result;
1592 }
1593 
1594 /// Estimate track momentum using range.
1595 ///
1596 /// Arguments:
1597 ///
1598 /// trg - Global track.
1599 /// prop - Propagator.
1600 /// tremom - Track with momentum estimate.
1601 ///
1602 /// Returns: True if success.
1603 ///
1604 /// This method generates a momentum-estimating track by extracting
1605 /// the last track from a global track, and setting its momentum to
1606 /// some small value.
1607 ///
1608 bool
1610 {
1611  if (!trg.isValid()) return false;
1612 
1613  // Extract track with lowest momentum.
1614 
1615  const KHitTrack& trh = trg.endTrack();
1616  if (trh.getStat() == KFitTrack::INVALID)
1617  throw cet::exception("KalmanFilterAlg")
1618  << "KalmanFilterAlg::fitMomentumRange(): invalid end track\n";
1619  tremom = trh;
1620 
1621  // Set track momentum to a small value.
1622 
1623  tremom.getVector()(4) = 100.;
1624  tremom.getError()(4, 0) = 0.;
1625  tremom.getError()(4, 1) = 0.;
1626  tremom.getError()(4, 2) = 0.;
1627  tremom.getError()(4, 3) = 0.;
1628  tremom.getError()(4, 4) = 10000.;
1629 
1630  // Done.
1631 
1632  return true;
1633 }
1634 
1635 /// Estimate track momentum using multiple scattering.
1636 ///
1637 /// Arguments:
1638 ///
1639 /// trg - Global track.
1640 /// prop - Propagator.
1641 /// tremom - Track containing momentum estimate.
1642 ///
1643 /// Returns: True if success.
1644 ///
1645 /// This method estimates the momentum of the specified track using
1646 /// multiple scattering. As a result of calling this method, the
1647 /// original global track is not updated, but a KETrack is produced
1648 /// near the starting surface that has the estimated momentum.
1649 ///
1650 /// The global track passed as argument should have been smoothed
1651 /// prior to calling this method so that all or most fits are optimal.
1652 /// If either the first or last fit is not optimal, return false.
1653 ///
1654 /// This method assumes that track parameter four is 1/p. This sort
1655 /// of momentum estimation only makes sense if the momentum track
1656 /// parameter is uncorrelated with any other track parameter. The
1657 /// error matrix of the first and last fit is checked for this. If it
1658 /// is found that either error matrix has correlated track parameter
1659 /// four with any other track parameter, this method returns without
1660 /// doing anything (return false).
1661 ///
1662 bool
1664  const Propagator& prop,
1665  KETrack& tremom) const
1666 {
1667  // Get iterators pointing to the first and last tracks.
1668 
1669  const std::multimap<double, KHitTrack>& trackmap = trg.getTrackMap();
1670  if (trackmap.size() < 2) return false;
1672  itend[0] = trackmap.begin();
1673  itend[1] = trackmap.end();
1674  --itend[1];
1675 
1676  // Check the fit status and error matrix of the first and last
1677  // track.
1678 
1679  bool result = true;
1680 
1681  for (int i = 0; result && i < 2; ++i) {
1682  const KHitTrack& trh = itend[i]->second;
1684  if (stat != KFitTrack::OPTIMAL) result = false;
1685  const TrackError& err = trh.getError();
1686  for (int j = 0; j < 4; ++j) {
1687  if (err(4, j) != 0.) result = false;
1688  }
1689  }
1690  if (!result) return result;
1691 
1692  // We will periodically sample the track trajectory. At each sample
1693  // point, collect the following information.
1694  //
1695  // 1. The path distance at the sample point.
1696  // 2. The original momentum of the track at the sample point and its error.
1697  // 3. One copy of the track (KETrack) that will be propagated without noise
1698  // (infinite momentum track).
1699  // 3. A second copy of the track (KETrack) that will be propagated with the
1700  // minimum allowed momentum (range out track).
1701  // 4. A third copy of the track (KETrack) that will propagated with noise
1702  // with some intermediate momentum (noise track).
1703  //
1704  // Collect the first sample from the maximum path distance track.
1705 
1706  double s_sample = itend[1]->first;
1707  const KETrack& tre = itend[1]->second;
1708  KETrack tre_inf(tre);
1709  KTrack trk_range(tre);
1710  KETrack tre_noise(tre);
1711  tre_inf.getVector()(4) = 0.;
1712  tre_inf.getError()(4, 4) = 0.;
1713  trk_range.getVector()(4) = 100.;
1714  tre_noise.getError()(4, 4) = 0.;
1715  tre_noise.getVector()(4) = 1.;
1716  tre_noise.getError()(4, 4) = 10.;
1717  double invp0 = tre_noise.getVector()(4);
1718  double var_invp0 = tre_noise.getError()(4, 4);
1719 
1720  // Loop over fits, starting at the high path distance (low momentum)
1721  // end.
1722 
1723  for (auto it = trackmap.rbegin(); it != trackmap.rend(); ++it) {
1724  double s = it->first;
1725  const KHitTrack& trh = it->second;
1726 
1727  // Ignore non-optimal fits.
1728 
1730  if (stat != KFitTrack::OPTIMAL) continue;
1731 
1732  // See if this track is far enough from the previous sample to
1733  // make a new sample point.
1734 
1735  if (std::abs(s - s_sample) > fMinSampleDist) {
1736 
1737  // Propagate tracks to the current track surface.
1738 
1739  std::shared_ptr<const Surface> psurf = trh.getSurface();
1740  std::optional<double> dist_inf = prop.err_prop(tre_inf, psurf, Propagator::UNKNOWN, false);
1741  std::optional<double> dist_range =
1742  prop.vec_prop(trk_range, psurf, Propagator::UNKNOWN, false);
1743  std::optional<double> dist_noise =
1744  prop.noise_prop(tre_noise, psurf, Propagator::UNKNOWN, true);
1745 
1746  // All propagations should normally succeed. If they don't,
1747  // ignore this sample for the purpose of updating the momentum.
1748 
1749  bool momentum_updated = false;
1750  if (dist_inf && dist_range && dist_noise) {
1751 
1752  // Extract the momentum at the new sample point.
1753 
1754  double invp1 = tre_noise.getVector()(4);
1755  double var_invp1 = tre_noise.getError()(4, 4);
1756 
1757  // Get the average momentum and error for this pair of
1758  // sample points, and other data.
1759 
1760  double invp = 0.5 * (invp0 + invp1);
1761  double var_invp = 0.5 * (var_invp0 + var_invp1);
1762  double mass = tre_inf.Mass();
1763  double beta = std::sqrt(1. + mass * mass * invp * invp);
1764  double invbp = invp / beta;
1765 
1766  // Extract slope subvectors and sub-error-matrices.
1767  // We have the following variables.
1768  //
1769  // slope0 - Predicted slope vector (from infinite momentum track).
1770  // slope1 - Measured slope vector (from new sample point).
1771  // defl - Deflection (slope residual = difference between measured
1772  // and predicted slope vector).
1773  // err0 - Slope error matrix of prediction.
1774  // err1 - Slope error matrix of measurement
1775  // errc - Slope residual error matrix = err0 + err1.
1776  // errn - Noise slope error matrix divided by (1/beta*momentum).
1777 
1778  KVector<2>::type slope0 = project(tre_inf.getVector(), ublas::range(2, 4));
1779  KVector<2>::type slope1 = project(trh.getVector(), ublas::range(2, 4));
1780  KVector<2>::type defl = slope1 - slope0;
1781  KSymMatrix<2>::type err0 =
1782  project(tre_inf.getError(), ublas::range(2, 4), ublas::range(2, 4));
1783  KSymMatrix<2>::type err1 = project(trh.getError(), ublas::range(2, 4), ublas::range(2, 4));
1784  KSymMatrix<2>::type errc = err0 + err1;
1785  KSymMatrix<2>::type errn =
1786  project(tre_noise.getError(), ublas::range(2, 4), ublas::range(2, 4));
1787  errn -= err0;
1788  errn /= invbp;
1789 
1790  // Calculate updated average momentum and error.
1791 
1792  double new_invp = invp;
1793  double new_var_invp = var_invp;
1794  update_momentum(defl, errc, errn, mass, new_invp, new_var_invp);
1795 
1796  // Calculate updated momentum and error at the second sample
1797  // point.
1798 
1799  double dp = 1. / new_invp - 1. / invp;
1800  invp0 = 1. / (1. / invp1 + dp);
1801  var_invp0 = new_var_invp;
1802  momentum_updated = true;
1803 
1804  // Make sure that updated momentum is not less than minimum
1805  // allowed momentum.
1806 
1807  double invp_range = trk_range.getVector()(4);
1808  if (invp0 > invp_range) invp0 = invp_range;
1809  }
1810 
1811  // Update sample.
1812 
1813  if (momentum_updated) {
1814  s_sample = s;
1815  tre_inf = trh;
1816  tre_inf.getVector()(4) = 0.;
1817  tre_inf.getError()(4, 4) = 0.;
1818  double invp_range = trk_range.getVector()(4);
1819  trk_range = trh;
1820  trk_range.getVector()(4) = invp_range;
1821  tre_noise = trh;
1822  tre_noise.getVector()(4) = invp0;
1823  tre_noise.getError()(4, 4) = var_invp0;
1824  }
1825  }
1826  }
1827 
1828  // Propagate noise track to starting (high momentum) track surface
1829  // to get final starting momentum. This propagation should normally
1830  // always succeed, but if it doesn't, don't update the track.
1831 
1832  const KHitTrack& trh0 = itend[0]->second;
1833  std::shared_ptr<const Surface> psurf = trh0.getSurface();
1834  std::optional<double> dist_noise = prop.noise_prop(tre_noise, psurf, Propagator::UNKNOWN, true);
1835  result = dist_noise.has_value();
1836 
1837  // Update momentum-estimating track.
1838 
1839  mf::LogInfo log("KalmanFilterAlg");
1840  if (result) tremom = tre_noise;
1841 
1842  // Done.
1843 
1844  return result;
1845 }
1846 
1847 /// Estimate track momentum using either range or multiple scattering.
1848 ///
1849 /// Arguments:
1850 ///
1851 /// trg - Global track whose momentum is to be updated.
1852 /// prop - Propagator.
1853 /// tremom - Track with momentum estimate.
1854 ///
1855 /// Returns: True if success.
1856 ///
1857 bool
1859  const Propagator& prop,
1860  KETrack& tremom) const
1861 {
1862  mf::LogInfo log("KalmanFilterAlg");
1863  double invp_range = 0.;
1864  double invp_ms = 0.;
1865 
1866  // Get multiple scattering momentum estimate.
1867 
1868  KETrack tremom_ms;
1869  bool ok_ms = false;
1870  if (fFitMomMS) {
1871  ok_ms = fitMomentumMS(trg, prop, tremom_ms);
1872  if (ok_ms) {
1873  KGTrack trg_ms(trg);
1874  ok_ms = updateMomentum(tremom_ms, prop, trg_ms);
1875  if (ok_ms) {
1876  invp_ms = trg_ms.startTrack().getVector()(4);
1877  double var_invp = trg_ms.startTrack().getError()(4, 4);
1878  double p = 0.;
1879  if (invp_ms != 0.) p = 1. / invp_ms;
1880  double err_p = p * p * std::sqrt(var_invp);
1881  log << "Multiple scattering momentum estimate = " << p << "+-" << err_p << "\n";
1882  }
1883  }
1884  }
1885 
1886  // Get range momentum estimate.
1887 
1888  KETrack tremom_range;
1889  bool ok_range = false;
1890  if (fFitMomRange) {
1891  ok_range = fitMomentumRange(trg, tremom_range);
1892  if (ok_range) {
1893  KGTrack trg_range(trg);
1894  ok_range = updateMomentum(tremom_range, prop, trg_range);
1895  if (ok_range) {
1896  invp_range = trg_range.startTrack().getVector()(4);
1897  double var_invp = trg_range.startTrack().getError()(4, 4);
1898  double p = 0.;
1899  if (invp_range != 0.) p = 1. / invp_range;
1900  double err_p = p * p * std::sqrt(var_invp);
1901  log << "Range momentum estimate = " << p << "+-" << err_p << "\n";
1902  }
1903  }
1904  }
1905 
1906  bool result = false;
1907  if (ok_range) {
1908  tremom = tremom_range;
1909  result = ok_range;
1910  }
1911  else if (ok_ms) {
1912  tremom = tremom_ms;
1913  result = ok_ms;
1914  }
1915  return result;
1916 }
1917 
1918 /// Set track momentum at each track surface.
1919 ///
1920 /// Arguments:
1921 ///
1922 /// tremom - Track containing momentum estimate.
1923 /// prop - Propagator.
1924 /// trg - Global track to be updated.
1925 ///
1926 /// The track containing the momentum estimate is propagated to the
1927 /// first or last track fit of the global track (whichever is closer),
1928 /// then the momentum estimate is transfered to that track fit. In
1929 /// similar fashion, the momentum estimate is successively transfered
1930 /// from that track fit to each track fit of the global track.
1931 ///
1932 /// Only momentum track parameters of the global track fits are
1933 /// updated. Other track parameters and their errors are unmodified.
1934 /// Unreachable track fits are deleted from the global track. Overall
1935 /// failure will occur if the momentum-estimating track can't be
1936 /// propagated to the initial track fit, or if the final global track
1937 /// has no valid track fits.
1938 ///
1939 bool
1941  const Propagator& prop,
1942  KGTrack& trg) const
1943 {
1944  // Get modifiable track map.
1945 
1946  std::multimap<double, KHitTrack>& trackmap = trg.getTrackMap();
1947 
1948  // If track map is empty, immediately return failure.
1949 
1950  if (trackmap.size() == 0) return false;
1951 
1952  // Make trial propagations to the first and last track fit to figure
1953  // out which track fit is closer to the momentum estimating track.
1954 
1955  // Find distance to first track fit.
1956 
1957  KETrack tre0(tremom);
1958  std::optional<double> dist0 =
1959  prop.vec_prop(tre0, trackmap.begin()->second.getSurface(), Propagator::UNKNOWN, false, 0, 0);
1960  // Find distance to last track fit.
1961 
1962  KETrack tre1(tremom);
1963  std::optional<double> dist1 =
1964  prop.vec_prop(tre1, trackmap.rbegin()->second.getSurface(), Propagator::UNKNOWN, false, 0, 0);
1965 
1966  // Based on distances, make starting iterator and direction flag.
1967 
1968  bool forward = true;
1969  std::multimap<double, KHitTrack>::iterator it = trackmap.begin();
1970  if (dist0) {
1971 
1972  // Propagation to first track succeeded.
1973 
1974  if (dist1) {
1975 
1976  // Propagation to both ends succeeded. If the last track is
1977  // closer, initialize iterator and direction flag for reverse
1978  // direction.
1979 
1980  if (std::abs(*dist0) > std::abs(*dist1)) {
1981  it = trackmap.end();
1982  --it;
1983  forward = false;
1984  }
1985  }
1986  }
1987  else {
1988 
1989  // Propagation to first track failed. Initialize iterator and
1990  // direction flag for reverse direction, provided that the
1991  // propagation to the last track succeeded.
1992 
1993  if (dist1) {
1994  it = trackmap.end();
1995  --it;
1996  forward = false;
1997  }
1998  else {
1999 
2000  // Propagation to both ends failed. Return failure.
2001 
2002  return false;
2003  }
2004  }
2005 
2006  // Loop over track fits in global track.
2007 
2008  KETrack tre(tremom);
2009  bool done = false;
2010  while (!done) {
2011  KHitTrack& trh = it->second;
2012  if (trh.getStat() == KFitTrack::INVALID)
2013  throw cet::exception("KalmanFilterAlg")
2014  << "KalmanFilterAlg::updateMomentum(): invalid track\n";
2015 
2016  // Propagate momentum-estimating track to current track surface
2017  // and update momentum.
2018 
2019  std::optional<double> dist = prop.noise_prop(tre, trh.getSurface(), Propagator::UNKNOWN, true);
2020 
2021  // Copy momentum to global track.
2022 
2023  std::multimap<double, KHitTrack>::iterator erase_it = trackmap.end();
2024  if (dist) {
2025  trh.getVector()(4) = tre.getVector()(4);
2026  trh.getError()(4, 0) = 0.;
2027  trh.getError()(4, 1) = 0.;
2028  trh.getError()(4, 2) = 0.;
2029  trh.getError()(4, 3) = 0.;
2030  trh.getError()(4, 4) = tre.getError()(4, 4);
2031  }
2032  else {
2033 
2034  // If the propagation failed, remember that we are supposed to
2035  // erase this track from the global track.
2036 
2037  erase_it = it;
2038  }
2039 
2040  // Advance the iterator and set the done flag.
2041 
2042  if (forward) {
2043  ++it;
2044  done = (it == trackmap.end());
2045  }
2046  else {
2047  done = (it == trackmap.begin());
2048  if (!done) --it;
2049  }
2050 
2051  // Update momentum-estimating track from just-updated global track
2052  // fit, or erase global track fit.
2053 
2054  if (erase_it == trackmap.end())
2055  tre = trh;
2056  else
2057  trackmap.erase(erase_it);
2058  }
2059 
2060  return not empty(trackmap);
2061 }
2062 
2063 /// Iteratively smooth a track.
2064 ///
2065 /// Arguments:
2066 ///
2067 /// nsmooth - Number of iterations.
2068 /// trg - Track to be smoothed.
2069 /// prop - Propagator.
2070 ///
2071 /// Returns: True if success.
2072 ///
2073 /// The initial track should have been unidirectionally fit.
2074 ///
2075 bool
2076 trkf::KalmanFilterAlg::smoothTrackIter(int nsmooth, KGTrack& trg, const Propagator& prop) const
2077 {
2078  bool ok = true;
2079 
2080  // Call smoothTrack method in a loop.
2081 
2082  for (int ismooth = 0; ok && ismooth < nsmooth - 1; ++ismooth) {
2083  KGTrack trg1{trg.getPrefPlane()};
2084  ok = smoothTrack(trg, &trg1, prop);
2085  if (ok) trg = trg1;
2086  }
2087 
2088  // Do one final smooth without generating a new unidirectional
2089  // track.
2090 
2091  if (ok) ok = smoothTrack(trg, 0, prop);
2092 
2093  // Done.
2094 
2095  return ok;
2096 }
2097 
2098 /// Clean track by removing noise hits near endpoints.
2099 ///
2100 /// Arguments:
2101 ///
2102 /// trg - Track.
2103 ///
2104 void
2106 {
2107  // Get hold of a modifiable track map from trg.
2108 
2109  std::multimap<double, KHitTrack>& trackmap = trg.getTrackMap();
2110 
2111  // Do an indefinite loop until we no longer find any dirt.
2112 
2113  bool done = false;
2114  while (!done) {
2115 
2116  // If track map has fewer than fMaxNoiseHits tracks, then this is a
2117  // noise track. Clear the map, making the whole track invalid.
2118 
2119  int ntrack = trackmap.size();
2120  if (ntrack <= fMaxNoiseHits) {
2121  for (auto const& trackmap_entry : trackmap) {
2122  const KHitTrack& trh = trackmap_entry.second;
2123  const KHitBase& hit = *(trh.getHit());
2124  auto marker_it = fMarkerMap.find(hit.getID());
2125  if (marker_it != fMarkerMap.end()) {
2126  TMarker* marker = marker_it->second;
2127  marker->SetMarkerColor(kGreen);
2128  }
2129  }
2130  trackmap.clear();
2131  done = true;
2132  break;
2133  }
2134 
2135  // Make sure the first two and last two tracks belong to different
2136  // views. If not, remove the first or last track.
2137 
2138  if (ntrack >= 2) {
2139 
2140  // Check start.
2141 
2142  std::multimap<double, KHitTrack>::iterator it = trackmap.begin();
2143  const KHitTrack& trh1a = (*it).second;
2144  const KHitBase& hit1a = *(trh1a.getHit());
2145  int plane1 = hit1a.getMeasPlane();
2146  double chisq1 = hit1a.getChisq();
2147  ++it;
2148  const KHitTrack& trh2a = (*it).second;
2149  const KHitBase& hit2a = *(trh2a.getHit());
2150  int plane2 = hit2a.getMeasPlane();
2151  double chisq2 = hit2a.getChisq();
2152  if ((plane1 >= 0 && plane2 >= 0 && plane1 == plane2) || chisq1 > fMaxEndChisq ||
2153  chisq2 > fMaxEndChisq) {
2154  auto marker_it = fMarkerMap.find(hit1a.getID());
2155  if (marker_it != fMarkerMap.end()) {
2156  TMarker* marker = marker_it->second;
2157  marker->SetMarkerColor(kGreen);
2158  }
2159  trackmap.erase(trackmap.begin(), it);
2160  done = false;
2161  continue;
2162  }
2163 
2164  // Check end.
2165 
2166  it = trackmap.end();
2167  --it;
2168  const KHitTrack& trh1b = (*it).second;
2169  const KHitBase& hit1b = *(trh1b.getHit());
2170  plane1 = hit1b.getMeasPlane();
2171  chisq1 = hit1b.getChisq();
2172  --it;
2173  const KHitTrack& trh2b = (*it).second;
2174  const KHitBase& hit2b = *(trh2b.getHit());
2175  plane2 = hit2b.getMeasPlane();
2176  chisq2 = hit2b.getChisq();
2177  if ((plane1 >= 0 && plane2 >= 0 && plane1 == plane2) || chisq1 > fMaxEndChisq ||
2178  chisq2 > fMaxEndChisq) {
2179  ++it;
2180  auto marker_it = fMarkerMap.find(hit1b.getID());
2181  if (marker_it != fMarkerMap.end()) {
2182  TMarker* marker = marker_it->second;
2183  marker->SetMarkerColor(kGreen);
2184  }
2185  trackmap.erase(it, trackmap.end());
2186  done = false;
2187  continue;
2188  }
2189  }
2190 
2191  // Loop over successive pairs of elements of track map. Look for
2192  // adjacent elements with distance separation greater than
2193  // fGapDist.
2194 
2195  std::multimap<double, KHitTrack>::iterator it = trackmap.begin();
2196  std::multimap<double, KHitTrack>::iterator jt = trackmap.end();
2197  int nb = 0; // Number of elements from begin to jt.
2198  int ne = ntrack; // Number of elements it to end.
2199  bool found_noise = false;
2200  for (; it != trackmap.end(); ++it, ++nb, --ne) {
2201  if (jt == trackmap.end())
2202  jt = trackmap.begin();
2203  else {
2204  if (nb < 1)
2205  throw cet::exception("KalmanFilterAlg")
2206  << "KalmanFilterAlg::cleanTrack(): nb not positive\n";
2207  if (ne < 1)
2208  throw cet::exception("KalmanFilterAlg")
2209  << "KalmanFilterAlg::cleanTrack(): ne not positive\n";
2210 
2211  double disti = (*it).first;
2212  double distj = (*jt).first;
2213  double sep = disti - distj;
2214  if (sep < 0.)
2215  throw cet::exception("KalmanFilterAlg")
2216  << "KalmanFilterAlg::fitMomentumRange(): negative separation\n";
2217 
2218  if (sep > fGapDist) {
2219 
2220  // Found a gap. See if we want to trim track.
2221 
2222  if (nb <= fMaxNoiseHits) {
2223 
2224  // Trim front.
2225 
2226  found_noise = true;
2227  for (auto jt = trackmap.begin(); jt != it; ++jt) {
2228  const KHitTrack& trh = jt->second;
2229  const KHitBase& hit = *(trh.getHit());
2230  auto marker_it = fMarkerMap.find(hit.getID());
2231  if (marker_it != fMarkerMap.end()) {
2232  TMarker* marker = marker_it->second;
2233  marker->SetMarkerColor(kGreen);
2234  }
2235  }
2236  trackmap.erase(trackmap.begin(), it);
2237  break;
2238  }
2239  if (ne <= fMaxNoiseHits) {
2240 
2241  // Trim back.
2242 
2243  found_noise = true;
2244  for (auto jt = it; jt != trackmap.end(); ++jt) {
2245  const KHitTrack& trh = jt->second;
2246  const KHitBase& hit = *(trh.getHit());
2247  auto marker_it = fMarkerMap.find(hit.getID());
2248  if (marker_it != fMarkerMap.end()) {
2249  TMarker* marker = marker_it->second;
2250  marker->SetMarkerColor(kGreen);
2251  }
2252  }
2253  trackmap.erase(it, trackmap.end());
2254  break;
2255  }
2256  }
2257  ++jt;
2258  }
2259  }
2260 
2261  // Decide if we are done.
2262 
2263  done = !found_noise;
2264  }
2265  if (fGTrace && fCanvases.size() > 0) fCanvases.back()->Update();
2266 }
bool fGTrace
Graphical trace flag.
ublas::symmetric_matrix< double, ublas::lower, ublas::row_major, ublas::bounded_array< double, N *(N+1)/2 > > type
intermediate_table::iterator iterator
int getPrefPlane() const
Definition: KGTrack.h:52
const TrackError & getError() const
Track error matrix.
Definition: KETrack.h:53
const std::shared_ptr< const KHitBase > & getHit() const
Measurement.
Definition: KHitTrack.h:53
bool fFitMomRange
Fit momentum using range.
double fMaxPropDist
Maximum propagation distance to candidate surface.
double fGapDist
Minimum gap distance.
static constexpr double nb
Definition: Units.h:81
QList< Entry > entry
double beta(double KE, const simb::MCParticle *part)
geo::Length_t DetHalfWidth(geo::TPCID const &tpcid) const
Returns the half width of the active volume of the specified TPC.
bool fFitMomMS
Fit momentum using multiple scattering.
double Mass() const
Based on pdg code.
Definition: KTrack.cxx:129
const std::shared_ptr< const Surface > & getSurface() const
Surface.
Definition: KTrack.h:55
static QCString result
double PointingError() const
Pointing error (radians).
Definition: KETrack.cxx:67
double getPredDistance() const
Prediction distance.
Definition: KHitBase.h:84
void recalibrate()
Recalibrate track map.
Definition: KGTrack.cxx:94
std::map< int, TMarker * > fMarkerMap
Markers in current canvas.
int getMeasPlane() const
Measurement plane index.
Definition: KHitBase.h:98
const std::list< KHitGroup > & getUnused() const
Definition: KHitContainer.h:80
ublas::matrix< double, ublas::row_major, ublas::bounded_array< double, N *M > > type
bool smoothTrackIter(int niter, KGTrack &trg, const Propagator &prop) const
Iteratively smooth a track.
const KHitTrack & endTrack() const
Track at end point.
Definition: KGTrack.cxx:42
KSymMatrix< 5 >::type TrackError
Track error matrix, dimension 5x5.
const std::vector< std::shared_ptr< const KHitBase > > & getHits() const
Measurement collection accessor.
Definition: KHitGroup.h:56
double fMaxEndChisq
Maximum incremental chisquare for endpoint hit.
bool fTrace
Trace flag.
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Planar surface parallel to x-axis.
KalmanFilterAlg(const fhicl::ParameterSet &pset)
Constructor.
struct vector vector
A collection of KHitGroups.
void cleanTrack(KGTrack &trg) const
Clean track by removing noise hits near endpoints.
void addTrack(const KHitTrack &trh)
Add track.
Definition: KGTrack.cxx:81
double getPath() const
Estimated path distance.
Definition: KHitGroup.h:62
intermediate_table::const_iterator const_iterator
FitStatus
Fit status enum.
Definition: KFitTrack.h:41
int getID() const
Unique id.
Definition: KHitBase.h:105
string dir
std::vector< double > fGTraceZMin
Graphical trace minimum z for each view.
const std::multimap< double, KHitTrack > & getTrackMap() const
KHitTrack collection, indexed by path distance.
Definition: KGTrack.h:59
art framework interface to geometry description
std::optional< double > vec_prop(KTrack &trk, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, TrackMatrix *prop_matrix=0, TrackError *noise_matrix=0) const
Propagate without error (long distance).
Definition: Propagator.cxx:53
double fMaxLDist
Maximum distance for linearized propagation.
double fMaxSortDist
Sort high distance threshold.
double fMaxSeedPredDist
Maximum prediciton distance to accept a hit in seed phase.
double fMaxPredDist
Maximum prediciton distance to accept a hit.
M::value_type trace(const M &m)
string project
Definition: conf.py:19
std::optional< double > noise_prop(KETrack &tre, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, KTrack *ref=0) const
Propagate with error and noise.
Definition: Propagator.cxx:386
T abs(T value)
bool fitMomentumMS(const KGTrack &trg, const Propagator &prop, KETrack &tremom) const
Estimate track momentum using multiple scattering.
void setChisq(double chisq)
Set chisquare.
Definition: KFitTrack.h:73
bool buildTrack(const KTrack &trk, KGTrack &trg, const Propagator &prop, const Propagator::PropDirection dir, KHitContainer &hits, bool linear) const
Make a new track.
double getPath() const
Propagation distance.
Definition: KFitTrack.h:66
Kalman filter measurement class template.
ublas::vector< double, ublas::bounded_array< double, N > > type
bool syminvert(ublas::symmetric_matrix< T, TRI, L, A > &m)
double fGTraceWH
Window height.
TVirtualPad * fInfoPad
Information pad.
T get(std::string const &key) const
Definition: ParameterSet.h:271
bool extendTrack(KGTrack &trg, const Propagator &prop, KHitContainer &hits) const
Add hits to existing track.
std::optional< double > err_prop(KETrack &tre, const std::shared_ptr< const Surface > &psurf, PropDirection dir, bool doDedx, KTrack *ref=0, TrackMatrix *prop_matrix=0) const
Propagate with error, but without noise.
Definition: Propagator.cxx:346
std::vector< std::unique_ptr< TCanvas > > fCanvases
Graphical trace canvases.
p
Definition: test.py:223
const std::list< KHitGroup > & getUnsorted() const
Definition: KHitContainer.h:75
double fGTraceXMin
Graphical trace minimum x.
double fGoodPErr
Pointing error threshold for switching to free propagation.
TPaveText * fMessages
Message box.
void setStat(FitStatus stat)
Set fit status.
Definition: KFitTrack.h:74
double fMaxSeedIncChisq
Maximum incremental chisquare to accept a hit in seed phase.
double fMaxSmoothIncChisq
Maximum incremental chisquare to accept a hit in smooth phase.
double fMaxIncChisq
Maximum incremental chisquare to accept a hit.
bool smoothTrack(KGTrack &trg, KGTrack *trg1, const Propagator &prop) const
Smooth track.
double fMinSampleDist
Minimum sample distance (for momentum measurement).
Kalman Filter.
void err(const char *fmt,...)
Definition: message.cpp:226
Detector simulation of raw signals on wires.
double fMinSortDist
Sort low distance threshold.
const TrackVector & getVector() const
Track state vector.
Definition: KTrack.h:56
virtual void update(KETrack &tre) const =0
Update track method.
virtual double getChisq() const =0
Return incremental chisquare.
const KVector< N >::type & getMeasVector() const
Measurement vector.
Definition: KHit.h:110
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
double fGTraceXMax
Graphical trace maximum x.
void sort(const KTrack &trk, bool addUnsorted, const Propagator &prop, Propagator::PropDirection dir=Propagator::UNKNOWN)
(Re)sort objects in unsorted and sorted lists.
std::vector< TVirtualPad * > fPads
View pads in current canvas.
size_t numHits() const
Number of measurements in track.
Definition: KGTrack.h:66
double fGTraceWW
Window width.
virtual bool predict(const KETrack &tre, const Propagator &prop, const KTrack *ref=0) const =0
Prediction method (return false if fail).
const std::shared_ptr< const Surface > & getSurface() const
Surface accessor.
Definition: KHitGroup.h:50
bool combineFit(const KFitTrack &trf)
Combine two tracks.
Definition: KFitTrack.cxx:66
void setPath(double path)
Set propagation distance.
Definition: KFitTrack.h:72
std::vector< double > fGTraceZMax
Graphical trace maximum z for each view.
int fMaxSamePlane
Maximum consecutive hits in same plane.
Line surface perpendicular to x-axis.
list x
Definition: train.py:276
int fMaxNoiseHits
Maximum number of hits in noise cluster.
double fMaxPErr
Maximum pointing error for free propagation.
bool fitMomentumRange(const KGTrack &trg, KETrack &tremom) const
Estimate track momentum using range.
bool isValid() const
Validity flag.
Definition: KGTrack.h:79
bool fitMomentum(const KGTrack &trg, const Propagator &prop, KETrack &tremom) const
Estimate track momentum using either range or multiple scattering.
bool updateMomentum(const KETrack &tremom, const Propagator &prop, KGTrack &trg) const
Set track momentum at each track surface.
double getChisq() const
Fit chisquare.
Definition: KFitTrack.h:67
static QCString * s
Definition: config.cpp:1042
const KHitTrack & startTrack() const
Track at start point.
Definition: KGTrack.cxx:31
const std::shared_ptr< const Surface > & getMeasSurface() const
Measurement surface.
Definition: KHitBase.h:91
int fMinLHits
Minimum number of hits to turn off linearized propagation.
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
decltype(auto) constexpr empty(T &&obj)
ADL-aware version of std::empty.
Definition: StdUtils.h:97
bool isValid() const
Test if track is valid.
Definition: KTrack.cxx:91
PropDirection
Propagation direction enum.
Definition: Propagator.h:94
FitStatus getStat() const
Fit status.
Definition: KFitTrack.h:68
int getPlane() const
Plane index.
Definition: KHitGroup.h:53
int fPlane
Preferred view plane.
const std::list< KHitGroup > & getSorted() const
Definition: KHitContainer.h:70