CMS 3D CMS Logo

KFFittingSmoother.h
Go to the documentation of this file.
1 
12 
15 
20 
21 
22 namespace {
23 
24  struct KFFittingSmootherParam {
25 
26  explicit KFFittingSmootherParam(const edm::ParameterSet & conf) :
27  theEstimateCut(conf.getParameter<double>("EstimateCut"))
28  ,theMaxFractionOutliers(conf.getParameter<double>("MaxFractionOutliers"))
29  ,theMaxNumberOfOutliers(conf.getParameter<int>("MaxNumberOfOutliers"))
30  ,theNoOutliersBeginEnd(conf.getParameter<bool>("NoOutliersBeginEnd"))
31  ,theMinDof(conf.getParameter<int>("MinDof"))
32  ,theMinNumberOfHits(conf.getParameter<int>("MinNumberOfHits"))
33  ,rejectTracksFlag(conf.getParameter<bool>("RejectTracks"))
34  ,breakTrajWith2ConsecutiveMissing(conf.getParameter<bool>("BreakTrajWith2ConsecutiveMissing"))
35  ,noInvalidHitsBeginEnd(conf.getParameter<bool>("NoInvalidHitsBeginEnd"))
36  {}
37 
38  double theEstimateCut;
39 
40  float theMaxFractionOutliers;
41  int theMaxNumberOfOutliers;
42  bool theNoOutliersBeginEnd;
43  int theMinDof;
44 
45  int theMinNumberOfHits;
46  bool rejectTracksFlag;
47  bool breakTrajWith2ConsecutiveMissing;
48  bool noInvalidHitsBeginEnd;
49 
50 
51  };
52 
53  class KFFittingSmoother final : public TrajectoryFitter, private KFFittingSmootherParam {
54 
55  public:
56  ~KFFittingSmoother() override {}
57 
58  KFFittingSmoother(const TrajectoryFitter& aFitter,
59  const TrajectorySmoother& aSmoother,
60  const edm::ParameterSet & conf) :
61  KFFittingSmootherParam(conf),
62  theFitter(aFitter.clone()),
63  theSmoother(aSmoother.clone())
64  {}
65 
66  private:
67 
68 
70  desc.add<double>("EstimateCut",-1);
71  desc.add<double>("MaxFractionOutliers",0.3);
72  desc.add<int>("MaxNumberOfOutliers",3);
73  desc.add<int>("MinDof",2);
74  desc.add<bool>("NoOutliersBeginEnd",false);
75  desc.add<int>("MinNumberOfHits",5);
76  desc.add<bool>("RejectTracks",true);
77  desc.add<bool>("BreakTrajWith2ConsecutiveMissing",true);
78  desc.add<bool>("NoInvalidHitsBeginEnd",true);
79  desc.add<double>("LogPixelProbabilityCut",0);
80 
81  }
82 
83 
84  Trajectory fitOne(const Trajectory& t, fitType type) const override;
85  Trajectory fitOne(const TrajectorySeed& aSeed,
86  const RecHitContainer& hits,
87  const TrajectoryStateOnSurface& firstPredTsos, fitType type) const override;
88  Trajectory fitOne(const TrajectorySeed& aSeed,
89  const RecHitContainer& hits, fitType type) const override;
90 
91  const TrajectoryFitter* fitter() const {return theFitter.get();}
92  const TrajectorySmoother* smoother() const {return theSmoother.get();}
93 
94  std::unique_ptr<TrajectoryFitter> clone() const override {
95  return std::unique_ptr<TrajectoryFitter>(
96  new KFFittingSmoother(*theFitter,
97  *theSmoother,
98  *this)
99  );
100  }
101 
102 
103  void setHitCloner(TkCloner const * hc) override {
104  theFitter->setHitCloner(hc);
105  theSmoother->setHitCloner(hc);
106  }
107 
108  KFFittingSmoother(const TrajectoryFitter& aFitter,
109  const TrajectorySmoother& aSmoother,
110  KFFittingSmootherParam const & other) :
111  KFFittingSmootherParam(other),
112  theFitter(aFitter.clone()),
113  theSmoother(aSmoother.clone())
114  {}
115 
116 
117 
118  Trajectory smoothingStep(Trajectory && fitted) const {
119  if (theEstimateCut>0) {
120  // remove "outlier" at the end of Traj
121  while (!fitted.empty() && fitted.foundHits()>=theMinNumberOfHits &&
122  ( !fitted.lastMeasurement().recHitR().isValid()
123  || ( fitted.lastMeasurement().recHitR().det()!=nullptr && fitted.lastMeasurement().estimate()>theEstimateCut)
124  )
125  )
126  fitted.pop();
127  if (fitted.foundHits() < theMinNumberOfHits ) return Trajectory();
128  }
129  return theSmoother->trajectory(fitted);
130  }
131 
132  private:
133 
134  const std::unique_ptr<TrajectoryFitter> theFitter;
135  const std::unique_ptr<TrajectorySmoother> theSmoother;
136 
137 
139  static bool checkForNans(const Trajectory &theTraj);
140 
141  friend class KFFittingSmootherESProducer;
142 
143  };
144 
145 
146 // #define VI_DEBUG
147 
148 #ifdef VI_DEBUG
149 #define DPRINT(x) std::cout << x << ": "
150 #define PRINT std::cout
151 #else
152 #define DPRINT(x) LogTrace(x)
153 #define PRINT LogTrace("")
154 #endif
155 
156 Trajectory KFFittingSmoother::fitOne(const Trajectory& t, fitType type) const {
157  if (!t.isValid() ) return Trajectory();
158  return smoothingStep(theFitter->fitOne(t,type));
159 }
160 
161 bool KFFittingSmoother::checkForNans(const Trajectory & theTraj)
162 {
163  if (edm::isNotFinite(theTraj.chiSquared ())) return false;
164  auto const & vtm = theTraj.measurements();
165  for (auto const & tm : vtm) {
166  if (edm::isNotFinite(tm.estimate())) return false;
167  auto const & v = tm.updatedState ().localParameters ().vector();
168  for (int i=0;i<5;++i) if (edm::isNotFinite(v[i])) return false;
169  if (!tm.updatedState ().curvilinearError ().posDef()) return false;//not a NaN by itself, but will lead to one
170  auto const & m = tm.updatedState ().curvilinearError ().matrix();
171  for (int i=0;i<5;++i)
172  for (int j=0;j<(i+1);++j) if (edm::isNotFinite(m(i,j))) return false;
173  }
174  return true;
175 }
176 
177 
178 namespace {
179  inline
180  void print(const std::string& header, const TrajectoryStateOnSurface & tsos) {
181  DPRINT("TrackFitters") << header << tsos.globalPosition().perp() << ' ' << tsos.globalPosition().z()
182  << ' ' << 1./tsos.signedInverseMomentum() << ' ' << 1./tsos.transverseCurvature()
183  << ' ' << tsos.globalMomentum().eta() << std::endl;
184  }
185 }
186 
187 
188 Trajectory KFFittingSmoother::fitOne(const TrajectorySeed& aSeed,
189  const RecHitContainer& hits,
190  const TrajectoryStateOnSurface& firstPredTsos,
191  fitType type) const
192 {
193  LogDebug("TrackFitters") << "In KFFittingSmoother::fit";
194 
195  print("firstPred ", firstPredTsos);
196 
197  if ( hits.empty() ) return Trajectory();
198 
199  RecHitContainer myHits = hits;
200  Trajectory tmp_first;
201 
202  //call the fitter
203  Trajectory smoothed = smoothingStep(theFitter->fitOne(aSeed, myHits, firstPredTsos));
204 
205 
206  do {
207 
208 #ifdef EDM_ML_DEBUG
209  //if no outliers the fit is done only once
210  for (unsigned int j=0;j<myHits.size();j++) {
211  if (myHits[j]->det())
212  LogTrace("TrackFitters") << "hit #:" << j+1 << " rawId=" << myHits[j]->det()->geographicalId().rawId()
213  << " validity=" << myHits[j]->isValid();
214  else
215  LogTrace("TrackFitters") << "hit #:" << j+1 << " Hit with no Det information";
216  }
217 #endif
218 
219 #if defined(VI_DEBUG) || defined(EDM_ML_DEBUG)
220 if (smoothed.isValid()) {
221  print("first state ",smoothed.firstMeasurement().updatedState());
222  print("last state ",smoothed.lastMeasurement().updatedState());
223 }
224 #endif
225 
226  bool hasNaN = false;
227  if ( !smoothed.isValid() || (hasNaN = !checkForNans(smoothed)) || ( smoothed.foundHits() < theMinNumberOfHits ) ) {
228  if(hasNaN) edm::LogWarning("TrackNaN")<<"Track has NaN or the cov is not pos-definite";
229  if ( smoothed.foundHits() < theMinNumberOfHits ) LogTrace("TrackFitters") << "smoothed.foundHits()<theMinNumberOfHits";
230  DPRINT("TrackFitters") << "smoothed invalid => trajectory rejected with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
231  if (rejectTracksFlag) {
232  return Trajectory();
233  } else {
234  std::swap(smoothed, tmp_first); // if first attempt, tmp_first would be invalid anyway
235  DPRINT("TrackFitters") << "smoothed invalid => returning orignal trajectory with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
236  }
237  break;
238  }
239 #ifdef EDM_ML_DEBUG
240  else {
241  LogTrace("TrackFitters") << "dump hits after smoothing";
242  Trajectory::DataContainer meas = smoothed.measurements();
243  for (Trajectory::DataContainer::iterator it=meas.begin();it!=meas.end();++it) {
244  LogTrace("TrackFitters") << "hit #" << meas.end()-it-1 << " validity=" << it->recHit()->isValid()
245  << " det=" << it->recHit()->geographicalId().rawId();
246  }
247  }
248 #endif
249 
250  if (myHits.size() !=smoothed.measurements().size())
251  DPRINT("TrackFitters") << "lost hits. before/after: " << myHits.size() <<'/' << smoothed.measurements().size()<< "\n";
252 
253  if ( theEstimateCut <= 0) break;
254 
255 
256  // Check if there are outliers
257 
258  auto msize = smoothed.measurements().size();
259  declareDynArray(unsigned int, msize, bad);
260  unsigned int nbad=0;
261  unsigned int ind=0;
262  unsigned int lastValid=smoothed.measurements().size();
263  for (auto const & tm : smoothed.measurements() ) {
264  if (tm.estimate()>theEstimateCut
265  && tm.recHitR().det()!=nullptr // do not consider outliers constraints and other special "hits"
266  ) bad[nbad++]=ind;
267  if (ind<lastValid && tm.recHitR().det()!=nullptr && tm.recHitR().isValid()) lastValid=ind;
268  ++ind;
269  }
270 
271  if (0==nbad) break;
272 
273  DPRINT("TrackFitters") << "size/found/outliers list " << smoothed.measurements().size() <<'/' << smoothed.foundHits() << ' ' << nbad << ": ";
274  for (auto i=0U; i<nbad; ++i) PRINT << bad[i] <<',';
275  PRINT << std::endl;
276 
277 
278 
279  if (
280  // (smoothed.foundHits() == theMinNumberOfHits) ||
281  int(nbad)>theMaxNumberOfOutliers ||
282  float(nbad) > theMaxFractionOutliers*float(smoothed.foundHits())
283  ) {
284  DPRINT("TrackFitters") << "smoothed low quality => trajectory with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
285  PRINT << "try to remove " << lastValid << std::endl;
286  nbad = 0; // try to short the traj... (below lastValid will be added)
287 
288  // do not perform outliers rejection if track is already low quality
289  /*
290  if ( rejectTracksFlag && (smoothed.chiSquared() > theEstimateCut*smoothed.ndof()) ) {
291  DPRINT("TrackFitters") << "smoothed low quality => trajectory rejected with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
292  return Trajectory();
293  } else {
294  DPRINT("TrackFitters") << "smoothed low quality => return original trajectory with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
295  }
296  break;
297  */
298  }
299 
300 
301  // always add last valid hit as outlier candidate
302  bad[nbad++]=lastValid;
303 
304  // if ( (smoothed.ndof()<theMinDof) | ) break;
305 
306 
307 
308  assert(smoothed.measurements().size() <= myHits.size());
309 
310  myHits.resize(smoothed.measurements().size()); // hits are only removed from the back...
311 
312  assert(smoothed.measurements().size() == myHits.size());
313 
314  declareDynArray(Trajectory,nbad,smoothedCand);
315 
316 
317  auto NHits = myHits.size();
319 
320  auto loc = nbad;
321  for (auto i=0U; i<nbad; ++i) {
322  auto j = NHits-bad[i]-1;
323  assert(myHits[j]->geographicalId() == smoothed.measurements()[bad[i]].recHitR().geographicalId());
324  auto removedHit = myHits[j];
325  myHits[j] = std::make_shared<InvalidTrackingRecHit>(*removedHit->det(), TrackingRecHit::missing);
326  smoothedCand[i] = smoothingStep(theFitter->fitOne(aSeed, myHits, firstPredTsos));
327  myHits[j] = removedHit;
328  if (smoothedCand[i].isValid() && smoothedCand[i].chiSquared()< minChi2) {
329  minChi2 = smoothedCand[i].chiSquared();
330  loc=i;
331  }
332  }
333 
334  if (loc == nbad) {
335  DPRINT("TrackFitters") <<"New trajectories all invalid"<< "\n";
336  return Trajectory();
337  }
338 
339  DPRINT("TrackFitters") << "outlier removed " << bad[loc] << '/' << minChi2 << " was " << smoothed.chiSquared()<< "\n";
340 
341  if (minChi2>smoothed.chiSquared()) {
342 
343  DPRINT("TrackFitters") << "removing outlier makes chi2 worse " << minChi2 << '/' << smoothed.chiSquared() << "\nOri: ";
344  for (auto const & tm : smoothed.measurements() )
345  PRINT << tm.recHitR().geographicalId() << '/' << tm.estimate() <<' ';
346  PRINT << "\nNew: ";
347  for (auto const & tm : smoothedCand[loc].measurements() )
348  PRINT << tm.recHitR().geographicalId() << '/' << tm.estimate() <<' ';
349  PRINT << "\n";
350 
351  // return Trajectory();
352  // break;
353  }
354 
355  std::swap(smoothed,tmp_first);
356  myHits[NHits-bad[loc]-1] = std::make_shared<InvalidTrackingRecHit>(*myHits[NHits-bad[loc]-1]->det(), TrackingRecHit::missing);
357  std::swap(smoothed,smoothedCand[loc]);
358  // firstTry=false;
359 
360  DPRINT("TrackFitters") << "new trajectory with nhits/chi2 " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
361 
362 
363 
364 
365  // Look if there are two consecutive invalid hits FIXME: take into account split matched hits!!!
366  if ( breakTrajWith2ConsecutiveMissing ) {
367  unsigned int firstinvalid = myHits.size();
368  for ( unsigned int j=0; j<myHits.size()-1; ++j )
369  {
370  if ( ((myHits[j ]->type() == TrackingRecHit::missing) && (myHits[j ]->geographicalId().rawId() != 0)) &&
371  ((myHits[j+1]->type() == TrackingRecHit::missing) && (myHits[j+1]->geographicalId().rawId() != 0))
372  && ( (myHits[j ]->geographicalId().rawId()&(~3)) != (myHits[j+1]->geographicalId().rawId()&(~3) ) ) // same gluedDet
373  )
374  {
375  firstinvalid = j;
376  DPRINT("TrackFitters") << "Found two consecutive missing hits. First invalid: " << firstinvalid << "\n";
377  break;
378  }
379  }
380 
381  //reject all the hits after the last valid before two consecutive invalid (missing) hits
382  //hits are sorted in the same order as in the track candidate FIXME??????
383  if (firstinvalid != myHits.size()) {
384  myHits.erase(myHits.begin()+firstinvalid,myHits.end());
385  smoothed = smoothingStep(theFitter->fitOne(aSeed, myHits, firstPredTsos));
386  DPRINT("TrackFitters") << "Trajectory shortened " << smoothed.foundHits() << '/' << smoothed.chiSquared() << "\n";
387  }
388  }
389 
390  } // do
391  while ( true );
392 
393  if ( smoothed.isValid() ) {
394  if ( noInvalidHitsBeginEnd
395  && !smoothed.empty() //should we send a warning ?
396  ) {
397  // discard latest dummy measurements
398  if (!smoothed.empty() && !smoothed.lastMeasurement().recHitR().isValid() )
399  LogTrace("TrackFitters") << "Last measurement is invalid";
400 
401  while (!smoothed.empty() && !smoothed.lastMeasurement().recHitR().isValid() )
402  smoothed.pop();
403 
404  //remove the invalid hits at the begin of the trajectory
405  if (!smoothed.empty() && !smoothed.firstMeasurement().recHitR().isValid() ) {
406  LogTrace("TrackFitters") << "First measurement is in`valid";
407  Trajectory tmpTraj(smoothed.seed(),smoothed.direction());
408  Trajectory::DataContainer & meas = smoothed.measurements();
409  auto it = meas.begin();
410  for ( ; it!=meas.end(); ++it )
411  if ( it->recHitR().isValid() ) break;
412  tmpTraj.push(std::move(*it),smoothed.chiSquared());//push the first valid measurement and set the same global chi2
413 
414  for (auto itt=it+1; itt!=meas.end();++itt)
415  tmpTraj.push(std::move(*itt),0);//add all the other measurements
416 
417  std::swap(smoothed,tmpTraj);
418 
419  } // if ( !smoothed.firstMeasurement().recHit()->isValid() )
420 
421  } // if ( noInvalidHitsBeginEnd )
422 
423  LogTrace("TrackFitters") << "end: returning smoothed trajectory with chi2="
424  << smoothed.chiSquared() ;
425 
426  //LogTrace("TrackFitters") << "dump hits before return";
427  //Trajectory::DataContainer meas = smoothed.measurements();
428  //for (Trajectory::DataContainer::iterator it=meas.begin();it!=meas.end();++it) {
429  //LogTrace("TrackFitters") << "hit #" << meas.end()-it-1 << " validity=" << it->recHit()->isValid()
430  //<< " det=" << it->recHit()->geographicalId().rawId();
431  //}
432 
433  }
434 
435  return smoothed;
436 
437 }
438 
439 
440 Trajectory KFFittingSmoother::fitOne(const TrajectorySeed& aSeed,
441  const RecHitContainer& hits,fitType type) const{
442 
443  throw cms::Exception("TrackFitters",
444  "KFFittingSmoother::fit(TrajectorySeed, <TransientTrackingRecHit>) not implemented");
445 
446  return Trajectory();
447 }
448 
449 } // namespace
#define LogDebug(id)
bool empty() const
True if trajectory has no measurements.
Definition: Trajectory.h:254
type
Definition: HCALResponse.h:21
int foundHits() const
Definition: Trajectory.h:225
T perp() const
Definition: PV3DBase.h:72
TrajectorySeed const & seed() const
Access to the seed used to reconstruct the Trajectory.
Definition: Trajectory.h:285
GlobalPoint globalPosition() const
S & print(S &os, JobReport::InputFile const &f)
Definition: JobReport.cc:66
PropagationDirection const & direction() const
Definition: Trajectory.cc:140
DataContainer const & measurements() const
Definition: Trajectory.h:196
std::vector< TrajectoryMeasurement > DataContainer
Definition: Trajectory.h:44
bool isNotFinite(T x)
Definition: isFinite.h:10
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:174
T z() const
Definition: PV3DBase.h:64
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define LogTrace(id)
bool isValid() const
Definition: Trajectory.h:279
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:187
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void pop()
Definition: Trajectory.cc:29
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
ConstRecHitPointer::element_type const & recHitR() const
#define DPRINT(x)
float chiSquared() const
Definition: Trajectory.h:262
T eta() const
Definition: PV3DBase.h:76
#define PRINT
GlobalVector globalMomentum() const
susybsm::HSCParticleCollection hc
Definition: classes.h:25
TrajectoryStateOnSurface const & updatedState() const
#define declareDynArray(T, n, x)
Definition: DynArray.h:59
def move(src, dest)
Definition: eostools.py:511